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 FacetCoupling | ||
10 | * \brief Base class for the finite volume grid geometry for box models in the | ||
11 | * context of models considering coupling of different domains across the | ||
12 | * bulk grid facets. This builds up the sub control volumes and sub control | ||
13 | * volume faces for each element of the grid partition. | ||
14 | */ | ||
15 | #ifndef DUMUX_FACETCOUPLING_BOX_GRID_FVGEOMETRY_HH | ||
16 | #define DUMUX_FACETCOUPLING_BOX_GRID_FVGEOMETRY_HH | ||
17 | |||
18 | #include <algorithm> | ||
19 | #include <utility> | ||
20 | |||
21 | #include <dune/grid/common/mcmgmapper.hh> | ||
22 | #include <dune/localfunctions/lagrange/lagrangelfecache.hh> | ||
23 | |||
24 | #include <dumux/common/indextraits.hh> | ||
25 | #include <dumux/discretization/method.hh> | ||
26 | #include <dumux/discretization/extrusion.hh> | ||
27 | #include <dumux/discretization/basegridgeometry.hh> | ||
28 | #include <dumux/discretization/box/boxgeometryhelper.hh> | ||
29 | #include <dumux/discretization/box/subcontrolvolume.hh> | ||
30 | |||
31 | #include <dumux/multidomain/facet/box/fvelementgeometry.hh> | ||
32 | #include <dumux/multidomain/facet/box/subcontrolvolumeface.hh> | ||
33 | #include <dumux/multidomain/facet/vertexmapper.hh> | ||
34 | |||
35 | namespace Dumux { | ||
36 | |||
37 | namespace Detail { | ||
38 | template<class GV, class T> | ||
39 | using BoxFacetCouplingGeometryHelper_t = Dune::Std::detected_or_t< | ||
40 | Dumux::BoxGeometryHelper<GV, GV::dimension, typename T::SubControlVolume, typename T::SubControlVolumeFace>, | ||
41 | SpecifiesGeometryHelper, | ||
42 | T | ||
43 | >; | ||
44 | } // end namespace Detail | ||
45 | |||
46 | /*! | ||
47 | * \ingroup FacetCoupling | ||
48 | * \brief The default traits for the finite volume grid geometry | ||
49 | * of the box scheme with coupling occurring across the element facets. | ||
50 | * Defines the scv and scvf types and the mapper types. | ||
51 | * \tparam the grid view type | ||
52 | */ | ||
53 | template<class GridView> | ||
54 | struct BoxFacetCouplingDefaultGridGeometryTraits | ||
55 | { | ||
56 | // use a specialized version of the box scvf | ||
57 | using SubControlVolume = BoxSubControlVolume<GridView>; | ||
58 | using SubControlVolumeFace = BoxFacetCouplingSubControlVolumeFace<GridView>; | ||
59 | |||
60 | template<class GridGeometry, bool enableCache> | ||
61 | using LocalView = BoxFacetCouplingFVElementGeometry<GridGeometry, enableCache>; | ||
62 | |||
63 | // per default we use an mcmg mapper for the elements | ||
64 | using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>; | ||
65 | // the default vertex mapper is the enriched vertex dof mapper | ||
66 | using VertexMapper = EnrichedVertexDofMapper<GridView>; | ||
67 | // Mapper type for mapping edges | ||
68 | using FacetMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>; | ||
69 | }; | ||
70 | |||
71 | /*! | ||
72 | * \ingroup FacetCoupling | ||
73 | * \brief Base class for the finite volume geometry vector for box schemes in the context | ||
74 | * of coupled models where the coupling occurs across the element facets. This builds | ||
75 | * up the sub control volumes and sub control volume faces. | ||
76 | * \note This class is specialized for versions with and without caching the fv geometries on the grid view | ||
77 | */ | ||
78 | template<class Scalar, | ||
79 | class GridView, | ||
80 | bool enableGridGeometryCache = false, | ||
81 | class Traits = BoxFacetCouplingDefaultGridGeometryTraits<GridView> > | ||
82 | class BoxFacetCouplingFVGridGeometry; | ||
83 | |||
84 | /*! | ||
85 | * \ingroup FacetCoupling | ||
86 | * \brief Base class for the finite volume geometry vector for box schemes in the context | ||
87 | * of coupled models where the coupling occurs across the element facets. This builds | ||
88 | * up the sub control volumes and sub control volume faces. | ||
89 | * \note This class is specialized for versions with and without caching the fv geometries on the grid view | ||
90 | */ | ||
91 | template<class Scalar, class GV, class Traits> | ||
92 | class BoxFacetCouplingFVGridGeometry<Scalar, GV, true, Traits> | ||
93 | : public BaseGridGeometry<GV, Traits> | ||
94 | { | ||
95 | using ThisType = BoxFacetCouplingFVGridGeometry<Scalar, GV, true, Traits>; | ||
96 | using ParentType = BaseGridGeometry<GV, Traits>; | ||
97 | using GridIndexType = typename IndexTraits<GV>::GridIndex; | ||
98 | using LocalIndexType = typename IndexTraits<GV>::LocalIndex; | ||
99 | |||
100 | using Element = typename GV::template Codim<0>::Entity; | ||
101 | using CoordScalar = typename GV::ctype; | ||
102 | static const int dim = GV::dimension; | ||
103 | static const int dimWorld = GV::dimensionworld; | ||
104 | |||
105 | public: | ||
106 | //! export the discretization method this geometry belongs to | ||
107 | using DiscretizationMethod = DiscretizationMethods::Box; | ||
108 | static constexpr DiscretizationMethod discMethod{}; | ||
109 | |||
110 | //! export the type of the fv element geometry (the local view type) | ||
111 | using LocalView = typename Traits::template LocalView<ThisType, true>; | ||
112 | //! export the type of sub control volume | ||
113 | using SubControlVolume = typename Traits::SubControlVolume; | ||
114 | //! export the type of sub control volume | ||
115 | using SubControlVolumeFace = typename Traits::SubControlVolumeFace; | ||
116 | //! export the type of extrusion | ||
117 | using Extrusion = Extrusion_t<Traits>; | ||
118 | //! export dof mapper type | ||
119 | using DofMapper = typename Traits::VertexMapper; | ||
120 | //! export the finite element cache type | ||
121 | using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>; | ||
122 | //! export the grid view type | ||
123 | using GridView = GV; | ||
124 | //! export the geometry helper type | ||
125 | using GeometryHelper = Detail::BoxFacetCouplingGeometryHelper_t<GV, Traits>; | ||
126 | |||
127 | //! Constructor | ||
128 | template<class FacetGridView, class CodimOneGridAdapter> | ||
129 | 2 | BoxFacetCouplingFVGridGeometry(const GridView& gridView, | |
130 | const FacetGridView& facetGridView, | ||
131 | const CodimOneGridAdapter& codimOneGridAdapter, | ||
132 | bool verbose = false) | ||
133 |
2/4✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
|
2 | : ParentType(gridView) |
134 | { | ||
135 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
2 | update_(facetGridView, codimOneGridAdapter, verbose); |
136 | 2 | } | |
137 | |||
138 | //! the vertex mapper is the dofMapper | ||
139 | const DofMapper& dofMapper() const | ||
140 | { return this->vertexMapper(); } | ||
141 | |||
142 | //! The total number of sub control volumes | ||
143 | std::size_t numScv() const | ||
144 | { return numScv_; } | ||
145 | |||
146 | //! The total number of sun control volume faces | ||
147 | std::size_t numScvf() const | ||
148 | { return numScvf_; } | ||
149 | |||
150 | //! The total number of boundary sub control volume faces | ||
151 | //! For compatibility reasons with cc methods | ||
152 | std::size_t numBoundaryScvf() const | ||
153 | { return numBoundaryScvf_; } | ||
154 | |||
155 | //! The total number of degrees of freedom | ||
156 | 2 | std::size_t numDofs() const | |
157 | 2 | { return this->vertexMapper().size(); } | |
158 | |||
159 | /*! | ||
160 | * \brief update all fvElementGeometries (call this after grid adaption) | ||
161 | * \note This assumes conforming grids! | ||
162 | * | ||
163 | * \param gridView The grid view of a dim-dimensional grid. | ||
164 | * | ||
165 | * \param facetGridView The grid view of a (dim-1)-dimensional grid conforming | ||
166 | * with the facets of this grid view, indicating on which facets | ||
167 | * nodal dofs should be enriched. | ||
168 | * \param codimOneGridAdapter Adapter class that allows access to information on the d- | ||
169 | * dimensional grid for entities of the (d-1)-dimensional grid | ||
170 | * \param verbose Verbosity level for vertex enrichment | ||
171 | */ | ||
172 | template<class FacetGridView, class CodimOneGridAdapter> | ||
173 | void update(const GridView& gridView, | ||
174 | const FacetGridView& facetGridView, | ||
175 | const CodimOneGridAdapter& codimOneGridAdapter, | ||
176 | bool verbose = false) | ||
177 | { | ||
178 | ParentType::update(gridView); | ||
179 | update_(facetGridView, codimOneGridAdapter, verbose); | ||
180 | } | ||
181 | |||
182 | //! update all fvElementGeometries (call this after grid adaption) | ||
183 | template<class FacetGridView, class CodimOneGridAdapter> | ||
184 | void update(GridView&& gridView, | ||
185 | const FacetGridView& facetGridView, | ||
186 | const CodimOneGridAdapter& codimOneGridAdapter, | ||
187 | bool verbose = false) | ||
188 | { | ||
189 | ParentType::update(std::move(gridView)); | ||
190 | update_(facetGridView, codimOneGridAdapter, verbose); | ||
191 | } | ||
192 | |||
193 | //! The finite element cache for creating local FE bases | ||
194 | const FeCache& feCache() const | ||
195 | { return feCache_; } | ||
196 | |||
197 | //! Get the local scvs for an element | ||
198 | 768 | const std::vector<SubControlVolume>& scvs(GridIndexType eIdx) const | |
199 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 768 times.
|
768 | { return scvs_[eIdx]; } |
200 | |||
201 | //! Get the local scvfs for an element | ||
202 | 2560 | const std::vector<SubControlVolumeFace>& scvfs(GridIndexType eIdx) const | |
203 |
4/8✓ Branch 1 taken 768 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 768 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 768 times.
✓ Branch 9 taken 768 times.
✗ Branch 10 not taken.
|
2560 | { return scvfs_[eIdx]; } |
204 | |||
205 | //! If a d.o.f. is on the boundary | ||
206 | bool dofOnBoundary(GridIndexType dofIdx) const | ||
207 | { return boundaryDofIndices_[dofIdx]; } | ||
208 | |||
209 | //! If a d.o.f. is on an interior boundary | ||
210 | bool dofOnInteriorBoundary(GridIndexType dofIdx) const | ||
211 | { return interiorBoundaryDofIndices_[dofIdx]; } | ||
212 | |||
213 | //! Periodic boundaries are not supported for the box facet coupling scheme | ||
214 | bool dofOnPeriodicBoundary(GridIndexType dofIdx) const | ||
215 | { return false; } | ||
216 | |||
217 | //! The index of the vertex / d.o.f. on the other side of the periodic boundary | ||
218 | GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const | ||
219 | { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme"); } | ||
220 | |||
221 | private: | ||
222 | |||
223 | template<class FacetGridView, class CodimOneGridAdapter> | ||
224 | 2 | void update_(const FacetGridView& facetGridView, | |
225 | const CodimOneGridAdapter& codimOneGridAdapter, | ||
226 | bool verbose = false) | ||
227 | { | ||
228 | // enrich the vertex mapper subject to the provided facet grid | ||
229 | 2 | this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose); | |
230 | |||
231 | // resize containers | ||
232 | 2 | const auto numDof = numDofs(); | |
233 | 2 | const auto numElements = this->gridView().size(0); | |
234 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2 | scvs_.clear(); |
235 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2 | scvfs_.clear(); |
236 | 2 | scvs_.resize(numElements); | |
237 | 2 | scvfs_.resize(numElements); | |
238 | 2 | boundaryDofIndices_.assign(numDof, false); | |
239 | 2 | interiorBoundaryDofIndices_.assign(numDof, false); | |
240 | |||
241 | // Build the SCV and SCV faces | ||
242 | 2 | numScv_ = 0; | |
243 | 2 | numScvf_ = 0; | |
244 | 2 | numBoundaryScvf_ = 0; | |
245 |
6/8✓ Branch 2 taken 492 times.
✓ Branch 3 taken 1 times.
✓ Branch 5 taken 1 times.
✗ Branch 6 not taken.
✓ Branch 9 taken 491 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 491 times.
✓ Branch 12 taken 1 times.
|
1969 | for (const auto& element : elements(this->gridView())) |
246 | { | ||
247 |
2/4✓ Branch 1 taken 491 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 491 times.
✗ Branch 5 not taken.
|
982 | auto eIdx = this->elementMapper().index(element); |
248 | |||
249 | // keep track of number of scvs and scvfs | ||
250 |
1/2✓ Branch 1 taken 491 times.
✗ Branch 2 not taken.
|
982 | numScv_ += element.subEntities(dim); |
251 |
1/2✓ Branch 1 taken 491 times.
✗ Branch 2 not taken.
|
982 | numScvf_ += element.subEntities(dim-1); |
252 | |||
253 | // get the element geometry | ||
254 |
1/2✓ Branch 1 taken 491 times.
✗ Branch 2 not taken.
|
982 | auto elementGeometry = element.geometry(); |
255 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 982 times.
|
982 | const auto refElement = referenceElement(elementGeometry); |
256 | |||
257 | // instantiate the geometry helper | ||
258 | 982 | GeometryHelper geometryHelper(elementGeometry); | |
259 | |||
260 | // construct the sub control volumes | ||
261 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 982 times.
|
982 | scvs_[eIdx].clear(); |
262 |
1/2✓ Branch 1 taken 491 times.
✗ Branch 2 not taken.
|
982 | scvs_[eIdx].reserve(elementGeometry.corners()); |
263 |
2/2✓ Branch 0 taken 3928 times.
✓ Branch 1 taken 982 times.
|
4910 | for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx) |
264 |
2/4✓ Branch 1 taken 1964 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1964 times.
✗ Branch 5 not taken.
|
7856 | scvs_[eIdx].emplace_back(geometryHelper.getScvCorners(scvLocalIdx), |
265 | scvLocalIdx, | ||
266 | eIdx, | ||
267 |
2/4✓ Branch 1 taken 1964 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1964 times.
✗ Branch 5 not taken.
|
7856 | this->vertexMapper().subIndex(element, scvLocalIdx, dim)); |
268 | |||
269 | // construct the sub control volume faces | ||
270 | 982 | LocalIndexType scvfLocalIdx = 0; | |
271 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 982 times.
|
982 | scvfs_[eIdx].clear(); |
272 |
2/4✓ Branch 1 taken 491 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 491 times.
✗ Branch 5 not taken.
|
1473 | scvfs_[eIdx].reserve(element.subEntities(dim-1)); |
273 |
4/4✓ Branch 1 taken 6383 times.
✓ Branch 2 taken 491 times.
✓ Branch 3 taken 2946 times.
✓ Branch 4 taken 491 times.
|
10311 | for (; scvfLocalIdx < element.subEntities(dim-1); ++scvfLocalIdx) |
274 | { | ||
275 | // find the global and local scv indices this scvf is belonging to | ||
276 |
1/2✓ Branch 2 taken 5892 times.
✗ Branch 3 not taken.
|
11784 | std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)), |
277 |
1/2✓ Branch 2 taken 2946 times.
✗ Branch 3 not taken.
|
5892 | static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))}); |
278 | |||
279 | // create the sub-control volume face | ||
280 |
1/2✓ Branch 1 taken 5892 times.
✗ Branch 2 not taken.
|
5892 | scvfs_[eIdx].emplace_back(geometryHelper, |
281 | element, | ||
282 | elementGeometry, | ||
283 | scvfLocalIdx, | ||
284 | std::move(localScvIndices)); | ||
285 | } | ||
286 | |||
287 | // construct the sub control volume faces on the domain/interior boundaries | ||
288 | // skip handled facets (necessary for e.g. Dune::FoamGrid) | ||
289 | 982 | std::vector<unsigned int> handledFacets; | |
290 |
12/18✓ Branch 1 taken 982 times.
✗ Branch 2 not taken.
✓ Branch 6 taken 1964 times.
✓ Branch 7 taken 491 times.
✓ Branch 8 taken 1964 times.
✓ Branch 9 taken 1964 times.
✓ Branch 12 taken 1964 times.
✓ Branch 13 taken 1964 times.
✓ Branch 14 taken 2455 times.
✓ Branch 15 taken 2455 times.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✓ Branch 21 taken 491 times.
✗ Branch 22 not taken.
✓ Branch 4 taken 491 times.
✗ Branch 5 not taken.
✗ Branch 10 not taken.
✓ Branch 20 taken 491 times.
|
20131 | for (const auto& intersection : intersections(this->gridView(), element)) |
291 | { | ||
292 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 3928 times.
|
7856 | if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside())) |
293 | ✗ | continue; | |
294 | |||
295 |
2/4✓ Branch 1 taken 3928 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 3928 times.
✗ Branch 5 not taken.
|
3928 | handledFacets.push_back(intersection.indexInInside()); |
296 | |||
297 | // determine if all corners live on the facet grid | ||
298 |
1/2✓ Branch 1 taken 1964 times.
✗ Branch 2 not taken.
|
3928 | const auto isGeometry = intersection.geometry(); |
299 |
2/3✓ Branch 1 taken 3160 times.
✗ Branch 2 not taken.
✓ Branch 0 taken 768 times.
|
3928 | const auto numFaceCorners = isGeometry.corners(); |
300 |
2/2✓ Branch 0 taken 768 times.
✓ Branch 1 taken 1196 times.
|
3928 | const auto idxInInside = intersection.indexInInside(); |
301 |
1/2✓ Branch 1 taken 3928 times.
✗ Branch 2 not taken.
|
3928 | const auto boundary = intersection.boundary(); |
302 | |||
303 |
1/2✓ Branch 1 taken 3928 times.
✗ Branch 2 not taken.
|
3928 | std::vector<LocalIndexType> vIndicesLocal(numFaceCorners); |
304 |
2/2✓ Branch 0 taken 11784 times.
✓ Branch 1 taken 3928 times.
|
15712 | for (int i = 0; i < numFaceCorners; ++i) |
305 | 11784 | vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim)); | |
306 | |||
307 |
1/4✓ Branch 1 taken 3928 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
3928 | std::vector<LocalIndexType> gridVertexIndices(numFaceCorners); |
308 |
2/2✓ Branch 0 taken 11784 times.
✓ Branch 1 taken 3928 times.
|
15712 | for (int i = 0; i < numFaceCorners; ++i) |
309 |
1/2✓ Branch 1 taken 11784 times.
✗ Branch 2 not taken.
|
11784 | gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim); |
310 | |||
311 | // if the vertices compose a facet element, the intersection is on facet grid | ||
312 |
1/2✓ Branch 1 taken 3928 times.
✗ Branch 2 not taken.
|
3928 | const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices); |
313 | |||
314 | // make sure there are no periodic boundaries | ||
315 |
3/4✓ Branch 0 taken 480 times.
✓ Branch 1 taken 3448 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 480 times.
|
4168 | if (boundary && intersection.neighbor()) |
316 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme"); | |
317 | |||
318 |
4/4✓ Branch 0 taken 3800 times.
✓ Branch 1 taken 128 times.
✓ Branch 2 taken 480 times.
✓ Branch 3 taken 3320 times.
|
3928 | if (isOnFacet || boundary) |
319 | { | ||
320 | // keep track of number of faces | ||
321 | 608 | numScvf_ += numFaceCorners; | |
322 | 608 | numBoundaryScvf_ += int(boundary)*numFaceCorners; | |
323 | |||
324 |
2/2✓ Branch 0 taken 1824 times.
✓ Branch 1 taken 608 times.
|
2432 | for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numFaceCorners; ++isScvfLocalIdx) |
325 | { | ||
326 | // find the inside scv this scvf is belonging to (localIdx = element local vertex index) | ||
327 |
1/4✓ Branch 1 taken 1824 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
1824 | std::vector<LocalIndexType> localScvIndices = {vIndicesLocal[isScvfLocalIdx], vIndicesLocal[isScvfLocalIdx]}; |
328 | |||
329 | // create the sub-control volume face | ||
330 |
1/2✓ Branch 1 taken 1824 times.
✗ Branch 2 not taken.
|
1824 | scvfs_[eIdx].emplace_back(geometryHelper, |
331 | intersection, | ||
332 | isGeometry, | ||
333 | isScvfLocalIdx, | ||
334 | scvfLocalIdx, | ||
335 | std::move(localScvIndices), | ||
336 | boundary, | ||
337 | isOnFacet); | ||
338 | |||
339 | // Mark vertices to be on domain and/or interior boundary | ||
340 |
1/2✓ Branch 1 taken 1824 times.
✗ Branch 2 not taken.
|
1824 | const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[isScvfLocalIdx], dim); |
341 |
4/6✓ Branch 0 taken 1440 times.
✓ Branch 1 taken 384 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1440 times.
✓ Branch 4 taken 1440 times.
✗ Branch 5 not taken.
|
1824 | if (boundary) boundaryDofIndices_[ dofIndex ] = boundary && !isOnFacet; |
342 |
2/2✓ Branch 0 taken 384 times.
✓ Branch 1 taken 1440 times.
|
1824 | if (isOnFacet) interiorBoundaryDofIndices_[ dofIndex ] = isOnFacet; |
343 | |||
344 | // increment local counter | ||
345 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1824 times.
|
1824 | scvfLocalIdx++; |
346 | } | ||
347 | } | ||
348 | } | ||
349 | } | ||
350 | 2 | } | |
351 | |||
352 | const FeCache feCache_; | ||
353 | |||
354 | std::vector<std::vector<SubControlVolume>> scvs_; | ||
355 | std::vector<std::vector<SubControlVolumeFace>> scvfs_; | ||
356 | |||
357 | // TODO do we need those? | ||
358 | std::size_t numScv_; | ||
359 | std::size_t numScvf_; | ||
360 | std::size_t numBoundaryScvf_; | ||
361 | |||
362 | // vertices on domain/interior boundaries | ||
363 | std::vector<bool> boundaryDofIndices_; | ||
364 | std::vector<bool> interiorBoundaryDofIndices_; | ||
365 | }; | ||
366 | |||
367 | /*! | ||
368 | * \ingroup FacetCoupling | ||
369 | * \brief Base class for the finite volume geometry vector for box schemes | ||
370 | * This builds up the sub control volumes and sub control volume faces | ||
371 | * \note For caching disabled we store only some essential index maps to build up local systems on-demand in | ||
372 | * the corresponding FVElementGeometry | ||
373 | */ | ||
374 | template<class Scalar, class GV, class Traits> | ||
375 | class BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits> | ||
376 | : public BaseGridGeometry<GV, Traits> | ||
377 | { | ||
378 | using ThisType = BoxFacetCouplingFVGridGeometry<Scalar, GV, false, Traits>; | ||
379 | using ParentType = BaseGridGeometry<GV, Traits>; | ||
380 | using GridIndexType = typename IndexTraits<GV>::GridIndex; | ||
381 | using LocalIndexType = typename IndexTraits<GV>::LocalIndex; | ||
382 | |||
383 | static const int dim = GV::dimension; | ||
384 | static const int dimWorld = GV::dimensionworld; | ||
385 | |||
386 | using Element = typename GV::template Codim<0>::Entity; | ||
387 | using Intersection = typename GV::Intersection; | ||
388 | using CoordScalar = typename GV::ctype; | ||
389 | |||
390 | public: | ||
391 | //! export the discretization method this geometry belongs to | ||
392 | using DiscretizationMethod = DiscretizationMethods::Box; | ||
393 | static constexpr DiscretizationMethod discMethod{}; | ||
394 | |||
395 | //! export the type of the fv element geometry (the local view type) | ||
396 | using LocalView = typename Traits::template LocalView<ThisType, false>; | ||
397 | //! export the type of sub control volume | ||
398 | using SubControlVolume = typename Traits::SubControlVolume; | ||
399 | //! export the type of sub control volume | ||
400 | using SubControlVolumeFace = typename Traits::SubControlVolumeFace; | ||
401 | //! export the type of extrusion | ||
402 | using Extrusion = Extrusion_t<Traits>; | ||
403 | //! export dof mapper type | ||
404 | using DofMapper = typename Traits::VertexMapper; | ||
405 | //! export the finite element cache type | ||
406 | using FeCache = Dune::LagrangeLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>; | ||
407 | //! export the grid view type | ||
408 | using GridView = GV; | ||
409 | //! export the geometry helper type | ||
410 | using GeometryHelper = Detail::BoxFacetCouplingGeometryHelper_t<GV, Traits>; | ||
411 | |||
412 | //! Constructor | ||
413 | template<class FacetGridView, class CodimOneGridAdapter> | ||
414 | 24 | BoxFacetCouplingFVGridGeometry(const GridView& gridView, | |
415 | const FacetGridView& facetGridView, | ||
416 | const CodimOneGridAdapter& codimOneGridAdapter, | ||
417 | bool verbose = false) | ||
418 | : ParentType(gridView) | ||
419 |
3/6✓ Branch 2 taken 22 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 22 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 22 times.
✗ Branch 9 not taken.
|
48 | , facetMapper_(gridView, Dune::mcmgLayout(Dune::template Codim<1>())) |
420 | { | ||
421 |
1/2✓ Branch 1 taken 22 times.
✗ Branch 2 not taken.
|
24 | update_(facetGridView, codimOneGridAdapter, verbose); |
422 | 24 | } | |
423 | |||
424 | //! the vertex mapper is the dofMapper | ||
425 | //! this is convenience to have better chance to have the same main files for box/tpfa/mpfa... | ||
426 | 12844044 | const DofMapper& dofMapper() const | |
427 |
6/22✓ Branch 0 taken 9 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 12 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 76 times.
✗ Branch 6 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 11 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 25 not taken.
✗ Branch 14 not taken.
|
12844046 | { return this->vertexMapper(); } |
428 | |||
429 | //! The total number of sub control volumes | ||
430 | std::size_t numScv() const | ||
431 | { return numScv_; } | ||
432 | |||
433 | //! The total number of sun control volume faces | ||
434 | std::size_t numScvf() const | ||
435 | { return numScvf_; } | ||
436 | |||
437 | //! The total number of boundary sub control volume faces | ||
438 | //! For compatibility reasons with cc methods | ||
439 | std::size_t numBoundaryScvf() const | ||
440 | { return numBoundaryScvf_; } | ||
441 | |||
442 | //! The total number of degrees of freedom | ||
443 | 135 | std::size_t numDofs() const | |
444 |
7/15✓ Branch 1 taken 21 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 21 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 11 taken 2 times.
✗ Branch 12 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✓ Branch 17 taken 1 times.
✗ Branch 18 not taken.
✓ Branch 9 taken 19 times.
✗ Branch 10 not taken.
✗ Branch 13 not taken.
|
112 | { return this->vertexMapper().size(); } |
445 | |||
446 | /*! | ||
447 | * \brief update all fvElementGeometries (call this after grid adaption) | ||
448 | * \note This assumes conforming grids! | ||
449 | * | ||
450 | * \param gridView The grid view of a dim-dimensional grid. | ||
451 | * | ||
452 | * \param facetGridView The grid view of a (dim-1)-dimensional grid conforming | ||
453 | * with the facets of this grid view, indicating on which facets | ||
454 | * nodal dofs should be enriched. | ||
455 | * \param codimOneGridAdapter Adapter class that allows access to information on the d- | ||
456 | * dimensional grid for entities of the (d-1)-dimensional grid | ||
457 | * \param verbose Verbosity level for vertex enrichment | ||
458 | */ | ||
459 | template<class FacetGridView, class CodimOneGridAdapter> | ||
460 | void update(const GridView& gridView, | ||
461 | const FacetGridView& facetGridView, | ||
462 | const CodimOneGridAdapter& codimOneGridAdapter, | ||
463 | bool verbose = false) | ||
464 | { | ||
465 | ParentType::update(gridView); | ||
466 | updateFacetMapper_(); | ||
467 | update_(facetGridView, codimOneGridAdapter, verbose); | ||
468 | } | ||
469 | |||
470 | //! update all fvElementGeometries (call this after grid adaption) | ||
471 | template<class FacetGridView, class CodimOneGridAdapter> | ||
472 | void update(GridView&& gridView, | ||
473 | const FacetGridView& facetGridView, | ||
474 | const CodimOneGridAdapter& codimOneGridAdapter, | ||
475 | bool verbose = false) | ||
476 | { | ||
477 | ParentType::update(std::move(gridView)); | ||
478 | updateFacetMapper_(); | ||
479 | update_(facetGridView, codimOneGridAdapter, verbose); | ||
480 | } | ||
481 | |||
482 | //! The finite element cache for creating local FE bases | ||
483 | 8936836 | const FeCache& feCache() const | |
484 |
3/6✓ Branch 1 taken 1151256 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 3733008 times.
✗ Branch 4 not taken.
✓ Branch 8 taken 36690 times.
✗ Branch 9 not taken.
|
8936836 | { return feCache_; } |
485 | |||
486 | //! If a d.o.f. is on the boundary | ||
487 | 4308984 | bool dofOnBoundary(unsigned int dofIdx) const | |
488 |
4/4✓ Branch 0 taken 174865 times.
✓ Branch 1 taken 4074335 times.
✓ Branch 2 taken 15209 times.
✓ Branch 3 taken 44575 times.
|
4308984 | { return boundaryDofIndices_[dofIdx]; } |
489 | |||
490 | //! If a d.o.f. is on an interior boundary | ||
491 | bool dofOnInteriorBoundary(unsigned int dofIdx) const | ||
492 | { return interiorBoundaryDofIndices_[dofIdx]; } | ||
493 | |||
494 | //! returns true if an intersection is on an interior boundary | ||
495 | 5478739 | bool isOnInteriorBoundary(const Element& element, const Intersection& intersection) const | |
496 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6216 times.
|
5478739 | { return facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, intersection.indexInInside(), 1) ]; } |
497 | |||
498 | //! Periodic boundaries are not supported for the box facet coupling scheme | ||
499 | bool dofOnPeriodicBoundary(GridIndexType dofIdx) const | ||
500 | { return false; } | ||
501 | |||
502 | //! The index of the vertex / d.o.f. on the other side of the periodic boundary | ||
503 | GridIndexType periodicallyMappedDof(GridIndexType dofIdx) const | ||
504 | { DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the facet coupling scheme"); } | ||
505 | |||
506 | private: | ||
507 | |||
508 | void updateFacetMapper_() | ||
509 | { | ||
510 | facetMapper_.update(this->gridView()); | ||
511 | } | ||
512 | |||
513 | template<class FacetGridView, class CodimOneGridAdapter> | ||
514 | 24 | void update_(const FacetGridView& facetGridView, | |
515 | const CodimOneGridAdapter& codimOneGridAdapter, | ||
516 | bool verbose) | ||
517 | { | ||
518 | // enrich the vertex mapper subject to the provided facet grid | ||
519 | 24 | this->vertexMapper().enrich(facetGridView, codimOneGridAdapter, verbose); | |
520 | |||
521 | // save global data on the grid's scvs and scvfs | ||
522 | // TODO do we need those information? | ||
523 | 24 | numScv_ = 0; | |
524 | 24 | numScvf_ = 0; | |
525 | 24 | numBoundaryScvf_ = 0; | |
526 | |||
527 | 24 | const auto numDof = numDofs(); | |
528 | 24 | boundaryDofIndices_.assign(numDof, false); | |
529 | 24 | interiorBoundaryDofIndices_.assign(numDof, false); | |
530 | 24 | facetIsOnInteriorBoundary_.assign(this->gridView().size(1), false); | |
531 |
7/10✓ Branch 2 taken 109 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 109 times.
✓ Branch 6 taken 1 times.
✓ Branch 8 taken 52427 times.
✗ Branch 9 not taken.
✓ Branch 12 taken 52427 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 52427 times.
✓ Branch 15 taken 21 times.
|
163504 | for (const auto& element : elements(this->gridView())) |
532 | { | ||
533 |
1/2✓ Branch 1 taken 52427 times.
✗ Branch 2 not taken.
|
54530 | numScv_ += element.subEntities(dim); |
534 |
1/2✓ Branch 1 taken 52515 times.
✗ Branch 2 not taken.
|
54530 | numScvf_ += element.subEntities(dim-1); |
535 | |||
536 |
1/2✓ Branch 1 taken 52515 times.
✗ Branch 2 not taken.
|
54530 | const auto elementGeometry = element.geometry(); |
537 | 54530 | const auto refElement = referenceElement(elementGeometry); | |
538 | |||
539 | // store the sub control volume face indices on the domain/interior boundary | ||
540 | // skip handled facets (necessary for e.g. Dune::FoamGrid) | ||
541 | 54530 | std::vector<unsigned int> handledFacets; | |
542 |
7/14✓ Branch 1 taken 52515 times.
✗ Branch 2 not taken.
✓ Branch 8 taken 296 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 384 times.
✗ Branch 11 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 296 times.
✓ Branch 6 taken 201348 times.
✗ Branch 7 not taken.
✓ Branch 12 taken 201348 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 253775 times.
✗ Branch 15 not taken.
|
847592 | for (const auto& intersection : intersections(this->gridView(), element)) |
543 | { | ||
544 |
2/2✓ Branch 0 taken 32 times.
✓ Branch 1 taken 201612 times.
|
419296 | if (std::count(handledFacets.begin(), handledFacets.end(), intersection.indexInInside())) |
545 | 64 | continue; | |
546 | |||
547 |
2/4✓ Branch 1 taken 201612 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 201612 times.
✗ Branch 5 not taken.
|
209584 | handledFacets.push_back(intersection.indexInInside()); |
548 | |||
549 | // determine if all corners live on the facet grid | ||
550 |
1/2✓ Branch 1 taken 201348 times.
✗ Branch 2 not taken.
|
209584 | const auto isGeometry = intersection.geometry(); |
551 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 201612 times.
✗ Branch 2 not taken.
|
209584 | const auto numFaceCorners = isGeometry.corners(); |
552 |
1/2✓ Branch 1 taken 264 times.
✗ Branch 2 not taken.
|
209584 | const auto idxInInside = intersection.indexInInside(); |
553 |
1/2✓ Branch 1 taken 201612 times.
✗ Branch 2 not taken.
|
209584 | const auto boundary = intersection.boundary(); |
554 | |||
555 |
1/2✓ Branch 1 taken 201612 times.
✗ Branch 2 not taken.
|
209584 | std::vector<LocalIndexType> vIndicesLocal(numFaceCorners); |
556 |
2/2✓ Branch 0 taken 410932 times.
✓ Branch 1 taken 201612 times.
|
644168 | for (int i = 0; i < numFaceCorners; ++i) |
557 | 434584 | vIndicesLocal[i] = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, i, dim)); | |
558 | |||
559 |
1/2✓ Branch 1 taken 201612 times.
✗ Branch 2 not taken.
|
209584 | std::vector<GridIndexType> gridVertexIndices(numFaceCorners); |
560 |
2/2✓ Branch 0 taken 410932 times.
✓ Branch 1 taken 201612 times.
|
644168 | for (int i = 0; i < numFaceCorners; ++i) |
561 |
1/2✓ Branch 1 taken 410932 times.
✗ Branch 2 not taken.
|
434584 | gridVertexIndices[i] = this->vertexMapper().vertexIndex(element, vIndicesLocal[i], dim); |
562 | |||
563 | // if all vertices are living on the facet grid, this is an interiour boundary | ||
564 |
1/2✓ Branch 1 taken 201612 times.
✗ Branch 2 not taken.
|
209584 | const bool isOnFacet = codimOneGridAdapter.composeFacetElement(gridVertexIndices); |
565 | |||
566 | // make sure there are no periodic boundaries | ||
567 |
3/4✓ Branch 0 taken 4064 times.
✓ Branch 1 taken 197548 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 4064 times.
|
214360 | if (boundary && intersection.neighbor()) |
568 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Periodic boundaries are not supported by the box facet coupling scheme"); | |
569 | |||
570 |
2/2✓ Branch 0 taken 5744 times.
✓ Branch 1 taken 195868 times.
|
209584 | if (isOnFacet || boundary) |
571 | { | ||
572 | 6712 | numScvf_ += numFaceCorners; | |
573 | 6712 | numBoundaryScvf_ += int(boundary)*numFaceCorners; | |
574 | |||
575 | // Mark vertices to be on domain and/or interior boundary | ||
576 |
2/2✓ Branch 0 taken 12408 times.
✓ Branch 1 taken 5744 times.
|
21976 | for (int i = 0; i < numFaceCorners; ++i) |
577 | { | ||
578 |
1/2✓ Branch 1 taken 12312 times.
✗ Branch 2 not taken.
|
15264 | const auto dofIndex = this->vertexMapper().subIndex(element, vIndicesLocal[i], dim); |
579 |
4/4✓ Branch 0 taken 8872 times.
✓ Branch 1 taken 3536 times.
✓ Branch 2 taken 8776 times.
✓ Branch 3 taken 96 times.
|
15264 | if (boundary) boundaryDofIndices_[ dofIndex ] = boundary && !isOnFacet; |
580 |
2/2✓ Branch 0 taken 3632 times.
✓ Branch 1 taken 8776 times.
|
15264 | if (isOnFacet) |
581 | { | ||
582 |
1/2✓ Branch 1 taken 3632 times.
✗ Branch 2 not taken.
|
4192 | interiorBoundaryDofIndices_[ dofIndex ] = true; |
583 |
1/2✓ Branch 1 taken 3632 times.
✗ Branch 2 not taken.
|
4192 | facetIsOnInteriorBoundary_[ facetMapper_.subIndex(element, idxInInside, 1) ] = true; |
584 | } | ||
585 | } | ||
586 | } | ||
587 | } | ||
588 | } | ||
589 | 24 | } | |
590 | |||
591 | const FeCache feCache_; | ||
592 | |||
593 | // Information on the global number of geometries | ||
594 | // TODO do we need those information? | ||
595 | std::size_t numScv_; | ||
596 | std::size_t numScvf_; | ||
597 | std::size_t numBoundaryScvf_; | ||
598 | |||
599 | // vertices on domain/interior boundaries | ||
600 | std::vector<bool> boundaryDofIndices_; | ||
601 | std::vector<bool> interiorBoundaryDofIndices_; | ||
602 | |||
603 | // facet mapper and markers which facets lie on interior boundaries | ||
604 | typename Traits::FacetMapper facetMapper_; | ||
605 | std::vector<bool> facetIsOnInteriorBoundary_; | ||
606 | }; | ||
607 | |||
608 | } // end namespace Dumux | ||
609 | |||
610 | #endif | ||
611 |