version 3.11-dev
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 <utility>
19#include <dune/common/exceptions.hh>
20#include <dune/geometry/type.hh>
21#include <dune/localfunctions/lagrange/pqkfactory.hh>
22
28
30
31namespace Dumux {
32
41template<class GG, bool enableGridGeometryCache>
43
45template<class GG>
47{
48 using GridView = typename GG::GridView;
49 static constexpr int dim = GridView::dimension;
50 static constexpr int dimWorld = GridView::dimensionworld;
51 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
52 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
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;
57
59 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
60 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
61 >;
62
63public:
65 using Element = typename GridView::template Codim<0>::Entity;
67 using SubControlVolume = typename GG::SubControlVolume;
69 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
71 using GridGeometry = GG;
73 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
75 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
77 using ElementQuadratureRule = typename GG::ElementQuadratureRule;
79 using IntersectionQuadratureRule = typename GG::IntersectionQuadratureRule;
81 static constexpr std::size_t maxNumElementDofs = GridGeometry::maxNumElementDofs;
82
84 PQ1BubbleFVElementGeometry(const GGCache& ggCache)
85 : ggCache_(&ggCache)
86 {}
87
89 const SubControlVolume& scv(LocalIndexType scvIdx) const
90 {
91 return ggCache_->scvs(eIdx_)[scvIdx];
92 }
93
95 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
96 {
97 return ggCache_->scvfs(eIdx_)[scvfIdx];
98 }
99
105 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
107 {
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());
111 }
112
114 friend inline auto cvLocalDofs(const PQ1BubbleFVElementGeometry& fvGeometry)
115 {
116 return Dune::transformedRangeView(
117 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type())),
118 [&](const auto i) { return CVFE::LocalDof
119 {
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())
124 }; }
125 );
126 }
127
129 template<bool enable = GridGeometry::enableHybridCVFE, std::enable_if_t<enable, int> = 0>
130 friend inline auto nonCVLocalDofs(const PQ1BubbleFVElementGeometry& fvGeometry)
131 {
132 return Dune::transformedRangeView(
133 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type()), fvGeometry.numLocalDofs()),
134 [&](const auto i) { return CVFE::LocalDof
135 {
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())
140 }; }
141 );
142 }
143
145 friend inline auto localDofs(const PQ1BubbleFVElementGeometry& fvGeometry)
146 {
147 return Dune::transformedRangeView(
148 Dune::range(fvGeometry.numLocalDofs()),
149 [&](const auto i) { return CVFE::LocalDof
150 {
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
160 template<class Intersection>
161 friend inline auto localDofs(const PQ1BubbleFVElementGeometry& fvGeometry, const Intersection& intersection)
162 {
163 return Dune::transformedRangeView(
164 Dune::range(GeometryHelper::numLocalDofsIntersection(fvGeometry.element().type(), intersection.indexInInside())),
165 [&](const auto i)
166 {
167 auto localDofIdx = GeometryHelper::localDofIndexIntersection(fvGeometry.element().type(), intersection.indexInInside(), i);
168 return CVFE::LocalDof
169 {
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())
174 };
175 }
176 );
177 }
178
184 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
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 GeometryHelper::numElementDofs(element().type());
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 if (scv.isOverlapping())
296 DUNE_THROW(Dune::NotImplemented, "Geometry of overlapping scv");
297
298 assert(isBound());
299 const GeometryHelper helper(*elementGeometry_);
300 return {
301 helper.getScvGeometryType(scv.indexInElement()),
302 helper.getScvCorners(scv.indexInElement())
303 };
304 }
305
307 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
308 {
309 assert(isBound());
310 if (scvf.boundary())
311 {
312 GeometryHelper helper(*elementGeometry_);
313 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
314 const auto [localFacetIndex, isScvfLocalIdx]
315 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
316 return {
317 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
318 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
319 };
320 }
321 else
322 {
323 GeometryHelper helper(*elementGeometry_);
324 return {
325 helper.getInteriorScvfGeometryType(scvf.index()),
326 helper.getScvfCorners(scvf.index())
327 };
328 }
329 }
330
332 friend inline auto ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const SubControlVolume& scv)
333 {
334 const auto type = fvGeometry.element().type();
335 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
336
337 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
338 }
339
341 template<class LocalDof>
342 friend inline auto ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const LocalDof& localDof)
343 {
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);
347
348 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
349 }
350
352 friend inline auto ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
353 {
354 // Create ipData that does not automatically calculate the local position but only if it is called
356 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
357 globalPos
358 };
359 }
360
362 friend inline auto ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
363 {
365 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
366 }
367
368private:
369 const GGCache* ggCache_;
370 GridIndexType eIdx_;
371
372 std::optional<Element> element_;
373 std::optional<typename Element::Geometry> elementGeometry_;
374};
375
376} // end namespace Dumux
377
378#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
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.
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