14#ifndef DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
20#include <unordered_map>
24#include <dune/geometry/type.hh>
25#include <dune/localfunctions/lagrange/pqkfactory.hh>
26#include <dune/common/reservedvector.hh>
45template<
class GG,
bool enableGr
idGeometryCache>
52 using GridView =
typename GG::GridView;
53 static constexpr int dim = GridView::dimension;
54 static constexpr int dimWorld = GridView::dimensionworld;
57 using CoordScalar =
typename GridView::ctype;
58 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
59 using GGCache =
typename GG::Cache;
60 using GeometryHelper =
typename GGCache::GeometryHelper;
63 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
64 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
69 using Element =
typename GridView::template Codim<0>::Entity;
90 : ggCache_(&ggCache) {}
95 return ggCache_->scvs(eIdx_)[scvIdx];
101 return ggCache_->scvfs(eIdx_)[scvfIdx];
109 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
112 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
113 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
114 return Dune::IteratorRange<Iter>(s.begin(), s.end());
122 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
125 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
126 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
127 return Dune::IteratorRange<Iter>(s.begin(), s.end());
132 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
135 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
136 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
137 const auto& range = fvGeometry.ggCache_->boundaryFaceScvfRanges(fvGeometry.eIdx_)[
boundaryFace.index()];
138 return Dune::IteratorRange<Iter>(s.begin() + range[0], s.begin() + range[0] + range[1]);
143 friend inline std::ranges::view
auto
146 const auto& v = fvGeometry.ggCache_->boundaryFaces(fvGeometry.eIdx_);
147 return std::ranges::views::all(v);
153 return Dune::transformedRangeView(
154 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(
boundaryFace.intersectionIndex(), 1, dim)),
157 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(boundaryFace.intersectionIndex(), 1, i, dim);
158 return CVFE::LocalDof
160 static_cast<LocalIndexType>(localDofIdx),
161 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
162 static_cast<GridIndexType>(fvGeometry.elementIndex())
171 return gridGeometry().feCache().get(element_->type()).localBasis();
183 return ggCache_->scvs(eIdx_).size();
189 return ggCache_->scvfs(eIdx_).size();
200 return std::move(*
this);
217 return std::move(*
this);
226 elementGeometry_.emplace(
element.geometry());
233 {
return static_cast<bool>(element_); }
237 {
return *element_; }
241 {
return *elementGeometry_; }
250 const auto localScvfIdx =
scvf.index() - GeometryHelper::numInteriorScvf(
element().type());
251 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
256 {
return ggCache_->gridGeometry(); }
260 {
return ggCache_->hasBoundaryScvf(eIdx_); }
268 {
return ggCache_->boundaryFaces(eIdx_)[bfIdx]; }
274 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(
scv.indexInElement()) };
281 const GeometryHelper geometryHelper(*elementGeometry_);
284 const auto localBoundaryIndex =
scvf.index() - geometryHelper.numInteriorScvf();
285 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
286 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
289 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(
scvf.index()) };
298 typename BoundaryFace::Traits::CornerStorage corners;
299 for (
int i = 0; i < faceGeoInRef.corners(); ++i)
300 corners.push_back(elemGeo.global(faceGeoInRef.corner(i)));
301 return { faceGeoInRef.type(), corners };
307 const auto type = fvGeometry.element().type();
308 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(
scv.localDofIndex());
314 template<
class LocalDof>
317 const auto type = fvGeometry.element().type();
318 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
319 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
329 [&] (
const typename Element::Geometry::GlobalCoordinate& pos)
330 {
return fvGeometry.elementGeometry().local(pos); },
339 {
scvf.unitOuterNormal(),
scvf.index(), fvGeometry.elementGeometry().local(
scvf.ipGlobal()),
scvf.ipGlobal() };
343 const GGCache* ggCache_;
346 std::optional<Element> element_;
347 std::optional<typename Element::Geometry> elementGeometry_;
354 using GridView =
typename GG::GridView;
355 static constexpr int dim = GridView::dimension;
356 static constexpr int dimWorld = GridView::dimensionworld;
359 using CoordScalar =
typename GridView::ctype;
360 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
361 using GGCache =
typename GG::Cache;
362 using GeometryHelper =
typename GGCache::GeometryHelper;
365 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
366 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
371 using Element =
typename GridView::template Codim<0>::Entity;
392 : ggCache_(&ggCache) {}
397 return scvs_[scvIdx];
403 return scvfs_[scvfIdx];
411 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
414 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
415 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
423 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
426 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
427 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
432 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
435 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
436 const auto& range = fvGeometry.boundaryFaceScvfRanges_[
boundaryFace.index()];
437 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin() + range[0], fvGeometry.scvfs_.begin() + range[0] + range[1]);
442 friend inline std::ranges::view
auto
444 {
return std::ranges::views::all(fvGeometry.boundaryFaces_); }
449 return Dune::transformedRangeView(
450 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(
boundaryFace.intersectionIndex(), 1, dim)),
453 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(boundaryFace.intersectionIndex(), 1, i, dim);
454 return CVFE::LocalDof
456 static_cast<LocalIndexType>(localDofIdx),
457 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
458 static_cast<GridIndexType>(fvGeometry.elementIndex())
468 return gridGeometry().feCache().get(element_->type()).localBasis();
486 return scvfs_.size();
497 return std::move(*
this);
514 return std::move(*
this);
524 elementGeometry_.emplace(
element.geometry());
525 makeElementGeometries_();
530 {
return static_cast<bool>(element_); }
534 {
return *element_; }
538 {
return *elementGeometry_; }
547 const auto localScvfIdx =
scvf.index() - GeometryHelper::numInteriorScvf(
element().type());
548 return scvfBoundaryGeometryKeys_[localScvfIdx][0];
553 {
return ggCache_->gridGeometry(); }
557 {
return hasBoundaryScvf_; }
565 {
return boundaryFaces_[bfIdx]; }
571 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(
scv.indexInElement()) };
578 const GeometryHelper geometryHelper(*elementGeometry_);
581 const auto localBoundaryIndex =
scvf.index() - geometryHelper.numInteriorScvf();
582 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
583 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
586 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(
scvf.index()) };
595 typename BoundaryFace::Traits::CornerStorage corners;
596 for (
int i = 0; i < faceGeoInRef.corners(); ++i)
597 corners.push_back(elemGeo.global(faceGeoInRef.corner(i)));
598 return { faceGeoInRef.type(), corners };
604 const auto type = fvGeometry.element().type();
605 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(
scv.localDofIndex());
611 template<
class LocalDof>
614 const auto type = fvGeometry.element().type();
615 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
616 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
626 [&] (
const typename Element::Geometry::GlobalCoordinate& pos) {
return fvGeometry.elementGeometry().local(pos); },
635 {
scvf.unitOuterNormal(),
scvf.index(), fvGeometry.elementGeometry().local(
scvf.ipGlobal()),
scvf.ipGlobal() };
639 void makeElementGeometries_()
641 hasBoundaryScvf_ =
false;
642 boundaryFaces_.clear();
643 boundaryFaceScvfRanges_.clear();
646 const auto& element = *element_;
647 const auto& elementGeometry = *elementGeometry_;
648 const auto refElement = referenceElement(elementGeometry);
651 GeometryHelper geometryHelper(elementGeometry);
654 scvs_.resize(elementGeometry.corners());
655 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
658 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
661 scvs_[scvLocalIdx] = SubControlVolume(
662 geometryHelper.getScvCorners(scvLocalIdx),
670 const auto numInnerScvf = geometryHelper.numInteriorScvf();
671 scvfs_.resize(numInnerScvf);
672 scvfBoundaryGeometryKeys_.clear();
674 LocalIndexType scvfLocalIdx = 0;
675 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
678 std::array<LocalIndexType, 2> localScvIndices{{
679 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
680 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
683 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
684 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
686 geometryHelper.normal(corners, localScvIndices),
689 std::move(localScvIndices)
694 LocalIndexType numBoundaryFaces = 0;
695 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
697 if (intersection.boundary() && !intersection.neighbor())
699 const auto isGeometry = intersection.geometry();
700 hasBoundaryScvf_ =
true;
703 boundaryFaces_.push_back(BoundaryFace{
706 intersection.centerUnitOuterNormal(),
708 static_cast<LocalIndexType
>(intersection.indexInInside()),
709 typename BoundaryFace::Traits::BoundaryFlag{intersection}
713 boundaryFaceScvfRanges_.push_back(std::array<LocalIndexType, 2>{{
715 static_cast<LocalIndexType
>(isGeometry.corners())
718 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
721 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
722 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
725 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
726 intersection.centerUnitOuterNormal(),
730 std::move(localScvIndices)
733 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
734 static_cast<LocalIndexType
>(intersection.indexInInside()),
735 static_cast<LocalIndexType
>(isScvfLocalIdx)
747 std::optional<Element> element_;
748 std::optional<typename Element::Geometry> elementGeometry_;
751 const GGCache* ggCache_;
754 std::vector<SubControlVolume> scvs_;
755 std::vector<SubControlVolumeFace> scvfs_;
756 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
757 Dune::ReservedVector<BoundaryFace, 2*dim> boundaryFaces_;
758 Dune::ReservedVector<std::array<LocalIndexType, 2>, 2*dim> boundaryFaceScvfRanges_;
760 bool hasBoundaryScvf_ =
false;
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition discretization/box/fvelementgeometry.hh:529
std::size_t numScvf() const
The total number of sub control volume faces.
Definition discretization/box/fvelementgeometry.hh:484
GG GridGeometry
export type of finite volume grid geometry
Definition discretization/box/fvelementgeometry.hh:377
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition discretization/box/fvelementgeometry.hh:575
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition discretization/box/fvelementgeometry.hh:556
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:412
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/box/fvelementgeometry.hh:478
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition discretization/box/fvelementgeometry.hh:541
void bindElement(const Element &element) &
Definition discretization/box/fvelementgeometry.hh:520
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition discretization/box/fvelementgeometry.hh:545
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition discretization/box/fvelementgeometry.hh:401
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition discretization/box/fvelementgeometry.hh:391
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/box/fvelementgeometry.hh:383
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition discretization/box/fvelementgeometry.hh:568
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition discretization/box/fvelementgeometry.hh:552
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
Definition discretization/box/fvelementgeometry.hh:433
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition discretization/box/fvelementgeometry.hh:371
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:424
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition discretization/box/fvelementgeometry.hh:602
BoundaryFace::Traits::Geometry geometry(const BoundaryFace &boundaryFace) const
Geometry of a boundary face.
Definition discretization/box/fvelementgeometry.hh:590
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition discretization/box/fvelementgeometry.hh:466
const Element & element() const
The bound element.
Definition discretization/box/fvelementgeometry.hh:533
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition discretization/box/fvelementgeometry.hh:395
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition discretization/box/fvelementgeometry.hh:622
BoxFVElementGeometry 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 discretization/box/fvelementgeometry.hh:494
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition discretization/box/fvelementgeometry.hh:632
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
an iterator over all local dofs related to a boundary face
Definition discretization/box/fvelementgeometry.hh:447
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition discretization/box/fvelementgeometry.hh:612
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition discretization/box/fvelementgeometry.hh:381
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition discretization/box/fvelementgeometry.hh:472
typename GG::BoundaryFace BoundaryFace
export the boundary face type
Definition discretization/box/fvelementgeometry.hh:379
friend std::ranges::view auto boundaryFaces(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:443
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition discretization/box/fvelementgeometry.hh:375
const BoundaryFace & boundaryFace(LocalIndexType bfIdx) const
Get a boundary face with a local boundary face index.
Definition discretization/box/fvelementgeometry.hh:564
BoxFVElementGeometry 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 discretization/box/fvelementgeometry.hh:511
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition discretization/box/fvelementgeometry.hh:537
bool hasBoundaryFaces() const
Returns whether the element has boundary faces.
Definition discretization/box/fvelementgeometry.hh:560
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition discretization/box/fvelementgeometry.hh:373
void bind(const Element &element) &
Definition discretization/box/fvelementgeometry.hh:503
static constexpr std::size_t maxNumElementScvs
the maximum number of scvs per element (2^dim for cubes)
Definition discretization/box/fvelementgeometry.hh:385
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition discretization/box/fvelementgeometry.hh:232
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition discretization/box/fvelementgeometry.hh:169
BoundaryFace::Traits::Geometry geometry(const BoundaryFace &boundaryFace) const
Geometry of a boundary face.
Definition discretization/box/fvelementgeometry.hh:293
BoxFVElementGeometry 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 discretization/box/fvelementgeometry.hh:214
const Element & element() const
The bound element.
Definition discretization/box/fvelementgeometry.hh:236
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition discretization/box/fvelementgeometry.hh:271
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:110
bool hasBoundaryFaces() const
Returns whether the element has boundary faces.
Definition discretization/box/fvelementgeometry.hh:263
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/box/fvelementgeometry.hh:181
typename GG::BoundaryFace BoundaryFace
export the boundary face type
Definition discretization/box/fvelementgeometry.hh:77
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition discretization/box/fvelementgeometry.hh:93
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition discretization/box/fvelementgeometry.hh:175
std::size_t numScvf() const
The total number of sub control volume faces.
Definition discretization/box/fvelementgeometry.hh:187
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition discretization/box/fvelementgeometry.hh:71
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition discretization/box/fvelementgeometry.hh:255
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
Definition discretization/box/fvelementgeometry.hh:133
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition discretization/box/fvelementgeometry.hh:259
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:123
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition discretization/box/fvelementgeometry.hh:305
typename GG::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition discretization/box/fvelementgeometry.hh:79
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition discretization/box/fvelementgeometry.hh:240
void bind(const Element &element) &
Definition discretization/box/fvelementgeometry.hh:206
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/box/fvelementgeometry.hh:81
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition discretization/box/fvelementgeometry.hh:248
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition discretization/box/fvelementgeometry.hh:278
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition discretization/box/fvelementgeometry.hh:325
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition discretization/box/fvelementgeometry.hh:244
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition discretization/box/fvelementgeometry.hh:336
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
an iterator over all local dofs related to a boundary face
Definition discretization/box/fvelementgeometry.hh:151
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition discretization/box/fvelementgeometry.hh:315
friend std::ranges::view auto boundaryFaces(const BoxFVElementGeometry &fvGeometry)
Definition discretization/box/fvelementgeometry.hh:144
const BoundaryFace & boundaryFace(LocalIndexType bfIdx) const
Get a boundary face with a local boundary face index.
Definition discretization/box/fvelementgeometry.hh:267
static constexpr std::size_t maxNumElementScvs
the maximum number of scvs per element (2^dim for cubes)
Definition discretization/box/fvelementgeometry.hh:83
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition discretization/box/fvelementgeometry.hh:73
void bindElement(const Element &element) &
Definition discretization/box/fvelementgeometry.hh:223
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition discretization/box/fvelementgeometry.hh:99
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition discretization/box/fvelementgeometry.hh:89
BoxFVElementGeometry 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 discretization/box/fvelementgeometry.hh:197
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition discretization/box/fvelementgeometry.hh:69
GG GridGeometry
export type of finite volume grid geometry
Definition discretization/box/fvelementgeometry.hh:75
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition discretization/box/fvelementgeometry.hh:46
An interpolation point related to a face of an element.
Definition cvfe/interpolationpointdata.hh:109
An interpolation point related to an element that includes global and local positions.
Definition cvfe/interpolationpointdata.hh:31
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition cvfe/interpolationpointdata.hh:82
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition cvfe/interpolationpointdata.hh:60
Classes representing interpolation point data for control-volume finite element schemes.
Class representing dofs on elements for control-volume finite element schemes.
Quadrature rules over sub-control volumes and sub-control volume faces.
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