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");
226 const auto& scvfIndices_()
const
228 return gridGeometry().scvfIndicesOfElement(eIdx_);
231 const Element* elementPtr_;
233 const GridGeometry* gridGeometry_;
246 using GridView =
typename GG::GridView;
251 using StaticInfo =
typename GG::StaticInformation;
252 static constexpr auto dim = StaticInfo::dim;
253 static constexpr auto numScvsPerElement = StaticInfo::numScvsPerElement;
254 static constexpr auto numLateralScvfsPerScv = StaticInfo::numLateralScvfsPerScv;
255 static constexpr auto numLateralScvfsPerElement = StaticInfo::numLateralScvfsPerElement;
256 static constexpr auto minNumScvfsPerElement = StaticInfo::minNumScvfsPerElement;
257 static constexpr auto maxNumScvfsPerElement = StaticInfo::maxNumScvfsPerElement;
265 using Element =
typename GridView::template Codim<0>::Entity;
269 static constexpr std::size_t maxNumElementScvs = numScvsPerElement;
272 : gridGeometry_(&gridGeometry)
273 , geometryHelper_(gridGeometry.gridView())
278 {
return scvfs_[findLocalIndex_(scvfIdx, scvfIndices_())]; }
288 using Iter =
typename decltype(fvGeometry.scvfs_)::const_iterator;
289 return Dune::IteratorRange<Iter>(fvGeometry.scvfs_.begin(), fvGeometry.scvfs_.end());
300 using IndexContainerType = std::decay_t<
decltype(fvGeometry.scvfIndices_())>;
302 auto begin = ScvfIterator::makeBegin(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
303 auto end = ScvfIterator::makeEnd(fvGeometry.scvfIndices_(), fvGeometry, scv.index());
304 return Dune::IteratorRange<ScvfIterator>(begin, end);
310 if (scvIdx >= scvIndicesOfElement_.front() && scvIdx <= scvIndicesOfElement_.back())
311 return scvs_[findLocalIndex_(scvIdx, scvIndicesOfElement_)];
313 return neighborScvs_[findLocalIndex_(scvIdx, neighborScvIndices_)];
324 using IteratorType =
typename std::array<SubControlVolume, 1>::const_iterator;
325 return Dune::IteratorRange<IteratorType>(g.scvs_.begin(), g.scvs_.end());
331 assert(scvf.isLateral());
332 const auto otherLocalIdx = GridGeometry::GeometryHelper::lateralOrthogonalScvfLocalIndex(scvf.localIndex());
333 return scvfs_[otherLocalIdx];
339 assert(scv.boundary());
342 auto scvfIter = scvfs(*
this, scv).begin();
343 const auto end = scvfs(*
this, scv).end();
344 while (!(scvfIter->isFrontal() && scvfIter->boundary()) && (scvfIter != end))
347 assert(scvfIter->isFrontal());
348 assert(scvfIter->boundary());
354 {
return gridGeometry().hasBoundaryScvf(eIdx_); }
358 {
return numScvsPerElement; }
362 {
return scvfIndices_().size(); }
371 this->bind_(element);
372 return std::move(*
this);
376 { this->bind_(element); }
385 typename GG::LocalIntersectionMapper localIsMapper;
386 localIsMapper.update(gridGeometry().gridView(), element);
387 this->bindElement_(element, localIsMapper);
388 return std::move(*
this);
393 typename GG::LocalIntersectionMapper localIsMapper;
394 localIsMapper.update(gridGeometry().gridView(), element);
395 this->bindElement_(element, localIsMapper);
403 {
return *elementPtr_; }
408 assert(gridGeometry_);
409 return *gridGeometry_;
414 {
return GG::GeometryHelper::scvfIntegrationPointInConcaveCorner(*
this, scvf); }
421 assert(!lateralOrthogonalScvf.boundary());
423 const auto otherLocalIdx = GG::GeometryHelper::localIndexOutsideScvfWithSameIntegrationPoint(scvf, scv(scvf.insideScvIdx()));
425 auto outsideFVGeometry =
localView(gridGeometry());
426 const auto outsideElementIdx = scv(lateralOrthogonalScvf.outsideScvIdx()).elementIndex();
427 outsideFVGeometry.bindElement(gridGeometry().element(outsideElementIdx));
429 for (
const auto& otherScvf : scvfs(outsideFVGeometry))
430 if (otherScvf.localIndex() == otherLocalIdx)
433 DUNE_THROW(Dune::InvalidStateException,
"No outside scvf found");
439 void bind_(
const Element& element)
441 typename GG::LocalIntersectionMapper localIsMapper;
442 localIsMapper.update(gridGeometry().gridView(), element);
444 bindElement_(element, localIsMapper);
445 neighborScvIndices_.clear();
446 neighborScvs_.clear();
448 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
450 if (intersection.neighbor())
452 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
453 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
454 const auto& neighborElement = intersection.outside();
455 const auto neighborElementIdx = gridGeometry().elementMapper().index(neighborElement);
456 const auto& neighborElementGeometry = neighborElement.geometry();
459 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
460 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), neighborElementIdx*numScvsPerElement);
462 typename GG::LocalIntersectionMapper localNeighborIsMapper;
463 localNeighborIsMapper.update(gridGeometry().gridView(), neighborElement);
465 for (
const auto& neighborIntersection : intersections(gridGeometry().gridView(), neighborElement))
467 const auto localNeighborScvIdx = localNeighborIsMapper.realToRefIdx(neighborIntersection.indexInInside());
468 if (localNeighborScvIdx != localScvIdx && localNeighborScvIdx != localOppositeScvIdx)
471 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(neighborElement, neighborIntersection.indexInInside());
472 neighborScvs_.push_back(SubControlVolume(
473 neighborElementGeometry,
474 neighborIntersection.geometry(),
475 globalScvIndicesOfNeighborElement[localNeighborScvIdx],
480 onDomainBoundary_(neighborIntersection)
483 neighborScvIndices_.push_back(globalScvIndicesOfNeighborElement[localNeighborScvIdx]);
491 void bindElement_(
const Element& element,
const typename GG::LocalIntersectionMapper& localIsMapper)
494 eIdx_ = gridGeometry().elementMapper().index(element);
495 std::iota(scvIndicesOfElement_.begin(), scvIndicesOfElement_.end(), eIdx_*numScvsPerElement);
497 scvfs_.resize(minNumScvfsPerElement);
499 const auto& elementGeometry =
element.geometry();
501 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
503 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
504 auto localScvfIdx = localScvIdx*(1 + numLateralScvfsPerScv);
505 const auto& intersectionGeometry = intersection.geometry();
506 const auto dofIndex = gridGeometry().intersectionMapper().globalIntersectionIndex(element, intersection.indexInInside());
508 scvs_[localScvIdx] = SubControlVolume(
510 intersectionGeometry,
511 scvIndicesOfElement_[localScvIdx],
516 onDomainBoundary_(intersection)
520 const auto localOppositeScvIdx = geometryHelper_.localOppositeIdx(localScvIdx);
521 scvfs_[localScvfIdx] = SubControlVolumeFace(
523 intersectionGeometry,
524 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localOppositeScvIdx]},
526 scvfIndices_()[localScvfIdx],
527 intersection.centerUnitOuterNormal(),
528 SubControlVolumeFace::FaceType::frontal,
529 SubControlVolumeFace::BoundaryType::interior
534 for (
const auto lateralFacetIndex : Dune::transformedRangeView(geometryHelper_.localLaterFaceIndices(localScvIdx),
535 [&](
auto&& idx) { return localIsMapper.refToRealIdx(idx) ;})
538 const auto& lateralIntersection = geometryHelper_.intersection(lateralFacetIndex, element);
539 const auto& lateralFacet = geometryHelper_.facet(lateralFacetIndex, element);
540 const auto& lateralFacetGeometry = lateralFacet.geometry();
543 const auto globalScvIndicesForLateralFace = [&]
545 const auto globalOutsideScvIdx = [&]
547 if (lateralIntersection.neighbor())
549 const auto parallelElemIdx = gridGeometry().elementMapper().index(lateralIntersection.outside());
551 std::array<GridIndexType, numScvsPerElement> globalScvIndicesOfNeighborElement;
552 std::iota(globalScvIndicesOfNeighborElement.begin(), globalScvIndicesOfNeighborElement.end(), parallelElemIdx*numScvsPerElement);
553 return globalScvIndicesOfNeighborElement[localScvIdx];
556 return gridGeometry().outsideVolVarIndex(scvfIndices_()[localScvfIdx]);
559 return std::array{scvIndicesOfElement_[localScvIdx], globalOutsideScvIdx};
562 const auto boundaryType = [&]
564 if (onProcessorBoundary_(lateralIntersection))
565 return SubControlVolumeFace::BoundaryType::processorBoundary;
566 else if (onDomainBoundary_(lateralIntersection))
567 return SubControlVolumeFace::BoundaryType::physicalBoundary;
569 return SubControlVolumeFace::BoundaryType::interior;
572 scvfs_[localScvfIdx] = SubControlVolumeFace(
574 intersectionGeometry,
575 lateralFacetGeometry,
576 globalScvIndicesForLateralFace,
578 scvfIndices_()[localScvfIdx],
579 lateralIntersection.centerUnitOuterNormal(),
580 SubControlVolumeFace::FaceType::lateral,
588 auto localScvfIdx = minNumScvfsPerElement;
589 for (
const auto& intersection : intersections(gridGeometry().gridView(), element))
592 if (onDomainBoundary_(intersection))
594 const auto localScvIdx = localIsMapper.realToRefIdx(intersection.indexInInside());
596 scvfs_.push_back(SubControlVolumeFace(
598 intersection.geometry(),
599 std::array{scvIndicesOfElement_[localScvIdx], scvIndicesOfElement_[localScvIdx]},
601 scvfIndices_()[localScvfIdx],
602 intersection.centerUnitOuterNormal(),
603 SubControlVolumeFace::FaceType::frontal,
604 SubControlVolumeFace::BoundaryType::physicalBoundary
610 if constexpr (!ConsistentlyOrientedGrid<typename GridView::Grid>{})
612 static const bool makeConsistentlyOriented = getParam<bool>(
"Grid.MakeConsistentlyOriented",
true);
613 if (makeConsistentlyOriented && scvfs_.size() > minNumScvfsPerElement)
616 std::sort(scvfs_.begin() + minNumScvfsPerElement, scvfs_.end(),
617 [](
const auto& scvfLeft,
const auto& scvfRight)
618 { return scvfLeft.insideScvIdx() < scvfRight.insideScvIdx(); }
624 const auto& scvfIndices_()
const
625 {
return gridGeometry().scvfIndicesOfElement(eIdx_); }
627 template<
class Entry,
class Container>
628 const LocalIndexType findLocalIndex_(
const Entry& entry,
629 const Container& container)
const
631 auto it = std::find(container.begin(), container.end(), entry);
632 assert(it != container.end() &&
"Could not find the entry! Make sure to properly bind this class!");
636 bool onDomainBoundary_(
const typename GridView::Intersection& intersection)
const
638 return !intersection.neighbor() && intersection.boundary();
641 bool onProcessorBoundary_(
const typename GridView::Intersection& intersection)
const
643 return !intersection.neighbor() && !intersection.boundary();
646 Dune::ReservedVector<SubControlVolumeFace, maxNumScvfsPerElement> scvfs_;
648 Dune::ReservedVector<SubControlVolume, numLateralScvfsPerElement> neighborScvs_;
649 Dune::ReservedVector<GridIndexType, numLateralScvfsPerElement> neighborScvIndices_;
651 std::array<SubControlVolume, numScvsPerElement> scvs_;
652 std::array<GridIndexType, numScvsPerElement> scvIndicesOfElement_;
654 const GridGeometry* gridGeometry_;
655 const Element* elementPtr_;
657 typename GridGeometry::GeometryHelper geometryHelper_;
Defines the index types used for grid and local indices.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Returns the normal axis index of a unit vector (0 = x, 1 = y, 2 = z)
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:292
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)
Definition: normalaxis.hh:39
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
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
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:243
std::size_t numScvf() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:361
bool hasBoundaryScvf() const
Returns whether one of the geometry's scvfs lies on a boundary.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:353
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:286
void bindElement(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:391
typename GG::SubControlVolume SubControlVolume
export type of subcontrol volume face
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:263
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:369
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:413
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:383
GG GridGeometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:266
GridIndexType elementIndex() const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:398
SubControlVolumeFace outsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf) const
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:418
typename GridView::template Codim< 0 >::Entity Element
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:265
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:329
const SubControlVolume & scv(const GridIndexType scvIdx) const
Get a half sub control volume with a global scv index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:308
const GridGeometry & gridGeometry() const
The grid geometry we are a restriction of.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:406
friend auto scvs(const FaceCenteredStaggeredFVElementGeometry &g)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:322
void bind(const Element &element) &
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:375
const SubControlVolumeFace & scvf(const GridIndexType scvfIdx) const
Get a sub control volume face with a global scvf index.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:277
const Element & element() const
The bound element.
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:402
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:337
friend auto scvfs(const FaceCenteredStaggeredFVElementGeometry &fvGeometry, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:298
FaceCenteredStaggeredFVElementGeometry(const GridGeometry &gridGeometry)
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:271
std::size_t numScv() const
number of sub control volumes in this fv element geometry
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:357
typename GG::SubControlVolumeFace SubControlVolumeFace
Definition: discretization/facecentered/staggered/fvelementgeometry.hh:264
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