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>
42template<
class GG,
bool enableGr
idGeometryCache>
49 using GridView =
typename GG::GridView;
50 static constexpr int dim = GridView::dimension;
51 static constexpr int dimWorld = GridView::dimensionworld;
54 using CoordScalar =
typename GridView::ctype;
55 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
56 using GGCache =
typename GG::Cache;
57 using GeometryHelper =
typename GGCache::GeometryHelper;
60 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
61 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
66 using Element =
typename GridView::template Codim<0>::Entity;
96 return ggCache_->scvs(eIdx_)[scvIdx];
102 return ggCache_->scvfs(eIdx_)[scvfIdx];
110 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
113 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
114 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
115 return Dune::IteratorRange<Iter>(s.begin(), s.end());
121 return Dune::transformedRangeView(
122 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type())),
123 [&](
const auto i) { return CVFE::LocalDof
125 static_cast<LocalIndexType>(i),
126 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
127 fvGeometry.feLocalCoefficients().localKey(i))),
128 static_cast<GridIndexType>(fvGeometry.elementIndex())
134 template<
bool enable = Gr
idGeometry::enableHybr
idCVFE, std::enable_if_t<enable,
int> = 0>
137 return Dune::transformedRangeView(
138 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type()), fvGeometry.numLocalDofs()),
139 [&](
const auto i) { return CVFE::LocalDof
141 static_cast<LocalIndexType>(i),
142 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
143 fvGeometry.feLocalCoefficients().localKey(i))),
144 static_cast<GridIndexType>(fvGeometry.elementIndex())
152 return Dune::transformedRangeView(
153 Dune::range(fvGeometry.numLocalDofs()),
154 [&](
const auto i) { return CVFE::LocalDof
156 static_cast<LocalIndexType>(i),
157 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
158 fvGeometry.feLocalCoefficients().localKey(i))),
159 static_cast<GridIndexType>(fvGeometry.elementIndex())
167 return Dune::transformedRangeView(
168 Dune::range(GeometryHelper::numLocalDofsIntersection(fvGeometry.element().type(),
boundaryFace.intersectionIndex())),
171 auto localDofIdx = GeometryHelper::localDofIndexIntersection(fvGeometry.element().type(), boundaryFace.intersectionIndex(), i);
172 return CVFE::LocalDof
174 static_cast<LocalIndexType>(localDofIdx),
175 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
176 fvGeometry.feLocalCoefficients().localKey(localDofIdx))),
177 static_cast<GridIndexType>(fvGeometry.elementIndex())
188 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
191 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
192 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
193 return Dune::IteratorRange<Iter>(s.begin(), s.end());
198 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
201 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
202 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
203 const auto& range = fvGeometry.ggCache_->boundaryFaceScvfRanges(fvGeometry.eIdx_)[
boundaryFace.index()];
204 return Dune::IteratorRange<Iter>(s.begin() + range[0], s.begin() + range[0] + range[1]);
209 friend inline std::ranges::view
auto
212 const auto& v = fvGeometry.ggCache_->boundaryFaces(fvGeometry.eIdx_);
213 return std::ranges::views::all(v);
219 return gridGeometry().feCache().get(element_->type()).localBasis();
225 return gridGeometry().feCache().get(element_->type()).localCoefficients();
231 return GeometryHelper::numElementDofs(
element().type());
237 return ggCache_->scvs(eIdx_).size();
243 return ggCache_->scvfs(eIdx_).size();
254 return std::move(*
this);
271 return std::move(*
this);
282 elementGeometry_.emplace(
element.geometry());
287 {
return static_cast<bool>(element_); }
291 {
return *element_; }
295 {
return *elementGeometry_; }
299 {
return ggCache_->gridGeometry(); }
303 {
return ggCache_->hasBoundaryScvf(eIdx_); }
311 {
return ggCache_->boundaryFaces(eIdx_)[bfIdx]; }
320 const auto localScvfIdx =
scvf.index() - GeometryHelper::numInteriorScvf(
element().type());
321 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
327 if (
scv.isOverlapping())
328 DUNE_THROW(Dune::NotImplemented,
"Geometry of overlapping scv");
331 const GeometryHelper helper(*elementGeometry_);
333 helper.getScvGeometryType(
scv.indexInElement()),
334 helper.getScvCorners(
scv.indexInElement())
344 GeometryHelper helper(*elementGeometry_);
345 const auto localScvfIdx =
scvf.index() - GeometryHelper::numInteriorScvf(
element().type());
346 const auto [localFacetIndex, isScvfLocalIdx]
347 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
349 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
350 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
355 GeometryHelper helper(*elementGeometry_);
357 helper.getInteriorScvfGeometryType(
scvf.index()),
358 helper.getScvfCorners(
scvf.index())
369 typename BoundaryFace::Traits::CornerStorage corners;
370 for (
int i = 0; i < faceGeoInRef.corners(); ++i)
371 corners.push_back(elemGeo.global(faceGeoInRef.corner(i)));
372 return { faceGeoInRef.type(), corners };
378 const auto type = fvGeometry.element().type();
379 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(
scv.localDofIndex());
385 template<
class LocalDof>
388 const auto type = fvGeometry.element().type();
389 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
390 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
400 [&] (
const typename Element::Geometry::GlobalCoordinate& pos) {
return fvGeometry.elementGeometry().local(pos); },
409 {
scvf.unitOuterNormal(),
scvf.index(), fvGeometry.elementGeometry().local(
scvf.ipGlobal()),
scvf.ipGlobal() };
413 const GGCache* ggCache_;
416 std::optional<Element> element_;
417 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
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition cvfe/interpolationpointdata.hh:60
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition discretization/pq1bubble/fvelementgeometry.hh:111
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition discretization/pq1bubble/fvelementgeometry.hh:189
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition discretization/pq1bubble/fvelementgeometry.hh:66
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition discretization/pq1bubble/fvelementgeometry.hh:302
friend auto nonCVLocalDofs(const PQ1BubbleFVElementGeometry &fvGeometry)
iterate over dof indices that are treated as hybrid dofs using the finite element method
Definition discretization/pq1bubble/fvelementgeometry.hh:135
typename GG::BoundaryFace BoundaryFace
export the boundary face type
Definition discretization/pq1bubble/fvelementgeometry.hh:74
friend std::ranges::view auto boundaryFaces(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition discretization/pq1bubble/fvelementgeometry.hh:210
bool hasBoundaryFaces() const
Returns whether the element has boundary faces.
Definition discretization/pq1bubble/fvelementgeometry.hh:306
PQ1BubbleFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition discretization/pq1bubble/fvelementgeometry.hh:89
typename GG::IntersectionQuadratureRule IntersectionQuadratureRule
the quadrature rule type for intersections
Definition discretization/pq1bubble/fvelementgeometry.hh:82
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition discretization/pq1bubble/fvelementgeometry.hh:94
friend auto localDofs(const PQ1BubbleFVElementGeometry &fvGeometry)
an iterator over all local dofs
Definition discretization/pq1bubble/fvelementgeometry.hh:150
GG GridGeometry
export type of finite volume grid geometry
Definition discretization/pq1bubble/fvelementgeometry.hh:72
friend auto ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition discretization/pq1bubble/fvelementgeometry.hh:386
const BoundaryFace & boundaryFace(LocalIndexType bfIdx) const
Get a boundary face with a local boundary face index.
Definition discretization/pq1bubble/fvelementgeometry.hh:310
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition discretization/pq1bubble/fvelementgeometry.hh:229
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition discretization/pq1bubble/fvelementgeometry.hh:217
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition discretization/pq1bubble/fvelementgeometry.hh:70
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition discretization/pq1bubble/fvelementgeometry.hh:68
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:268
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition discretization/pq1bubble/fvelementgeometry.hh:294
void bindElement(const Element &element) &
Definition discretization/pq1bubble/fvelementgeometry.hh:277
friend auto ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition discretization/pq1bubble/fvelementgeometry.hh:376
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ1BubbleFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
Definition discretization/pq1bubble/fvelementgeometry.hh:199
friend auto localDofs(const PQ1BubbleFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
an iterator over all local dofs related to a boundary face
Definition discretization/pq1bubble/fvelementgeometry.hh:165
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition discretization/pq1bubble/fvelementgeometry.hh:76
friend auto ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition discretization/pq1bubble/fvelementgeometry.hh:406
void bind(const Element &element) &
Definition discretization/pq1bubble/fvelementgeometry.hh:260
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:251
typename GG::BoundaryFaceQuadratureRule BoundaryFaceQuadratureRule
the quadrature rule type for boundary faces
Definition discretization/pq1bubble/fvelementgeometry.hh:84
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition discretization/pq1bubble/fvelementgeometry.hh:318
typename GG::ElementQuadratureRule ElementQuadratureRule
the quadrature rule type for elements
Definition discretization/pq1bubble/fvelementgeometry.hh:80
static constexpr std::size_t maxNumElementDofs
the maximum number of scvs per element (2^dim for cubes)
Definition discretization/pq1bubble/fvelementgeometry.hh:86
const auto & feLocalCoefficients() const
Get the local finite element coefficients.
Definition discretization/pq1bubble/fvelementgeometry.hh:223
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition discretization/pq1bubble/fvelementgeometry.hh:339
friend auto ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition discretization/pq1bubble/fvelementgeometry.hh:396
std::size_t numScv() const
The total number of sub control volumes.
Definition discretization/pq1bubble/fvelementgeometry.hh:235
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition discretization/pq1bubble/fvelementgeometry.hh:325
friend auto cvLocalDofs(const PQ1BubbleFVElementGeometry &fvGeometry)
iterate over dof indices that belong to dofs associated with control volumes
Definition discretization/pq1bubble/fvelementgeometry.hh:119
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition discretization/pq1bubble/fvelementgeometry.hh:100
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition discretization/pq1bubble/fvelementgeometry.hh:286
std::size_t elementIndex() const
The bound element index.
Definition discretization/pq1bubble/fvelementgeometry.hh:314
BoundaryFace::Traits::Geometry geometry(const BoundaryFace &boundaryFace) const
Geometry of a boundary face.
Definition discretization/pq1bubble/fvelementgeometry.hh:364
std::size_t numScvf() const
The total number of sub control volume faces.
Definition discretization/pq1bubble/fvelementgeometry.hh:241
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition discretization/pq1bubble/fvelementgeometry.hh:298
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/pq1bubble/fvelementgeometry.hh:78
const Element & element() const
The bound element.
Definition discretization/pq1bubble/fvelementgeometry.hh:290
Base class for the finite volume geometry vector for pq1bubble models This builds up the sub control ...
Definition discretization/pq1bubble/fvelementgeometry.hh:43
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