3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/cellcentered/mpfa/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 *****************************************************************************/
26#ifndef DUMUX_DISCRETIZATION_CCMPFA_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CCMPFA_FV_ELEMENT_GEOMETRY_HH
28
29#include <optional>
30#include <dune/common/exceptions.hh>
31#include <dune/common/iteratorrange.hh>
32
36
37namespace Dumux {
38
48template<class GG, bool enableGridGeometryCache>
50
57template<class GG>
59{
61 using GridView = typename GG::GridView;
62 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
63
64 static constexpr int dim = GridView::dimension;
65
66public:
68 using Element = typename GridView::template Codim<0>::Entity;
70 using SubControlVolume = typename GG::SubControlVolume;
72 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
74 using GridGeometry = GG;
76 static constexpr std::size_t maxNumElementScvs = 1;
78 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
79
82 : gridGeometryPtr_(&gridGeometry) {}
83
85 const SubControlVolume& scv(GridIndexType scvIdx) const
86 {
87 return gridGeometry().scv(scvIdx);
88 }
89
91 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
92 {
93 return gridGeometry().scvf(scvfIdx);
94 }
95
98 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
99 {
100 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
101 }
102
108 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
109 scvs(const CCMpfaFVElementGeometry& fvGeometry)
110 {
112 return Dune::IteratorRange<ScvIterator>(ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
113 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
114 }
115
121 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
123 {
124 const auto& g = fvGeometry.gridGeometry();
125 const auto scvIdx = fvGeometry.scvIndices_[0];
127 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
128 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
129 }
130
132 std::size_t numScv() const
133 {
134 return scvIndices_.size();
135 }
136
138 std::size_t numScvf() const
139 {
140 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
141 }
142
144 void bind(const Element& element)
145 {
146 this->bindElement(element);
147 }
148
150 void bindElement(const Element& element)
151 {
152 element_ = element;
153 scvIndices_[0] = gridGeometry().elementMapper().index(element);
154 }
155
157 bool isBound() const
158 { return static_cast<bool>(element_); }
159
161 const Element& element() const
162 { return *element_; }
163
166 { return *gridGeometryPtr_; }
167
169 bool hasBoundaryScvf() const
170 { return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
171
172private:
173
174 std::optional<Element> element_;
175 std::array<GridIndexType, 1> scvIndices_;
176 const GridGeometry* gridGeometryPtr_;
177};
178
184template<class GG>
186{
188 using GridView = typename GG::GridView;
189 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
190 using MpfaHelper = typename GG::MpfaHelper;
191
192 static const int dim = GridView::dimension;
193 static const int dimWorld = GridView::dimensionworld;
194 using CoordScalar = typename GridView::ctype;
195
196public:
198 using Element = typename GridView::template Codim<0>::Entity;
200 using SubControlVolume = typename GG::SubControlVolume;
202 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
204 using GridGeometry = GG;
206 static constexpr std::size_t maxNumElementScvs = 1;
208 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
209
212 : gridGeometryPtr_(&gridGeometry) {}
213
216 const SubControlVolume& scv(GridIndexType scvIdx) const
217 {
218 if (scvIdx == scvIndices_[0])
219 return scvs_[0];
220 else
221 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
222 }
223
226 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
227 {
228 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
229 if (it != scvfIndices_.end())
230 return scvfs_[std::distance(scvfIndices_.begin(), it)];
231 else
232 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
233 }
234
237 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
238 {
239 return scvf( gridGeometry().flipScvfIdx(scvfIdx, outsideScvIdx) );
240 }
241
247 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
248 scvs(const ThisType& g)
249 {
250 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
251 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
252 }
253
259 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
260 scvfs(const ThisType& g)
261 {
262 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
263 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
264 }
265
267 std::size_t numScv() const
268 { return scvs_.size(); }
269
271 std::size_t numScvf() const
272 { return scvfs_.size(); }
273
276 void bind(const Element& element)
277 {
278 // make inside geometries
279 bindElement(element);
280
281 // get some references for convenience
282 const auto globalI = gridGeometry().elementMapper().index(element);
283 const auto& assemblyMapI = gridGeometry().connectivityMap()[globalI];
284
285 // reserve memory
286 const auto numNeighbors = assemblyMapI.size();
287 const auto numNeighborScvfs = numNeighborScvfs_(assemblyMapI);
288 neighborScvs_.reserve(numNeighbors);
289 neighborScvIndices_.reserve(numNeighbors);
290 neighborScvfIndices_.reserve(numNeighborScvfs);
291 neighborScvfs_.reserve(numNeighborScvfs);
292
293 // make neighbor geometries
294 // use the assembly map to determine which faces are necessary
295 for (const auto& dataJ : assemblyMapI)
296 makeNeighborGeometries(gridGeometry().element(dataJ.globalJ),
297 dataJ.globalJ,
298 dataJ.scvfsJ,
299 dataJ.additionalScvfs);
300
301 // //! TODO Check if user added additional DOF dependencies, i.e. the residual of DOF globalI depends
302 // //! on additional DOFs not included in the discretization schemes' occupation pattern
303 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
304 // if (!additionalDofDependencies.empty())
305 // {
306 // const auto newNumNeighbors = neighborScvs_.size() + additionalDofDependencies.size();
307 // neighborScvs_.reserve(newNumNeighbors);
308 // neighborScvIndices_.reserve(newNumNeighbors);
309 // for (auto globalJ : additionalDofDependencies)
310 // {
311 // neighborScvs_.emplace_back(gridGeometry().element(globalJ).geometry(), globalJ);
312 // neighborScvIndices_.emplace_back(globalJ);
313 // }
314 // }
315 }
316
318 void bindElement(const Element& element)
319 {
320 clear();
321 element_ = element;
322 makeElementGeometries(element);
323 }
324
326 bool isBound() const
327 { return static_cast<bool>(element_); }
328
330 const Element& element() const
331 { return *element_; }
332
335 { return *gridGeometryPtr_; }
336
338 bool hasBoundaryScvf() const
339 { return hasBoundaryScvf_; }
340
341private:
342
344 template<class DataJContainer>
345 std::size_t numNeighborScvfs_(const DataJContainer& dataJContainer)
346 {
347 std::size_t numNeighborScvfs = 0;
348 for (const auto& dataJ : dataJContainer)
349 numNeighborScvfs += dataJ.scvfsJ.size() + dataJ.additionalScvfs.size();
350 return numNeighborScvfs;
351 }
352
354 void makeElementGeometries(const Element& element)
355 {
356 // make the scv
357 const auto eIdx = gridGeometry().elementMapper().index(element);
358 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
359 scvIndices_[0] = eIdx;
360
361 // get data on the scv faces
362 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
363 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
364
365 // the quadrature point parameterizaion to be used on scvfs
366 static const auto q = getParam<CoordScalar>("MPFA.Q");
367
368 // reserve memory for the scv faces
369 const auto numLocalScvf = scvFaceIndices.size();
370 scvfIndices_.reserve(numLocalScvf);
371 scvfs_.reserve(numLocalScvf);
372
373 // for network grids we only want to do one scvf per half facet
374 // this approach assumes conforming grids at branching facets
375 std::vector<bool> finishedFacets;
376 if (dim < dimWorld)
377 finishedFacets.resize(element.subEntities(1), false);
378
379 int scvfCounter = 0;
380 for (const auto& is : intersections(gridGeometry().gridView(), element))
381 {
382 // if we are dealing with a lower dimensional network
383 // only make a new scvf if we haven't handled it yet
384 if (dim < dimWorld)
385 {
386 const auto indexInInside = is.indexInInside();
387 if (finishedFacets[indexInInside])
388 continue;
389 else
390 finishedFacets[indexInInside] = true;
391 }
392
393 // if outside level > inside level, use the outside element in the following
394 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
395 const auto& e = useNeighbor ? is.outside() : element;
396 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
397 const auto eg = e.geometry();
398 const auto refElement = referenceElement(eg);
399
400 // Set up a container with all relevant positions for scvf corner computation
401 const auto numCorners = is.geometry().corners();
402 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
403 refElement,
404 indexInElement,
405 numCorners);
406
407 // make the scv faces belonging to each corner of the intersection
408 for (int c = 0; c < numCorners; ++c)
409 {
410 // get the global vertex index the scv face is connected to
411 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
412 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
413
414 // do not build scvfs connected to a processor boundary
415 if (gridGeometry().isGhostVertex(vIdxGlobal))
416 continue;
417
418 hasBoundaryScvf_ = (hasBoundaryScvf_ || is.boundary());
419
420 scvfs_.emplace_back(MpfaHelper(),
421 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
422 is,
423 vIdxGlobal,
424 vIdxLocal,
425 scvFaceIndices[scvfCounter],
426 eIdx,
427 neighborVolVarIndices[scvfCounter],
428 q,
429 is.boundary());
430
431 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
432 scvfCounter++;
433 }
434 }
435 }
436
438 template<typename IndexVector>
439 void makeNeighborGeometries(const Element& element,
440 GridIndexType eIdxGlobal,
441 const IndexVector& scvfIndices,
442 const IndexVector& additionalScvfs)
443 {
444 // create the neighbor scv if it doesn't exist yet
445 neighborScvs_.emplace_back(element.geometry(), eIdxGlobal);
446 neighborScvIndices_.push_back(eIdxGlobal);
447
448 // get data on the scv faces
449 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdxGlobal);
450 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdxGlobal);
451
452 // the quadrature point parameterizaion to be used on scvfs
453 static const auto q = getParam<CoordScalar>("MPFA.Q");
454
455 // for network grids we only want to do one scvf per half facet
456 // this approach assumes conforming grids at branching facets
457 std::vector<bool> finishedFacets;
458 if (dim < dimWorld)
459 finishedFacets.resize(element.subEntities(1), false);
460
461 int scvfCounter = 0;
462 for (const auto& is : intersections(gridGeometry().gridView(), element))
463 {
464 // if we are dealing with a lower dimensional network
465 // only make a new scvf if we haven't handled it yet
466 if (dim < dimWorld)
467 {
468 auto indexInInside = is.indexInInside();
469 if(finishedFacets[indexInInside])
470 continue;
471 else
472 finishedFacets[indexInInside] = true;
473 }
474
475 // if outside level > inside level, use the outside element in the following
476 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
477 const auto& e = useNeighbor ? is.outside() : element;
478 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
479 const auto eg = e.geometry();
480 const auto refElement = referenceElement(eg);
481
482 // Set up a container with all relevant positions for scvf corner computation
483 const auto numCorners = is.geometry().corners();
484 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
485 refElement,
486 indexInElement,
487 numCorners);
488
489 // make the scv faces belonging to each corner of the intersection
490 for (int c = 0; c < numCorners; ++c)
491 {
492 // get the global vertex index the scv face is connected to
493 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
494 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
495
496 // do not build scvfs connected to a processor boundary
497 if (gridGeometry().isGhostVertex(vIdxGlobal))
498 continue;
499
500 // only build the scvf if it is in the list of necessary indices
501 if (!MpfaHelper::vectorContainsValue(scvfIndices, scvFaceIndices[scvfCounter])
502 && !MpfaHelper::vectorContainsValue(additionalScvfs, scvFaceIndices[scvfCounter]))
503 {
504 // increment counter either way
505 scvfCounter++;
506 continue;
507 }
508
509 // build scvf
510 neighborScvfs_.emplace_back(MpfaHelper(),
511 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
512 is,
513 vIdxGlobal,
514 vIdxLocal,
515 scvFaceIndices[scvfCounter],
516 eIdxGlobal,
517 neighborVolVarIndices[scvfCounter],
518 q,
519 is.boundary());
520
521 neighborScvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
522
523 // increment counter
524 scvfCounter++;
525 }
526 }
527 }
528
530 unsigned int findLocalIndex(const GridIndexType idx,
531 const std::vector<GridIndexType>& indices) const
532 {
533 auto it = std::find(indices.begin(), indices.end(), idx);
534 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
535 return std::distance(indices.begin(), it);
536 }
537
539 void clear()
540 {
541 scvfIndices_.clear();
542 scvfs_.clear();
543
544 neighborScvIndices_.clear();
545 neighborScvfIndices_.clear();
546 neighborScvs_.clear();
547 neighborScvfs_.clear();
548
549 hasBoundaryScvf_ = false;
550 }
551
552 const GridGeometry* gridGeometryPtr_;
553 std::optional<Element> element_;
554
555 // local storage after binding an element
556 std::array<GridIndexType, 1> scvIndices_;
557 std::vector<GridIndexType> scvfIndices_;
558 std::array<SubControlVolume, 1> scvs_;
559 std::vector<SubControlVolumeFace> scvfs_;
560
561 std::vector<GridIndexType> neighborScvIndices_;
562 std::vector<GridIndexType> neighborScvfIndices_;
563 std::vector<SubControlVolume> neighborScvs_;
564 std::vector<SubControlVolumeFace> neighborScvfs_;
565
566 bool hasBoundaryScvf_ = false;
567};
568
569} // end namespace
570
571#endif
Defines the index types used for grid and local indices.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Class providing iterators over sub control volumes and sub control volume faces of an element.
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
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:215
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models This builds up th...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:49
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:59
const SubControlVolume & scv(GridIndexType scvIdx) const
Get an element sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:85
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:138
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:81
const Element & element() const
The bound element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:161
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:150
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:157
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:109
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:169
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:72
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:122
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get an element sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:91
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:98
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:165
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:132
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:68
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:70
void bind(const Element &element)
Binding of an element, called by the local assembler to prepare element assembly.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:144
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:74
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:186
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:226
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:216
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:260
void bind(const Element &element)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:276
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:334
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:318
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:237
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:200
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:198
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:267
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:202
const Element & element() const
The bound element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:330
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:326
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:271
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:338
GG GridGeometry
export type of finite volume grid geometrys
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:204
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:248
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:211
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82