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::FaceCenteredStaggeredFVGridGeometry | ||
11 | */ | ||
12 | #ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_GRID_GEOMETRY | ||
13 | #define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_GRID_GEOMETRY | ||
14 | |||
15 | #include <memory> | ||
16 | |||
17 | #include <dune/common/rangeutilities.hh> | ||
18 | #include <dune/grid/common/scsgmapper.hh> | ||
19 | |||
20 | #include <dumux/common/defaultmappertraits.hh> | ||
21 | #include <dumux/common/indextraits.hh> | ||
22 | #include <dumux/common/intersectionmapper.hh> | ||
23 | #include <dumux/common/math.hh> | ||
24 | |||
25 | #include <dumux/discretization/basegridgeometry.hh> | ||
26 | #include <dumux/discretization/checkoverlapsize.hh> | ||
27 | #include <dumux/discretization/method.hh> | ||
28 | #include <dumux/discretization/extrusion.hh> | ||
29 | |||
30 | #include <dumux/discretization/facecentered/staggered/subcontrolvolume.hh> | ||
31 | #include <dumux/discretization/facecentered/staggered/subcontrolvolumeface.hh> | ||
32 | #include <dumux/discretization/facecentered/staggered/fvelementgeometry.hh> | ||
33 | #include <dumux/discretization/facecentered/staggered/geometryhelper.hh> | ||
34 | #include <dumux/discretization/facecentered/staggered/connectivitymap.hh> | ||
35 | #include <dumux/discretization/facecentered/staggered/normalaxis.hh> | ||
36 | #include <dumux/discretization/facecentered/staggered/localintersectionindexmapper.hh> | ||
37 | |||
38 | #include <dumux/io/grid/periodicgridtraits.hh> | ||
39 | |||
40 | namespace Dumux { | ||
41 | |||
42 | /*! | ||
43 | * \ingroup FaceCenteredStaggeredDiscretization | ||
44 | * \brief The default traits for the face-center staggered finite volume grid geometry | ||
45 | * Defines the scv and scvf types and the mapper types | ||
46 | * \tparam GridView the grid view type | ||
47 | */ | ||
48 | template<class GridView> | ||
49 | struct FaceCenteredStaggeredDefaultGridGeometryTraits : public DefaultMapperTraits<GridView> | ||
50 | { | ||
51 | using SubControlVolume = FaceCenteredStaggeredSubControlVolume<GridView>; | ||
52 | using SubControlVolumeFace = FaceCenteredStaggeredSubControlVolumeFace<GridView>; | ||
53 | using IntersectionMapper = ConformingGridIntersectionMapper<GridView>; | ||
54 | using LocalIntersectionMapper = FaceCenteredStaggeredLocalIntersectionIndexMapper<GridView>; | ||
55 | using GeometryHelper = FaceCenteredStaggeredGeometryHelper<GridView>; | ||
56 | |||
57 | template<class GridGeometry> | ||
58 | using ConnectivityMap = FaceCenteredStaggeredConnectivityMap<GridGeometry>; | ||
59 | |||
60 | template<class GridGeometry, bool enableCache> | ||
61 | using LocalView = FaceCenteredStaggeredFVElementGeometry<GridGeometry, enableCache>; | ||
62 | |||
63 | struct StaticInfo | ||
64 | { | ||
65 | static constexpr auto dim = GridView::Grid::dimension; | ||
66 | static constexpr auto numFacesPerElement = dim * 2; | ||
67 | static constexpr auto numScvsPerElement = numFacesPerElement; | ||
68 | static constexpr auto numLateralScvfsPerScv = 2 * (dim - 1); | ||
69 | static constexpr auto numLateralScvfsPerElement = numFacesPerElement*numLateralScvfsPerScv; | ||
70 | static constexpr auto minNumScvfsPerElement = numLateralScvfsPerElement // number of lateral faces | ||
71 | + numFacesPerElement; // number of central frontal faces | ||
72 | static constexpr auto maxNumScvfsPerElement = minNumScvfsPerElement | ||
73 | + numFacesPerElement; // number of potential frontal faces on boundary | ||
74 | }; | ||
75 | }; | ||
76 | |||
77 | /*! | ||
78 | * \ingroup FaceCenteredStaggeredDiscretization | ||
79 | * \brief Base class for the finite volume geometry vector for face-centered staggered models | ||
80 | * This builds up the sub control volumes and sub control volume faces | ||
81 | * for each element. | ||
82 | */ | ||
83 | template<class GridView, | ||
84 | bool cachingEnabled = false, | ||
85 | class Traits = FaceCenteredStaggeredDefaultGridGeometryTraits<GridView>> | ||
86 | class FaceCenteredStaggeredFVGridGeometry; | ||
87 | |||
88 | /*! | ||
89 | * \ingroup FaceCenteredStaggeredDiscretization | ||
90 | * \brief Base class for the finite volume geometry vector for staggered models | ||
91 | * This builds up the sub control volumes and sub control volume faces | ||
92 | * for each element. Specialization in case the FVElementGeometries are stored. | ||
93 | */ | ||
94 | template<class GV, class Traits> | ||
95 | class FaceCenteredStaggeredFVGridGeometry<GV, true, Traits> | ||
96 | : public BaseGridGeometry<GV, Traits> | ||
97 | { | ||
98 | using ThisType = FaceCenteredStaggeredFVGridGeometry<GV, true, Traits>; | ||
99 | using ParentType = BaseGridGeometry<GV, Traits>; | ||
100 | using GridIndexType = typename IndexTraits<GV>::GridIndex; | ||
101 | using LocalIndexType = typename IndexTraits<GV>::LocalIndex; | ||
102 | using SmallLocalIndexType = typename IndexTraits<GV>::SmallLocalIndex; | ||
103 | using Element = typename GV::template Codim<0>::Entity; | ||
104 | |||
105 | using IntersectionMapper = typename Traits::IntersectionMapper; | ||
106 | using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>; | ||
107 | |||
108 | using Scalar = typename GV::ctype; | ||
109 | |||
110 | static constexpr auto dim = Traits::StaticInfo::dim; | ||
111 | static constexpr auto numScvsPerElement = Traits::StaticInfo::numScvsPerElement; | ||
112 | static constexpr auto numLateralScvfsPerScv = Traits::StaticInfo::numLateralScvfsPerScv; | ||
113 | static constexpr auto numLateralScvfsPerElement = Traits::StaticInfo::numLateralScvfsPerElement; | ||
114 | static constexpr auto minNumScvfsPerElement = Traits::StaticInfo::minNumScvfsPerElement; | ||
115 | static constexpr auto maxNumScvfsPerElement = Traits::StaticInfo::maxNumScvfsPerElement; | ||
116 | |||
117 | using ScvfCornerStorage = typename Traits::SubControlVolumeFace::Traits::CornerStorage; | ||
118 | using ScvCornerStorage = typename Traits::SubControlVolume::Traits::CornerStorage; | ||
119 | |||
120 | public: | ||
121 | //! export the discretization method this geometry belongs to | ||
122 | using DiscretizationMethod = DiscretizationMethods::FCStaggered; | ||
123 | static constexpr DiscretizationMethod discMethod{}; | ||
124 | |||
125 | static constexpr bool cachingEnabled = true; | ||
126 | |||
127 | //! export basic grid geometry type for the alternative constructor | ||
128 | using BasicGridGeometry = BasicGridGeometry_t<GV, Traits>; | ||
129 | //! export the type of the fv element geometry (the local view type) | ||
130 | using LocalView = typename Traits::template LocalView<ThisType, true>; | ||
131 | //! export the type of sub control volume | ||
132 | using SubControlVolume = typename Traits::SubControlVolume; | ||
133 | //! export the type of sub control volume | ||
134 | using SubControlVolumeFace = typename Traits::SubControlVolumeFace; | ||
135 | //! export the grid view type | ||
136 | using GridView = GV; | ||
137 | //! export the geometry helper type | ||
138 | using GeometryHelper = typename Traits::GeometryHelper; | ||
139 | //! export the local intersection mapper | ||
140 | using LocalIntersectionMapper = typename Traits::LocalIntersectionMapper; | ||
141 | //! export static information | ||
142 | using StaticInformation = typename Traits::StaticInfo; | ||
143 | //! export the type of extrusion | ||
144 | using Extrusion = Extrusion_t<Traits>; | ||
145 | //! export whether the grid(geometry) supports periodicity | ||
146 | using SupportsPeriodicity = typename PeriodicGridTraits<typename GV::Grid>::SupportsPeriodicity; | ||
147 | |||
148 | //! Constructor with basic grid geometry used to share state with another grid geometry on the same grid view | ||
149 | 70 | FaceCenteredStaggeredFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg, const std::string& paramGroup = "") | |
150 | : ParentType(std::move(gg)) | ||
151 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 68 times.
✗ Branch 2 not taken.
|
70 | , intersectionMapper_(this->gridView()) |
152 |
3/5✗ Branch 0 not taken.
✓ Branch 1 taken 68 times.
✓ Branch 3 taken 9 times.
✓ Branch 4 taken 3 times.
✗ Branch 2 not taken.
|
70 | , periodicGridTraits_(this->gridView().grid()) |
153 | { | ||
154 | // Check if the overlap size is what we expect | ||
155 |
2/2✓ Branch 1 taken 9 times.
✓ Branch 2 taken 61 times.
|
70 | if (!CheckOverlapSize<DiscretizationMethod>::isValid(this->gridView())) |
156 |
0/22✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 31 not taken.
✗ Branch 32 not taken.
|
9 | DUNE_THROW(Dune::InvalidStateException, "The staggered discretization method needs overlap of exactly 1 for parallel computations. " |
157 | << " Set the parameter \"Grid.Overlap\" in the input file."); | ||
158 | |||
159 |
1/2✓ Branch 1 taken 70 times.
✗ Branch 2 not taken.
|
70 | update_(); |
160 | 70 | } | |
161 | |||
162 | //! Constructor from gridView | ||
163 | 67 | FaceCenteredStaggeredFVGridGeometry(const GridView& gridView, const std::string& paramGroup = "") | |
164 |
1/2✓ Branch 1 taken 67 times.
✗ Branch 2 not taken.
|
67 | : FaceCenteredStaggeredFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView), paramGroup) |
165 | 67 | {} | |
166 | |||
167 | //! The total number of sub control volumes | ||
168 | 3101051 | std::size_t numScv() const | |
169 |
8/11✓ Branch 2 taken 2727226 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 161635 times.
✓ Branch 6 taken 195080 times.
✓ Branch 8 taken 1 times.
✓ Branch 9 taken 10 times.
✓ Branch 11 taken 1 times.
✗ Branch 12 not taken.
✗ Branch 7 not taken.
✓ Branch 1 taken 3536 times.
✓ Branch 10 taken 12441 times.
|
3101051 | { return scvs_.size(); } |
170 | |||
171 | //! The total number of sub control volume faces | ||
172 | 97 | std::size_t numScvf() const | |
173 | 97 | { return scvfs_.size(); } | |
174 | |||
175 | //! The total number of boundary sub control volumes | ||
176 | std::size_t numBoundaryScv() const | ||
177 | { return numBoundaryScv_; } | ||
178 | |||
179 | //! The total number of boundary sub control volume faces | ||
180 | std::size_t numBoundaryScvf() const | ||
181 | { return numBoundaryScvf_; } | ||
182 | |||
183 | //! The total number of intersections | ||
184 | std::size_t numIntersections() const | ||
185 | { return intersectionMapper_.numIntersections(); } | ||
186 | |||
187 | //! the total number of dofs | ||
188 | 693 | std::size_t numDofs() const | |
189 |
32/40✓ Branch 2 taken 62 times.
✗ Branch 3 not taken.
✓ Branch 13 taken 5 times.
✓ Branch 14 taken 3 times.
✓ Branch 16 taken 6 times.
✓ Branch 17 taken 2 times.
✓ Branch 21 taken 23 times.
✓ Branch 22 taken 6 times.
✓ Branch 24 taken 23 times.
✓ Branch 25 taken 4 times.
✓ Branch 27 taken 16 times.
✓ Branch 28 taken 4 times.
✓ Branch 30 taken 16 times.
✓ Branch 31 taken 9 times.
✓ Branch 11 taken 3 times.
✓ Branch 12 taken 13 times.
✓ Branch 15 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 4 taken 7 times.
✗ Branch 5 not taken.
✓ Branch 18 taken 8 times.
✓ Branch 19 taken 6 times.
✗ Branch 32 not taken.
✓ Branch 34 taken 9 times.
✗ Branch 35 not taken.
✓ Branch 37 taken 7 times.
✗ Branch 38 not taken.
✓ Branch 40 taken 7 times.
✗ Branch 41 not taken.
✓ Branch 20 taken 2 times.
✓ Branch 23 taken 2 times.
✓ Branch 9 taken 5 times.
✗ Branch 10 not taken.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 6 times.
✓ Branch 26 taken 2 times.
✗ Branch 29 not taken.
✓ Branch 6 taken 15 times.
✓ Branch 33 taken 12 times.
✓ Branch 36 taken 12 times.
|
661 | { return this->gridView().size(1); } |
190 | |||
191 | //! update all fvElementGeometries (call this after grid adaption) | ||
192 | void update(const GridView& gridView) | ||
193 | { | ||
194 | ParentType::update(gridView); | ||
195 | update_(); | ||
196 | } | ||
197 | |||
198 | //! update all fvElementGeometries (call this after grid adaption) | ||
199 | void update(GridView&& gridView) | ||
200 | { | ||
201 | ParentType::update(std::move(gridView)); | ||
202 | update_(); | ||
203 | } | ||
204 | |||
205 | //! Get a sub control volume with a global scv index | ||
206 | 2988115788 | const SubControlVolume& scv(GridIndexType scvIdx) const | |
207 |
55/58✓ Branch 0 taken 2477867 times.
✓ Branch 1 taken 229719267 times.
✓ Branch 3 taken 1358463 times.
✓ Branch 4 taken 80151511 times.
✓ Branch 6 taken 5671322 times.
✓ Branch 7 taken 196596861 times.
✓ Branch 8 taken 10326641 times.
✓ Branch 9 taken 145341678 times.
✓ Branch 10 taken 16709775 times.
✓ Branch 11 taken 33448219 times.
✓ Branch 12 taken 6379220 times.
✓ Branch 13 taken 158515514 times.
✓ Branch 15 taken 1792972 times.
✓ Branch 16 taken 108359150 times.
✓ Branch 18 taken 434830 times.
✓ Branch 19 taken 89792702 times.
✓ Branch 21 taken 21851656 times.
✓ Branch 22 taken 98029102 times.
✓ Branch 23 taken 149720 times.
✓ Branch 24 taken 2361752 times.
✓ Branch 26 taken 41793072 times.
✓ Branch 27 taken 1155832 times.
✓ Branch 28 taken 58243082 times.
✓ Branch 29 taken 20458072 times.
✓ Branch 30 taken 2140088 times.
✓ Branch 31 taken 3324788 times.
✓ Branch 32 taken 5554280 times.
✓ Branch 33 taken 5032090 times.
✓ Branch 36 taken 4829028 times.
✓ Branch 37 taken 1698710 times.
✓ Branch 38 taken 2745444 times.
✓ Branch 39 taken 2325636 times.
✓ Branch 40 taken 623226 times.
✓ Branch 41 taken 1266408 times.
✓ Branch 2 taken 924160 times.
✗ Branch 5 not taken.
✓ Branch 14 taken 3672 times.
✓ Branch 20 taken 3708656 times.
✓ Branch 34 taken 44283142 times.
✓ Branch 35 taken 2296484 times.
✓ Branch 42 taken 171864 times.
✓ Branch 17 taken 18365398 times.
✓ Branch 25 taken 50155974 times.
✓ Branch 43 taken 83874 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.
|
1974948960 | { return scvs_[scvIdx]; } |
208 | |||
209 | //! Iterator range for sub control volumes. Iterates over | ||
210 | //! all scvs of the element-local fvGeometry. | ||
211 | 86846386 | auto scvs(const LocalView& fvGeometry) const | |
212 | { | ||
213 | 86846386 | auto begin = scvs_.cbegin() + numScvsPerElement*fvGeometry.elementIndex(); | |
214 | 86846386 | const auto end = begin + numScvsPerElement; | |
215 | 86846386 | return Dune::IteratorRange<std::decay_t<decltype(begin)>>(begin, end); | |
216 | } | ||
217 | |||
218 | //! Get a sub control volume face with a global scvf index | ||
219 | 2005294246 | const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const | |
220 |
4/5✓ Branch 1 taken 92419950 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 77835960 times.
✓ Branch 4 taken 17018882 times.
✓ Branch 0 taken 432407460 times.
|
619682252 | { return scvfs_[scvfIdx]; } |
221 | |||
222 | //! Get the global sub control volume face indices of an element | ||
223 | 564482724 | const std::vector<GridIndexType>& scvfIndicesOfElement(GridIndexType eIdx) const | |
224 | 564482724 | { return scvfIndicesOfElement_[eIdx]; } | |
225 | |||
226 | /*! | ||
227 | * \brief Returns the connectivity map of which dofs have derivatives with respect | ||
228 | * to a given dof. | ||
229 | */ | ||
230 | const ConnectivityMap& connectivityMap() const | ||
231 |
1/2✓ Branch 1 taken 70 times.
✗ Branch 2 not taken.
|
8117112 | { return connectivityMap_; } |
232 | |||
233 | //! Returns whether one of the geometry's scvfs lies on a boundary | ||
234 | bool hasBoundaryScvf(GridIndexType eIdx) const | ||
235 | 5839695 | { return hasBoundaryScvf_[eIdx]; } | |
236 | |||
237 | //! Return a reference to the intersection mapper | ||
238 | const IntersectionMapper& intersectionMapper() const | ||
239 |
2/2✓ Branch 0 taken 336 times.
✓ Branch 1 taken 291656 times.
|
2151754 | { return intersectionMapper_; } |
240 | |||
241 | //! If a d.o.f. is on a periodic boundary | ||
242 | 4280752 | bool dofOnPeriodicBoundary(GridIndexType dofIdx) const | |
243 |
9/18✓ Branch 1 taken 3660598 times.
✓ Branch 2 taken 20358 times.
✓ Branch 5 taken 19364 times.
✗ Branch 6 not taken.
✓ Branch 13 taken 1650 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 30 times.
✓ Branch 17 taken 30 times.
✓ Branch 4 taken 3480124 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 3 not taken.
✓ Branch 10 taken 242 times.
✗ Branch 11 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 1136 times.
✗ Branch 9 not taken.
✗ Branch 12 not taken.
|
15559328 | { return periodicFaceMap_.count(dofIdx); } |
244 | |||
245 | //! The index of the d.o.f. on the other side of the periodic boundary | ||
246 | 1952 | GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const | |
247 |
6/15✓ Branch 2 taken 242 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 121 times.
✓ Branch 5 taken 121 times.
✓ Branch 7 taken 1650 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 1650 times.
✗ Branch 12 not taken.
✓ Branch 13 taken 60 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✗ Branch 11 not taken.
✗ Branch 14 not taken.
✗ Branch 6 not taken.
|
1952 | { return periodicFaceMap_.at(dofIdx); } |
248 | |||
249 | //! Returns the map between dofs across periodic boundaries | ||
250 | const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() const | ||
251 | 812 | { return periodicFaceMap_; } | |
252 | |||
253 | |||
254 | private: | ||
255 | |||
256 | 70 | void update_() | |
257 | { | ||
258 | // clear containers (necessary after grid refinement) | ||
259 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 70 times.
|
70 | scvs_.clear(); |
260 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 70 times.
|
70 | scvfs_.clear(); |
261 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 70 times.
|
70 | scvfIndicesOfElement_.clear(); |
262 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 65 times.
|
70 | intersectionMapper_.update(this->gridView()); |
263 | |||
264 | // determine size of containers | ||
265 | 70 | const auto numElements = this->gridView().size(0); | |
266 | 70 | scvfIndicesOfElement_.resize(numElements); | |
267 | 70 | hasBoundaryScvf_.resize(numElements, false); | |
268 | |||
269 | 70 | outSideBoundaryVolVarIdx_ = 0; | |
270 | 70 | numBoundaryScv_ = 0; | |
271 | 70 | numBoundaryScvf_ = 0; | |
272 | |||
273 | 70 | GeometryHelper geometryHelper(this->gridView()); | |
274 | |||
275 | // get the global scvf indices first | ||
276 | 70 | GridIndexType numScvfs = 0; | |
277 |
9/11✓ Branch 2 taken 1943 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 123745 times.
✓ Branch 6 taken 121696 times.
✓ Branch 8 taken 1768 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 1768 times.
✓ Branch 11 taken 3 times.
✓ Branch 4 taken 111 times.
✓ Branch 1 taken 340657 times.
✓ Branch 7 taken 4 times.
|
934460 | for (const auto& element : elements(this->gridView())) |
278 | { | ||
279 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 123571 times.
|
466107 | assert(numScvsPerElement == element.subEntities(1)); |
280 | |||
281 |
12/16✓ Branch 1 taken 1487846 times.
✓ Branch 2 taken 2 times.
✓ Branch 4 taken 123629 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 20330 times.
✓ Branch 7 taken 123629 times.
✓ Branch 11 taken 1935 times.
✓ Branch 12 taken 486776 times.
✓ Branch 0 taken 341037 times.
✓ Branch 3 taken 171 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 14 taken 1768 times.
✗ Branch 15 not taken.
✓ Branch 10 taken 486776 times.
✓ Branch 13 taken 121694 times.
|
3668817 | for (const auto& intersection : intersections(this->gridView(), element)) |
282 | { | ||
283 | // frontal scvf in element center | ||
284 | 1868114 | ++numScvfs; | |
285 | |||
286 | // lateral scvfs | ||
287 | 1868114 | numScvfs += numLateralScvfsPerScv; | |
288 | |||
289 | // handle physical domain boundary | ||
290 |
2/3✓ Branch 0 taken 47194 times.
✓ Branch 1 taken 1809430 times.
✗ Branch 2 not taken.
|
1870814 | if (onDomainBoundary_(intersection)) |
291 | { | ||
292 | 50170 | ++numBoundaryScv_; // frontal face | |
293 | 50168 | numBoundaryScv_ += numLateralScvfsPerScv; // boundary scvs for lateral faces | |
294 | |||
295 | // frontal scvf at boundary | ||
296 | 50170 | ++numScvfs; | |
297 | } | ||
298 | } | ||
299 | } | ||
300 | |||
301 | // allocate memory | ||
302 | 70 | const auto numScvs = numElements*numScvsPerElement; | |
303 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
70 | scvs_.resize(numScvs); |
304 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
70 | scvfs_.reserve(numScvfs); |
305 | |||
306 | // Build the scvs and scv faces | ||
307 | 70 | std::size_t globalScvfIdx = 0; | |
308 |
8/8✓ Branch 2 taken 1943 times.
✓ Branch 3 taken 111 times.
✓ Branch 5 taken 121703 times.
✓ Branch 6 taken 121694 times.
✓ Branch 7 taken 1772 times.
✓ Branch 8 taken 3 times.
✓ Branch 4 taken 109 times.
✓ Branch 1 taken 340657 times.
|
932583 | for (const auto& element : elements(this->gridView())) |
309 | { | ||
310 |
1/2✓ Branch 1 taken 1768 times.
✗ Branch 2 not taken.
|
466107 | const auto eIdx = this->elementMapper().index(element); |
311 |
1/2✓ Branch 1 taken 123571 times.
✗ Branch 2 not taken.
|
466107 | auto& globalScvfIndices = scvfIndicesOfElement_[eIdx]; |
312 |
1/2✓ Branch 1 taken 123571 times.
✗ Branch 2 not taken.
|
466107 | globalScvfIndices.resize(minNumScvfsPerElement); |
313 |
1/2✓ Branch 1 taken 123571 times.
✗ Branch 2 not taken.
|
466107 | globalScvfIndices.reserve(maxNumScvfsPerElement); |
314 | |||
315 | 7439018 | auto getGlobalScvIdx = [&](const auto elementIdx, const auto localScvIdx) | |
316 | 7439018 | { return numScvsPerElement*elementIdx + localScvIdx; }; | |
317 | |||
318 | 123571 | LocalIntersectionMapper localIsMapper; | |
319 |
1/2✓ Branch 1 taken 123571 times.
✗ Branch 2 not taken.
|
466107 | localIsMapper.update(this->gridView(), element); |
320 | |||
321 |
14/17✓ Branch 1 taken 1487737 times.
✓ Branch 2 taken 468 times.
✓ Branch 4 taken 814431 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 133184 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 11490 times.
✓ Branch 10 taken 9007 times.
✓ Branch 15 taken 1935 times.
✗ Branch 16 not taken.
✓ Branch 3 taken 681606 times.
✓ Branch 6 taken 7181 times.
✓ Branch 0 taken 340601 times.
✓ Branch 11 taken 486776 times.
✓ Branch 12 taken 8840 times.
✓ Branch 13 taken 486776 times.
✓ Branch 14 taken 121694 times.
|
3223752 | for (const auto& intersection : intersections(this->gridView(), element)) |
322 | { | ||
323 |
2/3✓ Branch 0 taken 681170 times.
✓ Branch 1 taken 1186508 times.
✗ Branch 2 not taken.
|
1868114 | const auto& intersectionUnitOuterNormal = intersection.centerUnitOuterNormal(); |
324 |
1/2✓ Branch 1 taken 493848 times.
✗ Branch 2 not taken.
|
1868114 | const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside()); |
325 | 1868114 | auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv); | |
326 | |||
327 | 1868114 | const auto globalScvIdx = getGlobalScvIdx(eIdx, localScvIdx); | |
328 |
4/5✓ Branch 1 taken 1186726 times.
✓ Branch 2 taken 218 times.
✓ Branch 4 taken 7072 times.
✗ Branch 5 not taken.
✓ Branch 0 taken 681170 times.
|
1868114 | const auto dofIndex = intersectionMapper().globalIntersectionIndex(element, intersection.indexInInside()); |
329 |
3/4✓ Branch 0 taken 252669 times.
✓ Branch 1 taken 252669 times.
✓ Branch 3 taken 7072 times.
✗ Branch 4 not taken.
|
2373452 | const auto localOppositeScvIdx = geometryHelper.localOppositeIdx(localScvIdx); |
330 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 18562 times.
✗ Branch 2 not taken.
|
1868114 | const auto& intersectionGeometry = intersection.geometry(); |
331 | 1868114 | const auto& elementGeometry = element.geometry(); | |
332 | |||
333 |
2/5✓ Branch 1 taken 505774 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 7072 times.
✗ Branch 0 not taken.
|
1868114 | assert(localIsMapper.refToRealIdx(localScvIdx) == intersection.indexInInside()); |
334 | |||
335 | // handle periodic boundaries | ||
336 |
2/2✓ Branch 0 taken 768 times.
✓ Branch 1 taken 497498 times.
|
499118 | if (onPeriodicBoundary_(intersection)) |
337 | { | ||
338 |
1/2✓ Branch 1 taken 768 times.
✗ Branch 2 not taken.
|
800 | this->setPeriodic(); |
339 | |||
340 | 800 | const auto& otherElement = intersection.outside(); | |
341 | |||
342 | 800 | SmallLocalIndexType otherIntersectionLocalIdx = 0; | |
343 | 800 | bool periodicFaceFound = false; | |
344 | |||
345 |
9/13✓ Branch 1 taken 868 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 896 times.
✓ Branch 5 taken 32 times.
✓ Branch 7 taken 800 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 3072 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 3072 times.
✓ Branch 13 taken 768 times.
✓ Branch 0 taken 28 times.
✓ Branch 3 taken 128 times.
✗ Branch 6 not taken.
|
10240 | for (const auto& otherIntersection : intersections(this->gridView(), otherElement)) |
346 | { | ||
347 |
2/2✓ Branch 0 taken 1148 times.
✓ Branch 1 taken 2052 times.
|
3200 | if (periodicFaceFound) |
348 | 1120 | continue; | |
349 | |||
350 |
5/5✓ Branch 0 taken 52 times.
✓ Branch 1 taken 848 times.
✓ Branch 2 taken 1184 times.
✓ Branch 3 taken 836 times.
✓ Branch 4 taken 1184 times.
|
4956 | if (Dune::FloatCmp::eq(intersectionUnitOuterNormal*otherIntersection.centerUnitOuterNormal(), -1.0, 1e-7)) |
351 | { | ||
352 |
1/2✓ Branch 1 taken 768 times.
✗ Branch 2 not taken.
|
800 | const auto periodicDofIdx = intersectionMapper().globalIntersectionIndex(otherElement, otherIntersectionLocalIdx); |
353 |
1/2✓ Branch 1 taken 800 times.
✗ Branch 2 not taken.
|
800 | periodicFaceMap_[dofIndex] = periodicDofIdx; |
354 | 800 | periodicFaceFound = true; | |
355 | } | ||
356 | |||
357 | 2052 | ++otherIntersectionLocalIdx; | |
358 | } | ||
359 | } | ||
360 | |||
361 | // the sub control volume | ||
362 |
3/5✓ Branch 1 taken 7104 times.
✓ Branch 2 taken 1861010 times.
✓ Branch 4 taken 7072 times.
✗ Branch 5 not taken.
✗ Branch 3 not taken.
|
1868210 | scvs_[globalScvIdx] = SubControlVolume( |
363 | elementGeometry, | ||
364 | intersectionGeometry, | ||
365 | globalScvIdx, | ||
366 | localScvIdx, | ||
367 | dofIndex, | ||
368 | 1868114 | Dumux::normalAxis(intersectionUnitOuterNormal), | |
369 | 1868114 | this->elementMapper().index(element), | |
370 |
1/2✓ Branch 1 taken 7072 times.
✗ Branch 2 not taken.
|
1868114 | onDomainBoundary_(intersection) |
371 | ); | ||
372 | |||
373 | // the frontal sub control volume face at the element center | ||
374 | 3736228 | scvfs_.emplace_back(elementGeometry, | |
375 | intersectionGeometry, | ||
376 | 1868114 | std::array{globalScvIdx, getGlobalScvIdx(eIdx, localOppositeScvIdx)}, | |
377 | localScvfIdx, | ||
378 | globalScvfIdx, | ||
379 | intersectionUnitOuterNormal, | ||
380 | 3736228 | SubControlVolumeFace::FaceType::frontal, | |
381 |
1/2✓ Branch 1 taken 1868050 times.
✗ Branch 2 not taken.
|
1868114 | SubControlVolumeFace::BoundaryType::interior |
382 | ); | ||
383 | |||
384 | 1868114 | globalScvfIndices[localScvfIdx] = globalScvfIdx++; | |
385 | 1868050 | ++localScvfIdx; | |
386 | |||
387 | // the lateral sub control volume faces | ||
388 |
2/2✓ Branch 0 taken 3758600 times.
✓ Branch 1 taken 1868050 times.
|
5626714 | for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper.localLaterFaceIndices(localScvIdx), |
389 |
1/2✓ Branch 1 taken 988568 times.
✗ Branch 2 not taken.
|
2856682 | [&](auto&& idx) { return localIsMapper.refToRealIdx(idx) ;}) |
390 | ) | ||
391 | { | ||
392 |
1/2✓ Branch 1 taken 3758600 times.
✗ Branch 2 not taken.
|
3758600 | const auto& lateralIntersection = geometryHelper.intersection(lateralFacetIndex, element); |
393 | |||
394 | // helper lambda to get the lateral scvf's global inside and outside scv indices | ||
395 | 11275800 | const auto globalScvIndicesForLateralFace = [&] | |
396 | { | ||
397 | 11275800 | const auto globalOutsideScvIdx = [&] | |
398 | { | ||
399 |
4/4✓ Branch 0 taken 3729600 times.
✓ Branch 1 taken 28448 times.
✓ Branch 2 taken 907352 times.
✓ Branch 3 taken 66752 times.
|
4732152 | if (lateralIntersection.neighbor()) |
400 | { | ||
401 |
1/2✓ Branch 2 taken 13632 times.
✗ Branch 3 not taken.
|
3652620 | const auto parallelElemIdx = this->elementMapper().index(lateralIntersection.outside()); |
402 | 3652620 | return getGlobalScvIdx(parallelElemIdx, localScvIdx); | |
403 | } | ||
404 |
2/2✓ Branch 0 taken 94384 times.
✓ Branch 1 taken 336 times.
|
106532 | else if (onDomainBoundary_(lateralIntersection)) |
405 | 105644 | return numScvs + outSideBoundaryVolVarIdx_++; | |
406 | else | ||
407 | 336 | return globalScvIdx; // fallback for parallel, won't be used anyway | |
408 | 3758600 | }(); | |
409 | |||
410 | 3758600 | return std::array{globalScvIdx, globalOutsideScvIdx}; | |
411 |
1/2✓ Branch 1 taken 1033176 times.
✗ Branch 2 not taken.
|
3758600 | }(); |
412 | |||
413 | 7532216 | const auto boundaryType = [&] | |
414 | { | ||
415 |
2/3✓ Branch 1 taken 60832 times.
✗ Branch 2 not taken.
✓ Branch 0 taken 3697768 times.
|
3758600 | if (onProcessorBoundary_(lateralIntersection)) |
416 | return SubControlVolumeFace::BoundaryType::processorBoundary; | ||
417 |
2/2✓ Branch 0 taken 3603384 times.
✓ Branch 1 taken 94384 times.
|
3818760 | else if (onDomainBoundary_(lateralIntersection)) |
418 | return SubControlVolumeFace::BoundaryType::physicalBoundary; | ||
419 | else | ||
420 | 3617848 | return SubControlVolumeFace::BoundaryType::interior; | |
421 | 3758600 | }(); | |
422 | |||
423 |
1/2✓ Branch 2 taken 988568 times.
✗ Branch 3 not taken.
|
3772744 | scvfs_.emplace_back( |
424 | elementGeometry, | ||
425 | intersectionGeometry, | ||
426 |
1/2✓ Branch 2 taken 2770032 times.
✗ Branch 3 not taken.
|
2770032 | geometryHelper.facet(lateralFacetIndex, element).geometry(), |
427 | globalScvIndicesForLateralFace, // TODO higher order | ||
428 | localScvfIdx, | ||
429 | globalScvfIdx, | ||
430 |
1/2✓ Branch 1 taken 14144 times.
✗ Branch 2 not taken.
|
3758600 | lateralIntersection.centerUnitOuterNormal(), |
431 |
2/3✓ Branch 0 taken 1362276 times.
✓ Branch 1 taken 1376420 times.
✗ Branch 2 not taken.
|
3713120 | SubControlVolumeFace::FaceType::lateral, |
432 | boundaryType | ||
433 | ); | ||
434 | |||
435 |
2/2✓ Branch 0 taken 139864 times.
✓ Branch 1 taken 3603720 times.
|
3758600 | globalScvfIndices[localScvfIdx] = globalScvfIdx++; |
436 | 3758600 | ++localScvfIdx; | |
437 | |||
438 |
2/2✓ Branch 0 taken 94384 times.
✓ Branch 1 taken 3603720 times.
|
3769308 | if (onDomainBoundary_(lateralIntersection)) |
439 | { | ||
440 | 105644 | ++numBoundaryScvf_; | |
441 | 105644 | hasBoundaryScvf_[eIdx] = true; | |
442 | } | ||
443 | } // end loop over lateral facets | ||
444 | |||
445 | } // end first loop over intersections | ||
446 | |||
447 | // do a second loop over all intersections to add frontal boundary faces | ||
448 | 466107 | int localScvfIdx = minNumScvfsPerElement; | |
449 |
13/16✓ Branch 0 taken 436 times.
✓ Branch 1 taken 464172 times.
✓ Branch 2 taken 1364275 times.
✓ Branch 3 taken 109 times.
✓ Branch 4 taken 121696 times.
✓ Branch 5 taken 1997 times.
✓ Branch 6 taken 8840 times.
✓ Branch 7 taken 133184 times.
✓ Branch 8 taken 1935 times.
✗ Branch 9 not taken.
✓ Branch 14 taken 1768 times.
✗ Branch 15 not taken.
✓ Branch 10 taken 486776 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 486776 times.
✓ Branch 13 taken 121694 times.
|
3668817 | for (const auto& intersection : intersections(this->gridView(), element)) |
450 | { | ||
451 | // the frontal sub control volume face at a domain boundary (coincides with element face) | ||
452 |
2/3✓ Branch 0 taken 47194 times.
✓ Branch 1 taken 1809430 times.
✗ Branch 2 not taken.
|
1870814 | if (onDomainBoundary_(intersection)) |
453 | { | ||
454 |
1/2✓ Branch 1 taken 256 times.
✗ Branch 2 not taken.
|
50170 | const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside()); |
455 | 50170 | const auto globalScvIdx = getGlobalScvIdx(eIdx, localScvIdx); | |
456 | 50170 | ++numBoundaryScvf_; | |
457 | |||
458 | // the frontal sub control volume face at the boundary | ||
459 |
3/5✓ Branch 1 taken 2956 times.
✓ Branch 2 taken 33396 times.
✓ Branch 4 taken 256 times.
✗ Branch 5 not taken.
✗ Branch 3 not taken.
|
50426 | scvfs_.emplace_back( |
460 | element.geometry(), | ||
461 |
1/2✓ Branch 1 taken 2700 times.
✗ Branch 2 not taken.
|
16518 | intersection.geometry(), |
462 | std::array{globalScvIdx, globalScvIdx}, // TODO outside boundary, periodic, parallel? | ||
463 | localScvfIdx, | ||
464 | globalScvfIdx, | ||
465 |
1/2✓ Branch 1 taken 256 times.
✗ Branch 2 not taken.
|
50170 | intersection.centerUnitOuterNormal(), |
466 | 100340 | SubControlVolumeFace::FaceType::frontal, | |
467 |
2/3✓ Branch 0 taken 6909 times.
✓ Branch 1 taken 7165 times.
✗ Branch 2 not taken.
|
50170 | SubControlVolumeFace::BoundaryType::physicalBoundary |
468 | ); | ||
469 | |||
470 |
1/2✓ Branch 1 taken 36352 times.
✗ Branch 2 not taken.
|
50170 | globalScvfIndices.push_back(globalScvfIdx); |
471 | 50170 | ++globalScvfIdx; | |
472 | 50170 | ++localScvfIdx; | |
473 | 50170 | hasBoundaryScvf_[eIdx] = true; | |
474 | } | ||
475 | } | ||
476 | } | ||
477 | |||
478 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
70 | connectivityMap_.update(*this); |
479 | 70 | } | |
480 | |||
481 |
6/10✓ Branch 0 taken 12986667 times.
✓ Branch 1 taken 62763 times.
✓ Branch 2 taken 11490 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 45480 times.
✓ Branch 6 taken 45480 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 11490 times.
✗ Branch 9 not taken.
|
13208418 | bool onDomainBoundary_(const typename GridView::Intersection& intersection) const |
482 | { | ||
483 |
15/20✓ Branch 0 taken 300660 times.
✓ Branch 1 taken 3192522 times.
✓ Branch 2 taken 276 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 11766 times.
✓ Branch 5 taken 7232 times.
✓ Branch 7 taken 552 times.
✓ Branch 8 taken 14464 times.
✓ Branch 9 taken 12042 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 276 times.
✓ Branch 12 taken 7232 times.
✓ Branch 13 taken 276 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 552 times.
✓ Branch 16 taken 14464 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 552 times.
✓ Branch 20 taken 552 times.
✗ Branch 21 not taken.
|
4093696 | return !intersection.neighbor() && intersection.boundary(); |
484 | } | ||
485 | |||
486 |
2/2✓ Branch 0 taken 3684120 times.
✓ Branch 1 taken 13984 times.
|
3713120 | bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const |
487 | { | ||
488 |
3/4✓ Branch 0 taken 94936 times.
✓ Branch 1 taken 921600 times.
✓ Branch 2 taken 67304 times.
✗ Branch 3 not taken.
|
1083840 | return !intersection.neighbor() && !intersection.boundary(); |
489 | } | ||
490 | |||
491 | 498702 | bool onPeriodicBoundary_(const typename GridView::Intersection& intersection) const | |
492 | { | ||
493 |
2/2✓ Branch 1 taken 32 times.
✓ Branch 2 taken 384 times.
|
1868582 | return periodicGridTraits_.isPeriodic(intersection); |
494 | } | ||
495 | |||
496 | // mappers | ||
497 | ConnectivityMap connectivityMap_; | ||
498 | IntersectionMapper intersectionMapper_; | ||
499 | |||
500 | std::vector<SubControlVolume> scvs_; | ||
501 | std::vector<SubControlVolumeFace> scvfs_; | ||
502 | GridIndexType numBoundaryScv_; | ||
503 | GridIndexType numBoundaryScvf_; | ||
504 | GridIndexType outSideBoundaryVolVarIdx_; | ||
505 | std::vector<bool> hasBoundaryScvf_; | ||
506 | |||
507 | std::vector<std::vector<GridIndexType>> scvfIndicesOfElement_; | ||
508 | |||
509 | // a map for periodic boundary vertices | ||
510 | std::unordered_map<GridIndexType, GridIndexType> periodicFaceMap_; | ||
511 | |||
512 | PeriodicGridTraits<typename GridView::Grid> periodicGridTraits_; | ||
513 | }; | ||
514 | |||
515 | /*! | ||
516 | * \ingroup FaceCenteredStaggeredDiscretization | ||
517 | * \brief Base class for the finite volume geometry vector for face-centered staggered models | ||
518 | * This builds up the sub control volumes and sub control volume faces | ||
519 | * for each element. Specialization in case the FVElementGeometries are stored. | ||
520 | */ | ||
521 | template<class GV, class Traits> | ||
522 | class FaceCenteredStaggeredFVGridGeometry<GV, false, Traits> | ||
523 | : public BaseGridGeometry<GV, Traits> | ||
524 | { | ||
525 | using ThisType = FaceCenteredStaggeredFVGridGeometry<GV, false, Traits>; | ||
526 | using ParentType = BaseGridGeometry<GV, Traits>; | ||
527 | using GridIndexType = typename IndexTraits<GV>::GridIndex; | ||
528 | using LocalIndexType = typename IndexTraits<GV>::LocalIndex; | ||
529 | using SmallLocalIndexType = typename IndexTraits<GV>::SmallLocalIndex; | ||
530 | using Element = typename GV::template Codim<0>::Entity; | ||
531 | |||
532 | using IntersectionMapper = typename Traits::IntersectionMapper; | ||
533 | using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>; | ||
534 | |||
535 | static constexpr auto dim = Traits::StaticInfo::dim; | ||
536 | static constexpr auto numScvsPerElement = Traits::StaticInfo::numScvsPerElement; | ||
537 | static constexpr auto numLateralScvfsPerScv = Traits::StaticInfo::numLateralScvfsPerScv; | ||
538 | static constexpr auto numLateralScvfsPerElement = Traits::StaticInfo::numLateralScvfsPerElement; | ||
539 | static constexpr auto minNumScvfsPerElement = Traits::StaticInfo::minNumScvfsPerElement; | ||
540 | static constexpr auto maxNumScvfsPerElement = Traits::StaticInfo::maxNumScvfsPerElement; | ||
541 | |||
542 | public: | ||
543 | //! export the discretization method this geometry belongs to | ||
544 | using DiscretizationMethod = DiscretizationMethods::FCStaggered; | ||
545 | static constexpr DiscretizationMethod discMethod{}; | ||
546 | |||
547 | static constexpr bool cachingEnabled = false; | ||
548 | |||
549 | //! export basic grid geometry type for the alternative constructor | ||
550 | using BasicGridGeometry = BasicGridGeometry_t<GV, Traits>; | ||
551 | //! export the type of the fv element geometry (the local view type) | ||
552 | using LocalView = typename Traits::template LocalView<ThisType, false>; | ||
553 | //! export the type of sub control volume | ||
554 | using SubControlVolume = typename Traits::SubControlVolume; | ||
555 | //! export the type of sub control volume | ||
556 | using SubControlVolumeFace = typename Traits::SubControlVolumeFace; | ||
557 | //! export the grid view type | ||
558 | using GridView = GV; | ||
559 | //! export the geometry helper type | ||
560 | using GeometryHelper = typename Traits::GeometryHelper; | ||
561 | //! export the local intersection mapper | ||
562 | using LocalIntersectionMapper = typename Traits::LocalIntersectionMapper; | ||
563 | //! export static information | ||
564 | using StaticInformation = typename Traits::StaticInfo; | ||
565 | //! export the type of extrusion | ||
566 | using Extrusion = Extrusion_t<Traits>; | ||
567 | //! export whether the grid(geometry) supports periodicity | ||
568 | using SupportsPeriodicity = typename PeriodicGridTraits<typename GV::Grid>::SupportsPeriodicity; | ||
569 | |||
570 | //! Constructor with basic grid geometry used to share state with another grid geometry on the same grid view | ||
571 | 4 | FaceCenteredStaggeredFVGridGeometry(std::shared_ptr<BasicGridGeometry> gg, const std::string& paramGroup = "") | |
572 | : ParentType(std::move(gg)) | ||
573 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | , intersectionMapper_(this->gridView()) |
574 |
2/5✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 3 times.
✗ Branch 2 not taken.
|
4 | , periodicGridTraits_(this->gridView().grid()) |
575 | { | ||
576 | // Check if the overlap size is what we expect | ||
577 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 4 times.
|
4 | if (!CheckOverlapSize<DiscretizationMethod>::isValid(this->gridView())) |
578 | ✗ | DUNE_THROW(Dune::InvalidStateException, "The staggered discretization method needs at least an overlap of 1 for parallel computations. " | |
579 | << " Set the parameter \"Grid.Overlap\" in the input file."); | ||
580 | |||
581 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | update_(); |
582 | 4 | } | |
583 | |||
584 | //! Constructor from gridView | ||
585 | 1 | FaceCenteredStaggeredFVGridGeometry(const GridView& gridView, const std::string& paramGroup = "") | |
586 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | : FaceCenteredStaggeredFVGridGeometry(std::make_shared<BasicGridGeometry>(gridView), paramGroup) |
587 | 1 | {} | |
588 | |||
589 | //! The total number of sub control volumes | ||
590 | 8 | std::size_t numScv() const | |
591 | 8 | { return numScvs_; } | |
592 | |||
593 | //! The total number of sub control volume faces | ||
594 | std::size_t numScvf() const | ||
595 | { return numScvf_; } | ||
596 | |||
597 | //! The total number of boundary sub control volumes | ||
598 | std::size_t numBoundaryScv() const | ||
599 | { return numBoundaryScv_; } | ||
600 | |||
601 | //! The total number of boundary sub control volume faces | ||
602 | std::size_t numBoundaryScvf() const | ||
603 | { return numBoundaryScvf_; } | ||
604 | |||
605 | //! The total number of intersections | ||
606 | std::size_t numIntersections() const | ||
607 | { return intersectionMapper_.numIntersections(); } | ||
608 | |||
609 | //! the total number of dofs | ||
610 | 34 | std::size_t numDofs() const | |
611 |
7/18✓ Branch 2 taken 4 times.
✗ Branch 3 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✓ Branch 21 taken 3 times.
✗ Branch 22 not taken.
✓ Branch 24 taken 3 times.
✗ Branch 25 not taken.
✓ Branch 27 taken 3 times.
✗ Branch 28 not taken.
✓ Branch 30 taken 3 times.
✗ Branch 31 not taken.
✗ Branch 12 not taken.
✓ Branch 15 taken 1 times.
✓ Branch 18 taken 1 times.
✗ Branch 19 not taken.
|
33 | { return this->gridView().size(1); } |
612 | |||
613 | /*! | ||
614 | * \brief Returns the connectivity map of which dofs have derivatives with respect | ||
615 | * to a given dof. | ||
616 | */ | ||
617 | const ConnectivityMap& connectivityMap() const | ||
618 | 248292 | { return connectivityMap_; } | |
619 | |||
620 | //! Returns whether one of the geometry's scvfs lies on a boundary | ||
621 | 62072 | bool hasBoundaryScvf(GridIndexType eIdx) const | |
622 | 62072 | { return hasBoundaryScvf_[eIdx]; } | |
623 | |||
624 | //! Return a reference to the intersection mapper | ||
625 | const IntersectionMapper& intersectionMapper() const | ||
626 |
5/6✓ Branch 1 taken 960578 times.
✓ Branch 2 taken 3553100 times.
✓ Branch 4 taken 68160 times.
✗ Branch 5 not taken.
✓ Branch 0 taken 121250 times.
✓ Branch 3 taken 3553100 times.
|
8256188 | { return intersectionMapper_; } |
627 | |||
628 | //! Get the global sub control volume face indices of an element | ||
629 | 65192472 | const std::vector<GridIndexType>& scvfIndicesOfElement(GridIndexType eIdx) const | |
630 | 65192472 | { return scvfIndicesOfElement_[eIdx]; } | |
631 | |||
632 | //! Get the global sub control volume face indices of an element | ||
633 | 422900 | GridIndexType outsideVolVarIndex(GridIndexType scvfIdx) const | |
634 | 422900 | { return outsideVolVarIndices_.at(scvfIdx); } | |
635 | |||
636 | //! update all fvElementGeometries (call this after grid adaption) | ||
637 | void update(const GridView& gridView) | ||
638 | { | ||
639 | ParentType::update(gridView); | ||
640 | update_(); | ||
641 | } | ||
642 | |||
643 | //! update all fvElementGeometries (call this after grid adaption) | ||
644 | void update(GridView&& gridView) | ||
645 | { | ||
646 | ParentType::update(std::move(gridView)); | ||
647 | update_(); | ||
648 | } | ||
649 | |||
650 | //! If a d.o.f. is on a periodic boundary | ||
651 | 3262 | bool dofOnPeriodicBoundary(GridIndexType dofIdx) const | |
652 |
0/4✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
6524 | { return periodicFaceMap_.count(dofIdx); } |
653 | |||
654 | //! The index of the d.o.f. on the other side of the periodic boundary | ||
655 | ✗ | GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const | |
656 | ✗ | { return periodicFaceMap_.at(dofIdx); } | |
657 | |||
658 | //! Returns the map between dofs across periodic boundaries | ||
659 | const std::unordered_map<GridIndexType, GridIndexType>& periodicDofMap() const | ||
660 | 28 | { return periodicFaceMap_; } | |
661 | |||
662 | private: | ||
663 | |||
664 | 4 | void update_() | |
665 | { | ||
666 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
4 | intersectionMapper_.update(this->gridView()); |
667 | |||
668 | // clear local data | ||
669 | 4 | numScvf_ = 0; | |
670 | 4 | numBoundaryScv_ = 0; | |
671 | 4 | numBoundaryScvf_ = 0; | |
672 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
4 | hasBoundaryScvf_.clear(); |
673 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4 times.
|
4 | scvfIndicesOfElement_.clear(); |
674 | 4 | outsideVolVarIndices_.clear(); | |
675 | |||
676 | // determine size of containers | ||
677 | 4 | const auto numElements = this->gridView().size(0); | |
678 | 4 | scvfIndicesOfElement_.resize(numElements); | |
679 | 4 | hasBoundaryScvf_.resize(numElements, false); | |
680 | 4 | numScvs_ = numElements*numScvsPerElement; | |
681 | |||
682 | 4 | GeometryHelper geometryHelper(this->gridView()); | |
683 | |||
684 | // get the global scv indices first | ||
685 | 4 | GridIndexType scvfIdx = 0; | |
686 | |||
687 | 4 | GridIndexType neighborVolVarIdx = numScvs_; | |
688 | |||
689 |
5/8✓ Branch 2 taken 1254 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 3 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 1768 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 1768 times.
✓ Branch 11 taken 3 times.
|
7299 | for (const auto& element : elements(this->gridView())) |
690 | { | ||
691 |
2/4✓ Branch 1 taken 1768 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1768 times.
✗ Branch 5 not taken.
|
3018 | const auto eIdx = this->elementMapper().index(element); |
692 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1768 times.
|
3018 | assert(numScvsPerElement == element.subEntities(1)); |
693 | |||
694 | // the element-wise index sets for finite volume geometry | ||
695 |
1/2✓ Branch 1 taken 1768 times.
✗ Branch 2 not taken.
|
3018 | auto& globalScvfIndices = scvfIndicesOfElement_[eIdx]; |
696 |
1/2✓ Branch 1 taken 1768 times.
✗ Branch 2 not taken.
|
3018 | globalScvfIndices.reserve(maxNumScvfsPerElement); |
697 |
1/2✓ Branch 1 taken 1768 times.
✗ Branch 2 not taken.
|
3018 | globalScvfIndices.resize(minNumScvfsPerElement); |
698 | |||
699 | // keep track of frontal boundary scvfs | ||
700 | 3018 | std::size_t numFrontalBoundaryScvfs = 0; | |
701 | |||
702 | using LocalIntersectionIndexMapper = FaceCenteredStaggeredLocalIntersectionIndexMapper<GridView>; | ||
703 | 1768 | LocalIntersectionIndexMapper localIsMapper; | |
704 |
1/2✓ Branch 1 taken 1768 times.
✗ Branch 2 not taken.
|
3018 | localIsMapper.update(this->gridView(), element); |
705 | |||
706 |
4/13✓ Branch 1 taken 6768 times.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 6 taken 7072 times.
✗ Branch 7 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✓ Branch 8 taken 8840 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✓ Branch 0 taken 1250 times.
|
25698 | for (const auto& intersection : intersections(this->gridView(), element)) |
707 | { | ||
708 |
1/2✓ Branch 1 taken 7072 times.
✗ Branch 2 not taken.
|
12072 | const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside()); |
709 | 12072 | auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv); | |
710 | |||
711 |
2/4✓ Branch 1 taken 7072 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 7072 times.
|
12072 | assert(localIsMapper.refToRealIdx(localScvIdx) == intersection.indexInInside()); |
712 | // the frontal sub control volume face at the element center | ||
713 | 12072 | globalScvfIndices[localScvfIdx] = scvfIdx++; | |
714 | 12072 | ++localScvfIdx; | |
715 | |||
716 | if constexpr(dim > 1) | ||
717 | { | ||
718 | // the lateral sub control volume faces | ||
719 |
2/2✓ Branch 0 taken 24144 times.
✓ Branch 1 taken 12072 times.
|
36216 | for (const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper.localLaterFaceIndices(localScvIdx), |
720 |
1/2✓ Branch 1 taken 14144 times.
✗ Branch 2 not taken.
|
26216 | [&](auto idx) { return localIsMapper.refToRealIdx(idx) ;}) |
721 | ) | ||
722 | { | ||
723 |
4/5✓ Branch 1 taken 24144 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 10212 times.
✓ Branch 5 taken 13632 times.
✓ Branch 3 taken 300 times.
|
38288 | if (onDomainBoundary_(geometryHelper.intersection(lateralFacetIndex, element))) |
724 | { | ||
725 |
1/2✓ Branch 1 taken 812 times.
✗ Branch 2 not taken.
|
812 | outsideVolVarIndices_[scvfIdx] = neighborVolVarIdx++; |
726 | 812 | ++numBoundaryScvf_; | |
727 | 812 | hasBoundaryScvf_[eIdx] = true; | |
728 | } | ||
729 | |||
730 | 24144 | globalScvfIndices[localScvfIdx] = scvfIdx++; | |
731 | 24144 | ++localScvfIdx; | |
732 | } | ||
733 | } | ||
734 | |||
735 | // handle physical domain boundary | ||
736 |
2/2✓ Branch 0 taken 150 times.
✓ Branch 1 taken 4850 times.
|
12328 | if (onDomainBoundary_(intersection)) |
737 | { | ||
738 | 406 | ++numBoundaryScv_; // frontal face | |
739 | 406 | numBoundaryScv_ += numLateralScvfsPerScv; // boundary scvs for lateral faces | |
740 | 406 | ++numFrontalBoundaryScvfs; | |
741 | 406 | ++numBoundaryScvf_; | |
742 | 406 | hasBoundaryScvf_[eIdx] = true; | |
743 | } | ||
744 | |||
745 | // handle periodic boundaries | ||
746 | ✗ | if (onPeriodicBoundary_(intersection)) | |
747 | { | ||
748 | ✗ | this->setPeriodic(); | |
749 | |||
750 | ✗ | const auto& otherElement = intersection.outside(); | |
751 | |||
752 | ✗ | SmallLocalIndexType otherIntersectionLocalIdx = 0; | |
753 | ✗ | bool periodicFaceFound = false; | |
754 | |||
755 | ✗ | for (const auto& otherIntersection : intersections(this->gridView(), otherElement)) | |
756 | { | ||
757 | ✗ | if (periodicFaceFound) | |
758 | ✗ | continue; | |
759 | |||
760 | ✗ | if (Dune::FloatCmp::eq(intersection.centerUnitOuterNormal()*otherIntersection.centerUnitOuterNormal(), -1.0, 1e-7)) | |
761 | { | ||
762 | ✗ | const auto periodicDofIdx = intersectionMapper().globalIntersectionIndex(otherElement, otherIntersectionLocalIdx); | |
763 | ✗ | const auto dofIndex = intersectionMapper().globalIntersectionIndex(element, localScvIdx); | |
764 | ✗ | periodicFaceMap_[dofIndex] = periodicDofIdx; | |
765 | ✗ | periodicFaceFound = true; | |
766 | } | ||
767 | |||
768 | ✗ | ++otherIntersectionLocalIdx; | |
769 | } | ||
770 | } | ||
771 | } | ||
772 | |||
773 | // add global indices of frontal boundary scvfs last | ||
774 |
2/2✓ Branch 0 taken 406 times.
✓ Branch 1 taken 3018 times.
|
3424 | for (std::size_t i = 0; i < numFrontalBoundaryScvfs; ++i) |
775 |
1/2✓ Branch 1 taken 256 times.
✗ Branch 2 not taken.
|
406 | globalScvfIndices.push_back(scvfIdx++); |
776 | } | ||
777 | |||
778 | // set number of subcontrolvolume faces | ||
779 | 4 | numScvf_ = scvfIdx; | |
780 | |||
781 | 4 | connectivityMap_.update(*this); | |
782 | 4 | } | |
783 | |||
784 |
2/4✓ Branch 0 taken 14775 times.
✓ Branch 1 taken 225 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
36216 | bool onDomainBoundary_(const typename GridView::Intersection& intersection) const |
785 | { | ||
786 |
5/7✓ Branch 0 taken 512 times.
✓ Branch 1 taken 13632 times.
✓ Branch 3 taken 256 times.
✓ Branch 4 taken 6816 times.
✓ Branch 5 taken 256 times.
✗ Branch 6 not taken.
✗ Branch 2 not taken.
|
21922 | return !intersection.neighbor() && intersection.boundary(); |
787 | } | ||
788 | |||
789 | bool onProcessorBoundary_(const typename GridView::Intersection& intersection) const | ||
790 | { | ||
791 | return !intersection.neighbor() && !intersection.boundary(); | ||
792 | } | ||
793 | |||
794 | ✗ | bool onPeriodicBoundary_(const typename GridView::Intersection& intersection) const | |
795 | { | ||
796 |
1/2✓ Branch 1 taken 7072 times.
✗ Branch 2 not taken.
|
12072 | return periodicGridTraits_.isPeriodic(intersection); |
797 | } | ||
798 | |||
799 | // mappers | ||
800 | ConnectivityMap connectivityMap_; | ||
801 | IntersectionMapper intersectionMapper_; | ||
802 | |||
803 | //! Information on the global number of geometries | ||
804 | std::size_t numScvs_; | ||
805 | std::size_t numScvf_; | ||
806 | std::size_t numBoundaryScv_; | ||
807 | std::size_t numBoundaryScvf_; | ||
808 | std::vector<bool> hasBoundaryScvf_; | ||
809 | |||
810 | std::vector<std::vector<GridIndexType>> scvfIndicesOfElement_; | ||
811 | |||
812 | // a map for periodic boundary vertices | ||
813 | std::unordered_map<GridIndexType, GridIndexType> periodicFaceMap_; | ||
814 | std::unordered_map<GridIndexType, GridIndexType> outsideVolVarIndices_; | ||
815 | |||
816 | PeriodicGridTraits<typename GridView::Grid> periodicGridTraits_; | ||
817 | }; | ||
818 | |||
819 | } // end namespace Dumux | ||
820 | |||
821 | #endif | ||
822 |