3.1-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 using FVGridGeometry [[deprecated ("Use GridGeometry instead. Will be removed after 3.1!")]] = GridGeometry;
103
106 template<class CellCenterOrFaceFVGridGeometry>
107 StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry& gridGeometry)
108 : gridGeometryPtr_(&gridGeometry.actualGridGeometry()) {}
109
112 : gridGeometryPtr_(&gridGeometry) {}
113
115 const SubControlVolumeFace& scvf(GridIndexType eIdx, LocalIndexType localScvfIdx) const
116 {
117 return scvf(this->gridGeometry().localToGlobalScvfIndex(eIdx, localScvfIdx));
118 }
119
122 const SubControlVolume& scv(GridIndexType scvIdx) const
123 {
124 if (scvIdx == scvIndices_[0])
125 return scvs_[0];
126 else
127 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
128 }
129
132 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
133 {
134 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
135 if (it != scvfIndices_.end())
136 return scvfs_[std::distance(scvfIndices_.begin(), it)];
137 else
138 return neighborScvfs_[findLocalIndex_(scvfIdx, neighborScvfIndices_)];
139 }
140
146 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
147 scvs(const ThisType& g)
148 {
149 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
150 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
151 }
152
158 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
159 scvfs(const ThisType& g)
160 {
161 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
162 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
163 }
164
166 std::size_t numScv() const
167 { return scvs_.size(); }
168
170 std::size_t numScvf() const
171 { return scvfs_.size(); }
172
175 void bind(const Element& element)
176 {
177 bindElement(element);
178
179 neighborScvs_.reserve(element.subEntities(1));
180 neighborScvfIndices_.reserve(element.subEntities(1));
181 neighborScvfs_.reserve(element.subEntities(1));
182
183 std::vector<GridIndexType> handledNeighbors;
184 handledNeighbors.reserve(element.subEntities(1));
185 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
186 {
187 if (intersection.neighbor())
188 {
189 const auto outside = intersection.outside();
190 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
191
192 // make outside geometries only if not done yet (could happen on non-conforming grids)
193 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
194 {
195 makeNeighborGeometries_(outside, outsideIdx);
196 handledNeighbors.push_back(outsideIdx);
197 }
198 }
199 }
200 }
201
203 void bindElement(const Element& element)
204 {
205 clear_();
206 elementPtr_ = &element;
207 scvfs_.reserve(element.subEntities(1));
208 scvfIndices_.reserve(element.subEntities(1));
209 makeElementGeometries_(element);
210 }
211
213 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
215 { return gridGeometry(); }
217 { return *gridGeometryPtr_; }
218
220 bool hasBoundaryScvf() const
221 { return hasBoundaryScvf_; }
222
223private:
224
226 void makeElementGeometries_(const Element& element)
227 {
228 const auto eIdx = gridGeometry().elementMapper().index(element);
229 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
230 scvIndices_[0] = eIdx;
231
232 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
233 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
234
235 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
236
237 int scvfCounter = 0;
238 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
239 {
240 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
241
242 if (intersection.neighbor() || intersection.boundary())
243 {
244 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
245 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
246 scvfs_.emplace_back(intersection,
247 intersection.geometry(),
248 scvFaceIndices[scvfCounter],
249 scvIndices,
250 geometryHelper);
251 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
252 scvfCounter++;
253
254 if (intersection.boundary())
255 hasBoundaryScvf_ = true;
256 }
257 }
258 }
259
261 void makeNeighborGeometries_(const Element& element, const GridIndexType eIdx)
262 {
263 // using ScvfGridIndexStorage = typename SubControlVolumeFace::Traits::GridIndexStorage;
264
265 // create the neighbor scv
266 neighborScvs_.emplace_back(element.geometry(), eIdx);
267 neighborScvIndices_.push_back(eIdx);
268
269 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
270 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
271
272 typename GridGeometry::GeometryHelper geometryHelper(element, gridGeometry().gridView());
273
274 int scvfCounter = 0;
275 for (const auto& intersection : intersections(gridGeometry().gridView(), element))
276 {
277 const auto& scvfNeighborVolVarIndex = neighborVolVarIndices[scvfCounter];
278 geometryHelper.updateLocalFace(gridGeometry().intersectionMapper(), intersection);
279
280 if (intersection.neighbor())
281 {
282 // only create subcontrol faces where the outside element is the bound element
283 if (intersection.outside() == *elementPtr_)
284 {
285 std::vector<GridIndexType> scvIndices{eIdx, scvfNeighborVolVarIndex};
286 neighborScvfs_.emplace_back(intersection,
287 intersection.geometry(),
288 scvFaceIndices[scvfCounter],
289 scvIndices,
290 geometryHelper);
291
292 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
293 }
294 scvfCounter++;
295 }
296 else if (intersection.boundary())
297 scvfCounter++;
298 }
299 }
300
301 const LocalIndexType findLocalIndex_(const GridIndexType idx,
302 const std::vector<GridIndexType>& indices) const
303 {
304 auto it = std::find(indices.begin(), indices.end(), idx);
305 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
306 return std::distance(indices.begin(), it);
307 }
308
310 void clear_()
311 {
312 scvfIndices_.clear();
313 scvfs_.clear();
314
315 neighborScvIndices_.clear();
316 neighborScvfIndices_.clear();
317 neighborScvs_.clear();
318 neighborScvfs_.clear();
319
320 hasBoundaryScvf_ = false;
321 }
322
323 const Element* elementPtr_;
324 const GridGeometry* gridGeometryPtr_;
325
326 // local storage after binding an element
327 std::array<GridIndexType, 1> scvIndices_;
328 std::array<SubControlVolume, 1> scvs_;
329
330 std::vector<GridIndexType> scvfIndices_;
331 std::vector<SubControlVolumeFace> scvfs_;
332
333 std::vector<GridIndexType> neighborScvIndices_;
334 std::vector<SubControlVolume> neighborScvs_;
335
336 std::vector<GridIndexType> neighborScvfIndices_;
337 std::vector<SubControlVolumeFace> neighborScvfs_;
338
339 bool hasBoundaryScvf_ = false;
340};
341
342
343} // end namespace
344
345#endif
Defines the index types used for grid and local indices.
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
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:220
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/staggered/fvelementgeometry.hh:159
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:122
StaggeredFVElementGeometry(const CellCenterOrFaceFVGridGeometry &gridGeometry)
Definition: discretization/staggered/fvelementgeometry.hh:107
const GridGeometry & fvGridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/staggered/fvelementgeometry.hh:214
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:166
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:115
const GridGeometry & gridGeometry() const
Definition: discretization/staggered/fvelementgeometry.hh:216
GridGeometry FVGridGeometry
Definition: discretization/staggered/fvelementgeometry.hh:102
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/staggered/fvelementgeometry.hh:170
StaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/staggered/fvelementgeometry.hh:111
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:147
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/staggered/fvelementgeometry.hh:203
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/staggered/fvelementgeometry.hh:132
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:175
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...