24#ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_ELEMENT_GEOMETRY_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_FV_ELEMENT_GEOMETRY_HH
29#include <dune/common/rangeutilities.hh>
30#include <dune/common/reservedvector.hh>
33#include <dune/common/iteratorrange.hh>
41template<
class GG,
bool cachingEnabled>
53 using GridView =
typename GG::GridView;
55 static constexpr auto numScvsPerElement = GG::StaticInformation::numScvsPerElement;
61 using Element =
typename GridView::template Codim<0>::Entity;
65 static constexpr std::size_t maxNumElementScvs = numScvsPerElement;
68 : gridGeometry_(&gridGeometry)
73 {
return gridGeometry().scv(scvIdx); }
77 {
return gridGeometry().scvf(scvfIdx); }
82 assert(scvf.isLateral());
83 const auto otherGlobalIdx = scvfIndices_()[GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex())];
84 return gridGeometry().scvf(otherGlobalIdx);
91 assert(scv.boundary());
94 auto scvfIter = scvfs(*
this, scv).begin();
95 const auto end = scvfs(*
this, scv).end();
96 while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end))
99 assert(scvfIter->isFrontal());
100 assert(scvfIter->boundary());
111 {
return fvGeometry.gridGeometry().scvs(fvGeometry); }
121 using IndexContainerType = std::decay_t<
decltype(fvGeometry.scvfIndices_())>;
123 return Dune::IteratorRange<ScvfIterator>(
ScvfIterator(fvGeometry.scvfIndices_().begin(), fvGeometry),
124 ScvfIterator(fvGeometry.scvfIndices_().end(), fvGeometry));
135 using IndexContainerType = std::decay_t<
decltype(fvGeometry.scvfIndices_())>;
137 auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
138 auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
139 return Dune::IteratorRange<ScvfIterator>(begin, end);
144 {
return numScvsPerElement; }
148 {
return scvfIndices_().size(); }
152 {
return gridGeometry().hasBoundaryScvf(eIdx_); }
162 return std::move(*
this);
167 { this->bindElement(element); }
176 this->bindElement(element);
177 return std::move(*
this);
183 elementPtr_ = &element;
184 eIdx_ = gridGeometry().elementMapper().index(element);
190 assert(gridGeometry_);
191 return *gridGeometry_;
199 {
return *elementPtr_; }
203 {
return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*
this, scvf); }
208 const auto& lateralOrthogonalScvf = this->lateralOrthogonalScvf(scvf);
209 assert(!lateralOrthogonalScvf.boundary());
211 const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx()));
213 auto outsideFVGeometry =
localView(gridGeometry());
214 const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex();
215 outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx));
217 for (
const auto& otherScvf : scvfs(outsideFVGeometry))
218 if (otherScvf.localIndex() == otherLocalIdx)
221 DUNE_THROW(Dune::InvalidStateException,
"No outside scvf found");
225 {
return scv.geometry(); }
228 {
return scvf.geometry(); }
232 const auto& scvfIndices_()
const
234 return gridGeometry().scvfIndicesOfElement(eIdx_);
237 const Element* elementPtr_;
239 const GridGeometry* gridGeometry_;
252 using GridView =
typename GG::GridView;
257 using StaticInfo =
typename GG::StaticInformation;
258 static constexpr auto dim = StaticInfo::dim;
259 static constexpr auto numScvsPerElement = StaticInfo::numScvsPerElement;
260 static constexpr auto numLateralScvfsPerScv = StaticInfo::numLateralScvfsPerScv;
261 static constexpr auto numLateralScvfsPerElement = StaticInfo::numLateralScvfsPerElement;
262 static constexpr auto minNumScvfsPerElement = StaticInfo::minNumScvfsPerElement;
263 static constexpr auto maxNumScvfsPerElement = StaticInfo::maxNumScvfsPerElement;
271 using Element =
typename GridView::template Codim<0>::Entity;
275 static constexpr std::size_t maxNumElementScvs = numScvsPerElement;
278 : gridGeometry_(&gridGeometry)
279 , geometryHelper_(gridGeometry.gridView())
284 {
return scvfs_[findLocalIndex_(scvfIdx, scvfIndices_())]; }
294 using Iter =
typename decltype(fvGeometry.scvfs_)::const_iterator;
295 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
306 using IndexContainerType = std::decay_t<
decltype(fvGeometry.scvfIndices_())>;
308 auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
309 auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
310 return Dune::IteratorRange<ScvfIterator>(begin, end);
316 if (scvIdx >= scvIndicesOfElement_.front() && scvIdx <= scvIndicesOfElement_.back())
317 return scvs_[findLocalIndex_(scvIdx, scvIndicesOfElement_)];
319 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
330 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
331 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
337 assert(scvf.isLateral());
338 const auto otherLocalIdx = GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex());
339 return scvfs_[otherLocalIdx];
345 assert(scv.boundary());
348 auto scvfIter = scvfs(*
this, scv).begin();
349 const auto end = scvfs(*
this, scv).end();
350 while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end))
353 assert(scvfIter->isFrontal());
354 assert(scvfIter->boundary());
360 {
return gridGeometry().hasBoundaryScvf(eIdx_); }
364 {
return numScvsPerElement; }
368 {
return scvfIndices_().size(); }
377 this->bind_(element);
378 return std::move(*
this);
382 { this->bind_(element); }
391 typename GG::LocalIntersectionMapper localIsMapper;
392 localIsMapper.update(gridGeometry().gridView(), element);
393 this->bindElement_(element, localIsMapper);
394 return std::move(*
this);
399 typename GG::LocalIntersectionMapper localIsMapper;
400 localIsMapper.update(gridGeometry().gridView(), element);
401 this->bindElement_(element, localIsMapper);
409 {
return *elementPtr_; }
414 assert(gridGeometry_);
415 return *gridGeometry_;
420 {
return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*
this, scvf); }
427 assert(!lateralOrthogonalScvf.boundary());
429 const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx()));
431 auto outsideFVGeometry =
localView(gridGeometry());
432 const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex();
433 outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx));
435 for (
const auto& otherScvf : scvfs(outsideFVGeometry))
436 if (otherScvf.localIndex() == otherLocalIdx)
439 DUNE_THROW(Dune::InvalidStateException,
"No outside scvf found");
445 void bind_(
const Element& element)
447 typename GG::LocalIntersectionMapper localIsMapper;
448 localIsMapper.update(gridGeometry().gridView(), element);
450 bindElement_(element, localIsMapper);
451 neighborScvIndices_.clear();
452 neighborScvs_.clear();
454 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
456 if (intersection.neighbor())
458 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
459 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
460 const auto& neighborElement = intersection.outside();
461 const auto neighborElementIdx = gridGeometry().elementMapper().index(neighborElement);
462 const auto& neighborElementGeometry = neighborElement.geometry();
465 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
466 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), neighborElementIdx*numScvsPerElement);
468 typename GG::LocalIntersectionMapper localNeighborIsMapper;
469 localNeighborIsMapper.update(gridGeometry().gridView(), neighborElement);
471 for (
const auto& neighborIntersection : intersections(gridGeometry().gridView(), neighborElement))
473 const auto localNeighborScvIdx = localNeighborIsMapper.realToRefIdx(neighborIntersection.indexInInside());
474 if (localNeighborScvIdx != localScvIdx && localNeighborScvIdx != localOppositeScvIdx)
477 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(neighborElement, neighborIntersection.indexInInside());
478 neighborScvs_.push_back(SubControlVolume(
479 neighborElementGeometry,
480 neighborIntersection.geometry(),
481 globalScvIndicesOfNeighborElement[localNeighborScvIdx],
486 onDomainBoundary_(neighborIntersection)
489 neighborScvIndices_.push_back(globalScvIndicesOfNeighborElement[localNeighborScvIdx]);
497 void bindElement_(
const Element& element,
const typename GG::LocalIntersectionMapper& localIsMapper)
500 eIdx_ = gridGeometry().elementMapper().index(element);
501 std::iota(scvIndicesOfElement_.begin(), scvIndicesOfElement_.end(), eIdx_*numScvsPerElement);
503 scvfs_.resize(minNumScvfsPerElement);
505 const auto& elementGeometry =
element.geometry();
507 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
509 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
510 auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv);
511 const auto& intersectionGeometry = intersection.geometry();
512 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(element, intersection.indexInInside());
514 scvs_[localScvIdx] = SubControlVolume(
516 intersectionGeometry,
517 scvIndicesOfElement_[localScvIdx],
522 onDomainBoundary_(intersection)
526 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
527 scvfs_[localScvfIdx] = SubControlVolumeFace(
529 intersectionGeometry,
530 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localOppositeScvIdx]},
532 scvfIndices_()[localScvfIdx],
533 intersection.centerUnitOuterNormal(),
534 SubControlVolumeFace::FaceType::frontal,
535 SubControlVolumeFace::BoundaryType::interior
540 for (
const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper_.localLaterFaceIndices(localScvIdx),
541 [&](
auto&& idx) { return localIsMapper.refToRealIdx(idx) ;})
544 const auto& lateralIntersection = geometryHelper_.intersection(lateralFacetIndex, element);
545 const auto& lateralFacet = geometryHelper_.facet(lateralFacetIndex, element);
546 const auto& lateralFacetGeometry = lateralFacet.geometry();
549 const auto globalScvIndicesForLateralFace = [&]
551 const auto globalOutsideScvIdx = [&]
553 if (lateralIntersection.neighbor())
555 const auto parallelElemIdx = gridGeometry().elementMapper().index(lateralIntersection.outside());
557 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
558 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), parallelElemIdx*numScvsPerElement);
559 return globalScvIndicesOfNeighborElement[localScvIdx];
562 return gridGeometry().outsideVolVarIndex(scvfIndices_()[localScvfIdx]);
565 return std::array{scvIndicesOfElement_[localScvIdx], globalOutsideScvIdx};
568 const auto boundaryType = [&]
570 if (onProcessorBoundary_(lateralIntersection))
571 return SubControlVolumeFace::BoundaryType::processorBoundary;
572 else if (onDomainBoundary_(lateralIntersection))
573 return SubControlVolumeFace::BoundaryType::physicalBoundary;
575 return SubControlVolumeFace::BoundaryType::interior;
578 scvfs_[localScvfIdx] = SubControlVolumeFace(
580 intersectionGeometry,
581 lateralFacetGeometry,
582 globalScvIndicesForLateralFace,
584 scvfIndices_()[localScvfIdx],
585 lateralIntersection.centerUnitOuterNormal(),
586 SubControlVolumeFace::FaceType::lateral,
594 auto localScvfIdx = minNumScvfsPerElement;
595 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
598 if (onDomainBoundary_(intersection))
600 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
602 scvfs_.push_back(SubControlVolumeFace(
604 intersection.geometry(),
605 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localScvIdx]},
607 scvfIndices_()[localScvfIdx],
608 intersection.centerUnitOuterNormal(),
609 SubControlVolumeFace::FaceType::frontal,
610 SubControlVolumeFace::BoundaryType::physicalBoundary
616 if constexpr (!ConsistentlyOrientedGrid<typename GridView::Grid>{})
618 static const bool makeConsistentlyOriented = getParam<bool>(
"Grid.MakeConsistentlyOriented",
true);
619 if (makeConsistentlyOriented && scvfs_.size() > minNumScvfsPerElement)
622 std::sort(scvfs_.begin() + minNumScvfsPerElement, scvfs_.end(),
623 [](
const auto& scvfLeft,
const auto& scvfRight)
624 { return scvfLeft.insideScvIdx() < scvfRight.insideScvIdx(); }
630 const auto& scvfIndices_()
const
631 {
return gridGeometry().scvfIndicesOfElement(eIdx_); }
633 template<
class Entry,
class Container>
634 const LocalIndexType findLocalIndex_(
const Entry& entry,
635 const Container& container)
const
637 auto it = std::find(container.begin(), container.end(), entry);
638 assert(it != container.end() &&
"Could not find the entry! Make sure to properly bind this class!");
642 bool onDomainBoundary_(
const typename GridView::Intersection& intersection)
const
644 return !intersection.neighbor() && intersection.boundary();
647 bool onProcessorBoundary_(
const typename GridView::Intersection& intersection)
const
649 return !intersection.neighbor() && !intersection.boundary();
652 Dune::ReservedVector<SubControlVolumeFace, maxNumScvfsPerElement> scvfs_;
654 Dune::ReservedVector<SubControlVolume, numLateralScvfsPerElement> neighborScvs_;
655 Dune::ReservedVector<GridIndexType, numLateralScvfsPerElement> neighborScvIndices_;
657 std::array<SubControlVolume, numScvsPerElement> scvs_;
658 std::array<GridIndexType, numScvsPerElement> scvIndicesOfElement_;
660 const GridGeometry* gridGeometry_;
661 const Element* elementPtr_;
663 typename GridGeometry::GeometryHelper geometryHelper_;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
static std::size_t normalAxis(const Vector &v)
Returns the normal axis index of a unit vector (0 = x, 1 = y, 2 = z)
Definition: normalaxis.hh:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
unsigned int LocalIndex
Definition: indextraits.hh:40
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:42
Stencil-local finite volume geometry (scvs and scvfs) for face-centered staggered models Specializati...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:51
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:59
FaceCenteredStaggeredFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:159
FaceCenteredStaggeredFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:174
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:143
const Element & element() const
The bound element.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:198
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:188
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:119
const SubControlVolumeFace & scvf(GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:76
SubControlVolumeFace::Traits::Geometry geometry(const SubControlVolumeFace &scvf) const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:227
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:151
void bind(const Element &element) &
Binding of an element, called by the local jacobian to prepare element assembly.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:166
const SubControlVolumeFace & outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf) const
Returns the the scvf of neighbor element with the same integration point and unit outer normal.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:206
friend auto scvs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:110
const SubControlVolume & scv(GridIndexType scvIdx) const
Get a sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:72
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:60
FaceCenteredStaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:67
GG GridGeometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:62
bool scvfIntegrationPointInConcaveCorner(const SubControlVolumeFace &scvf) const
Returns true if the IP of an scvf lies on a concave corner.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:202
SubControlVolume::Traits::Geometry geometry(const SubControlVolume &scv) const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:224
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:61
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:133
const SubControlVolumeFace & lateralOrthogonalScvf(const SubControlVolumeFace &scvf) const
Return a the lateral sub control volume face which is orthogonal to the given sub control volume face...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:80
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:147
const SubControlVolumeFace & frontalScvfOnBoundary(const SubControlVolume &scv) const
Return the frontal sub control volume face on a the boundary for a given sub control volume.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:89
void bindElement(const Element &element) &
Bind only element-local.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:181
std::size_t elementIndex() const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:194
Stencil-local finite volume geometry (scvs and scvfs) for face-centered staggered models Specializati...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:249
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:367
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:359
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:292
void bindElement(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:397
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:269
FaceCenteredStaggeredFVElementGeometry bind(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:375
bool scvfIntegrationPointInConcaveCorner(const SubControlVolumeFace &scvf) const
Returns true if the IP of an scvf lies on a concave corner.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:419
FaceCenteredStaggeredFVElementGeometry bindElement(const Element &element) &&
bind the local view (r-value overload) This overload is called when an instance of this class is a te...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:389
GG GridGeometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:272
GridIndexType elementIndex() const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:404
SubControlVolumeFace outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf) const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:424
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:271
const SubControlVolumeFace & lateralOrthogonalScvf(const SubControlVolumeFace &scvf) const
Return a the lateral sub control volume face which is orthogonal to the given sub control volume face...
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:335
const SubControlVolume & scv(const GridIndexType scvIdx) const
Get a half sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:314
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:412
friend auto scvs(const FaceCenteredStaggeredFVElementGeometry &g)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:328
void bind(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:381
const SubControlVolumeFace & scvf(const GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:283
const Element & element() const
The bound element.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:408
const SubControlVolumeFace & frontalScvfOnBoundary(const SubControlVolume &scv) const
Return the frontal sub control volume face on a the boundary for a given sub control volume.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:343
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:304
FaceCenteredStaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:277
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:363
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:270
Iterators over sub control volume faces of an fv geometry.
Definition: scvandscvfiterators.hh:82
Iterators over sub control volume faces of an fv geometry and a given sub control volume.
Definition: scvandscvfiterators.hh:122