12#ifndef DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
20#include <dune/common/reservedvector.hh>
21#include <dune/common/iteratorrange.hh>
37template<
class GG,
bool cachingEnabled>
48 using GridView =
typename GG::GridView;
51 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
52 using GGCache =
typename GG::Cache;
53 using GeometryHelper =
typename GGCache::GeometryHelper;
56 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
57 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
65 using Element =
typename GridView::template Codim<0>::Entity;
81 {
return ggCache_->scvs(eIdx_)[scvIdx]; }
85 {
return ggCache_->scvf(eIdx_)[scvfIdx]; }
95 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
96 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
97 return Dune::IteratorRange<Iter>(s.begin(), s.end());
108 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
109 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
110 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 const auto& range = fvGeometry.ggCache_->boundaryFaceScvfRanges(fvGeometry.eIdx_)[
boundaryFace.index()];
121 return Dune::IteratorRange<Iter>(s.begin() + range[0], s.begin() + range[0] + range[1]);
139 return ggCache_->scvs(eIdx_).size();
145 return ggCache_->scvfs(eIdx_).size();
150 {
return ggCache_->hasBoundaryScvf(eIdx_); }
158 {
return ggCache_->boundaryFaces(eIdx_)[bfIdx]; }
161 friend inline std::ranges::view
auto
164 const auto& v = fvGeometry.ggCache_->boundaryFaces(fvGeometry.eIdx_);
165 return std::ranges::views::all(v);
176 return std::move(*
this);
192 return std::move(*
this);
199 elementGeometry_.emplace(
element.geometry());
205 {
return static_cast<bool>(element_); }
209 {
return *element_; }
213 {
return *elementGeometry_; }
217 {
return ggCache_->gridGeometry(); }
227 return scvf.insideScvIdx();
235 SubControlVolume::Traits::geometryType((*elementGeometry_).type()),
236 GeometryHelper(*elementGeometry_).getScvCorners(
scv.indexInElement())
247 const auto localFacetIndex =
scvf.insideScvIdx();
249 referenceElement(*elementGeometry_).type(localFacetIndex, 1),
250 GeometryHelper(*elementGeometry_).getBoundaryScvfCorners(localFacetIndex)
256 SubControlVolumeFace::Traits::interiorGeometryType((*elementGeometry_).type()),
257 GeometryHelper(*elementGeometry_).getScvfCorners(
scvf.index())
268 typename BoundaryFace::Traits::CornerStorage corners;
269 for (
int i = 0; i < faceGeoInRef.corners(); ++i)
270 corners.push_back(elemGeo.global(faceGeoInRef.corner(i)));
271 return { faceGeoInRef.type(), corners };
278 const auto localDofIdx =
boundaryFace.intersectionIndex();
280 static_cast<LocalIndexType
>(localDofIdx),
281 static_cast<GridIndexType
>(fvGeometry.scv(localDofIdx).dofIndex()),
282 static_cast<GridIndexType
>(fvGeometry.elementIndex())
289 const auto type = fvGeometry.element().type();
290 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(
scv.localDofIndex());
296 template<
class LocalDof>
299 const auto type = fvGeometry.element().type();
300 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
301 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
311 [&] (
const typename Element::Geometry::GlobalCoordinate& pos) {
return fvGeometry.elementGeometry().local(pos); },
320 {
scvf.unitOuterNormal(),
scvf.index(), fvGeometry.elementGeometry().local(
scvf.ipGlobal()),
scvf.ipGlobal() };
324 std::optional<Element> element_;
325 std::optional<typename Element::Geometry> elementGeometry_;
327 const GGCache* ggCache_;
An interpolation point related to a face of an element.
Definition cvfe/interpolationpointdata.hh:109
An interpolation point related to an element that includes global and local positions.
Definition cvfe/interpolationpointdata.hh:31
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition cvfe/interpolationpointdata.hh:82
A local degree of freedom from an element perspective.
Definition localdof.hh:27
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition cvfe/interpolationpointdata.hh:60
friend auto scvs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition discretization/facecentered/diamond/fvelementgeometry.hh:93
void bindElement(const Element &element) &
Bind only element-local.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:196
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:216
static constexpr std::size_t maxNumElementScvs
the maximum number of scvs per element
Definition discretization/facecentered/diamond/fvelementgeometry.hh:73
const Element & element() const
The bound element.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:208
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition discretization/facecentered/diamond/fvelementgeometry.hh:68
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition discretization/facecentered/diamond/fvelementgeometry.hh:62
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:297
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:224
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:212
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:317
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:149
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition discretization/facecentered/diamond/fvelementgeometry.hh:63
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/facecentered/diamond/fvelementgeometry.hh:70
bool hasBoundaryFaces() const
Returns whether the element has boundary faces.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:153
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:287
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:204
const BoundaryFace & boundaryFace(LocalIndexType bfIdx) const
Get a boundary face with a local boundary face index.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:157
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition discretization/facecentered/diamond/fvelementgeometry.hh:143
friend auto scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition discretization/facecentered/diamond/fvelementgeometry.hh:106
std::size_t elementIndex() const
The bound element index.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:220
GG GridGeometry
Definition discretization/facecentered/diamond/fvelementgeometry.hh:66
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:307
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
Definition discretization/facecentered/diamond/fvelementgeometry.hh:116
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:84
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:80
typename GridView::template Codim< 0 >::Entity Element
Definition discretization/facecentered/diamond/fvelementgeometry.hh:65
friend auto localDofs(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
an iterator over all local dofs related to a boundary face
Definition discretization/facecentered/diamond/fvelementgeometry.hh:275
typename GG::BoundaryFace BoundaryFace
Definition discretization/facecentered/diamond/fvelementgeometry.hh:64
BoundaryFace::Traits::Geometry geometry(const BoundaryFace &boundaryFace) const
Geometry of a boundary face.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:263
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:173
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition discretization/facecentered/diamond/fvelementgeometry.hh:137
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:131
FaceCenteredDiamondFVElementGeometry(const GGCache &ggCache)
Definition discretization/facecentered/diamond/fvelementgeometry.hh:75
friend std::ranges::view auto boundaryFaces(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
iterator range for boundary faces
Definition discretization/facecentered/diamond/fvelementgeometry.hh:162
void bind(const Element &element) &
Definition discretization/facecentered/diamond/fvelementgeometry.hh:179
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:125
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:241
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:231
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:189
Element-wise grid geometry (local view).
Definition discretization/facecentered/diamond/fvelementgeometry.hh:38
Classes representing interpolation point data for control-volume finite element schemes.
Helper class to construct SCVs and SCVFs for the diamond scheme.
Quadrature rules over sub-control volumes and sub-control volume faces.
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