14#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
19#include <dune/common/exceptions.hh>
20#include <dune/geometry/type.hh>
21#include <dune/localfunctions/lagrange/pqkfactory.hh>
41template<
class GG,
bool enableGr
idGeometryCache>
48 using GridView =
typename GG::GridView;
49 static constexpr int dim = GridView::dimension;
50 static constexpr int dimWorld = GridView::dimensionworld;
53 using CoordScalar =
typename GridView::ctype;
54 using FeLocalBasis =
typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
55 using GGCache =
typename GG::Cache;
56 using GeometryHelper =
typename GGCache::GeometryHelper;
59 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
60 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
65 using Element =
typename GridView::template Codim<0>::Entity;
81 static constexpr std::size_t maxNumElementDofs = GridGeometry::maxNumElementDofs;
91 return ggCache_->scvs(eIdx_)[scvIdx];
97 return ggCache_->scvfs(eIdx_)[scvfIdx];
105 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
108 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
109 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
110 return Dune::IteratorRange<Iter>(s.begin(), s.end());
116 return Dune::transformedRangeView(
117 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type())),
118 [&](
const auto i) { return CVFE::LocalDof
120 static_cast<LocalIndexType>(i),
121 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
122 fvGeometry.feLocalCoefficients().localKey(i))),
123 static_cast<GridIndexType>(fvGeometry.elementIndex())
129 template<
bool enable = Gr
idGeometry::enableHybr
idCVFE, std::enable_if_t<enable,
int> = 0>
132 return Dune::transformedRangeView(
133 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type()), fvGeometry.numLocalDofs()),
134 [&](
const auto i) { return CVFE::LocalDof
136 static_cast<LocalIndexType>(i),
137 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
138 fvGeometry.feLocalCoefficients().localKey(i))),
139 static_cast<GridIndexType>(fvGeometry.elementIndex())
147 return Dune::transformedRangeView(
148 Dune::range(fvGeometry.numLocalDofs()),
149 [&](
const auto i) { return CVFE::LocalDof
151 static_cast<LocalIndexType>(i),
152 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
153 fvGeometry.feLocalCoefficients().localKey(i))),
154 static_cast<GridIndexType>(fvGeometry.elementIndex())
160 template<
class Intersection>
163 return Dune::transformedRangeView(
164 Dune::range(GeometryHelper::numLocalDofsIntersection(fvGeometry.element().type(), intersection.indexInInside())),
167 auto localDofIdx = GeometryHelper::localDofIndexIntersection(fvGeometry.element().type(), intersection.indexInInside(), i);
168 return CVFE::LocalDof
170 static_cast<LocalIndexType>(localDofIdx),
171 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
172 fvGeometry.feLocalCoefficients().localKey(localDofIdx))),
173 static_cast<GridIndexType>(fvGeometry.elementIndex())
184 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
187 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
188 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
189 return Dune::IteratorRange<Iter>(s.begin(), s.end());
195 return gridGeometry().feCache().get(element_->type()).localBasis();
201 return gridGeometry().feCache().get(element_->type()).localCoefficients();
207 return GeometryHelper::numElementDofs(element().type());
213 return ggCache_->scvs(eIdx_).size();
219 return ggCache_->scvfs(eIdx_).size();
229 this->bindElement(element);
230 return std::move(*
this);
237 { this->bindElement(element); }
246 this->bindElement(element);
247 return std::move(*
this);
257 eIdx_ = gridGeometry().elementMapper().index(element);
258 elementGeometry_.emplace(element.geometry());
263 {
return static_cast<bool>(element_); }
267 {
return *element_; }
271 {
return *elementGeometry_; }
275 {
return ggCache_->gridGeometry(); }
279 {
return ggCache_->hasBoundaryScvf(eIdx_); }
288 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
289 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
295 if (scv.isOverlapping())
296 DUNE_THROW(Dune::NotImplemented,
"Geometry of overlapping scv");
299 const GeometryHelper helper(*elementGeometry_);
301 helper.getScvGeometryType(scv.indexInElement()),
302 helper.getScvCorners(scv.indexInElement())
312 GeometryHelper helper(*elementGeometry_);
313 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
314 const auto [localFacetIndex, isScvfLocalIdx]
315 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
317 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
318 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
323 GeometryHelper helper(*elementGeometry_);
325 helper.getInteriorScvfGeometryType(scvf.index()),
326 helper.getScvfCorners(scvf.index())
334 const auto type = fvGeometry.element().type();
335 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
341 template<
class LocalDof>
344 const auto type = fvGeometry.element().type();
345 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
346 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
356 [&] (
const typename Element::Geometry::GlobalCoordinate& pos) {
return fvGeometry.elementGeometry().local(pos); },
365 { scvf.
unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
369 const GGCache* ggCache_;
372 std::optional<Element> element_;
373 std::optional<typename Element::Geometry> elementGeometry_;
An interpolation point related to a face of an element.
Definition: cvfe/interpolationpointdata.hh:103
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector at the quadrature point.
Definition: cvfe/interpolationpointdata.hh:117
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:25
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition: cvfe/interpolationpointdata.hh:76
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition: cvfe/interpolationpointdata.hh:54
LocalIndex localDofIndex() const
The local index of the corresponding dof.
Definition: cvfe/interpolationpointdata.hh:63
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:106
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:185
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/pq1bubble/fvelementgeometry.hh:65
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/pq1bubble/fvelementgeometry.hh:278
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:130
PQ1BubbleFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/pq1bubble/fvelementgeometry.hh:84
typename GG::IntersectionQuadratureRule IntersectionQuadratureRule
the quadrature rule type for intersections
Definition: discretization/pq1bubble/fvelementgeometry.hh:79
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:89
friend auto localDofs(const PQ1BubbleFVElementGeometry &fvGeometry)
an iterator over all local dofs
Definition: discretization/pq1bubble/fvelementgeometry.hh:145
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/pq1bubble/fvelementgeometry.hh:71
friend auto ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/pq1bubble/fvelementgeometry.hh:342
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/pq1bubble/fvelementgeometry.hh:205
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/pq1bubble/fvelementgeometry.hh:193
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/pq1bubble/fvelementgeometry.hh:69
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/pq1bubble/fvelementgeometry.hh:67
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:244
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/pq1bubble/fvelementgeometry.hh:270
void bindElement(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:253
friend auto ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/pq1bubble/fvelementgeometry.hh:332
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition: discretization/pq1bubble/fvelementgeometry.hh:73
friend auto ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition: discretization/pq1bubble/fvelementgeometry.hh:362
void bind(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:236
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:227
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/pq1bubble/fvelementgeometry.hh:286
typename GG::ElementQuadratureRule ElementQuadratureRule
the quadrature rule type for elements
Definition: discretization/pq1bubble/fvelementgeometry.hh:77
const auto & feLocalCoefficients() const
Get the local finite element coefficients.
Definition: discretization/pq1bubble/fvelementgeometry.hh:199
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/pq1bubble/fvelementgeometry.hh:307
friend auto ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/pq1bubble/fvelementgeometry.hh:352
friend auto localDofs(const PQ1BubbleFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/pq1bubble/fvelementgeometry.hh:161
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/pq1bubble/fvelementgeometry.hh:211
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/pq1bubble/fvelementgeometry.hh:293
friend auto cvLocalDofs(const PQ1BubbleFVElementGeometry &fvGeometry)
iterate over dof indices that belong to dofs associated with control volumes
Definition: discretization/pq1bubble/fvelementgeometry.hh:114
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:95
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/pq1bubble/fvelementgeometry.hh:262
std::size_t elementIndex() const
The bound element index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:282
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/pq1bubble/fvelementgeometry.hh:217
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/pq1bubble/fvelementgeometry.hh:274
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition: discretization/pq1bubble/fvelementgeometry.hh:75
const Element & element() const
The bound element.
Definition: discretization/pq1bubble/fvelementgeometry.hh:266
Base class for the finite volume geometry vector for pq1bubble models This builds up the sub control ...
Definition: discretization/pq1bubble/fvelementgeometry.hh:42
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