version 3.11-dev
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
18#include <dune/common/reservedvector.hh>
19#include <dune/common/iteratorrange.hh>
20
26
27namespace Dumux {
28
35template<class GG, bool cachingEnabled>
37
42template<class GG>
43class FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>
44{
45 using ThisType = FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>;
46 using GridView = typename GG::GridView;
47 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
48 using LocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
49 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
50 using GGCache = typename GG::Cache;
51 using GeometryHelper = typename GGCache::GeometryHelper;
52
54 typename GridView::template Codim<0>::Entity::Geometry::LocalCoordinate,
55 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate
56 >;
57
58public:
60 using SubControlVolume = typename GG::SubControlVolume;
61 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
62 using Element = typename GridView::template Codim<0>::Entity;
63 using GridGeometry = GG;
65 using ScvQuadratureRule = typename GG::ScvQuadratureRule;
67 using ScvfQuadratureRule = typename GG::ScvfQuadratureRule;
68
70 static constexpr std::size_t maxNumElementScvs = 2*GridView::dimension;
71
73 : ggCache_(&ggCache)
74 {}
75
77 const SubControlVolume& scv(LocalIndexType scvIdx) const
78 { return ggCache_->scvs(eIdx_)[scvIdx]; }
79
81 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
82 { return ggCache_->scvf(eIdx_)[scvfIdx]; }
83
89 friend inline auto
91 {
92 using Iter = typename std::vector<SubControlVolume>::const_iterator;
93 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
94 return Dune::IteratorRange<Iter>(s.begin(), s.end());
95 }
96
102 friend inline auto
104 {
105 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
106 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
107 return Dune::IteratorRange<Iter>(s.begin(), s.end());
108 }
109
111 const FeLocalBasis& feLocalBasis() const
112 {
113 return gridGeometry().feCache().get(element().type()).localBasis();
114 }
115
117 std::size_t numLocalDofs() const
118 {
119 return numScv();
120 }
121
123 std::size_t numScv() const
124 {
125 return ggCache_->scvs(eIdx_).size();
126 }
127
129 std::size_t numScvf() const
130 {
131 return ggCache_->scvfs(eIdx_).size();
132 }
133
135 bool hasBoundaryScvf() const
136 { return ggCache_->hasBoundaryScvf(eIdx_); }
137
144 {
145 this->bindElement(element);
146 return std::move(*this);
147 }
148
149 void bind(const Element& element) &
150 {
151 this->bindElement(element);
152 }
153
160 {
161 this->bindElement(element);
162 return std::move(*this);
163 }
164
166 void bindElement(const Element& element) &
167 {
168 element_ = element;
169 elementGeometry_.emplace(element.geometry());
170 eIdx_ = gridGeometry().elementMapper().index(element);
171 }
172
174 bool isBound() const
175 { return static_cast<bool>(element_); }
176
178 const Element& element() const
179 { return *element_; }
180
182 const typename Element::Geometry& elementGeometry() const
183 { return *elementGeometry_; }
184
187 { return ggCache_->gridGeometry(); }
188
190 std::size_t elementIndex() const
191 { return eIdx_; }
192
194 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
195 {
196 assert(isBound());
197 return {
198 SubControlVolume::Traits::geometryType((*elementGeometry_).type()),
199 GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement())
200 };
201 }
202
204 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
205 {
206 assert(isBound());
207 if (scvf.boundary())
208 {
209 // use the information that each boundary scvf corresponds to one scv constructed around the same facet
210 const auto localFacetIndex = scvf.insideScvIdx();
211 return {
212 referenceElement(*elementGeometry_).type(localFacetIndex, 1),
213 GeometryHelper(*elementGeometry_).getBoundaryScvfCorners(localFacetIndex)
214 };
215 }
216 else
217 {
218 return {
219 SubControlVolumeFace::Traits::interiorGeometryType((*elementGeometry_).type()),
220 GeometryHelper(*elementGeometry_).getScvfCorners(scvf.index())
221 };
222 }
223 }
224
226 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const SubControlVolume& scv)
227 {
228 const auto type = fvGeometry.element().type();
229 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
230
231 return CVFE::LocalDofInterpolationPointData{ GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex() };
232 }
233
235 template<class LocalDof>
236 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const LocalDof& localDof)
237 {
238 const auto type = fvGeometry.element().type();
239 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
240 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
241
242 return CVFE::LocalDofInterpolationPointData{ localPos, fvGeometry.elementGeometry().global(localPos), localDof.index() };
243 }
244
246 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
247 {
248 // Create ipData that does not automatically calculate the local position but only if it is called
250 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
251 globalPos
252 };
253 }
254
256 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
257 {
259 { scvf.unitOuterNormal(), scvf.index(), fvGeometry.elementGeometry().local(scvf.ipGlobal()), scvf.ipGlobal() };
260 }
261
262private:
263 std::optional<Element> element_;
264 std::optional<typename Element::Geometry> elementGeometry_;
265 GridIndexType eIdx_;
266 const GGCache* ggCache_;
267};
268
269} // end namespace Dumux
270
271#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
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:44
friend auto scvs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:90
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:166
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:186
const Element & element() const
The bound element.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:178
typename GG::ScvQuadratureRule ScvQuadratureRule
the quadrature rule type for scvs
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:65
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:60
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:236
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:182
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Interpolation point data for scvf.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:256
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:135
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:61
typename GG::ScvfQuadratureRule ScvfQuadratureRule
the quadrature rule type for scvfs
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:67
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:226
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:174
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:129
friend auto scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:103
std::size_t elementIndex() const
The bound element index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:190
GG GridGeometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:63
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:246
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:81
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:77
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:62
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:143
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:123
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:117
FaceCenteredDiamondFVElementGeometry(const GGCache &ggCache)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:72
void bind(const Element &element) &
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:149
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:111
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:204
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:194
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:159
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:36
Classes representing interpolation point data for control-volume finite element schemes.
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