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>
41template<
class GG,
bool enableGr
idGeometryCache>
48 using GridView =
typename GG::GridView;
49 static constexpr int dim = GridView::dimension;
50 static constexpr int dimWorld = GridView::dimensionworld;
53 using CoordScalar =
typename GridView::ctype;
54 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
55 using GGCache =
typename GG::Cache;
56 using GeometryHelper =
typename GGCache::GeometryHelper;
59 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
63 using Element =
typename GridView::template Codim<0>::Entity;
71 static constexpr std::size_t maxNumElementScvs = (1<<dim);
78 : ggCache_(&ggCache) {}
83 return ggCache_->scvs(eIdx_)[scvIdx];
89 return ggCache_->scvfs(eIdx_)[scvfIdx];
97 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
100 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
101 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
102 return Dune::IteratorRange<Iter>(s.begin(), s.end());
110 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
113 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
114 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
115 return Dune::IteratorRange<Iter>(s.begin(), s.end());
121 return gridGeometry().feCache().get(element_->type()).localBasis();
133 return ggCache_->scvs(eIdx_).size();
139 return ggCache_->scvfs(eIdx_).size();
149 this->bindElement(element);
150 return std::move(*
this);
157 { this->bindElement(element); }
166 this->bindElement(element);
167 return std::move(*
this);
176 elementGeometry_.emplace(element.geometry());
178 eIdx_ = gridGeometry().elementMapper().index(element);
183 {
return static_cast<bool>(element_); }
187 {
return *element_; }
191 {
return *elementGeometry_; }
199 {
return ggCache_->gridGeometry(); }
203 {
return ggCache_->hasBoundaryScvf(eIdx_); }
209 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
216 const GeometryHelper geometryHelper(*elementGeometry_);
219 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
220 const auto& key = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localBoundaryIndex];
221 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
224 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
230 const auto type = fvGeometry.element().type();
231 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
233 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
237 template<
class LocalDof>
240 const auto type = fvGeometry.element().type();
241 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
242 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
244 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
252 [&] (
const typename Element::Geometry::GlobalCoordinate& pos)
253 {
return fvGeometry.elementGeometry().local(pos); },
259 const GGCache* ggCache_;
262 std::optional<Element> element_;
263 std::optional<typename Element::Geometry> elementGeometry_;
270 using GridView =
typename GG::GridView;
271 static constexpr int dim = GridView::dimension;
272 static constexpr int dimWorld = GridView::dimensionworld;
275 using CoordScalar =
typename GridView::ctype;
276 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
277 using GGCache =
typename GG::Cache;
278 using GeometryHelper =
typename GGCache::GeometryHelper;
280 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
284 using Element =
typename GridView::template Codim<0>::Entity;
292 static constexpr std::size_t maxNumElementScvs = (1<<dim);
299 : ggCache_(&ggCache) {}
304 return scvs_[scvIdx];
310 return scvfs_[scvfIdx];
318 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
321 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
322 return Dune::IteratorRange<Iter>(fvGeometry.scvs_.begin(), fvGeometry.scvs_.end());
330 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
333 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
334 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
340 return gridGeometry().feCache().get(element_->type()).localBasis();
358 return scvfs_.size();
368 this->bindElement(element);
369 return std::move(*
this);
376 { this->bindElement(element); }
385 this->bindElement(element);
386 return std::move(*
this);
395 eIdx_ = gridGeometry().elementMapper().index(element);
396 elementGeometry_.emplace(element.geometry());
397 makeElementGeometries_();
402 {
return static_cast<bool>(element_); }
406 {
return *element_; }
410 {
return *elementGeometry_; }
418 {
return ggCache_->gridGeometry(); }
422 {
return hasBoundaryScvf_; }
428 return { Dune::GeometryTypes::cube(dim), GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement()) };
435 const GeometryHelper geometryHelper(*elementGeometry_);
438 const auto localBoundaryIndex = scvf.index() - geometryHelper.numInteriorScvf();
439 const auto& key = scvfBoundaryGeometryKeys_[localBoundaryIndex];
440 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getBoundaryScvfCorners(key[0], key[1]) };
443 return { Dune::GeometryTypes::cube(dim-1), geometryHelper.getScvfCorners(scvf.index()) };
449 const auto type = fvGeometry.element().type();
450 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
452 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
456 template<
class LocalDof>
459 const auto type = fvGeometry.element().type();
460 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
461 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
463 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
471 [&] (
const typename Element::Geometry::GlobalCoordinate& pos) {
return fvGeometry.elementGeometry().local(pos); },
477 void makeElementGeometries_()
479 hasBoundaryScvf_ =
false;
482 const auto& element = *element_;
483 const auto& elementGeometry = *elementGeometry_;
484 const auto refElement = referenceElement(elementGeometry);
487 GeometryHelper geometryHelper(elementGeometry);
490 scvs_.resize(elementGeometry.corners());
491 for (LocalIndexType scvLocalIdx = 0; scvLocalIdx < elementGeometry.corners(); ++scvLocalIdx)
494 const auto dofIdxGlobal = gridGeometry().vertexMapper().subIndex(element, scvLocalIdx, dim);
497 scvs_[scvLocalIdx] = SubControlVolume(
498 geometryHelper.getScvCorners(scvLocalIdx),
506 const auto numInnerScvf = geometryHelper.numInteriorScvf();
507 scvfs_.resize(numInnerScvf);
508 scvfBoundaryGeometryKeys_.clear();
510 LocalIndexType scvfLocalIdx = 0;
511 for (; scvfLocalIdx < numInnerScvf; ++scvfLocalIdx)
514 std::array<LocalIndexType, 2> localScvIndices{{
515 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 0, dim)),
516 static_cast<LocalIndexType
>(refElement.subEntity(scvfLocalIdx, dim-1, 1, dim))
519 const auto& corners = geometryHelper.getScvfCorners(scvfLocalIdx);
520 scvfs_[scvfLocalIdx] = SubControlVolumeFace(
522 geometryHelper.normal(corners, localScvIndices),
525 std::move(localScvIndices)
530 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
532 if (intersection.boundary() && !intersection.neighbor())
534 const auto isGeometry = intersection.geometry();
535 hasBoundaryScvf_ =
true;
537 for (
unsigned int isScvfLocalIdx = 0; isScvfLocalIdx < isGeometry.corners(); ++isScvfLocalIdx)
540 const LocalIndexType insideScvIdx =
static_cast<LocalIndexType
>(refElement.subEntity(intersection.indexInInside(), 1, isScvfLocalIdx, dim));
541 std::array<LocalIndexType, 2> localScvIndices{{insideScvIdx, insideScvIdx}};
544 geometryHelper.getBoundaryScvfCorners(intersection.indexInInside(), isScvfLocalIdx),
545 intersection.centerUnitOuterNormal(),
549 std::move(localScvIndices)
552 scvfBoundaryGeometryKeys_.emplace_back(std::array<LocalIndexType, 2>{{
553 static_cast<LocalIndexType
>(intersection.indexInInside()),
554 static_cast<LocalIndexType
>(isScvfLocalIdx)
566 std::optional<Element> element_;
567 std::optional<typename Element::Geometry> elementGeometry_;
570 const GGCache* ggCache_;
573 std::vector<SubControlVolume> scvs_;
574 std::vector<SubControlVolumeFace> scvfs_;
575 std::vector<std::array<LocalIndexType, 2>> scvfBoundaryGeometryKeys_;
577 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:401
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:356
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:290
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:432
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:421
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:319
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:350
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:413
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:392
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:308
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:298
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:425
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:417
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:284
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:331
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:447
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:338
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:405
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:302
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:467
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:366
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:344
friend IpData ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:457
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:288
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:383
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:409
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:286
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:375
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/box/fvelementgeometry.hh:182
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/box/fvelementgeometry.hh:119
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:164
const Element & element() const
The bound element.
Definition: discretization/box/fvelementgeometry.hh:186
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/box/fvelementgeometry.hh:206
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:98
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/box/fvelementgeometry.hh:131
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/box/fvelementgeometry.hh:81
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/box/fvelementgeometry.hh:125
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/box/fvelementgeometry.hh:137
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/box/fvelementgeometry.hh:65
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/box/fvelementgeometry.hh:198
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/box/fvelementgeometry.hh:202
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const BoxFVElementGeometry &fvGeometry)
Definition: discretization/box/fvelementgeometry.hh:111
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/box/fvelementgeometry.hh:228
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/box/fvelementgeometry.hh:190
void bind(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:156
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/box/fvelementgeometry.hh:213
friend auto ipData(const BoxFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/box/fvelementgeometry.hh:248
GridIndexType elementIndex() const
The bound element's index in the grid view.
Definition: discretization/box/fvelementgeometry.hh:194
friend IpData ipData(const BoxFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/box/fvelementgeometry.hh:238
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/box/fvelementgeometry.hh:67
void bindElement(const Element &element) &
Definition: discretization/box/fvelementgeometry.hh:173
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/box/fvelementgeometry.hh:87
BoxFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/box/fvelementgeometry.hh:77
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:147
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/box/fvelementgeometry.hh:63
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/box/fvelementgeometry.hh:69
Base class for the finite volume geometry vector for box models This builds up the sub control volume...
Definition: discretization/box/fvelementgeometry.hh:42
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 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