26#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
32#include <dune/common/exceptions.hh>
33#include <dune/geometry/type.hh>
34#include <dune/localfunctions/lagrange/pqkfactory.hh>
51template<
class GG,
bool enableGr
idGeometryCache>
58 using GridView =
typename GG::GridView;
59 static constexpr int dim = GridView::dimension;
60 static constexpr int dimWorld = GridView::dimensionworld;
63 using CoordScalar =
typename GridView::ctype;
64 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
65 using GGCache =
typename GG::Cache;
69 using Element =
typename GridView::template Codim<0>::Entity;
78 static constexpr std::size_t maxNumElementScvs = (1<<dim) + 1;
88 return ggCache_->scvs(eIdx_)[scvIdx];
94 return ggCache_->scvfs(eIdx_)[scvfIdx];
102 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
105 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
106 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
107 return Dune::IteratorRange<Iter>(s.begin(), s.end());
115 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
118 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
119 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
120 return Dune::IteratorRange<Iter>(s.begin(), s.end());
126 return gridGeometry().feCache().get(element_->type()).localBasis();
132 return ggCache_->scvs(eIdx_).size();
138 return ggCache_->scvfs(eIdx_).size();
148 this->bindElement(element);
149 return std::move(*
this);
156 { this->bindElement(element); }
165 this->bindElement(element);
166 return std::move(*
this);
176 eIdx_ = gridGeometry().elementMapper().index(element);
181 {
return static_cast<bool>(element_); }
185 {
return *element_; }
189 {
return ggCache_->gridGeometry(); }
193 {
return ggCache_->hasBoundaryScvf(eIdx_); }
202 if (scv.isOverlapping())
203 DUNE_THROW(Dune::NotImplemented,
"Geometry of overlapping scv");
206 const auto geo = element().geometry();
218 const auto geo = element().geometry();
223 const auto [localFacetIndex, isScvfLocalIdx]
224 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
226 Dune::GeometryTypes::cube(dim-1),
241 const GGCache* ggCache_;
244 std::optional<Element> element_;
Class providing iterators over sub control volumes and sub control volume faces of an element.
Defines the index types used for grid and local indices.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Base class for the finite volume geometry vector for pq1bubble models This builds up the sub control ...
Definition: discretization/pq1bubble/fvelementgeometry.hh:52
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:103
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:116
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/pq1bubble/fvelementgeometry.hh:69
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/pq1bubble/fvelementgeometry.hh:192
PQ1BubbleFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/pq1bubble/fvelementgeometry.hh:81
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:86
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/pq1bubble/fvelementgeometry.hh:75
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/pq1bubble/fvelementgeometry.hh:124
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/pq1bubble/fvelementgeometry.hh:73
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/pq1bubble/fvelementgeometry.hh:71
PQ1BubbleFVElementGeometry 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/pq1bubble/fvelementgeometry.hh:163
void bindElement(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:172
void bind(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:155
PQ1BubbleFVElementGeometry 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/pq1bubble/fvelementgeometry.hh:146
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/pq1bubble/fvelementgeometry.hh:215
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/pq1bubble/fvelementgeometry.hh:130
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/pq1bubble/fvelementgeometry.hh:200
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:92
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/pq1bubble/fvelementgeometry.hh:180
std::size_t elementIndex() const
The bound element index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:196
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/pq1bubble/fvelementgeometry.hh:136
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/pq1bubble/fvelementgeometry.hh:188
const Element & element() const
The bound element.
Definition: discretization/pq1bubble/fvelementgeometry.hh:184
A class to create sub control volume and sub control volume face geometries per element.
Definition: discretization/pq1bubble/geometryhelper.hh:178
std::size_t numInteriorScvf() const
number of interior sub control volume faces
Definition: discretization/pq1bubble/geometryhelper.hh:343
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition: discretization/pq1bubble/geometryhelper.hh:300
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:290
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition: discretization/pq1bubble/geometryhelper.hh:254
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition: discretization/pq1bubble/geometryhelper.hh:198
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:234
Helper class constructing the dual grid finite volume geometries for the cvfe discretizazion method.