3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
26#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_PQ1BUBBLE_FV_ELEMENT_GEOMETRY_HH
28
29#include <optional>
30#include <utility>
31
32#include <dune/common/exceptions.hh>
33#include <dune/geometry/type.hh>
34#include <dune/localfunctions/lagrange/pqkfactory.hh>
35
38
40
41namespace Dumux {
42
51template<class GG, bool enableGridGeometryCache>
53
55template<class GG>
57{
58 using GridView = typename GG::GridView;
59 static constexpr int dim = GridView::dimension;
60 static constexpr int dimWorld = GridView::dimensionworld;
61 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
62 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
63 using CoordScalar = typename GridView::ctype;
64 using FeLocalBasis = typename GG::FeCache::FiniteElementType::Traits::LocalBasisType;
65 using GGCache = typename GG::Cache;
67public:
69 using Element = typename GridView::template Codim<0>::Entity;
71 using SubControlVolume = typename GG::SubControlVolume;
73 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
75 using GridGeometry = GG;
77 // ToDo get this from GG
78 static constexpr std::size_t maxNumElementScvs = (1<<dim) + 1;
79
81 PQ1BubbleFVElementGeometry(const GGCache& ggCache)
82 : ggCache_(&ggCache)
83 {}
84
86 const SubControlVolume& scv(LocalIndexType scvIdx) const
87 {
88 return ggCache_->scvs(eIdx_)[scvIdx];
89 }
90
92 const SubControlVolumeFace& scvf(LocalIndexType scvfIdx) const
93 {
94 return ggCache_->scvfs(eIdx_)[scvfIdx];
95 }
96
102 friend inline Dune::IteratorRange<typename std::vector<SubControlVolume>::const_iterator>
104 {
105 using Iter = typename std::vector<SubControlVolume>::const_iterator;
106 const auto& s = fvGeometry.ggCache_->scvs(fvGeometry.eIdx_);
107 return Dune::IteratorRange<Iter>(s.begin(), s.end());
108 }
109
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 return Dune::IteratorRange<Iter>(s.begin(), s.end());
121 }
122
124 const FeLocalBasis& feLocalBasis() const
125 {
126 return gridGeometry().feCache().get(element_->type()).localBasis();
127 }
128
130 std::size_t numScv() const
131 {
132 return ggCache_->scvs(eIdx_).size();
133 }
134
136 std::size_t numScvf() const
137 {
138 return ggCache_->scvfs(eIdx_).size();
139 }
140
147 {
148 this->bindElement(element);
149 return std::move(*this);
150 }
151
155 void bind(const Element& element) &
156 { this->bindElement(element); }
157
164 {
165 this->bindElement(element);
166 return std::move(*this);
167 }
168
172 void bindElement(const Element& element) &
173 {
174 element_ = element;
175 // cache element index
176 eIdx_ = gridGeometry().elementMapper().index(element);
177 }
178
180 bool isBound() const
181 { return static_cast<bool>(element_); }
182
184 const Element& element() const
185 { return *element_; }
186
189 { return ggCache_->gridGeometry(); }
190
192 bool hasBoundaryScvf() const
193 { return ggCache_->hasBoundaryScvf(eIdx_); }
194
196 std::size_t elementIndex() const
197 { return eIdx_; }
198
200 typename SubControlVolume::Traits::Geometry geometry(const SubControlVolume& scv) const
201 {
202 if (scv.isOverlapping())
203 DUNE_THROW(Dune::NotImplemented, "Geometry of overlapping scv");
204
205 assert(isBound());
206 const auto geo = element().geometry();
207 const GeometryHelper helper(geo);
208 return {
209 helper.getScvGeometryType(scv.indexInElement()),
210 helper.getScvCorners(scv.indexInElement())
211 };
212 }
213
215 typename SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace& scvf) const
216 {
217 assert(isBound());
218 const auto geo = element().geometry();
219 if (scvf.boundary())
220 {
221 GeometryHelper helper(geo);
222 const auto localScvfIdx = scvf.index() - helper.numInteriorScvf();
223 const auto [localFacetIndex, isScvfLocalIdx]
224 = ggCache_->scvfBoundaryGeometryKeys(eIdx_)[localScvfIdx];
225 return {
226 Dune::GeometryTypes::cube(dim-1),
227 helper.getBoundaryScvfCorners(localFacetIndex, isScvfLocalIdx)
228 };
229 }
230 else
231 {
232 GeometryHelper helper(geo);
233 return {
234 helper.getInteriorScvfGeometryType(scvf.index()),
235 helper.getScvfCorners(scvf.index())
236 };
237 }
238 }
239
240private:
241 const GGCache* ggCache_;
242 GridIndexType eIdx_;
243
244 std::optional<Element> element_;
245};
246
247} // end namespace Dumux
248
249#endif
Class providing iterators over sub control volumes and sub control volume faces of an element.
Defines the index types used for grid and local indices.
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
unsigned int LocalIndex
Definition: indextraits.hh:40
Base class for the finite volume geometry vector for pq1bubble models This builds up the sub control ...
Definition: discretization/pq1bubble/fvelementgeometry.hh:52
friend Dune::IteratorRange< typename std::vector< SubControlVolume >::const_iterator > scvs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:103
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const PQ1BubbleFVElementGeometry &fvGeometry)
Definition: discretization/pq1bubble/fvelementgeometry.hh:116
typename GridView::template Codim< 0 >::Entity Element
export the element type
Definition: discretization/pq1bubble/fvelementgeometry.hh:69
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/pq1bubble/fvelementgeometry.hh:192
PQ1BubbleFVElementGeometry(const GGCache &ggCache)
Constructor.
Definition: discretization/pq1bubble/fvelementgeometry.hh:81
const SubControlVolume & scv(LocalIndexType scvIdx) const
Get a sub control volume with a local scv index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:86
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/pq1bubble/fvelementgeometry.hh:75
const FeLocalBasis & feLocalBasis() const
Get a local finite element basis.
Definition: discretization/pq1bubble/fvelementgeometry.hh:124
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/pq1bubble/fvelementgeometry.hh:73
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/pq1bubble/fvelementgeometry.hh:71
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:163
void bindElement(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:172
void bind(const Element &element) &
Definition: discretization/pq1bubble/fvelementgeometry.hh:155
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:146
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Geometry of a sub control volume face.
Definition: discretization/pq1bubble/fvelementgeometry.hh:215
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/pq1bubble/fvelementgeometry.hh:130
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Geometry of a sub control volume.
Definition: discretization/pq1bubble/fvelementgeometry.hh:200
const SubControlVolumeFace & scvf(LocalIndexType scvfIdx) const
Get a sub control volume face with a local scvf index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:92
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/pq1bubble/fvelementgeometry.hh:180
std::size_t elementIndex() const
The bound element index.
Definition: discretization/pq1bubble/fvelementgeometry.hh:196
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/pq1bubble/fvelementgeometry.hh:136
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/pq1bubble/fvelementgeometry.hh:188
const Element & element() const
The bound element.
Definition: discretization/pq1bubble/fvelementgeometry.hh:184
A class to create sub control volume and sub control volume face geometries per element.
Definition: discretization/pq1bubble/geometryhelper.hh:178
std::size_t numInteriorScvf() const
number of interior sub control volume faces
Definition: discretization/pq1bubble/geometryhelper.hh:343
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition: discretization/pq1bubble/geometryhelper.hh:300
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:290
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition: discretization/pq1bubble/geometryhelper.hh:254
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition: discretization/pq1bubble/geometryhelper.hh:198
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:234
Helper class constructing the dual grid finite volume geometries for the cvfe discretizazion method.