26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_FV_GRID_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_FV_GRID_GEOMETRY_HH
29#include <dune/geometry/referenceelements.hh>
50template<
class Gr
idView,
class Traits,
bool enableCache>
54template<
class Gr
idView>
59 DUNE_THROW(Dune::InvalidStateException,
"The ccmpfa discretization method needs at least an overlap of 1 for parallel computations. "
60 <<
" Set the parameter \"Grid.Overlap\" in the input file.");
69template<
class GV,
class Traits>
76 static constexpr int dim = GV::dimension;
77 static constexpr int dimWorld = GV::dimensionworld;
79 using Element =
typename GV::template Codim<0>::Entity;
80 using Vertex =
typename GV::template Codim<dim>::Entity;
81 using Intersection =
typename GV::Intersection;
83 using CoordScalar =
typename GV::ctype;
84 using ReferenceElements =
typename Dune::ReferenceElements<CoordScalar, dim>;
86 using ScvfOutsideGridIndexStorage =
typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
115 static constexpr int maxElementStencilSize = Traits::maxElementStencilSize;
118 static constexpr bool hasSingleInteractionVolumeType = !MpfaHelper::considerSecondaryIVs();
124 , secondaryIvIndicator_([] (const Element& e, const Intersection& is, bool isBranching)
125 {
return is.boundary() || isBranching; } )
133 , secondaryIvIndicator_(indicator)
141 {
return this->elementMapper(); }
145 {
return scvs_.size(); }
149 {
return scvfs_.size(); }
153 {
return numBoundaryScvf_; }
157 {
return this->gridView().size(0); }
161 template<
bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<useSecondary,
bool> = 0>
163 {
return secondaryInteractionVolumeVertices_[vIdxGlobal]; }
167 template<
bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<!useSecondary,
bool> = 0>
174 ParentType::update();
183 const auto numVert = this->gridView().size(dim);
184 const auto numScvs = numDofs();
185 std::size_t numScvf = MpfaHelper::getGlobalNumScvf(this->gridView());
188 scvs_.resize(numScvs);
189 scvfs_.reserve(numScvf);
190 scvfIndicesOfScv_.resize(numScvs);
191 hasBoundaryScvf_.resize(numScvs,
false);
195 secondaryInteractionVolumeVertices_.resize(numVert,
false);
198 const auto isGhostVertex = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
201 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
204 GridIndexType scvfIdx = 0;
205 numBoundaryScvf_ = 0;
206 for (
const auto& element : elements(this->gridView()))
208 const auto eIdx = this->elementMapper().index(element);
211 auto elementGeometry = element.geometry();
214 std::vector<GridIndexType> scvfIndexSet;
215 scvfIndexSet.reserve(MpfaHelper::getNumLocalScvfs(elementGeometry.type()));
219 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
222 outsideIndices.resize(element.subEntities(1));
223 for (
const auto& intersection :
intersections(this->gridView(), element))
225 if (intersection.neighbor())
227 const auto nIdx = this->elementMapper().index( intersection.outside() );
228 outsideIndices[intersection.indexInInside()].push_back(nIdx);
234 for (
const auto& is :
intersections(this->gridView(), element))
236 const auto indexInInside = is.indexInInside();
237 const bool boundary = is.boundary();
238 const bool neighbor = is.neighbor();
241 hasBoundaryScvf_[eIdx] =
true;
244 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
248 const bool useNeighbor = neighbor && is.outside().level() > element.level();
249 const auto& e = useNeighbor ? is.outside() : element;
250 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
251 const auto eg = e.geometry();
252 const auto refElement = ReferenceElements::general(eg.type());
255 const auto numCorners = is.geometry().corners();
256 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
262 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 :
false;
263 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
266 for (std::size_t c = 0; c < numCorners; ++c)
269 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
270 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
273 if (isGhostVertex[vIdxGlobal])
278 secondaryInteractionVolumeVertices_[vIdxGlobal] =
true;
281 static const auto q = getParam<CoordScalar>(
"Mpfa.Q");
284 const auto& outsideScvIndices = [&] ()
287 return dim == dimWorld ?
288 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
289 outsideIndices[indexInInside];
291 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs) + numBoundaryScvf_++});
294 scvfIndexSet.push_back(scvfIdx);
296 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
297 is.centerUnitOuterNormal(),
307 dualIdSet[vIdxGlobal].insert(scvfs_.back());
315 outsideIndices[indexInInside].clear();
322 scvfIndicesOfScv_[eIdx] = scvfIndexSet;
326 flipScvfIndices_.resize(scvfs_.size());
327 for (
const auto& scvf : scvfs_)
332 const auto numOutsideScvs = scvf.numOutsideScvs();
333 const auto vIdxGlobal = scvf.vertexIndex();
334 const auto insideScvIdx = scvf.insideScvIdx();
336 flipScvfIndices_[scvf.index()].resize(numOutsideScvs);
337 for (std::size_t i = 0; i < numOutsideScvs; ++i)
339 const auto outsideScvIdx = scvf.outsideScvIdx(i);
340 for (
auto outsideScvfIndex : scvfIndicesOfScv_[outsideScvIdx])
342 const auto& outsideScvf = this->scvf(outsideScvfIndex);
343 if (outsideScvf.vertexIndex() == vIdxGlobal &&
344 MpfaHelper::vectorContainsValue(outsideScvf.outsideScvIndices(), insideScvIdx))
346 flipScvfIndices_[scvf.index()][i] = outsideScvfIndex;
355 std::cout <<
"Initializing of the grid finite volume geometry took " << timer.elapsed() <<
" seconds." << std::endl;
359 ivIndexSets_.update(*
this, std::move(dualIdSet));
360 std::cout <<
"Initializing of the grid interaction volume index sets took " << timer.elapsed() <<
" seconds." << std::endl;
364 connectivityMap_.update(*
this);
365 std::cout <<
"Initializing of the connectivity map took " << timer.elapsed() <<
" seconds." << std::endl;
374 {
return scvs_[scvIdx]; }
378 {
return scvfs_[scvfIdx]; }
383 {
return connectivityMap_; }
387 {
return ivIndexSets_; }
391 {
return scvfIndicesOfScv_[scvIdx]; }
395 {
return flipScvfIndices_; }
400 {
return scvfs_[flipScvfIndices_[scvfIdx][outsideScvfIdx]]; }
404 {
return hasBoundaryScvf_[eIdx]; }
408 ConnectivityMap connectivityMap_;
411 std::vector<SubControlVolume> scvs_;
412 std::vector<SubControlVolumeFace> scvfs_;
415 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
416 std::vector<bool> secondaryInteractionVolumeVertices_;
417 GridIndexType numBoundaryScvf_;
418 std::vector<bool> hasBoundaryScvf_;
421 FlipScvfIndexSet flipScvfIndices_;
424 GridIVIndexSets ivIndexSets_;
427 SecondaryIvIndicatorType secondaryIvIndicator_;
437template<
class GV,
class Traits>
444 static constexpr int dim = GV::dimension;
445 static constexpr int dimWorld = GV::dimensionworld;
447 using Element =
typename GV::template Codim<0>::Entity;
448 using Vertex =
typename GV::template Codim<dim>::Entity;
449 using Intersection =
typename GV::Intersection;
451 using CoordScalar =
typename GV::ctype;
452 using ReferenceElements =
typename Dune::ReferenceElements<CoordScalar, dim>;
454 using ScvfOutsideGridIndexStorage =
typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
483 static constexpr int maxElementStencilSize = Traits::maxElementStencilSize;
486 static constexpr bool hasSingleInteractionVolumeType = !MpfaHelper::considerSecondaryIVs();
492 , secondaryIvIndicator_([] (const Element& e, const Intersection& is, bool isBranching)
493 {
return is.boundary() || isBranching; } )
501 , secondaryIvIndicator_(indicator)
509 {
return this->elementMapper(); }
521 {
return numBoundaryScvf_; }
525 {
return this->gridView().size(0); }
529 template<
bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<useSecondary,
bool> = 0>
531 {
return secondaryInteractionVolumeVertices_[vIdxGlobal]; }
535 template<
bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<!useSecondary,
bool> = 0>
541 {
return isGhostVertex_[this->vertexMapper().index(v)]; }
545 {
return isGhostVertex_[vIdxGlobal]; }
550 ParentType::update();
556 numScvs_ = numDofs();
557 scvfIndicesOfScv_.resize(numScvs_);
558 neighborVolVarIndices_.resize(numScvs_);
562 const auto numVert = this->gridView().size(dim);
563 secondaryInteractionVolumeVertices_.resize(numVert,
false);
566 isGhostVertex_ = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
569 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
572 const auto maxNumScvfs = numScvs_*LocalView::maxNumElementScvfs;
573 std::vector<bool> scvfIsOnBoundary;
574 std::vector<GridIndexType> scvfVertexIndex;
575 scvfIsOnBoundary.reserve(maxNumScvfs);
576 scvfVertexIndex.reserve(maxNumScvfs);
580 numBoundaryScvf_ = 0;
581 for (
const auto& element : elements(this->gridView()))
583 const auto eIdx = this->elementMapper().index(element);
586 auto elementGeometry = element.geometry();
589 const auto numLocalFaces = MpfaHelper::getNumLocalScvfs(elementGeometry.type());
590 std::vector<GridIndexType> scvfsIndexSet;
591 std::vector<ScvfOutsideGridIndexStorage> neighborVolVarIndexSet;
592 scvfsIndexSet.reserve(numLocalFaces);
593 neighborVolVarIndexSet.reserve(numLocalFaces);
597 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
600 outsideIndices.resize(element.subEntities(1));
601 for (
const auto& intersection :
intersections(this->gridView(), element))
603 if (intersection.neighbor())
605 auto nIdx = this->elementMapper().index( intersection.outside() );
606 outsideIndices[intersection.indexInInside()].push_back(nIdx);
612 for (
const auto& is :
intersections(this->gridView(), element))
614 const auto indexInInside = is.indexInInside();
615 const bool boundary = is.boundary();
616 const bool neighbor = is.neighbor();
619 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
623 const bool useNeighbor = neighbor && is.outside().level() > element.level();
624 const auto& e = useNeighbor ? is.outside() : element;
625 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
626 const auto eg = e.geometry();
627 const auto refElement = ReferenceElements::general(eg.type());
630 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 :
false;
631 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
634 for (std::size_t c = 0; c < is.geometry().corners(); ++c)
637 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
638 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
641 if (isGhostVertex_[vIdxGlobal])
646 secondaryInteractionVolumeVertices_[vIdxGlobal] =
true;
649 const auto& outsideScvIndices = [&] ()
652 return dim == dimWorld ?
653 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
654 outsideIndices[indexInInside];
656 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs_) + numBoundaryScvf_++});
660 dualIdSet[vIdxGlobal].insert(numScvf_, eIdx, boundary);
663 scvfsIndexSet.push_back(numScvf_++);
664 scvfIsOnBoundary.push_back(boundary);
665 scvfVertexIndex.push_back(vIdxGlobal);
666 neighborVolVarIndexSet.emplace_back(std::move(outsideScvIndices));
671 outsideIndices[indexInInside].clear();
675 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
676 neighborVolVarIndices_[eIdx] = neighborVolVarIndexSet;
680 flipScvfIndices_.resize(numScvf_);
681 for (std::size_t scvIdx = 0; scvIdx < numScvs_; ++scvIdx)
683 const auto& scvfIndices = scvfIndicesOfScv_[scvIdx];
684 for (
unsigned int i = 0; i < scvfIndices.size(); ++i)
687 if (scvfIsOnBoundary[ scvfIndices[i] ])
690 const auto scvfIdx = scvfIndices[i];
691 const auto vIdxGlobal = scvfVertexIndex[scvfIdx];
692 const auto numOutsideScvs = neighborVolVarIndices_[scvIdx][i].size();
694 flipScvfIndices_[scvfIdx].resize(numOutsideScvs);
695 for (
unsigned int j = 0; j < numOutsideScvs; ++j)
697 const auto outsideScvIdx = neighborVolVarIndices_[scvIdx][i][j];
698 const auto& outsideScvfIndices = scvfIndicesOfScv_[outsideScvIdx];
699 for (
unsigned int k = 0; k < outsideScvfIndices.size(); ++k)
701 const auto outsideScvfIndex = outsideScvfIndices[k];
702 const auto outsideScvfVertexIndex = scvfVertexIndex[outsideScvfIndex];
703 const auto& outsideScvfNeighborIndices = neighborVolVarIndices_[outsideScvIdx][k];
704 if (outsideScvfVertexIndex == vIdxGlobal &&
705 MpfaHelper::vectorContainsValue(outsideScvfNeighborIndices, scvIdx))
707 flipScvfIndices_[scvfIdx][j] = outsideScvfIndex;
717 std::cout <<
"Initializing of the grid finite volume geometry took " << timer.elapsed() <<
" seconds." << std::endl;
721 ivIndexSets_.update(*
this, std::move(dualIdSet));
722 std::cout <<
"Initializing of the grid interaction volume index sets took " << timer.elapsed() <<
" seconds." << std::endl;
726 connectivityMap_.update(*
this);
727 std::cout <<
"Initializing of the connectivity map took " << timer.elapsed() <<
" seconds." << std::endl;
736 {
return scvfIndicesOfScv_[scvIdx]; }
740 {
return neighborVolVarIndices_[scvIdx]; }
744 const GridIndexType
flipScvfIdx(GridIndexType scvfIdx,
unsigned int outsideScvfIdx = 0)
const
745 {
return flipScvfIndices_[scvfIdx][outsideScvfIdx]; }
749 {
return flipScvfIndices_; }
754 {
return connectivityMap_; }
758 {
return ivIndexSets_; }
762 ConnectivityMap connectivityMap_;
765 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
766 std::vector<std::vector<ScvfOutsideGridIndexStorage>> neighborVolVarIndices_;
767 std::vector<bool> secondaryInteractionVolumeVertices_;
768 std::vector<bool> isGhostVertex_;
769 GridIndexType numScvs_;
770 GridIndexType numScvf_;
771 GridIndexType numBoundaryScvf_;
774 FlipScvfIndexSet flipScvfIndices_;
777 GridIVIndexSets ivIndexSets_;
780 SecondaryIvIndicatorType secondaryIvIndicator_;
Defines the index types used for grid and local indices.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Base class for grid geometries.
Check the overlap size for different discretization methods.
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
void checkOverlapSizeCCMpfa(const GridView &gridView)
check the overlap size for parallel computations
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:55
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
Base class for all finite volume grid geometries.
Definition: basegridgeometry.hh:49
GV GridView
export the grid view type
Definition: basegridgeometry.hh:64
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:51
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:72
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:105
void update()
update all fvElementGeometries (do this again after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:172
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:101
std::size_t numScvf() const
The total number of sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:148
typename Traits::template MpfaHelper< ThisType > MpfaHelper
export the mpfa helper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:109
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:399
const std::vector< GridIndexType > & scvfIndicesOfScv(GridIndexType scvIdx) const
Get the sub control volume face indices of an scv by global index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:390
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:162
std::vector< ScvfOutsideGridIndexStorage > FlipScvfIndexSet
export the flip scvf index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:90
const ConnectivityMap & connectivityMap() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:382
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:152
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:99
typename Traits::template GridIvIndexSets< ThisType > GridIVIndexSets
export the grid interaction volume index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:92
std::size_t numDofs() const
The total number of degrees of freedom.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:156
typename Traits::template LocalView< ThisType, true > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:97
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:377
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:144
const GridIVIndexSets & gridInteractionVolumeIndexSets() const
Returns the grid interaction volume index set class.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:386
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:394
CCMpfaFVGridGeometry(const GridView &gridView, const SecondaryIvIndicatorType &indicator)
Constructor with user-defined indicator function for secondary interaction volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:131
std::function< bool(const Element &, const Intersection &, bool)> SecondaryIvIndicatorType
export the type to be used for indicators where to use the secondary ivs
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:94
typename Traits::template ConnectivityMap< ThisType > ConnectivityMap
export the connectivity map type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:103
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:140
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:403
CCMpfaFVGridGeometry(const GridView &gridView)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:122
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:369
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:373
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:168
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:440
std::function< bool(const Element &, const Intersection &, bool)> SecondaryIvIndicatorType
export the type to be used for indicators where to use the secondary ivs
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:462
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:469
typename Traits::template ConnectivityMap< ThisType > ConnectivityMap
export the connectivity map type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:471
bool isGhostVertex(GridIndexType vIdxGlobal) const
Returns true if the vertex (index) lies on a processor boundary inside a ghost element.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:544
const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:744
const GridIVIndexSets & gridInteractionVolumeIndexSets() const
Returns the grid interaction volume seeds class.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:757
typename Traits::template GridIvIndexSets< ThisType > GridIVIndexSets
export the grid interaction volume index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:460
const std::vector< ScvfOutsideGridIndexStorage > & neighborVolVarIndices(GridIndexType scvIdx) const
Returns the neighboring vol var indices for each scvf contained in an scv.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:739
CCMpfaFVGridGeometry(const GridView &gridView, const SecondaryIvIndicatorType &indicator)
Constructor with user-defined indicator function for secondary interaction volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:499
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:748
CCMpfaFVGridGeometry(const GridView &gridView)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:490
std::size_t numScv() const
Returns the total number of sub control volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:512
const std::vector< GridIndexType > & scvfIndicesOfScv(GridIndexType scvIdx) const
Returns the sub control volume face indices of an scv by global index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:735
std::vector< ScvfOutsideGridIndexStorage > FlipScvfIndexSet
export the flip scvf index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:458
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:467
std::size_t numBoundaryScvf() const
Returns the number of scvfs on the domain boundary.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:520
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:473
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:731
void update()
Updates all finite volume geometries of the grid. Has to be called again after grid adaption.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:548
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:508
const ConnectivityMap & connectivityMap() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:753
bool isGhostVertex(const Vertex &v) const
Returns true if a given vertex lies on a processor boundary inside a ghost element.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:540
typename Traits::template MpfaHelper< ThisType > MpfaHelper
export the mpfa helper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:477
typename Traits::template LocalView< ThisType, false > LocalView
export the type of the fv element geometry (the local view type)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:465
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:530
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:536
std::size_t numScvf() const
Returns the total number of sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:516
std::size_t numDofs() const
Returns the total number of degrees of freedom.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:524
Check if the overlap size is valid for a given discretization method.
Definition: checkoverlapsize.hh:40