3.2-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;
75 static constexpr std::size_t maxNumElementScvs = 1;
77 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
78
81 : gridGeometryPtr_(&gridGeometry) {}
82
84 const SubControlVolume& scv(GridIndexType scvIdx) const
85 {
86 return gridGeometry().scv(scvIdx);
87 }
88
90 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
91 {
92 return gridGeometry().scvf(scvfIdx);
93 }
94
97 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
98 {
99 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
100 }
101
107 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
108 scvs(const CCMpfaFVElementGeometry& fvGeometry)
109 {
111 return Dune::IteratorRange<ScvIterator>(ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
112 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
113 }
114
120 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
122 {
123 const auto& g = fvGeometry.gridGeometry();
124 const auto scvIdx = fvGeometry.scvIndices_[0];
126 return Dune::IteratorRange<ScvfIterator>(ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
127 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
128 }
129
131 std::size_t numScv() const
132 {
133 return scvIndices_.size();
134 }
135
137 std::size_t numScvf() const
138 {
139 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
140 }
141
143 void bind(const Element& element)
144 {
145 this->bindElement(element);
146 }
147
149 void bindElement(const Element& element)
150 {
151 scvIndices_[0] = gridGeometry().elementMapper().index(element);
152 }
153
156 { return *gridGeometryPtr_; }
157
159 bool hasBoundaryScvf() const
160 { return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
161
162private:
163
164 std::array<GridIndexType, 1> scvIndices_;
165 const GridGeometry* gridGeometryPtr_;
166};
167
173template<class GG>
175{
177 using GridView = typename GG::GridView;
178 using Element = typename GridView::template Codim<0>::Entity;
179 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
180 using MpfaHelper = typename GG::MpfaHelper;
181
182 static const int dim = GridView::dimension;
183 static const int dimWorld = GridView::dimensionworld;
184 using CoordScalar = typename GridView::ctype;
185 using ReferenceElements = typename Dune::ReferenceElements<CoordScalar, dim>;
186
187public:
189 using SubControlVolume = typename GG::SubControlVolume;
191 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
193 using GridGeometry = GG;
195 static constexpr std::size_t maxNumElementScvs = 1;
197 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
198
201 : gridGeometryPtr_(&gridGeometry) {}
202
205 const SubControlVolume& scv(GridIndexType scvIdx) const
206 {
207 if (scvIdx == scvIndices_[0])
208 return scvs_[0];
209 else
210 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
211 }
212
215 const SubControlVolumeFace& scvf(GridIndexType scvfIdx) const
216 {
217 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
218 if (it != scvfIndices_.end())
219 return scvfs_[std::distance(scvfIndices_.begin(), it)];
220 else
221 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
222 }
223
226 const SubControlVolumeFace& flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx = 0) const
227 {
228 return scvf( gridGeometry().flipScvfIdx(scvfIdx, outsideScvIdx) );
229 }
230
236 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
237 scvs(const ThisType& g)
238 {
239 using IteratorType = typename std::array<SubControlVolume, 1>::const_iterator;
240 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
241 }
242
248 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
249 scvfs(const ThisType& g)
250 {
251 using IteratorType = typename std::vector<SubControlVolumeFace>::const_iterator;
252 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
253 }
254
256 std::size_t numScv() const
257 { return scvs_.size(); }
258
260 std::size_t numScvf() const
261 { return scvfs_.size(); }
262
265 void bind(const Element& element)
266 {
267 // make inside geometries
268 bindElement(element);
269
270 // get some references for convenience
271 const auto globalI = gridGeometry().elementMapper().index(element);
272 const auto& assemblyMapI = gridGeometry().connectivityMap()[globalI];
273
274 // reserve memory
275 const auto numNeighbors = assemblyMapI.size();
276 const auto numNeighborScvfs = numNeighborScvfs_(assemblyMapI);
277 neighborScvs_.reserve(numNeighbors);
278 neighborScvIndices_.reserve(numNeighbors);
279 neighborScvfIndices_.reserve(numNeighborScvfs);
280 neighborScvfs_.reserve(numNeighborScvfs);
281
282 // make neighbor geometries
283 // use the assembly map to determine which faces are necessary
284 for (const auto& dataJ : assemblyMapI)
285 makeNeighborGeometries(gridGeometry().element(dataJ.globalJ),
286 dataJ.globalJ,
287 dataJ.scvfsJ,
288 dataJ.additionalScvfs);
289
290 // //! TODO Check if user added additional DOF dependencies, i.e. the residual of DOF globalI depends
291 // //! on additional DOFs not included in the discretization schemes' occupation pattern
292 // const auto& additionalDofDependencies = problem.getAdditionalDofDependencies(globalI);
293 // if (!additionalDofDependencies.empty())
294 // {
295 // const auto newNumNeighbors = neighborScvs_.size() + additionalDofDependencies.size();
296 // neighborScvs_.reserve(newNumNeighbors);
297 // neighborScvIndices_.reserve(newNumNeighbors);
298 // for (auto globalJ : additionalDofDependencies)
299 // {
300 // neighborScvs_.emplace_back(gridGeometry().element(globalJ).geometry(), globalJ);
301 // neighborScvIndices_.emplace_back(globalJ);
302 // }
303 // }
304 }
305
307 void bindElement(const Element& element)
308 {
309 clear();
310 makeElementGeometries(element);
311 }
312
315 { return *gridGeometryPtr_; }
316
318 bool hasBoundaryScvf() const
319 { return hasBoundaryScvf_; }
320
321private:
322
324 template<class DataJContainer>
325 std::size_t numNeighborScvfs_(const DataJContainer& dataJContainer)
326 {
327 std::size_t numNeighborScvfs = 0;
328 for (const auto& dataJ : dataJContainer)
329 numNeighborScvfs += dataJ.scvfsJ.size() + dataJ.additionalScvfs.size();
330 return numNeighborScvfs;
331 }
332
334 void makeElementGeometries(const Element& element)
335 {
336 // make the scv
337 const auto eIdx = gridGeometry().elementMapper().index(element);
338 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
339 scvIndices_[0] = eIdx;
340
341 // get data on the scv faces
342 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
343 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
344
345 // the quadrature point parameterizaion to be used on scvfs
346 static const auto q = []{ // REMOVE deprecated version after release 3.2
347 CoordScalar q = 0.0;
348 bool hasOldParamName = hasParam("Mpfa.Q");
349 if (hasOldParamName) {
350 std::cerr << "Deprecation warning: Parameter Mpfa.Q is deprecated, use MPFA.Q (uppercase MPFA)" << std::endl;
351 q = getParam<CoordScalar>("Mpfa.Q");
352 }
353 bool hasNewParamName = hasParam("MPFA.Q");
354 if (hasNewParamName) {
355 q = getParam<CoordScalar>("MPFA.Q");
356 }
357 if (hasOldParamName | hasNewParamName)
358 return q;
359 else
360 return getParam<CoordScalar>("MPFA.Q");
361 }();
362
363 // reserve memory for the scv faces
364 const auto numLocalScvf = scvFaceIndices.size();
365 scvfIndices_.reserve(numLocalScvf);
366 scvfs_.reserve(numLocalScvf);
367
368 // for network grids we only want to do one scvf per half facet
369 // this approach assumes conforming grids at branching facets
370 std::vector<bool> finishedFacets;
371 if (dim < dimWorld)
372 finishedFacets.resize(element.subEntities(1), false);
373
374 int scvfCounter = 0;
375 for (const auto& is : intersections(gridGeometry().gridView(), element))
376 {
377 // if we are dealing with a lower dimensional network
378 // only make a new scvf if we haven't handled it yet
379 if (dim < dimWorld)
380 {
381 const auto indexInInside = is.indexInInside();
382 if (finishedFacets[indexInInside])
383 continue;
384 else
385 finishedFacets[indexInInside] = true;
386 }
387
388 // if outside level > inside level, use the outside element in the following
389 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
390 const auto& e = useNeighbor ? is.outside() : element;
391 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
392 const auto eg = e.geometry();
393 const auto refElement = ReferenceElements::general(eg.type());
394
395 // Set up a container with all relevant positions for scvf corner computation
396 const auto numCorners = is.geometry().corners();
397 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
398 refElement,
399 indexInElement,
400 numCorners);
401
402 // make the scv faces belonging to each corner of the intersection
403 for (int c = 0; c < numCorners; ++c)
404 {
405 // get the global vertex index the scv face is connected to
406 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
407 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
408
409 // do not build scvfs connected to a processor boundary
410 if (gridGeometry().isGhostVertex(vIdxGlobal))
411 continue;
412
413 hasBoundaryScvf_ = (hasBoundaryScvf_ || is.boundary());
414
415 scvfs_.emplace_back(MpfaHelper(),
416 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
417 is.centerUnitOuterNormal(),
418 vIdxGlobal,
419 vIdxLocal,
420 scvFaceIndices[scvfCounter],
421 eIdx,
422 neighborVolVarIndices[scvfCounter],
423 q,
424 is.boundary());
425
426 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
427 scvfCounter++;
428 }
429 }
430 }
431
433 template<typename IndexVector>
434 void makeNeighborGeometries(const Element& element,
435 GridIndexType eIdxGlobal,
436 const IndexVector& scvfIndices,
437 const IndexVector& additionalScvfs)
438 {
439 // create the neighbor scv if it doesn't exist yet
440 neighborScvs_.emplace_back(element.geometry(), eIdxGlobal);
441 neighborScvIndices_.push_back(eIdxGlobal);
442
443 // get data on the scv faces
444 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdxGlobal);
445 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdxGlobal);
446
447 // the quadrature point parameterizaion to be used on scvfs
448 static const auto q = []{ // REMOVE deprecated version after release 3.2
449 CoordScalar q = 0.0;
450 bool hasOldParamName = hasParam("Mpfa.Q");
451 if (hasOldParamName) {
452 std::cerr << "Deprecation warning: Parameter Mpfa.Q is deprecated, use MPFA.Q (uppercase MPFA)" << std::endl;
453 q = getParam<CoordScalar>("Mpfa.Q");
454 }
455 bool hasNewParamName = hasParam("MPFA.Q");
456 if (hasNewParamName) {
457 q = getParam<CoordScalar>("MPFA.Q");
458 }
459 if (hasOldParamName | hasNewParamName)
460 return q;
461 else
462 return getParam<CoordScalar>("MPFA.Q");
463 }();
464
465 // for network grids we only want to do one scvf per half facet
466 // this approach assumes conforming grids at branching facets
467 std::vector<bool> finishedFacets;
468 if (dim < dimWorld)
469 finishedFacets.resize(element.subEntities(1), false);
470
471 int scvfCounter = 0;
472 for (const auto& is : intersections(gridGeometry().gridView(), element))
473 {
474 // if we are dealing with a lower dimensional network
475 // only make a new scvf if we haven't handled it yet
476 if (dim < dimWorld)
477 {
478 auto indexInInside = is.indexInInside();
479 if(finishedFacets[indexInInside])
480 continue;
481 else
482 finishedFacets[indexInInside] = true;
483 }
484
485 // if outside level > inside level, use the outside element in the following
486 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
487 const auto& e = useNeighbor ? is.outside() : element;
488 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
489 const auto eg = e.geometry();
490 const auto refElement = ReferenceElements::general(eg.type());
491
492 // Set up a container with all relevant positions for scvf corner computation
493 const auto numCorners = is.geometry().corners();
494 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
495 refElement,
496 indexInElement,
497 numCorners);
498
499 // make the scv faces belonging to each corner of the intersection
500 for (int c = 0; c < numCorners; ++c)
501 {
502 // get the global vertex index the scv face is connected to
503 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
504 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
505
506 // do not build scvfs connected to a processor boundary
507 if (gridGeometry().isGhostVertex(vIdxGlobal))
508 continue;
509
510 // only build the scvf if it is in the list of necessary indices
511 if (!MpfaHelper::vectorContainsValue(scvfIndices, scvFaceIndices[scvfCounter])
512 && !MpfaHelper::vectorContainsValue(additionalScvfs, scvFaceIndices[scvfCounter]))
513 {
514 // increment counter either way
515 scvfCounter++;
516 continue;
517 }
518
519 // build scvf
520 neighborScvfs_.emplace_back(MpfaHelper(),
521 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
522 is.centerUnitOuterNormal(),
523 vIdxGlobal,
524 vIdxLocal,
525 scvFaceIndices[scvfCounter],
526 eIdxGlobal,
527 neighborVolVarIndices[scvfCounter],
528 q,
529 is.boundary());
530
531 neighborScvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
532
533 // increment counter
534 scvfCounter++;
535 }
536 }
537 }
538
540 unsigned int findLocalIndex(const GridIndexType idx,
541 const std::vector<GridIndexType>& indices) const
542 {
543 auto it = std::find(indices.begin(), indices.end(), idx);
544 assert(it != indices.end() && "Could not find the scv/scvf! Make sure to properly bind this class!");
545 return std::distance(indices.begin(), it);
546 }
547
549 void clear()
550 {
551 scvfIndices_.clear();
552 scvfs_.clear();
553
554 neighborScvIndices_.clear();
555 neighborScvfIndices_.clear();
556 neighborScvs_.clear();
557 neighborScvfs_.clear();
558
559 hasBoundaryScvf_ = false;
560 }
561
562 const GridGeometry* gridGeometryPtr_;
563
564 // local storage after binding an element
565 std::array<GridIndexType, 1> scvIndices_;
566 std::vector<GridIndexType> scvfIndices_;
567 std::array<SubControlVolume, 1> scvs_;
568 std::vector<SubControlVolumeFace> scvfs_;
569
570 std::vector<GridIndexType> neighborScvIndices_;
571 std::vector<GridIndexType> neighborScvfIndices_;
572 std::vector<SubControlVolume> neighborScvs_;
573 std::vector<SubControlVolumeFace> neighborScvfs_;
574
575 bool hasBoundaryScvf_ = false;
576};
577
578} // end namespace
579
580#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.
bool hasParam(const std::string &param)
Check whether a key exists in the parameter tree.
Definition: parameters.hh:383
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:84
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:137
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:80
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:149
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:108
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:159
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:121
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get an element sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:90
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:97
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:155
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:131
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:143
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:73
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:175
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:215
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:205
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:249
void bind(const Element &element)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:265
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:314
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:307
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:226
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:189
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:256
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:191
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:260
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:318
GG GridGeometry
export type of finite volume grid geometrys
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:193
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:237
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:200
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82