3.5-git
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
182 {
183 this->bind_(element);
184 return std::move(*this);
185 }
186
188 void bind(const Element& element) &
189 { this->bind_(element); }
190
197 {
198 this->bindElement_(element);
199 return std::move(*this);
200 }
201
203 void bindElement(const Element& element) &
204 { this->bindElement_(element); }
205
207 bool isBound() const
208 { return static_cast<bool>(element_); }
209
211 const Element& element() const
212 { return *element_; }
213
216 { return *gridGeometryPtr_; }
217
219 bool hasBoundaryScvf() const
220 { return hasBoundaryScvf_; }
221
222private:
224 void bindElement_(const Element& element)
225 {
226 clear_();
227 element_ = element;
228 scvfs_.reserve(element.subEntities(1));
229 scvfIndices_.reserve(element.subEntities(1));
230 makeElementGeometries_();
231 }
232
235 void bind_(const Element& element)
236 {
237 bindElement_(element);
238
239 neighborScvs_.reserve(element.subEntities(1));
240 neighborScvfIndices_.reserve(element.subEntities(1));
241 neighborScvfs_.reserve(element.subEntities(1));
242
243 std::vector<GridIndexType> handledNeighbors;
244 handledNeighbors.reserve(element.subEntities(1));
245 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
246 {
247 if (intersection.neighbor())
248 {
249 const auto outside = intersection.outside();
250 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
251
252 // make outside geometries only if not done yet (could happen on non-conforming grids)
253 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
254 {
255 makeNeighborGeometries_(outside, outsideIdx);
256 handledNeighbors.push_back(outsideIdx);
257 }
258 }
259 }
260 }
261
263 void makeElementGeometries_()
264 {
265 const auto& element = *element_;
266 const auto eIdx = gridGeometry().elementMapper().index(element);
267 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
268 scvIndices_[0] = eIdx;
269
270 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
271 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
272
273 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
274
275 int scvfCounter = 0;
276 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
277 {
278 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
279
280 if (intersection.neighbor() || intersection.boundary())
281 {
282 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
283 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
284 scvfs_.emplace_back(intersection,
285 intersection.geometry(),
286 scvFaceIndices[scvfCounter],
287 scvIndices,
288 geometryHelper);
289 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
290 scvfCounter++;
291
292 if (intersection.boundary())
293 hasBoundaryScvf_ = true;
294 }
295 }
296 }
297
299 void makeNeighborGeometries_(const Element& element, const GridIndexType eIdx)
300 {
301 // using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
302
303 // create the neighbor scv
304 neighborScvs_.emplace_back(element.geometry(), eIdx);
305 neighborScvIndices_.push_back(eIdx);
306
307 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
308 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
309
310 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
311
312 int scvfCounter = 0;
313 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
314 {
315 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
316 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
317
318 if (intersection.neighbor())
319 {
320 // only create subcontrol faces where the outside element is the bound element
321 if (intersection.outside() == *element_)
322 {
323 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
324 neighborScvfs_.emplace_back(intersection,
325 intersection.geometry(),
326 scvFaceIndices[scvfCounter],
327 scvIndices,
328 geometryHelper);
329
330 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
331 }
332 scvfCounter++;
333 }
334 else if (intersection.boundary())
335 scvfCounter++;
336 }
337 }
338
339 const LocalIndexType findLocalIndex_(const GridIndexType idx,
340 const std::vector<GridIndexType>& indices) const
341 {
342 auto it = std::find(indices.begin(), indices.end(), idx);
343 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
344 return std::distance(indices.begin(), it);
345 }
346
348 void clear_()
349 {
350 scvfIndices_.clear();
351 scvfs_.clear();
352
353 neighborScvIndices_.clear();
354 neighborScvfIndices_.clear();
355 neighborScvs_.clear();
356 neighborScvfs_.clear();
357
358 hasBoundaryScvf_ = false;
359 }
360
361 std::optional<Element> element_;
362 const GridGeometry* gridGeometryPtr_;
363
364 // local storage after binding an element
365 std::array<GridIndexType, 1> scvIndices_;
366 std::array<SubControlVolume, 1> scvs_;
367
368 std::vector<GridIndexType> scvfIndices_;
369 std::vector<SubControlVolumeFace> scvfs_;
370
371 std::vector<GridIndexType> neighborScvIndices_;
372 std::vector<SubControlVolume> neighborScvs_;
373
374 std::vector<GridIndexType> neighborScvfIndices_;
375 std::vector<SubControlVolumeFace> neighborScvfs_;
376
377 bool hasBoundaryScvf_ = false;
378};
379
380
381} // end namespace
382
383#endif
Defines the index types used for grid and local indices.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:292
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:52
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:62
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:69
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:73
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
StaggeredFVElementGeometry 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/staggered/fvelementgeometry.hh:181
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:219
void bindElement(const Element &element) &
bind this local view to a specific element
Definition: discretization/staggered/fvelementgeometry.hh:203
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:215
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:173
StaggeredFVElementGeometry 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/staggered/fvelementgeometry.hh:196
StaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/staggered/fvelementgeometry.hh:114
const Element & element() const
The bound element.
Definition: discretization/staggered/fvelementgeometry.hh:211
void bind(const Element &element) &
bind this local view to a specific element (full stencil)
Definition: discretization/staggered/fvelementgeometry.hh:188
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/staggered/fvelementgeometry.hh:207
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
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
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...