14#ifndef DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
15#define DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
23#include <dune/common/exceptions.hh>
25#include <dune/common/iteratorrange.hh>
30namespace Detail::Tpfa {
32template<
class Gr
idIndexType>
34 const std::vector<GridIndexType>& indices)
36 auto it = std::find(indices.begin(), indices.end(), idx);
37 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
52template<
class GG,
bool enableGr
idGeometryCache>
65 using GridView =
typename GG::GridView;
71 using Element =
typename GridView::template Codim<0>::Entity;
80 static constexpr std::size_t maxNumElementScvs = 1;
82 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
86 : gridGeometryPtr_(&gridGeometry) {}
92 return gridGeometry().scv(scvIdx);
99 return gridGeometry().scvf(scvfIdx);
106 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
114 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
118 return Dune::IteratorRange<ScvIterator>(
ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
119 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
127 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
130 const auto& g = fvGeometry.gridGeometry();
131 const auto scvIdx = fvGeometry.scvIndices_[0];
133 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
134 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
140 return scvIndices_.size();
146 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
156 this->bindElement(element);
157 return std::move(*
this);
162 this->bindElement(element);
172 this->bindElement(element);
173 return std::move(*
this);
180 scvIndices_[0] = gridGeometry().elementMapper().index(*element_);
185 {
return static_cast<bool>(element_); }
189 {
return *element_; }
193 {
return *gridGeometryPtr_; }
197 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
200 {
return gridGeometryPtr_->element(scv.dofIndex()).geometry(); }
204 const auto element = gridGeometryPtr_->element(scvf.insideScvIdx());
205 const auto& scvfIndices = gridGeometryPtr_->scvfIndicesOfScv(scvf.insideScvIdx());
207 LocalIndexType localIdx = 0;
208 for (
const auto& intersection : intersections(gridGeometryPtr_->gridView(), element))
210 if (intersection.neighbor() || intersection.boundary())
212 if (localIdx == localScvfIdx)
213 return intersection.geometry();
219 DUNE_THROW(Dune::InvalidStateException,
"Could not find scvf geometry");
224 std::optional<Element> element_;
225 std::array<GridIndexType, 1> scvIndices_;
226 const GridGeometry* gridGeometryPtr_;
238 using GridView =
typename GG::GridView;
242 static const int dim = GridView::dimension;
243 static const int dimWorld = GridView::dimensionworld;
247 using Element =
typename GridView::template Codim<0>::Entity;
255 static constexpr std::size_t maxNumElementScvs = 1;
257 static constexpr std::size_t maxNumElementScvfs = 2*dim;
261 : gridGeometryPtr_(&gridGeometry) {}
267 if (scvIdx == scvIndices_[0])
277 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
278 if (it != scvfIndices_.end())
288 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
289 if (it != scvfIndices_.end())
291 const auto localScvfIdx =
std::distance(scvfIndices_.begin(), it);
292 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
297 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
298 if (localFlippedIndex < scvfs_.size())
299 return scvfs_[localFlippedIndex];
301 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
310 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
313 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
314 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
322 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
325 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
326 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
331 {
return scvs_.size(); }
335 {
return scvfs_.size(); }
344 this->bind_(element);
345 return std::move(*
this);
350 this->bind_(element);
360 this->bindElement_(element);
361 return std::move(*
this);
366 this->bindElement_(element);
371 {
return static_cast<bool>(element_); }
375 {
return *element_; }
379 {
return *gridGeometryPtr_; }
383 {
return hasBoundaryScvf_; }
386 {
return gridGeometryPtr_->element(scv.dofIndex()).geometry(); }
390 const auto element = gridGeometryPtr_->element(scvf.insideScvIdx());
391 const auto& scvfIndices = gridGeometryPtr_->scvfIndicesOfScv(scvf.insideScvIdx());
393 LocalIndexType localIdx = 0;
394 for (
const auto& intersection : intersections(gridGeometryPtr_->gridView(), element))
396 if (intersection.neighbor() || intersection.boundary())
398 if (localIdx == localScvfIdx)
399 return intersection.geometry();
405 DUNE_THROW(Dune::InvalidStateException,
"Could not find scvf geometry");
411 void bind_(
const Element& element)
413 bindElement_(element);
415 neighborScvs_.reserve(element.subEntities(1));
416 neighborScvfIndices_.reserve(element.subEntities(1));
417 neighborScvfs_.reserve(element.subEntities(1));
419 std::vector<GridIndexType> handledNeighbors;
420 handledNeighbors.reserve(element.subEntities(1));
422 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
425 if (intersection.neighbor())
427 const auto outside = intersection.outside();
428 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
431 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
433 makeNeighborGeometries(outside, outsideIdx);
434 handledNeighbors.push_back(outsideIdx);
440 if (dim < dimWorld || gridGeometry().isPeriodic())
442 flippedScvfIndices_.resize(scvfs_.size());
443 for (
unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
445 const auto& scvf = scvfs_[localScvfIdx];
449 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
450 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
452 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
453 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
455 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
457 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
464 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
465 for (
unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
467 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
468 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
469 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
471 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
478 void bindElement_(
const Element& element)
482 scvfs_.reserve(
element.subEntities(1));
483 scvfIndices_.reserve(
element.subEntities(1));
484 makeElementGeometries(element);
487 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
489 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
491 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
493 return scvfs_.size() + localNeighborScvfIdx;
498 for (
unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
500 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
501 for (
unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
503 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
505 return localOutsideScvfIdx;
510 DUNE_THROW(Dune::InvalidStateException,
"No flipped version of this scvf found!");
514 void makeElementGeometries(
const Element& element)
516 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
518 const auto eIdx = gridGeometry().elementMapper().index(element);
519 scvs_[0] = SubControlVolume(
element.geometry(), eIdx);
520 scvIndices_[0] = eIdx;
522 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
523 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
528 std::vector<bool> handledScvf;
530 handledScvf.resize(
element.subEntities(1),
false);
533 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
536 if (handledScvf[intersection.indexInInside()])
539 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
540 if (intersection.neighbor() || intersection.boundary())
542 ScvfGridIndexStorage scvIndices;
543 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
544 scvIndices[0] = eIdx;
545 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
547 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
548 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
550 scvfs_.emplace_back(intersection,
551 intersection.geometry(),
552 scvFaceIndices[scvfCounter],
555 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
560 handledScvf[intersection.indexInInside()] =
true;
566 void makeNeighborGeometries(
const Element& element,
const GridIndexType eIdx)
568 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
571 neighborScvs_.emplace_back(
element.geometry(), eIdx);
572 neighborScvIndices_.push_back(eIdx);
574 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
575 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
580 std::vector<bool> handledScvf;
582 handledScvf.resize(
element.subEntities(1),
false);
585 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
588 if (handledScvf[intersection.indexInInside()])
591 if (intersection.neighbor())
594 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
595 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
600 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*element_))
602 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
603 neighborScvfs_.emplace_back(intersection,
604 intersection.geometry(),
605 scvFaceIndices[scvfCounter],
609 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
618 for (
unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
620 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*element_))
622 ScvfGridIndexStorage scvIndices;
623 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
624 scvIndices[0] = eIdx;
625 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
626 neighborScvfs_.emplace_back(intersection,
627 intersection.geometry(),
628 scvFaceIndices[scvfCounter],
632 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
640 handledScvf[intersection.indexInInside()] =
true;
647 else if (intersection.boundary())
650 handledScvf[intersection.indexInInside()] =
true;
659 scvfIndices_.clear();
661 flippedScvfIndices_.clear();
663 neighborScvIndices_.clear();
664 neighborScvfIndices_.clear();
665 neighborScvs_.clear();
666 neighborScvfs_.clear();
667 flippedNeighborScvfIndices_.clear();
669 hasBoundaryScvf_ =
false;
672 std::optional<Element> element_;
673 const GridGeometry* gridGeometryPtr_;
676 std::array<GridIndexType, 1> scvIndices_;
677 std::array<SubControlVolume, 1> scvs_;
679 std::vector<GridIndexType> scvfIndices_;
680 std::vector<SubControlVolumeFace> scvfs_;
681 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
683 std::vector<GridIndexType> neighborScvIndices_;
684 std::vector<SubControlVolume> neighborScvs_;
686 std::vector<GridIndexType> neighborScvfIndices_;
687 std::vector<SubControlVolumeFace> neighborScvfs_;
688 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
690 bool hasBoundaryScvf_ =
false;
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:236
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:370
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:323
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:374
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:251
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:265
Element::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:385
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:330
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:260
GridView::Intersection::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:388
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:334
void bindElement(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:364
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:247
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:342
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:378
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:348
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:286
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:275
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:249
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:311
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:358
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:253
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:382
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:63
GridView::Intersection::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:202
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:170
Element::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:199
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:71
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:75
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:188
void bind(const Element &element) &
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:160
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:85
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:128
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:138
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:196
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:77
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:192
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:144
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:184
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:154
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:90
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:97
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:104
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:73
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:177
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:115
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models This builds up th...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:53
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
Class providing iterators over sub control volumes and sub control volume faces of an element.
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
unsigned int LocalIndex
Definition: indextraits.hh:28