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]); }
199 std::optional<Element> element_;
200 std::array<GridIndexType, 1> scvIndices_;
201 const GridGeometry* gridGeometryPtr_;
213 using GridView =
typename GG::GridView;
217 static const int dim = GridView::dimension;
218 static const int dimWorld = GridView::dimensionworld;
222 using Element =
typename GridView::template Codim<0>::Entity;
230 static constexpr std::size_t maxNumElementScvs = 1;
232 static constexpr std::size_t maxNumElementScvfs = 2*dim;
236 : gridGeometryPtr_(&gridGeometry) {}
242 if (scvIdx == scvIndices_[0])
245 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
252 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
253 if (it != scvfIndices_.end())
256 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
263 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
264 if (it != scvfIndices_.end())
266 const auto localScvfIdx =
std::distance(scvfIndices_.begin(), it);
267 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
271 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
272 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
273 if (localFlippedIndex < scvfs_.size())
274 return scvfs_[localFlippedIndex];
276 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
285 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
288 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
289 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
297 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
300 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
301 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
306 {
return scvs_.size(); }
310 {
return scvfs_.size(); }
319 this->bind_(element);
320 return std::move(*
this);
325 this->bind_(element);
335 this->bindElement_(element);
336 return std::move(*
this);
341 this->bindElement_(element);
346 {
return static_cast<bool>(element_); }
350 {
return *element_; }
354 {
return *gridGeometryPtr_; }
358 {
return hasBoundaryScvf_; }
363 void bind_(
const Element& element)
365 bindElement_(element);
367 neighborScvs_.reserve(element.subEntities(1));
368 neighborScvfIndices_.reserve(element.subEntities(1));
369 neighborScvfs_.reserve(element.subEntities(1));
371 std::vector<GridIndexType> handledNeighbors;
372 handledNeighbors.reserve(element.subEntities(1));
374 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
377 if (intersection.neighbor())
379 const auto outside = intersection.outside();
380 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
383 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
385 makeNeighborGeometries(outside, outsideIdx);
386 handledNeighbors.push_back(outsideIdx);
392 if (dim < dimWorld || gridGeometry().isPeriodic())
394 flippedScvfIndices_.resize(scvfs_.size());
395 for (
unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
397 const auto& scvf = scvfs_[localScvfIdx];
401 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
402 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
404 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
405 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
407 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
409 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
416 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
417 for (
unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
419 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
420 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
421 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
423 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
430 void bindElement_(
const Element& element)
434 scvfs_.reserve(
element.subEntities(1));
435 scvfIndices_.reserve(
element.subEntities(1));
436 makeElementGeometries(element);
439 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
441 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
443 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
445 return scvfs_.size() + localNeighborScvfIdx;
450 for (
unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
452 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
453 for (
unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
455 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
457 return localOutsideScvfIdx;
462 DUNE_THROW(Dune::InvalidStateException,
"No flipped version of this scvf found!");
466 void makeElementGeometries(
const Element& element)
468 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
470 const auto eIdx = gridGeometry().elementMapper().index(element);
471 scvs_[0] = SubControlVolume(
element.geometry(), eIdx);
472 scvIndices_[0] = eIdx;
474 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
475 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
480 std::vector<bool> handledScvf;
482 handledScvf.resize(
element.subEntities(1),
false);
485 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
488 if (handledScvf[intersection.indexInInside()])
491 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
492 if (intersection.neighbor() || intersection.boundary())
494 ScvfGridIndexStorage scvIndices;
495 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
496 scvIndices[0] = eIdx;
497 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
499 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
500 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
502 scvfs_.emplace_back(intersection,
503 intersection.geometry(),
504 scvFaceIndices[scvfCounter],
507 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
512 handledScvf[intersection.indexInInside()] =
true;
518 void makeNeighborGeometries(
const Element& element,
const GridIndexType eIdx)
520 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
523 neighborScvs_.emplace_back(
element.geometry(), eIdx);
524 neighborScvIndices_.push_back(eIdx);
526 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
527 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
532 std::vector<bool> handledScvf;
534 handledScvf.resize(
element.subEntities(1),
false);
537 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
540 if (handledScvf[intersection.indexInInside()])
543 if (intersection.neighbor())
546 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
547 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
552 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*element_))
554 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
555 neighborScvfs_.emplace_back(intersection,
556 intersection.geometry(),
557 scvFaceIndices[scvfCounter],
561 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
570 for (
unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
572 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*element_))
574 ScvfGridIndexStorage scvIndices;
575 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
576 scvIndices[0] = eIdx;
577 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
578 neighborScvfs_.emplace_back(intersection,
579 intersection.geometry(),
580 scvFaceIndices[scvfCounter],
584 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
592 handledScvf[intersection.indexInInside()] =
true;
599 else if (intersection.boundary())
602 handledScvf[intersection.indexInInside()] =
true;
608 const LocalIndexType findLocalIndex(
const GridIndexType idx,
609 const std::vector<GridIndexType>& indices)
const
611 auto it = std::find(indices.begin(), indices.end(), idx);
612 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
620 scvfIndices_.clear();
622 flippedScvfIndices_.clear();
624 neighborScvIndices_.clear();
625 neighborScvfIndices_.clear();
626 neighborScvs_.clear();
627 neighborScvfs_.clear();
628 flippedNeighborScvfIndices_.clear();
630 hasBoundaryScvf_ =
false;
633 std::optional<Element> element_;
634 const GridGeometry* gridGeometryPtr_;
637 std::array<GridIndexType, 1> scvIndices_;
638 std::array<SubControlVolume, 1> scvs_;
640 std::vector<GridIndexType> scvfIndices_;
641 std::vector<SubControlVolumeFace> scvfs_;
642 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
644 std::vector<GridIndexType> neighborScvIndices_;
645 std::vector<SubControlVolume> neighborScvs_;
647 std::vector<GridIndexType> neighborScvfIndices_;
648 std::vector<SubControlVolumeFace> neighborScvfs_;
649 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
651 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:292
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
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
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:211
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:345
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:298
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:349
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:226
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:240
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:305
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:235
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:309
void bindElement(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:339
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:222
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:317
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:353
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:323
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:261
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:250
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:224
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:286
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:333
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:228
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:357
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82