24#ifndef DUMUX_DISCRETIZATION_STAGGERED_GEOMETRY_HELPER_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_GEOMETRY_HELPER_HH
27#include <dune/geometry/multilineargeometry.hh>
46template<
class Gr
idView,
int upwindSchemeOrder>
49 using Scalar =
typename GridView::ctype;
50 using GlobalPosition =
typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate;
69template<
class Gr
idView,
int upwindSchemeOrder>
75template<
class Gr
idView,
int upwindSchemeOrder>
78 using Scalar =
typename GridView::ctype;
96template<
class Gr
idView>
99 using Scalar =
typename GridView::ctype;
113template<
class Vector>
116 const auto eps = 1e-8;
117 return std::find_if(vector.begin(), vector.end(), [eps](
const auto& x) { return std::abs(x) > eps; } ) - vector.begin();
125template<
class Gr
idView,
int upwindSchemeOrder>
128 using Scalar =
typename GridView::ctype;
129 static constexpr auto dim = GridView::dimension;
131 using Element =
typename GridView::template Codim<0>::Entity;
132 using Intersection =
typename GridView::Intersection;
133 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
139 static constexpr auto codimIntersection = 1;
140 static constexpr auto numfacets = dim * 2;
141 static constexpr auto numPairs = 2 * (dim - 1);
143 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
153 template<
class IntersectionMapper>
156 intersection_ = intersection;
166 return intersection_.indexInInside();
209 fillOuterAxisData_(std::integral_constant<bool, useHigherOrder>{});
211 const auto inIdx = intersection_.indexInInside();
212 const auto oppIdx = localOppositeIdx_(inIdx);
215 axisData_.
selfDof = gridView_.indexSet().subIndex(intersection_.inside(), inIdx, codimIntersection);
218 axisData_.
oppositeDof = gridView_.indexSet().subIndex(intersection_.inside(), oppIdx, codimIntersection);
221 const auto self = getFacet_(inIdx, element_);
222 const auto opposite = getFacet_(oppIdx, element_);
231 void fillOuterAxisData_(std::false_type) {}
237 void fillOuterAxisData_(std::true_type)
245 std::stack<Element> inAxisForwardElementStack;
246 const auto inIdx = intersection_.indexInInside();
247 auto selfFace = getFacet_(inIdx, element_);
248 if(intersection_.neighbor())
250 inAxisForwardElementStack.push(intersection_.outside());
251 bool keepStackingForward = (inAxisForwardElementStack.size() < numForwardBackwardAxisDofs);
252 while(keepStackingForward)
254 auto e = inAxisForwardElementStack.top();
255 for(
const auto& intersection : intersections(gridView_,e))
257 if( (intersection.indexInInside() == inIdx ) )
259 if( intersection.neighbor())
261 inAxisForwardElementStack.push(intersection.outside());
262 keepStackingForward = (inAxisForwardElementStack.size() < numForwardBackwardAxisDofs);
266 keepStackingForward =
false;
273 std::vector<GlobalPosition> forwardFaceCoordinates(inAxisForwardElementStack.size(), selfFace.geometry().center());
274 while(!inAxisForwardElementStack.empty())
276 int forwardIdx = inAxisForwardElementStack.size()-1;
278 axisData_.
inAxisForwardDofs[forwardIdx] = gridView_.indexSet().subIndex(inAxisForwardElementStack.top(), inIdx, codimIntersection);
279 selfFace = getFacet_(inIdx, inAxisForwardElementStack.top());
280 forwardFaceCoordinates[forwardIdx] = selfFace.geometry().center();
281 inAxisForwardElementStack.pop();
284 const auto self = getFacet_(inIdx, element_);
285 for (
int i = 0; i< forwardFaceCoordinates.size(); i++)
290 axisData_.
inAxisForwardDistances[i] = (forwardFaceCoordinates[i]- forwardFaceCoordinates[i-1]).two_norm();
294 std::stack<Element> inAxisBackwardElementStack;
295 const auto oppIdx = localOppositeIdx_(inIdx);
296 auto oppFace = getFacet_(oppIdx, element_);
297 for(
const auto& intersection : intersections(gridView_, element_))
299 if(intersection.indexInInside() == oppIdx)
301 if(intersection.neighbor())
303 inAxisBackwardElementStack.push(intersection.outside());
304 bool keepStackingBackward = (inAxisBackwardElementStack.size() < numForwardBackwardAxisDofs);
305 while(keepStackingBackward)
307 auto e = inAxisBackwardElementStack.top();
308 for(
const auto& intersectionOut : intersections(gridView_,e))
311 if( (intersectionOut.indexInInside() == oppIdx ) )
313 if( intersectionOut.neighbor())
315 inAxisBackwardElementStack.push(intersectionOut.outside());
316 keepStackingBackward = (inAxisBackwardElementStack.size() < numForwardBackwardAxisDofs);
320 keepStackingBackward =
false;
329 std::vector<GlobalPosition> backwardFaceCoordinates(inAxisBackwardElementStack.size(), oppFace.geometry().center());
330 while(!inAxisBackwardElementStack.empty())
332 int backwardIdx = inAxisBackwardElementStack.size()-1;
334 axisData_.
inAxisBackwardDofs[backwardIdx] = gridView_.indexSet().subIndex(inAxisBackwardElementStack.top(), oppIdx, codimIntersection);
335 oppFace = getFacet_(oppIdx, inAxisBackwardElementStack.top());
336 backwardFaceCoordinates[backwardIdx] = oppFace.geometry().center();
337 inAxisBackwardElementStack.pop();
340 const auto opposite = getFacet_(oppIdx, element_);
341 for (
int i = 0; i< backwardFaceCoordinates.size(); i++)
360 const auto& selfElementCenter = element_.geometry().center();
361 const auto& selfFacetCenter = intersection_.geometry().center();
364 SmallLocalIndexType numPairInnerLateralIdx = 0;
365 for (
const auto& innerElementIntersection : intersections(gridView_, element_))
367 if (facetIsNormal_(innerElementIntersection.indexInInside(), intersection_.indexInInside()))
369 const auto innerElementIntersectionIdx = innerElementIntersection.indexInInside();
370 setLateralPairFirstInfo_(innerElementIntersectionIdx, element_, numPairInnerLateralIdx);
372 const auto distance = innerElementIntersection.geometry().center() - selfElementCenter;
373 pairData_[numPairInnerLateralIdx].lateralStaggeredFaceCenter = std::move(selfFacetCenter + distance);
375 numPairInnerLateralIdx++;
380 SmallLocalIndexType numPairOuterLateralIdx = 0;
381 if (intersection_.neighbor())
384 const auto& directNeighborElement = intersection_.outside();
386 for (
const auto& directNeighborElementIntersection : intersections(gridView_, directNeighborElement))
389 if (facetIsNormal_(directNeighborElementIntersection.indexInInside(), intersection_.indexInOutside()))
391 const auto directNeighborElemIsIdx = directNeighborElementIntersection.indexInInside();
392 setLateralPairSecondInfo_(directNeighborElemIsIdx, directNeighborElement, numPairOuterLateralIdx);
393 numPairOuterLateralIdx++;
399 for (
const auto& intersection : intersections(gridView_, element_))
401 if (facetIsNormal_(intersection.indexInInside(), intersection_.indexInInside()))
403 assert(!pairData_[numPairOuterLateralIdx].hasOuterLateral);
405 const auto normalDistanceoffset = selfFacetCenter - selfElementCenter;
406 pairData_[numPairOuterLateralIdx].lateralDistance = std::move(normalDistanceoffset.two_norm());
407 numPairOuterLateralIdx++;
412 fillParallelPairData_(std::integral_constant<bool, useHigherOrder>{});
419 void fillParallelPairData_(std::false_type)
423 const auto parallelLocalIdx = intersection_.indexInInside();
424 SmallLocalIndexType numPairParallelIdx = 0;
426 for (
const auto& intersection : intersections(gridView_, element_))
428 if (facetIsNormal_(intersection.indexInInside(), parallelLocalIdx))
431 if (intersection.neighbor())
434 const auto& outerElement = intersection.outside();
435 pairData_[numPairParallelIdx].hasParallelNeighbor.set(0,
true);
436 pairData_[numPairParallelIdx].parallelDofs[0] = gridView_.indexSet().subIndex(outerElement, parallelLocalIdx, codimIntersection);
437 pairData_[numPairParallelIdx].parallelCellWidths[0] = setParallelPairCellWidths_(outerElement, parallelAxisIdx);
440 numPairParallelIdx++;
449 void fillParallelPairData_(std::true_type)
452 const auto numParallelFaces = pairData_[0].parallelCellWidths.size();
455 const auto parallelLocalIdx = intersection_.indexInInside();
456 SmallLocalIndexType numPairParallelIdx = 0;
457 std::stack<Element> parallelElementStack;
458 for(
const auto& intersection : intersections(gridView_, element_))
460 if( facetIsNormal_(intersection.indexInInside(), parallelLocalIdx) )
462 if( intersection.neighbor() )
465 auto localLateralIntersectionIndex = intersection.indexInInside();
468 bool keepStacking = (parallelElementStack.size() < numParallelFaces);
471 for(
const auto& lateralIntersection : intersections(gridView_, e))
473 if( lateralIntersection.indexInInside() == localLateralIntersectionIndex )
475 if( lateralIntersection.neighbor() )
477 parallelElementStack.push(lateralIntersection.outside());
478 keepStacking = (parallelElementStack.size() < numParallelFaces);
482 keepStacking =
false;
486 e = parallelElementStack.top();
489 while(!parallelElementStack.empty())
491 pairData_[numPairParallelIdx].hasParallelNeighbor.set(parallelElementStack.size()-1,
true);
492 pairData_[numPairParallelIdx].parallelDofs[parallelElementStack.size()-1] = gridView_.indexSet().subIndex(parallelElementStack.top(), parallelLocalIdx, codimIntersection);
493 pairData_[numPairParallelIdx].parallelCellWidths[parallelElementStack.size()-1] = setParallelPairCellWidths_(parallelElementStack.top(), parallelAxisIdx);
494 parallelElementStack.pop();
499 numPairParallelIdx++;
509 int localOppositeIdx_(
const int idx)
const
511 return (idx % 2) ? (idx - 1) : (idx + 1);
520 bool facetIsNormal_(
const int selfIdx,
const int otherIdx)
const
522 return !(selfIdx == otherIdx || localOppositeIdx_(selfIdx) == otherIdx);
525 auto getFacet_(
const int localFacetIdx,
const Element& element)
const
527 return element.template subEntity <1> (localFacetIdx);
531 void setLateralPairFirstInfo_(
const int isIdx,
const Element& element,
const int numPairsIdx)
534 pairData_[numPairsIdx].lateralPair.first = gridView_.indexSet().subIndex(element, isIdx, codimIntersection);
537 pairData_[numPairsIdx].localLateralFaceIdx = isIdx;
541 void setLateralPairSecondInfo_(
const int isIdx,
const Element& element,
const int numPairsIdx)
544 pairData_[numPairsIdx].hasOuterLateral =
true;
545 pairData_[numPairsIdx].lateralPair.second = gridView_.indexSet().subIndex(element, isIdx, codimIntersection);
548 const auto& selfFacetCenter = intersection_.geometry().center();
549 const auto& selfElementCenter = element_.geometry().center();
551 const auto& neighborElement = intersection_.outside();
552 const auto& neighborElementCenter = neighborElement.geometry().center();
553 const auto& neighborFacetCenter = getFacet_(intersection_.indexInOutside(), neighborElement).geometry().center();
555 const Scalar insideLateralDistance = (selfFacetCenter - selfElementCenter).two_norm();
556 const Scalar outsideLateralDistance = (neighborFacetCenter - neighborElementCenter).two_norm();
558 pairData_[numPairsIdx].lateralDistance = insideLateralDistance + outsideLateralDistance;
562 Scalar setParallelPairCellWidths_(
const Element& element,
const int parallelAxisIdx)
const
564 std::vector<GlobalPosition> faces;
565 faces.reserve(numfacets);
566 for (
const auto& intElement : intersections(gridView_, element))
568 faces.push_back(intElement.geometry().center());
570 switch(parallelAxisIdx)
573 return (faces[1] - faces[0]).two_norm();
575 return (faces[3] - faces[2]).two_norm();
579 return (faces[5] - faces[4]).two_norm();
582 DUNE_THROW(Dune::InvalidStateException,
"Something went terribly wrong");
586 Intersection intersection_;
587 const Element element_;
588 const typename Element::Geometry elementGeometry_;
589 const GridView gridView_;
591 std::array<PairData, numPairs> pairData_;
Defines the index types used for grid and local indices.
Define some often used mathematical functions.
static unsigned int directionIndex(Vector &&vector)
Returns the dirction index of the facet (0 = x, 1 = y, 2 = z)
Definition: staggeredgeometryhelper.hh:114
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
defines an intersection mapper for mapping of global DOFs assigned to faces which also works for adap...
Definition: intersectionmapper.hh:185
Parallel Data stored per sub face.
Definition: staggeredgeometryhelper.hh:48
std::pair< GridIndexType, GridIndexType > lateralPair
Definition: staggeredgeometryhelper.hh:58
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: staggeredgeometryhelper.hh:51
bool hasOuterLateral
Definition: staggeredgeometryhelper.hh:57
typename GridView::ctype Scalar
Definition: staggeredgeometryhelper.hh:49
GlobalPosition lateralStaggeredFaceCenter
Definition: staggeredgeometryhelper.hh:61
std::bitset< upwindSchemeOrder > hasParallelNeighbor
Definition: staggeredgeometryhelper.hh:54
SmallLocalIndexType localLateralFaceIdx
Definition: staggeredgeometryhelper.hh:59
typename IndexTraits< GridView >::SmallLocalIndex SmallLocalIndexType
Definition: staggeredgeometryhelper.hh:52
std::array< Scalar, upwindSchemeOrder > parallelCellWidths
Definition: staggeredgeometryhelper.hh:56
typename GridView::template Codim< 0 >::Entity::Geometry::GlobalCoordinate GlobalPosition
Definition: staggeredgeometryhelper.hh:50
std::array< GridIndexType, upwindSchemeOrder > parallelDofs
Definition: staggeredgeometryhelper.hh:55
Scalar lateralDistance
Definition: staggeredgeometryhelper.hh:60
In Axis Data stored per sub face.
Definition: staggeredgeometryhelper.hh:77
Scalar selfToOppositeDistance
Definition: staggeredgeometryhelper.hh:87
GridIndexType selfDof
Definition: staggeredgeometryhelper.hh:81
std::array< Scalar, upwindSchemeOrder-1 > inAxisForwardDistances
Definition: staggeredgeometryhelper.hh:88
std::bitset< upwindSchemeOrder-1 > hasBackwardNeighbor
Definition: staggeredgeometryhelper.hh:84
std::array< GridIndexType, upwindSchemeOrder-1 > inAxisForwardDofs
Definition: staggeredgeometryhelper.hh:85
std::array< Scalar, upwindSchemeOrder-1 > inAxisBackwardDistances
Definition: staggeredgeometryhelper.hh:89
std::bitset< upwindSchemeOrder-1 > hasForwardNeighbor
Definition: staggeredgeometryhelper.hh:83
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: staggeredgeometryhelper.hh:79
typename GridView::ctype Scalar
Definition: staggeredgeometryhelper.hh:78
GridIndexType oppositeDof
Definition: staggeredgeometryhelper.hh:82
std::array< GridIndexType, upwindSchemeOrder-1 > inAxisBackwardDofs
Definition: staggeredgeometryhelper.hh:86
GridIndexType oppositeDof
Definition: staggeredgeometryhelper.hh:103
typename IndexTraits< GridView >::GridIndex GridIndexType
Definition: staggeredgeometryhelper.hh:100
Scalar selfToOppositeDistance
Definition: staggeredgeometryhelper.hh:104
typename GridView::ctype Scalar
Definition: staggeredgeometryhelper.hh:99
GridIndexType selfDof
Definition: staggeredgeometryhelper.hh:102
Helper class constructing the dual grid finite volume geometries for the free flow staggered discreti...
Definition: staggeredgeometryhelper.hh:127
Detail::AxisData< GridView, upwindSchemeOrder > AxisData
Definition: staggeredgeometryhelper.hh:147
SmallLocalIndexType localFaceIndex() const
Returns the local index of the face (i.e. the intersection)
Definition: staggeredgeometryhelper.hh:164
FreeFlowStaggeredGeometryHelper(const Element &element, const GridView &gridView)
Definition: staggeredgeometryhelper.hh:149
AxisData axisData() const
Returns a copy of the axis data.
Definition: staggeredgeometryhelper.hh:172
unsigned int directionIndex() const
Returns the dirction index of the primary facet (0 = x, 1 = y, 2 = z)
Definition: staggeredgeometryhelper.hh:188
void updateLocalFace(const IntersectionMapper &intersectionMapper_, const Intersection &intersection)
update the local face
Definition: staggeredgeometryhelper.hh:154
std::array< PairData, numPairs > pairData() const
Returns a copy of the pair data.
Definition: staggeredgeometryhelper.hh:180
unsigned int directionIndex(const Intersection &intersection) const
Returns the dirction index of the facet passed as an argument (0 = x, 1 = y, 2 = z)
Definition: staggeredgeometryhelper.hh:196