14#ifndef DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_BOX_FV_ELEMENT_GEOMETRY_HH
19#include <unordered_map>
23#include <dune/geometry/type.hh>
24#include <dune/localfunctions/lagrange/pqkfactory.hh>
43template<
class GG,
bool enableGr
idGeometryCache>
50 using GridView =
typename GG::GridView;
51 static constexpr int dim = GridView::dimension;
52 static constexpr int dimWorld = GridView::dimensionworld;
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;
61 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
62 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
67 using Element =
typename GridView::template Codim<0>::Entity;
79 static constexpr std::size_t maxNumElementScvs = (1<<dim);
86 : ggCache_(&ggCache) {}
91 return ggCache_->scvs(eIdx_)[scvIdx];
97 return ggCache_->scvfs(eIdx_)[scvfIdx];
105 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
108 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
109 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
110 return Dune::IteratorRange<Iter>(s.begin(), s.end());
118 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
121 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
122 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
123 return Dune::IteratorRange<Iter>(s.begin(), s.end());
127 template<
class Intersection>
130 return Dune::transformedRangeView(
131 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(intersection.indexInInside(), 1, dim)),
134 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(intersection.indexInInside(), 1, i, dim);
135 return CVFE::LocalDof
137 static_cast<LocalIndexType>(localDofIdx),
138 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
139 static_cast<GridIndexType>(fvGeometry.elementIndex())
148 return gridGeometry().feCache().get(element_->type()).localBasis();
160 return ggCache_->scvs(eIdx_).size();
166 return ggCache_->scvfs(eIdx_).size();
176 this->bindElement(element);
177 return std::move(*
this);
184 { this->bindElement(element); }
193 this->bindElement(element);
194 return std::move(*
this);
203 elementGeometry_.emplace(element.geometry());
205 eIdx_ = gridGeometry().elementMapper().index(element);
210 {
return static_cast<bool>(element_); }
214 {
return *element_; }
218 {
return *elementGeometry_; }
227 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
228 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
233 {
return ggCache_->gridGeometry(); }
237 {
return ggCache_->hasBoundaryScvf(eIdx_); }
243 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
250 const GeometryHelper geometryHelper(*elementGeometry_);
253 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
254 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
255 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
258 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
264 const auto type = fvGeometry.element().type();
265 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
271 template<
class LocalDof>
274 const auto type = fvGeometry.element().type();
275 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
276 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
286 [&] (
const typename Element::Geometry::GlobalCoordinate& pos)
287 {
return fvGeometry.elementGeometry().local(pos); },
296 { scvf.
unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
300 const GGCache* ggCache_;
303 std::optional<Element> element_;
304 std::optional<typename Element::Geometry> elementGeometry_;
311 using GridView =
typename GG::GridView;
312 static constexpr int dim = GridView::dimension;
313 static constexpr int dimWorld = GridView::dimensionworld;
316 using CoordScalar =
typename GridView::ctype;
317 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
318 using GGCache =
typename GG::Cache;
319 using GeometryHelper =
typename GGCache::GeometryHelper;
322 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
323 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
328 using Element =
typename GridView::template Codim<0>::Entity;
340 static constexpr std::size_t maxNumElementScvs = (1<<dim);
347 : ggCache_(&ggCache) {}
352 return scvs_[scvIdx];
358 return scvfs_[scvfIdx];
366 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
369 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
370 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
378 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
381 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
382 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
386 template<
class Intersection>
389 return Dune::transformedRangeView(
390 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(intersection.indexInInside(), 1, dim)),
393 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(intersection.indexInInside(), 1, i, dim);
394 return CVFE::LocalDof
396 static_cast<LocalIndexType>(localDofIdx),
397 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
398 static_cast<GridIndexType>(fvGeometry.elementIndex())
408 return gridGeometry().feCache().get(element_->type()).localBasis();
426 return scvfs_.size();
436 this->bindElement(element);
437 return std::move(*
this);
444 { this->bindElement(element); }
453 this->bindElement(element);
454 return std::move(*
this);
463 eIdx_ = gridGeometry().elementMapper().index(element);
464 elementGeometry_.emplace(element.geometry());
465 makeElementGeometries_();
470 {
return static_cast<bool>(element_); }
474 {
return *element_; }
478 {
return *elementGeometry_; }
487 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
488 return scvfBoundaryGeometryKeys_[localScvfIdx][0];
493 {
return ggCache_->gridGeometry(); }
497 {
return hasBoundaryScvf_; }
503 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
510 const GeometryHelper geometryHelper(*elementGeometry_);
513 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
514 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
515 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
518 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
524 const auto type = fvGeometry.element().type();
525 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
531 template<
class LocalDof>
534 const auto type = fvGeometry.element().type();
535 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
536 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
546 [&] (
const typename Element::Geometry::GlobalCoordinate& pos) {
return fvGeometry.elementGeometry().local(pos); },
555 { scvf.
unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
559 void makeElementGeometries_()
561 hasBoundaryScvf_ =
false;
564 const auto& element = *element_;
565 const auto& elementGeometry = *elementGeometry_;
566 const auto refElement = referenceElement(elementGeometry);
569 GeometryHelper geometryHelper(elementGeometry);
572 scvs_.resize(elementGeometry.corners());
573 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
576 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
579 scvs_[scvLocalIdx] = SubControlVolume(
580 geometryHelper.getScvCorners(scvLocalIdx),
588 const auto numInnerScvf = geometryHelper.numInteriorScvf();
589 scvfs_.resize(numInnerScvf);
590 scvfBoundaryGeometryKeys_.clear();
592 LocalIndexType scvfLocalIdx = 0;
593 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
596 std::array<LocalIndexType, 2> localScvIndices{{
597 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
598 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
601 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
602 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
604 geometryHelper.normal(corners, localScvIndices),
607 std::move(localScvIndices)
612 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
614 if (intersection.boundary() && !intersection.neighbor())
616 const auto isGeometry = intersection.geometry();
617 hasBoundaryScvf_ =
true;
619 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
622 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
623 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
626 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
627 intersection.centerUnitOuterNormal(),
631 std::move(localScvIndices)
634 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
635 static_cast<LocalIndexType
>(intersection.indexInInside()),
636 static_cast<LocalIndexType
>(isScvfLocalIdx)
648 std::optional<Element> element_;
649 std::optional<typename Element::Geometry> elementGeometry_;
652 const GGCache* ggCache_;
655 std::vector<SubControlVolume> scvs_;
656 std::vector<SubControlVolumeFace> scvfs_;
657 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
659 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:469
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/box/fvelementgeometry.hh:387
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:424
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:334
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:507
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:496
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:367
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:418
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:481
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:460
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/box/fvelementgeometry.hh:485
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:356
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:346
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition: discretization/box/fvelementgeometry.hh:338
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:500
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:492
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:328
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:379
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:522
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:406
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:473
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:350
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:542
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:434
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition: discretization/box/fvelementgeometry.hh:552
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:532
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition: discretization/box/fvelementgeometry.hh:336
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:412
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:332
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:451
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:477
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:330
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:443
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:209
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:146
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/box/fvelementgeometry.hh:128
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:191
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:213
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:240
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:106
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:158
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:89
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:152
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:164
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:69
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:232
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:236
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:119
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:262
typename GG::ScvQuadratureRule ScvQuadratureRule
export the scv interpolation point data type
Definition: discretization/box/fvelementgeometry.hh:75
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:217
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:183
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition: discretization/box/fvelementgeometry.hh:77
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/box/fvelementgeometry.hh:225
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:247
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:282
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:221
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition: discretization/box/fvelementgeometry.hh:293
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:272
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:71
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:200
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:95
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:85
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:174
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:67
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:73
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:44
An interpolation point related to a face of an element.
Definition: cvfe/interpolationpointdata.hh:103
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector at the quadrature point.
Definition: cvfe/interpolationpointdata.hh:117
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:25
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition: cvfe/interpolationpointdata.hh:76
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition: cvfe/interpolationpointdata.hh:54
LocalIndex localDofIndex() const
The local index of the corresponding dof.
Definition: cvfe/interpolationpointdata.hh:63
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