26#ifndef DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
34#include <dune/common/exceptions.hh>
36#include <dune/common/iteratorrange.hh>
50template<
class GG,
bool enableGr
idGeometryCache>
63 using GridView =
typename GG::GridView;
68 using Element =
typename GridView::template Codim<0>::Entity;
77 static constexpr std::size_t maxNumElementScvs = 1;
79 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
83 : gridGeometryPtr_(&gridGeometry) {}
89 return gridGeometry().scv(scvIdx);
96 return gridGeometry().scvf(scvfIdx);
103 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
111 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
115 return Dune::IteratorRange<ScvIterator>(
ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
116 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
124 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
127 const auto& g = fvGeometry.gridGeometry();
128 const auto scvIdx = fvGeometry.scvIndices_[0];
130 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
131 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
137 return scvIndices_.size();
143 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
149 this->bindElement(element);
156 scvIndices_[0] = gridGeometry().elementMapper().index(*element_);
161 {
return static_cast<bool>(element_); }
165 {
return *element_; }
169 {
return *gridGeometryPtr_; }
173 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
177 std::optional<Element> element_;
178 std::array<GridIndexType, 1> scvIndices_;
179 const GridGeometry* gridGeometryPtr_;
191 using GridView =
typename GG::GridView;
195 static const int dim = GridView::dimension;
196 static const int dimWorld = GridView::dimensionworld;
200 using Element =
typename GridView::template Codim<0>::Entity;
208 static constexpr std::size_t maxNumElementScvs = 1;
210 static constexpr std::size_t maxNumElementScvfs = 2*dim;
214 : gridGeometryPtr_(&gridGeometry) {}
220 if (scvIdx == scvIndices_[0])
223 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
230 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
231 if (it != scvfIndices_.end())
234 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
241 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
242 if (it != scvfIndices_.end())
244 const auto localScvfIdx =
std::distance(scvfIndices_.begin(), it);
245 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
249 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
250 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
251 if (localFlippedIndex < scvfs_.size())
252 return scvfs_[localFlippedIndex];
254 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
263 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
266 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
267 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
275 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
278 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
279 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
284 {
return scvs_.size(); }
288 {
return scvfs_.size(); }
294 bindElement(element);
296 neighborScvs_.reserve(element.subEntities(1));
297 neighborScvfIndices_.reserve(element.subEntities(1));
298 neighborScvfs_.reserve(element.subEntities(1));
300 std::vector<GridIndexType> handledNeighbors;
301 handledNeighbors.reserve(element.subEntities(1));
303 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
306 if (intersection.neighbor())
308 const auto outside = intersection.outside();
309 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
312 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
314 makeNeighborGeometries(outside, outsideIdx);
315 handledNeighbors.push_back(outsideIdx);
321 if (dim < dimWorld || gridGeometry().isPeriodic())
323 flippedScvfIndices_.resize(scvfs_.size());
324 for (
unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
326 const auto& scvf = scvfs_[localScvfIdx];
330 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
331 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
333 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
334 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
336 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
338 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
345 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
346 for (
unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
348 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
349 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
350 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
352 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
363 scvfs_.reserve(element.subEntities(1));
364 scvfIndices_.reserve(element.subEntities(1));
365 makeElementGeometries(element);
370 {
return static_cast<bool>(element_); }
374 {
return *element_; }
378 {
return *gridGeometryPtr_; }
382 {
return hasBoundaryScvf_; }
386 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
388 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
390 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
392 return scvfs_.size() + localNeighborScvfIdx;
397 for (
unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
399 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
400 for (
unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
402 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
404 return localOutsideScvfIdx;
409 DUNE_THROW(Dune::InvalidStateException,
"No flipped version of this scvf found!");
413 void makeElementGeometries(
const Element& element)
415 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
417 const auto eIdx = gridGeometry().elementMapper().index(element);
418 scvs_[0] = SubControlVolume(
element.geometry(), eIdx);
419 scvIndices_[0] = eIdx;
421 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
422 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
427 std::vector<bool> handledScvf;
429 handledScvf.resize(
element.subEntities(1),
false);
432 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
435 if (handledScvf[intersection.indexInInside()])
438 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
439 if (intersection.neighbor() || intersection.boundary())
441 ScvfGridIndexStorage scvIndices;
442 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
443 scvIndices[0] = eIdx;
444 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
446 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
447 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
449 scvfs_.emplace_back(intersection,
450 intersection.geometry(),
451 scvFaceIndices[scvfCounter],
454 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
459 handledScvf[intersection.indexInInside()] =
true;
465 void makeNeighborGeometries(
const Element& element,
const GridIndexType eIdx)
467 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
470 neighborScvs_.emplace_back(
element.geometry(), eIdx);
471 neighborScvIndices_.push_back(eIdx);
473 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
474 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
479 std::vector<bool> handledScvf;
481 handledScvf.resize(
element.subEntities(1),
false);
484 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
487 if (handledScvf[intersection.indexInInside()])
490 if (intersection.neighbor())
493 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
494 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
499 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*element_))
501 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
502 neighborScvfs_.emplace_back(intersection,
503 intersection.geometry(),
504 scvFaceIndices[scvfCounter],
508 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
517 for (
unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
519 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*element_))
521 ScvfGridIndexStorage scvIndices;
522 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
523 scvIndices[0] = eIdx;
524 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
525 neighborScvfs_.emplace_back(intersection,
526 intersection.geometry(),
527 scvFaceIndices[scvfCounter],
531 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
539 handledScvf[intersection.indexInInside()] =
true;
546 else if (intersection.boundary())
549 handledScvf[intersection.indexInInside()] =
true;
555 const LocalIndexType findLocalIndex(
const GridIndexType idx,
556 const std::vector<GridIndexType>& indices)
const
558 auto it = std::find(indices.begin(), indices.end(), idx);
559 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
567 scvfIndices_.clear();
569 flippedScvfIndices_.clear();
571 neighborScvIndices_.clear();
572 neighborScvfIndices_.clear();
573 neighborScvs_.clear();
574 neighborScvfs_.clear();
575 flippedNeighborScvfIndices_.clear();
577 hasBoundaryScvf_ =
false;
580 std::optional<Element> element_;
581 const GridGeometry* gridGeometryPtr_;
584 std::array<GridIndexType, 1> scvIndices_;
585 std::array<SubControlVolume, 1> scvs_;
587 std::vector<GridIndexType> scvfIndices_;
588 std::vector<SubControlVolumeFace> scvfs_;
589 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
591 std::vector<GridIndexType> neighborScvIndices_;
592 std::vector<SubControlVolume> neighborScvs_;
594 std::vector<GridIndexType> neighborScvfIndices_;
595 std::vector<SubControlVolumeFace> neighborScvfs_;
596 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
598 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.
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:138
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:51
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:61
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:68
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:72
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:164
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:82
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:125
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:135
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:172
void bind(const Element &element)
Binding of an element, called by the local jacobian to prepare element assembly.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:147
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:74
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:168
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:141
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:160
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:87
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:153
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:94
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:101
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:70
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:112
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:189
bool isBound() const
Returns true if bind/bindElement has already been called.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:369
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:276
const Element & element() const
The bound element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:373
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:204
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:218
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:283
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:213
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:359
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:287
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:200
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:377
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:239
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:228
void bind(const Element &element)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:292
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:202
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:264
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:206
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:381
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82