GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: dumux/dumux/porousmediumflow/boxdfm/fvelementgeometry.hh
Date: 2025-04-12 19:19:20
Exec Total Coverage
Lines: 79 79 100.0%
Functions: 12 12 100.0%
Branches: 68 122 55.7%

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
23 #include <dune/geometry/type.hh>
24 #include <dune/localfunctions/lagrange/pqkfactory.hh>
25
26 #include <dumux/discretization/scvandscvfiterators.hh>
27 #include "geometryhelper.hh"
28
29 namespace Dumux {
30
31 /*!
32 * \ingroup BoxDFMModel
33 * \brief Base class for the finite volume geometry vector for box discrete fracture model.
34 *
35 * This builds up the sub control volumes and sub control volume faces for each element.
36 *
37 * \tparam GG the finite volume grid geometry type
38 * \tparam enableGridGeometryCache if the grid geometry is cached or not
39 */
40 template<class GG, bool enableGridGeometryCache>
41 class BoxDfmFVElementGeometry;
42
43 //! Specialization in case the FVElementGeometries are stored
44 template<class GG>
45 class BoxDfmFVElementGeometry<GG, true>
46 {
47 using GridView = typename GG::GridView;
48 static constexpr int dim = GridView::dimension;
49 static constexpr int dimWorld = GridView::dimensionworld;
50 using GridIndexType = typename GridView::IndexSet::IndexType;
51 using CoordScalar = typename GridView::ctype;
52 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
53 public:
54 //! export type of the element
55 using Element = typename GridView::template Codim<0>::Entity;
56 //! Export type of subcontrol volume
57 using SubControlVolume = typename GG::SubControlVolume;
58 //! Export type of subcontrol volume face
59 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
60 //! Export type of finite volume grid geometry
61 using GridGeometry = GG;
62
63 //! The maximum number of scvs per element (2^dim for cubes)
64 //! multiplied by 3 for the maximum number of fracture scvs per vertex
65 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
66
67 //! Constructor
68 BoxDfmFVElementGeometry(const GridGeometry& gridGeometry)
69 : gridGeometryPtr_(&gridGeometry) {}
70
71 //! Get a sub control volume with a local scv index
72 const SubControlVolume& scv(std::size_t scvIdx) const
73 { return gridGeometry().scvs(eIdx_)[scvIdx]; }
74
75 //! Get a sub control volume face with a local scvf index
76 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
77 { return gridGeometry().scvfs(eIdx_)[scvfIdx]; }
78
79 /*!
80 * \brief Iterator range for sub control volumes.
81 *
82 * Iterates over all scvs of the bound element.
83 * This is a free function found by means of ADL.
84 * To iterate over all sub control volumes of this FVElementGeometry use
85 * for (auto&& scv : scvs(fvGeometry)).
86 */
87 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
88 scvs(const BoxDfmFVElementGeometry& fvGeometry)
89 {
90 const auto& g = fvGeometry.gridGeometry();
91 using Iter = typename std::vector<SubControlVolume>::const_iterator;
92 return Dune::IteratorRange<Iter>(g.scvs(fvGeometry.eIdx_).begin(), g.scvs(fvGeometry.eIdx_).end());
93 }
94
95 /*!
96 * \brief Iterator range for sub control volumes faces.
97 *
98 * Iterates over all scvfs of the bound element.
99 * This is a free function found by means of ADL.
100 * To iterate over all sub control volume faces of this FVElementGeometry use
101 * for (auto&& scvf : scvfs(fvGeometry)).
102 */
103 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
104 scvfs(const BoxDfmFVElementGeometry& fvGeometry)
105 {
106 const auto& g = fvGeometry.gridGeometry();
107 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
108 return Dune::IteratorRange<Iter>(g.scvfs(fvGeometry.eIdx_).begin(), g.scvfs(fvGeometry.eIdx_).end());
109 }
110
111 //! Get a local finite element basis
112 const FeLocalBasis& feLocalBasis() const
113 { return gridGeometry().feCache().get(element_->type()).localBasis(); }
114
115 //! The total number of sub control volumes
116 std::size_t numScv() const
117 { return gridGeometry().scvs(eIdx_).size(); }
118
119 //! The total number of sub control volume faces
120 std::size_t numScvf() const
121 { return gridGeometry().scvfs(eIdx_).size(); }
122
123 /*!
124 * \brief bind the local view (r-value overload)
125 * This overload is called when an instance of this class is a temporary in the usage context
126 * This allows a usage like this: `const auto view = localView(...).bind(element);`
127 */
128 BoxDfmFVElementGeometry bind(const Element& element) &&
129 {
130 this->bindElement(element);
131 return std::move(*this);
132 }
133
134 //! This function is for compatibility reasons with cc methods
135 //! The box stencil is always element-local so bind and bindElement are identical.
136 void bind(const Element& element) &
137 { this->bindElement(element); }
138
139 /*!
140 * \brief bind the local view (r-value overload)
141 * This overload is called when an instance of this class is a temporary in the usage context
142 * This allows a usage like this: `const auto view = localView(...).bindElement(element);`
143 */
144 BoxDfmFVElementGeometry bindElement(const Element& element) &&
145 {
146 this->bindElement(element);
147 return std::move(*this);
148 }
149
150 /*!
151 * \brief Binding of an element, has to be called before using the fvgeometries
152 *
153 * Prepares all the volume variables within the element.
154 * For compatibility reasons with the FVGeometry cache being disabled.
155 */
156 void bindElement(const Element& element) &
157 {
158 element_ = element;
159 eIdx_ = gridGeometry().elementMapper().index(element);
160 }
161
162 //! The global finite volume geometry we are a restriction of
163 const GridGeometry& gridGeometry() const
164 { return *gridGeometryPtr_; }
165
166 //! Returns true if bind/bindElement has already been called
167 bool isBound() const
168 { return static_cast<bool>(element_); }
169
170 //! The bound element
171 const Element& element() const
172 { return *element_; }
173
174 //! Create the geometry of a given sub control volume
175 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
176 {
177 if (scv.isOnFracture())
178 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
179 "because the number of known corners is insufficient. "
180 "You can do this manually by extracting the corners from this scv "
181 "and extruding them by the corresponding aperture. ");
182
183 const typename GG::GeometryHelper geometryHelper(element().geometry());
184 const auto corners = geometryHelper.getScvCorners(scv.index());
185 using ScvGeometry = typename SubControlVolume::Traits::Geometry;
186 return { Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners };
187 }
188
189 //! Create the geometry of a given sub control volume face
190 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
191 {
192 if (scvf.isOnFracture())
193 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
194 "because the number of known corners is insufficient. "
195 "You can do this manually by extracting the corners from this scv "
196 "and extruding them by the corresponding aperture. ");
197 const typename GG::GeometryHelper geometryHelper(element().geometry());
198 const auto corners = geometryHelper.getScvfCorners(scvf.indexInElement());
199 using ScvfGeometry = typename SubControlVolumeFace::Traits::Geometry;
200 return { Dune::GeometryTypes::cube(ScvfGeometry::mydimension), corners };
201 }
202
203 private:
204 const GridGeometry* gridGeometryPtr_;
205
206 std::optional<Element> element_;
207 GridIndexType eIdx_;
208 };
209
210 //! Specialization in case the FVElementGeometries are not stored
211 template<class GG>
212 class BoxDfmFVElementGeometry<GG, false>
213 {
214 using GridView = typename GG::GridView;
215 static constexpr int dim = GridView::dimension;
216 static constexpr int dimWorld = GridView::dimensionworld;
217
218 using GridIndexType = typename GridView::IndexSet::IndexType;
219
220 using CoordScalar = typename GridView::ctype;
221 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
222 using GeometryHelper = typename GG::GeometryHelper;
223 public:
224 //! export type of the element
225 using Element = typename GridView::template Codim<0>::Entity;
226 //! Export type of subcontrol volume
227 using SubControlVolume = typename GG::SubControlVolume;
228 //! Export type of subcontrol volume face
229 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
230 //! Export type of finite volume grid geometry
231 using GridGeometry = GG;
232 //! The maximum number of scvs per element (2^dim for cubes)
233 //! multiplied by 3 for the maximum number of fracture scvs per vertex
234 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
235
236 //! Constructor
237 2377024 BoxDfmFVElementGeometry(const GridGeometry& gridGeometry)
238
3/12
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✓ Branch 7 taken 2119356 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.
2377024 : gridGeometryPtr_(&gridGeometry) {}
239
240 //! Get a sub control volume with a local scv index
241 298479260 const SubControlVolume& scv(std::size_t scvIdx) const
242
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 116749288 times.
✓ Branch 2 taken 116749288 times.
✗ Branch 3 not taken.
181729972 { return scvs_[scvIdx]; }
243
244 //! Get a sub control volume face with a local scvf index
245 const SubControlVolumeFace& scvf(std::size_t scvfIdx) const
246 { return scvfs_[scvfIdx]; }
247
248 /*!
249 * \brief Iterator range for sub control volumes.
250 *
251 * Iterates over all scvs of the bound element.
252 * This is a free function found by means of ADL.
253 * To iterate over all sub control volumes of this FVElementGeometry use
254 * for (auto&& scv : scvs(fvGeometry)).
255 */
256 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
257 182946758 scvs(const BoxDfmFVElementGeometry& fvGeometry)
258 {
259 using Iter = typename std::vector<SubControlVolume>::const_iterator;
260 182946758 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
261 }
262
263 /*!
264 * \brief Iterator range for sub control volumes faces.
265 *
266 * Iterates over all scvfs of the bound element.
267 * This is a free function found by means of ADL.
268 * To iterate over all sub control volume faces of this FVElementGeometry use
269 * for (auto&& scvf : scvfs(fvGeometry)).
270 */
271 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
272 18391160 scvfs(const BoxDfmFVElementGeometry& fvGeometry)
273 {
274 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
275 18391160 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
276 }
277
278 //! Get a local finite element basis
279 7675224 const FeLocalBasis& feLocalBasis() const
280 7675224 { return gridGeometry().feCache().get(element_->type()).localBasis(); }
281
282 //! The total number of sub control volumes
283 56629136 std::size_t numScv() const
284
2/4
✓ Branch 3 taken 7675224 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 7675224 times.
✗ Branch 7 not taken.
56629136 { return scvs_.size(); }
285
286 //! The total number of sub control volume faces
287 2119356 std::size_t numScvf() const
288 2119356 { return scvfs_.size(); }
289
290 /*!
291 * \brief bind the local view (r-value overload)
292 * This overload is called when an instance of this class is a temporary in the usage context
293 * This allows a usage like this: `const auto view = localView(...).bind(element);`
294 */
295 BoxDfmFVElementGeometry bind(const Element& element) &&
296 {
297 this->bindElement(element);
298 return std::move(*this);
299 }
300
301 /*!
302 * \brief Binding of an element, has to be called before using the fvgeometries
303 * Prepares all the volume variables within the element.
304 * \note For the box scheme, bind() and bindElement() are identical, but the
305 * distinction is here for the sake of compatibility with cc schemes.
306 */
307 2128208 void bind(const Element& element) &
308
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.
2128208 { this->bindElement(element); }
309
310 /*!
311 * \brief bind the local view (r-value overload)
312 * This overload is called when an instance of this class is a temporary in the usage context
313 * This allows a usage like this: `const auto view = localView(...).bindElement(element);`
314 */
315 BoxDfmFVElementGeometry bindElement(const Element& element) &&
316 {
317 this->bindElement(element);
318 return std::move(*this);
319 }
320
321 /*!
322 * \brief Binding of an element, has to be called before using the fvgeometries
323 * Prepares all the volume variables within the element.
324 */
325 2394716 void bindElement(const Element& element) &
326 {
327
2/2
✓ Branch 0 taken 8846 times.
✓ Branch 1 taken 1188512 times.
2394716 element_ = element;
328 4789432 eIdx_ = gridGeometry().elementMapper().index(element);
329 2394716 makeElementGeometries_();
330 2394716 }
331
332 //! The global finite volume geometry we are a restriction of
333 57393036 const GridGeometry& gridGeometry() const
334
3/4
✓ Branch 0 taken 643678 times.
✓ Branch 1 taken 6432546 times.
✓ Branch 4 taken 3837612 times.
✗ Branch 5 not taken.
23761888 { return *gridGeometryPtr_; }
335
336 //! Returns true if bind/bindElement has already been called
337 bool isBound() const
338 { return static_cast<bool>(element_); }
339
340 //! The bound element
341 14152448 const Element& element() const
342 14152448 { return *element_; }
343
344 //! Create the geometry of a given sub control volume
345 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
346 {
347 if (scv.isOnFracture())
348 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
349 "because the number of known corners is insufficient. "
350 "You can do this manually by extracting the corners from this scv "
351 "and extruding them by the corresponding aperture. ");
352
353 const GeometryHelper geometryHelper(element().geometry());
354 const auto corners = geometryHelper.getScvCorners(scv.index());
355 using ScvGeometry = typename SubControlVolume::Traits::Geometry;
356 return { Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners };
357 }
358
359 //! Create the geometry of a given sub control volume face
360 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
361 {
362 if (scvf.isOnFracture())
363 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
364 "because the number of known corners is insufficient. "
365 "You can do this manually by extracting the corners from this scv "
366 "and extruding them by the corresponding aperture. ");
367 const GeometryHelper geometryHelper(element().geometry());
368 const auto corners = geometryHelper.getScvfCorners(scvf.indexInElement());
369 using ScvfGeometry = typename SubControlVolumeFace::Traits::Geometry;
370 return { Dune::GeometryTypes::cube(ScvfGeometry::mydimension), corners };
371 }
372
373 private:
374
375 2394716 void makeElementGeometries_()
376 {
377 // get the element geometry
378 2394716 const auto& element = *element_;
379
1/2
✓ Branch 1 taken 1197358 times.
✗ Branch 2 not taken.
2394716 const auto elementGeometry = element.geometry();
380
1/2
✓ Branch 1 taken 1197358 times.
✗ Branch 2 not taken.
2394716 const auto refElement = referenceElement(element);
381
382 // get the sub control volume geometries of this element
383 2394716 GeometryHelper geometryHelper(elementGeometry);
384
385 // construct the sub control volumes
386
1/2
✓ Branch 1 taken 1197358 times.
✗ Branch 2 not taken.
2394716 scvs_.resize(elementGeometry.corners());
387 using LocalIndexType = typename SubControlVolumeFace::Traits::LocalIndexType;
388
2/2
✓ Branch 0 taken 7535184 times.
✓ Branch 1 taken 2394716 times.
9929900 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
389 {
390 // get associated dof index
391
1/2
✓ Branch 1 taken 3767592 times.
✗ Branch 2 not taken.
7535184 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
392
393 // add scv to the local container
394
1/2
✓ Branch 1 taken 3767592 times.
✗ Branch 2 not taken.
7535184 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
395 scvLocalIdx,
396 eIdx_,
397 dofIdxGlobal);
398 }
399
400 // construct the sub control volume faces
401 2394716 const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1);
402
1/2
✓ Branch 1 taken 1197358 times.
✗ Branch 2 not taken.
2394716 scvfs_.resize(numInnerScvf);
403
404 2394716 unsigned int scvfLocalIdx = 0;
405
2/2
✓ Branch 0 taken 8083032 times.
✓ Branch 1 taken 2394716 times.
10477748 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
406 {
407 // find the local scv indices this scvf is connected to
408
1/2
✓ Branch 2 taken 8083032 times.
✗ Branch 3 not taken.
16166064 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
409
1/2
✓ Branch 2 taken 4041516 times.
✗ Branch 3 not taken.
8083032 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))});
410
411
1/2
✓ Branch 1 taken 8083032 times.
✗ Branch 2 not taken.
8083032 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
412 element,
413 elementGeometry,
414 scvfLocalIdx,
415 std::move(localScvIndices));
416 }
417
418 // construct the ...
419 // ... sub-control volume faces on the domain boundary
420 // ... sub-control volumes on fracture facets
421 // ... sub-control volume faces on fracture facets
422 // NOTE We do not construct fracture scvfs on boundaries here!
423 // That means specifying Neumann fluxes on only the fractures is not possible
424 // However, it is difficult to interpret this here as a fracture ending on the boundary
425 // could also be connected to a facet which is both boundary and fracture at the same time!
426 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
427 // we would have to find only those fractures that are at the boundary and aren't connected
428 // to a fracture which is a boundary.
429 2394716 LocalIndexType scvLocalIdx = element.subEntities(dim);
430
11/16
✓ Branch 2 taken 1197358 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 1197358 times.
✓ Branch 6 taken 3767592 times.
✓ Branch 8 taken 3219744 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 6987336 times.
✓ Branch 12 taken 684810 times.
✓ Branch 13 taken 9245090 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 3219744 times.
✗ Branch 17 not taken.
✓ Branch 1 taken 1197358 times.
✗ Branch 7 not taken.
✓ Branch 15 taken 547848 times.
✓ Branch 10 taken 547848 times.
35404160 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
431 {
432 // first, obtain all vertex indices on this intersection
433
1/2
✓ Branch 1 taken 3767592 times.
✗ Branch 2 not taken.
7535184 const auto& isGeometry = intersection.geometry();
434
1/2
✓ Branch 1 taken 7535184 times.
✗ Branch 2 not taken.
7535184 const auto numCorners = isGeometry.corners();
435
1/2
✓ Branch 1 taken 7535184 times.
✗ Branch 2 not taken.
7535184 const auto idxInInside = intersection.indexInInside();
436
437
1/4
✓ Branch 1 taken 7535184 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
7535184 std::vector<GridIndexType> isVertexIndices(numCorners);
438
2/2
✓ Branch 0 taken 16166064 times.
✓ Branch 1 taken 7535184 times.
23701248 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
439
1/2
✓ Branch 2 taken 16166064 times.
✗ Branch 3 not taken.
16166064 isVertexIndices[vIdxLocal] = gridGeometry().vertexMapper().subIndex(element,
440 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
441 dim);
442
443
2/2
✓ Branch 0 taken 95638 times.
✓ Branch 1 taken 3671954 times.
7630822 if (intersection.boundary())
444 {
445
2/2
✓ Branch 0 taken 484612 times.
✓ Branch 1 taken 191276 times.
675888 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
446 {
447 // find the scv this scvf is connected to
448 484612 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
449
1/2
✓ Branch 1 taken 484612 times.
✗ Branch 2 not taken.
484612 std::vector<LocalIndexType> localScvIndices = {insideScvIdx, insideScvIdx};
450
451
1/2
✓ Branch 1 taken 484612 times.
✗ Branch 2 not taken.
484612 scvfs_.emplace_back(geometryHelper,
452 intersection,
453 isGeometry,
454 isScvfLocalIdx,
455 scvfLocalIdx,
456 std::move(localScvIndices));
457
458 // increment local counter
459 484612 scvfLocalIdx++;
460 }
461 }
462
463 // maybe add fracture scvs & scvfs
464
3/4
✓ Branch 1 taken 7535184 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 244004 times.
✓ Branch 4 taken 7291180 times.
7535184 if (this->gridGeometry().isOnFracture(element, intersection))
465 {
466 // add fracture scv for each vertex of intersection
467
1/2
✓ Branch 1 taken 244004 times.
✗ Branch 2 not taken.
244004 const auto curNumScvs = scvs_.size();
468
1/2
✓ Branch 1 taken 244004 times.
✗ Branch 2 not taken.
244004 scvs_.reserve(curNumScvs+numCorners);
469
2/2
✓ Branch 0 taken 515728 times.
✓ Branch 1 taken 244004 times.
759732 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
470 1547184 scvs_.emplace_back(geometryHelper,
471 intersection,
472 isGeometry,
473 vIdxLocal,
474
1/2
✓ Branch 1 taken 515728 times.
✗ Branch 2 not taken.
515728 static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
475 scvLocalIdx++,
476 idxInInside,
477 515728 eIdx_,
478 515728 isVertexIndices[vIdxLocal]);
479
480 // add fracture scvf for each edge of the intersection in 3d
481 if (dim == 3)
482 {
483 27720 const auto& faceRefElement = referenceElement(isGeometry);
484
2/2
✓ Branch 1 taken 83160 times.
✓ Branch 2 taken 27720 times.
110880 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
485 {
486 // inside/outside scv indices in face local node numbering
487 166320 std::vector<LocalIndexType> localScvIndices({static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
488
1/2
✓ Branch 2 taken 83160 times.
✗ Branch 3 not taken.
83160 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))});
489
490 // add offset to get the right scv indices
491 83160 std::for_each( localScvIndices.begin(),
492 localScvIndices.end(),
493 166320 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
494
495 // add scvf
496
2/2
✓ Branch 0 taken 7182 times.
✓ Branch 1 taken 34398 times.
83160 scvfs_.emplace_back(geometryHelper,
497 intersection,
498 isGeometry,
499 edgeIdx,
500 scvfLocalIdx++,
501 std::move(localScvIndices),
502
1/2
✓ Branch 1 taken 83160 times.
✗ Branch 2 not taken.
83160 intersection.boundary());
503 }
504 }
505
506 // dim == 2, intersection is an edge, make 1 scvf
507 else
508 {
509 // inside/outside scv indices in face local node numbering
510
1/2
✓ Branch 1 taken 216284 times.
✗ Branch 2 not taken.
216284 std::vector<LocalIndexType> localScvIndices({0, 1});
511
512 // add offset such that the fracture scvs above are addressed
513 216284 std::for_each( localScvIndices.begin(),
514 localScvIndices.end(),
515 432568 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
516
517 // add scvf
518
2/2
✓ Branch 0 taken 1828 times.
✓ Branch 1 taken 106314 times.
216284 scvfs_.emplace_back(geometryHelper,
519 intersection,
520 isGeometry,
521 432568 /*idxOnIntersection*/0,
522 scvfLocalIdx++,
523 std::move(localScvIndices),
524
1/2
✓ Branch 1 taken 216284 times.
✗ Branch 2 not taken.
216284 intersection.boundary());
525 216284 }
526 }
527 }
528 2394716 }
529
530 //! The bound element
531 std::optional<Element> element_;
532 GridIndexType eIdx_;
533
534 //! The global geometry this is a restriction of
535 const GridGeometry* gridGeometryPtr_;
536
537 //! The element-local scvs & scvfs
538 std::vector<SubControlVolume> scvs_;
539 std::vector<SubControlVolumeFace> scvfs_;
540 };
541
542 } // end namespace Dumux
543
544 #endif
545