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;
64 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();
146 void bind(
const Element& element)
148 this->bindElement(element);
154 elementPtr_ = &element;
155 scvIndices_[0] = gridGeometry().elementMapper().index(*elementPtr_);
160 {
return *gridGeometryPtr_; }
161 [[deprecated(
"Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
163 {
return gridGeometry(); }
167 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
171 const Element* elementPtr_;
172 std::array<GridIndexType, 1> scvIndices_;
173 const GridGeometry* gridGeometryPtr_;
185 using GridView =
typename GG::GridView;
188 using Element =
typename GridView::template Codim<0>::Entity;
190 static const int dim = GridView::dimension;
191 static const int dimWorld = GridView::dimensionworld;
202 static constexpr std::size_t maxNumElementScvs = 1;
204 static constexpr std::size_t maxNumElementScvfs = 2*dim;
208 : gridGeometryPtr_(&gridGeometry) {}
214 if (scvIdx == scvIndices_[0])
217 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
224 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
225 if (it != scvfIndices_.end())
226 return scvfs_[std::distance(scvfIndices_.begin(), it)];
228 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
235 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
236 if (it != scvfIndices_.end())
238 const auto localScvfIdx = std::distance(scvfIndices_.begin(), it);
239 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
243 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
244 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
245 if (localFlippedIndex < scvfs_.size())
246 return scvfs_[localFlippedIndex];
248 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
257 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
260 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
261 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
269 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
272 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
273 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
278 {
return scvs_.size(); }
282 {
return scvfs_.size(); }
286 void bind(
const Element& element)
288 bindElement(element);
290 neighborScvs_.reserve(element.subEntities(1));
291 neighborScvfIndices_.reserve(element.subEntities(1));
292 neighborScvfs_.reserve(element.subEntities(1));
294 std::vector<GridIndexType> handledNeighbors;
295 handledNeighbors.reserve(element.subEntities(1));
297 for (
const auto& intersection :
intersections(gridGeometry().gridView(), element))
300 if (intersection.neighbor())
302 const auto outside = intersection.outside();
303 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
306 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
308 makeNeighborGeometries(outside, outsideIdx);
309 handledNeighbors.push_back(outsideIdx);
315 if (dim < dimWorld || gridGeometry().isPeriodic())
317 flippedScvfIndices_.resize(scvfs_.size());
318 for (
unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
320 const auto& scvf = scvfs_[localScvfIdx];
324 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
325 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
327 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
328 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
330 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
332 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
339 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
340 for (
unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
342 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
343 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
344 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
346 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
356 elementPtr_ = &element;
357 scvfs_.reserve(element.subEntities(1));
358 scvfIndices_.reserve(element.subEntities(1));
359 makeElementGeometries(element);
364 {
return *gridGeometryPtr_; }
365 [[deprecated(
"Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
367 {
return gridGeometry(); }
371 {
return hasBoundaryScvf_; }
375 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
377 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
379 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
381 return scvfs_.size() + localNeighborScvfIdx;
386 for (
unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
388 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
389 for (
unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
391 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
393 return localOutsideScvfIdx;
398 DUNE_THROW(Dune::InvalidStateException,
"No flipped version of this scvf found!");
402 void makeElementGeometries(
const Element& element)
404 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
406 const auto eIdx = gridGeometry().elementMapper().index(element);
407 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
408 scvIndices_[0] = eIdx;
410 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
411 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
416 std::vector<bool> handledScvf;
418 handledScvf.resize(element.subEntities(1),
false);
421 for (
const auto& intersection :
intersections(gridGeometry().gridView(), element))
424 if (handledScvf[intersection.indexInInside()])
427 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
428 if (intersection.neighbor() || intersection.boundary())
430 ScvfGridIndexStorage scvIndices;
431 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
432 scvIndices[0] = eIdx;
433 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
435 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
436 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
438 scvfs_.emplace_back(intersection,
439 intersection.geometry(),
440 scvFaceIndices[scvfCounter],
443 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
448 handledScvf[intersection.indexInInside()] =
true;
454 void makeNeighborGeometries(
const Element& element,
const GridIndexType eIdx)
456 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
459 neighborScvs_.emplace_back(element.geometry(), eIdx);
460 neighborScvIndices_.push_back(eIdx);
462 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
463 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
468 std::vector<bool> handledScvf;
470 handledScvf.resize(element.subEntities(1),
false);
473 for (
const auto& intersection :
intersections(gridGeometry().gridView(), element))
476 if (handledScvf[intersection.indexInInside()])
479 if (intersection.neighbor())
482 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
483 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
488 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*elementPtr_))
490 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
491 neighborScvfs_.emplace_back(intersection,
492 intersection.geometry(),
493 scvFaceIndices[scvfCounter],
497 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
506 for (
unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
508 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*elementPtr_))
510 ScvfGridIndexStorage scvIndices;
511 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
512 scvIndices[0] = eIdx;
513 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
514 neighborScvfs_.emplace_back(intersection,
515 intersection.geometry(),
516 scvFaceIndices[scvfCounter],
520 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
528 handledScvf[intersection.indexInInside()] =
true;
535 else if (intersection.boundary())
538 handledScvf[intersection.indexInInside()] =
true;
544 const LocalIndexType findLocalIndex(
const GridIndexType idx,
545 const std::vector<GridIndexType>& indices)
const
547 auto it = std::find(indices.begin(), indices.end(), idx);
548 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
549 return std::distance(indices.begin(), it);
556 scvfIndices_.clear();
558 flippedScvfIndices_.clear();
560 neighborScvIndices_.clear();
561 neighborScvfIndices_.clear();
562 neighborScvs_.clear();
563 neighborScvfs_.clear();
564 flippedNeighborScvfIndices_.clear();
566 hasBoundaryScvf_ =
false;
569 const Element* elementPtr_;
570 const GridGeometry* gridGeometryPtr_;
573 std::array<GridIndexType, 1> scvIndices_;
574 std::array<SubControlVolume, 1> scvs_;
576 std::vector<GridIndexType> scvfIndices_;
577 std::vector<SubControlVolumeFace> scvfs_;
578 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
580 std::vector<GridIndexType> neighborScvIndices_;
581 std::vector<SubControlVolume> neighborScvs_;
583 std::vector<GridIndexType> neighborScvfIndices_;
584 std::vector<SubControlVolumeFace> neighborScvfs_;
585 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
587 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.
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
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:50
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:60
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:70
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:166
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:72
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:68
GridGeometry FVGridGeometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:73
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:111
const GridGeometry & fvGridGeometry() const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:162
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:183
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:270
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:197
const GridGeometry & fvGridGeometry() const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:366
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:212
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:277
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:207
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:353
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:281
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:363
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:233
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:222
void bind(const Element &element)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:286
GridGeometry FVGridGeometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:200
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:195
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:258
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:199
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:370
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82