14#ifndef DUMUX_DISCRETIZATION_PQ2_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_PQ2_FV_ELEMENT_GEOMETRY_HH
25#include <dune/common/rangeutilities.hh>
26#include <dune/geometry/type.hh>
27#include <dune/localfunctions/lagrange/pqkfactory.hh>
47template<
class GG,
bool enableGr
idGeometryCache>
54 using GridView =
typename GG::GridView;
55 static constexpr int dim = GridView::dimension;
56 static constexpr int dimWorld = GridView::dimensionworld;
59 using CoordScalar =
typename GridView::ctype;
60 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
61 using GGCache =
typename GG::Cache;
62 using GeometryHelper =
typename GGCache::GeometryHelper;
65 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
66 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
71 using Element =
typename GridView::template Codim<0>::Entity;
101 return ggCache_->scvs(eIdx_)[scvIdx];
107 return ggCache_->scvfs(eIdx_)[scvfIdx];
115 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
118 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
119 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
120 return Dune::IteratorRange<Iter>(s.begin(), s.end());
125 return std::views::iota(std::size_t(0), fvGeometry.numLocalDofs())
126 | std::views::filter([&](
size_t i) {
return fvGeometry.feLocalCoefficients().localKey(i).codim() == dim; })
127 | std::views::transform([&](
size_t i) {
129 static_cast<LocalIndexType
>(i),
130 [&]() -> GridIndexType {
131 const auto& lk = fvGeometry.feLocalCoefficients().localKey(i);
132 if constexpr (
requires { fvGeometry.gridGeometry().dofIndex(fvGeometry.element(), lk); })
133 return static_cast<GridIndexType
>(fvGeometry.gridGeometry().dofIndex(fvGeometry.element(), lk));
135 return static_cast<GridIndexType
>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), lk));
137 static_cast<GridIndexType
>(fvGeometry.elementIndex())
144 return std::views::iota(std::size_t(0), fvGeometry.numLocalDofs())
145 | std::views::filter([&](
size_t i) {
return !(fvGeometry.feLocalCoefficients().localKey(i).codim() == dim); })
146 | std::views::transform([&](
size_t i) {
148 static_cast<LocalIndexType
>(i),
149 [&]() -> GridIndexType {
150 const auto& lk = fvGeometry.feLocalCoefficients().localKey(i);
151 if constexpr (
requires { fvGeometry.gridGeometry().dofIndex(fvGeometry.element(), lk); })
152 return static_cast<GridIndexType
>(fvGeometry.gridGeometry().dofIndex(fvGeometry.element(), lk));
154 return static_cast<GridIndexType
>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), lk));
156 static_cast<GridIndexType
>(fvGeometry.elementIndex())
163 return Dune::transformedRangeView(
164 Dune::range(std::size_t(0), fvGeometry.numLocalDofs()), [&](
const auto i) {
165 return CVFE::LocalDof{
166 static_cast<LocalIndexType>(i),
167 [&]() -> GridIndexType {
168 const auto& lk = fvGeometry.feLocalCoefficients().localKey(i);
169 if constexpr (requires { fvGeometry.gridGeometry().dofIndex(fvGeometry.element(), lk); })
170 return static_cast<GridIndexType>(fvGeometry.gridGeometry().dofIndex(fvGeometry.element(), lk));
172 return static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), lk));
174 static_cast<GridIndexType>(fvGeometry.elementIndex())
183 return std::views::iota(std::size_t(0), fvGeometry.numLocalDofs())
184 | std::views::filter([&](
size_t i) {
185 return GeometryHelper::localDofOnIntersection(fvGeometry.element().type(),
187 fvGeometry.feLocalCoefficients().localKey(i)); })
188 | std::views::transform([&](
size_t i) {
190 static_cast<LocalIndexType
>(i),
191 [&]() -> GridIndexType {
192 const auto& lk = fvGeometry.feLocalCoefficients().localKey(i);
193 if constexpr (
requires { fvGeometry.gridGeometry().dofIndex(fvGeometry.element(), lk); })
194 return static_cast<GridIndexType
>(fvGeometry.gridGeometry().dofIndex(fvGeometry.element(), lk));
196 return static_cast<GridIndexType
>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), lk));
198 static_cast<GridIndexType
>(fvGeometry.elementIndex())
208 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
211 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
212 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
213 return Dune::IteratorRange<Iter>(s.begin(), s.end());
218 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
221 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
222 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
223 const auto& range = fvGeometry.ggCache_->boundaryFaceScvfRanges(fvGeometry.eIdx_)[
boundaryFace.index()];
224 return Dune::IteratorRange<Iter>(s.begin() + range[0], s.begin() + range[0] + range[1]);
229 friend inline std::ranges::view
auto
232 const auto& v = fvGeometry.ggCache_->boundaryFaces(fvGeometry.eIdx_);
233 return std::ranges::views::all(v);
239 return gridGeometry().feCache().get(element_->type()).localBasis();
245 return gridGeometry().feCache().get(element_->type()).localCoefficients();
257 return ggCache_->scvs(eIdx_).size();
263 return ggCache_->scvfs(eIdx_).size();
274 return std::move(*
this);
291 return std::move(*
this);
302 elementGeometry_.emplace(
element.geometry());
307 {
return static_cast<bool>(element_); }
311 {
return *element_; }
315 {
return *elementGeometry_; }
319 {
return ggCache_->gridGeometry(); }
323 {
return ggCache_->hasBoundaryScvf(eIdx_); }
331 {
return ggCache_->boundaryFaces(eIdx_)[bfIdx]; }
340 const auto localScvfIdx =
scvf.index() - GeometryHelper::numInteriorScvf(
element().type());
341 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
349 const GeometryHelper helper(geo);
352 helper.getScvGeometryType(scvIdx),
353 helper.getScvCorners(scvIdx)
363 GeometryHelper helper(*elementGeometry_);
364 const auto localScvfIdx =
scvf.index() - GeometryHelper::numInteriorScvf(
element().type());
365 const auto [localFacetIndex, isScvfLocalIdx]
366 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
368 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
369 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
374 GeometryHelper helper(*elementGeometry_);
376 helper.getInteriorScvfGeometryType(
scvf.index()),
377 helper.getScvfCorners(
scvf.index())
388 typename BoundaryFace::Traits::CornerStorage corners;
389 for (
int i = 0; i < faceGeoInRef.corners(); ++i)
390 corners.push_back(elemGeo.global(faceGeoInRef.corner(i)));
391 return { faceGeoInRef.type(), corners };
397 const auto type = fvGeometry.element().type();
398 const auto& localKey = fvGeometry.feLocalCoefficients().localKey(
scv.localDofIndex());
404 template<
class LocalDof>
407 const auto type = fvGeometry.element().type();
408 const auto& localKey = fvGeometry.feLocalCoefficients().localKey(localDof.index());
409 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
419 [&] (
const typename Element::Geometry::GlobalCoordinate& pos) {
return fvGeometry.elementGeometry().local(pos); },
428 {
scvf.unitOuterNormal(),
scvf.index(), fvGeometry.elementGeometry().local(
scvf.ipGlobal()),
scvf.ipGlobal() };
432 const GGCache* ggCache_;
435 std::optional<Element> element_;
436 std::optional<typename Element::Geometry> elementGeometry_;
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
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition discretization/pq2/fvelementgeometry.hh:306
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ2FVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
Definition discretization/pq2/fvelementgeometry.hh:219
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition discretization/pq2/fvelementgeometry.hh:338
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition discretization/pq2/fvelementgeometry.hh:314
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/pq2/fvelementgeometry.hh:83
friend auto ipData(const PQ2FVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition discretization/pq2/fvelementgeometry.hh:395
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ2FVElementGeometry &fvGeometry)
Definition discretization/pq2/fvelementgeometry.hh:209
friend auto ipData(const PQ2FVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition discretization/pq2/fvelementgeometry.hh:405
void bind(const Element &element) &
Definition discretization/pq2/fvelementgeometry.hh:280
BoundaryFace::Traits::Geometry geometry(const BoundaryFace &boundaryFace) const
Geometry of a boundary face.
Definition discretization/pq2/fvelementgeometry.hh:383
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition discretization/pq2/fvelementgeometry.hh:358
PQ2FVElementGeometry(const GGCache &ggCache)
Constructor.
Definition discretization/pq2/fvelementgeometry.hh:94
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const PQ2FVElementGeometry &fvGeometry)
Definition discretization/pq2/fvelementgeometry.hh:116
friend auto ipData(const PQ2FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition discretization/pq2/fvelementgeometry.hh:425
const auto & feLocalCoefficients() const
Get a local finite element basis.
Definition discretization/pq2/fvelementgeometry.hh:243
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition discretization/pq2/fvelementgeometry.hh:75
friend auto nonCVLocalDofs(const PQ2FVElementGeometry &fvGeometry)
Definition discretization/pq2/fvelementgeometry.hh:142
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition discretization/pq2/fvelementgeometry.hh:345
friend auto ipData(const PQ2FVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition discretization/pq2/fvelementgeometry.hh:415
PQ2FVElementGeometry 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/pq2/fvelementgeometry.hh:271
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/pq2/fvelementgeometry.hh:255
std::size_t numScvf() const
The total number of sub control volume faces.
Definition discretization/pq2/fvelementgeometry.hh:261
void bindElement(const Element &element) &
Definition discretization/pq2/fvelementgeometry.hh:297
static constexpr std::size_t maxNumElementDofs
the maximum number of scvs per element for hypercubes
Definition discretization/pq2/fvelementgeometry.hh:91
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition discretization/pq2/fvelementgeometry.hh:322
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition discretization/pq2/fvelementgeometry.hh:81
bool hasBoundaryFaces() const
Returns whether the element has boundary faces.
Definition discretization/pq2/fvelementgeometry.hh:326
std::size_t elementIndex() const
The bound element index.
Definition discretization/pq2/fvelementgeometry.hh:334
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition discretization/pq2/fvelementgeometry.hh:105
friend auto localDofs(const PQ2FVElementGeometry &fvGeometry)
Definition discretization/pq2/fvelementgeometry.hh:161
const BoundaryFace & boundaryFace(LocalIndexType bfIdx) const
Get a boundary face with a local boundary face index.
Definition discretization/pq2/fvelementgeometry.hh:330
friend auto cvLocalDofs(const PQ2FVElementGeometry &fvGeometry)
Definition discretization/pq2/fvelementgeometry.hh:123
typename GG::ElementQuadratureRule ElementQuadratureRule
the quadrature rule type for elements
Definition discretization/pq2/fvelementgeometry.hh:85
typename GG::IntersectionQuadratureRule IntersectionQuadratureRule
the quadrature rule type for intersections
Definition discretization/pq2/fvelementgeometry.hh:87
friend std::ranges::view auto boundaryFaces(const PQ2FVElementGeometry &fvGeometry)
Definition discretization/pq2/fvelementgeometry.hh:230
typename GG::BoundaryFace BoundaryFace
export the boundary face type
Definition discretization/pq2/fvelementgeometry.hh:79
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition discretization/pq2/fvelementgeometry.hh:99
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition discretization/pq2/fvelementgeometry.hh:249
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition discretization/pq2/fvelementgeometry.hh:73
typename GG::BoundaryFaceQuadratureRule BoundaryFaceQuadratureRule
the quadrature rule type for boundary faces
Definition discretization/pq2/fvelementgeometry.hh:89
GG GridGeometry
export type of finite volume grid geometry
Definition discretization/pq2/fvelementgeometry.hh:77
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition discretization/pq2/fvelementgeometry.hh:71
const Element & element() const
The bound element.
Definition discretization/pq2/fvelementgeometry.hh:310
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition discretization/pq2/fvelementgeometry.hh:318
friend auto localDofs(const PQ2FVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
an iterator over all local dofs related to a boundary face
Definition discretization/pq2/fvelementgeometry.hh:181
PQ2FVElementGeometry 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/pq2/fvelementgeometry.hh:288
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition discretization/pq2/fvelementgeometry.hh:237
Base class for the finite volume geometry vector for pq2 models This builds up the sub control volume...
Definition discretization/pq2/fvelementgeometry.hh:48
Classes representing interpolation point data for control-volume finite element schemes.
Helper class constructing the dual grid finite volume geometries for the cvfe discretizazion method.
Class representing dofs on elements for control-volume finite element schemes.
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
unsigned int LocalIndex
Definition indextraits.hh:28