version 3.11-dev
Loading...
Searching...
No Matches
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 <span>
22#include <utility>
23#include <vector>
24
25#include <dune/common/rangeutilities.hh>
26#include <dune/geometry/type.hh>
27#include <dune/localfunctions/lagrange/pqkfactory.hh>
28
34
36
37namespace Dumux {
38
47template<class GG, bool enableGridGeometryCache>
49
51template<class GG>
52class PQ2FVElementGeometry<GG, true>
53{
54 using GridView = typename GG::GridView;
55 static constexpr int dim = GridView::dimension;
56 static constexpr int dimWorld = GridView::dimensionworld;
57 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
58 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
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;
63
64 using BaseIpData = CVFE::InterpolationPointData<
65 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
66 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
67 >;
68
69public:
71 using Element = typename GridView::template Codim<0>::Entity;
73 using SubControlVolume = typename GG::SubControlVolume;
75 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
77 using GridGeometry = GG;
79 using BoundaryFace = typename GG::BoundaryFace;
81 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
83 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
85 using ElementQuadratureRule = typename GG::ElementQuadratureRule;
87 using IntersectionQuadratureRule = typename GG::IntersectionQuadratureRule;
89 using BoundaryFaceQuadratureRule = typename GG::BoundaryFaceQuadratureRule;
91 static constexpr std::size_t maxNumElementDofs = GridGeometry::maxNumElementDofs;
92
94 PQ2FVElementGeometry(const GGCache& ggCache)
95 : ggCache_(&ggCache)
96 {}
97
99 const SubControlVolume& scv(LocalIndexType scvIdx) const
100 {
101 return ggCache_->scvs(eIdx_)[scvIdx];
102 }
103
105 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
106 {
107 return ggCache_->scvfs(eIdx_)[scvfIdx];
108 }
109
115 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
116 scvs(const PQ2FVElementGeometry& fvGeometry)
117 {
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());
121 }
122
123 friend inline auto cvLocalDofs(const PQ2FVElementGeometry& fvGeometry)
124 {
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) {
128 return CVFE::LocalDof{
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));
134 else
135 return static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), lk));
136 }(),
137 static_cast<GridIndexType>(fvGeometry.elementIndex())
138 };
139 });
140 }
141
142 friend inline auto nonCVLocalDofs(const PQ2FVElementGeometry& fvGeometry)
143 {
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) {
147 return CVFE::LocalDof{
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));
153 else
154 return static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), lk));
155 }(),
156 static_cast<GridIndexType>(fvGeometry.elementIndex())
157 };
158 });
159 }
160
161 friend inline auto localDofs(const PQ2FVElementGeometry& fvGeometry)
162 {
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));
171 else
172 return static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), lk));
173 }(),
174 static_cast<GridIndexType>(fvGeometry.elementIndex())
175 };
176 }
177 );
178 }
179
181 friend inline auto localDofs(const PQ2FVElementGeometry& fvGeometry, const BoundaryFace& boundaryFace)
182 {
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(),
186 boundaryFace.intersectionIndex(),
187 fvGeometry.feLocalCoefficients().localKey(i)); })
188 | std::views::transform([&](size_t i) {
189 return CVFE::LocalDof{
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));
195 else
196 return static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), lk));
197 }(),
198 static_cast<GridIndexType>(fvGeometry.elementIndex())
199 };
200 });
201 }
202
208 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
209 scvfs(const PQ2FVElementGeometry& fvGeometry)
210 {
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());
214 }
215
218 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
220 {
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]);
225 }
226
229 friend inline std::ranges::view auto
231 {
232 const auto& v = fvGeometry.ggCache_->boundaryFaces(fvGeometry.eIdx_);
233 return std::ranges::views::all(v);
234 }
235
237 const FeLocalBasis& feLocalBasis() const
238 {
239 return gridGeometry().feCache().get(element_->type()).localBasis();
240 }
241
243 const auto& feLocalCoefficients() const
244 {
245 return gridGeometry().feCache().get(element_->type()).localCoefficients();
246 }
247
249 std::size_t numLocalDofs() const
250 {
251 return feLocalCoefficients().size();
252 }
253
255 std::size_t numScv() const
256 {
257 return ggCache_->scvs(eIdx_).size();
258 }
259
261 std::size_t numScvf() const
262 {
263 return ggCache_->scvfs(eIdx_).size();
264 }
265
272 {
273 this->bindElement(element);
274 return std::move(*this);
275 }
276
280 void bind(const Element& element) &
281 { this->bindElement(element); }
282
289 {
290 this->bindElement(element);
291 return std::move(*this);
292 }
293
298 {
299 element_ = element;
300 // cache element index
301 eIdx_ = gridGeometry().elementMapper().index(element);
302 elementGeometry_.emplace(element.geometry());
303 }
304
306 bool isBound() const
307 { return static_cast<bool>(element_); }
308
310 const Element& element() const
311 { return *element_; }
312
314 const typename Element::Geometry& elementGeometry() const
315 { return *elementGeometry_; }
316
319 { return ggCache_->gridGeometry(); }
320
322 bool hasBoundaryScvf() const
323 { return ggCache_->hasBoundaryScvf(eIdx_); }
324
326 bool hasBoundaryFaces() const
327 { return hasBoundaryScvf(); }
328
330 const BoundaryFace& boundaryFace(LocalIndexType bfIdx) const
331 { return ggCache_->boundaryFaces(eIdx_)[bfIdx]; }
332
334 std::size_t elementIndex() const
335 { return eIdx_; }
336
339 {
340 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
341 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
342 }
343
345 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
346 {
347 assert(isBound());
348 const auto& geo = elementGeometry();
349 const GeometryHelper helper(geo);
350 const auto scvIdx = this->feLocalCoefficients().localKey(scv.localDofIndex()).subEntity();
351 return {
352 helper.getScvGeometryType(scvIdx),
353 helper.getScvCorners(scvIdx)
354 };
355 }
356
358 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
359 {
360 assert(isBound());
361 if (scvf.boundary())
362 {
363 GeometryHelper helper(*elementGeometry_);
364 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
365 const auto [localFacetIndex, isScvfLocalIdx]
366 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
367 return {
368 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
369 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
370 };
371 }
372 else
373 {
374 GeometryHelper helper(*elementGeometry_);
375 return {
376 helper.getInteriorScvfGeometryType(scvf.index()),
377 helper.getScvfCorners(scvf.index())
378 };
379 }
380 }
381
383 typename BoundaryFace::Traits::Geometry geometry(const BoundaryFace& boundaryFace) const
384 {
385 assert(isBound());
386 const auto& elemGeo = elementGeometry();
387 const auto faceGeoInRef = referenceElement(elemGeo).template geometry<1>(boundaryFace.intersectionIndex());
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 };
392 }
393
395 friend inline auto ipData(const PQ2FVElementGeometry& fvGeometry, const SubControlVolume& scv)
396 {
397 const auto type = fvGeometry.element().type();
398 const auto& localKey = fvGeometry.feLocalCoefficients().localKey(scv.localDofIndex());
399
400 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
401 }
402
404 template<class LocalDof>
405 friend inline auto ipData(const PQ2FVElementGeometry& fvGeometry, const LocalDof& localDof)
406 {
407 const auto type = fvGeometry.element().type();
408 const auto& localKey = fvGeometry.feLocalCoefficients().localKey(localDof.index());
409 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
410
411 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
412 }
413
415 friend inline auto ipData(const PQ2FVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
416 {
417 // Create ipData that does not automatically calculate the local position but only if it is called
419 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
420 globalPos
421 };
422 }
423
425 friend inline auto ipData(const PQ2FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
426 {
428 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
429 }
430
431private:
432 const GGCache* ggCache_;
433 GridIndexType eIdx_;
434
435 std::optional<Element> element_;
436 std::optional<typename Element::Geometry> elementGeometry_;
437};
438
439} // end namespace Dumux
440
441#endif
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.
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