26#ifndef DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CCTPFA_FV_ELEMENT_GEOMETRY_HH
33#include <dune/common/exceptions.hh>
35#include <dune/common/iteratorrange.hh>
49template<
class GG,
bool enableGr
idGeometryCache>
62 using GridView =
typename GG::GridView;
67 using Element =
typename GridView::template Codim<0>::Entity;
76 static constexpr std::size_t maxNumElementScvs = 1;
78 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
82 : gridGeometryPtr_(&gridGeometry) {}
88 return gridGeometry().scv(scvIdx);
95 return gridGeometry().scvf(scvfIdx);
102 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
110 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
114 return Dune::IteratorRange<ScvIterator>(
ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
115 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
123 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
126 const auto& g = fvGeometry.gridGeometry();
127 const auto scvIdx = fvGeometry.scvIndices_[0];
129 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
130 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
136 return scvIndices_.size();
142 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
148 this->bindElement(element);
154 elementPtr_ = &element;
155 scvIndices_[0] = gridGeometry().elementMapper().index(*elementPtr_);
160 {
return *gridGeometryPtr_; }
164 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
168 const Element* elementPtr_;
169 std::array<GridIndexType, 1> scvIndices_;
170 const GridGeometry* gridGeometryPtr_;
182 using GridView =
typename GG::GridView;
186 static const int dim = GridView::dimension;
187 static const int dimWorld = GridView::dimensionworld;
191 using Element =
typename GridView::template Codim<0>::Entity;
199 static constexpr std::size_t maxNumElementScvs = 1;
201 static constexpr std::size_t maxNumElementScvfs = 2*dim;
205 : gridGeometryPtr_(&gridGeometry) {}
211 if (scvIdx == scvIndices_[0])
214 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
221 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
222 if (it != scvfIndices_.end())
225 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
232 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
233 if (it != scvfIndices_.end())
235 const auto localScvfIdx =
std::distance(scvfIndices_.begin(), it);
236 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
240 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
241 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
242 if (localFlippedIndex < scvfs_.size())
243 return scvfs_[localFlippedIndex];
245 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
254 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
257 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
258 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
266 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
269 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
270 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
275 {
return scvs_.size(); }
279 {
return scvfs_.size(); }
285 bindElement(element);
287 neighborScvs_.reserve(element.subEntities(1));
288 neighborScvfIndices_.reserve(element.subEntities(1));
289 neighborScvfs_.reserve(element.subEntities(1));
291 std::vector<GridIndexType> handledNeighbors;
292 handledNeighbors.reserve(element.subEntities(1));
294 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
297 if (intersection.neighbor())
299 const auto outside = intersection.outside();
300 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
303 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
305 makeNeighborGeometries(outside, outsideIdx);
306 handledNeighbors.push_back(outsideIdx);
312 if (dim < dimWorld || gridGeometry().isPeriodic())
314 flippedScvfIndices_.resize(scvfs_.size());
315 for (
unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
317 const auto& scvf = scvfs_[localScvfIdx];
321 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
322 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
324 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
325 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
327 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
329 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
336 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
337 for (
unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
339 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
340 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
341 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
343 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
353 elementPtr_ = &element;
354 scvfs_.reserve(element.subEntities(1));
355 scvfIndices_.reserve(element.subEntities(1));
356 makeElementGeometries(element);
361 {
return *gridGeometryPtr_; }
365 {
return hasBoundaryScvf_; }
369 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
371 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
373 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
375 return scvfs_.size() + localNeighborScvfIdx;
380 for (
unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
382 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
383 for (
unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
385 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
387 return localOutsideScvfIdx;
392 DUNE_THROW(Dune::InvalidStateException,
"No flipped version of this scvf found!");
396 void makeElementGeometries(
const Element& element)
398 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
400 const auto eIdx = gridGeometry().elementMapper().index(element);
401 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
402 scvIndices_[0] = eIdx;
404 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
405 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
410 std::vector<bool> handledScvf;
412 handledScvf.resize(element.subEntities(1),
false);
415 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
418 if (handledScvf[intersection.indexInInside()])
421 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
422 if (intersection.neighbor() || intersection.boundary())
424 ScvfGridIndexStorage scvIndices;
425 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
426 scvIndices[0] = eIdx;
427 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
429 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
430 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
432 scvfs_.emplace_back(intersection,
433 intersection.geometry(),
434 scvFaceIndices[scvfCounter],
437 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
442 handledScvf[intersection.indexInInside()] =
true;
448 void makeNeighborGeometries(
const Element& element,
const GridIndexType eIdx)
450 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
453 neighborScvs_.emplace_back(element.geometry(), eIdx);
454 neighborScvIndices_.push_back(eIdx);
456 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
457 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
462 std::vector<bool> handledScvf;
464 handledScvf.resize(element.subEntities(1),
false);
467 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
470 if (handledScvf[intersection.indexInInside()])
473 if (intersection.neighbor())
476 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
477 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
482 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*elementPtr_))
484 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
485 neighborScvfs_.emplace_back(intersection,
486 intersection.geometry(),
487 scvFaceIndices[scvfCounter],
491 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
500 for (
unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
502 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*elementPtr_))
504 ScvfGridIndexStorage scvIndices;
505 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
506 scvIndices[0] = eIdx;
507 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
508 neighborScvfs_.emplace_back(intersection,
509 intersection.geometry(),
510 scvFaceIndices[scvfCounter],
514 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
522 handledScvf[intersection.indexInInside()] =
true;
529 else if (intersection.boundary())
532 handledScvf[intersection.indexInInside()] =
true;
538 const LocalIndexType findLocalIndex(
const GridIndexType idx,
539 const std::vector<GridIndexType>& indices)
const
541 auto it = std::find(indices.begin(), indices.end(), idx);
542 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
550 scvfIndices_.clear();
552 flippedScvfIndices_.clear();
554 neighborScvIndices_.clear();
555 neighborScvfIndices_.clear();
556 neighborScvs_.clear();
557 neighborScvfs_.clear();
558 flippedNeighborScvfIndices_.clear();
560 hasBoundaryScvf_ =
false;
563 const Element* elementPtr_;
564 const GridGeometry* gridGeometryPtr_;
567 std::array<GridIndexType, 1> scvIndices_;
568 std::array<SubControlVolume, 1> scvs_;
570 std::vector<GridIndexType> scvfIndices_;
571 std::vector<SubControlVolumeFace> scvfs_;
572 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
574 std::vector<GridIndexType> neighborScvIndices_;
575 std::vector<SubControlVolume> neighborScvs_;
577 std::vector<GridIndexType> neighborScvfIndices_;
578 std::vector<SubControlVolumeFace> neighborScvfs_;
579 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
581 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: geometry/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:50
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:60
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:67
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:71
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:81
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:124
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:134
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:163
void bind(const Element &element)
Binding of an element, called by the local jacobian to prepare element assembly.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:146
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:73
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:159
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:140
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:86
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:152
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:93
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:100
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:69
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:111
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:180
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:267
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:195
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:209
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:274
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:204
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:350
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:278
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:191
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:360
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:230
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:219
void bind(const Element &element)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:283
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:193
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:255
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:197
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:364
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82