24#ifndef DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
30#include <dune/common/reservedvector.hh>
31#include <dune/common/iteratorrange.hh>
43template<
class GG,
bool cachingEnabled>
54 using GridView =
typename GG::GridView;
57 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
58 using GGCache =
typename GG::Cache;
65 using Element =
typename GridView::template Codim<0>::Entity;
69 static constexpr std::size_t maxNumElementScvs = 2*GridView::dimension;
77 {
return ggCache_->scvs(eIdx_)[scvIdx]; }
81 {
return ggCache_->scvf(eIdx_)[scvfIdx]; }
91 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
92 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
93 return Dune::IteratorRange<Iter>(s.begin(), s.end());
104 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
105 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
106 return Dune::IteratorRange<Iter>(s.begin(), s.end());
112 return gridGeometry().feCache().get(element().type()).localBasis();
118 return ggCache_->scvs(eIdx_).size();
124 return ggCache_->scvfs(eIdx_).size();
129 {
return ggCache_->hasBoundaryScvf(eIdx_); }
138 this->bindElement(element);
139 return std::move(*
this);
144 this->bindElement(element);
154 this->bindElement(element);
155 return std::move(*
this);
162 eIdx_ = gridGeometry().elementMapper().index(element);
167 {
return static_cast<bool>(element_); }
171 {
return *element_; }
175 {
return ggCache_->gridGeometry(); }
185 const auto geo = element().geometry();
187 SubControlVolume::Traits::geometryType(geo.type()),
196 const auto geo = element().geometry();
200 const auto localFacetIndex = scvf.insideScvIdx();
202 referenceElement(geo).type(localFacetIndex, 1),
209 SubControlVolumeFace::Traits::interiorGeometryType(geo.type()),
216 std::optional<Element> element_;
218 const GGCache* ggCache_;
Defines the index types used for grid and local indices.
Class providing iterators over sub control volumes and sub control volume faces of an element.
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
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:44
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:52
friend auto scvs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:89
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:159
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:174
const Element & element() const
The bound element.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:170
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:63
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:128
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:64
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:166
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:122
friend auto scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:102
std::size_t elementIndex() const
The bound element index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:178
GG GridGeometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:66
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:80
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:76
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:65
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:136
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:116
FaceCenteredDiamondFVElementGeometry(const GGCache &ggCache)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:71
void bind(const Element &element) &
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:142
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:110
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:193
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:182
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:152
Helper class to construct SCVs and SCVFs for the diamond scheme.
Definition: discretization/facecentered/diamond/geometryhelper.hh:267
ScvfCornerStorage getScvfCorners(unsigned int localEdgeIndex) const
Create a corner storage with the scvf corners for a given edge (codim-2) index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:315
ScvCornerStorage getScvCorners(unsigned int localFacetIndex) const
Create a corner storage with the scv corners for a given face (codim-1) index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:287
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex) const
Create the sub control volume face geometries on the boundary for a given face index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:343