version 3.11-dev
porousmediumflow/boxdfm/fvelementgeometry.hh
Go to the documentation of this file.
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//
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
30#include "geometryhelper.hh"
31
32namespace Dumux {
33
43template<class GG, bool enableGridGeometryCache>
45
47template<class GG>
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 using GGCache = typename GG::Cache;
58 using GeometryHelper = typename GGCache::GeometryHelper;
59public:
61 using Element = typename GridView::template Codim<0>::Entity;
63 using SubControlVolume = typename GG::SubControlVolume;
65 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
67 using GridGeometry = GG;
68
71 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
72
74 [[deprecated("This Constructor is deprecated and will be removed after release 3.11. Always use localView(gridGeometry).")]]
76 : BoxDfmFVElementGeometry(gridGeometry.cache_) {}
77
79 BoxDfmFVElementGeometry(const GGCache& ggCache)
80 : ggCache_(&ggCache) {}
81
83 const SubControlVolume& scv(LocalIndexType scvIdx) const
84 { return ggCache_->scvs(eIdx_)[scvIdx]; }
85
87 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
88 { return ggCache_->scvfs(eIdx_)[scvfIdx]; }
89
98 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
99 scvs(const BoxDfmFVElementGeometry& fvGeometry)
100 {
101 using Iter = typename std::vector<SubControlVolume>::const_iterator;
102 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
103 return Dune::IteratorRange<Iter>(s.begin(), s.end());
104 }
105
107 template<class LocalDof>
108 friend inline std::ranges::range auto
109 scvs(const BoxDfmFVElementGeometry& fvGeometry, const LocalDof& localDof)
110 {
111 return std::views::iota(0u, fvGeometry.numScv())
112 | std::views::filter([&](std::size_t i) { return fvGeometry.scv(i).localDofIndex() == localDof.index(); })
113 | std::views::transform([&](std::size_t i) -> const SubControlVolume& { return fvGeometry.scv(i); });
114 }
115
117 friend inline auto localDofs(const BoxDfmFVElementGeometry& fvGeometry)
118 {
119 return Dune::transformedRangeView(
120 Dune::range(fvGeometry.numLocalDofs()),
121 [&](const auto i) { return CVFE::LocalDof
122 {
123 static_cast<LocalIndexType>(fvGeometry.scv(i).localDofIndex()),
124 static_cast<GridIndexType>(fvGeometry.scv(i).dofIndex()),
125 static_cast<GridIndexType>(fvGeometry.scv(i).elementIndex())
126 }; }
127 );
128 }
129
131 friend inline auto cvLocalDofs(const BoxDfmFVElementGeometry& fvGeometry)
132 {
133 // Default it that all dofs are cv dofs
134 return localDofs(fvGeometry);
135 }
136
145 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
147 {
148 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
149 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
150 return Dune::IteratorRange<Iter>(s.begin(), s.end());
151 }
152
154 const FeLocalBasis& feLocalBasis() const
155 { return gridGeometry().feCache().get(element_->type()).localBasis(); }
156
158 std::size_t numScv() const
159 { return ggCache_->scvs(eIdx_).size(); }
160
162 std::size_t numLocalDofs() const
163 { return element().geometry().corners(); }
164
166 std::size_t numScvf() const
167 { return ggCache_->scvfs(eIdx_).size(); }
168
175 {
176 this->bindElement(element);
177 return std::move(*this);
178 }
179
182 void bind(const Element& element) &
183 { this->bindElement(element); }
184
191 {
192 this->bindElement(element);
193 return std::move(*this);
194 }
195
202 void bindElement(const Element& element) &
203 {
204 element_ = element;
205 eIdx_ = gridGeometry().elementMapper().index(element);
206 }
207
210 { return ggCache_->gridGeometry(); }
211
213 bool isBound() const
214 { return static_cast<bool>(element_); }
215
217 const Element& element() const
218 { return *element_; }
219
221 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
222 {
223 if (scv.isOnFracture())
224 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
225 "because the number of known corners is insufficient. "
226 "You can do this manually by extracting the corners from this scv "
227 "and extruding them by the corresponding aperture. ");
228
229 const GeometryHelper geometryHelper(element().geometry());
230 const auto corners = geometryHelper.getScvCorners(scv.index());
231 using ScvGeometry = typename SubControlVolume::Traits::Geometry;
232 return { Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners };
233 }
234
236 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
237 {
238 if (scvf.isOnFracture())
239 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
240 "because the number of known corners is insufficient. "
241 "You can do this manually by extracting the corners from this scv "
242 "and extruding them by the corresponding aperture. ");
243 const GeometryHelper geometryHelper(element().geometry());
244 const auto corners = geometryHelper.getScvfCorners(scvf.indexInElement());
245 using ScvfGeometry = typename SubControlVolumeFace::Traits::Geometry;
246 return { Dune::GeometryTypes::cube(ScvfGeometry::mydimension), corners };
247 }
248
249private:
250 const GGCache* ggCache_;
251
252 std::optional<Element> element_;
253 GridIndexType eIdx_;
254};
255
257template<class GG>
259{
260 using GridView = typename GG::GridView;
261 static constexpr int dim = GridView::dimension;
262 static constexpr int dimWorld = GridView::dimensionworld;
263
264 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
265 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
266
267 using CoordScalar = typename GridView::ctype;
268 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
269 using GGCache = typename GG::Cache;
270 using GeometryHelper = typename GGCache::GeometryHelper;
271public:
273 using Element = typename GridView::template Codim<0>::Entity;
275 using SubControlVolume = typename GG::SubControlVolume;
277 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
279 using GridGeometry = GG;
282 static constexpr std::size_t maxNumElementScvs = (1<<dim)*3;
283
285 [[deprecated("This Constructor is deprecated and will be removed after release 3.11. Always use localView(gridGeometry).")]]
287 : BoxDfmFVElementGeometry(gridGeometry.cache_) {}
288
290 BoxDfmFVElementGeometry(const GGCache& ggCache)
291 : ggCache_(&ggCache) {}
292
294 const SubControlVolume& scv(LocalIndexType scvIdx) const
295 { return scvs_[scvIdx]; }
296
298 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
299 { return scvfs_[scvfIdx]; }
300
309 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
310 scvs(const BoxDfmFVElementGeometry& fvGeometry)
311 {
312 using Iter = typename std::vector<SubControlVolume>::const_iterator;
313 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
314 }
315
317 template<class LocalDof>
318 friend inline std::ranges::range auto
319 scvs(const BoxDfmFVElementGeometry& fvGeometry, const LocalDof& localDof)
320 {
321 return std::views::iota(0u, fvGeometry.numScv())
322 | std::views::filter([&](std::size_t i) { return fvGeometry.scv(i).localDofIndex() == localDof.index(); })
323 | std::views::transform([&](std::size_t i) -> const SubControlVolume& { return fvGeometry.scv(i); });
324 }
325
327 friend inline auto localDofs(const BoxDfmFVElementGeometry& fvGeometry)
328 {
329 return Dune::transformedRangeView(
330 Dune::range(fvGeometry.numLocalDofs()),
331 [&](const auto i) { return CVFE::LocalDof
332 {
333 static_cast<LocalIndexType>(fvGeometry.scv(i).localDofIndex()),
334 static_cast<GridIndexType>(fvGeometry.scv(i).dofIndex()),
335 static_cast<GridIndexType>(fvGeometry.scv(i).elementIndex())
336 }; }
337 );
338 }
339
341 friend inline auto cvLocalDofs(const BoxDfmFVElementGeometry& fvGeometry)
342 {
343 // Default it that all dofs are cv dofs
344 return localDofs(fvGeometry);
345 }
346
355 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
357 {
358 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
359 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
360 }
361
363 const FeLocalBasis& feLocalBasis() const
364 { return gridGeometry().feCache().get(element_->type()).localBasis(); }
365
367 std::size_t numScv() const
368 { return scvs_.size(); }
369
371 std::size_t numLocalDofs() const
372 { return element().geometry().corners(); }
373
375 std::size_t numScvf() const
376 { return scvfs_.size(); }
377
384 {
385 this->bindElement(element);
386 return std::move(*this);
387 }
388
395 void bind(const Element& element) &
396 { this->bindElement(element); }
397
404 {
405 this->bindElement(element);
406 return std::move(*this);
407 }
408
413 void bindElement(const Element& element) &
414 {
415 element_ = element;
416 eIdx_ = gridGeometry().elementMapper().index(element);
417 makeElementGeometries_();
418 }
419
422 { return ggCache_->gridGeometry(); }
423
425 bool isBound() const
426 { return static_cast<bool>(element_); }
427
429 const Element& element() const
430 { return *element_; }
431
433 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
434 {
435 if (scv.isOnFracture())
436 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
437 "because the number of known corners is insufficient. "
438 "You can do this manually by extracting the corners from this scv "
439 "and extruding them by the corresponding aperture. ");
440
441 const GeometryHelper geometryHelper(element().geometry());
442 const auto corners = geometryHelper.getScvCorners(scv.index());
443 using ScvGeometry = typename SubControlVolume::Traits::Geometry;
444 return { Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners };
445 }
446
448 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
449 {
450 if (scvf.isOnFracture())
451 DUNE_THROW(Dune::InvalidStateException, "The geometry object cannot be defined for fracture scvs "
452 "because the number of known corners is insufficient. "
453 "You can do this manually by extracting the corners from this scv "
454 "and extruding them by the corresponding aperture. ");
455 const GeometryHelper geometryHelper(element().geometry());
456 const auto corners = geometryHelper.getScvfCorners(scvf.indexInElement());
457 using ScvfGeometry = typename SubControlVolumeFace::Traits::Geometry;
458 return { Dune::GeometryTypes::cube(ScvfGeometry::mydimension), corners };
459 }
460
461private:
462
463 void makeElementGeometries_()
464 {
465 // get the element geometry
466 const auto& element = *element_;
467 const auto elementGeometry = element.geometry();
468 const auto refElement = referenceElement(element);
469
470 // get the sub control volume geometries of this element
471 GeometryHelper geometryHelper(elementGeometry);
472
473 // construct the sub control volumes
474 scvs_.resize(elementGeometry.corners());
475 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
476 {
477 // get associated dof index
478 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
479
480 // add scv to the local container
481 scvs_[scvLocalIdx] = SubControlVolume(geometryHelper,
482 scvLocalIdx,
483 eIdx_,
484 dofIdxGlobal);
485 }
486
487 // construct the sub control volume faces
488 const auto numInnerScvf = (dim==1) ? 1 : element.subEntities(dim-1);
489 scvfs_.resize(numInnerScvf);
490
491 LocalIndexType scvfLocalIdx = 0;
492 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
493 {
494 // find the local scv indices this scvf is connected to
495 std::array<LocalIndexType, 2> localScvIndices{{
496 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
497 static_cast<LocalIndexType>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
498 }};
499
500 scvfs_[scvfLocalIdx] = SubControlVolumeFace(geometryHelper,
501 element,
502 scvfLocalIdx,
503 std::move(localScvIndices));
504 }
505
506 // construct the ...
507 // ... sub-control volume faces on the domain boundary
508 // ... sub-control volumes on fracture facets
509 // ... sub-control volume faces on fracture facets
510 // NOTE We do not construct fracture scvfs on boundaries here!
511 // That means specifying Neumann fluxes on only the fractures is not possible
512 // However, it is difficult to interpret this here as a fracture ending on the boundary
513 // could also be connected to a facet which is both boundary and fracture at the same time!
514 // In that case, the fracture boundary scvf wouldn't make sense. In order to do it properly
515 // we would have to find only those fractures that are at the boundary and aren't connected
516 // to a fracture which is a boundary.
517 LocalIndexType scvLocalIdx = element.subEntities(dim);
518 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
519 {
520 // first, obtain all vertex indices on this intersection
521 const auto& isGeometry = intersection.geometry();
522 const auto numCorners = isGeometry.corners();
523 const auto idxInInside = intersection.indexInInside();
524
525 std::vector<GridIndexType> isVertexIndices(numCorners);
526 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
527 isVertexIndices[vIdxLocal] = gridGeometry().vertexMapper().subIndex(element,
528 refElement.subEntity(idxInInside, 1, vIdxLocal, dim),
529 dim);
530
531 if (intersection.boundary())
532 {
533 for (unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < numCorners; ++isScvfLocalIdx)
534 {
535 // find the scv this scvf is connected to
536 const LocalIndexType insideScvIdx = static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, isScvfLocalIdx, dim));
537 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
538
539 scvfs_.emplace_back(geometryHelper,
540 intersection,
541 isScvfLocalIdx,
542 scvfLocalIdx,
543 std::move(localScvIndices));
544
545 // increment local counter
546 scvfLocalIdx++;
547 }
548 }
549
550 // maybe add fracture scvs & scvfs
551 if (this->gridGeometry().isOnFracture(element, intersection))
552 {
553 // add fracture scv for each vertex of intersection
554 const auto curNumScvs = scvs_.size();
555 scvs_.reserve(curNumScvs+numCorners);
556 for (unsigned int vIdxLocal = 0; vIdxLocal < numCorners; ++vIdxLocal)
557 scvs_.emplace_back(geometryHelper,
558 intersection,
559 isGeometry,
560 vIdxLocal,
561 static_cast<LocalIndexType>(refElement.subEntity(idxInInside, 1, vIdxLocal, dim)),
562 scvLocalIdx++,
563 idxInInside,
564 eIdx_,
565 isVertexIndices[vIdxLocal]);
566
567 // add fracture scvf for each edge of the intersection in 3d
568 if (dim == 3)
569 {
570 const auto& faceRefElement = referenceElement(isGeometry);
571 for (unsigned int edgeIdx = 0; edgeIdx < faceRefElement.size(1); ++edgeIdx)
572 {
573 // inside/outside scv indices in face local node numbering
574 std::array<LocalIndexType, 2> localScvIndices{{
575 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 0, dim-1)),
576 static_cast<LocalIndexType>(faceRefElement.subEntity(edgeIdx, 1, 1, dim-1))
577 }};
578
579 // add offset to get the right scv indices
580 std::for_each( localScvIndices.begin(),
581 localScvIndices.end(),
582 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
583
584 // add scvf
585 scvfs_.emplace_back(geometryHelper,
586 intersection,
587 edgeIdx,
588 scvfLocalIdx++,
589 std::move(localScvIndices),
590 intersection.boundary());
591 }
592 }
593
594 // dim == 2, intersection is an edge, make 1 scvf
595 else
596 {
597 // inside/outside scv indices in face local node numbering
598 std::array<LocalIndexType, 2> localScvIndices{{0, 1}};
599
600 // add offset such that the fracture scvs above are addressed
601 std::for_each( localScvIndices.begin(),
602 localScvIndices.end(),
603 [curNumScvs] (auto& elemLocalIdx) { elemLocalIdx += curNumScvs; } );
604
605 // add scvf
606 scvfs_.emplace_back(geometryHelper,
607 intersection,
608 /*idxOnIntersection*/0,
609 scvfLocalIdx++,
610 std::move(localScvIndices),
611 intersection.boundary());
612 }
613 }
614 }
615 }
616
618 std::optional<Element> element_;
619 GridIndexType eIdx_;
620
622 const GGCache* ggCache_;
623
625 std::vector<SubControlVolume> scvs_;
626 std::vector<SubControlVolumeFace> scvfs_;
627};
628
629} // end namespace Dumux
630
631#endif
friend auto localDofs(const BoxDfmFVElementGeometry &fvGeometry)
range over local dofs
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:327
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Create the geometry of a given sub control volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:433
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:363
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:421
BoxDfmFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:290
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:371
BoxDfmFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:403
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:273
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:375
friend auto cvLocalDofs(const BoxDfmFVElementGeometry &fvGeometry)
range over control-volume local dofs
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:341
void bind(const Element &element) &
Binding of an element, has to be called before using the fvgeometries Prepares all the volume variabl...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:395
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:294
BoxDfmFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:383
const Element & element() const
The bound element.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:429
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Create the geometry of a given sub control volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:448
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:286
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:367
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:425
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:310
friend std::ranges::range auto scvs(const BoxDfmFVElementGeometry &fvGeometry, const LocalDof &localDof)
range over sub control volumes related to a local dof.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:319
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:279
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:277
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:356
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:298
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:275
void bindElement(const Element &element) &
Binding of an element, has to be called before using the fvgeometries Prepares all the volume variabl...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:413
typename GG::SubControlVolumeFace SubControlVolumeFace
Export type of subcontrol volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:65
void bindElement(const Element &element) &
Binding of an element, has to be called before using the fvgeometries.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:202
friend auto localDofs(const BoxDfmFVElementGeometry &fvGeometry)
range over local dofs
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:117
std::size_t numScv() const
The total number of sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:158
GG GridGeometry
Export type of finite volume grid geometry.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:67
const Element & element() const
The bound element.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:217
void bind(const Element &element) &
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:182
BoxDfmFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:75
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:162
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:166
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:154
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:61
typename GG::SubControlVolume SubControlVolume
Export type of subcontrol volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:63
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Create the geometry of a given sub control volume face.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:236
friend auto cvLocalDofs(const BoxDfmFVElementGeometry &fvGeometry)
range over control-volume local dofs
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:131
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:87
BoxDfmFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:79
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:83
BoxDfmFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:174
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:99
friend std::ranges::range auto scvs(const BoxDfmFVElementGeometry &fvGeometry, const LocalDof &localDof)
range over sub control volumes related to a local dof.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:109
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:209
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:213
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Create the geometry of a given sub control volume.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:221
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxDfmFVElementGeometry &fvGeometry)
Iterator range for sub control volumes faces.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:146
BoxDfmFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:190
Base class for the finite volume geometry vector for box discrete fracture model.
Definition: porousmediumflow/boxdfm/fvelementgeometry.hh:44
Defines the index types used for grid and local indices.
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
Definition: adapt.hh:17
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50
Helper class constructing the dual grid finite volume geometries for the box discrete fracture model.
Class providing iterators over sub control volumes and sub control volume faces of an element.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
unsigned int LocalIndex
Definition: indextraits.hh:28