14#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
20#include <dune/common/exceptions.hh>
21#include <dune/geometry/type.hh>
22#include <dune/localfunctions/lagrange/pqkfactory.hh>
39template<
class GG,
bool enableGr
idGeometryCache>
46 using GridView =
typename GG::GridView;
47 static constexpr int dim = GridView::dimension;
48 static constexpr int dimWorld = GridView::dimensionworld;
51 using CoordScalar =
typename GridView::ctype;
52 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
53 using GGCache =
typename GG::Cache;
54 using GeometryHelper =
typename GGCache::GeometryHelper;
57 using Element =
typename GridView::template Codim<0>::Entity;
66 static constexpr std::size_t maxNumElementScvs = (1<<dim) + 1;
76 return ggCache_->scvs(eIdx_)[scvIdx];
82 return ggCache_->scvfs(eIdx_)[scvfIdx];
90 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
93 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
94 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
95 return Dune::IteratorRange<Iter>(s.begin(), s.end());
103 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
106 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
107 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
108 return Dune::IteratorRange<Iter>(s.begin(), s.end());
114 return gridGeometry().feCache().get(element_->type()).localBasis();
120 return ggCache_->scvs(eIdx_).size();
126 return ggCache_->scvfs(eIdx_).size();
136 this->bindElement(element);
137 return std::move(*
this);
144 { this->bindElement(element); }
153 this->bindElement(element);
154 return std::move(*
this);
164 eIdx_ = gridGeometry().elementMapper().index(element);
169 {
return static_cast<bool>(element_); }
173 {
return *element_; }
177 {
return ggCache_->gridGeometry(); }
181 {
return ggCache_->hasBoundaryScvf(eIdx_); }
190 if (scv.isOverlapping())
191 DUNE_THROW(Dune::NotImplemented,
"Geometry of overlapping scv");
194 const auto geo = element().geometry();
195 const GeometryHelper helper(geo);
197 helper.getScvGeometryType(scv.indexInElement()),
198 helper.getScvCorners(scv.indexInElement())
206 const auto geo = element().geometry();
209 GeometryHelper helper(geo);
210 const auto localScvfIdx = scvf.index() - helper.numInteriorScvf();
211 const auto [localFacetIndex, isScvfLocalIdx]
212 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
214 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
215 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
220 GeometryHelper helper(geo);
222 helper.getInteriorScvfGeometryType(scvf.index()),
223 helper.getScvfCorners(scvf.index())
229 const GGCache* ggCache_;
232 std::optional<Element> element_;
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:91
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:104
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/pq1bubble/fvelementgeometry.hh:57
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/pq1bubble/fvelementgeometry.hh:180
PQ1BubbleFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/pq1bubble/fvelementgeometry.hh:69
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:74
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/pq1bubble/fvelementgeometry.hh:63
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/pq1bubble/fvelementgeometry.hh:112
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/pq1bubble/fvelementgeometry.hh:61
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/pq1bubble/fvelementgeometry.hh:59
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:151
void bindElement(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:160
void bind(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:143
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:134
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/pq1bubble/fvelementgeometry.hh:203
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/pq1bubble/fvelementgeometry.hh:118
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/pq1bubble/fvelementgeometry.hh:188
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:80
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/pq1bubble/fvelementgeometry.hh:168
std::size_t elementIndex() const
The bound element index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:184
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/pq1bubble/fvelementgeometry.hh:124
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/pq1bubble/fvelementgeometry.hh:176
const Element & element() const
The bound element.
Definition: discretization/pq1bubble/fvelementgeometry.hh:172
Base class for the finite volume geometry vector for pq1bubble models This builds up the sub control ...
Definition: discretization/pq1bubble/fvelementgeometry.hh:40
Helper class constructing the dual grid finite volume geometries for the cvfe discretizazion method.
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