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>
42template<
class GG,
bool enableGr
idGeometryCache>
49 using GridView =
typename GG::GridView;
50 static constexpr int dim = GridView::dimension;
51 static constexpr int dimWorld = GridView::dimensionworld;
54 using CoordScalar =
typename GridView::ctype;
55 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
56 using GGCache =
typename GG::Cache;
57 using GeometryHelper =
typename GGCache::GeometryHelper;
60 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
64 using Element =
typename GridView::template Codim<0>::Entity;
72 static constexpr std::size_t maxNumElementScvs = (1<<dim);
79 : ggCache_(&ggCache) {}
84 return ggCache_->scvs(eIdx_)[scvIdx];
90 return ggCache_->scvfs(eIdx_)[scvfIdx];
98 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
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());
111 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
114 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
115 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
116 return Dune::IteratorRange<Iter>(s.begin(), s.end());
120 template<
class Intersection>
123 return Dune::transformedRangeView(
124 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(intersection.indexInInside(), 1, dim)),
127 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(intersection.indexInInside(), 1, i, dim);
128 return CVFE::LocalDof
130 static_cast<LocalIndexType>(localDofIdx),
131 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
132 static_cast<GridIndexType>(fvGeometry.elementIndex())
141 return gridGeometry().feCache().get(element_->type()).localBasis();
153 return ggCache_->scvs(eIdx_).size();
159 return ggCache_->scvfs(eIdx_).size();
169 this->bindElement(element);
170 return std::move(*
this);
177 { this->bindElement(element); }
186 this->bindElement(element);
187 return std::move(*
this);
196 elementGeometry_.emplace(element.geometry());
198 eIdx_ = gridGeometry().elementMapper().index(element);
203 {
return static_cast<bool>(element_); }
207 {
return *element_; }
211 {
return *elementGeometry_; }
220 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
221 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
226 {
return ggCache_->gridGeometry(); }
230 {
return ggCache_->hasBoundaryScvf(eIdx_); }
236 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
243 const GeometryHelper geometryHelper(*elementGeometry_);
246 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
247 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
248 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
251 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
257 const auto type = fvGeometry.element().type();
258 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
260 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
264 template<
class LocalDof>
267 const auto type = fvGeometry.element().type();
268 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
269 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
271 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
279 [&] (
const typename Element::Geometry::GlobalCoordinate& pos)
280 {
return fvGeometry.elementGeometry().local(pos); },
286 const GGCache* ggCache_;
289 std::optional<Element> element_;
290 std::optional<typename Element::Geometry> elementGeometry_;
297 using GridView =
typename GG::GridView;
298 static constexpr int dim = GridView::dimension;
299 static constexpr int dimWorld = GridView::dimensionworld;
302 using CoordScalar =
typename GridView::ctype;
303 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
304 using GGCache =
typename GG::Cache;
305 using GeometryHelper =
typename GGCache::GeometryHelper;
307 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
311 using Element =
typename GridView::template Codim<0>::Entity;
319 static constexpr std::size_t maxNumElementScvs = (1<<dim);
326 : ggCache_(&ggCache) {}
331 return scvs_[scvIdx];
337 return scvfs_[scvfIdx];
345 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
348 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
349 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
357 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
360 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
361 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
365 template<
class Intersection>
368 return Dune::transformedRangeView(
369 Dune::range(Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).size(intersection.indexInInside(), 1, dim)),
372 auto localDofIdx = Dune::referenceElement<CoordScalar, dim>(fvGeometry.element().type()).subEntity(intersection.indexInInside(), 1, i, dim);
373 return CVFE::LocalDof
375 static_cast<LocalIndexType>(localDofIdx),
376 fvGeometry.gridGeometry().dofMapper().subIndex(fvGeometry.element(), localDofIdx, dim),
377 static_cast<GridIndexType>(fvGeometry.elementIndex())
387 return gridGeometry().feCache().get(element_->type()).localBasis();
405 return scvfs_.size();
415 this->bindElement(element);
416 return std::move(*
this);
423 { this->bindElement(element); }
432 this->bindElement(element);
433 return std::move(*
this);
442 eIdx_ = gridGeometry().elementMapper().index(element);
443 elementGeometry_.emplace(element.geometry());
444 makeElementGeometries_();
449 {
return static_cast<bool>(element_); }
453 {
return *element_; }
457 {
return *elementGeometry_; }
466 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
467 return scvfBoundaryGeometryKeys_[localScvfIdx][0];
472 {
return ggCache_->gridGeometry(); }
476 {
return hasBoundaryScvf_; }
482 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
489 const GeometryHelper geometryHelper(*elementGeometry_);
492 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
493 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
494 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
497 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
503 const auto type = fvGeometry.element().type();
504 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
506 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
510 template<
class LocalDof>
513 const auto type = fvGeometry.element().type();
514 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
515 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
517 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
525 [&] (
const typename Element::Geometry::GlobalCoordinate& pos) {
return fvGeometry.elementGeometry().local(pos); },
531 void makeElementGeometries_()
533 hasBoundaryScvf_ =
false;
536 const auto& element = *element_;
537 const auto& elementGeometry = *elementGeometry_;
538 const auto refElement = referenceElement(elementGeometry);
541 GeometryHelper geometryHelper(elementGeometry);
544 scvs_.resize(elementGeometry.corners());
545 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
548 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
551 scvs_[scvLocalIdx] = SubControlVolume(
552 geometryHelper.getScvCorners(scvLocalIdx),
560 const auto numInnerScvf = geometryHelper.numInteriorScvf();
561 scvfs_.resize(numInnerScvf);
562 scvfBoundaryGeometryKeys_.clear();
564 LocalIndexType scvfLocalIdx = 0;
565 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
568 std::array<LocalIndexType, 2> localScvIndices{{
569 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
570 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
573 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
574 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
576 geometryHelper.normal(corners, localScvIndices),
579 std::move(localScvIndices)
584 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
586 if (intersection.boundary() && !intersection.neighbor())
588 const auto isGeometry = intersection.geometry();
589 hasBoundaryScvf_ =
true;
591 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
594 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
595 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
598 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
599 intersection.centerUnitOuterNormal(),
603 std::move(localScvIndices)
606 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
607 static_cast<LocalIndexType
>(intersection.indexInInside()),
608 static_cast<LocalIndexType
>(isScvfLocalIdx)
620 std::optional<Element> element_;
621 std::optional<typename Element::Geometry> elementGeometry_;
624 const GGCache* ggCache_;
627 std::vector<SubControlVolume> scvs_;
628 std::vector<SubControlVolumeFace> scvfs_;
629 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
631 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:448
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/box/fvelementgeometry.hh:366
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:403
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:317
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:486
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:475
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:346
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:397
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:460
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:439
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/box/fvelementgeometry.hh:464
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:335
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:325
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:479
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:471
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:311
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:358
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:501
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:385
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:452
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:329
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:521
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:413
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:391
friend IpData ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:511
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:315
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:430
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:456
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:313
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:422
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:202
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:139
friend auto localDofs(const BoxFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/box/fvelementgeometry.hh:121
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:184
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:206
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:233
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:99
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:151
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:82
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:145
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:157
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:66
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:225
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:229
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:112
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:255
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:210
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:176
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/box/fvelementgeometry.hh:218
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:240
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:275
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:214
friend IpData ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:265
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:68
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:193
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:88
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:78
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:167
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:64
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:70
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:43
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition: cvfe/interpolationpointdata.hh:72
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition: cvfe/interpolationpointdata.hh:50
Classes representing interpolation point data for control-volume finite element schemes.
Class representing dofs on elements for control-volume finite element schemes.
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