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 BoxDFMModel | ||
10 | * \brief Base class for the local finite volume geometry for the box discrete | ||
11 | * fracture model. | ||
12 | * | ||
13 | * This builds up the sub control volumes and sub control | ||
14 | * volume faces for an element. | ||
15 | */ | ||
16 | |||
17 | #ifndef DUMUX_POROUSMEDIUMFLOW_BOXDFM_FV_ELEMENT_GEOMETRY_HH | ||
18 | #define DUMUX_POROUSMEDIUMFLOW_BOXDFM_FV_ELEMENT_GEOMETRY_HH | ||
19 | |||
20 | #include <optional> | ||
21 | #include <utility> | ||
22 | #include <ranges> | ||
23 | |||
24 | #include <dune/geometry/type.hh> | ||
25 | #include <dune/localfunctions/lagrange/pqkfactory.hh> | ||
26 | |||
27 | |||
28 | #include <dumux/common/indextraits.hh> | ||
29 | #include <dumux/discretization/scvandscvfiterators.hh> | ||
30 | #include "geometryhelper.hh" | ||
31 | |||
32 | namespace Dumux { | ||
33 | |||
34 | /*! | ||
35 | * \ingroup BoxDFMModel | ||
36 | * \brief Base class for the finite volume geometry vector for box discrete fracture model. | ||
37 | * | ||
38 | * This builds up the sub control volumes and sub control volume faces for each element. | ||
39 | * | ||
40 | * \tparam GG the finite volume grid geometry type | ||
41 | * \tparam enableGridGeometryCache if the grid geometry is cached or not | ||
42 | */ | ||
43 | template<class GG, bool enableGridGeometryCache> | ||
44 | class BoxDfmFVElementGeometry; | ||
45 | |||
46 | //! Specialization in case the FVElementGeometries are stored | ||
47 | template<class GG> | ||
48 | class BoxDfmFVElementGeometry<GG, true> | ||
49 | { | ||
50 | using GridView = typename GG::GridView; | ||
51 | static constexpr int dim = GridView::dimension; | ||
52 | static constexpr int dimWorld = GridView::dimensionworld; | ||
53 | using GridIndexType = typename IndexTraits<GridView>::GridIndex; | ||
54 | using LocalIndexType = typename IndexTraits<GridView>::LocalIndex; | ||
55 | using CoordScalar = typename GridView::ctype; | ||
56 | using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType; | ||
57 | public: | ||
58 | //! export type of the element | ||
59 | using Element = typename GridView::template Codim<0>::Entity; | ||
60 | //! Export type of subcontrol volume | ||
61 | using SubControlVolume = typename GG::SubControlVolume; | ||
62 | //! Export type of subcontrol volume face | ||
63 | using SubControlVolumeFace = typename GG::SubControlVolumeFace; | ||
64 | //! Export type of finite volume grid geometry | ||
65 | using GridGeometry = GG; | ||
66 | |||
67 | //! The maximum number of scvs per element (2^dim for cubes) | ||
68 | //! multiplied by 3 for the maximum number of fracture scvs per vertex | ||
69 | static constexpr std::size_t maxNumElementScvs = (1<<dim)*3; | ||
70 | |||
71 | //! Constructor | ||
72 | BoxDfmFVElementGeometry(const GridGeometry& gridGeometry) | ||
73 | : gridGeometryPtr_(&gridGeometry) {} | ||
74 | |||
75 | //! Get a sub control volume with a local scv index | ||
76 | const SubControlVolume& scv(std::size_t scvIdx) const | ||
77 | { return gridGeometry().scvs(eIdx_)[scvIdx]; } | ||
78 | |||
79 | //! Get a sub control volume face with a local scvf index | ||
80 | const SubControlVolumeFace& scvf(std::size_t scvfIdx) const | ||
81 | { return gridGeometry().scvfs(eIdx_)[scvfIdx]; } | ||
82 | |||
83 | /*! | ||
84 | * \brief Iterator range for sub control volumes. | ||
85 | * | ||
86 | * Iterates over all scvs of the bound element. | ||
87 | * This is a free function found by means of ADL. | ||
88 | * To iterate over all sub control volumes of this FVElementGeometry use | ||
89 | * for (auto&& scv : scvs(fvGeometry)). | ||
90 | */ | ||
91 | friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator> | ||
92 | scvs(const BoxDfmFVElementGeometry& fvGeometry) | ||
93 | { | ||
94 | const auto& g = fvGeometry.gridGeometry(); | ||
95 | using Iter = typename std::vector<SubControlVolume>::const_iterator; | ||
96 | return Dune::IteratorRange<Iter>(g.scvs(fvGeometry.eIdx_).begin(), g.scvs(fvGeometry.eIdx_).end()); | ||
97 | } | ||
98 | |||
99 | //! range over sub control volumes related to a local dof. | ||
100 | template<class LocalDof> | ||
101 | friend inline std::ranges::range auto | ||
102 | scvs(const BoxDfmFVElementGeometry& fvGeometry, const LocalDof& localDof) | ||
103 | { | ||
104 | return std::views::iota(0u, fvGeometry.numScv()) | ||
105 | | std::views::filter([&](std::size_t i) { return fvGeometry.scv(i).localDofIndex() == localDof.index(); }) | ||
106 | | std::views::transform([&](std::size_t i) -> const SubControlVolume& { return fvGeometry.scv(i); }); | ||
107 | } | ||
108 | |||
109 | //! range over local dofs | ||
110 | friend inline auto localDofs(const BoxDfmFVElementGeometry& fvGeometry) | ||
111 | { | ||
112 | return Dune::transformedRangeView( | ||
113 | Dune::range(fvGeometry.numLocalDofs()), | ||
114 | [&](const auto i) { return CVFE::LocalDof | ||
115 | { | ||
116 | static_cast<LocalIndexType>(fvGeometry.scv(i).localDofIndex()), | ||
117 | static_cast<GridIndexType>(fvGeometry.scv(i).dofIndex()), | ||
118 | static_cast<GridIndexType>(fvGeometry.scv(i).elementIndex()) | ||
119 | }; } | ||
120 | ); | ||
121 | } | ||
122 | |||
123 | //! range over control-volume local dofs | ||
124 | friend inline auto cvLocalDofs(const BoxDfmFVElementGeometry& fvGeometry) | ||
125 | { | ||
126 | // Default it that all dofs are cv dofs | ||
127 | return localDofs(fvGeometry); | ||
128 | } | ||
129 | |||
130 | /*! | ||
131 | * \brief Iterator range for sub control volumes faces. | ||
132 | * | ||
133 | * Iterates over all scvfs of the bound element. | ||
134 | * This is a free function found by means of ADL. | ||
135 | * To iterate over all sub control volume faces of this FVElementGeometry use | ||
136 | * for (auto&& scvf : scvfs(fvGeometry)). | ||
137 | */ | ||
138 | friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator> | ||
139 | scvfs(const BoxDfmFVElementGeometry& fvGeometry) | ||
140 | { | ||
141 | const auto& g = fvGeometry.gridGeometry(); | ||
142 | using Iter = typename std::vector<SubControlVolumeFace>::const_iterator; | ||
143 | return Dune::IteratorRange<Iter>(g.scvfs(fvGeometry.eIdx_).begin(), g.scvfs(fvGeometry.eIdx_).end()); | ||
144 | } | ||
145 | |||
146 | //! Get a local finite element basis | ||
147 | const FeLocalBasis& feLocalBasis() const | ||
148 | { return gridGeometry().feCache().get(element_->type()).localBasis(); } | ||
149 | |||
150 | //! The total number of sub control volumes | ||
151 | std::size_t numScv() const | ||
152 | { return gridGeometry().scvs(eIdx_).size(); } | ||
153 | |||
154 | //! The total number of element-local dofs | ||
155 | std::size_t numLocalDofs() const | ||
156 | { return element().geometry().corners(); } | ||
157 | |||
158 | //! The total number of sub control volume faces | ||
159 | std::size_t numScvf() const | ||
160 | { return gridGeometry().scvfs(eIdx_).size(); } | ||
161 | |||
162 | /*! | ||
163 | * \brief bind the local view (r-value overload) | ||
164 | * This overload is called when an instance of this class is a temporary in the usage context | ||
165 | * This allows a usage like this: `const auto view = localView(...).bind(element);` | ||
166 | */ | ||
167 | BoxDfmFVElementGeometry bind(const Element& element) && | ||
168 | { | ||
169 | this->bindElement(element); | ||
170 | return std::move(*this); | ||
171 | } | ||
172 | |||
173 | //! This function is for compatibility reasons with cc methods | ||
174 | //! The box stencil is always element-local so bind and bindElement are identical. | ||
175 | void bind(const Element& element) & | ||
176 | { this->bindElement(element); } | ||
177 | |||
178 | /*! | ||
179 | * \brief bind the local view (r-value overload) | ||
180 | * This overload is called when an instance of this class is a temporary in the usage context | ||
181 | * This allows a usage like this: `const auto view = localView(...).bindElement(element);` | ||
182 | */ | ||
183 | BoxDfmFVElementGeometry bindElement(const Element& element) && | ||
184 | { | ||
185 | this->bindElement(element); | ||
186 | return std::move(*this); | ||
187 | } | ||
188 | |||
189 | /*! | ||
190 | * \brief Binding of an element, has to be called before using the fvgeometries | ||
191 | * | ||
192 | * Prepares all the volume variables within the element. | ||
193 | * For compatibility reasons with the FVGeometry cache being disabled. | ||
194 | */ | ||
195 | void bindElement(const Element& element) & | ||
196 | { | ||
197 | element_ = element; | ||
198 | eIdx_ = gridGeometry().elementMapper().index(element); | ||
199 | } | ||
200 | |||
201 | //! The global finite volume geometry we are a restriction of | ||
202 | const GridGeometry& gridGeometry() const | ||
203 | { return *gridGeometryPtr_; } | ||
204 | |||
205 | //! Returns true if bind/bindElement has already been called | ||
206 | bool isBound() const | ||
207 | { return static_cast<bool>(element_); } | ||
208 | |||
209 | //! The bound element | ||
210 | const Element& element() const | ||
211 | { return *element_; } | ||
212 | |||
213 | //! Create the geometry of a given sub control volume | ||
214 | typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const | ||
215 | { | ||
216 | if (scv.isOnFracture()) | ||
217 | DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs " | ||
218 | "because the number of known corners is insufficient. " | ||
219 | "You can do this manually by extracting the corners from this scv " | ||
220 | "and extruding them by the corresponding aperture. "); | ||
221 | |||
222 | const typename GG::GeometryHelper geometryHelper(element().geometry()); | ||
223 | const auto corners = geometryHelper.getScvCorners(scv.index()); | ||
224 | using ScvGeometry = typename SubControlVolume::Traits::Geometry; | ||
225 | return { Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners }; | ||
226 | } | ||
227 | |||
228 | //! Create the geometry of a given sub control volume face | ||
229 | typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const | ||
230 | { | ||
231 | if (scvf.isOnFracture()) | ||
232 | DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs " | ||
233 | "because the number of known corners is insufficient. " | ||
234 | "You can do this manually by extracting the corners from this scv " | ||
235 | "and extruding them by the corresponding aperture. "); | ||
236 | const typename GG::GeometryHelper geometryHelper(element().geometry()); | ||
237 | const auto corners = geometryHelper.getScvfCorners(scvf.indexInElement()); | ||
238 | using ScvfGeometry = typename SubControlVolumeFace::Traits::Geometry; | ||
239 | return { Dune::GeometryTypes::cube(ScvfGeometry::mydimension), corners }; | ||
240 | } | ||
241 | |||
242 | private: | ||
243 | const GridGeometry* gridGeometryPtr_; | ||
244 | |||
245 | std::optional<Element> element_; | ||
246 | GridIndexType eIdx_; | ||
247 | }; | ||
248 | |||
249 | //! Specialization in case the FVElementGeometries are not stored | ||
250 | template<class GG> | ||
251 | class BoxDfmFVElementGeometry<GG, false> | ||
252 | { | ||
253 | using GridView = typename GG::GridView; | ||
254 | static constexpr int dim = GridView::dimension; | ||
255 | static constexpr int dimWorld = GridView::dimensionworld; | ||
256 | |||
257 | using GridIndexType = typename IndexTraits<GridView>::GridIndex; | ||
258 | using LocalIndexType = typename IndexTraits<GridView>::LocalIndex; | ||
259 | |||
260 | using CoordScalar = typename GridView::ctype; | ||
261 | using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType; | ||
262 | using GeometryHelper = typename GG::GeometryHelper; | ||
263 | public: | ||
264 | //! export type of the element | ||
265 | using Element = typename GridView::template Codim<0>::Entity; | ||
266 | //! Export type of subcontrol volume | ||
267 | using SubControlVolume = typename GG::SubControlVolume; | ||
268 | //! Export type of subcontrol volume face | ||
269 | using SubControlVolumeFace = typename GG::SubControlVolumeFace; | ||
270 | //! Export type of finite volume grid geometry | ||
271 | using GridGeometry = GG; | ||
272 | //! The maximum number of scvs per element (2^dim for cubes) | ||
273 | //! multiplied by 3 for the maximum number of fracture scvs per vertex | ||
274 | static constexpr std::size_t maxNumElementScvs = (1<<dim)*3; | ||
275 | |||
276 | //! Constructor | ||
277 | 2631222 | BoxDfmFVElementGeometry(const GridGeometry& gridGeometry) | |
278 |
3/12✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 7 taken 2416380 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 7 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 7 times.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
|
2631222 | : gridGeometryPtr_(&gridGeometry) {} |
279 | |||
280 | //! Get a sub control volume with a local scv index | ||
281 | 330250212 | const SubControlVolume& scv(std::size_t scvIdx) const | |
282 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 80159952 times.
✓ Branch 2 taken 80159952 times.
✗ Branch 3 not taken.
|
128215536 | { return scvs_[scvIdx]; } |
283 | |||
284 | //! Get a sub control volume face with a local scvf index | ||
285 | const SubControlVolumeFace& scvf(std::size_t scvfIdx) const | ||
286 | { return scvfs_[scvfIdx]; } | ||
287 | |||
288 | /*! | ||
289 | * \brief Iterator range for sub control volumes. | ||
290 | * | ||
291 | * Iterates over all scvs of the bound element. | ||
292 | * This is a free function found by means of ADL. | ||
293 | * To iterate over all sub control volumes of this FVElementGeometry use | ||
294 | * for (auto&& scv : scvs(fvGeometry)). | ||
295 | */ | ||
296 | friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator> | ||
297 | 118583478 | scvs(const BoxDfmFVElementGeometry& fvGeometry) | |
298 | { | ||
299 | using Iter = typename std::vector<SubControlVolume>::const_iterator; | ||
300 | 118583478 | return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end()); | |
301 | } | ||
302 | |||
303 | //! range over sub control volumes related to a local dof. | ||
304 | template<class LocalDof> | ||
305 | friend inline std::ranges::range auto | ||
306 | 18206280 | scvs(const BoxDfmFVElementGeometry& fvGeometry, const LocalDof& localDof) | |
307 | { | ||
308 | 18206280 | return std::views::iota(0u, fvGeometry.numScv()) | |
309 |
8/8✓ Branch 0 taken 21756780 times.
✓ Branch 1 taken 18206280 times.
✓ Branch 2 taken 9892752 times.
✓ Branch 3 taken 499080 times.
✓ Branch 4 taken 9892752 times.
✓ Branch 5 taken 499080 times.
✓ Branch 6 taken 4946376 times.
✓ Branch 7 taken 249540 times.
|
65942640 | | std::views::filter([&](std::size_t i) { return fvGeometry.scv(i).localDofIndex() == localDof.index(); }) |
310 | 37660260 | | std::views::transform([&](std::size_t i) -> const SubControlVolume& { return fvGeometry.scv(i); }); | |
311 | } | ||
312 | |||
313 | //! range over local dofs | ||
314 | 10845568 | friend inline auto localDofs(const BoxDfmFVElementGeometry& fvGeometry) | |
315 | { | ||
316 | return Dune::transformedRangeView( | ||
317 | 10845568 | Dune::range(fvGeometry.numLocalDofs()), | |
318 | 47283424 | [&](const auto i) { return CVFE::LocalDof | |
319 | { | ||
320 |
3/6✗ Branch 1 not taken.
✓ Branch 2 taken 24687936 times.
✓ Branch 4 taken 107880 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 57320 times.
✗ Branch 8 not taken.
|
36478104 | static_cast<LocalIndexType>(fvGeometry.scv(i).localDofIndex()), |
321 |
3/6✗ Branch 1 not taken.
✓ Branch 2 taken 24687936 times.
✓ Branch 4 taken 107880 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 57320 times.
✗ Branch 8 not taken.
|
36478104 | static_cast<GridIndexType>(fvGeometry.scv(i).dofIndex()), |
322 |
3/6✗ Branch 1 not taken.
✓ Branch 2 taken 24687936 times.
✓ Branch 4 taken 107880 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 57320 times.
✗ Branch 8 not taken.
|
36478104 | static_cast<GridIndexType>(fvGeometry.scv(i).elementIndex()) |
323 |
3/6✗ Branch 1 not taken.
✓ Branch 2 taken 24687936 times.
✓ Branch 4 taken 107880 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 57320 times.
✗ Branch 8 not taken.
|
36478104 | }; } |
324 | ); | ||
325 | } | ||
326 | |||
327 | //! range over control-volume local dofs | ||
328 | friend inline auto cvLocalDofs(const BoxDfmFVElementGeometry& fvGeometry) | ||
329 | { | ||
330 | // Default it that all dofs are cv dofs | ||
331 | return localDofs(fvGeometry); | ||
332 | } | ||
333 | |||
334 | /*! | ||
335 | * \brief Iterator range for sub control volumes faces. | ||
336 | * | ||
337 | * Iterates over all scvfs of the bound element. | ||
338 | * This is a free function found by means of ADL. | ||
339 | * To iterate over all sub control volume faces of this FVElementGeometry use | ||
340 | * for (auto&& scvf : scvfs(fvGeometry)). | ||
341 | */ | ||
342 | friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator> | ||
343 | 12115272 | scvfs(const BoxDfmFVElementGeometry& fvGeometry) | |
344 | { | ||
345 | using Iter = typename std::vector<SubControlVolumeFace>::const_iterator; | ||
346 | 12115272 | return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end()); | |
347 | } | ||
348 | |||
349 | //! Get a local finite element basis | ||
350 | 10326420 | const FeLocalBasis& feLocalBasis() const | |
351 | 10326420 | { return gridGeometry().feCache().get(element_->type()).localBasis(); } | |
352 | |||
353 | //! The total number of sub control volumes | ||
354 | 65720872 | std::size_t numScv() const | |
355 |
2/4✓ Branch 5 taken 10326420 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 10326420 times.
✗ Branch 9 not taken.
|
65720872 | { return scvs_.size(); } |
356 | |||
357 | //! The total number of element-local dofs | ||
358 | 11940970 | std::size_t numLocalDofs() const | |
359 |
3/6✓ Branch 6 taken 3249936 times.
✗ Branch 7 not taken.
✓ Branch 16 taken 4426 times.
✗ Branch 17 not taken.
✓ Branch 20 taken 15776 times.
✗ Branch 21 not taken.
|
10845568 | { return element().geometry().corners(); } |
360 | |||
361 | //! The total number of sub control volume faces | ||
362 | 2416380 | std::size_t numScvf() const | |
363 | 2416380 | { return scvfs_.size(); } | |
364 | |||
365 | /*! | ||
366 | * \brief bind the local view (r-value overload) | ||
367 | * This overload is called when an instance of this class is a temporary in the usage context | ||
368 | * This allows a usage like this: `const auto view = localView(...).bind(element);` | ||
369 | */ | ||
370 | BoxDfmFVElementGeometry bind(const Element& element) && | ||
371 | { | ||
372 | this->bindElement(element); | ||
373 | return std::move(*this); | ||
374 | } | ||
375 | |||
376 | /*! | ||
377 | * \brief Binding of an element, has to be called before using the fvgeometries | ||
378 | * Prepares all the volume variables within the element. | ||
379 | * \note For the box scheme, bind() and bindElement() are identical, but the | ||
380 | * distinction is here for the sake of compatibility with cc schemes. | ||
381 | */ | ||
382 | 2427406 | void bind(const Element& element) & | |
383 |
1/10✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 8 taken 11026 times.
✗ Branch 9 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
2427406 | { this->bindElement(element); } |
384 | |||
385 | /*! | ||
386 | * \brief bind the local view (r-value overload) | ||
387 | * This overload is called when an instance of this class is a temporary in the usage context | ||
388 | * This allows a usage like this: `const auto view = localView(...).bindElement(element);` | ||
389 | */ | ||
390 | BoxDfmFVElementGeometry bindElement(const Element& element) && | ||
391 | { | ||
392 | this->bindElement(element); | ||
393 | return std::move(*this); | ||
394 | } | ||
395 | |||
396 | /*! | ||
397 | * \brief Binding of an element, has to be called before using the fvgeometries | ||
398 | * Prepares all the volume variables within the element. | ||
399 | */ | ||
400 | 2653260 | void bindElement(const Element& element) & | |
401 | { | ||
402 |
2/2✓ Branch 0 taken 13192 times.
✓ Branch 1 taken 1433480 times.
|
2653260 | element_ = element; |
403 | 5306520 | eIdx_ = gridGeometry().elementMapper().index(element); | |
404 | 2653260 | makeElementGeometries_(); | |
405 | 2653260 | } | |
406 | |||
407 | //! The global finite volume geometry we are a restriction of | ||
408 | 68210698 | const GridGeometry& gridGeometry() const | |
409 |
3/4✓ Branch 0 taken 1404365 times.
✓ Branch 1 taken 7100755 times.
✓ Branch 4 taken 4421880 times.
✗ Branch 5 not taken.
|
27627790 | { return *gridGeometryPtr_; } |
410 | |||
411 | //! Returns true if bind/bindElement has already been called | ||
412 | bool isBound() const | ||
413 | { return static_cast<bool>(element_); } | ||
414 | |||
415 | //! The bound element | ||
416 | 7781592 | const Element& element() const | |
417 |
3/6✓ Branch 5 taken 3249936 times.
✗ Branch 6 not taken.
✓ Branch 11 taken 4426 times.
✗ Branch 12 not taken.
✓ Branch 14 taken 15776 times.
✗ Branch 15 not taken.
|
19125946 | { return *element_; } |
418 | |||
419 | //! Create the geometry of a given sub control volume | ||
420 | typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const | ||
421 | { | ||
422 | if (scv.isOnFracture()) | ||
423 | DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs " | ||
424 | "because the number of known corners is insufficient. " | ||
425 | "You can do this manually by extracting the corners from this scv " | ||
426 | "and extruding them by the corresponding aperture. "); | ||
427 | |||
428 | const GeometryHelper geometryHelper(element().geometry()); | ||
429 | const auto corners = geometryHelper.getScvCorners(scv.index()); | ||
430 | using ScvGeometry = typename SubControlVolume::Traits::Geometry; | ||
431 | return { Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners }; | ||
432 | } | ||
433 | |||
434 | //! Create the geometry of a given sub control volume face | ||
435 | typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const | ||
436 | { | ||
437 | if (scvf.isOnFracture()) | ||
438 | DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs " | ||
439 | "because the number of known corners is insufficient. " | ||
440 | "You can do this manually by extracting the corners from this scv " | ||
441 | "and extruding them by the corresponding aperture. "); | ||
442 | const GeometryHelper geometryHelper(element().geometry()); | ||
443 | const auto corners = geometryHelper.getScvfCorners(scvf.indexInElement()); | ||
444 | using ScvfGeometry = typename SubControlVolumeFace::Traits::Geometry; | ||
445 | return { Dune::GeometryTypes::cube(ScvfGeometry::mydimension), corners }; | ||
446 | } | ||
447 | |||
448 | private: | ||
449 | |||
450 | 2653260 | void makeElementGeometries_() | |
451 | { | ||
452 | // get the element geometry | ||
453 | 2653260 | const auto& element = *element_; | |
454 |
1/2✓ Branch 1 taken 1206588 times.
✗ Branch 2 not taken.
|
2653260 | const auto elementGeometry = element.geometry(); |
455 |
1/2✓ Branch 1 taken 1206588 times.
✗ Branch 2 not taken.
|
2653260 | const auto refElement = referenceElement(element); |
456 | |||
457 | // get the sub control volume geometries of this element | ||
458 | 2653260 | GeometryHelper geometryHelper(elementGeometry); | |
459 | |||
460 | // construct the sub control volumes | ||
461 |
1/2✓ Branch 1 taken 1206588 times.
✗ Branch 2 not taken.
|
2653260 | scvs_.resize(elementGeometry.corners()); |
462 | using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType; | ||
463 |
2/2✓ Branch 0 taken 8791080 times.
✓ Branch 1 taken 2653260 times.
|
11444340 | for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx) |
464 | { | ||
465 | // get associated dof index | ||
466 |
1/2✓ Branch 1 taken 3920192 times.
✗ Branch 2 not taken.
|
8791080 | const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim); |
467 | |||
468 | // add scv to the local container | ||
469 |
1/2✓ Branch 1 taken 3920192 times.
✗ Branch 2 not taken.
|
8791080 | scvs_[scvLocalIdx] = SubControlVolume(geometryHelper, |
470 | scvLocalIdx, | ||
471 | eIdx_, | ||
472 | dofIdxGlobal); | ||
473 | } | ||
474 | |||
475 | // construct the sub control volume faces | ||
476 | 2653260 | const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1); | |
477 |
1/2✓ Branch 1 taken 1206588 times.
✗ Branch 2 not taken.
|
2653260 | scvfs_.resize(numInnerScvf); |
478 | |||
479 | 2653260 | unsigned int scvfLocalIdx = 0; | |
480 |
2/2✓ Branch 0 taken 10173744 times.
✓ Branch 1 taken 2653260 times.
|
12827004 | for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx) |
481 | { | ||
482 | // find the local scv indices this scvf is connected to | ||
483 |
1/2✓ Branch 2 taken 10173744 times.
✗ Branch 3 not taken.
|
20347488 | std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)), |
484 |
1/2✓ Branch 2 taken 4381080 times.
✗ Branch 3 not taken.
|
10173744 | static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))}); |
485 | |||
486 |
1/2✓ Branch 1 taken 10173744 times.
✗ Branch 2 not taken.
|
10173744 | scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper, |
487 | element, | ||
488 | elementGeometry, | ||
489 | scvfLocalIdx, | ||
490 | std::move(localScvIndices)); | ||
491 | } | ||
492 | |||
493 | // construct the ... | ||
494 | // ... sub-control volume faces on the domain boundary | ||
495 | // ... sub-control volumes on fracture facets | ||
496 | // ... sub-control volume faces on fracture facets | ||
497 | // NOTE We do not construct fracture scvfs on boundaries here! | ||
498 | // That means specifying Neumann fluxes on only the fractures is not possible | ||
499 | // However, it is difficult to interpret this here as a fracture ending on the boundary | ||
500 | // could also be connected to a facet which is both boundary and fracture at the same time! | ||
501 | // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly | ||
502 | // we would have to find only those fractures that are at the boundary and aren't connected | ||
503 | // to a fracture which is a boundary. | ||
504 | 2653260 | LocalIndexType scvLocalIdx = element.subEntities(dim); | |
505 |
11/16✓ Branch 2 taken 1446672 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 1446672 times.
✓ Branch 6 taken 3920192 times.
✓ Branch 8 taken 3027336 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 6947528 times.
✓ Branch 12 taken 2304440 times.
✓ Branch 13 taken 9139900 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 3027336 times.
✗ Branch 17 not taken.
✓ Branch 1 taken 1206588 times.
✗ Branch 7 not taken.
✓ Branch 15 taken 1843552 times.
✓ Branch 10 taken 1843552 times.
|
40013616 | for (const auto& intersection : intersections(gridGeometry().gridView(), element)) |
506 | { | ||
507 | // first, obtain all vertex indices on this intersection | ||
508 |
1/2✓ Branch 1 taken 3920192 times.
✗ Branch 2 not taken.
|
8791080 | const auto& isGeometry = intersection.geometry(); |
509 |
1/2✓ Branch 1 taken 8791080 times.
✗ Branch 2 not taken.
|
8791080 | const auto numCorners = isGeometry.corners(); |
510 |
1/2✓ Branch 1 taken 8791080 times.
✗ Branch 2 not taken.
|
8791080 | const auto idxInInside = intersection.indexInInside(); |
511 | |||
512 |
1/4✓ Branch 1 taken 8791080 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
8791080 | std::vector<GridIndexType> isVertexIndices(numCorners); |
513 |
2/2✓ Branch 0 taken 20347488 times.
✓ Branch 1 taken 8791080 times.
|
29138568 | for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal) |
514 |
1/2✓ Branch 2 taken 20347488 times.
✗ Branch 3 not taken.
|
20347488 | isVertexIndices[vIdxLocal] = gridGeometry().vertexMapper().subIndex(element, |
515 | refElement.subEntity(idxInInside, 1, vIdxLocal, dim), | ||
516 | dim); | ||
517 | |||
518 |
2/2✓ Branch 0 taken 133372 times.
✓ Branch 1 taken 3786820 times.
|
9010652 | if (intersection.boundary()) |
519 | { | ||
520 |
2/2✓ Branch 0 taken 963468 times.
✓ Branch 1 taken 352944 times.
|
1316412 | for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx) |
521 | { | ||
522 | // find the scv this scvf is connected to | ||
523 | 963468 | const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim)); | |
524 |
1/2✓ Branch 1 taken 963468 times.
✗ Branch 2 not taken.
|
963468 | std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx}; |
525 | |||
526 |
1/2✓ Branch 1 taken 963468 times.
✗ Branch 2 not taken.
|
963468 | scvfs_.emplace_back(geometryHelper, |
527 | intersection, | ||
528 | isGeometry, | ||
529 | isScvfLocalIdx, | ||
530 | scvfLocalIdx, | ||
531 | std::move(localScvIndices)); | ||
532 | |||
533 | // increment local counter | ||
534 | 963468 | scvfLocalIdx++; | |
535 | } | ||
536 | } | ||
537 | |||
538 | // maybe add fracture scvs & scvfs | ||
539 |
3/4✓ Branch 1 taken 8791080 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 266736 times.
✓ Branch 4 taken 8524344 times.
|
8791080 | if (this->gridGeometry().isOnFracture(element, intersection)) |
540 | { | ||
541 | // add fracture scv for each vertex of intersection | ||
542 |
1/2✓ Branch 1 taken 266736 times.
✗ Branch 2 not taken.
|
266736 | const auto curNumScvs = scvs_.size(); |
543 |
1/2✓ Branch 1 taken 266736 times.
✗ Branch 2 not taken.
|
266736 | scvs_.reserve(curNumScvs+numCorners); |
544 |
2/2✓ Branch 0 taken 603432 times.
✓ Branch 1 taken 266736 times.
|
870168 | for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal) |
545 | 1810296 | scvs_.emplace_back(geometryHelper, | |
546 | intersection, | ||
547 | isGeometry, | ||
548 | vIdxLocal, | ||
549 |
1/2✓ Branch 1 taken 603432 times.
✗ Branch 2 not taken.
|
603432 | static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)), |
550 | scvLocalIdx++, | ||
551 | idxInInside, | ||
552 | 603432 | eIdx_, | |
553 | 603432 | isVertexIndices[vIdxLocal]); | |
554 | |||
555 | // add fracture scvf for each edge of the intersection in 3d | ||
556 | if (dim == 3) | ||
557 | { | ||
558 | 69960 | const auto& faceRefElement = referenceElement(isGeometry); | |
559 |
2/2✓ Branch 1 taken 209880 times.
✓ Branch 2 taken 69960 times.
|
279840 | for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx) |
560 | { | ||
561 | // inside/outside scv indices in face local node numbering | ||
562 | 419760 | std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)), | |
563 |
1/2✓ Branch 2 taken 209880 times.
✗ Branch 3 not taken.
|
209880 | static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))}); |
564 | |||
565 | // add offset to get the right scv indices | ||
566 | 209880 | std::for_each( localScvIndices.begin(), | |
567 | localScvIndices.end(), | ||
568 | 419760 | [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } ); | |
569 | |||
570 | // add scvf | ||
571 |
2/2✓ Branch 0 taken 24168 times.
✓ Branch 1 taken 115752 times.
|
209880 | scvfs_.emplace_back(geometryHelper, |
572 | intersection, | ||
573 | isGeometry, | ||
574 | edgeIdx, | ||
575 | scvfLocalIdx++, | ||
576 | std::move(localScvIndices), | ||
577 |
1/2✓ Branch 1 taken 209880 times.
✗ Branch 2 not taken.
|
209880 | intersection.boundary()); |
578 | } | ||
579 | } | ||
580 | |||
581 | // dim == 2, intersection is an edge, make 1 scvf | ||
582 | else | ||
583 | { | ||
584 | // inside/outside scv indices in face local node numbering | ||
585 |
1/2✓ Branch 1 taken 196776 times.
✗ Branch 2 not taken.
|
196776 | std::vector<LocalIndexType> localScvIndices({0, 1}); |
586 | |||
587 | // add offset such that the fracture scvs above are addressed | ||
588 | 196776 | std::for_each( localScvIndices.begin(), | |
589 | localScvIndices.end(), | ||
590 | 393552 | [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } ); | |
591 | |||
592 | // add scvf | ||
593 |
2/2✓ Branch 0 taken 1857 times.
✓ Branch 1 taken 97031 times.
|
196776 | scvfs_.emplace_back(geometryHelper, |
594 | intersection, | ||
595 | isGeometry, | ||
596 | 393552 | /*idxOnIntersection*/0, | |
597 | scvfLocalIdx++, | ||
598 | std::move(localScvIndices), | ||
599 |
1/2✓ Branch 1 taken 196776 times.
✗ Branch 2 not taken.
|
196776 | intersection.boundary()); |
600 | 196776 | } | |
601 | } | ||
602 | } | ||
603 | 2653260 | } | |
604 | |||
605 | //! The bound element | ||
606 | std::optional<Element> element_; | ||
607 | GridIndexType eIdx_; | ||
608 | |||
609 | //! The global geometry this is a restriction of | ||
610 | const GridGeometry* gridGeometryPtr_; | ||
611 | |||
612 | //! The element-local scvs & scvfs | ||
613 | std::vector<SubControlVolume> scvs_; | ||
614 | std::vector<SubControlVolumeFace> scvfs_; | ||
615 | }; | ||
616 | |||
617 | } // end namespace Dumux | ||
618 | |||
619 | #endif | ||
620 |