3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_FV_ELEMENT_GEOMETRY_HH
26
27#include <type_traits>
28#include <optional>
29
30#include <dune/common/reservedvector.hh>
31#include <dune/common/iteratorrange.hh>
32
36
37namespace Dumux {
38
43template<class GG, bool cachingEnabled>
45
50template<class GG>
51class FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>
52{
53 using ThisType = FaceCenteredDiamondFVElementGeometry<GG, /*cachingEnabled*/true>;
54 using GridView = typename GG::GridView;
55 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
56 using LocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
57 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
58 using GGCache = typename GG::Cache;
60
61public:
63 using SubControlVolume = typename GG::SubControlVolume;
64 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
65 using Element = typename GridView::template Codim<0>::Entity;
66 using GridGeometry = GG;
67
69 static constexpr std::size_t maxNumElementScvs = 2*GridView::dimension;
70
72 : ggCache_(&ggCache)
73 {}
74
76 const SubControlVolume& scv(LocalIndexType scvIdx) const
77 { return ggCache_->scvs(eIdx_)[scvIdx]; }
78
80 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
81 { return ggCache_->scvf(eIdx_)[scvfIdx]; }
82
88 friend inline auto
90 {
91 using Iter = typename std::vector<SubControlVolume>::const_iterator;
92 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
93 return Dune::IteratorRange<Iter>(s.begin(), s.end());
94 }
95
101 friend inline auto
103 {
104 using Iter = typename std::vector<SubControlVolumeFace>::const_iterator;
105 const auto& s = fvGeometry.ggCache_->scvfs(fvGeometry.eIdx_);
106 return Dune::IteratorRange<Iter>(s.begin(), s.end());
107 }
108
110 const FeLocalBasis& feLocalBasis() const
111 {
112 return gridGeometry().feCache().get(element().type()).localBasis();
113 }
114
116 std::size_t numScv() const
117 {
118 return ggCache_->scvs(eIdx_).size();
119 }
120
122 std::size_t numScvf() const
123 {
124 return ggCache_->scvfs(eIdx_).size();
125 }
126
128 bool hasBoundaryScvf() const
129 { return ggCache_->hasBoundaryScvf(eIdx_); }
130
137 {
138 this->bindElement(element);
139 return std::move(*this);
140 }
141
142 void bind(const Element& element) &
143 {
144 this->bindElement(element);
145 }
146
153 {
154 this->bindElement(element);
155 return std::move(*this);
156 }
157
159 void bindElement(const Element& element) &
160 {
161 element_ = element;
162 eIdx_ = gridGeometry().elementMapper().index(element);
163 }
164
166 bool isBound() const
167 { return static_cast<bool>(element_); }
168
170 const Element& element() const
171 { return *element_; }
172
175 { return ggCache_->gridGeometry(); }
176
178 std::size_t elementIndex() const
179 { return eIdx_; }
180
182 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
183 {
184 assert(isBound());
185 const auto geo = element().geometry();
186 return {
187 SubControlVolume::Traits::geometryType(geo.type()),
188 GeometryHelper(geo).getScvCorners(scv.indexInElement())
189 };
190 }
191
193 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
194 {
195 assert(isBound());
196 const auto geo = element().geometry();
197 if (scvf.boundary())
198 {
199 // use the information that each boundary scvf corresponds to one scv constructed around the same facet
200 const auto localFacetIndex = scvf.insideScvIdx();
201 return {
202 referenceElement(geo).type(localFacetIndex, 1),
203 GeometryHelper(geo).getBoundaryScvfCorners(localFacetIndex)
204 };
205 }
206 else
207 {
208 return {
209 SubControlVolumeFace::Traits::interiorGeometryType(geo.type()),
210 GeometryHelper(geo).getScvfCorners(scvf.index())
211 };
212 }
213 }
214
215private:
216 std::optional<Element> element_;
217 GridIndexType eIdx_;
218 const GGCache* ggCache_;
219};
220
221} // end namespace Dumux
222
223#endif
Defines the index types used for grid and local indices.
Class providing iterators over sub control volumes and sub control volume faces of an element.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:44
Element-wise grid geometry (local view)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:52
friend auto scvs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:89
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:159
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:174
const Element & element() const
The bound element.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:170
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:63
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:128
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:64
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:166
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:122
friend auto scvfs(const FaceCenteredDiamondFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:102
std::size_t elementIndex() const
The bound element index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:178
GG GridGeometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:66
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:80
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:76
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:65
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:136
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:116
FaceCenteredDiamondFVElementGeometry(const GGCache &ggCache)
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:71
void bind(const Element &element) &
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:142
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:110
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:193
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/facecentered/diamond/fvelementgeometry.hh:182
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:152
Helper class to construct SCVs and SCVFs for the diamond scheme.
Definition: discretization/facecentered/diamond/geometryhelper.hh:267
ScvfCornerStorage getScvfCorners(unsigned int localEdgeIndex) const
Create a corner storage with the scvf corners for a given edge (codim-2) index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:315
ScvCornerStorage getScvCorners(unsigned int localFacetIndex) const
Create a corner storage with the scv corners for a given face (codim-1) index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:287
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex) const
Create the sub control volume face geometries on the boundary for a given face index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:343