26#ifndef DUMUX_DISCRETIZATION_CCMPFA_FV_ELEMENT_GEOMETRY_HH
27#define DUMUX_DISCRETIZATION_CCMPFA_FV_ELEMENT_GEOMETRY_HH
29#include <dune/common/exceptions.hh>
30#include <dune/common/iteratorrange.hh>
31#include <dune/geometry/referenceelements.hh>
48template<
class GG,
bool enableGr
idGeometryCache>
61 using GridView =
typename GG::GridView;
62 using Element =
typename GridView::template Codim<0>::Entity;
65 static constexpr int dim = GridView::dimension;
76 static constexpr std::size_t maxNumElementScvs = 1;
78 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
82 : gridGeometryPtr_(&gridGeometry) {}
87 return gridGeometry().scv(scvIdx);
93 return gridGeometry().scvf(scvfIdx);
100 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
108 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
112 return Dune::IteratorRange<ScvIterator>(
ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
113 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
121 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
124 const auto& g = fvGeometry.gridGeometry();
125 const auto scvIdx = fvGeometry.scvIndices_[0];
127 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
128 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
134 return scvIndices_.size();
140 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
144 void bind(
const Element& element)
146 this->bindElement(element);
152 scvIndices_[0] = gridGeometry().elementMapper().index(element);
156 [[deprecated(
"Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
158 {
return gridGeometry(); }
160 {
return *gridGeometryPtr_; }
163 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
167 std::array<GridIndexType, 1> scvIndices_;
168 const GridGeometry* gridGeometryPtr_;
180 using GridView =
typename GG::GridView;
181 using Element =
typename GridView::template Codim<0>::Entity;
183 using MpfaHelper =
typename GG::MpfaHelper;
185 static const int dim = GridView::dimension;
186 static const int dimWorld = GridView::dimensionworld;
187 using CoordScalar =
typename GridView::ctype;
188 using ReferenceElements =
typename Dune::ReferenceElements<CoordScalar, dim>;
199 static constexpr std::size_t maxNumElementScvs = 1;
201 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
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())
223 return scvfs_[std::distance(scvfIndices_.begin(), it)];
225 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
232 return scvf( gridGeometry().flipScvfIdx(scvfIdx, outsideScvIdx) );
240 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
243 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
244 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
252 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
255 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
256 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
261 {
return scvs_.size(); }
265 {
return scvfs_.size(); }
269 void bind(
const Element& element)
272 bindElement(element);
275 const auto globalI = gridGeometry().elementMapper().index(element);
276 const auto& assemblyMapI = gridGeometry().connectivityMap()[globalI];
279 const auto numNeighbors = assemblyMapI.size();
280 const auto numNeighborScvfs = numNeighborScvfs_(assemblyMapI);
281 neighborScvs_.reserve(numNeighbors);
282 neighborScvIndices_.reserve(numNeighbors);
283 neighborScvfIndices_.reserve(numNeighborScvfs);
284 neighborScvfs_.reserve(numNeighborScvfs);
288 for (
const auto& dataJ : assemblyMapI)
289 makeNeighborGeometries(gridGeometry().element(dataJ.globalJ),
292 dataJ.additionalScvfs);
314 makeElementGeometries(element);
318 [[deprecated(
"Use gridGeometry() instead. fvGridGeometry() will be removed after 3.1!")]]
320 {
return gridGeometry(); }
322 {
return *gridGeometryPtr_; }
325 {
return hasBoundaryScvf_; }
330 template<
class DataJContainer>
331 std::size_t numNeighborScvfs_(
const DataJContainer& dataJContainer)
333 std::size_t numNeighborScvfs = 0;
334 for (
const auto& dataJ : dataJContainer)
335 numNeighborScvfs += dataJ.scvfsJ.size() + dataJ.additionalScvfs.size();
336 return numNeighborScvfs;
340 void makeElementGeometries(
const Element& element)
343 const auto eIdx = gridGeometry().elementMapper().index(element);
344 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
345 scvIndices_[0] = eIdx;
348 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
349 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
352 static const auto q = getParam<CoordScalar>(
"Mpfa.Q");
355 const auto numLocalScvf = scvFaceIndices.size();
356 scvfIndices_.reserve(numLocalScvf);
357 scvfs_.reserve(numLocalScvf);
361 std::vector<bool> finishedFacets;
363 finishedFacets.resize(element.subEntities(1),
false);
366 for (
const auto& is :
intersections(gridGeometry().gridView(), element))
372 const auto indexInInside = is.indexInInside();
373 if (finishedFacets[indexInInside])
376 finishedFacets[indexInInside] =
true;
380 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
381 const auto& e = useNeighbor ? is.outside() : element;
382 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
383 const auto eg = e.geometry();
384 const auto refElement = ReferenceElements::general(eg.type());
387 const auto numCorners = is.geometry().corners();
388 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
394 for (
int c = 0; c < numCorners; ++c)
397 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
398 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
401 if (gridGeometry().isGhostVertex(vIdxGlobal))
404 hasBoundaryScvf_ = (hasBoundaryScvf_ || is.boundary());
406 scvfs_.emplace_back(MpfaHelper(),
407 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
408 is.centerUnitOuterNormal(),
411 scvFaceIndices[scvfCounter],
413 neighborVolVarIndices[scvfCounter],
417 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
424 template<
typename IndexVector>
425 void makeNeighborGeometries(
const Element& element,
426 GridIndexType eIdxGlobal,
427 const IndexVector& scvfIndices,
428 const IndexVector& additionalScvfs)
431 neighborScvs_.emplace_back(element.geometry(), eIdxGlobal);
432 neighborScvIndices_.push_back(eIdxGlobal);
435 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdxGlobal);
436 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdxGlobal);
439 static const auto q = getParam<CoordScalar>(
"Mpfa.Q");
443 std::vector<bool> finishedFacets;
445 finishedFacets.resize(element.subEntities(1),
false);
448 for (
const auto& is :
intersections(gridGeometry().gridView(), element))
454 auto indexInInside = is.indexInInside();
455 if(finishedFacets[indexInInside])
458 finishedFacets[indexInInside] =
true;
462 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
463 const auto& e = useNeighbor ? is.outside() : element;
464 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
465 const auto eg = e.geometry();
466 const auto refElement = ReferenceElements::general(eg.type());
469 const auto numCorners = is.geometry().corners();
470 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
476 for (
int c = 0; c < numCorners; ++c)
479 if (!MpfaHelper::vectorContainsValue(scvfIndices, scvFaceIndices[scvfCounter])
480 && !MpfaHelper::vectorContainsValue(additionalScvfs, scvFaceIndices[scvfCounter]))
488 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
489 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
492 if (gridGeometry().isGhostVertex(vIdxGlobal))
496 neighborScvfs_.emplace_back(MpfaHelper(),
497 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
498 is.centerUnitOuterNormal(),
501 scvFaceIndices[scvfCounter],
503 neighborVolVarIndices[scvfCounter],
507 neighborScvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
516 unsigned int findLocalIndex(
const GridIndexType idx,
517 const std::vector<GridIndexType>& indices)
const
519 auto it = std::find(indices.begin(), indices.end(), idx);
520 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
521 return std::distance(indices.begin(), it);
527 scvfIndices_.clear();
530 neighborScvIndices_.clear();
531 neighborScvfIndices_.clear();
532 neighborScvs_.clear();
533 neighborScvfs_.clear();
535 hasBoundaryScvf_ =
false;
538 const GridGeometry* gridGeometryPtr_;
541 std::array<GridIndexType, 1> scvIndices_;
542 std::vector<GridIndexType> scvfIndices_;
543 std::array<SubControlVolume, 1> scvs_;
544 std::vector<SubControlVolumeFace> scvfs_;
546 std::vector<GridIndexType> neighborScvIndices_;
547 std::vector<GridIndexType> neighborScvfIndices_;
548 std::vector<SubControlVolume> neighborScvs_;
549 std::vector<SubControlVolumeFace> neighborScvfs_;
551 bool hasBoundaryScvf_ =
false;
Defines the index types used for grid and local indices.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models This builds up th...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:49
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:59
const SubControlVolume & scv(GridIndexType scvIdx) const
Get an element sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:85
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:138
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:81
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:150
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:109
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:162
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:71
friend Dune::IteratorRange< ScvfIterator< SubControlVolumeFace, std::vector< GridIndexType >, ThisType > > scvfs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:122
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get an element sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:91
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:98
const GridGeometry & gridGeometry() const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:159
GridGeometry FVGridGeometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:74
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:132
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:69
void bind(const Element &element)
Binding of an element, called by the local assembler to prepare element assembly.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:144
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:73
const GridGeometry & fvGridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:157
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:178
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:219
const GridGeometry & fvGridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:319
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:209
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:253
GridGeometry FVGridGeometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:197
void bind(const Element &element)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:269
const GridGeometry & gridGeometry() const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:321
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:311
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:230
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:192
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:260
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:194
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:264
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:324
GG GridGeometry
export type of finite volume grid geometrys
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:196
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:241
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:204
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82