24#ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_GEOMETRY_HELPER_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_GEOMETRY_HELPER_HH
29#include <dune/common/exceptions.hh>
38template<
class Gr
idView>
43 using Element =
typename GridView::template Codim<0>::Entity;
44 using Facet =
typename GridView::template Codim<1>::Entity;
47 static constexpr auto dim = GridView::Grid::dimension;
54 static constexpr SmallLocalIndexType
localOppositeIdx(
const SmallLocalIndexType ownLocalFaceIndex)
56 return isOdd_(ownLocalFaceIndex) ? (ownLocalFaceIndex - 1) : (ownLocalFaceIndex + 1);
62 if constexpr(GridView::Grid::dimension == 1)
64 assert(
false &&
"There are no lateral scvfs in 1D");
68 if constexpr (GridView::Grid::dimension == 2)
70 switch (ownLocalScvfIndex)
82 assert(
false &&
"No lateral orthogonal scvf found");
89 switch (ownLocalScvfIndex)
118 assert(
false &&
"No lateral orthogonal scvf found");
128 constexpr auto table = []
130 using Table = std::array<std::array<SmallLocalIndexType, numLateralFacesPerScv>,
numElementFaces>;
131 if constexpr (
dim == 1)
133 else if constexpr (
dim == 2)
134 return Table {{ {2,3}, {2,3}, {0,1}, {0,1} }};
136 return Table {{ {2,3,4,5}, {2,3,4,5}, {0,1,4,5}, {0,1,4,5}, {0,1,2,3}, {0,1,2,3} }};
139 return table[ownLocalFaceIndex];
143 static Facet
facet(
const SmallLocalIndexType localFacetIdx,
const Element& element)
145 return element.template subEntity <1> (localFacetIdx);
149 auto intersection(
const SmallLocalIndexType localFacetIdx,
const Element& element)
const
156 DUNE_THROW(Dune::InvalidStateException,
"localFacetIdx " << localFacetIdx <<
" out of range");
160 template<
class FVElementGeometry,
class SubControlVolumeFace>
163 assert (scvf.isLateral());
165 using Grid =
typename FVElementGeometry::GridGeometry::Grid;
173 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
174 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
175 int onBoundaryCounter = 0;
176 onBoundaryCounter +=
static_cast<int>(insideScv.boundary());
177 onBoundaryCounter +=
static_cast<int>(outsideScv.boundary());
178 return onBoundaryCounter == 1;
182 template<
class SubControlVolumeFace,
class SubControlVolume>
184 const SubControlVolume& scv)
188 const SmallLocalIndexType offset = (
dim == 2) ? 3 : 5;
191 return isOdd_(scv.indexInElement()) ? scvf.localIndex() - offset : scvf.localIndex() + offset;
195 {
return gridView_; }
199 static constexpr bool isOdd_(
int number)
200 {
return number % 2; }
Defines the index types used for grid and local indices.
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
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
Face centered staggered geometry helper.
Definition: discretization/facecentered/staggered/geometryhelper.hh:40
static Facet facet(const SmallLocalIndexType localFacetIdx, const Element &element)
Returns an element's facet based on the local facet index.
Definition: discretization/facecentered/staggered/geometryhelper.hh:143
static constexpr auto localLaterFaceIndices(const SmallLocalIndexType ownLocalFaceIndex)
Returns the local indices of the faces lateral to the own one.
Definition: discretization/facecentered/staggered/geometryhelper.hh:126
static bool scvfIntegrationPointInConcaveCorner(const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Returns true if the IP of an scvf lies on a concave corner.
Definition: discretization/facecentered/staggered/geometryhelper.hh:161
static constexpr SmallLocalIndexType localOppositeIdx(const SmallLocalIndexType ownLocalFaceIndex)
Returns the local index of the opposing face.
Definition: discretization/facecentered/staggered/geometryhelper.hh:54
static constexpr auto numLateralFacesPerScv
Definition: discretization/facecentered/staggered/geometryhelper.hh:49
const GridView & gridView() const
Definition: discretization/facecentered/staggered/geometryhelper.hh:194
static constexpr auto numElementFaces
Definition: discretization/facecentered/staggered/geometryhelper.hh:48
FaceCenteredStaggeredGeometryHelper(const GridView &gridView)
Definition: discretization/facecentered/staggered/geometryhelper.hh:51
static SmallLocalIndexType localIndexOutsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/geometryhelper.hh:183
auto intersection(const SmallLocalIndexType localFacetIdx, const Element &element) const
Returns an element's intersection based on the local facet index.
Definition: discretization/facecentered/staggered/geometryhelper.hh:149
static constexpr int lateralOrthogonalScvfLocalIndex(const SmallLocalIndexType ownLocalScvfIndex)
Return the local index of a lateral orthogonal scvf.
Definition: discretization/facecentered/staggered/geometryhelper.hh:60
static constexpr auto dim
Definition: discretization/facecentered/staggered/geometryhelper.hh:47
Type trait to determine if a grid supports concave corners (e.g. by cutting out a hole from the domai...
Definition: gridsupportsconcavecorners.hh:43