26#ifndef DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
35#include <dune/common/exceptions.hh>
37#include <dune/common/iteratorrange.hh>
51template<
class GG,
bool enableGr
idGeometryCache>
64 using GridView =
typename GG::GridView;
69 using Element =
typename GridView::template Codim<0>::Entity;
78 static constexpr std::size_t maxNumElementScvs = 1;
80 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
84 : gridGeometryPtr_(&gridGeometry) {}
90 return gridGeometry().scv(scvIdx);
97 return gridGeometry().scvf(scvfIdx);
104 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
112 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
116 return Dune::IteratorRange<ScvIterator>(
ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
117 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
125 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
128 const auto& g = fvGeometry.gridGeometry();
129 const auto scvIdx = fvGeometry.scvIndices_[0];
131 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
132 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
138 return scvIndices_.size();
144 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
154 this->bindElement(element);
155 return std::move(*
this);
160 this->bindElement(element);
170 this->bindElement(element);
171 return std::move(*
this);
178 scvIndices_[0] = gridGeometry().elementMapper().index(*element_);
183 {
return static_cast<bool>(element_); }
187 {
return *element_; }
191 {
return *gridGeometryPtr_; }
195 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
200 return scv.geometry();
206 return scvf.geometry();
211 std::optional<Element> element_;
212 std::array<GridIndexType, 1> scvIndices_;
213 const GridGeometry* gridGeometryPtr_;
225 using GridView =
typename GG::GridView;
229 static const int dim = GridView::dimension;
230 static const int dimWorld = GridView::dimensionworld;
234 using Element =
typename GridView::template Codim<0>::Entity;
242 static constexpr std::size_t maxNumElementScvs = 1;
244 static constexpr std::size_t maxNumElementScvfs = 2*dim;
248 : gridGeometryPtr_(&gridGeometry) {}
254 if (scvIdx == scvIndices_[0])
257 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
264 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
265 if (it != scvfIndices_.end())
268 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
275 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
276 if (it != scvfIndices_.end())
278 const auto localScvfIdx =
std::distance(scvfIndices_.begin(), it);
279 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
283 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
284 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
285 if (localFlippedIndex < scvfs_.size())
286 return scvfs_[localFlippedIndex];
288 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
297 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
300 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
301 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
309 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
312 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
313 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
318 {
return scvs_.size(); }
322 {
return scvfs_.size(); }
331 this->bind_(element);
332 return std::move(*
this);
337 this->bind_(element);
347 this->bindElement_(element);
348 return std::move(*
this);
353 this->bindElement_(element);
358 {
return static_cast<bool>(element_); }
362 {
return *element_; }
366 {
return *gridGeometryPtr_; }
370 {
return hasBoundaryScvf_; }
375 return scv.geometry();
381 return scvf.geometry();
387 void bind_(
const Element& element)
389 bindElement_(element);
391 neighborScvs_.reserve(element.subEntities(1));
392 neighborScvfIndices_.reserve(element.subEntities(1));
393 neighborScvfs_.reserve(element.subEntities(1));
395 std::vector<GridIndexType> handledNeighbors;
396 handledNeighbors.reserve(element.subEntities(1));
398 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
401 if (intersection.neighbor())
403 const auto outside = intersection.outside();
404 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
407 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
409 makeNeighborGeometries(outside, outsideIdx);
410 handledNeighbors.push_back(outsideIdx);
416 if (dim < dimWorld || gridGeometry().isPeriodic())
418 flippedScvfIndices_.resize(scvfs_.size());
419 for (
unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
421 const auto& scvf = scvfs_[localScvfIdx];
425 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
426 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
428 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
429 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
431 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
433 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
440 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
441 for (
unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
443 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
444 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
445 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
447 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
454 void bindElement_(
const Element& element)
458 scvfs_.reserve(
element.subEntities(1));
459 scvfIndices_.reserve(
element.subEntities(1));
460 makeElementGeometries(element);
463 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
465 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
467 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
469 return scvfs_.size() + localNeighborScvfIdx;
474 for (
unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
476 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
477 for (
unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
479 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
481 return localOutsideScvfIdx;
486 DUNE_THROW(Dune::InvalidStateException,
"No flipped version of this scvf found!");
490 void makeElementGeometries(
const Element& element)
492 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
494 const auto eIdx = gridGeometry().elementMapper().index(element);
495 scvs_[0] = SubControlVolume(
element.geometry(), eIdx);
496 scvIndices_[0] = eIdx;
498 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
499 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
504 std::vector<bool> handledScvf;
506 handledScvf.resize(
element.subEntities(1),
false);
509 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
512 if (handledScvf[intersection.indexInInside()])
515 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
516 if (intersection.neighbor() || intersection.boundary())
518 ScvfGridIndexStorage scvIndices;
519 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
520 scvIndices[0] = eIdx;
521 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
523 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
524 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
526 scvfs_.emplace_back(intersection,
527 intersection.geometry(),
528 scvFaceIndices[scvfCounter],
531 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
536 handledScvf[intersection.indexInInside()] =
true;
542 void makeNeighborGeometries(
const Element& element,
const GridIndexType eIdx)
544 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
547 neighborScvs_.emplace_back(
element.geometry(), eIdx);
548 neighborScvIndices_.push_back(eIdx);
550 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
551 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
556 std::vector<bool> handledScvf;
558 handledScvf.resize(
element.subEntities(1),
false);
561 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
564 if (handledScvf[intersection.indexInInside()])
567 if (intersection.neighbor())
570 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
571 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
576 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*element_))
578 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
579 neighborScvfs_.emplace_back(intersection,
580 intersection.geometry(),
581 scvFaceIndices[scvfCounter],
585 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
594 for (
unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
596 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*element_))
598 ScvfGridIndexStorage scvIndices;
599 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
600 scvIndices[0] = eIdx;
601 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
602 neighborScvfs_.emplace_back(intersection,
603 intersection.geometry(),
604 scvFaceIndices[scvfCounter],
608 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
616 handledScvf[intersection.indexInInside()] =
true;
623 else if (intersection.boundary())
626 handledScvf[intersection.indexInInside()] =
true;
632 const LocalIndexType findLocalIndex(
const GridIndexType idx,
633 const std::vector<GridIndexType>& indices)
const
635 auto it = std::find(indices.begin(), indices.end(), idx);
636 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
644 scvfIndices_.clear();
646 flippedScvfIndices_.clear();
648 neighborScvIndices_.clear();
649 neighborScvfIndices_.clear();
650 neighborScvs_.clear();
651 neighborScvfs_.clear();
652 flippedNeighborScvfIndices_.clear();
654 hasBoundaryScvf_ =
false;
657 std::optional<Element> element_;
658 const GridGeometry* gridGeometryPtr_;
661 std::array<GridIndexType, 1> scvIndices_;
662 std::array<SubControlVolume, 1> scvs_;
664 std::vector<GridIndexType> scvfIndices_;
665 std::vector<SubControlVolumeFace> scvfs_;
666 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
668 std::vector<GridIndexType> neighborScvIndices_;
669 std::vector<SubControlVolume> neighborScvs_;
671 std::vector<GridIndexType> neighborScvfIndices_;
672 std::vector<SubControlVolumeFace> neighborScvfs_;
673 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
675 bool hasBoundaryScvf_ =
false;
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
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:52
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:62
CCTpfaFVElementGeometry 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/tpfa/fvelementgeometry.hh:168
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:69
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:73
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:186
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:158
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:83
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:197
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:126
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:136
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:194
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:75
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:190
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:142
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:182
CCTpfaFVElementGeometry 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/tpfa/fvelementgeometry.hh:152
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:88
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:95
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:102
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:71
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:175
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:113
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:203
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:223
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:357
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:310
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:361
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:238
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:252
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:317
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:247
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:321
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:372
void bindElement(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:351
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:234
CCTpfaFVElementGeometry 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/tpfa/fvelementgeometry.hh:329
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:378
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:365
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:335
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:273
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:262
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:236
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:298
CCTpfaFVElementGeometry 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/tpfa/fvelementgeometry.hh:345
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:240
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:369
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82