26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_FV_GRID_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_FV_GRID_GEOMETRY_HH
51template<
class Gr
idView,
class Traits,
bool enableCache>
55template<
class Gr
idView>
60 DUNE_THROW(Dune::InvalidStateException,
"The ccmpfa discretization method needs at least an overlap of 1 for parallel computations. "
61 <<
" Set the parameter \"Grid.Overlap\" in the input file.");
70template<
class GV,
class Traits>
77 static constexpr int dim = GV::dimension;
78 static constexpr int dimWorld = GV::dimensionworld;
80 using Element =
typename GV::template Codim<0>::Entity;
81 using Vertex =
typename GV::template Codim<dim>::Entity;
82 using Intersection =
typename GV::Intersection;
84 using CoordScalar =
typename GV::ctype;
86 using ScvfOutsideGridIndexStorage =
typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
118 static constexpr int maxElementStencilSize = Traits::maxElementStencilSize;
121 static constexpr bool hasSingleInteractionVolumeType = !MpfaHelper::considerSecondaryIVs();
127 , secondaryIvIndicator_([] (const Element& e, const Intersection& is, bool isBranching)
128 {
return is.boundary() || isBranching; } )
137 , secondaryIvIndicator_(indicator)
146 {
return this->elementMapper(); }
150 {
return scvs_.size(); }
154 {
return scvfs_.size(); }
158 {
return numBoundaryScvf_; }
162 {
return this->gridView().size(0); }
166 template<
bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<useSecondary,
bool> = 0>
168 {
return secondaryInteractionVolumeVertices_[vIdxGlobal]; }
172 template<
bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<!useSecondary,
bool> = 0>
180 ParentType::update(gridView);
187 ParentType::update(std::move(gridView));
197 {
return scvs_[scvIdx]; }
201 {
return scvfs_[scvfIdx]; }
206 {
return connectivityMap_; }
210 {
return ivIndexSets_; }
214 {
return scvfIndicesOfScv_[scvIdx]; }
218 {
return flipScvfIndices_; }
223 {
return scvfs_[flipScvfIndices_[scvfIdx][outsideScvfIdx]]; }
227 {
return hasBoundaryScvf_[eIdx]; }
240 const auto numVert = this->gridView().size(dim);
241 const auto numScvs = numDofs();
242 std::size_t numScvf = MpfaHelper::getGlobalNumScvf(this->gridView());
245 scvs_.resize(numScvs);
246 scvfs_.reserve(numScvf);
247 scvfIndicesOfScv_.resize(numScvs);
248 hasBoundaryScvf_.assign(numScvs,
false);
252 secondaryInteractionVolumeVertices_.assign(numVert,
false);
255 const auto isGhostVertex = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
258 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
261 GridIndexType scvfIdx = 0;
262 numBoundaryScvf_ = 0;
263 for (
const auto& element : elements(this->gridView()))
265 const auto eIdx = this->elementMapper().index(element);
268 auto elementGeometry = element.geometry();
271 std::vector<GridIndexType> scvfIndexSet;
272 scvfIndexSet.reserve(MpfaHelper::getNumLocalScvfs(elementGeometry.type()));
276 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
279 outsideIndices.resize(element.subEntities(1));
280 for (
const auto& intersection : intersections(this->gridView(), element))
282 if (intersection.neighbor())
284 const auto nIdx = this->elementMapper().index( intersection.outside() );
285 outsideIndices[intersection.indexInInside()].push_back(nIdx);
291 for (
const auto& is : intersections(this->gridView(), element))
293 const auto indexInInside = is.indexInInside();
294 const bool boundary = is.boundary();
295 const bool neighbor = is.neighbor();
298 hasBoundaryScvf_[eIdx] =
true;
301 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
305 const bool useNeighbor = neighbor && is.outside().level() >
element.level();
306 const auto& e = useNeighbor ? is.outside() :
element;
307 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
308 const auto eg = e.geometry();
309 const auto refElement = referenceElement(eg);
312 const auto numCorners = is.geometry().corners();
313 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
319 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 :
false;
320 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
326 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
327 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
330 if (isGhostVertex[vIdxGlobal])
335 secondaryInteractionVolumeVertices_[vIdxGlobal] =
true;
338 static const auto q = getParam<CoordScalar>(
"MPFA.Q");
341 const auto& outsideScvIndices = [&] ()
344 return dim == dimWorld ?
345 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
346 outsideIndices[indexInInside];
348 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs) + numBoundaryScvf_++});
351 scvfIndexSet.push_back(scvfIdx);
352 scvfs_.emplace_back(MpfaHelper(),
353 MpfaHelper::getScvfCorners(isPositions,
numCorners, c),
364 dualIdSet[vIdxGlobal].insert(scvfs_.back());
372 outsideIndices[indexInInside].clear();
376 scvs_[eIdx] = SubControlVolume(std::move(elementGeometry), eIdx);
379 scvfIndicesOfScv_[eIdx] = scvfIndexSet;
383 flipScvfIndices_.resize(scvfs_.size());
384 for (
const auto& scvf : scvfs_)
389 const auto numOutsideScvs = scvf.numOutsideScvs();
390 const auto vIdxGlobal = scvf.vertexIndex();
391 const auto insideScvIdx = scvf.insideScvIdx();
393 flipScvfIndices_[scvf.index()].resize(numOutsideScvs);
394 for (std::size_t i = 0; i < numOutsideScvs; ++i)
396 const auto outsideScvIdx = scvf.outsideScvIdx(i);
397 for (
auto outsideScvfIndex : scvfIndicesOfScv_[outsideScvIdx])
399 const auto& outsideScvf = this->scvf(outsideScvfIndex);
400 if (outsideScvf.vertexIndex() == vIdxGlobal &&
401 MpfaHelper::vectorContainsValue(outsideScvf.outsideScvIndices(), insideScvIdx))
403 flipScvfIndices_[scvf.index()][i] = outsideScvfIndex;
412 std::cout <<
"Initializing of the grid finite volume geometry took " << timer.elapsed() <<
" seconds." << std::endl;
416 ivIndexSets_.update(*
this, std::move(dualIdSet));
417 std::cout <<
"Initializing of the grid interaction volume index sets took " << timer.elapsed() <<
" seconds." << std::endl;
421 connectivityMap_.update(*
this);
422 std::cout <<
"Initializing of the connectivity map took " << timer.elapsed() <<
" seconds." << std::endl;
426 ConnectivityMap connectivityMap_;
429 std::vector<SubControlVolume> scvs_;
430 std::vector<SubControlVolumeFace> scvfs_;
433 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
434 std::vector<bool> secondaryInteractionVolumeVertices_;
435 GridIndexType numBoundaryScvf_;
436 std::vector<bool> hasBoundaryScvf_;
439 FlipScvfIndexSet flipScvfIndices_;
442 GridIVIndexSets ivIndexSets_;
445 SecondaryIvIndicatorType secondaryIvIndicator_;
455template<
class GV,
class Traits>
462 static constexpr int dim = GV::dimension;
463 static constexpr int dimWorld = GV::dimensionworld;
465 using Element =
typename GV::template Codim<0>::Entity;
466 using Vertex =
typename GV::template Codim<dim>::Entity;
467 using Intersection =
typename GV::Intersection;
469 using CoordScalar =
typename GV::ctype;
471 using ScvfOutsideGridIndexStorage =
typename Traits::SubControlVolumeFace::Traits::OutsideGridIndexStorage;
503 static constexpr int maxElementStencilSize = Traits::maxElementStencilSize;
506 static constexpr bool hasSingleInteractionVolumeType = !MpfaHelper::considerSecondaryIVs();
512 , secondaryIvIndicator_([] (const Element& e, const Intersection& is, bool isBranching)
513 {
return is.boundary() || isBranching; } )
522 , secondaryIvIndicator_(indicator)
531 {
return this->elementMapper(); }
543 {
return numBoundaryScvf_; }
547 {
return this->gridView().size(0); }
551 template<
bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<useSecondary,
bool> = 0>
553 {
return secondaryInteractionVolumeVertices_[vIdxGlobal]; }
557 template<
bool useSecondary = !hasSingleInteractionVolumeType, std::enable_if_t<!useSecondary,
bool> = 0>
563 {
return isGhostVertex_[this->vertexMapper().index(v)]; }
567 {
return isGhostVertex_[vIdxGlobal]; }
573 ParentType::update(gridView);
580 ParentType::update(std::move(gridView));
590 {
return scvfIndicesOfScv_[scvIdx]; }
594 {
return neighborVolVarIndices_[scvIdx]; }
598 const GridIndexType
flipScvfIdx(GridIndexType scvfIdx,
unsigned int outsideScvfIdx = 0)
const
599 {
return flipScvfIndices_[scvfIdx][outsideScvfIdx]; }
603 {
return flipScvfIndices_; }
608 {
return connectivityMap_; }
612 {
return ivIndexSets_; }
623 numScvs_ = numDofs();
624 scvfIndicesOfScv_.resize(numScvs_);
625 neighborVolVarIndices_.resize(numScvs_);
629 const auto numVert = this->gridView().size(dim);
630 secondaryInteractionVolumeVertices_.assign(numVert,
false);
633 isGhostVertex_ = MpfaHelper::findGhostVertices(this->gridView(), this->vertexMapper());
636 typename GridIVIndexSets::DualGridIndexSet dualIdSet(this->gridView());
639 const auto maxNumScvfs = numScvs_*LocalView::maxNumElementScvfs;
640 std::vector<bool> scvfIsOnBoundary;
641 std::vector<GridIndexType> scvfVertexIndex;
642 scvfIsOnBoundary.reserve(maxNumScvfs);
643 scvfVertexIndex.reserve(maxNumScvfs);
647 numBoundaryScvf_ = 0;
648 for (
const auto& element : elements(this->gridView()))
650 const auto eIdx = this->elementMapper().index(element);
653 auto elementGeometry = element.geometry();
656 const auto numLocalFaces = MpfaHelper::getNumLocalScvfs(elementGeometry.type());
657 std::vector<GridIndexType> scvfsIndexSet;
658 std::vector<ScvfOutsideGridIndexStorage> neighborVolVarIndexSet;
659 scvfsIndexSet.reserve(numLocalFaces);
660 neighborVolVarIndexSet.reserve(numLocalFaces);
664 std::vector<ScvfOutsideGridIndexStorage> outsideIndices;
667 outsideIndices.resize(element.subEntities(1));
668 for (
const auto& intersection : intersections(this->gridView(), element))
670 if (intersection.neighbor())
672 auto nIdx = this->elementMapper().index( intersection.outside() );
673 outsideIndices[intersection.indexInInside()].push_back(nIdx);
679 for (
const auto& is : intersections(this->gridView(), element))
681 const auto indexInInside = is.indexInInside();
682 const bool boundary = is.boundary();
683 const bool neighbor = is.neighbor();
686 if (dim < dimWorld && neighbor && outsideIndices[indexInInside].empty())
690 const bool useNeighbor = neighbor && is.outside().level() >
element.level();
691 const auto& e = useNeighbor ? is.outside() :
element;
692 const auto indexInElement = useNeighbor ? is.indexInOutside() : indexInInside;
693 const auto eg = e.geometry();
694 const auto refElement = referenceElement(eg);
697 const bool isBranchingPoint = dim < dimWorld ? outsideIndices[indexInInside].size() > 1 :
false;
698 const bool usesSecondaryIV = secondaryIvIndicator_(element, is, isBranchingPoint);
701 for (std::size_t c = 0; c < is.geometry().corners(); ++c)
704 const auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
705 const auto vIdxGlobal = this->vertexMapper().subIndex(e, vIdxLocal, dim);
708 if (isGhostVertex_[vIdxGlobal])
713 secondaryInteractionVolumeVertices_[vIdxGlobal] =
true;
716 const auto& outsideScvIndices = [&] ()
719 return dim == dimWorld ?
720 ScvfOutsideGridIndexStorage({this->elementMapper().index(is.outside())}) :
721 outsideIndices[indexInInside];
723 return ScvfOutsideGridIndexStorage({GridIndexType(numScvs_) + numBoundaryScvf_++});
727 dualIdSet[vIdxGlobal].insert(numScvf_, eIdx, boundary);
730 scvfsIndexSet.push_back(numScvf_++);
731 scvfIsOnBoundary.push_back(boundary);
732 scvfVertexIndex.push_back(vIdxGlobal);
733 neighborVolVarIndexSet.emplace_back(std::move(outsideScvIndices));
738 outsideIndices[indexInInside].clear();
742 scvfIndicesOfScv_[eIdx] = scvfsIndexSet;
743 neighborVolVarIndices_[eIdx] = neighborVolVarIndexSet;
747 flipScvfIndices_.resize(numScvf_);
748 for (std::size_t scvIdx = 0; scvIdx < numScvs_; ++scvIdx)
750 const auto& scvfIndices = scvfIndicesOfScv_[scvIdx];
751 for (
unsigned int i = 0; i < scvfIndices.size(); ++i)
754 if (scvfIsOnBoundary[ scvfIndices[i] ])
757 const auto scvfIdx = scvfIndices[i];
758 const auto vIdxGlobal = scvfVertexIndex[scvfIdx];
759 const auto numOutsideScvs = neighborVolVarIndices_[scvIdx][i].size();
761 flipScvfIndices_[scvfIdx].resize(numOutsideScvs);
762 for (
unsigned int j = 0; j < numOutsideScvs; ++j)
764 const auto outsideScvIdx = neighborVolVarIndices_[scvIdx][i][j];
765 const auto& outsideScvfIndices = scvfIndicesOfScv_[outsideScvIdx];
766 for (
unsigned int k = 0; k < outsideScvfIndices.size(); ++k)
768 const auto outsideScvfIndex = outsideScvfIndices[k];
769 const auto outsideScvfVertexIndex = scvfVertexIndex[outsideScvfIndex];
770 const auto& outsideScvfNeighborIndices = neighborVolVarIndices_[outsideScvIdx][k];
771 if (outsideScvfVertexIndex == vIdxGlobal &&
772 MpfaHelper::vectorContainsValue(outsideScvfNeighborIndices, scvIdx))
774 flipScvfIndices_[scvfIdx][j] = outsideScvfIndex;
784 std::cout <<
"Initializing of the grid finite volume geometry took " << timer.elapsed() <<
" seconds." << std::endl;
788 ivIndexSets_.update(*
this, std::move(dualIdSet));
789 std::cout <<
"Initializing of the grid interaction volume index sets took " << timer.elapsed() <<
" seconds." << std::endl;
793 connectivityMap_.update(*
this);
794 std::cout <<
"Initializing of the connectivity map took " << timer.elapsed() <<
" seconds." << std::endl;
798 ConnectivityMap connectivityMap_;
801 std::vector<std::vector<GridIndexType>> scvfIndicesOfScv_;
802 std::vector<std::vector<ScvfOutsideGridIndexStorage>> neighborVolVarIndices_;
803 std::vector<bool> secondaryInteractionVolumeVertices_;
804 std::vector<bool> isGhostVertex_;
805 GridIndexType numScvs_;
806 GridIndexType numScvf_;
807 GridIndexType numBoundaryScvf_;
810 FlipScvfIndexSet flipScvfIndices_;
813 GridIVIndexSets ivIndexSets_;
816 SecondaryIvIndicatorType secondaryIvIndicator_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Defines the index types used for grid and local indices.
Check the overlap size for different discretization methods.
Base class for grid geometries.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
void checkOverlapSizeCCMpfa(const GridView &gridView)
check the overlap size for parallel computations
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:56
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:232
Structure to define the index types used for grid and local indices.
Definition: indextraits.hh:38
Base class for all grid geometries.
Definition: basegridgeometry.hh:61
typename BaseImplementation::GridView GridView
export the grid view type
Definition: basegridgeometry.hh:69
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:52
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:73
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:107
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:153
typename Traits::template MpfaHelper< ThisType > MpfaHelper
export the mpfa helper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:111
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:222
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:213
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:103
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:167
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:205
std::size_t numBoundaryScvf() const
The total number of boundary sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:157
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:161
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
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:178
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:185
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:200
std::size_t numScv() const
The total number of sub control volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:149
const GridIVIndexSets & gridInteractionVolumeIndexSets() const
Returns the grid interaction volume index set class.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:209
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:217
CCMpfaFVGridGeometry(const GridView &gridView, const SecondaryIvIndicatorType &indicator)
Constructor with user-defined indicator function for secondary interaction volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:135
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:105
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:145
bool hasBoundaryScvf(GridIndexType eIdx) const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:226
CCMpfaFVGridGeometry(const GridView &gridView)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:125
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:192
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:196
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:173
The finite volume geometry (scvs and scvfs) for cell-centered mpfa models on a grid view This builds ...
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:458
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:479
typename Traits::SubControlVolumeFace SubControlVolumeFace
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:486
typename Traits::template ConnectivityMap< ThisType > ConnectivityMap
export the connectivity map type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:490
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:566
const GridIndexType flipScvfIdx(GridIndexType scvfIdx, unsigned int outsideScvfIdx=0) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:598
void update(GridView &&gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:578
const GridIVIndexSets & gridInteractionVolumeIndexSets() const
Returns the grid interaction volume seeds class.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:611
Extrusion_t< Traits > Extrusion
export the type of extrusion
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:488
typename Traits::template GridIvIndexSets< ThisType > GridIVIndexSets
export the grid interaction volume index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:477
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:593
CCMpfaFVGridGeometry(const GridView &gridView, const SecondaryIvIndicatorType &indicator)
Constructor with user-defined indicator function for secondary interaction volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:520
const FlipScvfIndexSet & flipScvfIndexSet() const
Returns the flip scvf index set.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:602
CCMpfaFVGridGeometry(const GridView &gridView)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:510
std::size_t numScv() const
Returns the total number of sub control volumes.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:534
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:589
std::vector< ScvfOutsideGridIndexStorage > FlipScvfIndexSet
export the flip scvf index set type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:475
typename Traits::SubControlVolume SubControlVolume
export the type of sub control volume
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:484
std::size_t numBoundaryScvf() const
Returns the number of scvfs on the domain boundary.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:542
typename Traits::ElementMapper DofMapper
export dof mapper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:492
MpfaHelper mpfaHelper() const
Returns instance of the mpfa helper type.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:585
const DofMapper & dofMapper() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:530
const ConnectivityMap & connectivityMap() const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:607
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:562
typename Traits::template MpfaHelper< ThisType > MpfaHelper
export the mpfa helper type
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:496
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:482
bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:552
constexpr bool vertexUsesSecondaryInteractionVolume(GridIndexType vIdxGlobal) const
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:558
std::size_t numScvf() const
Returns the total number of sub control volume faces.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:538
void update(const GridView &gridView)
update all fvElementGeometries (call this after grid adaption)
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:571
std::size_t numDofs() const
Returns the total number of degrees of freedom.
Definition: discretization/cellcentered/mpfa/fvgridgeometry.hh:546
Check if the overlap size is valid for a given discretization method.
Definition: checkoverlapsize.hh:40