3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/staggered/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_STAGGERED_FV_ELEMENT_GEOMETRY_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_FV_ELEMENT_GEOMETRY_HH
26
27#include <optional>
30
31namespace Dumux {
32
42template<class GG, bool enableGridGeometryCache>
44
53template<class GG>
54class StaggeredFVElementGeometry<GG, true> : public CCTpfaFVElementGeometry<GG, true>
55{
57 using GridView = typename GG::GridView;
58 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
59 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
60public:
62 using Element = typename GridView::template Codim<0>::Entity;
64 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
65
66 using ParentType::ParentType;
67
70 template<class CellCenterOrFaceFVGridGeometry>
71 StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry& gridGeometry)
72 : ParentType(gridGeometry.actualGridGeometry()) {}
73
75 using ParentType::scvf;
76 const SubControlVolumeFace& scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
77 {
78 return this->gridGeometry().scvf(eIdx, localScvfIdx);
79 }
80};
81
90template<class GG>
92{
94 using GridView = typename GG::GridView;
95 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
96 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
97public:
99 using Element = typename GridView::template Codim<0>::Entity;
101 using SubControlVolume = typename GG::SubControlVolume;
103 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
105 using GridGeometry = GG;
106
109 template<class CellCenterOrFaceFVGridGeometry>
110 StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry& gridGeometry)
111 : gridGeometryPtr_(&gridGeometry.actualGridGeometry()) {}
112
115 : gridGeometryPtr_(&gridGeometry) {}
116
118 const SubControlVolumeFace& scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
119 {
120 return scvf(this->gridGeometry().localToGlobalScvfIndex(eIdx, localScvfIdx));
121 }
122
125 const SubControlVolume& scv(GridIndexType scvIdx) const
126 {
127 if (scvIdx == scvIndices_[0])
128 return scvs_[0];
129 else
130 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
131 }
132
135 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
136 {
137 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
138 if (it != scvfIndices_.end())
139 return scvfs_[std::distance(scvfIndices_.begin(), it)];
140 else
141 return neighborScvfs_[findLocalIndex_(scvfIdx, neighborScvfIndices_)];
142 }
143
149 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
150 scvs(const ThisType& g)
151 {
152 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
153 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
154 }
155
161 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
162 scvfs(const ThisType& g)
163 {
164 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
165 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
166 }
167
169 std::size_t numScv() const
170 { return scvs_.size(); }
171
173 std::size_t numScvf() const
174 { return scvfs_.size(); }
175
178 void bind(const Element& element)
179 {
180 bindElement(element);
181
182 neighborScvs_.reserve(element.subEntities(1));
183 neighborScvfIndices_.reserve(element.subEntities(1));
184 neighborScvfs_.reserve(element.subEntities(1));
185
186 std::vector<GridIndexType> handledNeighbors;
187 handledNeighbors.reserve(element.subEntities(1));
188 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
189 {
190 if (intersection.neighbor())
191 {
192 const auto outside = intersection.outside();
193 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
194
195 // make outside geometries only if not done yet (could happen on non-conforming grids)
196 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
197 {
198 makeNeighborGeometries_(outside, outsideIdx);
199 handledNeighbors.push_back(outsideIdx);
200 }
201 }
202 }
203 }
204
206 void bindElement(const Element& element)
207 {
208 clear_();
209 element_ = element;
210 scvfs_.reserve(element.subEntities(1));
211 scvfIndices_.reserve(element.subEntities(1));
212 makeElementGeometries_();
213 }
214
216 bool isBound() const
217 { return static_cast<bool>(element_); }
218
220 const Element& element() const
221 { return *element_; }
222
225 { return *gridGeometryPtr_; }
226
228 bool hasBoundaryScvf() const
229 { return hasBoundaryScvf_; }
230
231private:
232
234 void makeElementGeometries_()
235 {
236 const auto& element = *element_;
237 const auto eIdx = gridGeometry().elementMapper().index(element);
238 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
239 scvIndices_[0] = eIdx;
240
241 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
242 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
243
244 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
245
246 int scvfCounter = 0;
247 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
248 {
249 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
250
251 if (intersection.neighbor() || intersection.boundary())
252 {
253 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
254 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
255 scvfs_.emplace_back(intersection,
256 intersection.geometry(),
257 scvFaceIndices[scvfCounter],
258 scvIndices,
259 geometryHelper);
260 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
261 scvfCounter++;
262
263 if (intersection.boundary())
264 hasBoundaryScvf_ = true;
265 }
266 }
267 }
268
270 void makeNeighborGeometries_(const Element& element, const GridIndexType eIdx)
271 {
272 // using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
273
274 // create the neighbor scv
275 neighborScvs_.emplace_back(element.geometry(), eIdx);
276 neighborScvIndices_.push_back(eIdx);
277
278 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
279 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
280
281 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
282
283 int scvfCounter = 0;
284 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
285 {
286 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
287 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
288
289 if (intersection.neighbor())
290 {
291 // only create subcontrol faces where the outside element is the bound element
292 if (intersection.outside() == *element_)
293 {
294 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
295 neighborScvfs_.emplace_back(intersection,
296 intersection.geometry(),
297 scvFaceIndices[scvfCounter],
298 scvIndices,
299 geometryHelper);
300
301 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
302 }
303 scvfCounter++;
304 }
305 else if (intersection.boundary())
306 scvfCounter++;
307 }
308 }
309
310 const LocalIndexType findLocalIndex_(const GridIndexType idx,
311 const std::vector<GridIndexType>& indices) const
312 {
313 auto it = std::find(indices.begin(), indices.end(), idx);
314 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
315 return std::distance(indices.begin(), it);
316 }
317
319 void clear_()
320 {
321 scvfIndices_.clear();
322 scvfs_.clear();
323
324 neighborScvIndices_.clear();
325 neighborScvfIndices_.clear();
326 neighborScvs_.clear();
327 neighborScvfs_.clear();
328
329 hasBoundaryScvf_ = false;
330 }
331
332 std::optional<Element> element_;
333 const GridGeometry* gridGeometryPtr_;
334
335 // local storage after binding an element
336 std::array<GridIndexType, 1> scvIndices_;
337 std::array<SubControlVolume, 1> scvs_;
338
339 std::vector<GridIndexType> scvfIndices_;
340 std::vector<SubControlVolumeFace> scvfs_;
341
342 std::vector<GridIndexType> neighborScvIndices_;
343 std::vector<SubControlVolume> neighborScvs_;
344
345 std::vector<GridIndexType> neighborScvfIndices_;
346 std::vector<SubControlVolumeFace> neighborScvfs_;
347
348 bool hasBoundaryScvf_ = false;
349};
350
351
352} // end namespace
353
354#endif
Defines the index types used for grid and local indices.
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:138
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:51
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:61
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:68
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:72
Stencil-local finite volume geometry (scvs and scvfs) for staggered models This builds up the sub con...
Definition: discretization/staggered/fvelementgeometry.hh:43
const SubControlVolumeFace & scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:76
StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry &gridGeometry)
Definition: discretization/staggered/fvelementgeometry.hh:71
Base class for the finite volume geometry vector for staggered models This locally builds up the sub ...
Definition: discretization/staggered/fvelementgeometry.hh:92
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/staggered/fvelementgeometry.hh:99
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/staggered/fvelementgeometry.hh:228
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/staggered/fvelementgeometry.hh:162
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/staggered/fvelementgeometry.hh:101
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:125
StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry &gridGeometry)
Definition: discretization/staggered/fvelementgeometry.hh:110
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:169
const SubControlVolumeFace & scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
Get a sub control volume face with an element index and a local scvf index.
Definition: discretization/staggered/fvelementgeometry.hh:118
const GridGeometry & gridGeometry() const
The grid finite volume geometry we are a restriction of.
Definition: discretization/staggered/fvelementgeometry.hh:224
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:173
StaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/staggered/fvelementgeometry.hh:114
const Element & element() const
The bound element.
Definition: discretization/staggered/fvelementgeometry.hh:220
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/staggered/fvelementgeometry.hh:216
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/staggered/fvelementgeometry.hh:103
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/staggered/fvelementgeometry.hh:150
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/staggered/fvelementgeometry.hh:206
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:135
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/staggered/fvelementgeometry.hh:105
void bind(const Element &element)
Definition: discretization/staggered/fvelementgeometry.hh:178
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...