12#ifndef DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
18#include <dune/common/reservedvector.hh>
19#include <dune/common/iteratorrange.hh>
31template<
class GG,
bool cachingEnabled>
42 using GridView =
typename GG::GridView;
45 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
46 using GGCache =
typename GG::Cache;
47 using GeometryHelper =
typename GGCache::GeometryHelper;
53 using Element =
typename GridView::template Codim<0>::Entity;
57 static constexpr std::size_t maxNumElementScvs = 2*GridView::dimension;
65 {
return ggCache_->scvs(eIdx_)[scvIdx]; }
69 {
return ggCache_->scvf(eIdx_)[scvfIdx]; }
79 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
80 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
81 return Dune::IteratorRange<Iter>(s.begin(), s.end());
92 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
93 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
94 return Dune::IteratorRange<Iter>(s.begin(), s.end());
100 return gridGeometry().feCache().get(element().type()).localBasis();
106 return ggCache_->scvs(eIdx_).size();
112 return ggCache_->scvfs(eIdx_).size();
117 {
return ggCache_->hasBoundaryScvf(eIdx_); }
126 this->bindElement(element);
127 return std::move(*
this);
132 this->bindElement(element);
142 this->bindElement(element);
143 return std::move(*
this);
150 eIdx_ = gridGeometry().elementMapper().index(element);
155 {
return static_cast<bool>(element_); }
159 {
return *element_; }
163 {
return ggCache_->gridGeometry(); }
173 const auto geo = element().geometry();
175 SubControlVolume::Traits::geometryType(geo.type()),
176 GeometryHelper(geo).getScvCorners(scv.indexInElement())
184 const auto geo = element().geometry();
188 const auto localFacetIndex = scvf.insideScvIdx();
190 referenceElement(geo).type(localFacetIndex, 1),
191 GeometryHelper(geo).getBoundaryScvfCorners(localFacetIndex)
197 SubControlVolumeFace::Traits::interiorGeometryType(geo.type()),
198 GeometryHelper(geo).getScvfCorners(scvf.index())
204 std::optional<Element> element_;
206 const GGCache* ggCache_;
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:40
friend auto scvs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:77
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:147
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:162
const Element & element() const
The bound element.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:158
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:51
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:116
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:52
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:154
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:110
friend auto scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:90
std::size_t elementIndex() const
The bound element index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:166
GG GridGeometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:54
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:68
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:64
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:53
FaceCenteredDiamondFVElementGeometry 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/facecentered/diamond/fvelementgeometry.hh:124
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:104
FaceCenteredDiamondFVElementGeometry(const GGCache &ggCache)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:59
void bind(const Element &element) &
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:130
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:98
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:181
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:170
FaceCenteredDiamondFVElementGeometry 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/facecentered/diamond/fvelementgeometry.hh:140
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:32
Class providing iterators over sub control volumes and sub control volume faces of an element.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:29