Line | Branch | Exec | Source |
---|---|---|---|
1 | // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- | ||
2 | // vi: set et ts=4 sw=4 sts=4: | ||
3 | // | ||
4 | // SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup PQ1BubbleDiscretization | ||
10 | * \brief Helper class constructing the dual grid finite volume geometries | ||
11 | * for the cvfe discretizazion method | ||
12 | */ | ||
13 | #ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH | ||
14 | #define DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH | ||
15 | |||
16 | #include <array> | ||
17 | |||
18 | #include <dune/common/exceptions.hh> | ||
19 | |||
20 | #include <dune/geometry/type.hh> | ||
21 | #include <dune/geometry/referenceelements.hh> | ||
22 | #include <dune/geometry/multilineargeometry.hh> | ||
23 | #include <dune/common/reservedvector.hh> | ||
24 | |||
25 | #include <dumux/common/math.hh> | ||
26 | #include <dumux/geometry/volume.hh> | ||
27 | #include <dumux/discretization/box/boxgeometryhelper.hh> | ||
28 | |||
29 | namespace Dumux { | ||
30 | |||
31 | //! Traits for an efficient corner storage for the PQ1Bubble method | ||
32 | template <class ct> | ||
33 | struct PQ1BubbleMLGeometryTraits : public Dune::MultiLinearGeometryTraits<ct> | ||
34 | { | ||
35 | // we use static vectors to store the corners as we know | ||
36 | // the maximum number of corners in advance (2^dim) | ||
37 | template< int mydim, int cdim > | ||
38 | struct CornerStorage | ||
39 | { | ||
40 | using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)+1>; | ||
41 | }; | ||
42 | }; | ||
43 | |||
44 | namespace Detail::PQ1Bubble { | ||
45 | |||
46 | template<Dune::GeometryType::Id gt> | ||
47 | struct OverlappingScvCorners; | ||
48 | |||
49 | template<> | ||
50 | struct OverlappingScvCorners<Dune::GeometryTypes::line> | ||
51 | { | ||
52 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
53 | static constexpr std::array<std::array<Key, 2>, 1> keys = {{ | ||
54 | { Key{0, 1}, Key{1, 1} } | ||
55 | }}; | ||
56 | }; | ||
57 | |||
58 | template<> | ||
59 | struct OverlappingScvCorners<Dune::GeometryTypes::triangle> | ||
60 | { | ||
61 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
62 | static constexpr std::array<std::array<Key, 3>, 1> keys = {{ | ||
63 | { Key{0, 1}, Key{1, 1}, Key{2, 1} } | ||
64 | }}; | ||
65 | }; | ||
66 | |||
67 | template<> | ||
68 | struct OverlappingScvCorners<Dune::GeometryTypes::quadrilateral> | ||
69 | { | ||
70 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
71 | static constexpr std::array<std::array<Key, 4>, 1> keys = {{ | ||
72 | { Key{2, 1}, Key{1, 1}, Key{0, 1}, Key{3, 1} } | ||
73 | }}; | ||
74 | }; | ||
75 | |||
76 | template<> | ||
77 | struct OverlappingScvCorners<Dune::GeometryTypes::tetrahedron> | ||
78 | { | ||
79 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
80 | static constexpr std::array<std::array<Key, 4>, 1> keys = {{ | ||
81 | { Key{0, 1}, Key{1, 1}, Key{2, 1}, Key{3, 1} } | ||
82 | }}; | ||
83 | }; | ||
84 | |||
85 | template<> | ||
86 | struct OverlappingScvCorners<Dune::GeometryTypes::hexahedron> | ||
87 | { | ||
88 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
89 | static constexpr std::array<std::array<Key, 6>, 1> keys = {{ | ||
90 | { Key{0, 1}, Key{2, 1}, Key{3, 1}, Key{1, 1}, Key{4, 1}, Key{5, 1} } | ||
91 | }}; | ||
92 | }; | ||
93 | |||
94 | |||
95 | template<Dune::GeometryType::Id gt> | ||
96 | struct OverlappingScvfCorners; | ||
97 | |||
98 | template<> | ||
99 | struct OverlappingScvfCorners<Dune::GeometryTypes::line> | ||
100 | { | ||
101 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
102 | static constexpr std::array<std::array<Key, 1>, 1> keys = {{ | ||
103 | { Key{0, 0} } | ||
104 | }}; | ||
105 | }; | ||
106 | |||
107 | template<> | ||
108 | struct OverlappingScvfCorners<Dune::GeometryTypes::triangle> | ||
109 | { | ||
110 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
111 | static constexpr std::array<std::array<Key, 2>, 3> keys = {{ | ||
112 | { Key{0, 1}, Key{1, 1} }, | ||
113 | { Key{0, 1}, Key{2, 1} }, | ||
114 | { Key{1, 1}, Key{2, 1} } | ||
115 | }}; | ||
116 | }; | ||
117 | |||
118 | template<> | ||
119 | struct OverlappingScvfCorners<Dune::GeometryTypes::quadrilateral> | ||
120 | { | ||
121 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
122 | static constexpr std::array<std::array<Key, 2>, 4> keys = {{ | ||
123 | { Key{0, 1}, Key{2, 1} }, | ||
124 | { Key{2, 1}, Key{1, 1} }, | ||
125 | { Key{0, 1}, Key{3, 1} }, | ||
126 | { Key{1, 1}, Key{3, 1} } | ||
127 | }}; | ||
128 | }; | ||
129 | |||
130 | template<> | ||
131 | struct OverlappingScvfCorners<Dune::GeometryTypes::tetrahedron> | ||
132 | { | ||
133 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
134 | static constexpr std::array<std::array<Key, 3>, 10> keys = {{ | ||
135 | { Key{0, 1}, Key{1, 1}, Key{2, 1} }, | ||
136 | { Key{0, 1}, Key{1, 1}, Key{3, 1} }, | ||
137 | { Key{0, 1}, Key{2, 1}, Key{3, 1} }, | ||
138 | { Key{1, 1}, Key{2, 1}, Key{3, 1} } | ||
139 | }}; | ||
140 | }; | ||
141 | |||
142 | template<> | ||
143 | struct OverlappingScvfCorners<Dune::GeometryTypes::hexahedron> | ||
144 | { | ||
145 | using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim) | ||
146 | static constexpr std::array<std::array<Key, 3>, 8> keys = {{ | ||
147 | { Key{4, 1}, Key{0, 1}, Key{2, 1} }, | ||
148 | { Key{4, 1}, Key{2, 1}, Key{1, 1} }, | ||
149 | { Key{4, 1}, Key{0, 1}, Key{3, 1} }, | ||
150 | { Key{4, 1}, Key{1, 1}, Key{3, 1} }, | ||
151 | { Key{5, 1}, Key{0, 1}, Key{2, 1} }, | ||
152 | { Key{5, 1}, Key{2, 1}, Key{1, 1} }, | ||
153 | { Key{5, 1}, Key{0, 1}, Key{3, 1} }, | ||
154 | { Key{5, 1}, Key{1, 1}, Key{3, 1} } | ||
155 | }}; | ||
156 | }; | ||
157 | |||
158 | } // end namespace Detail::PQ1Bubble | ||
159 | |||
160 | /*! | ||
161 | * \ingroup PQ1BubbleDiscretization | ||
162 | * \brief A class to create sub control volume and sub control volume face geometries per element | ||
163 | */ | ||
164 | template <class GridView, class ScvType, class ScvfType> | ||
165 | class PQ1BubbleGeometryHelper | ||
166 | { | ||
167 | using Scalar = typename GridView::ctype; | ||
168 | using GlobalPosition = typename Dune::FieldVector<Scalar, GridView::dimensionworld>; | ||
169 | using ScvCornerStorage = typename ScvType::Traits::CornerStorage; | ||
170 | using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage; | ||
171 | using LocalIndexType = typename ScvType::Traits::LocalIndexType; | ||
172 | |||
173 | using Element = typename GridView::template Codim<0>::Entity; | ||
174 | using Intersection = typename GridView::Intersection; | ||
175 | |||
176 | static constexpr auto dim = GridView::dimension; | ||
177 | static constexpr auto dimWorld = GridView::dimensionworld; | ||
178 | public: | ||
179 | |||
180 | PQ1BubbleGeometryHelper(const typename Element::Geometry& geometry) | ||
181 | : geo_(geometry) | ||
182 |
2/4✓ Branch 1 taken 12000 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12000 times.
✗ Branch 5 not taken.
|
247906 | , boxHelper_(geometry) |
183 | {} | ||
184 | |||
185 | //! Create a vector with the scv corners | ||
186 | 547821 | ScvCornerStorage getScvCorners(unsigned int localScvIdx) const | |
187 | { | ||
188 | // proceed according to number of corners of the element | ||
189 | 547821 | const auto type = geo_.type(); | |
190 | 547821 | const auto numBoxScv = boxHelper_.numScv(); | |
191 | // reuse box geometry helper for the corner scvs | ||
192 |
2/2✓ Branch 0 taken 424560 times.
✓ Branch 1 taken 123261 times.
|
547821 | if (localScvIdx < numBoxScv) |
193 | 424560 | return boxHelper_.getScvCorners(localScvIdx); | |
194 | |||
195 | 123261 | const auto localOverlappingScvIdx = localScvIdx-numBoxScv; | |
196 |
1/2✓ Branch 0 taken 123261 times.
✗ Branch 1 not taken.
|
123261 | if (type == Dune::GeometryTypes::triangle) |
197 | { | ||
198 | using Corners = Detail::PQ1Bubble::OverlappingScvCorners<Dune::GeometryTypes::triangle>; | ||
199 | 136992 | return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]); | |
200 | } | ||
201 |
1/2✓ Branch 0 taken 54765 times.
✗ Branch 1 not taken.
|
54765 | else if (type == Dune::GeometryTypes::quadrilateral) |
202 | { | ||
203 | using Corners = Detail::PQ1Bubble::OverlappingScvCorners<Dune::GeometryTypes::quadrilateral>; | ||
204 | 100704 | return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]); | |
205 | } | ||
206 |
1/2✓ Branch 0 taken 4413 times.
✗ Branch 1 not taken.
|
4413 | else if (type == Dune::GeometryTypes::tetrahedron) |
207 | { | ||
208 | using Corners = Detail::PQ1Bubble::OverlappingScvCorners<Dune::GeometryTypes::tetrahedron>; | ||
209 | 8824 | return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]); | |
210 | } | ||
211 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
1 | else if (type == Dune::GeometryTypes::hexahedron) |
212 | { | ||
213 | using Corners = Detail::PQ1Bubble::OverlappingScvCorners<Dune::GeometryTypes::hexahedron>; | ||
214 | 2 | return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]); | |
215 | } | ||
216 | else | ||
217 | ✗ | DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv geometries for dim=" << dim | |
218 | << " dimWorld=" << dimWorld | ||
219 | << " type=" << type); | ||
220 | } | ||
221 | |||
222 | 547821 | Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const | |
223 | { | ||
224 | // proceed according to number of corners of the element | ||
225 | 547821 | const auto type = geo_.type(); | |
226 | 547821 | const auto numBoxScv = boxHelper_.numScv(); | |
227 |
2/2✓ Branch 0 taken 424560 times.
✓ Branch 1 taken 123261 times.
|
547821 | if (localScvIdx < numBoxScv) |
228 | 424560 | return Dune::GeometryTypes::cube(dim); | |
229 |
1/2✓ Branch 0 taken 123261 times.
✗ Branch 1 not taken.
|
123261 | else if (type == Dune::GeometryTypes::simplex(dim)) |
230 | 72908 | return Dune::GeometryTypes::simplex(dim); | |
231 |
1/2✓ Branch 0 taken 50353 times.
✗ Branch 1 not taken.
|
50353 | else if (type == Dune::GeometryTypes::quadrilateral) |
232 | 50352 | return Dune::GeometryTypes::quadrilateral; | |
233 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
1 | else if (type == Dune::GeometryTypes::hexahedron) |
234 | 1 | return Dune::GeometryTypes::none(dim); // octahedron | |
235 | else | ||
236 | ✗ | DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv geometries for dim=" << dim | |
237 | << " dimWorld=" << dimWorld | ||
238 | << " type=" << type); | ||
239 | ✗ | } | |
240 | |||
241 | //! Create a vector with the corners of sub control volume faces | ||
242 | 857952 | ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const | |
243 | { | ||
244 | // proceed according to number of corners | ||
245 | 857952 | const auto type = geo_.type(); | |
246 | 857952 | const auto numBoxScvf = boxHelper_.numInteriorScvf(); | |
247 | // reuse box geometry helper for the corner scvs | ||
248 |
2/2✓ Branch 0 taken 433392 times.
✓ Branch 1 taken 424560 times.
|
857952 | if (localScvfIdx < numBoxScvf) |
249 | 433392 | return boxHelper_.getScvfCorners(localScvfIdx); | |
250 | |||
251 | 424560 | const auto localOverlappingScvfIdx = localScvfIdx-numBoxScvf; | |
252 |
1/2✓ Branch 0 taken 424560 times.
✗ Branch 1 not taken.
|
424560 | if (type == Dune::GeometryTypes::triangle) |
253 | { | ||
254 | using Corners = Detail::PQ1Bubble::OverlappingScvfCorners<Dune::GeometryTypes::triangle>; | ||
255 | 410976 | return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]); | |
256 | } | ||
257 |
1/2✓ Branch 0 taken 219072 times.
✗ Branch 1 not taken.
|
219072 | else if (type == Dune::GeometryTypes::quadrilateral) |
258 | { | ||
259 | using Corners = Detail::PQ1Bubble::OverlappingScvfCorners<Dune::GeometryTypes::quadrilateral>; | ||
260 | 402816 | return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]); | |
261 | } | ||
262 |
1/2✓ Branch 0 taken 17664 times.
✗ Branch 1 not taken.
|
17664 | else if (type == Dune::GeometryTypes::tetrahedron) |
263 | { | ||
264 | using Corners = Detail::PQ1Bubble::OverlappingScvfCorners<Dune::GeometryTypes::tetrahedron>; | ||
265 | 35296 | return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]); | |
266 | } | ||
267 |
1/2✓ Branch 0 taken 16 times.
✗ Branch 1 not taken.
|
16 | else if (type == Dune::GeometryTypes::hexahedron) |
268 | { | ||
269 | using Corners = Detail::PQ1Bubble::OverlappingScvfCorners<Dune::GeometryTypes::hexahedron>; | ||
270 | 32 | return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]); | |
271 | } | ||
272 | else | ||
273 | ✗ | DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scvf geometries for dim=" << dim | |
274 | << " dimWorld=" << dimWorld | ||
275 | << " type=" << type); | ||
276 | } | ||
277 | |||
278 | 613912 | Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const | |
279 | { | ||
280 |
1/2✓ Branch 2 taken 75200 times.
✗ Branch 3 not taken.
|
857952 | const auto numBoxScvf = boxHelper_.numInteriorScvf(); |
281 |
6/6✓ Branch 0 taken 311368 times.
✓ Branch 1 taken 302544 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 12 times.
✓ Branch 4 taken 8 times.
✓ Branch 5 taken 12 times.
|
613952 | if (localScvfIdx < numBoxScvf) |
282 | 311368 | return Dune::GeometryTypes::cube(dim-1); | |
283 | else | ||
284 | 302544 | return Dune::GeometryTypes::simplex(dim-1); | |
285 | } | ||
286 | |||
287 | //! Create the sub control volume face geometries on the boundary | ||
288 | ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, | ||
289 | unsigned int indexInFacet) const | ||
290 | { | ||
291 |
1/2✓ Branch 1 taken 14246 times.
✗ Branch 2 not taken.
|
18094 | return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet); |
292 | } | ||
293 | |||
294 | template<int d = dimWorld, std::enable_if_t<(d==3), int> = 0> | ||
295 | 44140 | GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair) | |
296 | { | ||
297 | 308980 | auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]); | |
298 | 88280 | normal /= normal.two_norm(); | |
299 | |||
300 | 88280 | GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]); | |
301 | |||
302 | 44140 | const auto s = v*normal; | |
303 |
4/4✓ Branch 0 taken 8828 times.
✓ Branch 1 taken 35312 times.
✓ Branch 2 taken 8828 times.
✓ Branch 3 taken 35312 times.
|
88280 | if (std::signbit(s)) |
304 | normal *= -1; | ||
305 | |||
306 | 44140 | return normal; | |
307 | } | ||
308 | |||
309 | template<int d = dimWorld, std::enable_if_t<(d==2), int> = 0> | ||
310 | 813792 | GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair) | |
311 | { | ||
312 | //! obtain normal vector by 90° counter-clockwise rotation of t | ||
313 | 2441376 | const auto t = p[1] - p[0]; | |
314 | 2441376 | GlobalPosition normal({-t[1], t[0]}); | |
315 | 1627584 | normal /= normal.two_norm(); | |
316 | |||
317 | 1627584 | GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]); | |
318 | |||
319 | 813792 | const auto s = v*normal; | |
320 |
4/4✓ Branch 0 taken 261152 times.
✓ Branch 1 taken 552640 times.
✓ Branch 2 taken 261152 times.
✓ Branch 3 taken 552640 times.
|
1627584 | if (std::signbit(s)) |
321 | normal *= -1; | ||
322 | |||
323 | 813792 | return normal; | |
324 | } | ||
325 | |||
326 | //! the wrapped element geometry | ||
327 | const typename Element::Geometry& elementGeometry() const | ||
328 | { return geo_; } | ||
329 | |||
330 | //! number of interior sub control volume faces | ||
331 | 885032 | std::size_t numInteriorScvf() const | |
332 | { | ||
333 | 885032 | return boxHelper_.numInteriorScvf() + referenceElement(geo_).size(dim); | |
334 | } | ||
335 | |||
336 | //! number of boundary sub control volume faces for face localFacetIndex | ||
337 | ✗ | std::size_t numBoundaryScvf(unsigned int localFacetIndex) const | |
338 | { | ||
339 | ✗ | return referenceElement(geo_).size(localFacetIndex, 1, dim); | |
340 | } | ||
341 | |||
342 | //! number of sub control volumes (number of codim-1 entities) | ||
343 | ✗ | std::size_t numScv() const | |
344 | { | ||
345 |
9/13✓ Branch 2 taken 12000 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 404704 times.
✓ Branch 6 taken 195669 times.
✓ Branch 7 taken 21101 times.
✓ Branch 8 taken 12000 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 12000 times.
✗ Branch 12 not taken.
✓ Branch 14 taken 61600 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 49600 times.
✓ Branch 17 taken 12000 times.
|
671074 | return boxHelper_.numScv() + 1; |
346 | } | ||
347 | |||
348 | //! get scv volume | ||
349 | 22069 | Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage& p) const | |
350 | { | ||
351 |
1/2✓ Branch 1 taken 49600 times.
✗ Branch 2 not taken.
|
547813 | const auto scvType = getScvGeometryType(localScvIdx); |
352 | if constexpr (dim == 3) | ||
353 |
2/2✓ Branch 0 taken 1 times.
✓ Branch 1 taken 22068 times.
|
22069 | if (scvType == Dune::GeometryTypes::none(dim)) |
354 | 1 | return octahedronVolume_(p); | |
355 | |||
356 | 547812 | return Dumux::convexPolytopeVolume<dim>( | |
357 | scvType, | ||
358 | 4753264 | [&](unsigned int i){ return p[i]; } | |
359 | ); | ||
360 | } | ||
361 | |||
362 | template<class DofMapper> | ||
363 | 392704 | auto dofIndex(const DofMapper& dofMapper, const Element& element, unsigned int localScvIdx) const | |
364 | { | ||
365 |
2/2✓ Branch 1 taken 302544 times.
✓ Branch 2 taken 90160 times.
|
392704 | if (localScvIdx < numScv()-1) |
366 | 302544 | return dofMapper.subIndex(element, localScvIdx, dim); | |
367 | else | ||
368 | 90160 | return dofMapper.index(element); | |
369 | } | ||
370 | |||
371 | 2063677 | GlobalPosition dofPosition(unsigned int localScvIdx) const | |
372 | { | ||
373 |
2/2✓ Branch 1 taken 1565464 times.
✓ Branch 2 taken 498213 times.
|
2063677 | if (localScvIdx < numScv()-1) |
374 | 1565464 | return geo_.corner(localScvIdx); | |
375 | else | ||
376 | 498213 | return geo_.center(); | |
377 | } | ||
378 | |||
379 | 857932 | std::array<LocalIndexType, 2> getScvPairForScvf(unsigned int localScvfIndex) const | |
380 | { | ||
381 |
2/2✓ Branch 1 taken 433380 times.
✓ Branch 2 taken 424552 times.
|
857932 | const auto numEdges = referenceElement(geo_).size(dim-1); |
382 |
2/2✓ Branch 0 taken 433380 times.
✓ Branch 1 taken 424552 times.
|
857932 | if (localScvfIndex < numEdges) |
383 | return { | ||
384 | 433380 | static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)), | |
385 | 433380 | static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim)) | |
386 | 866760 | }; | |
387 | else | ||
388 | return { | ||
389 | 424552 | static_cast<LocalIndexType>(numScv()-1), | |
390 | 424552 | static_cast<LocalIndexType>(localScvfIndex-numEdges) | |
391 | 424552 | }; | |
392 | } | ||
393 | |||
394 | ✗ | std::array<LocalIndexType, 2> getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const | |
395 | { | ||
396 | ✗ | const LocalIndexType insideScvIdx | |
397 | ✗ | = static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim)); | |
398 | ✗ | return { insideScvIdx, insideScvIdx }; | |
399 | } | ||
400 | |||
401 | ✗ | bool isOverlappingScvf(unsigned int localScvfIndex) const | |
402 | { | ||
403 |
6/6✓ Branch 1 taken 302544 times.
✓ Branch 2 taken 311368 times.
✓ Branch 4 taken 159608 times.
✓ Branch 5 taken 84412 times.
✓ Branch 6 taken 37600 times.
✓ Branch 7 taken 37600 times.
|
857932 | if (localScvfIndex < boxHelper_.numInteriorScvf()) |
404 | ✗ | return false; | |
405 | else | ||
406 | ✗ | return true; | |
407 | } | ||
408 | |||
409 | ✗ | bool isOverlappingScv(unsigned int localScvIndex) const | |
410 | { | ||
411 |
6/6✓ Branch 1 taken 90160 times.
✓ Branch 2 taken 302544 times.
✓ Branch 4 taken 70701 times.
✓ Branch 5 taken 84408 times.
✓ Branch 6 taken 12000 times.
✓ Branch 7 taken 37600 times.
|
547813 | if (localScvIndex < boxHelper_.numScv()) |
412 | ✗ | return false; | |
413 | else | ||
414 | ✗ | return true; | |
415 | } | ||
416 | |||
417 | private: | ||
418 | 1 | Scalar octahedronVolume_(const ScvCornerStorage& p) const | |
419 | { | ||
420 | using std::abs; | ||
421 | return 1.0/6.0 * ( | ||
422 | 11 | abs(Dumux::tripleProduct(p[4]-p[0], p[1]-p[0], p[2]-p[0])) | |
423 | 10 | + abs(Dumux::tripleProduct(p[5]-p[0], p[1]-p[0], p[2]-p[0])) | |
424 | 2 | ); | |
425 | } | ||
426 | |||
427 | const typename Element::Geometry& geo_; //!< Reference to the element geometry | ||
428 | Dumux::BoxGeometryHelper<GridView, dim, ScvType, ScvfType> boxHelper_; | ||
429 | }; | ||
430 | |||
431 | } // end namespace Dumux | ||
432 | |||
433 | #endif | ||
434 |