12#ifndef DUMUX_DISCRETIZATION_STAGGERED_GEOMETRY_HELPER_HH
13#define DUMUX_DISCRETIZATION_STAGGERED_GEOMETRY_HELPER_HH
15#include <dune/geometry/multilineargeometry.hh>
58template<
class Gr
idView,
int upwindSchemeOrder>
61 using Scalar =
typename GridView::ctype;
62 using GlobalPosition =
typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate;
82template<
class Gr
idView,
int upwindSchemeOrder>
85 using Scalar =
typename GridView::ctype;
103template<
class Gr
idView>
120template<
class Vector>
123 const auto eps = 1e-8;
124 return std::find_if(vector.begin(), vector.end(), [eps](
const auto& x) { return std::abs(x) > eps; } ) - vector.begin();
132template<
class Gr
idView,
int upwindSchemeOrder>
135 using Scalar =
typename GridView::ctype;
136 static constexpr auto dim = GridView::dimension;
138 using Element =
typename GridView::template Codim<0>::Entity;
139 using Intersection =
typename GridView::Intersection;
140 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
146 static constexpr auto codimIntersection = 1;
147 static constexpr auto numFacets = dim * 2;
148 static constexpr auto numPairs = 2 * (dim - 1);
150 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
158 , gridView_(gridView)
162 template<
class IntersectionMapper>
165 intersection_ = intersection;
175 return intersection_.indexInInside();
218 if constexpr (useHigherOrder)
219 fillOuterAxisData_();
221 const auto inIdx = intersection_.indexInInside();
222 const auto oppIdx = localOppositeIdx_(inIdx);
225 axisData_.
selfDof = gridView_.indexSet().subIndex(intersection_.inside(), inIdx, codimIntersection);
228 axisData_.
oppositeDof = gridView_.indexSet().subIndex(intersection_.inside(), oppIdx, codimIntersection);
231 const auto self = getFacet_(inIdx, element_);
232 const auto opposite = getFacet_(oppIdx, element_);
242 void fillOuterAxisData_()
248 const auto thisFaceLocalIdx = intersection_.indexInInside();
249 const auto& faceCenterPos = getFacet_(thisFaceLocalIdx, element_).geometry().center();
250 if (intersection_.neighbor())
251 addForwardNeighborAxisData_(intersection_.outside(), 0, thisFaceLocalIdx, faceCenterPos);
254 const auto oppositeFaceLocalIdx = localOppositeIdx_(thisFaceLocalIdx);
255 const auto& oppositeFaceCenterPos = getFacet_(oppositeFaceLocalIdx, element_).geometry().center();
256 for (
const auto& intersection : intersections(gridView_, element_))
258 if (intersection.indexInInside() == oppositeFaceLocalIdx && intersection.neighbor())
260 addBackwardNeighborAxisData_(intersection.outside(), 0, oppositeFaceLocalIdx, oppositeFaceCenterPos);
269 void addForwardNeighborAxisData_(
const Element& neighbor,
271 const unsigned int localIntersectionIndex,
272 const GlobalPosition& faceCenterPos)
276 axisData_.
inAxisForwardDofs[
distance] = gridView_.indexSet().subIndex(neighbor, localIntersectionIndex, codimIntersection);
277 const auto& forwardFacePos = getFacet_(localIntersectionIndex, neighbor).geometry().center();
282 if (
distance >= numForwardBackwardAxisDofs-1)
285 for (
const auto& intersection : intersections(gridView_, neighbor))
287 if (intersection.indexInInside() == localIntersectionIndex && intersection.neighbor())
289 addForwardNeighborAxisData_(intersection.outside(),
distance+1, localIntersectionIndex, forwardFacePos);
298 void addBackwardNeighborAxisData_(
const Element& neighbor,
300 const unsigned int localIntersectionIndex,
301 const GlobalPosition& faceCenterPos)
306 const auto& backwardFacePos = getFacet_(localIntersectionIndex, neighbor).geometry().center();
311 if (
distance >= numForwardBackwardAxisDofs-1)
314 for (
const auto& intersection : intersections(gridView_, neighbor))
316 if (intersection.indexInInside() == localIntersectionIndex && intersection.neighbor())
318 addBackwardNeighborAxisData_(intersection.outside(),
distance+1, localIntersectionIndex, backwardFacePos);
334 const auto& selfElementCenter = element_.geometry().center();
335 const auto& selfFacetCenter = intersection_.geometry().center();
338 SmallLocalIndexType numPairInnerLateralIdx = 0;
339 for (
const auto& innerElementIntersection : intersections(gridView_, element_))
341 if (facetIsNormal_(innerElementIntersection.indexInInside(), intersection_.indexInInside()))
343 const auto innerElementIntersectionIdx = innerElementIntersection.indexInInside();
344 setLateralPairFirstInfo_(innerElementIntersectionIdx, element_, numPairInnerLateralIdx);
346 const auto distance = innerElementIntersection.geometry().center() - selfElementCenter;
347 pairData_[numPairInnerLateralIdx].lateralStaggeredFaceCenter = selfFacetCenter+
distance;
349 numPairInnerLateralIdx++;
354 SmallLocalIndexType numPairOuterLateralIdx = 0;
355 if (intersection_.neighbor())
358 const auto& directNeighborElement = intersection_.outside();
360 for (
const auto& directNeighborElementIntersection : intersections(gridView_, directNeighborElement))
363 if (facetIsNormal_(directNeighborElementIntersection.indexInInside(), intersection_.indexInOutside()))
365 const auto directNeighborElemIsIdx = directNeighborElementIntersection.indexInInside();
366 setLateralPairSecondInfo_(directNeighborElemIsIdx, directNeighborElement, numPairOuterLateralIdx);
367 numPairOuterLateralIdx++;
373 for (
const auto& intersection : intersections(gridView_, element_))
375 if (facetIsNormal_(intersection.indexInInside(), intersection_.indexInInside()))
377 assert(!pairData_[numPairOuterLateralIdx].hasOuterLateral);
379 const auto normalDistanceoffset = selfFacetCenter - selfElementCenter;
380 pairData_[numPairOuterLateralIdx].lateralDistance = normalDistanceoffset.two_norm();
381 numPairOuterLateralIdx++;
387 const auto parallelLocalIdx = intersection_.indexInInside();
388 SmallLocalIndexType numPairParallelIdx = 0;
389 for (
const auto& intersection : intersections(gridView_, element_))
391 if (facetIsNormal_(intersection.indexInInside(), parallelLocalIdx))
393 if (intersection.neighbor())
397 const auto localLateralIntersectionIndex = intersection.indexInInside();
412 if (intersection_.neighbor())
413 for (
const auto& outerIntersection : intersections(gridView_, intersection_.outside()))
414 if (intersection.indexInInside() == outerIntersection.indexInInside())
415 if (!outerIntersection.neighbor())
416 pairData_[numPairParallelIdx].hasHalfParallelNeighbor =
true;
430 if (!intersection_.neighbor())
431 for (
const auto& outerIntersection : intersections(gridView_, intersection.outside()))
432 if (intersection_.indexInInside() == outerIntersection.indexInInside())
433 if (outerIntersection.neighbor())
434 pairData_[numPairParallelIdx].hasCornerParallelNeighbor =
true;
437 addParallelNeighborPairData_(intersection.outside(), 0, localLateralIntersectionIndex, parallelLocalIdx, parallelAxisIdx, numPairParallelIdx);
440 ++numPairParallelIdx;
446 void addParallelNeighborPairData_(
const Element& neighbor,
448 const unsigned int localLateralIntersectionIndex,
449 const unsigned int parallelLocalIdx,
450 const unsigned int parallelAxisIdx,
451 const SmallLocalIndexType dataIdx)
454 pairData_[dataIdx].hasParallelNeighbor.set(
distance,
true);
455 pairData_[dataIdx].parallelDofs[
distance] = gridView_.indexSet().subIndex(neighbor, parallelLocalIdx, codimIntersection);
456 pairData_[dataIdx].parallelCellWidths[
distance] = setParallelPairCellWidths_(neighbor, parallelAxisIdx);
459 const auto numParallelFaces = pairData_[0].parallelCellWidths.size();
464 for (
const auto& lateralIntersection : intersections(gridView_, neighbor))
467 if (lateralIntersection.indexInInside() == localLateralIntersectionIndex)
470 if (lateralIntersection.neighbor())
471 addParallelNeighborPairData_(lateralIntersection.outside(),
distance+1,
472 localLateralIntersectionIndex, parallelLocalIdx, parallelAxisIdx, dataIdx);
483 int localOppositeIdx_(
const int idx)
const
485 return (idx % 2) ? (idx - 1) : (idx + 1);
494 bool facetIsNormal_(
const int selfIdx,
const int otherIdx)
const
496 return !(selfIdx == otherIdx || localOppositeIdx_(selfIdx) == otherIdx);
499 auto getFacet_(
const int localFacetIdx,
const Element& element)
const
501 return element.template subEntity <1> (localFacetIdx);
505 void setLateralPairFirstInfo_(
const int isIdx,
const Element& element,
const int numPairsIdx)
508 pairData_[numPairsIdx].lateralPair.first = gridView_.indexSet().subIndex(element, isIdx, codimIntersection);
511 pairData_[numPairsIdx].localLateralFaceIdx = isIdx;
515 void setLateralPairSecondInfo_(
const int isIdx,
const Element& element,
const int numPairsIdx)
518 pairData_[numPairsIdx].hasOuterLateral =
true;
519 pairData_[numPairsIdx].lateralPair.second = gridView_.indexSet().subIndex(element, isIdx, codimIntersection);
522 const auto& selfFacetCenter = intersection_.geometry().center();
523 const auto& selfElementCenter = element_.geometry().center();
525 const auto& neighborElement = intersection_.outside();
526 const auto& neighborElementCenter = neighborElement.geometry().center();
527 const auto& neighborFacetCenter = getFacet_(intersection_.indexInOutside(), neighborElement).geometry().center();
529 const Scalar insideLateralDistance = (selfFacetCenter - selfElementCenter).two_norm();
530 const Scalar outsideLateralDistance = (neighborFacetCenter - neighborElementCenter).two_norm();
532 pairData_[numPairsIdx].lateralDistance = insideLateralDistance + outsideLateralDistance;
536 Scalar setParallelPairCellWidths_(
const Element& element,
const int parallelAxisIdx)
const
538 std::vector<GlobalPosition> faces;
539 faces.reserve(numFacets);
540 for (
const auto& intersection : intersections(gridView_, element))
541 faces.push_back(intersection.geometry().center());
543 switch (parallelAxisIdx)
546 return (faces[1] - faces[0]).two_norm();
548 return (faces[3] - faces[2]).two_norm();
552 return (faces[5] - faces[4]).two_norm();
555 DUNE_THROW(Dune::InvalidStateException,
"Something went terribly wrong");
559 Intersection intersection_;
560 const Element& element_;
561 const GridView gridView_;
563 std::array<PairData, numPairs> pairData_;
Helper class constructing the dual grid finite volume geometries for the free flow staggered discreti...
Definition: staggeredgeometryhelper.hh:134
Detail::AxisData< GridView, upwindSchemeOrder > AxisData
Definition: staggeredgeometryhelper.hh:154
SmallLocalIndexType localFaceIndex() const
Returns the local index of the face (i.e. the intersection)
Definition: staggeredgeometryhelper.hh:173
FreeFlowStaggeredGeometryHelper(const Element &element, const GridView &gridView)
Definition: staggeredgeometryhelper.hh:156
AxisData axisData() const
Returns a copy of the axis data.
Definition: staggeredgeometryhelper.hh:181
unsigned int directionIndex() const
Returns the direction index of the primary facet (0 = x, 1 = y, 2 = z)
Definition: staggeredgeometryhelper.hh:197
void updateLocalFace(const IntersectionMapper &, const Intersection &intersection)
update the local face
Definition: staggeredgeometryhelper.hh:163
std::array< PairData, numPairs > pairData() const
Returns a copy of the pair data.
Definition: staggeredgeometryhelper.hh:189
unsigned int directionIndex(const Intersection &intersection) const
Returns the direction index of the facet passed as an argument (0 = x, 1 = y, 2 = z)
Definition: staggeredgeometryhelper.hh:205
defines an intersection mapper for mapping of global DOFs assigned to faces which also works for adap...
Definition: intersectionmapper.hh:200
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:282
static unsigned int directionIndex(Vector &&vector)
Returns the direction index of the facet (0 = x, 1 = y, 2 = z)
Definition: staggeredgeometryhelper.hh:121
Define some often used mathematical functions.
GridIndexType oppositeDof
Definition: staggeredgeometryhelper.hh:110
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: staggeredgeometryhelper.hh:107
Scalar selfToOppositeDistance
Definition: staggeredgeometryhelper.hh:111
typename GridView::ctype Scalar
Definition: staggeredgeometryhelper.hh:106
GridIndexType selfDof
Definition: staggeredgeometryhelper.hh:109
In Axis Data stored per sub face.
Definition: staggeredgeometryhelper.hh:84
Scalar selfToOppositeDistance
Definition: staggeredgeometryhelper.hh:94
GridIndexType selfDof
Definition: staggeredgeometryhelper.hh:88
std::array< Scalar, upwindSchemeOrder-1 > inAxisForwardDistances
Definition: staggeredgeometryhelper.hh:95
std::bitset< upwindSchemeOrder-1 > hasBackwardNeighbor
Definition: staggeredgeometryhelper.hh:91
std::array< GridIndexType, upwindSchemeOrder-1 > inAxisForwardDofs
Definition: staggeredgeometryhelper.hh:92
std::array< Scalar, upwindSchemeOrder-1 > inAxisBackwardDistances
Definition: staggeredgeometryhelper.hh:96
std::bitset< upwindSchemeOrder-1 > hasForwardNeighbor
Definition: staggeredgeometryhelper.hh:90
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: staggeredgeometryhelper.hh:86
typename GridView::ctype Scalar
Definition: staggeredgeometryhelper.hh:85
GridIndexType oppositeDof
Definition: staggeredgeometryhelper.hh:89
std::array< GridIndexType, upwindSchemeOrder-1 > inAxisBackwardDofs
Definition: staggeredgeometryhelper.hh:93
Parallel Data stored per sub face.
Definition: staggeredgeometryhelper.hh:60
std::pair< GridIndexType, GridIndexType > lateralPair
Definition: staggeredgeometryhelper.hh:72
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: staggeredgeometryhelper.hh:63
bool hasOuterLateral
Definition: staggeredgeometryhelper.hh:71
typename GridView::ctype Scalar
Definition: staggeredgeometryhelper.hh:61
GlobalPosition lateralStaggeredFaceCenter
Definition: staggeredgeometryhelper.hh:75
std::bitset< upwindSchemeOrder > hasParallelNeighbor
Definition: staggeredgeometryhelper.hh:66
bool hasHalfParallelNeighbor
Definition: staggeredgeometryhelper.hh:67
bool hasCornerParallelNeighbor
Definition: staggeredgeometryhelper.hh:68
SmallLocalIndexType localLateralFaceIdx
Definition: staggeredgeometryhelper.hh:73
typename IndexTraits< GridView >::SmallLocalIndex SmallLocalIndexType
Definition: staggeredgeometryhelper.hh:64
std::array< Scalar, upwindSchemeOrder > parallelCellWidths
Definition: staggeredgeometryhelper.hh:70
typename GridView::template Codim< 0 >::Entity::Geometry::GlobalCoordinate GlobalPosition
Definition: staggeredgeometryhelper.hh:62
std::array< GridIndexType, upwindSchemeOrder > parallelDofs
Definition: staggeredgeometryhelper.hh:69
Scalar lateralDistance
Definition: staggeredgeometryhelper.hh:74
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:29