26#ifndef DUMUX_DISCRETIZATION_CCMPFA_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CCMPFA_FV_ELEMENT_GEOMETRY_HH
32#include <dune/common/exceptions.hh>
33#include <dune/common/iteratorrange.hh>
50template<
class GG,
bool enableGr
idGeometryCache>
63 using GridView =
typename GG::GridView;
66 static constexpr int dim = GridView::dimension;
70 using Element =
typename GridView::template Codim<0>::Entity;
78 static constexpr std::size_t maxNumElementScvs = 1;
80 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
84 : gridGeometryPtr_(&gridGeometry) {}
89 return gridGeometry().scv(scvIdx);
95 return gridGeometry().scvf(scvfIdx);
102 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
110 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
114 return Dune::IteratorRange<ScvIterator>(
ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
115 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
123 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
126 const auto& g = fvGeometry.gridGeometry();
127 const auto scvIdx = fvGeometry.scvIndices_[0];
129 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
130 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
136 return scvIndices_.size();
142 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
152 this->bindElement(element);
153 return std::move(*
this);
158 this->bindElement(element);
168 this->bindElement(element);
169 return std::move(*
this);
176 scvIndices_[0] = gridGeometry().elementMapper().index(element);
181 {
return static_cast<bool>(element_); }
185 {
return *element_; }
189 {
return *gridGeometryPtr_; }
193 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
197 std::optional<Element> element_;
198 std::array<GridIndexType, 1> scvIndices_;
199 const GridGeometry* gridGeometryPtr_;
211 using GridView =
typename GG::GridView;
213 using MpfaHelper =
typename GG::MpfaHelper;
215 static const int dim = GridView::dimension;
216 static const int dimWorld = GridView::dimensionworld;
217 using CoordScalar =
typename GridView::ctype;
221 using Element =
typename GridView::template Codim<0>::Entity;
229 static constexpr std::size_t maxNumElementScvs = 1;
231 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
235 : gridGeometryPtr_(&gridGeometry) {}
241 if (scvIdx == scvIndices_[0])
244 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
251 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
252 if (it != scvfIndices_.end())
255 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
262 return scvf( gridGeometry().flipScvfIdx(scvfIdx, outsideScvIdx) );
270 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
273 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
274 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
282 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
285 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
286 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
291 {
return scvs_.size(); }
295 {
return scvfs_.size(); }
304 this->bindElement_(element);
305 return std::move(*
this);
310 this->bindElement_(element);
320 this->bind_(element);
321 return std::move(*
this);
326 this->bind_(element);
331 {
return static_cast<bool>(element_); }
335 {
return *element_; }
339 {
return *gridGeometryPtr_; }
343 {
return hasBoundaryScvf_; }
349 void bind_(
const Element& element)
352 bindElement_(element);
355 const auto globalI = gridGeometry().elementMapper().index(element);
356 const auto& assemblyMapI = gridGeometry().connectivityMap()[globalI];
359 const auto numNeighbors = assemblyMapI.size();
360 const auto numNeighborScvfs = numNeighborScvfs_(assemblyMapI);
361 neighborScvs_.reserve(numNeighbors);
362 neighborScvIndices_.reserve(numNeighbors);
363 neighborScvfIndices_.reserve(numNeighborScvfs);
364 neighborScvfs_.reserve(numNeighborScvfs);
368 for (
const auto& dataJ : assemblyMapI)
369 makeNeighborGeometries(gridGeometry().element(dataJ.globalJ),
372 dataJ.additionalScvfs);
391 void bindElement_(
const Element& element)
395 makeElementGeometries(element);
399 template<
class DataJContainer>
400 std::size_t numNeighborScvfs_(
const DataJContainer& dataJContainer)
402 std::size_t numNeighborScvfs = 0;
403 for (
const auto& dataJ : dataJContainer)
404 numNeighborScvfs += dataJ.scvfsJ.size() + dataJ.additionalScvfs.size();
405 return numNeighborScvfs;
409 void makeElementGeometries(
const Element& element)
412 const auto eIdx = gridGeometry().elementMapper().index(element);
413 scvs_[0] = SubControlVolume(
element.geometry(), eIdx);
414 scvIndices_[0] = eIdx;
417 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
418 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
421 static const auto q = getParam<CoordScalar>(
"MPFA.Q");
424 const auto numLocalScvf = scvFaceIndices.size();
425 scvfIndices_.reserve(numLocalScvf);
426 scvfs_.reserve(numLocalScvf);
430 std::vector<bool> finishedFacets;
432 finishedFacets.resize(
element.subEntities(1),
false);
435 for (
const auto& is : intersections(gridGeometry().gridView(), element))
441 const auto indexInInside = is.indexInInside();
442 if (finishedFacets[indexInInside])
445 finishedFacets[indexInInside] =
true;
449 bool useNeighbor = is.neighbor() && is.outside().level() >
element.level();
450 const auto& e = useNeighbor ? is.outside() :
element;
451 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
452 const auto eg = e.geometry();
453 const auto refElement = referenceElement(eg);
456 const auto numCorners = is.geometry().corners();
457 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
466 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
467 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
470 if (gridGeometry().isGhostVertex(vIdxGlobal))
473 hasBoundaryScvf_ = (hasBoundaryScvf_ || is.boundary());
475 scvfs_.emplace_back(MpfaHelper(),
476 MpfaHelper::getScvfCorners(isPositions,
numCorners, c),
480 scvFaceIndices[scvfCounter],
482 neighborVolVarIndices[scvfCounter],
486 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
493 template<
typename IndexVector>
494 void makeNeighborGeometries(
const Element& element,
495 GridIndexType eIdxGlobal,
496 const IndexVector& scvfIndices,
497 const IndexVector& additionalScvfs)
500 neighborScvs_.emplace_back(
element.geometry(), eIdxGlobal);
501 neighborScvIndices_.push_back(eIdxGlobal);
504 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdxGlobal);
505 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdxGlobal);
508 static const auto q = getParam<CoordScalar>(
"MPFA.Q");
512 std::vector<bool> finishedFacets;
514 finishedFacets.resize(
element.subEntities(1),
false);
517 for (
const auto& is : intersections(gridGeometry().gridView(), element))
523 auto indexInInside = is.indexInInside();
524 if(finishedFacets[indexInInside])
527 finishedFacets[indexInInside] =
true;
531 bool useNeighbor = is.neighbor() && is.outside().level() >
element.level();
532 const auto& e = useNeighbor ? is.outside() :
element;
533 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
534 const auto eg = e.geometry();
535 const auto refElement = referenceElement(eg);
538 const auto numCorners = is.geometry().corners();
539 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
548 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
549 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
552 if (gridGeometry().isGhostVertex(vIdxGlobal))
556 if (!MpfaHelper::vectorContainsValue(scvfIndices, scvFaceIndices[scvfCounter])
557 && !MpfaHelper::vectorContainsValue(additionalScvfs, scvFaceIndices[scvfCounter]))
565 neighborScvfs_.emplace_back(MpfaHelper(),
566 MpfaHelper::getScvfCorners(isPositions,
numCorners, c),
570 scvFaceIndices[scvfCounter],
572 neighborVolVarIndices[scvfCounter],
576 neighborScvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
585 unsigned int findLocalIndex(
const GridIndexType idx,
586 const std::vector<GridIndexType>& indices)
const
588 auto it = std::find(indices.begin(), indices.end(), idx);
589 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
596 scvfIndices_.clear();
599 neighborScvIndices_.clear();
600 neighborScvfIndices_.clear();
601 neighborScvs_.clear();
602 neighborScvfs_.clear();
604 hasBoundaryScvf_ =
false;
607 const GridGeometry* gridGeometryPtr_;
608 std::optional<Element> element_;
611 std::array<GridIndexType, 1> scvIndices_;
612 std::vector<GridIndexType> scvfIndices_;
613 std::array<SubControlVolume, 1> scvs_;
614 std::vector<SubControlVolumeFace> scvfs_;
616 std::vector<GridIndexType> neighborScvIndices_;
617 std::vector<GridIndexType> neighborScvfIndices_;
618 std::vector<SubControlVolume> neighborScvs_;
619 std::vector<SubControlVolumeFace> neighborScvfs_;
621 bool hasBoundaryScvf_ =
false;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Defines the index types used for grid and local indices.
Class providing iterators over sub control volumes and sub control volume faces of an element.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:232
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:51
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:61
const SubControlVolume & scv(GridIndexType scvIdx) const
Get an element sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:87
CCMpfaFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:166
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:173
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:140
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:83
const Element & element() const
The bound element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:184
CCMpfaFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:150
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:180
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:111
void bind(const Element &element) &
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:156
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:192
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:74
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:124
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get an element sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:93
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:100
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:188
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:134
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:70
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:72
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:76
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:209
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:249
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:239
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:283
void bindElement(const Element &element) &
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:308
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:338
CCMpfaFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:318
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:260
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:223
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:221
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:290
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:225
const Element & element() const
The bound element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:334
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:330
void bind(const Element &element) &
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:324
CCMpfaFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:302
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:294
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:342
GG GridGeometry
export type of finite volume grid geometries
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:227
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:271
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:234
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82