3.3.0
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
29
30namespace Dumux {
31
41template<class GG, bool enableGridGeometryCache>
43
52template<class GG>
53class StaggeredFVElementGeometry<GG, true> : public CCTpfaFVElementGeometry<GG, true>
54{
56 using GridView = typename GG::GridView;
57 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
58 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
59public:
61 using Element = typename GridView::template Codim<0>::Entity;
63 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
64
65 using ParentType::ParentType;
66
69 template<class CellCenterOrFaceFVGridGeometry>
70 StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry& gridGeometry)
71 : ParentType(gridGeometry.actualGridGeometry()) {}
72
74 using ParentType::scvf;
75 const SubControlVolumeFace& scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
76 {
77 return this->gridGeometry().scvf(eIdx, localScvfIdx);
78 }
79};
80
89template<class GG>
91{
93 using GridView = typename GG::GridView;
94 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
95 using LocalIndexType = typename IndexTraits<GridView>::LocalIndex;
96public:
98 using Element = typename GridView::template Codim<0>::Entity;
100 using SubControlVolume = typename GG::SubControlVolume;
102 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
104 using GridGeometry = GG;
105
108 template<class CellCenterOrFaceFVGridGeometry>
109 StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry& gridGeometry)
110 : gridGeometryPtr_(&gridGeometry.actualGridGeometry()) {}
111
114 : gridGeometryPtr_(&gridGeometry) {}
115
117 const SubControlVolumeFace& scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
118 {
119 return scvf(this->gridGeometry().localToGlobalScvfIndex(eIdx, localScvfIdx));
120 }
121
124 const SubControlVolume& scv(GridIndexType scvIdx) const
125 {
126 if (scvIdx == scvIndices_[0])
127 return scvs_[0];
128 else
129 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
130 }
131
134 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
135 {
136 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
137 if (it != scvfIndices_.end())
138 return scvfs_[std::distance(scvfIndices_.begin(), it)];
139 else
140 return neighborScvfs_[findLocalIndex_(scvfIdx, neighborScvfIndices_)];
141 }
142
148 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
149 scvs(const ThisType& g)
150 {
151 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
152 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
153 }
154
160 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
161 scvfs(const ThisType& g)
162 {
163 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
164 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
165 }
166
168 std::size_t numScv() const
169 { return scvs_.size(); }
170
172 std::size_t numScvf() const
173 { return scvfs_.size(); }
174
177 void bind(const Element& element)
178 {
179 bindElement(element);
180
181 neighborScvs_.reserve(element.subEntities(1));
182 neighborScvfIndices_.reserve(element.subEntities(1));
183 neighborScvfs_.reserve(element.subEntities(1));
184
185 std::vector<GridIndexType> handledNeighbors;
186 handledNeighbors.reserve(element.subEntities(1));
187 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
188 {
189 if (intersection.neighbor())
190 {
191 const auto outside = intersection.outside();
192 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
193
194 // make outside geometries only if not done yet (could happen on non-conforming grids)
195 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
196 {
197 makeNeighborGeometries_(outside, outsideIdx);
198 handledNeighbors.push_back(outsideIdx);
199 }
200 }
201 }
202 }
203
205 void bindElement(const Element& element)
206 {
207 clear_();
208 elementPtr_ = &element;
209 scvfs_.reserve(element.subEntities(1));
210 scvfIndices_.reserve(element.subEntities(1));
211 makeElementGeometries_(element);
212 }
213
216 { return *gridGeometryPtr_; }
217
219 bool hasBoundaryScvf() const
220 { return hasBoundaryScvf_; }
221
222private:
223
225 void makeElementGeometries_(const Element& element)
226 {
227 const auto eIdx = gridGeometry().elementMapper().index(element);
228 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
229 scvIndices_[0] = eIdx;
230
231 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
232 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
233
234 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
235
236 int scvfCounter = 0;
237 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
238 {
239 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
240
241 if (intersection.neighbor() || intersection.boundary())
242 {
243 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
244 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
245 scvfs_.emplace_back(intersection,
246 intersection.geometry(),
247 scvFaceIndices[scvfCounter],
248 scvIndices,
249 geometryHelper);
250 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
251 scvfCounter++;
252
253 if (intersection.boundary())
254 hasBoundaryScvf_ = true;
255 }
256 }
257 }
258
260 void makeNeighborGeometries_(const Element& element, const GridIndexType eIdx)
261 {
262 // using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
263
264 // create the neighbor scv
265 neighborScvs_.emplace_back(element.geometry(), eIdx);
266 neighborScvIndices_.push_back(eIdx);
267
268 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
269 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
270
271 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
272
273 int scvfCounter = 0;
274 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
275 {
276 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
277 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
278
279 if (intersection.neighbor())
280 {
281 // only create subcontrol faces where the outside element is the bound element
282 if (intersection.outside() == *elementPtr_)
283 {
284 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
285 neighborScvfs_.emplace_back(intersection,
286 intersection.geometry(),
287 scvFaceIndices[scvfCounter],
288 scvIndices,
289 geometryHelper);
290
291 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
292 }
293 scvfCounter++;
294 }
295 else if (intersection.boundary())
296 scvfCounter++;
297 }
298 }
299
300 const LocalIndexType findLocalIndex_(const GridIndexType idx,
301 const std::vector<GridIndexType>& indices) const
302 {
303 auto it = std::find(indices.begin(), indices.end(), idx);
304 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
305 return std::distance(indices.begin(), it);
306 }
307
309 void clear_()
310 {
311 scvfIndices_.clear();
312 scvfs_.clear();
313
314 neighborScvIndices_.clear();
315 neighborScvfIndices_.clear();
316 neighborScvs_.clear();
317 neighborScvfs_.clear();
318
319 hasBoundaryScvf_ = false;
320 }
321
322 const Element* elementPtr_;
323 const GridGeometry* gridGeometryPtr_;
324
325 // local storage after binding an element
326 std::array<GridIndexType, 1> scvIndices_;
327 std::array<SubControlVolume, 1> scvs_;
328
329 std::vector<GridIndexType> scvfIndices_;
330 std::vector<SubControlVolumeFace> scvfs_;
331
332 std::vector<GridIndexType> neighborScvIndices_;
333 std::vector<SubControlVolume> neighborScvs_;
334
335 std::vector<GridIndexType> neighborScvfIndices_;
336 std::vector<SubControlVolumeFace> neighborScvfs_;
337
338 bool hasBoundaryScvf_ = false;
339};
340
341
342} // end namespace
343
344#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: geometry/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:50
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:60
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:67
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:71
Stencil-local finite volume geometry (scvs and scvfs) for staggered models This builds up the sub con...
Definition: discretization/staggered/fvelementgeometry.hh:42
const SubControlVolumeFace & scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:75
StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry &gridGeometry)
Definition: discretization/staggered/fvelementgeometry.hh:70
Base class for the finite volume geometry vector for staggered models This locally builds up the sub ...
Definition: discretization/staggered/fvelementgeometry.hh:91
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/staggered/fvelementgeometry.hh:98
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/staggered/fvelementgeometry.hh:219
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/staggered/fvelementgeometry.hh:161
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/staggered/fvelementgeometry.hh:100
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:124
StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry &gridGeometry)
Definition: discretization/staggered/fvelementgeometry.hh:109
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:168
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:117
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:172
StaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/staggered/fvelementgeometry.hh:113
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/staggered/fvelementgeometry.hh:102
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/staggered/fvelementgeometry.hh:149
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/staggered/fvelementgeometry.hh:205
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:134
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/staggered/fvelementgeometry.hh:104
void bind(const Element &element)
Definition: discretization/staggered/fvelementgeometry.hh:177
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...