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>
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;
58 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
63 using Element =
typename GridView::template Codim<0>::Entity;
71 static constexpr std::size_t maxNumElementDofs = GridGeometry::maxNumElementDofs;
81 return ggCache_->scvs(eIdx_)[scvIdx];
87 return ggCache_->scvfs(eIdx_)[scvfIdx];
95 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
98 using Iter =
typename std::vector<SubControlVolume>::const_iterator;
99 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
100 return Dune::IteratorRange<Iter>(s.begin(), s.end());
106 return Dune::transformedRangeView(
107 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type())),
108 [&](
const auto i) { return CVFE::LocalDof
110 static_cast<LocalIndexType>(i),
111 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), i)),
112 static_cast<GridIndexType>(fvGeometry.elementIndex())
118 template<
bool enable = Gr
idGeometry::enableHybr
idCVFE, std::enable_if_t<enable,
int> = 0>
121 return Dune::transformedRangeView(
122 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type()), fvGeometry.numLocalDofs()),
123 [&](
const auto i) { return CVFE::LocalDof
125 static_cast<LocalIndexType>(i),
126 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), i)),
127 static_cast<GridIndexType>(fvGeometry.elementIndex())
135 return Dune::transformedRangeView(
136 Dune::range(fvGeometry.numLocalDofs()),
137 [&](
const auto i) { return CVFE::LocalDof
139 static_cast<LocalIndexType>(i),
140 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), i)),
141 static_cast<GridIndexType>(fvGeometry.elementIndex())
147 template<
class Intersection>
150 return Dune::transformedRangeView(
151 Dune::range(GeometryHelper::numLocalDofsIntersection(fvGeometry.element().type(), intersection.indexInInside())),
154 auto localDofIdx = GeometryHelper::localDofIndexIntersection(fvGeometry.element().type(), intersection.indexInInside(), i);
155 return CVFE::LocalDof
157 static_cast<LocalIndexType>(localDofIdx),
158 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), localDofIdx)),
159 static_cast<GridIndexType>(fvGeometry.elementIndex())
170 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
173 using Iter =
typename std::vector<SubControlVolumeFace>::const_iterator;
174 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
175 return Dune::IteratorRange<Iter>(s.begin(), s.end());
181 return gridGeometry().feCache().get(element_->type()).localBasis();
187 return GeometryHelper::numElementDofs(element().type());
193 return ggCache_->scvs(eIdx_).size();
199 return ggCache_->scvfs(eIdx_).size();
209 this->bindElement(element);
210 return std::move(*
this);
217 { this->bindElement(element); }
226 this->bindElement(element);
227 return std::move(*
this);
237 eIdx_ = gridGeometry().elementMapper().index(element);
238 elementGeometry_.emplace(element.geometry());
243 {
return static_cast<bool>(element_); }
247 {
return *element_; }
251 {
return *elementGeometry_; }
255 {
return ggCache_->gridGeometry(); }
259 {
return ggCache_->hasBoundaryScvf(eIdx_); }
268 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
269 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
275 if (scv.isOverlapping())
276 DUNE_THROW(Dune::NotImplemented,
"Geometry of overlapping scv");
279 const GeometryHelper helper(*elementGeometry_);
281 helper.getScvGeometryType(scv.indexInElement()),
282 helper.getScvCorners(scv.indexInElement())
292 GeometryHelper helper(*elementGeometry_);
293 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
294 const auto [localFacetIndex, isScvfLocalIdx]
295 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
297 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
298 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
303 GeometryHelper helper(*elementGeometry_);
305 helper.getInteriorScvfGeometryType(scvf.index()),
306 helper.getScvfCorners(scvf.index())
314 const auto type = fvGeometry.element().type();
315 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
317 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
321 template<
class LocalDof>
324 const auto type = fvGeometry.element().type();
325 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
326 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
328 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
336 [&] (
const typename Element::Geometry::GlobalCoordinate& pos) {
return fvGeometry.elementGeometry().local(pos); },
342 const GGCache* ggCache_;
345 std::optional<Element> element_;
346 std::optional<typename Element::Geometry> elementGeometry_;
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition: cvfe/interpolationpointdata.hh:72
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition: cvfe/interpolationpointdata.hh:50
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:96
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:171
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/pq1bubble/fvelementgeometry.hh:63
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/pq1bubble/fvelementgeometry.hh:258
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:119
PQ1BubbleFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/pq1bubble/fvelementgeometry.hh:74
friend IpData ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/pq1bubble/fvelementgeometry.hh:322
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:79
friend auto localDofs(const PQ1BubbleFVElementGeometry &fvGeometry)
an iterator over all local dofs
Definition: discretization/pq1bubble/fvelementgeometry.hh:133
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/pq1bubble/fvelementgeometry.hh:69
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/pq1bubble/fvelementgeometry.hh:185
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/pq1bubble/fvelementgeometry.hh:179
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/pq1bubble/fvelementgeometry.hh:67
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/pq1bubble/fvelementgeometry.hh:65
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:224
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/pq1bubble/fvelementgeometry.hh:250
void bindElement(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:233
void bind(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:216
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:207
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/pq1bubble/fvelementgeometry.hh:266
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/pq1bubble/fvelementgeometry.hh:287
friend auto ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/pq1bubble/fvelementgeometry.hh:332
friend auto localDofs(const PQ1BubbleFVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/pq1bubble/fvelementgeometry.hh:148
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/pq1bubble/fvelementgeometry.hh:191
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/pq1bubble/fvelementgeometry.hh:273
friend auto cvLocalDofs(const PQ1BubbleFVElementGeometry &fvGeometry)
iterate over dof indices that belong to dofs associated with control volumes
Definition: discretization/pq1bubble/fvelementgeometry.hh:104
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:85
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/pq1bubble/fvelementgeometry.hh:242
friend IpData ipData(const PQ1BubbleFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/pq1bubble/fvelementgeometry.hh:312
std::size_t elementIndex() const
The bound element index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:262
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/pq1bubble/fvelementgeometry.hh:197
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/pq1bubble/fvelementgeometry.hh:254
const Element & element() const
The bound element.
Definition: discretization/pq1bubble/fvelementgeometry.hh:246
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.
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