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;
75 static constexpr std::size_t maxNumElementScvs = 1;
77 static constexpr std::size_t maxNumElementScvfs = 2*GridView::dimension;
81 : gridGeometryPtr_(&gridGeometry) {}
87 return gridGeometry().scv(scvIdx);
94 return gridGeometry().scvf(scvfIdx);
101 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
109 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
113 return Dune::IteratorRange<ScvIterator>(
ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
114 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
122 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
125 const auto& g = fvGeometry.gridGeometry();
126 const auto scvIdx = fvGeometry.scvIndices_[0];
128 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
129 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
135 return scvIndices_.size();
141 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
145 void bind(
const Element& element)
147 this->bindElement(element);
153 elementPtr_ = &element;
154 scvIndices_[0] = gridGeometry().elementMapper().index(*elementPtr_);
159 {
return *gridGeometryPtr_; }
163 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
167 const Element* elementPtr_;
168 std::array<GridIndexType, 1> scvIndices_;
169 const GridGeometry* gridGeometryPtr_;
181 using GridView =
typename GG::GridView;
184 using Element =
typename GridView::template Codim<0>::Entity;
186 static const int dim = GridView::dimension;
187 static const int dimWorld = GridView::dimensionworld;
197 static constexpr std::size_t maxNumElementScvs = 1;
199 static constexpr std::size_t maxNumElementScvfs = 2*dim;
203 : gridGeometryPtr_(&gridGeometry) {}
209 if (scvIdx == scvIndices_[0])
212 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
219 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
220 if (it != scvfIndices_.end())
221 return scvfs_[std::distance(scvfIndices_.begin(), it)];
223 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
230 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
231 if (it != scvfIndices_.end())
233 const auto localScvfIdx = std::distance(scvfIndices_.begin(), it);
234 return neighborScvfs_[flippedScvfIndices_[localScvfIdx][outsideScvIdx]];
238 const auto localScvfIdx = findLocalIndex(scvfIdx, neighborScvfIndices_);
239 const auto localFlippedIndex = flippedNeighborScvfIndices_[localScvfIdx][outsideScvIdx];
240 if (localFlippedIndex < scvfs_.size())
241 return scvfs_[localFlippedIndex];
243 return neighborScvfs_[localFlippedIndex - scvfs_.size()];
252 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
255 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
256 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
264 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
267 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
268 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
273 {
return scvs_.size(); }
277 {
return scvfs_.size(); }
281 void bind(
const Element& element)
283 bindElement(element);
285 neighborScvs_.reserve(element.subEntities(1));
286 neighborScvfIndices_.reserve(element.subEntities(1));
287 neighborScvfs_.reserve(element.subEntities(1));
289 std::vector<GridIndexType> handledNeighbors;
290 handledNeighbors.reserve(element.subEntities(1));
292 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
295 if (intersection.neighbor())
297 const auto outside = intersection.outside();
298 const auto outsideIdx = gridGeometry().elementMapper().index(outside);
301 if ( std::find(handledNeighbors.begin(), handledNeighbors.end(), outsideIdx) == handledNeighbors.end() )
303 makeNeighborGeometries(outside, outsideIdx);
304 handledNeighbors.push_back(outsideIdx);
310 if (dim < dimWorld || gridGeometry().isPeriodic())
312 flippedScvfIndices_.resize(scvfs_.size());
313 for (
unsigned int localScvfIdx = 0; localScvfIdx < scvfs_.size(); ++localScvfIdx)
315 const auto& scvf = scvfs_[localScvfIdx];
319 flippedScvfIndices_[localScvfIdx].resize(scvf.numOutsideScvs());
320 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < scvf.numOutsideScvs(); ++localOutsideScvIdx)
322 const auto globalOutsideScvIdx = scvf.outsideScvIdx(localOutsideScvIdx);
323 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
325 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
327 flippedScvfIndices_[localScvfIdx][localOutsideScvIdx] = localNeighborScvfIdx;
334 flippedNeighborScvfIndices_.resize(neighborScvfs_.size());
335 for (
unsigned int localScvfIdx = 0; localScvfIdx < neighborScvfs_.size(); ++localScvfIdx)
337 const auto& neighborScvf = neighborScvfs_[localScvfIdx];
338 flippedNeighborScvfIndices_[localScvfIdx].resize(neighborScvf.numOutsideScvs());
339 for (
unsigned int localOutsideScvIdx = 0; localOutsideScvIdx < neighborScvf.numOutsideScvs(); ++localOutsideScvIdx)
341 flippedNeighborScvfIndices_[localScvfIdx][localOutsideScvIdx] = findFlippedScvfIndex_(neighborScvf.insideScvIdx(), neighborScvf.outsideScvIdx(localOutsideScvIdx));
351 elementPtr_ = &element;
352 scvfs_.reserve(element.subEntities(1));
353 scvfIndices_.reserve(element.subEntities(1));
354 makeElementGeometries(element);
359 {
return *gridGeometryPtr_; }
363 {
return hasBoundaryScvf_; }
367 GridIndexType findFlippedScvfIndex_(GridIndexType insideScvIdx, GridIndexType globalOutsideScvIdx)
369 for (
unsigned int localNeighborScvfIdx = 0; localNeighborScvfIdx < neighborScvfs_.size(); ++localNeighborScvfIdx)
371 if (neighborScvfs_[localNeighborScvfIdx].insideScvIdx() == globalOutsideScvIdx)
373 return scvfs_.size() + localNeighborScvfIdx;
378 for (
unsigned int localOutsideScvfIdx = 0; localOutsideScvfIdx < scvfs_.size(); ++localOutsideScvfIdx)
380 const auto& outsideScvf = scvfs_[localOutsideScvfIdx];
381 for (
unsigned int j = 0; j < outsideScvf.numOutsideScvs(); ++j)
383 if (outsideScvf.outsideScvIdx(j) == insideScvIdx)
385 return localOutsideScvfIdx;
390 DUNE_THROW(Dune::InvalidStateException,
"No flipped version of this scvf found!");
394 void makeElementGeometries(
const Element& element)
396 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
398 const auto eIdx = gridGeometry().elementMapper().index(element);
399 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
400 scvIndices_[0] = eIdx;
402 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
403 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
408 std::vector<bool> handledScvf;
410 handledScvf.resize(element.subEntities(1),
false);
413 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
416 if (handledScvf[intersection.indexInInside()])
419 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
420 if (intersection.neighbor() || intersection.boundary())
422 ScvfGridIndexStorage scvIndices;
423 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
424 scvIndices[0] = eIdx;
425 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
427 const bool onBoundary = intersection.boundary() && !intersection.neighbor();
428 hasBoundaryScvf_ = (hasBoundaryScvf_ || onBoundary);
430 scvfs_.emplace_back(intersection,
431 intersection.geometry(),
432 scvFaceIndices[scvfCounter],
435 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
440 handledScvf[intersection.indexInInside()] =
true;
446 void makeNeighborGeometries(
const Element& element,
const GridIndexType eIdx)
448 using ScvfGridIndexStorage =
typename SubControlVolumeFace::Traits::GridIndexStorage;
451 neighborScvs_.emplace_back(element.geometry(), eIdx);
452 neighborScvIndices_.push_back(eIdx);
454 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
455 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
460 std::vector<bool> handledScvf;
462 handledScvf.resize(element.subEntities(1),
false);
465 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
468 if (handledScvf[intersection.indexInInside()])
471 if (intersection.neighbor())
474 const auto& scvfNeighborVolVarIndices = neighborVolVarIndices[scvfCounter];
475 if (scvfNeighborVolVarIndices[0] < gridGeometry().gridView().size(0))
480 if (scvfNeighborVolVarIndices[0] == gridGeometry().elementMapper().index(*elementPtr_))
482 ScvfGridIndexStorage scvIndices({eIdx, scvfNeighborVolVarIndices[0]});
483 neighborScvfs_.emplace_back(intersection,
484 intersection.geometry(),
485 scvFaceIndices[scvfCounter],
489 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
498 for (
unsigned outsideScvIdx = 0; outsideScvIdx < scvfNeighborVolVarIndices.size(); ++outsideScvIdx)
500 if (scvfNeighborVolVarIndices[outsideScvIdx] == gridGeometry().elementMapper().index(*elementPtr_))
502 ScvfGridIndexStorage scvIndices;
503 scvIndices.resize(scvfNeighborVolVarIndices.size() + 1);
504 scvIndices[0] = eIdx;
505 std::copy(scvfNeighborVolVarIndices.begin(), scvfNeighborVolVarIndices.end(), scvIndices.begin()+1);
506 neighborScvfs_.emplace_back(intersection,
507 intersection.geometry(),
508 scvFaceIndices[scvfCounter],
512 neighborScvfIndices_.push_back(scvFaceIndices[scvfCounter]);
520 handledScvf[intersection.indexInInside()] =
true;
527 else if (intersection.boundary())
530 handledScvf[intersection.indexInInside()] =
true;
536 const LocalIndexType findLocalIndex(
const GridIndexType idx,
537 const std::vector<GridIndexType>& indices)
const
539 auto it = std::find(indices.begin(), indices.end(), idx);
540 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
541 return std::distance(indices.begin(), it);
548 scvfIndices_.clear();
550 flippedScvfIndices_.clear();
552 neighborScvIndices_.clear();
553 neighborScvfIndices_.clear();
554 neighborScvs_.clear();
555 neighborScvfs_.clear();
556 flippedNeighborScvfIndices_.clear();
558 hasBoundaryScvf_ =
false;
561 const Element* elementPtr_;
562 const GridGeometry* gridGeometryPtr_;
565 std::array<GridIndexType, 1> scvIndices_;
566 std::array<SubControlVolume, 1> scvs_;
568 std::vector<GridIndexType> scvfIndices_;
569 std::vector<SubControlVolumeFace> scvfs_;
570 std::vector<std::vector<GridIndexType>> flippedScvfIndices_;
572 std::vector<GridIndexType> neighborScvIndices_;
573 std::vector<SubControlVolume> neighborScvs_;
575 std::vector<GridIndexType> neighborScvfIndices_;
576 std::vector<SubControlVolumeFace> neighborScvfs_;
577 std::vector<std::vector<GridIndexType>> flippedNeighborScvfIndices_;
579 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.
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:80
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:123
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:133
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:162
void bind(const Element &element)
Binding of an element, called by the local jacobian to prepare element assembly.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:145
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:158
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:139
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:85
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:151
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:92
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:99
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:68
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCTpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:110
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:179
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:265
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:193
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:207
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:272
CCTpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:202
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:348
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:276
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:358
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:228
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:217
void bind(const Element &element)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:281
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:191
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:253
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:195
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/tpfa/fvelementgeometry.hh:362
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82