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
20#include <dune/common/exceptions.hh>
21#include <dune/geometry/type.hh>
22#include <dune/localfunctions/lagrange/pqkfactory.hh>
23
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;
58 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
59 LocalIndexType>;
60
61public:
63 using Element = typename GridView::template Codim<0>::Entity;
65 using SubControlVolume = typename GG::SubControlVolume;
67 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
69 using GridGeometry = GG;
71 static constexpr std::size_t maxNumElementDofs = GridGeometry::maxNumElementDofs;
72
74 PQ1BubbleFVElementGeometry(const GGCache& ggCache)
75 : ggCache_(&ggCache)
76 {}
77
79 const SubControlVolume& scv(LocalIndexType scvIdx) const
80 {
81 return ggCache_->scvs(eIdx_)[scvIdx];
82 }
83
85 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
86 {
87 return ggCache_->scvfs(eIdx_)[scvfIdx];
88 }
89
95 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
97 {
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());
101 }
102
104 friend inline auto cvLocalDofs(const PQ1BubbleFVElementGeometry& fvGeometry)
105 {
106 return Dune::transformedRangeView(
107 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type())),
108 [&](const auto i) { return CVFE::LocalDof
109 {
110 static_cast<LocalIndexType>(i),
111 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), i)),
112 static_cast<GridIndexType>(fvGeometry.elementIndex())
113 }; }
114 );
115 }
116
118 template<bool enable = GridGeometry::enableHybridCVFE, std::enable_if_t<enable, int> = 0>
119 friend inline auto nonCVLocalDofs(const PQ1BubbleFVElementGeometry& fvGeometry)
120 {
121 return Dune::transformedRangeView(
122 Dune::range(fvGeometry.numLocalDofs()-GeometryHelper::numNonCVLocalDofs(fvGeometry.element().type()), fvGeometry.numLocalDofs()),
123 [&](const auto i) { return CVFE::LocalDof
124 {
125 static_cast<LocalIndexType>(i),
126 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), i)),
127 static_cast<GridIndexType>(fvGeometry.elementIndex())
128 }; }
129 );
130 }
131
133 friend inline auto localDofs(const PQ1BubbleFVElementGeometry& fvGeometry)
134 {
135 return Dune::transformedRangeView(
136 Dune::range(fvGeometry.numLocalDofs()),
137 [&](const auto i) { return CVFE::LocalDof
138 {
139 static_cast<LocalIndexType>(i),
140 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), i)),
141 static_cast<GridIndexType>(fvGeometry.elementIndex())
142 }; }
143 );
144 }
145
147 template<class Intersection>
148 friend inline auto localDofs(const PQ1BubbleFVElementGeometry& fvGeometry, const Intersection& intersection)
149 {
150 return Dune::transformedRangeView(
151 Dune::range(GeometryHelper::numLocalDofsIntersection(fvGeometry.element().type(), intersection.indexInInside())),
152 [&](const auto i)
153 {
154 auto localDofIdx = GeometryHelper::localDofIndexIntersection(fvGeometry.element().type(), intersection.indexInInside(), i);
155 return CVFE::LocalDof
156 {
157 static_cast<LocalIndexType>(localDofIdx),
158 static_cast<GridIndexType>(GeometryHelper::dofIndex(fvGeometry.gridGeometry().dofMapper(), fvGeometry.element(), localDofIdx)),
159 static_cast<GridIndexType>(fvGeometry.elementIndex())
160 };
161 }
162 );
163 }
164
170 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
172 {
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());
176 }
177
179 const FeLocalBasis& feLocalBasis() const
180 {
181 return gridGeometry().feCache().get(element_->type()).localBasis();
182 }
183
185 std::size_t numLocalDofs() const
186 {
187 return GeometryHelper::numElementDofs(element().type());
188 }
189
191 std::size_t numScv() const
192 {
193 return ggCache_->scvs(eIdx_).size();
194 }
195
197 std::size_t numScvf() const
198 {
199 return ggCache_->scvfs(eIdx_).size();
200 }
201
208 {
209 this->bindElement(element);
210 return std::move(*this);
211 }
212
216 void bind(const Element& element) &
217 { this->bindElement(element); }
218
225 {
226 this->bindElement(element);
227 return std::move(*this);
228 }
229
233 void bindElement(const Element& element) &
234 {
235 element_ = element;
236 // cache element index
237 eIdx_ = gridGeometry().elementMapper().index(element);
238 elementGeometry_.emplace(element.geometry());
239 }
240
242 bool isBound() const
243 { return static_cast<bool>(element_); }
244
246 const Element& element() const
247 { return *element_; }
248
250 const typename Element::Geometry& elementGeometry() const
251 { return *elementGeometry_; }
252
255 { return ggCache_->gridGeometry(); }
256
258 bool hasBoundaryScvf() const
259 { return ggCache_->hasBoundaryScvf(eIdx_); }
260
262 std::size_t elementIndex() const
263 { return eIdx_; }
264
266 std::size_t intersectionIndex(const SubControlVolumeFace& scvf) const
267 {
268 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
269 return ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx][0];
270 }
271
273 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
274 {
275 if (scv.isOverlapping())
276 DUNE_THROW(Dune::NotImplemented, "Geometry of overlapping scv");
277
278 assert(isBound());
279 const GeometryHelper helper(*elementGeometry_);
280 return {
281 helper.getScvGeometryType(scv.indexInElement()),
282 helper.getScvCorners(scv.indexInElement())
283 };
284 }
285
287 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
288 {
289 assert(isBound());
290 if (scvf.boundary())
291 {
292 GeometryHelper helper(*elementGeometry_);
293 const auto localScvfIdx = scvf.index() - GeometryHelper::numInteriorScvf(element().type());
294 const auto [localFacetIndex, isScvfLocalIdx]
295 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
296 return {
297 helper.getBoundaryScvfGeometryType(isScvfLocalIdx),
298 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
299 };
300 }
301 else
302 {
303 GeometryHelper helper(*elementGeometry_);
304 return {
305 helper.getInteriorScvfGeometryType(scvf.index()),
306 helper.getScvfCorners(scvf.index())
307 };
308 }
309 }
310
312 friend inline IpData ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const SubControlVolume& scv)
313 {
314 const auto type = fvGeometry.element().type();
315 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
316
317 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
318 }
319
321 template<class LocalDof>
322 friend inline IpData ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const LocalDof& localDof)
323 {
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);
327
328 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
329 }
330
332 friend inline auto ipData(const PQ1BubbleFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
333 {
334 // Create ipData that does not automatically calculate the local position but only if it is called
336 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
337 globalPos
338 };
339 }
340
341private:
342 const GGCache* ggCache_;
343 GridIndexType eIdx_;
344
345 std::optional<Element> element_;
346 std::optional<typename Element::Geometry> elementGeometry_;
347};
348
349} // end namespace Dumux
350
351#endif
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.
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
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