version 3.11-dev
Loading...
Searching...
No Matches
discretization/facecentered/diamond/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//
12#ifndef DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
14
15#include <type_traits>
16#include <optional>
17#include <ranges>
18#include <span>
19
20#include <dune/common/reservedvector.hh>
21#include <dune/common/iteratorrange.hh>
22
28
29namespace Dumux {
30
37template<class GG, bool cachingEnabled>
39
44template<class GG>
45class FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>
46{
47 using ThisType = FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>;
48 using GridView = typename GG::GridView;
49 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
50 using LocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
51 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
52 using GGCache = typename GG::Cache;
53 using GeometryHelper = typename GGCache::GeometryHelper;
54
55 using BaseIpData = CVFE::InterpolationPointData<
56 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
57 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
58 >;
59
60public:
62 using SubControlVolume = typename GG::SubControlVolume;
63 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
64 using BoundaryFace = typename GG::BoundaryFace;
65 using Element = typename GridView::template Codim<0>::Entity;
66 using GridGeometry = GG;
68 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
70 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
71
73 static constexpr std::size_t maxNumElementScvs = 2*GridView::dimension;
74
76 : ggCache_(&ggCache)
77 {}
78
80 const SubControlVolume& scv(LocalIndexType scvIdx) const
81 { return ggCache_->scvs(eIdx_)[scvIdx]; }
82
84 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
85 { return ggCache_->scvf(eIdx_)[scvfIdx]; }
86
92 friend inline auto
94 {
95 using Iter = typename std::vector<SubControlVolume>::const_iterator;
96 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
97 return Dune::IteratorRange<Iter>(s.begin(), s.end());
98 }
99
105 friend inline auto
107 {
108 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
109 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
110 return Dune::IteratorRange<Iter>(s.begin(), s.end());
111 }
112
115 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
117 {
118 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
119 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
120 const auto& range = fvGeometry.ggCache_->boundaryFaceScvfRanges(fvGeometry.eIdx_)[boundaryFace.index()];
121 return Dune::IteratorRange<Iter>(s.begin() + range[0], s.begin() + range[0] + range[1]);
122 }
123
125 const FeLocalBasis& feLocalBasis() const
126 {
127 return gridGeometry().feCache().get(element().type()).localBasis();
128 }
129
131 std::size_t numLocalDofs() const
132 {
133 return numScv();
134 }
135
137 std::size_t numScv() const
138 {
139 return ggCache_->scvs(eIdx_).size();
140 }
141
143 std::size_t numScvf() const
144 {
145 return ggCache_->scvfs(eIdx_).size();
146 }
147
149 bool hasBoundaryScvf() const
150 { return ggCache_->hasBoundaryScvf(eIdx_); }
151
153 bool hasBoundaryFaces() const
154 { return hasBoundaryScvf(); }
155
157 const BoundaryFace& boundaryFace(LocalIndexType bfIdx) const
158 { return ggCache_->boundaryFaces(eIdx_)[bfIdx]; }
159
161 friend inline std::ranges::view auto
163 {
164 const auto& v = fvGeometry.ggCache_->boundaryFaces(fvGeometry.eIdx_);
165 return std::ranges::views::all(v);
166 }
167
174 {
175 this->bindElement(element);
176 return std::move(*this);
177 }
178
179 void bind(const Element& element) &
180 {
181 this->bindElement(element);
182 }
183
190 {
191 this->bindElement(element);
192 return std::move(*this);
193 }
194
197 {
198 element_ = element;
199 elementGeometry_.emplace(element.geometry());
200 eIdx_ = gridGeometry().elementMapper().index(element);
201 }
202
204 bool isBound() const
205 { return static_cast<bool>(element_); }
206
208 const Element& element() const
209 { return *element_; }
210
212 const typename Element::Geometry& elementGeometry() const
213 { return *elementGeometry_; }
214
217 { return ggCache_->gridGeometry(); }
218
220 std::size_t elementIndex() const
221 { return eIdx_; }
222
225 {
226 // Since scvs are constructed around intersections, they have the same index
227 return scvf.insideScvIdx();
228 }
229
231 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
232 {
233 assert(isBound());
234 return {
235 SubControlVolume::Traits::geometryType((*elementGeometry_).type()),
236 GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement())
237 };
238 }
239
241 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
242 {
243 assert(isBound());
244 if (scvf.boundary())
245 {
246 // use the information that each boundary scvf corresponds to one scv constructed around the same facet
247 const auto localFacetIndex = scvf.insideScvIdx();
248 return {
249 referenceElement(*elementGeometry_).type(localFacetIndex, 1),
250 GeometryHelper(*elementGeometry_).getBoundaryScvfCorners(localFacetIndex)
251 };
252 }
253 else
254 {
255 return {
256 SubControlVolumeFace::Traits::interiorGeometryType((*elementGeometry_).type()),
257 GeometryHelper(*elementGeometry_).getScvfCorners(scvf.index())
258 };
259 }
260 }
261
263 typename BoundaryFace::Traits::Geometry geometry(const BoundaryFace& boundaryFace) const
264 {
265 assert(isBound());
266 const auto& elemGeo = elementGeometry();
267 const auto faceGeoInRef = referenceElement(elemGeo).template geometry<1>(boundaryFace.intersectionIndex());
268 typename BoundaryFace::Traits::CornerStorage corners;
269 for (int i = 0; i < faceGeoInRef.corners(); ++i)
270 corners.push_back(elemGeo.global(faceGeoInRef.corner(i)));
271 return { faceGeoInRef.type(), corners };
272 }
273
275 friend inline auto localDofs(const FaceCenteredDiamondFVElementGeometry& fvGeometry,
277 {
278 const auto localDofIdx = boundaryFace.intersectionIndex();
279 return std::views::single(CVFE::LocalDof{
280 static_cast<LocalIndexType>(localDofIdx),
281 static_cast<GridIndexType>(fvGeometry.scv(localDofIdx).dofIndex()),
282 static_cast<GridIndexType>(fvGeometry.elementIndex())
283 });
284 }
285
287 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const SubControlVolume& scv)
288 {
289 const auto type = fvGeometry.element().type();
290 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
291
292 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
293 }
294
296 template<class LocalDof>
297 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const LocalDof& localDof)
298 {
299 const auto type = fvGeometry.element().type();
300 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
301 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
302
303 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
304 }
305
307 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
308 {
309 // Create ipData that does not automatically calculate the local position but only if it is called
311 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
312 globalPos
313 };
314 }
315
317 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
318 {
320 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
321 }
322
323private:
324 std::optional<Element> element_;
325 std::optional<typename Element::Geometry> elementGeometry_;
326 GridIndexType eIdx_;
327 const GGCache* ggCache_;
328};
329
330} // end namespace Dumux
331
332#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
friend auto scvs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition discretization/facecentered/diamond/fvelementgeometry.hh:93
void bindElement(const Element &element) &
Bind only element-local.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:196
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:216
static constexpr std::size_t maxNumElementScvs
the maximum number of scvs per element
Definition discretization/facecentered/diamond/fvelementgeometry.hh:73
const Element & element() const
The bound element.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:208
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition discretization/facecentered/diamond/fvelementgeometry.hh:68
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition discretization/facecentered/diamond/fvelementgeometry.hh:62
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:297
std::size_t intersectionIndex(const SubControlVolumeFace &scvf) const
The intersection index the scvf belongs to.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:224
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:212
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:317
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:149
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition discretization/facecentered/diamond/fvelementgeometry.hh:63
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition discretization/facecentered/diamond/fvelementgeometry.hh:70
bool hasBoundaryFaces() const
Returns whether the element has boundary faces.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:153
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:287
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:204
const BoundaryFace & boundaryFace(LocalIndexType bfIdx) const
Get a boundary face with a local boundary face index.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:157
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition discretization/facecentered/diamond/fvelementgeometry.hh:143
friend auto scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition discretization/facecentered/diamond/fvelementgeometry.hh:106
std::size_t elementIndex() const
The bound element index.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:220
GG GridGeometry
Definition discretization/facecentered/diamond/fvelementgeometry.hh:66
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const typename Element::Geometry::GlobalCoordinate &globalPos)
Interpolation point data for a global position.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:307
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
Definition discretization/facecentered/diamond/fvelementgeometry.hh:116
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:84
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:80
typename GridView::template Codim< 0 >::Entity Element
Definition discretization/facecentered/diamond/fvelementgeometry.hh:65
friend auto localDofs(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const BoundaryFace &boundaryFace)
an iterator over all local dofs related to a boundary face
Definition discretization/facecentered/diamond/fvelementgeometry.hh:275
typename GG::BoundaryFace BoundaryFace
Definition discretization/facecentered/diamond/fvelementgeometry.hh:64
BoundaryFace::Traits::Geometry geometry(const BoundaryFace &boundaryFace) const
Geometry of a boundary face.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:263
FaceCenteredDiamondFVElementGeometry 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/facecentered/diamond/fvelementgeometry.hh:173
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition discretization/facecentered/diamond/fvelementgeometry.hh:137
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:131
FaceCenteredDiamondFVElementGeometry(const GGCache &ggCache)
Definition discretization/facecentered/diamond/fvelementgeometry.hh:75
friend std::ranges::view auto boundaryFaces(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
iterator range for boundary faces
Definition discretization/facecentered/diamond/fvelementgeometry.hh:162
void bind(const Element &element) &
Definition discretization/facecentered/diamond/fvelementgeometry.hh:179
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:125
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:241
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition discretization/facecentered/diamond/fvelementgeometry.hh:231
FaceCenteredDiamondFVElementGeometry 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/facecentered/diamond/fvelementgeometry.hh:189
Element-wise grid geometry (local view).
Definition discretization/facecentered/diamond/fvelementgeometry.hh:38
Classes representing interpolation point data for control-volume finite element schemes.
Helper class to construct SCVs and SCVFs for the diamond scheme.
Defines the index types used for grid and local indices.
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
std::uint_least8_t SmallLocalIndex
Definition indextraits.hh:29