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-FileCopyrightText: 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 FaceCenteredStaggeredDiscretization | ||
10 | * \copydoc Dumux::FaceCenteredStaggeredFVElementGeometry | ||
11 | */ | ||
12 | #ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_ELEMENT_GEOMETRY_HH | ||
13 | #define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_ELEMENT_GEOMETRY_HH | ||
14 | |||
15 | #include <utility> | ||
16 | #include <bitset> | ||
17 | |||
18 | #include <dune/common/rangeutilities.hh> | ||
19 | #include <dune/common/reservedvector.hh> | ||
20 | #include <dune/common/iteratorrange.hh> | ||
21 | #include <dune/common/exceptions.hh> | ||
22 | |||
23 | #include <dumux/common/indextraits.hh> | ||
24 | #include <dumux/common/parameters.hh> | ||
25 | #include <dumux/discretization/scvandscvfiterators.hh> | ||
26 | #include <dumux/discretization/facecentered/staggered/normalaxis.hh> | ||
27 | #include <dumux/discretization/facecentered/staggered/consistentlyorientedgrid.hh> | ||
28 | |||
29 | namespace Dumux { | ||
30 | |||
31 | #ifndef DOXYGEN | ||
32 | namespace Detail::FCStaggered { | ||
33 | |||
34 | template<class FVElementGeometry, class SubControlVolume> | ||
35 | 112 | typename SubControlVolume::Traits::Geometry scvGeometry(const FVElementGeometry& fvGeometry, | |
36 | const SubControlVolume& scv) | ||
37 | { | ||
38 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 112 times.
|
112 | typename SubControlVolume::Traits::CornerStorage corners{}; |
39 | // select the containing element | ||
40 | 112 | const auto elementGeometry = (scv.elementIndex() != fvGeometry.elementIndex()) ? | |
41 | ✗ | fvGeometry.element().geometry() : | |
42 | 112 | fvGeometry.gridGeometry().element(scv.elementIndex()).geometry(); | |
43 | |||
44 | 112 | const auto center = elementGeometry.center(); | |
45 | 112 | const auto dofAxis = scv.dofAxis(); | |
46 |
2/2✓ Branch 0 taken 448 times.
✓ Branch 1 taken 112 times.
|
560 | for (int i = 0; i < corners.size(); ++i) |
47 | { | ||
48 | 448 | auto& corner = corners[i]; | |
49 | |||
50 | // copy the corner of the corresponding element | ||
51 | 448 | corner = elementGeometry.corner(i); | |
52 | |||
53 | // shift the corner such that the scv covers half of the element | ||
54 | // (keep the corner positions at the face with the staggered dof) | ||
55 |
2/2✓ Branch 0 taken 224 times.
✓ Branch 1 taken 224 times.
|
448 | if ((corner[dofAxis] - center[dofAxis]) * scv.directionSign() < 0.0) |
56 | 224 | corner[dofAxis] = center[dofAxis]; | |
57 | } | ||
58 | |||
59 | 112 | return {corners.front(), corners.back()}; | |
60 | } | ||
61 | |||
62 | template<class FVElementGeometry, class SubControlVolumeFace> | ||
63 | 9470 | typename SubControlVolumeFace::Traits::Geometry scvfGeometry(const FVElementGeometry& fvGeometry, | |
64 | const SubControlVolumeFace& scvf) | ||
65 | { | ||
66 | 9470 | const auto normalAxis = scvf.normalAxis(); | |
67 | 9470 | const auto center = scvf.center(); | |
68 |
2/2✓ Branch 0 taken 9070 times.
✓ Branch 1 taken 400 times.
|
18940 | const auto shift = scvf.ipGlobal() - center; |
69 |
3/4✓ Branch 0 taken 9070 times.
✓ Branch 1 taken 400 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 9470 times.
|
9470 | const auto dofAxis = scvf.isLateral() ? Dumux::normalAxis(shift) : normalAxis; |
70 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 9470 times.
|
9470 | const auto insideElementIndex = (fvGeometry.scv(scvf.insideScvIdx())).elementIndex(); |
71 | 9470 | const auto elementGeometry = (insideElementIndex != fvGeometry.elementIndex()) ? | |
72 | ✗ | fvGeometry.element().geometry() : | |
73 | 9470 | fvGeometry.gridGeometry().element(insideElementIndex).geometry(); | |
74 | |||
75 | auto corners = std::array{ | ||
76 | 9470 | elementGeometry.corner(0), | |
77 | 9470 | elementGeometry.corner(elementGeometry.corners() - 1) | |
78 | }; | ||
79 | |||
80 | // shift corners to scvf plane and halve lateral faces | ||
81 |
2/2✓ Branch 0 taken 18940 times.
✓ Branch 1 taken 9470 times.
|
28410 | for (int i = 0; i < corners.size(); ++i) |
82 | { | ||
83 |
2/2✓ Branch 0 taken 18140 times.
✓ Branch 1 taken 800 times.
|
18940 | auto& corner = corners[i]; |
84 |
2/2✓ Branch 0 taken 18140 times.
✓ Branch 1 taken 800 times.
|
18940 | corner[normalAxis] = center[normalAxis]; |
85 |
4/4✓ Branch 0 taken 18140 times.
✓ Branch 1 taken 800 times.
✓ Branch 2 taken 9070 times.
✓ Branch 3 taken 9070 times.
|
55220 | if (scvf.isLateral() && (corner - center)*shift < 0.0) |
86 | 9070 | corner[dofAxis] = elementGeometry.center()[dofAxis]; | |
87 | } | ||
88 | |||
89 | 9470 | auto inPlaneAxes = std::move(std::bitset<SubControlVolumeFace::Traits::dimWorld>{}.set()); | |
90 | 9470 | inPlaneAxes.set(normalAxis, false); | |
91 | |||
92 | 9470 | return {corners[0], corners[1], inPlaneAxes}; | |
93 | } | ||
94 | |||
95 | //! Get the scv on the outside side of a periodic boundary | ||
96 | template<class FVElementGeometry, class SubControlVolume> | ||
97 | 72 | const SubControlVolume& outsidePeriodicScv(const FVElementGeometry& fvGeometry, | |
98 | const SubControlVolume& selfScv) | ||
99 | { | ||
100 |
2/2✓ Branch 1 taken 36 times.
✓ Branch 2 taken 36 times.
|
144 | assert(fvGeometry.gridGeometry().dofOnPeriodicBoundary(selfScv.dofIndex())); |
101 | |||
102 |
2/2✓ Branch 0 taken 36 times.
✓ Branch 1 taken 36 times.
|
72 | auto localOppositeIndex = FVElementGeometry::GridGeometry::GeometryHelper::localOppositeIdx(selfScv.localDofIndex()); |
103 | 144 | const auto& normalScvf = std::next(scvfs(fvGeometry, selfScv).begin()); | |
104 | |||
105 | // at a lateral velocity, find the inner and outer normal velocities | ||
106 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 72 times.
|
72 | const auto& orthogonalScvf = fvGeometry.lateralOrthogonalScvf(*normalScvf); |
107 | 144 | const auto orthogonalOutsideScv = fvGeometry.scv(orthogonalScvf.outsideScvIdx()); | |
108 | |||
109 | 72 | auto outsidePeriodicFVGeometry = localView(fvGeometry.gridGeometry()); | |
110 | 72 | const auto& periodicElement = fvGeometry.gridGeometry().element(orthogonalOutsideScv.elementIndex()); | |
111 | 72 | outsidePeriodicFVGeometry.bindElement(periodicElement); | |
112 | |||
113 |
3/4✓ Branch 1 taken 72 times.
✓ Branch 2 taken 168 times.
✓ Branch 3 taken 240 times.
✗ Branch 4 not taken.
|
240 | for (const auto& outerPeriodicScv : scvs(outsidePeriodicFVGeometry)) |
114 |
2/2✓ Branch 0 taken 72 times.
✓ Branch 1 taken 168 times.
|
240 | if (outerPeriodicScv.localDofIndex() == localOppositeIndex) |
115 | 72 | return outerPeriodicScv; | |
116 | |||
117 | ✗ | DUNE_THROW(Dune::InvalidStateException, "No outside periodic scv found"); | |
118 | } | ||
119 | |||
120 | } // end namespace Detail::FCStaggered | ||
121 | #endif // DOXYGEN | ||
122 | |||
123 | template<class GG, bool cachingEnabled> | ||
124 | class FaceCenteredStaggeredFVElementGeometry; | ||
125 | |||
126 | /*! | ||
127 | * \ingroup FaceCenteredStaggeredDiscretization | ||
128 | * \brief Stencil-local finite volume geometry (scvs and scvfs) for face-centered staggered models | ||
129 | * Specialization for grid caching enabled | ||
130 | */ | ||
131 | template<class GG> | ||
132 | class FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/true> | ||
133 | { | ||
134 | using ThisType = FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/true>; | ||
135 | using GridView = typename GG::GridView; | ||
136 | using GridIndexType = typename IndexTraits<GridView>::GridIndex; | ||
137 | static constexpr auto numScvsPerElement = GG::StaticInformation::numScvsPerElement; | ||
138 | |||
139 | public: | ||
140 | //! export type of subcontrol volume face | ||
141 | using SubControlVolume = typename GG::SubControlVolume; | ||
142 | using SubControlVolumeFace = typename GG::SubControlVolumeFace; | ||
143 | using Element = typename GridView::template Codim<0>::Entity; | ||
144 | using GridGeometry = GG; | ||
145 | |||
146 | //! the maximum number of scvs per element | ||
147 | static constexpr std::size_t maxNumElementScvs = numScvsPerElement; | ||
148 | |||
149 | 74445102 | FaceCenteredStaggeredFVElementGeometry(const GridGeometry& gridGeometry) | |
150 |
23/26✓ Branch 2 taken 6768 times.
✓ Branch 3 taken 8 times.
✓ Branch 5 taken 45 times.
✓ Branch 6 taken 8 times.
✓ Branch 9 taken 2168341 times.
✓ Branch 10 taken 7215 times.
✓ Branch 12 taken 1384 times.
✓ Branch 13 taken 72189 times.
✓ Branch 18 taken 1777 times.
✓ Branch 19 taken 25 times.
✓ Branch 4 taken 18 times.
✓ Branch 7 taken 15 times.
✗ Branch 11 not taken.
✓ Branch 15 taken 17 times.
✓ Branch 16 taken 13 times.
✓ Branch 14 taken 386401 times.
✓ Branch 20 taken 3 times.
✓ Branch 21 taken 5312 times.
✓ Branch 17 taken 14 times.
✓ Branch 22 taken 19 times.
✗ Branch 23 not taken.
✓ Branch 1 taken 3 times.
✓ Branch 8 taken 88544 times.
✓ Branch 24 taken 19 times.
✓ Branch 25 taken 1 times.
✗ Branch 26 not taken.
|
74445102 | : gridGeometry_(&gridGeometry) |
151 | {} | ||
152 | |||
153 | //! Get a sub control volume with a global scv index | ||
154 | 105690480 | const SubControlVolume& scv(GridIndexType scvIdx) const | |
155 |
55/58✓ Branch 0 taken 2477867 times.
✓ Branch 1 taken 229719267 times.
✓ Branch 3 taken 3316797 times.
✓ Branch 4 taken 77705433 times.
✓ Branch 6 taken 6092028 times.
✓ Branch 7 taken 196176155 times.
✓ Branch 8 taken 10326653 times.
✓ Branch 9 taken 151273450 times.
✓ Branch 10 taken 10777991 times.
✓ Branch 11 taken 39380015 times.
✓ Branch 12 taken 447424 times.
✓ Branch 13 taken 158515514 times.
✓ Branch 15 taken 6796840 times.
✓ Branch 16 taken 106945690 times.
✓ Branch 18 taken 806166 times.
✓ Branch 19 taken 93003622 times.
✓ Branch 21 taken 25148064 times.
✓ Branch 22 taken 94530190 times.
✓ Branch 23 taken 1397092 times.
✓ Branch 24 taken 4533852 times.
✓ Branch 26 taken 41975104 times.
✓ Branch 27 taken 6627988 times.
✓ Branch 28 taken 52588894 times.
✓ Branch 29 taken 20925344 times.
✓ Branch 30 taken 3541332 times.
✓ Branch 31 taken 1496832 times.
✓ Branch 32 taken 5789288 times.
✓ Branch 33 taken 6455270 times.
✓ Branch 36 taken 3641088 times.
✓ Branch 37 taken 2040074 times.
✓ Branch 38 taken 2523266 times.
✓ Branch 39 taken 2719694 times.
✓ Branch 40 taken 1091196 times.
✓ Branch 41 taken 323346 times.
✓ Branch 2 taken 1411904 times.
✗ Branch 5 not taken.
✓ Branch 14 taken 94208 times.
✓ Branch 20 taken 229248 times.
✓ Branch 34 taken 43907694 times.
✓ Branch 35 taken 2168452 times.
✓ Branch 42 taken 127926 times.
✓ Branch 17 taken 14784110 times.
✓ Branch 25 taken 46736502 times.
✓ Branch 43 taken 82332 times.
✓ Branch 44 taken 9648 times.
✓ Branch 45 taken 19224 times.
✓ Branch 46 taken 3745280 times.
✓ Branch 50 taken 32640 times.
✓ Branch 51 taken 24 times.
✓ Branch 52 taken 5040 times.
✓ Branch 53 taken 9600 times.
✓ Branch 54 taken 1120 times.
✓ Branch 55 taken 6 times.
✓ Branch 58 taken 60 times.
✗ Branch 59 not taken.
✗ Branch 47 not taken.
✓ Branch 48 taken 542080 times.
✓ Branch 49 taken 1075200 times.
|
1974854032 | { return gridGeometry().scv(scvIdx); } |
156 | |||
157 | //! Get a sub control volume face with a global scvf index | ||
158 | 1664658728 | const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const | |
159 |
2/2✓ Branch 0 taken 510243420 times.
✓ Branch 1 taken 109434912 times.
|
619678332 | { return gridGeometry().scvf(scvfIdx); } |
160 | |||
161 | //! Return a the lateral sub control volume face which is orthogonal to the given sub control volume face | ||
162 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 340635518 times.
|
340635518 | const SubControlVolumeFace& lateralOrthogonalScvf(const SubControlVolumeFace& scvf) const |
163 | { | ||
164 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 340635518 times.
|
340635518 | assert(scvf.isLateral()); |
165 |
1/2✗ Branch 2 not taken.
✓ Branch 3 taken 340635518 times.
|
340635518 | const auto otherGlobalIdx = scvfIndices_()[GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex())]; |
166 | 340635518 | return gridGeometry().scvf(otherGlobalIdx); | |
167 | |||
168 | } | ||
169 | |||
170 | //! Return the frontal sub control volume face on a the boundary for a given sub control volume | ||
171 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 8392389 times.
|
8392389 | const SubControlVolumeFace& frontalScvfOnBoundary(const SubControlVolume& scv) const |
172 | { | ||
173 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 8392389 times.
|
8392389 | assert(scv.boundary()); |
174 | |||
175 | // frontal boundary faces are always stored after the lateral faces | ||
176 | 8392389 | auto scvfIter = scvfs(*this, scv).begin(); | |
177 | 8392389 | const auto end = scvfs(*this, scv).end(); | |
178 |
5/6✓ Branch 1 taken 16784778 times.
✓ Branch 2 taken 18594290 times.
✓ Branch 3 taken 8392389 times.
✓ Branch 4 taken 8392389 times.
✓ Branch 5 taken 26986679 times.
✗ Branch 6 not taken.
|
35379068 | while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end)) |
179 | 26986679 | ++scvfIter; | |
180 | |||
181 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 8392389 times.
|
8392389 | assert(scvfIter->isFrontal()); |
182 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 8392389 times.
|
8392389 | assert(scvfIter->boundary()); |
183 | 8392389 | return *scvfIter; | |
184 | } | ||
185 | |||
186 | //! iterator range for sub control volumes. Iterates over | ||
187 | //! all scvs of the bound element (not including neighbor scvs) | ||
188 | //! This is a free function found by means of ADL | ||
189 | //! To iterate over all sub control volumes of this FVElementGeometry use | ||
190 | //! for (auto&& scv : scvs(fvGeometry)) | ||
191 | friend inline auto | ||
192 | 86846386 | scvs(const FaceCenteredStaggeredFVElementGeometry& fvGeometry) | |
193 | 86846386 | { return fvGeometry.gridGeometry().scvs(fvGeometry); } | |
194 | |||
195 | //! iterator range for sub control volumes faces. Iterates over | ||
196 | //! all scvfs of the bound element. | ||
197 | //! This is a free function found by means of ADL | ||
198 | //! To iterate over all sub control volume faces of this FVElementGeometry use | ||
199 | //! for (auto&& scvf : scvfs(fvGeometry)) | ||
200 | friend inline auto | ||
201 | 4431248 | scvfs(const FaceCenteredStaggeredFVElementGeometry& fvGeometry) | |
202 | { | ||
203 | using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>; | ||
204 | using ScvfIterator = Dumux::ScvfIterator<SubControlVolumeFace, IndexContainerType, ThisType>; | ||
205 | 4431248 | return Dune::IteratorRange<ScvfIterator>(ScvfIterator(fvGeometry.scvfIndices_().begin(), fvGeometry), | |
206 | 4431248 | ScvfIterator(fvGeometry.scvfIndices_().end(), fvGeometry)); | |
207 | } | ||
208 | |||
209 | //! iterator range for sub control volumes faces. Iterates over | ||
210 | //! all scvfs of the bound element belonging to the given sub control volume. | ||
211 | //! This is a free function found by means of ADL | ||
212 | //! To iterate over all sub control volume faces of this FVElementGeometry use | ||
213 | //! for (auto&& scvf : scvfs(fvGeometry, scv)) | ||
214 | friend inline auto | ||
215 | 109434912 | scvfs(const FaceCenteredStaggeredFVElementGeometry& fvGeometry, const SubControlVolume& scv) | |
216 | { | ||
217 | using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>; | ||
218 | using ScvfIterator = Dumux::SkippingScvfIterator<SubControlVolumeFace, IndexContainerType, ThisType>; | ||
219 | 109434912 | auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index()); | |
220 | 109434912 | auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index()); | |
221 | 109434912 | return Dune::IteratorRange<ScvfIterator>(begin, end); | |
222 | } | ||
223 | |||
224 | //! number of sub control volumes in this fv element geometry | ||
225 | std::size_t numScv() const | ||
226 | { return numScvsPerElement; } | ||
227 | |||
228 | //! number of sub control volumes in this fv element geometry | ||
229 | 546134 | std::size_t numScvf() const | |
230 | 367816 | { return scvfIndices_().size(); } | |
231 | |||
232 | //! Returns whether one of the geometry's scvfs lies on a boundary | ||
233 | 5839695 | bool hasBoundaryScvf() const | |
234 | 5839695 | { return gridGeometry().hasBoundaryScvf(eIdx_); } | |
235 | |||
236 | /*! | ||
237 | * \brief bind the local view (r-value overload) | ||
238 | * This overload is called when an instance of this class is a temporary in the usage context | ||
239 | * This allows a usage like this: `const auto view = localView(...).bind(element);` | ||
240 | */ | ||
241 | 368865 | FaceCenteredStaggeredFVElementGeometry bind(const Element& element) && | |
242 | { | ||
243 | 368865 | this->bind(element); | |
244 |
2/4✓ Branch 1 taken 356714 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 12151 times.
✗ Branch 5 not taken.
|
368865 | return std::move(*this); |
245 | } | ||
246 | |||
247 | //! Binding of an element, called by the local jacobian to prepare element assembly | ||
248 | 4349278 | void bind(const Element& element) & | |
249 |
4/8✓ Branch 2 taken 3703 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 1768 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 2711 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 1768 times.
✗ Branch 12 not taken.
|
4349278 | { this->bindElement(element); } |
250 | |||
251 | /*! | ||
252 | * \brief bind the local view (r-value overload) | ||
253 | * This overload is called when an instance of this class is a temporary in the usage context | ||
254 | * This allows a usage like this: `const auto view = localView(...).bind(element);` | ||
255 | */ | ||
256 | 71202228 | FaceCenteredStaggeredFVElementGeometry bindElement(const Element& element) && | |
257 | { | ||
258 |
1/2✓ Branch 2 taken 5304 times.
✗ Branch 3 not taken.
|
71202228 | this->bindElement(element); |
259 |
2/3✓ Branch 1 taken 13360 times.
✓ Branch 2 taken 5304 times.
✗ Branch 3 not taken.
|
71202228 | return std::move(*this); |
260 | } | ||
261 | |||
262 | //! Bind only element-local | ||
263 | 104043931 | void bindElement(const Element& element) & | |
264 | { | ||
265 | 104043931 | elementPtr_ = &element; | |
266 | 104043931 | eIdx_ = gridGeometry().elementMapper().index(element); | |
267 | 104043931 | } | |
268 | |||
269 | //! The grid geometry we are a restriction of | ||
270 |
1/2✓ Branch 1 taken 3920 times.
✗ Branch 2 not taken.
|
4377556106 | const GridGeometry& gridGeometry() const |
271 | { | ||
272 |
93/135✗ Branch 0 not taken.
✓ Branch 1 taken 230047994 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 71295925 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 131723206 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 156122432 times.
✗ Branch 12 not taken.
✓ Branch 13 taken 183597871 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 159165552 times.
✗ Branch 18 not taken.
✓ Branch 19 taken 143067778 times.
✗ Branch 21 not taken.
✓ Branch 22 taken 68371212 times.
✗ Branch 24 not taken.
✓ Branch 25 taken 99666054 times.
✗ Branch 27 not taken.
✓ Branch 28 taken 137329992 times.
✗ Branch 30 not taken.
✓ Branch 31 taken 95585159 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 3401928 times.
✗ Branch 36 not taken.
✓ Branch 37 taken 20843570 times.
✓ Branch 38 taken 4666844 times.
✓ Branch 39 taken 1152 times.
✗ Branch 42 not taken.
✓ Branch 43 taken 21337018 times.
✗ Branch 45 not taken.
✓ Branch 46 taken 73675664 times.
✗ Branch 48 not taken.
✓ Branch 49 taken 441549744 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 3577930 times.
✗ Branch 54 not taken.
✓ Branch 55 taken 86466598 times.
✗ Branch 57 not taken.
✓ Branch 58 taken 366722644 times.
✗ Branch 60 not taken.
✓ Branch 61 taken 141059702 times.
✗ Branch 63 not taken.
✓ Branch 64 taken 37660344 times.
✗ Branch 66 not taken.
✓ Branch 67 taken 87460000 times.
✗ Branch 69 not taken.
✓ Branch 70 taken 71891244 times.
✗ Branch 72 not taken.
✓ Branch 73 taken 107030148 times.
✗ Branch 75 not taken.
✓ Branch 76 taken 50597959 times.
✗ Branch 78 not taken.
✓ Branch 79 taken 253755474 times.
✗ Branch 81 not taken.
✓ Branch 82 taken 73939401 times.
✗ Branch 84 not taken.
✓ Branch 85 taken 49339503 times.
✗ Branch 87 not taken.
✓ Branch 88 taken 201751165 times.
✗ Branch 90 not taken.
✓ Branch 91 taken 24861464 times.
✗ Branch 93 not taken.
✓ Branch 94 taken 149678154 times.
✗ Branch 96 not taken.
✓ Branch 97 taken 88318816 times.
✗ Branch 99 not taken.
✓ Branch 100 taken 56621630 times.
✗ Branch 102 not taken.
✓ Branch 103 taken 105152630 times.
✓ Branch 105 taken 329568 times.
✓ Branch 106 taken 13714200 times.
✓ Branch 108 taken 2699542 times.
✓ Branch 109 taken 60092042 times.
✗ Branch 111 not taken.
✓ Branch 112 taken 6091872 times.
✓ Branch 114 taken 450862 times.
✓ Branch 115 taken 21678264 times.
✗ Branch 117 not taken.
✓ Branch 118 taken 4943589 times.
✓ Branch 121 taken 35824871 times.
✓ Branch 122 taken 289780 times.
✓ Branch 40 taken 120925336 times.
✓ Branch 92 taken 47343726 times.
✓ Branch 95 taken 60847096 times.
✓ Branch 20 taken 4567279 times.
✓ Branch 26 taken 3920520 times.
✓ Branch 35 taken 3584878 times.
✓ Branch 98 taken 9322583 times.
✓ Branch 74 taken 233978 times.
✓ Branch 101 taken 7425148 times.
✓ Branch 104 taken 15599464 times.
✓ Branch 5 taken 2933822 times.
✓ Branch 8 taken 420706 times.
✓ Branch 11 taken 383302 times.
✓ Branch 14 taken 2250852 times.
✓ Branch 17 taken 5931796 times.
✓ Branch 23 taken 1503996 times.
✓ Branch 29 taken 3680944 times.
✓ Branch 32 taken 2225976 times.
✓ Branch 41 taken 4914236 times.
✓ Branch 44 taken 530256 times.
✓ Branch 47 taken 2155276 times.
✓ Branch 50 taken 2194316 times.
✓ Branch 53 taken 4273276 times.
✓ Branch 56 taken 24070928 times.
✓ Branch 59 taken 39075648 times.
✓ Branch 62 taken 1146473 times.
✓ Branch 65 taken 15088676 times.
✓ Branch 68 taken 31930396 times.
✓ Branch 71 taken 3451374 times.
✓ Branch 77 taken 1182026 times.
✓ Branch 80 taken 1349998 times.
✓ Branch 83 taken 1593638 times.
✓ Branch 86 taken 2049948 times.
✓ Branch 89 taken 3540042 times.
✓ Branch 107 taken 8418544 times.
✓ Branch 110 taken 676729 times.
✓ Branch 113 taken 223859 times.
✓ Branch 116 taken 1542 times.
✗ Branch 120 not taken.
✗ Branch 123 not taken.
✓ Branch 124 taken 8904480 times.
✗ Branch 126 not taken.
✓ Branch 127 taken 195135 times.
✗ Branch 129 not taken.
✓ Branch 130 taken 13672 times.
✗ Branch 133 not taken.
✓ Branch 134 taken 3756900 times.
✗ Branch 136 not taken.
✓ Branch 137 taken 2439524 times.
✓ Branch 131 taken 575250 times.
✓ Branch 128 taken 668000 times.
✓ Branch 119 taken 15585360 times.
✓ Branch 125 taken 709700 times.
|
4617265095 | assert(gridGeometry_); |
273 | return *gridGeometry_; | ||
274 | } | ||
275 | |||
276 | 686863530 | std::size_t elementIndex() const | |
277 |
6/10✗ Branch 2 not taken.
✓ Branch 3 taken 93628 times.
✓ Branch 4 taken 82974 times.
✓ Branch 5 taken 50 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 93628 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 23798 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 112 times.
|
498514503 | { return eIdx_; } |
278 | |||
279 | //! The bound element | ||
280 | 6951974 | const Element& element() const | |
281 |
6/7✓ Branch 1 taken 669946 times.
✓ Branch 2 taken 83198 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 72998 times.
✓ Branch 7 taken 94744 times.
✓ Branch 8 taken 351480 times.
✓ Branch 0 taken 118568 times.
|
7018742 | { return *elementPtr_; } |
282 | |||
283 | //! Returns true if the IP of an scvf lies on a concave corner | ||
284 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 124241256 times.
|
128855548 | bool scvfIntegrationPointInConcaveCorner(const SubControlVolumeFace& scvf) const |
285 | 128855548 | { return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*this, scvf); } | |
286 | |||
287 | //! Returns the the scvf of neighbor element with the same integration point and unit outer normal | ||
288 | 94928 | const SubControlVolumeFace& outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace& scvf) const | |
289 | { | ||
290 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 94928 times.
|
94928 | const auto& lateralOrthogonalScvf = this->lateralOrthogonalScvf(scvf); |
291 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 94928 times.
|
94928 | assert(!lateralOrthogonalScvf.boundary()); |
292 | |||
293 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 94928 times.
✓ Branch 2 taken 49784 times.
✓ Branch 3 taken 45144 times.
|
94928 | const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx())); |
294 | |||
295 | 94928 | auto outsideFVGeometry = localView(gridGeometry()); | |
296 | 94928 | const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex(); | |
297 |
0/2✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
94928 | outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx)); |
298 | |||
299 |
3/4✓ Branch 0 taken 94928 times.
✓ Branch 1 taken 707472 times.
✓ Branch 2 taken 802400 times.
✗ Branch 3 not taken.
|
897328 | for (const auto& otherScvf : scvfs(outsideFVGeometry)) |
300 |
2/2✓ Branch 0 taken 94928 times.
✓ Branch 1 taken 707472 times.
|
802400 | if (otherScvf.localIndex() == otherLocalIdx) |
301 | 94928 | return otherScvf; | |
302 | |||
303 | ✗ | DUNE_THROW(Dune::InvalidStateException, "No outside scvf found"); | |
304 | } | ||
305 | |||
306 | //! Get the scv on the outside side of a periodic boundary | ||
307 | 72 | const SubControlVolume& outsidePeriodicScv(const SubControlVolume& selfScv) const | |
308 | { | ||
309 |
1/2✓ Branch 1 taken 60 times.
✗ Branch 2 not taken.
|
72 | return Detail::FCStaggered::outsidePeriodicScv(*this, selfScv); |
310 | } | ||
311 | |||
312 | //! Create the geometry of a given sub control volume | ||
313 | 112 | typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const | |
314 | { | ||
315 | 112 | return Detail::FCStaggered::scvGeometry(*this, scv); | |
316 | } | ||
317 | |||
318 | //! Create the geometry of a given sub control volume face | ||
319 | 9470 | typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const | |
320 | { | ||
321 |
1/2✓ Branch 1 taken 9470 times.
✗ Branch 2 not taken.
|
9470 | return Detail::FCStaggered::scvfGeometry(*this, scvf); |
322 | } | ||
323 | |||
324 | private: | ||
325 | |||
326 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 564482724 times.
|
564482724 | const auto& scvfIndices_() const |
327 | { | ||
328 | 564482724 | return gridGeometry().scvfIndicesOfElement(eIdx_); | |
329 | } | ||
330 | |||
331 | const Element* elementPtr_; | ||
332 | GridIndexType eIdx_; | ||
333 | const GridGeometry* gridGeometry_; | ||
334 | }; | ||
335 | |||
336 | /*! | ||
337 | * \ingroup FaceCenteredStaggeredDiscretization | ||
338 | * \brief Stencil-local finite volume geometry (scvs and scvfs) for face-centered staggered models | ||
339 | * Specialization for grid caching disabled | ||
340 | */ | ||
341 | template<class GG> | ||
342 | class FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/false> | ||
343 | { | ||
344 | using ThisType = FaceCenteredStaggeredFVElementGeometry<GG, /*cachingEnabled*/false>; | ||
345 | |||
346 | using GridView = typename GG::GridView; | ||
347 | using GridIndexType = typename IndexTraits<GridView>::GridIndex; | ||
348 | |||
349 | //TODO include assert that checks for quad geometry | ||
350 | |||
351 | using StaticInfo = typename GG::StaticInformation; | ||
352 | static constexpr auto dim = StaticInfo::dim; | ||
353 | static constexpr auto numScvsPerElement = StaticInfo::numScvsPerElement; | ||
354 | static constexpr auto numLateralScvfsPerScv = StaticInfo::numLateralScvfsPerScv; | ||
355 | static constexpr auto numLateralScvfsPerElement = StaticInfo::numLateralScvfsPerElement; | ||
356 | static constexpr auto minNumScvfsPerElement = StaticInfo::minNumScvfsPerElement; | ||
357 | static constexpr auto maxNumScvfsPerElement = StaticInfo::maxNumScvfsPerElement; | ||
358 | |||
359 | using LocalIndexType = typename IndexTraits<GridView>::LocalIndex; | ||
360 | |||
361 | public: | ||
362 | //! export type of subcontrol volume face | ||
363 | using SubControlVolume = typename GG::SubControlVolume; | ||
364 | using SubControlVolumeFace = typename GG::SubControlVolumeFace; | ||
365 | using Element = typename GridView::template Codim<0>::Entity; | ||
366 | using GridGeometry = GG; | ||
367 | |||
368 | //! the maximum number of scvs per element | ||
369 | static constexpr std::size_t maxNumElementScvs = numScvsPerElement; | ||
370 | |||
371 | 1483778 | FaceCenteredStaggeredFVElementGeometry(const GridGeometry& gridGeometry) | |
372 | 1483778 | : gridGeometry_(&gridGeometry) | |
373 | 4451334 | , geometryHelper_(gridGeometry.gridView()) | |
374 | 1483778 | {} | |
375 | |||
376 | //! Get a sub control volume face with a global scvf index | ||
377 | 35742434 | const SubControlVolumeFace& scvf(const GridIndexType scvfIdx) const | |
378 |
2/2✓ Branch 4 taken 11609472 times.
✓ Branch 5 taken 2472200 times.
|
35742434 | { return scvfs_[findLocalIndex_(scvfIdx, scvfIndices_())]; } |
379 | |||
380 | //! iterator range for sub control volumes faces. Iterates over | ||
381 | //! all scvfs of the bound element. | ||
382 | //! This is a free function found by means of ADL | ||
383 | //! To iterate over all sub control volume faces of this FVElementGeometry use | ||
384 | //! for (auto&& scvf : scvfs(fvGeometry)) | ||
385 | friend inline auto | ||
386 | 44552 | scvfs(const FaceCenteredStaggeredFVElementGeometry& fvGeometry) | |
387 | { | ||
388 | using Iter = typename decltype(fvGeometry.scvfs_)::const_iterator; | ||
389 | 44552 | return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end()); | |
390 | } | ||
391 | |||
392 | //! iterator range for sub control volumes faces. Iterates over | ||
393 | //! all scvfs of the bound element belonging to the given sub control volume. | ||
394 | //! This is a free function found by means of ADL | ||
395 | //! To iterate over all sub control volume faces of this FVElementGeometry use | ||
396 | //! for (auto&& scvf : scvfs(fvGeometry, scv)) | ||
397 | friend inline auto | ||
398 | 2472200 | scvfs(const FaceCenteredStaggeredFVElementGeometry& fvGeometry, const SubControlVolume& scv) | |
399 | { | ||
400 | using IndexContainerType = std::decay_t<decltype(fvGeometry.scvfIndices_())>; | ||
401 | using ScvfIterator = Dumux::SkippingScvfIterator<SubControlVolumeFace, IndexContainerType, ThisType>; | ||
402 | 2472200 | auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index()); | |
403 | 2472200 | auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index()); | |
404 | 2472200 | return Dune::IteratorRange<ScvfIterator>(begin, end); | |
405 | } | ||
406 | |||
407 | //! Get a half sub control volume with a global scv index | ||
408 | 67610688 | const SubControlVolume& scv(const GridIndexType scvIdx) const | |
409 | { | ||
410 |
4/4✓ Branch 0 taken 59306326 times.
✓ Branch 1 taken 8304362 times.
✓ Branch 2 taken 51001964 times.
✓ Branch 3 taken 8304362 times.
|
67610688 | if (scvIdx >= scvIndicesOfElement_.front() && scvIdx <= scvIndicesOfElement_.back()) |
411 | 51001964 | return scvs_[findLocalIndex_(scvIdx, scvIndicesOfElement_)]; | |
412 | else | ||
413 | 16608724 | return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)]; | |
414 | } | ||
415 | |||
416 | //! iterator range for sub control volumes. Iterates over | ||
417 | //! all scvs of the bound element (not including neighbor scvs) | ||
418 | //! This is a free function found by means of ADL | ||
419 | //! To iterate over all sub control volumes of this FVElementGeometry use | ||
420 | //! for (auto&& scv : scvs(fvGeometry)) | ||
421 | friend inline auto | ||
422 | 1845740 | scvs(const FaceCenteredStaggeredFVElementGeometry& g) | |
423 | { | ||
424 | using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator; | ||
425 | 1845740 | return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end()); | |
426 | } | ||
427 | |||
428 | //! Return a the lateral sub control volume face which is orthogonal to the given sub control volume face | ||
429 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 7920188 times.
|
7920188 | const SubControlVolumeFace& lateralOrthogonalScvf(const SubControlVolumeFace& scvf) const |
430 | { | ||
431 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 7920188 times.
|
7920188 | assert(scvf.isLateral()); |
432 | 7920188 | const auto otherLocalIdx = GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex()); | |
433 | 7920188 | return scvfs_[otherLocalIdx]; | |
434 | } | ||
435 | |||
436 | //! Return the frontal sub control volume face on a the boundary for a given sub control volume | ||
437 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 245404 times.
|
245404 | const SubControlVolumeFace& frontalScvfOnBoundary(const SubControlVolume& scv) const |
438 | { | ||
439 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 245404 times.
|
245404 | assert(scv.boundary()); |
440 | |||
441 | // frontal boundary faces are always stored after the lateral faces | ||
442 | 245404 | auto scvfIter = scvfs(*this, scv).begin(); | |
443 | 245404 | const auto end = scvfs(*this, scv).end(); | |
444 |
5/6✓ Branch 1 taken 490808 times.
✓ Branch 2 taken 490808 times.
✓ Branch 3 taken 245404 times.
✓ Branch 4 taken 245404 times.
✓ Branch 5 taken 736212 times.
✗ Branch 6 not taken.
|
981616 | while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end)) |
445 | 736212 | ++scvfIter; | |
446 | |||
447 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 245404 times.
|
245404 | assert(scvfIter->isFrontal()); |
448 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 245404 times.
|
245404 | assert(scvfIter->boundary()); |
449 | 245404 | return *scvfIter; | |
450 | } | ||
451 | |||
452 | //! Returns whether one of the geometry's scvfs lies on a boundary | ||
453 | 62072 | bool hasBoundaryScvf() const | |
454 | 62072 | { return gridGeometry().hasBoundaryScvf(eIdx_); } | |
455 | |||
456 | //! number of sub control volumes in this fv element geometry | ||
457 | std::size_t numScv() const | ||
458 | { return numScvsPerElement; } | ||
459 | |||
460 | //! number of sub control volumes in this fv element geometry | ||
461 | 34704 | std::size_t numScvf() const | |
462 | 34704 | { return scvfIndices_().size(); } | |
463 | |||
464 | /*! | ||
465 | * \brief bind the local view (r-value overload) | ||
466 | * This overload is called when an instance of this class is a temporary in the usage context | ||
467 | * This allows a usage like this: `const auto view = localView(...).bind(element);` | ||
468 | */ | ||
469 | FaceCenteredStaggeredFVElementGeometry bind(const Element& element) && | ||
470 | { | ||
471 | this->bind_(element); | ||
472 | return std::move(*this); | ||
473 | } | ||
474 | |||
475 | 40090 | void bind(const Element& element) & | |
476 |
3/6✓ Branch 2 taken 3018 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 1768 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 1768 times.
✗ Branch 9 not taken.
|
40090 | { this->bind_(element); } |
477 | |||
478 | /*! | ||
479 | * \brief bind the local view (r-value overload) | ||
480 | * This overload is called when an instance of this class is a temporary in the usage context | ||
481 | * This allows a usage like this: `const auto view = localView(...).bind(element);` | ||
482 | */ | ||
483 | 1452712 | FaceCenteredStaggeredFVElementGeometry bindElement(const Element& element) && | |
484 | { | ||
485 | 152512 | typename GG::LocalIntersectionMapper localIsMapper; | |
486 | 1452712 | localIsMapper.update(gridGeometry().gridView(), element); | |
487 | 1452712 | this->bindElement_(element, localIsMapper); | |
488 | 1452712 | return std::move(*this); | |
489 | } | ||
490 | |||
491 | 493580 | void bindElement(const Element& element) & | |
492 | { | ||
493 | 48480 | typename GG::LocalIntersectionMapper localIsMapper; | |
494 | 493580 | localIsMapper.update(gridGeometry().gridView(), element); | |
495 | 493580 | this->bindElement_(element, localIsMapper); | |
496 | 493580 | } | |
497 | |||
498 | 13839544 | GridIndexType elementIndex() const | |
499 | 11786976 | { return eIdx_; } | |
500 | |||
501 | //! The bound element | ||
502 | 40656 | const Element& element() const | |
503 |
2/2✓ Branch 0 taken 8624 times.
✓ Branch 1 taken 32032 times.
|
40656 | { return *elementPtr_; } |
504 | |||
505 | //! The grid geometry we are a restriction of | ||
506 | 97609084 | const GridGeometry& gridGeometry() const | |
507 | { | ||
508 |
17/31✓ Branch 1 taken 1300200 times.
✓ Branch 2 taken 152512 times.
✓ Branch 4 taken 55000 times.
✓ Branch 5 taken 7072 times.
✓ Branch 7 taken 27500 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 27500 times.
✓ Branch 11 taken 3536 times.
✓ Branch 13 taken 110000 times.
✓ Branch 14 taken 14144 times.
✓ Branch 16 taken 866800 times.
✓ Branch 17 taken 111104 times.
✓ Branch 19 taken 866800 times.
✓ Branch 20 taken 111104 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 34080 times.
✗ Branch 36 not taken.
✓ Branch 37 taken 68160 times.
✗ Branch 39 not taken.
✓ Branch 40 taken 48480 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 0 not taken.
✗ Branch 3 not taken.
✗ Branch 6 not taken.
✗ Branch 9 not taken.
✗ Branch 12 not taken.
✗ Branch 15 not taken.
✗ Branch 18 not taken.
✗ Branch 29 not taken.
✓ Branch 30 taken 445100 times.
|
4249092 | assert(gridGeometry_); |
509 | return *gridGeometry_; | ||
510 | } | ||
511 | |||
512 | //! Returns true if the IP of an scvf lies on a concave corner | ||
513 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 3699872 times.
|
3699872 | bool scvfIntegrationPointInConcaveCorner(const SubControlVolumeFace& scvf) const |
514 | 3699872 | { return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*this, scvf); } | |
515 | |||
516 | //! Returns the the scvf of neighbor element with the same integration point and unit outer normal | ||
517 | //! Todo: this code can likely be improved, we don't need to build all of the outside geometry if we know how to build scvf. | ||
518 | ✗ | SubControlVolumeFace outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace& scvf) const | |
519 | { | ||
520 | ✗ | const SubControlVolumeFace& lateralOrthogonalScvf = this->lateralOrthogonalScvf(scvf); | |
521 | ✗ | assert(!lateralOrthogonalScvf.boundary()); | |
522 | |||
523 | ✗ | const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx())); | |
524 | |||
525 | ✗ | auto outsideFVGeometry = localView(gridGeometry()); | |
526 | ✗ | const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex(); | |
527 | ✗ | outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx)); | |
528 | |||
529 | ✗ | for (const auto& otherScvf : scvfs(outsideFVGeometry)) | |
530 | ✗ | if (otherScvf.localIndex() == otherLocalIdx) | |
531 | ✗ | return otherScvf; | |
532 | |||
533 | ✗ | DUNE_THROW(Dune::InvalidStateException, "No outside scvf found"); | |
534 | } | ||
535 | |||
536 | //! Get the scv on the outside side of a periodic boundary | ||
537 | SubControlVolume outsidePeriodicScv(const SubControlVolume& selfScv) const | ||
538 | { | ||
539 | return Detail::FCStaggered::outsidePeriodicScv(*this, selfScv); | ||
540 | } | ||
541 | |||
542 | //! Create the geometry of a given sub control volume | ||
543 | typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const | ||
544 | { | ||
545 | return Detail::FCStaggered::scvGeometry(*this, scv); | ||
546 | } | ||
547 | |||
548 | //! Create the geometry of a given sub control volume face | ||
549 | typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const | ||
550 | { | ||
551 | return Detail::FCStaggered::scvfGeometry(*this, scvf); | ||
552 | } | ||
553 | |||
554 | private: | ||
555 | //! Binding of an element preparing the geometries of the whole stencil | ||
556 | //! called by the local jacobian to prepare element assembly | ||
557 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 31250 times.
|
40090 | void bind_(const Element& element) |
558 | { | ||
559 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 8840 times.
|
8840 | typename GG::LocalIntersectionMapper localIsMapper; |
560 | 40090 | localIsMapper.update(gridGeometry().gridView(), element); | |
561 | |||
562 | 40090 | bindElement_(element, localIsMapper); | |
563 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 40090 times.
|
40090 | neighborScvIndices_.clear(); |
564 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 40090 times.
|
40090 | neighborScvs_.clear(); |
565 | |||
566 |
4/9✓ Branch 2 taken 125000 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 35360 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 44200 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 1 taken 31250 times.
|
288850 | for (const auto& intersection : intersections(gridGeometry().gridView(), element)) |
567 | { | ||
568 |
2/4✓ Branch 0 taken 94705 times.
✓ Branch 1 taken 61905 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
281610 | if (intersection.neighbor()) |
569 | { | ||
570 |
2/2✓ Branch 0 taken 77665 times.
✓ Branch 1 taken 77665 times.
|
155330 | const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside()); |
571 |
1/2✓ Branch 1 taken 34080 times.
✗ Branch 2 not taken.
|
155330 | const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx); |
572 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 34080 times.
|
155330 | const auto& neighborElement = intersection.outside(); |
573 |
2/4✓ Branch 1 taken 34080 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 34080 times.
✗ Branch 5 not taken.
|
155330 | const auto neighborElementIdx = gridGeometry().elementMapper().index(neighborElement); |
574 | 155330 | const auto& neighborElementGeometry = neighborElement.geometry(); | |
575 | |||
576 | // todo: could be done easier? | ||
577 | std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement; | ||
578 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
155330 | std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), neighborElementIdx*numScvsPerElement); |
579 | |||
580 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 34080 times.
|
34080 | typename GG::LocalIntersectionMapper localNeighborIsMapper; |
581 |
1/2✓ Branch 1 taken 34080 times.
✗ Branch 2 not taken.
|
155330 | localNeighborIsMapper.update(gridGeometry().gridView(), neighborElement); |
582 | |||
583 |
7/14✓ Branch 1 taken 155330 times.
✓ Branch 2 taken 485000 times.
✓ Branch 4 taken 242500 times.
✓ Branch 5 taken 242500 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 136320 times.
✗ Branch 10 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✓ Branch 6 taken 136320 times.
✓ Branch 11 taken 170400 times.
|
1272780 | for (const auto& neighborIntersection : intersections(gridGeometry().gridView(), neighborElement)) |
584 | { | ||
585 |
2/2✓ Branch 0 taken 68160 times.
✓ Branch 1 taken 68160 times.
|
621320 | const auto localNeighborScvIdx = localNeighborIsMapper.realToRefIdx(neighborIntersection.indexInInside()); |
586 |
2/2✓ Branch 0 taken 310660 times.
✓ Branch 1 taken 310660 times.
|
621320 | if (localNeighborScvIdx != localScvIdx && localNeighborScvIdx != localOppositeScvIdx) |
587 | { | ||
588 | |||
589 |
3/5✓ Branch 1 taken 189410 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 68160 times.
✗ Branch 5 not taken.
✓ Branch 0 taken 121250 times.
|
310660 | const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(neighborElement, neighborIntersection.indexInInside()); |
590 |
1/2✓ Branch 1 taken 68160 times.
✗ Branch 2 not taken.
|
310660 | neighborScvs_.push_back(SubControlVolume( |
591 | neighborElementGeometry, | ||
592 | 242500 | neighborIntersection.geometry(), | |
593 |
1/2✓ Branch 1 taken 68160 times.
✗ Branch 2 not taken.
|
310660 | globalScvIndicesOfNeighborElement[localNeighborScvIdx], |
594 | localNeighborScvIdx, | ||
595 | dofIndex, | ||
596 | 310660 | Dumux::normalAxis(neighborIntersection.centerUnitOuterNormal()), | |
597 | neighborElementIdx, | ||
598 |
2/3✓ Branch 1 taken 189410 times.
✗ Branch 2 not taken.
✓ Branch 0 taken 121250 times.
|
310660 | onDomainBoundary_(neighborIntersection) |
599 | )); | ||
600 | |||
601 | 310660 | neighborScvIndices_.push_back(globalScvIndicesOfNeighborElement[localNeighborScvIdx]); | |
602 | } | ||
603 | } | ||
604 | 34080 | } | |
605 | } | ||
606 | 40090 | } | |
607 | |||
608 | //! Bind only element-local | ||
609 | 1986382 | void bindElement_(const Element& element, const typename GG::LocalIntersectionMapper& localIsMapper) | |
610 | { | ||
611 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1986382 times.
|
1986382 | elementPtr_ = &element; |
612 | 1986382 | eIdx_ = gridGeometry().elementMapper().index(element); | |
613 | 1986382 | std::iota(scvIndicesOfElement_.begin(), scvIndicesOfElement_.end(), eIdx_*numScvsPerElement); | |
614 | 1986382 | scvfs_.clear(); | |
615 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
1986382 | scvfs_.resize(minNumScvfsPerElement); |
616 | |||
617 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 209832 times.
|
1986382 | const auto& elementGeometry = element.geometry(); |
618 | |||
619 |
5/12✓ Branch 1 taken 1986382 times.
✓ Branch 2 taken 7106200 times.
✓ Branch 6 taken 839328 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 839328 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 1049160 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 3 not taken.
✗ Branch 5 not taken.
|
12030230 | for (const auto& intersection : intersections(gridGeometry().gridView(), element)) |
620 | { | ||
621 |
1/2✓ Branch 1 taken 839328 times.
✗ Branch 2 not taken.
|
7945528 | const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside()); |
622 |
1/2✓ Branch 1 taken 839328 times.
✗ Branch 2 not taken.
|
7945528 | auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv); |
623 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 7945528 times.
|
7945528 | const auto& intersectionGeometry = intersection.geometry(); |
624 |
3/5✓ Branch 1 taken 4392428 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 839328 times.
✗ Branch 5 not taken.
✓ Branch 0 taken 3553100 times.
|
7945528 | const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(element, intersection.indexInInside()); |
625 | |||
626 |
2/2✓ Branch 0 taken 3972764 times.
✓ Branch 1 taken 3972764 times.
|
7945528 | scvs_[localScvIdx] = SubControlVolume( |
627 | elementGeometry, | ||
628 | intersectionGeometry, | ||
629 |
1/2✓ Branch 1 taken 839328 times.
✗ Branch 2 not taken.
|
7945528 | scvIndicesOfElement_[localScvIdx], |
630 | localScvIdx, | ||
631 | dofIndex, | ||
632 |
2/2✓ Branch 1 taken 3972764 times.
✓ Branch 2 taken 3972764 times.
|
15891056 | Dumux::normalAxis(intersection.centerUnitOuterNormal()), |
633 | eIdx_, | ||
634 |
2/3✓ Branch 1 taken 4392428 times.
✗ Branch 2 not taken.
✓ Branch 0 taken 3553100 times.
|
7945528 | onDomainBoundary_(intersection) |
635 | ); | ||
636 | |||
637 | // the frontal sub control volume face at the element center | ||
638 |
1/2✓ Branch 1 taken 839328 times.
✗ Branch 2 not taken.
|
7945528 | const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx); |
639 | 7945528 | scvfs_[localScvfIdx] = SubControlVolumeFace( | |
640 | elementGeometry, | ||
641 | intersectionGeometry, | ||
642 |
1/2✓ Branch 1 taken 839328 times.
✗ Branch 2 not taken.
|
7945528 | std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localOppositeScvIdx]}, |
643 | localScvfIdx, | ||
644 |
1/2✓ Branch 1 taken 839328 times.
✗ Branch 2 not taken.
|
7945528 | scvfIndices_()[localScvfIdx], |
645 | 7945528 | intersection.centerUnitOuterNormal(), | |
646 | SubControlVolumeFace::FaceType::frontal, | ||
647 | SubControlVolumeFace::BoundaryType::interior | ||
648 | ); | ||
649 | 7945528 | ++localScvfIdx; | |
650 | |||
651 | // the lateral sub control volume faces | ||
652 |
4/4✓ Branch 0 taken 14212400 times.
✓ Branch 1 taken 7106200 times.
✓ Branch 2 taken 1678656 times.
✓ Branch 3 taken 839328 times.
|
23836584 | for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper_.localLaterFaceIndices(localScvIdx), |
653 |
1/2✓ Branch 1 taken 1678656 times.
✗ Branch 2 not taken.
|
9624184 | [&](auto&& idx) { return localIsMapper.refToRealIdx(idx) ;}) |
654 | ) | ||
655 | { | ||
656 |
2/4✓ Branch 1 taken 15891056 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1678656 times.
✗ Branch 5 not taken.
|
15891056 | const auto& lateralIntersection = geometryHelper_.intersection(lateralFacetIndex, element); |
657 | 15891056 | const auto& lateralFacet = geometryHelper_.facet(lateralFacetIndex, element); | |
658 | 15891056 | const auto& lateralFacetGeometry = lateralFacet.geometry(); | |
659 | |||
660 | // helper lambda to get the lateral scvf's global inside and outside scv indices | ||
661 | 47673168 | const auto globalScvIndicesForLateralFace = [&] | |
662 | { | ||
663 | 47673168 | const auto globalOutsideScvIdx = [&] | |
664 | { | ||
665 |
3/3✓ Branch 0 taken 14026934 times.
✓ Branch 1 taken 1812154 times.
✓ Branch 2 taken 51968 times.
|
15891056 | if (lateralIntersection.neighbor()) |
666 | { | ||
667 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 15468156 times.
✓ Branch 3 taken 1626688 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 1626688 times.
✗ Branch 7 not taken.
|
15468156 | const auto parallelElemIdx = gridGeometry().elementMapper().index(lateralIntersection.outside()); |
668 | // todo: could be done easier? | ||
669 | std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement; | ||
670 | 15468156 | std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), parallelElemIdx*numScvsPerElement); | |
671 | 15468156 | return globalScvIndicesOfNeighborElement[localScvIdx]; | |
672 | } | ||
673 | else | ||
674 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 422900 times.
|
422900 | return gridGeometry().outsideVolVarIndex(scvfIndices_()[localScvfIdx]); |
675 | 15891056 | }(); | |
676 | |||
677 | 15891056 | return std::array{scvIndicesOfElement_[localScvIdx], globalOutsideScvIdx}; | |
678 |
1/2✓ Branch 1 taken 15891056 times.
✗ Branch 2 not taken.
|
15891056 | }(); |
679 | |||
680 | 33460768 | const auto boundaryType = [&] | |
681 | { | ||
682 |
2/3✓ Branch 0 taken 14212400 times.
✓ Branch 1 taken 1678656 times.
✗ Branch 2 not taken.
|
15891056 | if (onProcessorBoundary_(lateralIntersection)) |
683 | return SubControlVolumeFace::BoundaryType::processorBoundary; | ||
684 |
2/2✓ Branch 0 taken 13841468 times.
✓ Branch 1 taken 370932 times.
|
17569712 | else if (onDomainBoundary_(lateralIntersection)) |
685 | return SubControlVolumeFace::BoundaryType::physicalBoundary; | ||
686 | else | ||
687 | 15468156 | return SubControlVolumeFace::BoundaryType::interior; | |
688 |
3/4✓ Branch 2 taken 1678656 times.
✗ Branch 3 not taken.
✓ Branch 0 taken 7106200 times.
✓ Branch 1 taken 7106200 times.
|
15891056 | }(); |
689 | |||
690 |
1/2✓ Branch 1 taken 1678656 times.
✗ Branch 2 not taken.
|
15891056 | scvfs_[localScvfIdx] = SubControlVolumeFace( |
691 | elementGeometry, | ||
692 | intersectionGeometry, | ||
693 | lateralFacetGeometry, | ||
694 | globalScvIndicesForLateralFace, // TODO higher order | ||
695 | localScvfIdx, | ||
696 |
1/2✓ Branch 1 taken 1678656 times.
✗ Branch 2 not taken.
|
15891056 | scvfIndices_()[localScvfIdx], |
697 | 15891056 | lateralIntersection.centerUnitOuterNormal(), | |
698 | SubControlVolumeFace::FaceType::lateral, | ||
699 | boundaryType | ||
700 | ); | ||
701 | 15891056 | ++localScvfIdx; | |
702 | } | ||
703 | } | ||
704 | |||
705 | // do a second loop over all intersections to add frontal boundary faces | ||
706 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1986382 times.
|
1986382 | auto localScvfIdx = minNumScvfsPerElement; |
707 |
3/8✓ Branch 2 taken 7106200 times.
✗ Branch 3 not taken.
✗ Branch 5 not taken.
✓ Branch 6 taken 1049160 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 1 taken 1986382 times.
✗ Branch 9 not taken.
|
11190902 | for (const auto& intersection : intersections(gridGeometry().gridView(), element)) |
708 | { | ||
709 | // the frontal sub control volume face at a domain boundary (coincides with element face) | ||
710 |
2/3✓ Branch 1 taken 7760062 times.
✗ Branch 2 not taken.
✓ Branch 0 taken 185466 times.
|
7945528 | if (onDomainBoundary_(intersection)) |
711 | { | ||
712 |
1/2✓ Branch 1 taken 25984 times.
✗ Branch 2 not taken.
|
211450 | const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside()); |
713 | // the frontal sub control volume face at the boundary | ||
714 |
4/7✓ Branch 1 taken 118717 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 25984 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 25984 times.
✗ Branch 8 not taken.
✓ Branch 0 taken 92733 times.
|
422900 | scvfs_.push_back(SubControlVolumeFace( |
715 | element.geometry(), | ||
716 | 185466 | intersection.geometry(), | |
717 |
1/2✓ Branch 1 taken 25984 times.
✗ Branch 2 not taken.
|
211450 | std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localScvIdx]}, |
718 | localScvfIdx, | ||
719 |
1/2✓ Branch 1 taken 25984 times.
✗ Branch 2 not taken.
|
211450 | scvfIndices_()[localScvfIdx], |
720 | 211450 | intersection.centerUnitOuterNormal(), | |
721 | SubControlVolumeFace::FaceType::frontal, | ||
722 | SubControlVolumeFace::BoundaryType::physicalBoundary | ||
723 | )); | ||
724 | 211450 | ++localScvfIdx; | |
725 | } | ||
726 | } | ||
727 | |||
728 | if constexpr (!ConsistentlyOrientedGrid<typename GridView::Grid>{}) | ||
729 | { | ||
730 |
4/6✓ Branch 0 taken 3 times.
✓ Branch 1 taken 209829 times.
✓ Branch 3 taken 3 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 3 times.
✗ Branch 7 not taken.
|
209832 | static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true); |
731 |
3/4✓ Branch 0 taken 209832 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 23492 times.
✓ Branch 3 taken 186340 times.
|
209832 | if (makeConsistentlyOriented && scvfs_.size() > minNumScvfsPerElement) |
732 | { | ||
733 | // make sure frontal boundary scvfs are sorted correctly | ||
734 | 23492 | std::sort(scvfs_.begin() + minNumScvfsPerElement, scvfs_.end(), | |
735 |
2/20✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2492 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2492 times.
|
4984 | [](const auto& scvfLeft, const auto& scvfRight) |
736 |
2/24✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2492 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2492 times.
|
4984 | { return scvfLeft.insideScvIdx() < scvfRight.insideScvIdx(); } |
737 | ); | ||
738 | } | ||
739 | } | ||
740 | 1986382 | } | |
741 | |||
742 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 65192472 times.
|
65192472 | const auto& scvfIndices_() const |
743 | 65192472 | { return gridGeometry().scvfIndicesOfElement(eIdx_); } | |
744 | |||
745 | template<class Entry, class Container> | ||
746 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 103353122 times.
|
206706244 | const LocalIndexType findLocalIndex_(const Entry& entry, |
747 | const Container& container) const | ||
748 | { | ||
749 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 103353122 times.
|
206706244 | auto it = std::find(container.begin(), container.end(), entry); |
750 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 103353122 times.
|
206706244 | assert(it != container.end() && "Could not find the entry! Make sure to properly bind this class!"); |
751 | 206706244 | return std::distance(container.begin(), it); | |
752 | } | ||
753 | |||
754 |
2/8✓ Branch 0 taken 28292718 times.
✓ Branch 1 taken 374582 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
32092772 | bool onDomainBoundary_(const typename GridView::Intersection& intersection) const |
755 | { | ||
756 |
10/12✓ Branch 0 taken 25984 times.
✓ Branch 1 taken 813344 times.
✓ Branch 3 taken 25984 times.
✓ Branch 4 taken 813344 times.
✓ Branch 5 taken 25984 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 51968 times.
✓ Branch 8 taken 1626688 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 51968 times.
✓ Branch 11 taken 2280 times.
✓ Branch 12 taken 65880 times.
|
4252588 | return !intersection.neighbor() && intersection.boundary(); |
757 | } | ||
758 | |||
759 |
2/2✓ Branch 0 taken 14026934 times.
✓ Branch 1 taken 185466 times.
|
15891056 | bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const |
760 | { | ||
761 |
3/4✓ Branch 0 taken 422900 times.
✓ Branch 1 taken 1626688 times.
✓ Branch 2 taken 51968 times.
✗ Branch 3 not taken.
|
2101556 | return !intersection.neighbor() && !intersection.boundary(); |
762 | } | ||
763 | |||
764 | Dune::ReservedVector<SubControlVolumeFace, maxNumScvfsPerElement> scvfs_; | ||
765 | |||
766 | Dune::ReservedVector<SubControlVolume, numLateralScvfsPerElement> neighborScvs_; | ||
767 | Dune::ReservedVector<GridIndexType, numLateralScvfsPerElement> neighborScvIndices_; | ||
768 | |||
769 | std::array<SubControlVolume, numScvsPerElement> scvs_; | ||
770 | std::array<GridIndexType, numScvsPerElement> scvIndicesOfElement_; | ||
771 | |||
772 | const GridGeometry* gridGeometry_; | ||
773 | const Element* elementPtr_; | ||
774 | GridIndexType eIdx_; | ||
775 | typename GridGeometry::GeometryHelper geometryHelper_; | ||
776 | }; | ||
777 | |||
778 | } // end namespace Dumux | ||
779 | |||
780 | #endif | ||
781 |