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