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>
47template<
class GG,
bool enableGr
idGeometryCache>
60 using GridView =
typename GG::GridView;
63 static constexpr int dim = GridView::dimension;
67 using Element =
typename GridView::template Codim<0>::Entity;
75 static constexpr std::size_t maxNumElementScvs = 1;
77 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
81 : gridGeometryPtr_(&gridGeometry) {}
86 return gridGeometry().scv(scvIdx);
92 return gridGeometry().scvf(scvfIdx);
99 return gridGeometry().flipScvf(scvfIdx, outsideScvIdx);
107 friend inline Dune::IteratorRange< ScvIterator<SubControlVolume, std::array<GridIndexType, 1>, ThisType> >
111 return Dune::IteratorRange<ScvIterator>(
ScvIterator(fvGeometry.scvIndices_.begin(), fvGeometry),
112 ScvIterator(fvGeometry.scvIndices_.end(), fvGeometry));
120 friend inline Dune::IteratorRange< ScvfIterator<SubControlVolumeFace, std::vector<GridIndexType>, ThisType> >
123 const auto& g = fvGeometry.gridGeometry();
124 const auto scvIdx = fvGeometry.scvIndices_[0];
126 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(g.scvfIndicesOfScv(scvIdx).begin(), fvGeometry),
127 ScvfIterator(g.scvfIndicesOfScv(scvIdx).end(), fvGeometry));
133 return scvIndices_.size();
139 return gridGeometry().scvfIndicesOfScv(scvIndices_[0]).size();
145 this->bindElement(element);
151 scvIndices_[0] = gridGeometry().elementMapper().index(element);
156 {
return *gridGeometryPtr_; }
160 {
return gridGeometry().hasBoundaryScvf(scvIndices_[0]); }
164 std::array<GridIndexType, 1> scvIndices_;
165 const GridGeometry* gridGeometryPtr_;
177 using GridView =
typename GG::GridView;
179 using MpfaHelper =
typename GG::MpfaHelper;
181 static const int dim = GridView::dimension;
182 static const int dimWorld = GridView::dimensionworld;
183 using CoordScalar =
typename GridView::ctype;
187 using Element =
typename GridView::template Codim<0>::Entity;
195 static constexpr std::size_t maxNumElementScvs = 1;
197 static constexpr std::size_t maxNumElementScvfs = dim == 3 ? 24 : 8;
201 : gridGeometryPtr_(&gridGeometry) {}
207 if (scvIdx == scvIndices_[0])
210 return neighborScvs_[findLocalIndex(scvIdx, neighborScvIndices_)];
217 auto it = std::find(scvfIndices_.begin(), scvfIndices_.end(), scvfIdx);
218 if (it != scvfIndices_.end())
221 return neighborScvfs_[findLocalIndex(scvfIdx, neighborScvfIndices_)];
228 return scvf( gridGeometry().flipScvfIdx(scvfIdx, outsideScvIdx) );
236 friend inline Dune::IteratorRange<typename std::array<SubControlVolume, 1>::const_iterator>
239 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
240 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
248 friend inline Dune::IteratorRange<typename std::vector<SubControlVolumeFace>::const_iterator>
251 using IteratorType =
typename std::vector<SubControlVolumeFace>::const_iterator;
252 return Dune::IteratorRange<IteratorType>(g.scvfs_.begin(), g.scvfs_.end());
257 {
return scvs_.size(); }
261 {
return scvfs_.size(); }
268 bindElement(element);
271 const auto globalI = gridGeometry().elementMapper().index(element);
272 const auto& assemblyMapI = gridGeometry().connectivityMap()[globalI];
275 const auto numNeighbors = assemblyMapI.size();
276 const auto numNeighborScvfs = numNeighborScvfs_(assemblyMapI);
277 neighborScvs_.reserve(numNeighbors);
278 neighborScvIndices_.reserve(numNeighbors);
279 neighborScvfIndices_.reserve(numNeighborScvfs);
280 neighborScvfs_.reserve(numNeighborScvfs);
284 for (
const auto& dataJ : assemblyMapI)
285 makeNeighborGeometries(gridGeometry().element(dataJ.globalJ),
288 dataJ.additionalScvfs);
310 makeElementGeometries(element);
315 {
return *gridGeometryPtr_; }
319 {
return hasBoundaryScvf_; }
324 template<
class DataJContainer>
325 std::size_t numNeighborScvfs_(
const DataJContainer& dataJContainer)
327 std::size_t numNeighborScvfs = 0;
328 for (
const auto& dataJ : dataJContainer)
329 numNeighborScvfs += dataJ.scvfsJ.size() + dataJ.additionalScvfs.size();
330 return numNeighborScvfs;
334 void makeElementGeometries(
const Element& element)
337 const auto eIdx = gridGeometry().elementMapper().index(element);
338 scvs_[0] = SubControlVolume(element.geometry(), eIdx);
339 scvIndices_[0] = eIdx;
342 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdx);
343 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdx);
346 static const auto q = getParam<CoordScalar>(
"MPFA.Q");
349 const auto numLocalScvf = scvFaceIndices.size();
350 scvfIndices_.reserve(numLocalScvf);
351 scvfs_.reserve(numLocalScvf);
355 std::vector<bool> finishedFacets;
357 finishedFacets.resize(element.subEntities(1),
false);
360 for (
const auto& is : intersections(gridGeometry().gridView(), element))
366 const auto indexInInside = is.indexInInside();
367 if (finishedFacets[indexInInside])
370 finishedFacets[indexInInside] =
true;
374 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
375 const auto& e = useNeighbor ? is.outside() : element;
376 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
377 const auto eg = e.geometry();
378 const auto refElement = referenceElement(eg);
381 const auto numCorners = is.geometry().corners();
382 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
388 for (
int c = 0; c < numCorners; ++c)
391 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
392 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
395 if (gridGeometry().isGhostVertex(vIdxGlobal))
398 hasBoundaryScvf_ = (hasBoundaryScvf_ || is.boundary());
400 scvfs_.emplace_back(MpfaHelper(),
401 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
402 is.centerUnitOuterNormal(),
405 scvFaceIndices[scvfCounter],
407 neighborVolVarIndices[scvfCounter],
411 scvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
418 template<
typename IndexVector>
419 void makeNeighborGeometries(
const Element& element,
420 GridIndexType eIdxGlobal,
421 const IndexVector& scvfIndices,
422 const IndexVector& additionalScvfs)
425 neighborScvs_.emplace_back(element.geometry(), eIdxGlobal);
426 neighborScvIndices_.push_back(eIdxGlobal);
429 const auto& scvFaceIndices = gridGeometry().scvfIndicesOfScv(eIdxGlobal);
430 const auto& neighborVolVarIndices = gridGeometry().neighborVolVarIndices(eIdxGlobal);
433 static const auto q = getParam<CoordScalar>(
"MPFA.Q");
437 std::vector<bool> finishedFacets;
439 finishedFacets.resize(element.subEntities(1),
false);
442 for (
const auto& is : intersections(gridGeometry().gridView(), element))
448 auto indexInInside = is.indexInInside();
449 if(finishedFacets[indexInInside])
452 finishedFacets[indexInInside] =
true;
456 bool useNeighbor = is.neighbor() && is.outside().level() > element.level();
457 const auto& e = useNeighbor ? is.outside() : element;
458 const auto indexInElement = useNeighbor ? is.indexInOutside() : is.indexInInside();
459 const auto eg = e.geometry();
460 const auto refElement = referenceElement(eg);
463 const auto numCorners = is.geometry().corners();
464 const auto isPositions = MpfaHelper::computeScvfCornersOnIntersection(eg,
470 for (
int c = 0; c < numCorners; ++c)
473 auto vIdxLocal = refElement.subEntity(indexInElement, 1, c, dim);
474 auto vIdxGlobal = gridGeometry().vertexMapper().subIndex(e, vIdxLocal, dim);
477 if (gridGeometry().isGhostVertex(vIdxGlobal))
481 if (!MpfaHelper::vectorContainsValue(scvfIndices, scvFaceIndices[scvfCounter])
482 && !MpfaHelper::vectorContainsValue(additionalScvfs, scvFaceIndices[scvfCounter]))
490 neighborScvfs_.emplace_back(MpfaHelper(),
491 MpfaHelper::getScvfCorners(isPositions, numCorners, c),
492 is.centerUnitOuterNormal(),
495 scvFaceIndices[scvfCounter],
497 neighborVolVarIndices[scvfCounter],
501 neighborScvfIndices_.emplace_back(scvFaceIndices[scvfCounter]);
510 unsigned int findLocalIndex(
const GridIndexType idx,
511 const std::vector<GridIndexType>& indices)
const
513 auto it = std::find(indices.begin(), indices.end(), idx);
514 assert(it != indices.end() &&
"Could not find the scv/scvf! Make sure to properly bind this class!");
521 scvfIndices_.clear();
524 neighborScvIndices_.clear();
525 neighborScvfIndices_.clear();
526 neighborScvs_.clear();
527 neighborScvfs_.clear();
529 hasBoundaryScvf_ =
false;
532 const GridGeometry* gridGeometryPtr_;
535 std::array<GridIndexType, 1> scvIndices_;
536 std::vector<GridIndexType> scvfIndices_;
537 std::array<SubControlVolume, 1> scvs_;
538 std::vector<SubControlVolumeFace> scvfs_;
540 std::vector<GridIndexType> neighborScvIndices_;
541 std::vector<GridIndexType> neighborScvfIndices_;
542 std::vector<SubControlVolume> neighborScvs_;
543 std::vector<SubControlVolumeFace> neighborScvfs_;
545 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.
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
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models This builds up th...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:48
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered mpfa models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:58
const SubControlVolume & scv(GridIndexType scvIdx) const
Get an element sub control volume with a global scv index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:84
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:137
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:80
void bindElement(const Element &element)
Bind only element-local.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:149
friend Dune::IteratorRange< ScvIterator< SubControlVolume, std::array< GridIndexType, 1 >, ThisType > > scvs(const CCMpfaFVElementGeometry &fvGeometry)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:108
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:159
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:121
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get an element sub control volume face with a global scvf index.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:90
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:97
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:155
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:131
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:67
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:143
GG GridGeometry
export type of finite volume grid geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:73
Stencil-local finite volume geometry (scvs and scvfs) for cell-centered TPFA models Specialization fo...
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:175
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:215
const SubControlVolume & scv(GridIndexType scvIdx) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:205
friend Dune::IteratorRange< typename std::vector< SubControlVolumeFace >::const_iterator > scvfs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:249
void bind(const Element &element)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:265
const GridGeometry & gridGeometry() const
The global finite volume geometry we are a restriction of.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:314
void bindElement(const Element &element)
Binding of an element preparing the geometries only inside the element.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:307
const SubControlVolumeFace & flipScvf(GridIndexType scvfIdx, unsigned int outsideScvIdx=0) const
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:226
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:189
typename GridView::template Codim< 0 >::Entity Element
export type of the element
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:187
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:256
typename GG::SubControlVolumeFace SubControlVolumeFace
export type of subcontrol volume face
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:191
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:260
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:318
GG GridGeometry
export type of finite volume grid geometrys
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:193
friend Dune::IteratorRange< typename std::array< SubControlVolume, 1 >::const_iterator > scvs(const ThisType &g)
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:237
CCMpfaFVElementGeometry(const GridGeometry &gridGeometry)
Constructor.
Definition: discretization/cellcentered/mpfa/fvelementgeometry.hh:200
Iterators over sub control volumes.
Definition: scvandscvfiterators.hh:42
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82