version 3.11-dev
Loading...
Searching...
No Matches
discretization/pq1bubble/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_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
16
17#include <optional>
18#include <span>
19#include <utility>
20#include <dune/common/exceptions.hh>
21#include <dune/geometry/type.hh>
22#include <dune/localfunctions/lagrange/pqkfactory.hh>
23
29
31
32namespace Dumux {
33
42template<class GG, bool enableGridGeometryCache>
44
46template<class GG>
48{
49 using GridView = typename GG::GridView;
50 static constexpr int dim = GridView::dimension;
51 static constexpr int dimWorld = GridView::dimensionworld;
52 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
53 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
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;
58
59 using BaseIpData = CVFE::InterpolationPointData<
60 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
61 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
62 >;
63
64public:
66 using Element = typename GridView::template Codim<0>::Entity;
68 using SubControlVolume = typename GG::SubControlVolume;
70 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
72 using GridGeometry = GG;
74 using BoundaryFace = typename GG::BoundaryFace;
76 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
78 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
80 using ElementQuadratureRule = typename GG::ElementQuadratureRule;
82 using IntersectionQuadratureRule = typename GG::IntersectionQuadratureRule;
84 using BoundaryFaceQuadratureRule = typename GG::BoundaryFaceQuadratureRule;
86 static constexpr std::size_t maxNumElementDofs = GridGeometry::maxNumElementDofs;
87
89 PQ1BubbleFVElementGeometry(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>
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
119 friend inline auto cvLocalDofs(const PQ1BubbleFVElementGeometry& fvGeometry)
120 {
121 return Dune::transformedRangeView(
122 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type())),
123 [&](const auto i) { return CVFE::LocalDof
124 {
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())
129 }; }
130 );
131 }
132
134 template<bool enable = GridGeometry::enableHybridCVFE, std::enable_if_t<enable, int> = 0>
135 friend inline auto nonCVLocalDofs(const PQ1BubbleFVElementGeometry& fvGeometry)
136 {
137 return Dune::transformedRangeView(
138 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type()), fvGeometry.numLocalDofs()),
139 [&](const auto i) { return CVFE::LocalDof
140 {
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())
145 }; }
146 );
147 }
148
150 friend inline auto localDofs(const PQ1BubbleFVElementGeometry& fvGeometry)
151 {
152 return Dune::transformedRangeView(
153 Dune::range(fvGeometry.numLocalDofs()),
154 [&](const auto i) { return CVFE::LocalDof
155 {
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())
160 }; }
161 );
162 }
163
165 friend inline auto localDofs(const PQ1BubbleFVElementGeometry& fvGeometry, const BoundaryFace& boundaryFace)
166 {
167 return Dune::transformedRangeView(
168 Dune::range(GeometryHelper::numLocalDofsIntersection(fvGeometry.element().type(), boundaryFace.intersectionIndex())),
169 [&](const auto i)
170 {
171 auto localDofIdx = GeometryHelper::localDofIndexIntersection(fvGeometry.element().type(), boundaryFace.intersectionIndex(), i);
172 return CVFE::LocalDof
173 {
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())
178 };
179 }
180 );
181 }
182
188 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
190 {
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());
194 }
195
198 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
200 {
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]);
205 }
206
209 friend inline std::ranges::view auto
211 {
212 const auto& v = fvGeometry.ggCache_->boundaryFaces(fvGeometry.eIdx_);
213 return std::ranges::views::all(v);
214 }
215
217 const FeLocalBasis& feLocalBasis() const
218 {
219 return gridGeometry().feCache().get(element_->type()).localBasis();
220 }
221
223 const auto& feLocalCoefficients() const
224 {
225 return gridGeometry().feCache().get(element_->type()).localCoefficients();
226 }
227
229 std::size_t numLocalDofs() const
230 {
231 return GeometryHelper::numElementDofs(element().type());
232 }
233
235 std::size_t numScv() const
236 {
237 return ggCache_->scvs(eIdx_).size();
238 }
239
241 std::size_t numScvf() const
242 {
243 return ggCache_->scvfs(eIdx_).size();
244 }
245
252 {
253 this->bindElement(element);
254 return std::move(*this);
255 }
256
260 void bind(const Element& element) &
261 { this->bindElement(element); }
262
269 {
270 this->bindElement(element);
271 return std::move(*this);
272 }
273
278 {
279 element_ = element;
280 // cache element index
281 eIdx_ = gridGeometry().elementMapper().index(element);
282 elementGeometry_.emplace(element.geometry());
283 }
284
286 bool isBound() const
287 { return static_cast<bool>(element_); }
288
290 const Element& element() const
291 { return *element_; }
292
294 const typename Element::Geometry& elementGeometry() const
295 { return *elementGeometry_; }
296
299 { return ggCache_->gridGeometry(); }
300
302 bool hasBoundaryScvf() const
303 { return ggCache_->hasBoundaryScvf(eIdx_); }
304
306 bool hasBoundaryFaces() const
307 { return hasBoundaryScvf(); }
308
310 const BoundaryFace& boundaryFace(LocalIndexType bfIdx) const
311 { return ggCache_->boundaryFaces(eIdx_)[bfIdx]; }
312
314 std::size_t elementIndex() const
315 { return eIdx_; }
316
319 {
320 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
321 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
322 }
323
325 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
326 {
327 if (scv.isOverlapping())
328 DUNE_THROW(Dune::NotImplemented, "Geometry of overlapping scv");
329
330 assert(isBound());
331 const GeometryHelper helper(*elementGeometry_);
332 return {
333 helper.getScvGeometryType(scv.indexInElement()),
334 helper.getScvCorners(scv.indexInElement())
335 };
336 }
337
339 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
340 {
341 assert(isBound());
342 if (scvf.boundary())
343 {
344 GeometryHelper helper(*elementGeometry_);
345 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
346 const auto [localFacetIndex, isScvfLocalIdx]
347 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
348 return {
349 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
350 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
351 };
352 }
353 else
354 {
355 GeometryHelper helper(*elementGeometry_);
356 return {
357 helper.getInteriorScvfGeometryType(scvf.index()),
358 helper.getScvfCorners(scvf.index())
359 };
360 }
361 }
362
364 typename BoundaryFace::Traits::Geometry geometry(const BoundaryFace& boundaryFace) const
365 {
366 assert(isBound());
367 const auto& elemGeo = elementGeometry();
368 const auto faceGeoInRef = referenceElement(elemGeo).template geometry<1>(boundaryFace.intersectionIndex());
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 };
373 }
374
376 friend inline auto ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const SubControlVolume& scv)
377 {
378 const auto type = fvGeometry.element().type();
379 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
380
381 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
382 }
383
385 template<class LocalDof>
386 friend inline auto ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const LocalDof& localDof)
387 {
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);
391
392 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
393 }
394
396 friend inline auto ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
397 {
398 // Create ipData that does not automatically calculate the local position but only if it is called
400 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
401 globalPos
402 };
403 }
404
406 friend inline auto ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
407 {
409 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
410 }
411
412private:
413 const GGCache* ggCache_;
414 GridIndexType eIdx_;
415
416 std::optional<Element> element_;
417 std::optional<typename Element::Geometry> elementGeometry_;
418};
419
420} // end namespace Dumux
421
422#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
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.
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