3.1-git
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 <dune/common/exceptions.hh>
30#include <dune/common/iteratorrange.hh>
31#include <dune/geometry/referenceelements.hh>
32
36
37namespace Dumux {
38
48template<class GG, bool enableGridGeometryCache>
50
57template<class GG>
59{
61 using GridView = typename GG::GridView;
62 using Element = typename GridView::template Codim<0>::Entity;
63 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
64
65 static constexpr int dim = GridView::dimension;
66
67public:
69 using SubControlVolume = typename GG::SubControlVolume;
71 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
73 using GridGeometry = GG;
74 using FVGridGeometry [[deprecated("Use GridGeometry instead. Will be removed after 3.1!" )]] = GridGeometry;
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 scvIndices_[0] = gridGeometry().elementMapper().index(element);
153 }
154
156 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
158 { return gridGeometry(); }
160 { return *gridGeometryPtr_; }
162 bool hasBoundaryScvf() const
163 { return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
164
165private:
166
167 std::array<GridIndexType, 1> scvIndices_;
168 const GridGeometry* gridGeometryPtr_;
169};
170
176template<class GG>
178{
180 using GridView = typename GG::GridView;
181 using Element = typename GridView::template Codim<0>::Entity;
182 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
183 using MpfaHelper = typename GG::MpfaHelper;
184
185 static const int dim = GridView::dimension;
186 static const int dimWorld = GridView::dimensionworld;
187 using CoordScalar = typename GridView::ctype;
188 using ReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim>;
189
190public:
192 using SubControlVolume = typename GG::SubControlVolume;
194 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
196 using GridGeometry = GG;
197 using FVGridGeometry [[deprecated("Use GridGeometry instead. Will be removed after 3.1!")]] = GridGeometry;
199 static constexpr std::size_t maxNumElementScvs = 1;
201 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
202
205 : gridGeometryPtr_(&gridGeometry) {}
206
209 const SubControlVolume& scv(GridIndexType scvIdx) const
210 {
211 if (scvIdx == scvIndices_[0])
212 return scvs_[0];
213 else
214 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
215 }
216
219 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
220 {
221 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
222 if (it != scvfIndices_.end())
223 return scvfs_[std::distance(scvfIndices_.begin(), it)];
224 else
225 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
226 }
227
230 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
231 {
232 return scvf( gridGeometry().flipScvfIdx(scvfIdx, outsideScvIdx) );
233 }
234
240 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
241 scvs(const ThisType& g)
242 {
243 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
244 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
245 }
246
252 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
253 scvfs(const ThisType& g)
254 {
255 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
256 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
257 }
258
260 std::size_t numScv() const
261 { return scvs_.size(); }
262
264 std::size_t numScvf() const
265 { return scvfs_.size(); }
266
269 void bind(const Element& element)
270 {
271 // make inside geometries
272 bindElement(element);
273
274 // get some references for convenience
275 const auto globalI = gridGeometry().elementMapper().index(element);
276 const auto& assemblyMapI = gridGeometry().connectivityMap()[globalI];
277
278 // reserve memory
279 const auto numNeighbors = assemblyMapI.size();
280 const auto numNeighborScvfs = numNeighborScvfs_(assemblyMapI);
281 neighborScvs_.reserve(numNeighbors);
282 neighborScvIndices_.reserve(numNeighbors);
283 neighborScvfIndices_.reserve(numNeighborScvfs);
284 neighborScvfs_.reserve(numNeighborScvfs);
285
286 // make neighbor geometries
287 // use the assembly map to determine which faces are necessary
288 for (const auto& dataJ : assemblyMapI)
289 makeNeighborGeometries(gridGeometry().element(dataJ.globalJ),
290 dataJ.globalJ,
291 dataJ.scvfsJ,
292 dataJ.additionalScvfs);
293
294 // //! TODO Check if user added additional DOF dependencies, i.e. the residual of DOF globalI depends
295 // //! on additional DOFs not included in the discretization schemes' occupation pattern
296 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
297 // if (!additionalDofDependencies.empty())
298 // {
299 // const auto newNumNeighbors = neighborScvs_.size() + additionalDofDependencies.size();
300 // neighborScvs_.reserve(newNumNeighbors);
301 // neighborScvIndices_.reserve(newNumNeighbors);
302 // for (auto globalJ : additionalDofDependencies)
303 // {
304 // neighborScvs_.emplace_back(gridGeometry().element(globalJ).geometry(), globalJ);
305 // neighborScvIndices_.emplace_back(globalJ);
306 // }
307 // }
308 }
309
311 void bindElement(const Element& element)
312 {
313 clear();
314 makeElementGeometries(element);
315 }
316
318 [[deprecated("Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
320 { return gridGeometry(); }
322 { return *gridGeometryPtr_; }
324 bool hasBoundaryScvf() const
325 { return hasBoundaryScvf_; }
326
327private:
328
330 template<class DataJContainer>
331 std::size_t numNeighborScvfs_(const DataJContainer& dataJContainer)
332 {
333 std::size_t numNeighborScvfs = 0;
334 for (const auto& dataJ : dataJContainer)
335 numNeighborScvfs += dataJ.scvfsJ.size() + dataJ.additionalScvfs.size();
336 return numNeighborScvfs;
337 }
338
340 void makeElementGeometries(const Element& element)
341 {
342 // make the scv
343 const auto eIdx = gridGeometry().elementMapper().index(element);
344 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
345 scvIndices_[0] = eIdx;
346
347 // get data on the scv faces
348 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
349 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
350
351 // the quadrature point parameterizaion to be used on scvfs
352 static const auto q = getParam<CoordScalar>("Mpfa.Q");
353
354 // reserve memory for the scv faces
355 const auto numLocalScvf = scvFaceIndices.size();
356 scvfIndices_.reserve(numLocalScvf);
357 scvfs_.reserve(numLocalScvf);
358
359 // for network grids we only want to do one scvf per half facet
360 // this approach assumes conforming grids at branching facets
361 std::vector<bool> finishedFacets;
362 if (dim < dimWorld)
363 finishedFacets.resize(element.subEntities(1), false);
364
365 int scvfCounter = 0;
366 for (const auto& is : intersections(gridGeometry().gridView(), element))
367 {
368 // if we are dealing with a lower dimensional network
369 // only make a new scvf if we haven't handled it yet
370 if (dim < dimWorld)
371 {
372 const auto indexInInside = is.indexInInside();
373 if (finishedFacets[indexInInside])
374 continue;
375 else
376 finishedFacets[indexInInside] = true;
377 }
378
379 // if outside level > inside level, use the outside element in the following
380 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
381 const auto& e = useNeighbor ? is.outside() : element;
382 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
383 const auto eg = e.geometry();
384 const auto refElement = ReferenceElements::general(eg.type());
385
386 // Set up a container with all relevant positions for scvf corner computation
387 const auto numCorners = is.geometry().corners();
388 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
389 refElement,
390 indexInElement,
391 numCorners);
392
393 // make the scv faces belonging to each corner of the intersection
394 for (int c = 0; c < numCorners; ++c)
395 {
396 // get the global vertex index the scv face is connected to
397 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
398 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
399
400 // do not build scvfs connected to a processor boundary
401 if (gridGeometry().isGhostVertex(vIdxGlobal))
402 continue;
403
404 hasBoundaryScvf_ = (hasBoundaryScvf_ || is.boundary());
405
406 scvfs_.emplace_back(MpfaHelper(),
407 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
408 is.centerUnitOuterNormal(),
409 vIdxGlobal,
410 vIdxLocal,
411 scvFaceIndices[scvfCounter],
412 eIdx,
413 neighborVolVarIndices[scvfCounter],
414 q,
415 is.boundary());
416
417 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
418 scvfCounter++;
419 }
420 }
421 }
422
424 template<typename IndexVector>
425 void makeNeighborGeometries(const Element& element,
426 GridIndexType eIdxGlobal,
427 const IndexVector& scvfIndices,
428 const IndexVector& additionalScvfs)
429 {
430 // create the neighbor scv if it doesn't exist yet
431 neighborScvs_.emplace_back(element.geometry(), eIdxGlobal);
432 neighborScvIndices_.push_back(eIdxGlobal);
433
434 // get data on the scv faces
435 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdxGlobal);
436 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdxGlobal);
437
438 // the quadrature point parameterizaion to be used on scvfs
439 static const auto q = getParam<CoordScalar>("Mpfa.Q");
440
441 // for network grids we only want to do one scvf per half facet
442 // this approach assumes conforming grids at branching facets
443 std::vector<bool> finishedFacets;
444 if (dim < dimWorld)
445 finishedFacets.resize(element.subEntities(1), false);
446
447 int scvfCounter = 0;
448 for (const auto& is : intersections(gridGeometry().gridView(), element))
449 {
450 // if we are dealing with a lower dimensional network
451 // only make a new scvf if we haven't handled it yet
452 if (dim < dimWorld)
453 {
454 auto indexInInside = is.indexInInside();
455 if(finishedFacets[indexInInside])
456 continue;
457 else
458 finishedFacets[indexInInside] = true;
459 }
460
461 // if outside level > inside level, use the outside element in the following
462 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
463 const auto& e = useNeighbor ? is.outside() : element;
464 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
465 const auto eg = e.geometry();
466 const auto refElement = ReferenceElements::general(eg.type());
467
468 // Set up a container with all relevant positions for scvf corner computation
469 const auto numCorners = is.geometry().corners();
470 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
471 refElement,
472 indexInElement,
473 numCorners);
474
475 // make the scv faces belonging to each corner of the intersection
476 for (int c = 0; c < numCorners; ++c)
477 {
478 // only build the scvf if it is in the list of necessary indices
479 if (!MpfaHelper::vectorContainsValue(scvfIndices, scvFaceIndices[scvfCounter])
480 && !MpfaHelper::vectorContainsValue(additionalScvfs, scvFaceIndices[scvfCounter]))
481 {
482 // increment counter either way
483 scvfCounter++;
484 continue;
485 }
486
487 // get the global vertex index the scv face is connected to
488 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
489 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
490
491 // do not build scvfs connected to a processor boundary
492 if (gridGeometry().isGhostVertex(vIdxGlobal))
493 continue;
494
495 // build scvf
496 neighborScvfs_.emplace_back(MpfaHelper(),
497 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
498 is.centerUnitOuterNormal(),
499 vIdxGlobal,
500 vIdxLocal,
501 scvFaceIndices[scvfCounter],
502 eIdxGlobal,
503 neighborVolVarIndices[scvfCounter],
504 q,
505 is.boundary());
506
507 neighborScvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
508
509 // increment counter
510 scvfCounter++;
511 }
512 }
513 }
514
516 unsigned int findLocalIndex(const GridIndexType idx,
517 const std::vector<GridIndexType>& indices) const
518 {
519 auto it = std::find(indices.begin(), indices.end(), idx);
520 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
521 return std::distance(indices.begin(), it);
522 }
523
525 void clear()
526 {
527 scvfIndices_.clear();
528 scvfs_.clear();
529
530 neighborScvIndices_.clear();
531 neighborScvfIndices_.clear();
532 neighborScvs_.clear();
533 neighborScvfs_.clear();
534
535 hasBoundaryScvf_ = false;
536 }
537
538 const GridGeometry* gridGeometryPtr_;
539
540 // local storage after binding an element
541 std::array<GridIndexType, 1> scvIndices_;
542 std::vector<GridIndexType> scvfIndices_;
543 std::array<SubControlVolume, 1> scvs_;
544 std::vector<SubControlVolumeFace> scvfs_;
545
546 std::vector<GridIndexType> neighborScvIndices_;
547 std::vector<GridIndexType> neighborScvfIndices_;
548 std::vector<SubControlVolume> neighborScvs_;
549 std::vector<SubControlVolumeFace> neighborScvfs_;
550
551 bool hasBoundaryScvf_ = false;
552};
553
554} // end namespace
555
556#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.
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
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
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:150
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:162
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:71
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
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:159
GridGeometry FVGridGeometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:74
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:132
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:69
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:73
const GridGeometry & fvGridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:157
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:178
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:219
const GridGeometry & fvGridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:319
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:209
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:253
GridGeometry FVGridGeometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:197
void bind(const Element &element)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:269
const GridGeometry & gridGeometry() const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:321
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:311
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:230
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:192
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:260
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:194
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:264
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:324
GG GridGeometry
export type of finite volume grid geometrys
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:196
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:241
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:204
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82