GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: dumux/dumux/porousmediumflow/boxdfm/fvelementgeometry.hh
Date: 2025-05-03 19:19:02
Exec Total Coverage
Lines: 89 89 100.0%
Functions: 17 17 100.0%
Branches: 79 138 57.2%

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 auto
102 scvs(const BoxDfmFVElementGeometry& fvGeometry, const LocalDof& localDof)
103 {
104 std::vector<std::size_t> indices;
105 for(std::size_t i = 0; i<fvGeometry.numScv(); i++)
106 if(fvGeometry.scv(i).localDofIndex() == localDof.index())
107 indices.push_back(fvGeometry.scv(i).indexInElement());
108
109 return Dune::transformedRangeView(
110 std::move(indices),
111 [&](const auto idx) { return fvGeometry.scv(idx); }
112 );
113 }
114
115 //! range over local dofs
116 template<class FVElementGeometry>
117 inline auto localDofs(const FVElementGeometry& fvGeometry)
118 {
119 return Dune::transformedRangeView(
120 Dune::range(numLocalDofs()),
121 [&](const auto i) { return CVFE::LocalDof
122 {
123 static_cast<LocalIndexType>(i),
124 static_cast<GridIndexType>(fvGeometry.scv(i).dofIndex()),
125 static_cast<GridIndexType>(fvGeometry.scv(i).elementIndex())
126 }; }
127 );
128 }
129
130 //! range over control-volume local dofs
131 template<class FVElementGeometry>
132 inline auto cvLocalDofs(const FVElementGeometry& fvGeometry)
133 {
134 // Default it that all dofs are cv dofs
135 return localDofs(fvGeometry);
136 }
137
138 /*!
139 * \brief Iterator range for sub control volumes faces.
140 *
141 * Iterates over all scvfs of the bound element.
142 * This is a free function found by means of ADL.
143 * To iterate over all sub control volume faces of this FVElementGeometry use
144 * for (auto&& scvf : scvfs(fvGeometry)).
145 */
146 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
147 scvfs(const BoxDfmFVElementGeometry& fvGeometry)
148 {
149 const auto& g = fvGeometry.gridGeometry();
150 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
151 return Dune::IteratorRange<Iter>(g.scvfs(fvGeometry.eIdx_).begin(), g.scvfs(fvGeometry.eIdx_).end());
152 }
153
154 //! Get a local finite element basis
155 const FeLocalBasis& feLocalBasis() const
156 { return gridGeometry().feCache().get(element_->type()).localBasis(); }
157
158 //! The total number of sub control volumes
159 std::size_t numScv() const
160 { return gridGeometry().scvs(eIdx_).size(); }
161
162 //! The total number of element-local dofs
163 std::size_t numLocalDofs() const
164 { return element().geometry().corners(); }
165
166 //! The total number of sub control volume faces
167 std::size_t numScvf() const
168 { return gridGeometry().scvfs(eIdx_).size(); }
169
170 /*!
171 * \brief bind the local view (r-value overload)
172 * This overload is called when an instance of this class is a temporary in the usage context
173 * This allows a usage like this: `const auto view = localView(...).bind(element);`
174 */
175 BoxDfmFVElementGeometry bind(const Element& element) &&
176 {
177 this->bindElement(element);
178 return std::move(*this);
179 }
180
181 //! This function is for compatibility reasons with cc methods
182 //! The box stencil is always element-local so bind and bindElement are identical.
183 void bind(const Element& element) &
184 { this->bindElement(element); }
185
186 /*!
187 * \brief bind the local view (r-value overload)
188 * This overload is called when an instance of this class is a temporary in the usage context
189 * This allows a usage like this: `const auto view = localView(...).bindElement(element);`
190 */
191 BoxDfmFVElementGeometry bindElement(const Element& element) &&
192 {
193 this->bindElement(element);
194 return std::move(*this);
195 }
196
197 /*!
198 * \brief Binding of an element, has to be called before using the fvgeometries
199 *
200 * Prepares all the volume variables within the element.
201 * For compatibility reasons with the FVGeometry cache being disabled.
202 */
203 void bindElement(const Element& element) &
204 {
205 element_ = element;
206 eIdx_ = gridGeometry().elementMapper().index(element);
207 }
208
209 //! The global finite volume geometry we are a restriction of
210 const GridGeometry& gridGeometry() const
211 { return *gridGeometryPtr_; }
212
213 //! Returns true if bind/bindElement has already been called
214 bool isBound() const
215 { return static_cast<bool>(element_); }
216
217 //! The bound element
218 const Element& element() const
219 { return *element_; }
220
221 //! Create the geometry of a given sub control volume
222 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
223 {
224 if (scv.isOnFracture())
225 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
226 "because the number of known corners is insufficient. "
227 "You can do this manually by extracting the corners from this scv "
228 "and extruding them by the corresponding aperture. ");
229
230 const typename GG::GeometryHelper geometryHelper(element().geometry());
231 const auto corners = geometryHelper.getScvCorners(scv.index());
232 using ScvGeometry = typename SubControlVolume::Traits::Geometry;
233 return { Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners };
234 }
235
236 //! Create the geometry of a given sub control volume face
237 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
238 {
239 if (scvf.isOnFracture())
240 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
241 "because the number of known corners is insufficient. "
242 "You can do this manually by extracting the corners from this scv "
243 "and extruding them by the corresponding aperture. ");
244 const typename GG::GeometryHelper geometryHelper(element().geometry());
245 const auto corners = geometryHelper.getScvfCorners(scvf.indexInElement());
246 using ScvfGeometry = typename SubControlVolumeFace::Traits::Geometry;
247 return { Dune::GeometryTypes::cube(ScvfGeometry::mydimension), corners };
248 }
249
250 private:
251 const GridGeometry* gridGeometryPtr_;
252
253 std::optional<Element> element_;
254 GridIndexType eIdx_;
255 };
256
257 //! Specialization in case the FVElementGeometries are not stored
258 template<class GG>
259 class BoxDfmFVElementGeometry<GG, false>
260 {
261 using GridView = typename GG::GridView;
262 static constexpr int dim = GridView::dimension;
263 static constexpr int dimWorld = GridView::dimensionworld;
264
265 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
266 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
267
268 using CoordScalar = typename GridView::ctype;
269 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
270 using GeometryHelper = typename GG::GeometryHelper;
271 public:
272 //! export type of the element
273 using Element = typename GridView::template Codim<0>::Entity;
274 //! Export type of subcontrol volume
275 using SubControlVolume = typename GG::SubControlVolume;
276 //! Export type of subcontrol volume face
277 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
278 //! Export type of finite volume grid geometry
279 using GridGeometry = GG;
280 //! The maximum number of scvs per element (2^dim for cubes)
281 //! multiplied by 3 for the maximum number of fracture scvs per vertex
282 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
283
284 //! Constructor
285 2405124 BoxDfmFVElementGeometry(const GridGeometry& gridGeometry)
286
3/12
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 7 taken 2220720 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 6 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 6 times.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
2405124 : gridGeometryPtr_(&gridGeometry) {}
287
288 //! Get a sub control volume with a local scv index
289 305201536 const SubControlVolume& scv(std::size_t scvIdx) const
290
6/10
✗ Branch 0 not taken.
✓ Branch 1 taken 72577584 times.
✓ Branch 2 taken 72577584 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 17348880 times.
✓ Branch 6 taken 46343640 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 25477008 times.
✓ Branch 11 taken 136328 times.
✗ Branch 12 not taken.
215275072 { return scvs_[scvIdx]; }
291
292 //! Get a sub control volume face with a local scvf index
293 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
294 { return scvfs_[scvfIdx]; }
295
296 /*!
297 * \brief Iterator range for sub control volumes.
298 *
299 * Iterates over all scvs of the bound element.
300 * This is a free function found by means of ADL.
301 * To iterate over all sub control volumes of this FVElementGeometry use
302 * for (auto&& scv : scvs(fvGeometry)).
303 */
304 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
305 107633054 scvs(const BoxDfmFVElementGeometry& fvGeometry)
306 {
307 using Iter = typename std::vector<SubControlVolume>::const_iterator;
308 107633054 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
309 }
310
311 //! range over sub control volumes related to a local dof.
312 template<class LocalDof>
313 friend inline auto
314 17348880 scvs(const BoxDfmFVElementGeometry& fvGeometry, const LocalDof& localDof)
315 {
316 17348880 std::vector<std::size_t> indices;
317 81041400 for(std::size_t i = 0; i<fvGeometry.numScv(); i++)
318
2/2
✓ Branch 0 taken 17348880 times.
✓ Branch 1 taken 46343640 times.
63692520 if(fvGeometry.scv(i).localDofIndex() == localDof.index())
319
1/2
✓ Branch 1 taken 17348880 times.
✗ Branch 2 not taken.
17348880 indices.push_back(fvGeometry.scv(i).indexInElement());
320
321 return Dune::transformedRangeView(
322 std::move(indices),
323
1/2
✓ Branch 1 taken 6939552 times.
✗ Branch 2 not taken.
34697760 [&](const auto idx) { return fvGeometry.scv(idx); }
324 17348880 );
325 17348880 }
326
327 //! range over local dofs
328 template<class FVElementGeometry>
329 inline auto localDofs(const FVElementGeometry& fvGeometry)
330 {
331 return Dune::transformedRangeView(
332 Dune::range(numLocalDofs()),
333 [&](const auto i) { return CVFE::LocalDof
334 {
335 static_cast<LocalIndexType>(i),
336 static_cast<GridIndexType>(fvGeometry.scv(i).dofIndex()),
337 static_cast<GridIndexType>(fvGeometry.scv(i).elementIndex())
338 }; }
339 );
340 }
341
342 //! range over control-volume local dofs
343 template<class FVElementGeometry>
344 inline auto cvLocalDofs(const FVElementGeometry& fvGeometry)
345 {
346 // Default it that all dofs are cv dofs
347 return localDofs(fvGeometry);
348 }
349
350 /*!
351 * \brief Iterator range for sub control volumes faces.
352 *
353 * Iterates over all scvfs of the bound element.
354 * This is a free function found by means of ADL.
355 * To iterate over all sub control volume faces of this FVElementGeometry use
356 * for (auto&& scvf : scvfs(fvGeometry)).
357 */
358 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
359 11380992 scvfs(const BoxDfmFVElementGeometry& fvGeometry)
360 {
361 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
362 11380992 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
363 }
364
365 //! Get a local finite element basis
366 8874360 const FeLocalBasis& feLocalBasis() const
367 8874360 { return gridGeometry().feCache().get(element_->type()).localBasis(); }
368
369 //! The total number of sub control volumes
370 134157656 std::size_t numScv() const
371
4/6
✓ Branch 0 taken 63692520 times.
✓ Branch 1 taken 17348880 times.
✓ Branch 4 taken 8874360 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 8874360 times.
✗ Branch 8 not taken.
134157656 { return scvs_.size(); }
372
373 //! The total number of element-local dofs
374 997572 std::size_t numLocalDofs() const
375 997572 { return element().geometry().corners(); }
376
377 //! The total number of sub control volume faces
378 2220720 std::size_t numScvf() const
379 2220720 { return scvfs_.size(); }
380
381 /*!
382 * \brief bind the local view (r-value overload)
383 * This overload is called when an instance of this class is a temporary in the usage context
384 * This allows a usage like this: `const auto view = localView(...).bind(element);`
385 */
386 BoxDfmFVElementGeometry bind(const Element& element) &&
387 {
388 this->bindElement(element);
389 return std::move(*this);
390 }
391
392 /*!
393 * \brief Binding of an element, has to be called before using the fvgeometries
394 * Prepares all the volume variables within the element.
395 * \note For the box scheme, bind() and bindElement() are identical, but the
396 * distinction is here for the sake of compatibility with cc schemes.
397 */
398 2229572 void bind(const Element& element) &
399
1/10
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 8 taken 8852 times.
✗ Branch 9 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
2229572 { this->bindElement(element); }
400
401 /*!
402 * \brief bind the local view (r-value overload)
403 * This overload is called when an instance of this class is a temporary in the usage context
404 * This allows a usage like this: `const auto view = localView(...).bindElement(element);`
405 */
406 BoxDfmFVElementGeometry bindElement(const Element& element) &&
407 {
408 this->bindElement(element);
409 return std::move(*this);
410 }
411
412 /*!
413 * \brief Binding of an element, has to be called before using the fvgeometries
414 * Prepares all the volume variables within the element.
415 */
416 2422816 void bindElement(const Element& element) &
417 {
418
2/2
✓ Branch 0 taken 8846 times.
✓ Branch 1 taken 1207382 times.
2422816 element_ = element;
419 4845632 eIdx_ = gridGeometry().elementMapper().index(element);
420 2422816 makeElementGeometries_();
421 2422816 }
422
423 //! The global finite volume geometry we are a restriction of
424 60327244 const GridGeometry& gridGeometry() const
425
3/4
✓ Branch 0 taken 1020605 times.
✓ Branch 1 taken 6642475 times.
✓ Branch 4 taken 4421880 times.
✗ Branch 5 not taken.
24583660 { return *gridGeometryPtr_; }
426
427 //! Returns true if bind/bindElement has already been called
428 bool isBound() const
429 { return static_cast<bool>(element_); }
430
431 //! The bound element
432 6939552 const Element& element() const
433
1/2
✓ Branch 2 taken 6939552 times.
✗ Branch 3 not taken.
7937124 { return *element_; }
434
435 //! Create the geometry of a given sub control volume
436 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
437 {
438 if (scv.isOnFracture())
439 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
440 "because the number of known corners is insufficient. "
441 "You can do this manually by extracting the corners from this scv "
442 "and extruding them by the corresponding aperture. ");
443
444 const GeometryHelper geometryHelper(element().geometry());
445 const auto corners = geometryHelper.getScvCorners(scv.index());
446 using ScvGeometry = typename SubControlVolume::Traits::Geometry;
447 return { Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners };
448 }
449
450 //! Create the geometry of a given sub control volume face
451 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
452 {
453 if (scvf.isOnFracture())
454 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
455 "because the number of known corners is insufficient. "
456 "You can do this manually by extracting the corners from this scv "
457 "and extruding them by the corresponding aperture. ");
458 const GeometryHelper geometryHelper(element().geometry());
459 const auto corners = geometryHelper.getScvfCorners(scvf.indexInElement());
460 using ScvfGeometry = typename SubControlVolumeFace::Traits::Geometry;
461 return { Dune::GeometryTypes::cube(ScvfGeometry::mydimension), corners };
462 }
463
464 private:
465
466 2422816 void makeElementGeometries_()
467 {
468 // get the element geometry
469 2422816 const auto& element = *element_;
470
1/2
✓ Branch 1 taken 1206588 times.
✗ Branch 2 not taken.
2422816 const auto elementGeometry = element.geometry();
471
1/2
✓ Branch 1 taken 1206588 times.
✗ Branch 2 not taken.
2422816 const auto refElement = referenceElement(element);
472
473 // get the sub control volume geometries of this element
474 2422816 GeometryHelper geometryHelper(elementGeometry);
475
476 // construct the sub control volumes
477
1/2
✓ Branch 1 taken 1206588 times.
✗ Branch 2 not taken.
2422816 scvs_.resize(elementGeometry.corners());
478 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
479
2/2
✓ Branch 0 taken 7869304 times.
✓ Branch 1 taken 2422816 times.
10292120 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
480 {
481 // get associated dof index
482
1/2
✓ Branch 1 taken 3920192 times.
✗ Branch 2 not taken.
7869304 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
483
484 // add scv to the local container
485
1/2
✓ Branch 1 taken 3920192 times.
✗ Branch 2 not taken.
7869304 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
486 scvLocalIdx,
487 eIdx_,
488 dofIdxGlobal);
489 }
490
491 // construct the sub control volume faces
492 2422816 const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1);
493
1/2
✓ Branch 1 taken 1206588 times.
✗ Branch 2 not taken.
2422816 scvfs_.resize(numInnerScvf);
494
495 2422816 unsigned int scvfLocalIdx = 0;
496
2/2
✓ Branch 0 taken 8791080 times.
✓ Branch 1 taken 2422816 times.
11213896 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
497 {
498 // find the local scv indices this scvf is connected to
499
1/2
✓ Branch 2 taken 8791080 times.
✗ Branch 3 not taken.
17582160 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
500
1/2
✓ Branch 2 taken 4381080 times.
✗ Branch 3 not taken.
8791080 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
501
502
1/2
✓ Branch 1 taken 8791080 times.
✗ Branch 2 not taken.
8791080 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
503 element,
504 elementGeometry,
505 scvfLocalIdx,
506 std::move(localScvIndices));
507 }
508
509 // construct the ...
510 // ... sub-control volume faces on the domain boundary
511 // ... sub-control volumes on fracture facets
512 // ... sub-control volume faces on fracture facets
513 // NOTE We do not construct fracture scvfs on boundaries here!
514 // That means specifying Neumann fluxes on only the fractures is not possible
515 // However, it is difficult to interpret this here as a fracture ending on the boundary
516 // could also be connected to a facet which is both boundary and fracture at the same time!
517 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
518 // we would have to find only those fractures that are at the boundary and aren't connected
519 // to a fracture which is a boundary.
520 2422816 LocalIndexType scvLocalIdx = element.subEntities(dim);
521
11/16
✓ Branch 2 taken 1216228 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 1216228 times.
✓ Branch 6 taken 3920192 times.
✓ Branch 8 taken 3027336 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 6947528 times.
✓ Branch 12 taken 1152220 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 921776 times.
✓ Branch 10 taken 921776 times.
36326512 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
522 {
523 // first, obtain all vertex indices on this intersection
524
1/2
✓ Branch 1 taken 3920192 times.
✗ Branch 2 not taken.
7869304 const auto& isGeometry = intersection.geometry();
525
1/2
✓ Branch 1 taken 7869304 times.
✗ Branch 2 not taken.
7869304 const auto numCorners = isGeometry.corners();
526
1/2
✓ Branch 1 taken 7869304 times.
✗ Branch 2 not taken.
7869304 const auto idxInInside = intersection.indexInInside();
527
528
1/4
✓ Branch 1 taken 7869304 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
7869304 std::vector<GridIndexType> isVertexIndices(numCorners);
529
2/2
✓ Branch 0 taken 17582160 times.
✓ Branch 1 taken 7869304 times.
25451464 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
530
1/2
✓ Branch 2 taken 17582160 times.
✗ Branch 3 not taken.
17582160 isVertexIndices[vIdxLocal] = gridGeometry().vertexMapper().subIndex(element,
531 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
532 dim);
533
534
2/2
✓ Branch 0 taken 133372 times.
✓ Branch 1 taken 3786820 times.
8003016 if (intersection.boundary())
535 {
536
2/2
✓ Branch 0 taken 705888 times.
✓ Branch 1 taken 267084 times.
972972 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
537 {
538 // find the scv this scvf is connected to
539 705888 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
540
1/2
✓ Branch 1 taken 705888 times.
✗ Branch 2 not taken.
705888 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
541
542
1/2
✓ Branch 1 taken 705888 times.
✗ Branch 2 not taken.
705888 scvfs_.emplace_back(geometryHelper,
543 intersection,
544 isGeometry,
545 isScvfLocalIdx,
546 scvfLocalIdx,
547 std::move(localScvIndices));
548
549 // increment local counter
550 705888 scvfLocalIdx++;
551 }
552 }
553
554 // maybe add fracture scvs & scvfs
555
3/4
✓ Branch 1 taken 7869304 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 243416 times.
✓ Branch 4 taken 7625888 times.
7869304 if (this->gridGeometry().isOnFracture(element, intersection))
556 {
557 // add fracture scv for each vertex of intersection
558
1/2
✓ Branch 1 taken 243416 times.
✗ Branch 2 not taken.
243416 const auto curNumScvs = scvs_.size();
559
1/2
✓ Branch 1 taken 243416 times.
✗ Branch 2 not taken.
243416 scvs_.reserve(curNumScvs+numCorners);
560
2/2
✓ Branch 0 taken 533472 times.
✓ Branch 1 taken 243416 times.
776888 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
561 1600416 scvs_.emplace_back(geometryHelper,
562 intersection,
563 isGeometry,
564 vIdxLocal,
565
1/2
✓ Branch 1 taken 533472 times.
✗ Branch 2 not taken.
533472 static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
566 scvLocalIdx++,
567 idxInInside,
568 533472 eIdx_,
569 533472 isVertexIndices[vIdxLocal]);
570
571 // add fracture scvf for each edge of the intersection in 3d
572 if (dim == 3)
573 {
574 46640 const auto& faceRefElement = referenceElement(isGeometry);
575
2/2
✓ Branch 1 taken 139920 times.
✓ Branch 2 taken 46640 times.
186560 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
576 {
577 // inside/outside scv indices in face local node numbering
578 279840 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
579
1/2
✓ Branch 2 taken 139920 times.
✗ Branch 3 not taken.
139920 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
580
581 // add offset to get the right scv indices
582 139920 std::for_each( localScvIndices.begin(),
583 localScvIndices.end(),
584 279840 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
585
586 // add scvf
587
2/2
✓ Branch 0 taken 12084 times.
✓ Branch 1 taken 57876 times.
139920 scvfs_.emplace_back(geometryHelper,
588 intersection,
589 isGeometry,
590 edgeIdx,
591 scvfLocalIdx++,
592 std::move(localScvIndices),
593
1/2
✓ Branch 1 taken 139920 times.
✗ Branch 2 not taken.
139920 intersection.boundary());
594 }
595 }
596
597 // dim == 2, intersection is an edge, make 1 scvf
598 else
599 {
600 // inside/outside scv indices in face local node numbering
601
1/2
✓ Branch 1 taken 196776 times.
✗ Branch 2 not taken.
196776 std::vector<LocalIndexType> localScvIndices({0, 1});
602
603 // add offset such that the fracture scvs above are addressed
604 196776 std::for_each( localScvIndices.begin(),
605 localScvIndices.end(),
606 393552 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
607
608 // add scvf
609
2/2
✓ Branch 0 taken 1857 times.
✓ Branch 1 taken 97031 times.
196776 scvfs_.emplace_back(geometryHelper,
610 intersection,
611 isGeometry,
612 393552 /*idxOnIntersection*/0,
613 scvfLocalIdx++,
614 std::move(localScvIndices),
615
1/2
✓ Branch 1 taken 196776 times.
✗ Branch 2 not taken.
196776 intersection.boundary());
616 196776 }
617 }
618 }
619 2422816 }
620
621 //! The bound element
622 std::optional<Element> element_;
623 GridIndexType eIdx_;
624
625 //! The global geometry this is a restriction of
626 const GridGeometry* gridGeometryPtr_;
627
628 //! The element-local scvs & scvfs
629 std::vector<SubControlVolume> scvs_;
630 std::vector<SubControlVolumeFace> scvfs_;
631 };
632
633 } // end namespace Dumux
634
635 #endif
636