14#ifndef DUMUX_DISCRETIZATION_CCMPFA_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_CCMPFA_FV_ELEMENT_GEOMETRY_HH
20#include <dune/common/exceptions.hh>
21#include <dune/common/iteratorrange.hh>
38template<
class GG,
bool enableGr
idGeometryCache>
51 using GridView =
typename GG::GridView;
54 static constexpr int dim = GridView::dimension;
58 using Element =
typename GridView::template Codim<0>::Entity;
66 static constexpr std::size_t maxNumElementScvs = 1;
68 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
72 : gridGeometryPtr_(&gridGeometry) {}
77 return gridGeometry().scv(scvIdx);
83 return gridGeometry().scvf(scvfIdx);
90 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
98 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
102 return Dune::IteratorRange<ScvIterator>(
ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
103 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
111 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
114 const auto& g = fvGeometry.gridGeometry();
115 const auto scvIdx = fvGeometry.scvIndices_[0];
117 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
118 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
124 return scvIndices_.size();
130 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
140 this->bindElement(element);
141 return std::move(*
this);
146 this->bindElement(element);
156 this->bindElement(element);
157 return std::move(*
this);
164 scvIndices_[0] = gridGeometry().elementMapper().index(element);
169 {
return static_cast<bool>(element_); }
173 {
return *element_; }
177 {
return *gridGeometryPtr_; }
181 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
185 #pragma GCC diagnostic push
186 #pragma GCC diagnostic ignored "-Wdeprecated-declarations"
190 {
return scv.geometry(); }
194 {
return scvf.geometry(); }
196 #pragma GCC diagnostic pop
200 std::optional<Element> element_;
201 std::array<GridIndexType, 1> scvIndices_;
202 const GridGeometry* gridGeometryPtr_;
214 using GridView =
typename GG::GridView;
216 using MpfaHelper =
typename GG::MpfaHelper;
218 static const int dim = GridView::dimension;
219 static const int dimWorld = GridView::dimensionworld;
220 using CoordScalar =
typename GridView::ctype;
224 using Element =
typename GridView::template Codim<0>::Entity;
232 static constexpr std::size_t maxNumElementScvs = 1;
234 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
238 : gridGeometryPtr_(&gridGeometry) {}
244 if (scvIdx == scvIndices_[0])
247 return neighborScvs_[
findLocalIndex(scvIdx, neighborScvIndices_)];
254 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
255 if (it != scvfIndices_.end())
258 return neighborScvfs_[
findLocalIndex(scvfIdx, neighborScvfIndices_)];
265 return scvf( gridGeometry().flipScvfIdx(scvfIdx, outsideScvIdx) );
273 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
276 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
277 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
285 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
288 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
289 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
294 {
return scvs_.size(); }
298 {
return scvfs_.size(); }
307 this->bindElement_(element);
308 return std::move(*
this);
313 this->bindElement_(element);
323 this->bind_(element);
324 return std::move(*
this);
329 this->bind_(element);
334 {
return static_cast<bool>(element_); }
338 {
return *element_; }
342 {
return *gridGeometryPtr_; }
346 {
return hasBoundaryScvf_; }
350 #pragma GCC diagnostic push
351 #pragma GCC diagnostic ignored "-Wdeprecated-declarations"
355 {
return scv.geometry(); }
359 {
return scvf.geometry(); }
361 #pragma GCC diagnostic pop
367 void bind_(
const Element& element)
370 bindElement_(element);
373 const auto globalI = gridGeometry().elementMapper().index(element);
374 const auto& assemblyMapI = gridGeometry().connectivityMap()[globalI];
377 const auto numNeighbors = assemblyMapI.size();
378 const auto numNeighborScvfs = numNeighborScvfs_(assemblyMapI);
379 neighborScvs_.reserve(numNeighbors);
380 neighborScvIndices_.reserve(numNeighbors);
381 neighborScvfIndices_.reserve(numNeighborScvfs);
382 neighborScvfs_.reserve(numNeighborScvfs);
386 for (
const auto& dataJ : assemblyMapI)
387 makeNeighborGeometries(gridGeometry().element(dataJ.globalJ),
390 dataJ.additionalScvfs);
409 void bindElement_(
const Element& element)
413 makeElementGeometries(element);
417 template<
class DataJContainer>
418 std::size_t numNeighborScvfs_(
const DataJContainer& dataJContainer)
420 std::size_t numNeighborScvfs = 0;
421 for (
const auto& dataJ : dataJContainer)
422 numNeighborScvfs += dataJ.scvfsJ.size() + dataJ.additionalScvfs.size();
423 return numNeighborScvfs;
427 void makeElementGeometries(
const Element& element)
430 const auto eIdx = gridGeometry().elementMapper().index(element);
431 scvs_[0] = SubControlVolume(
element.geometry(), eIdx);
432 scvIndices_[0] = eIdx;
435 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
436 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
439 static const auto q = getParam<CoordScalar>(
"MPFA.Q");
442 const auto numLocalScvf = scvFaceIndices.size();
443 scvfIndices_.reserve(numLocalScvf);
444 scvfs_.reserve(numLocalScvf);
448 std::vector<bool> finishedFacets;
450 finishedFacets.resize(
element.subEntities(1),
false);
453 for (
const auto& is : intersections(gridGeometry().gridView(), element))
459 const auto indexInInside = is.indexInInside();
460 if (finishedFacets[indexInInside])
463 finishedFacets[indexInInside] =
true;
467 bool useNeighbor = is.neighbor() && is.outside().level() >
element.level();
468 const auto& e = useNeighbor ? is.outside() :
element;
469 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
470 const auto eg = e.geometry();
471 const auto refElement = referenceElement(eg);
474 const auto numCorners = is.geometry().corners();
475 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
484 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
485 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
488 if (gridGeometry().isGhostVertex(vIdxGlobal))
491 hasBoundaryScvf_ = (hasBoundaryScvf_ || is.boundary());
493 scvfs_.emplace_back(MpfaHelper(),
494 MpfaHelper::getScvfCorners(isPositions,
numCorners, c),
498 scvFaceIndices[scvfCounter],
500 neighborVolVarIndices[scvfCounter],
504 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
511 template<
typename IndexVector>
512 void makeNeighborGeometries(
const Element& element,
513 GridIndexType eIdxGlobal,
514 const IndexVector& scvfIndices,
515 const IndexVector& additionalScvfs)
518 neighborScvs_.emplace_back(
element.geometry(), eIdxGlobal);
519 neighborScvIndices_.push_back(eIdxGlobal);
522 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdxGlobal);
523 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdxGlobal);
526 static const auto q = getParam<CoordScalar>(
"MPFA.Q");
530 std::vector<bool> finishedFacets;
532 finishedFacets.resize(
element.subEntities(1),
false);
535 for (
const auto& is : intersections(gridGeometry().gridView(), element))
541 auto indexInInside = is.indexInInside();
542 if(finishedFacets[indexInInside])
545 finishedFacets[indexInInside] =
true;
549 bool useNeighbor = is.neighbor() && is.outside().level() >
element.level();
550 const auto& e = useNeighbor ? is.outside() :
element;
551 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
552 const auto eg = e.geometry();
553 const auto refElement = referenceElement(eg);
556 const auto numCorners = is.geometry().corners();
557 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
566 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
567 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
570 if (gridGeometry().isGhostVertex(vIdxGlobal))
574 if (!MpfaHelper::vectorContainsValue(scvfIndices, scvFaceIndices[scvfCounter])
575 && !MpfaHelper::vectorContainsValue(additionalScvfs, scvFaceIndices[scvfCounter]))
583 neighborScvfs_.emplace_back(MpfaHelper(),
584 MpfaHelper::getScvfCorners(isPositions,
numCorners, c),
588 scvFaceIndices[scvfCounter],
590 neighborVolVarIndices[scvfCounter],
594 neighborScvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
604 const std::vector<GridIndexType>& indices)
const
606 auto it = std::find(indices.begin(), indices.end(), idx);
607 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
614 scvfIndices_.clear();
617 neighborScvIndices_.clear();
618 neighborScvfIndices_.clear();
619 neighborScvs_.clear();
620 neighborScvfs_.clear();
622 hasBoundaryScvf_ =
false;
625 const GridGeometry* gridGeometryPtr_;
626 std::optional<Element> element_;
629 std::array<GridIndexType, 1> scvIndices_;
630 std::vector<GridIndexType> scvfIndices_;
631 std::array<SubControlVolume, 1> scvs_;
632 std::vector<SubControlVolumeFace> scvfs_;
634 std::vector<GridIndexType> neighborScvIndices_;
635 std::vector<GridIndexType> neighborScvfIndices_;
636 std::vector<SubControlVolume> neighborScvs_;
637 std::vector<SubControlVolumeFace> neighborScvfs_;
639 bool hasBoundaryScvf_ =
false;
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:212
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:252
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:242
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:286
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Create the geometry of a given sub control volume.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:354
void bindElement(const Element &element) &
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:311
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:341
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:321
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:263
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:226
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:224
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:293
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:228
const Element & element() const
The bound element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:337
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:333
void bind(const Element &element) &
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:327
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:305
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:297
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:345
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Create the geometry of a given sub control volume face.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:358
GG GridGeometry
export type of finite volume grid geometries
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:230
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:274
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:237
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:49
const SubControlVolume & scv(GridIndexType scvIdx) const
Get an element sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:75
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:154
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Create the geometry of a given sub control volume face.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:193
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:161
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:128
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:71
const Element & element() const
The bound element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:172
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:138
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:168
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:99
void bind(const Element &element) &
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:144
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:180
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:62
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:112
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get an element sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:81
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:88
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:176
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:122
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:58
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:60
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Create the geometry of a given sub control volume.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:189
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:64
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models This builds up th...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:39
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:30
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:70
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:282
auto findLocalIndex(const GridIndexType idx, const std::vector< GridIndexType > &indices)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:33
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:220
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.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27