version 3.11-dev
discretization/pq2/fvelementgeometry.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_DISCRETIZATION_PQ2_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_PQ2_FV_ELEMENT_GEOMETRY_HH
16
17#include <cstddef>
18#include <cassert>
19#include <optional>
20#include <ranges>
21#include <utility>
22#include <vector>
23
24#include <dune/common/rangeutilities.hh>
25#include <dune/geometry/type.hh>
26#include <dune/localfunctions/lagrange/pqkfactory.hh>
27
33
35
36namespace Dumux {
37
46template<class GG, bool enableGridGeometryCache>
48
50template<class GG>
51class PQ2FVElementGeometry<GG, true>
52{
53 using GridView = typename GG::GridView;
54 static constexpr int dim = GridView::dimension;
55 static constexpr int dimWorld = GridView::dimensionworld;
56 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
57 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
58 using CoordScalar = typename GridView::ctype;
59 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
60 using GGCache = typename GG::Cache;
61 using GeometryHelper = typename GGCache::GeometryHelper;
62
64 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
65 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
66 >;
67
68public:
70 using Element = typename GridView::template Codim<0>::Entity;
72 using SubControlVolume = typename GG::SubControlVolume;
74 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
76 using GridGeometry = GG;
78 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
80 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
82 using ElementQuadratureRule = typename GG::ElementQuadratureRule;
84 using IntersectionQuadratureRule = typename GG::IntersectionQuadratureRule;
86 static constexpr std::size_t maxNumElementDofs = GridGeometry::maxNumElementDofs;
87
89 PQ2FVElementGeometry(const GGCache& ggCache)
90 : ggCache_(&ggCache)
91 {}
92
94 const SubControlVolume& scv(LocalIndexType scvIdx) const
95 {
96 return ggCache_->scvs(eIdx_)[scvIdx];
97 }
98
100 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
101 {
102 return ggCache_->scvfs(eIdx_)[scvfIdx];
103 }
104
110 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
111 scvs(const PQ2FVElementGeometry& fvGeometry)
112 {
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());
116 }
117
118 friend inline auto cvLocalDofs(const PQ2FVElementGeometry& fvGeometry)
119 {
120 return std::views::iota(std::size_t(0), fvGeometry.numLocalDofs())
121 | std::views::filter([&](size_t i) { return fvGeometry.feLocalCoefficients().localKey(i).codim() == dim; })
122 | std::views::transform([&](size_t i) {
123 return CVFE::LocalDof{
124 static_cast<LocalIndexType>(i),
125 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
126 fvGeometry.feLocalCoefficients().localKey(i))),
127 static_cast<GridIndexType>(fvGeometry.elementIndex())
128 };
129 });
130 }
131
132 friend inline auto nonCVLocalDofs(const PQ2FVElementGeometry& fvGeometry)
133 {
134 return std::views::iota(std::size_t(0), fvGeometry.numLocalDofs())
135 | std::views::filter([&](size_t i) { return !(fvGeometry.feLocalCoefficients().localKey(i).codim() == dim); })
136 | std::views::transform([&](size_t i) {
137 return CVFE::LocalDof{
138 static_cast<LocalIndexType>(i),
139 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
140 fvGeometry.feLocalCoefficients().localKey(i))),
141 static_cast<GridIndexType>(fvGeometry.elementIndex())
142 };
143 });
144 }
145
146 friend inline auto localDofs(const PQ2FVElementGeometry& fvGeometry)
147 {
148 return Dune::transformedRangeView(
149 Dune::range(std::size_t(0), fvGeometry.numLocalDofs()), [&](const auto i) {
150 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())
155 };
156 }
157 );
158 }
159
161 template<class Intersection>
162 friend inline auto localDofs(const PQ2FVElementGeometry& fvGeometry, const Intersection& intersection)
163 {
164 return std::views::iota(std::size_t(0), fvGeometry.numLocalDofs())
165 | std::views::filter([&](size_t i) {
166 return GeometryHelper::localDofOnIntersection(fvGeometry.element().type(),
167 intersection.indexInInside(),
168 fvGeometry.feLocalCoefficients().localKey(i)); })
169 | std::views::transform([&](size_t i) {
170 return CVFE::LocalDof{
171 static_cast<LocalIndexType>(i),
172 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(),
173 fvGeometry.feLocalCoefficients().localKey(i))),
174 static_cast<GridIndexType>(fvGeometry.elementIndex())
175 };
176 });
177 }
178
184 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
185 scvfs(const PQ2FVElementGeometry& fvGeometry)
186 {
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());
190 }
191
193 const FeLocalBasis& feLocalBasis() const
194 {
195 return gridGeometry().feCache().get(element_->type()).localBasis();
196 }
197
199 const auto& feLocalCoefficients() const
200 {
201 return gridGeometry().feCache().get(element_->type()).localCoefficients();
202 }
203
205 std::size_t numLocalDofs() const
206 {
207 return feLocalCoefficients().size();
208 }
209
211 std::size_t numScv() const
212 {
213 return ggCache_->scvs(eIdx_).size();
214 }
215
217 std::size_t numScvf() const
218 {
219 return ggCache_->scvfs(eIdx_).size();
220 }
221
228 {
229 this->bindElement(element);
230 return std::move(*this);
231 }
232
236 void bind(const Element& element) &
237 { this->bindElement(element); }
238
245 {
246 this->bindElement(element);
247 return std::move(*this);
248 }
249
253 void bindElement(const Element& element) &
254 {
255 element_ = element;
256 // cache element index
257 eIdx_ = gridGeometry().elementMapper().index(element);
258 elementGeometry_.emplace(element.geometry());
259 }
260
262 bool isBound() const
263 { return static_cast<bool>(element_); }
264
266 const Element& element() const
267 { return *element_; }
268
270 const typename Element::Geometry& elementGeometry() const
271 { return *elementGeometry_; }
272
275 { return ggCache_->gridGeometry(); }
276
278 bool hasBoundaryScvf() const
279 { return ggCache_->hasBoundaryScvf(eIdx_); }
280
282 std::size_t elementIndex() const
283 { return eIdx_; }
284
286 std::size_t intersectionIndex(const SubControlVolumeFace& scvf) const
287 {
288 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
289 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
290 }
291
293 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
294 {
295 assert(isBound());
296 const auto& geo = elementGeometry();
297 const GeometryHelper helper(geo);
298 const auto scvIdx = this->feLocalCoefficients().localKey(scv.localDofIndex()).subEntity();
299 return {
300 helper.getScvGeometryType(scvIdx),
301 helper.getScvCorners(scvIdx)
302 };
303 }
304
306 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
307 {
308 assert(isBound());
309 if (scvf.boundary())
310 {
311 GeometryHelper helper(*elementGeometry_);
312 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
313 const auto [localFacetIndex, isScvfLocalIdx]
314 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
315 return {
316 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
317 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
318 };
319 }
320 else
321 {
322 GeometryHelper helper(*elementGeometry_);
323 return {
324 helper.getInteriorScvfGeometryType(scvf.index()),
325 helper.getScvfCorners(scvf.index())
326 };
327 }
328 }
329
331 friend inline auto ipData(const PQ2FVElementGeometry& fvGeometry, const SubControlVolume& scv)
332 {
333 const auto type = fvGeometry.element().type();
334 const auto& localKey = fvGeometry.feLocalCoefficients().localKey(scv.localDofIndex());
335
336 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
337 }
338
340 template<class LocalDof>
341 friend inline auto ipData(const PQ2FVElementGeometry& fvGeometry, const LocalDof& localDof)
342 {
343 const auto type = fvGeometry.element().type();
344 const auto& localKey = fvGeometry.feLocalCoefficients().localKey(localDof.index());
345 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
346
347 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
348 }
349
351 friend inline auto ipData(const PQ2FVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
352 {
353 // Create ipData that does not automatically calculate the local position but only if it is called
355 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
356 globalPos
357 };
358 }
359
361 friend inline auto ipData(const PQ2FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
362 {
364 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
365 }
366
367private:
368 const GGCache* ggCache_;
369 GridIndexType eIdx_;
370
371 std::optional<Element> element_;
372 std::optional<typename Element::Geometry> elementGeometry_;
373};
374
375} // end namespace Dumux
376
377#endif
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
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:54
LocalIndex localDofIndex() const
The local index of the corresponding dof.
Definition: cvfe/interpolationpointdata.hh:63
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/pq2/fvelementgeometry.hh:262
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition: discretization/pq2/fvelementgeometry.hh:286
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/pq2/fvelementgeometry.hh:270
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition: discretization/pq2/fvelementgeometry.hh:80
friend auto ipData(const PQ2FVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/pq2/fvelementgeometry.hh:331
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ2FVElementGeometry &fvGeometry)
Definition: discretization/pq2/fvelementgeometry.hh:185
friend auto ipData(const PQ2FVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/pq2/fvelementgeometry.hh:341
void bind(const Element &element) &
Definition: discretization/pq2/fvelementgeometry.hh:236
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/pq2/fvelementgeometry.hh:306
PQ2FVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/pq2/fvelementgeometry.hh:89
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const PQ2FVElementGeometry &fvGeometry)
Definition: discretization/pq2/fvelementgeometry.hh:111
friend auto ipData(const PQ2FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition: discretization/pq2/fvelementgeometry.hh:361
const auto & feLocalCoefficients() const
Get a local finite element basis.
Definition: discretization/pq2/fvelementgeometry.hh:199
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/pq2/fvelementgeometry.hh:74
friend auto nonCVLocalDofs(const PQ2FVElementGeometry &fvGeometry)
Definition: discretization/pq2/fvelementgeometry.hh:132
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/pq2/fvelementgeometry.hh:293
friend auto ipData(const PQ2FVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition: discretization/pq2/fvelementgeometry.hh:351
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:227
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/pq2/fvelementgeometry.hh:211
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/pq2/fvelementgeometry.hh:217
friend auto localDofs(const PQ2FVElementGeometry &fvGeometry, const Intersection &intersection)
an iterator over all local dofs related to an intersection
Definition: discretization/pq2/fvelementgeometry.hh:162
void bindElement(const Element &element) &
Definition: discretization/pq2/fvelementgeometry.hh:253
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/pq2/fvelementgeometry.hh:278
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition: discretization/pq2/fvelementgeometry.hh:78
std::size_t elementIndex() const
The bound element index.
Definition: discretization/pq2/fvelementgeometry.hh:282
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/pq2/fvelementgeometry.hh:100
friend auto localDofs(const PQ2FVElementGeometry &fvGeometry)
Definition: discretization/pq2/fvelementgeometry.hh:146
friend auto cvLocalDofs(const PQ2FVElementGeometry &fvGeometry)
Definition: discretization/pq2/fvelementgeometry.hh:118
typename GG::ElementQuadratureRule ElementQuadratureRule
the quadrature rule type for elements
Definition: discretization/pq2/fvelementgeometry.hh:82
typename GG::IntersectionQuadratureRule IntersectionQuadratureRule
the quadrature rule type for intersections
Definition: discretization/pq2/fvelementgeometry.hh:84
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/pq2/fvelementgeometry.hh:94
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/pq2/fvelementgeometry.hh:205
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/pq2/fvelementgeometry.hh:72
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/pq2/fvelementgeometry.hh:76
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/pq2/fvelementgeometry.hh:70
const Element & element() const
The bound element.
Definition: discretization/pq2/fvelementgeometry.hh:266
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/pq2/fvelementgeometry.hh:274
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:244
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/pq2/fvelementgeometry.hh:193
Base class for the finite volume geometry vector for pq2 models This builds up the sub control volume...
Definition: discretization/pq2/fvelementgeometry.hh:47
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.
Defines the index types used for grid and local indices.
Class representing dofs on elements for control-volume finite element schemes.
Definition: adapt.hh:17
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