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
25
26namespace Dumux {
27
32template<class GG, bool cachingEnabled>
34
39template<class GG>
40class FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>
41{
42 using ThisType = FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>;
43 using GridView = typename GG::GridView;
44 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
45 using LocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
46 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
47 using GGCache = typename GG::Cache;
48 using GeometryHelper = typename GGCache::GeometryHelper;
50 typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate,
51 LocalIndexType>;
52
53public:
55 using SubControlVolume = typename GG::SubControlVolume;
56 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
57 using Element = typename GridView::template Codim<0>::Entity;
58 using GridGeometry = GG;
59
61 static constexpr std::size_t maxNumElementScvs = 2*GridView::dimension;
62
64 : ggCache_(&ggCache)
65 {}
66
68 const SubControlVolume& scv(LocalIndexType scvIdx) const
69 { return ggCache_->scvs(eIdx_)[scvIdx]; }
70
72 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
73 { return ggCache_->scvf(eIdx_)[scvfIdx]; }
74
80 friend inline auto
82 {
83 using Iter = typename std::vector<SubControlVolume>::const_iterator;
84 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
85 return Dune::IteratorRange<Iter>(s.begin(), s.end());
86 }
87
93 friend inline auto
95 {
96 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
97 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
98 return Dune::IteratorRange<Iter>(s.begin(), s.end());
99 }
100
102 const FeLocalBasis& feLocalBasis() const
103 {
104 return gridGeometry().feCache().get(element().type()).localBasis();
105 }
106
108 std::size_t numLocalDofs() const
109 {
110 return numScv();
111 }
112
114 std::size_t numScv() const
115 {
116 return ggCache_->scvs(eIdx_).size();
117 }
118
120 std::size_t numScvf() const
121 {
122 return ggCache_->scvfs(eIdx_).size();
123 }
124
126 bool hasBoundaryScvf() const
127 { return ggCache_->hasBoundaryScvf(eIdx_); }
128
135 {
136 this->bindElement(element);
137 return std::move(*this);
138 }
139
140 void bind(const Element& element) &
141 {
142 this->bindElement(element);
143 }
144
151 {
152 this->bindElement(element);
153 return std::move(*this);
154 }
155
157 void bindElement(const Element& element) &
158 {
159 element_ = element;
160 elementGeometry_.emplace(element.geometry());
161 eIdx_ = gridGeometry().elementMapper().index(element);
162 }
163
165 bool isBound() const
166 { return static_cast<bool>(element_); }
167
169 const Element& element() const
170 { return *element_; }
171
173 const typename Element::Geometry& elementGeometry() const
174 { return *elementGeometry_; }
175
178 { return ggCache_->gridGeometry(); }
179
181 std::size_t elementIndex() const
182 { return eIdx_; }
183
185 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
186 {
187 assert(isBound());
188 return {
189 SubControlVolume::Traits::geometryType((*elementGeometry_).type()),
190 GeometryHelper(*elementGeometry_).getScvCorners(scv.indexInElement())
191 };
192 }
193
195 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
196 {
197 assert(isBound());
198 if (scvf.boundary())
199 {
200 // use the information that each boundary scvf corresponds to one scv constructed around the same facet
201 const auto localFacetIndex = scvf.insideScvIdx();
202 return {
203 referenceElement(*elementGeometry_).type(localFacetIndex, 1),
204 GeometryHelper(*elementGeometry_).getBoundaryScvfCorners(localFacetIndex)
205 };
206 }
207 else
208 {
209 return {
210 SubControlVolumeFace::Traits::interiorGeometryType((*elementGeometry_).type()),
211 GeometryHelper(*elementGeometry_).getScvfCorners(scvf.index())
212 };
213 }
214 }
215
217 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const SubControlVolume& scv)
218 {
219 const auto type = fvGeometry.element().type();
220 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(scv.localDofIndex());
221
222 return IpData(GeometryHelper::localDofPosition(type, localKey), scv.dofPosition(), scv.localDofIndex());
223 }
224
226 template<class LocalDof>
227 friend inline IpData ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const LocalDof& localDof)
228 {
229 const auto type = fvGeometry.element().type();
230 const auto& localKey = fvGeometry.gridGeometry().feCache().get(type).localCoefficients().localKey(localDof.index());
231 const auto& localPos = GeometryHelper::localDofPosition(type, localKey);
232
233 return IpData(localPos, fvGeometry.elementGeometry().global(localPos), localDof.index());
234 }
235
237 friend inline auto ipData(const FaceCenteredDiamondFVElementGeometry& fvGeometry, const typename Element::Geometry::GlobalCoordinate& globalPos)
238 {
239 // Create ipData that does not automatically calculate the local position but only if it is called
241 [&] (const typename Element::Geometry::GlobalCoordinate& pos) { return fvGeometry.elementGeometry().local(pos); },
242 globalPos
243 };
244 }
245
246private:
247 std::optional<Element> element_;
248 std::optional<typename Element::Geometry> elementGeometry_;
249 GridIndexType eIdx_;
250 const GGCache* ggCache_;
251};
252
253} // end namespace Dumux
254
255#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
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:41
friend auto scvs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:81
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:157
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:177
friend IpData ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const LocalDof &localDof)
Interpolation point data for a localDof.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:227
const Element & element() const
The bound element.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:169
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:55
const Element::Geometry & elementGeometry() const
The bound element geometry.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:173
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:126
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:56
friend auto ipData(const FaceCenteredDiamondFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Interpolation point data for an scv.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:217
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:165
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:120
friend auto scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:94
std::size_t elementIndex() const
The bound element index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:181
GG GridGeometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:58
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:237
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:72
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:68
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:57
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:134
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:114
std::size_t numLocalDofs() const
The total number of element-local dofs.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:108
FaceCenteredDiamondFVElementGeometry(const GGCache &ggCache)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:63
void bind(const Element &element) &
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:140
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:102
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:195
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:185
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:150
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:33
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
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