version 3.8
staggeredgeometryhelper.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_STAGGERED_GEOMETRY_HELPER_HH
13#define DUMUX_DISCRETIZATION_STAGGERED_GEOMETRY_HELPER_HH
14
15#include <dune/geometry/multilineargeometry.hh>
16
17#include <dumux/common/math.hh>
19#include <type_traits>
20#include <algorithm>
21#include <array>
22#include <bitset>
23
24namespace Dumux {
25namespace Detail {
26
58template<class GridView, int upwindSchemeOrder>
60{
61 using Scalar = typename GridView::ctype;
62 using GlobalPosition = typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate;
65
66 std::bitset<upwindSchemeOrder> hasParallelNeighbor;
69 std::array<GridIndexType, upwindSchemeOrder> parallelDofs;
70 std::array<Scalar, upwindSchemeOrder> parallelCellWidths;
71 bool hasOuterLateral = false;
72 std::pair<GridIndexType, GridIndexType> lateralPair;
76};
77
82template<class GridView, int upwindSchemeOrder>
84{
85 using Scalar = typename GridView::ctype;
87
90 std::bitset<upwindSchemeOrder-1> hasForwardNeighbor;
91 std::bitset<upwindSchemeOrder-1> hasBackwardNeighbor;
92 std::array<GridIndexType, upwindSchemeOrder-1> inAxisForwardDofs;
93 std::array<GridIndexType, upwindSchemeOrder-1> inAxisBackwardDofs;
95 std::array<Scalar, upwindSchemeOrder-1> inAxisForwardDistances;
96 std::array<Scalar, upwindSchemeOrder-1> inAxisBackwardDistances;
97};
98
103template<class GridView>
104struct AxisData<GridView, 1>
105{
106 using Scalar = typename GridView::ctype;
108
112};
113
114} // namespace Detail
115
120template<class Vector>
121inline static unsigned int directionIndex(Vector&& vector)
122{
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();
125}
126
132template<class GridView, int upwindSchemeOrder>
134{
135 using Scalar = typename GridView::ctype;
136 static constexpr auto dim = GridView::dimension;
137
138 using Element = typename GridView::template Codim<0>::Entity;
139 using Intersection = typename GridView::Intersection;
140 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
141
142 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
143 using SmallLocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
144
145 //TODO include assert that checks for quad geometry
146 static constexpr auto codimIntersection = 1;
147 static constexpr auto numFacets = dim * 2;
148 static constexpr auto numPairs = 2 * (dim - 1);
149
150 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
151
152public:
155
156 FreeFlowStaggeredGeometryHelper(const Element& element, const GridView& gridView)
157 : element_(element)
158 , gridView_(gridView)
159 { }
160
162 template<class IntersectionMapper>
163 void updateLocalFace(const IntersectionMapper&, const Intersection& intersection)
164 {
165 intersection_ = intersection;
166 fillAxisData_();
167 fillPairData_();
168 }
169
173 SmallLocalIndexType localFaceIndex() const
174 {
175 return intersection_.indexInInside();
176 }
177
182 {
183 return axisData_;
184 }
185
189 std::array<PairData, numPairs> pairData() const
190 {
191 return pairData_;
192 }
193
197 unsigned int directionIndex() const
198 {
199 return Dumux::directionIndex(intersection_.centerUnitOuterNormal());
200 }
201
205 unsigned int directionIndex(const Intersection& intersection) const
206 {
207 return Dumux::directionIndex(std::move(intersection.centerUnitOuterNormal()));
208 }
209
210private:
211
216 void fillAxisData_()
217 {
218 if constexpr (useHigherOrder)
219 fillOuterAxisData_();
220
221 const auto inIdx = intersection_.indexInInside();
222 const auto oppIdx = localOppositeIdx_(inIdx);
223
224 // Set the self Dof
225 axisData_.selfDof = gridView_.indexSet().subIndex(intersection_.inside(), inIdx, codimIntersection);
226
227 // Set the opposite Dof
228 axisData_.oppositeDof = gridView_.indexSet().subIndex(intersection_.inside(), oppIdx, codimIntersection);
229
230 // Set the Self to Opposite Distance
231 const auto self = getFacet_(inIdx, element_);
232 const auto opposite = getFacet_(oppIdx, element_);
233 axisData_.selfToOppositeDistance = (self.geometry().center() - opposite.geometry().center()).two_norm();
234
235 }
236
242 void fillOuterAxisData_()
243 {
244 // reset the axis data struct
245 axisData_ = {};
246
247 // Set the forward dofs
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);
252
253 // Set the backward dofs
254 const auto oppositeFaceLocalIdx = localOppositeIdx_(thisFaceLocalIdx);
255 const auto& oppositeFaceCenterPos = getFacet_(oppositeFaceLocalIdx, element_).geometry().center();
256 for (const auto& intersection : intersections(gridView_, element_))
257 {
258 if (intersection.indexInInside() == oppositeFaceLocalIdx && intersection.neighbor())
259 {
260 addBackwardNeighborAxisData_(intersection.outside(), 0, oppositeFaceLocalIdx, oppositeFaceCenterPos);
261 break;
262 }
263 }
264 }
265
269 void addForwardNeighborAxisData_(const Element& neighbor,
270 const std::size_t distance,
271 const unsigned int localIntersectionIndex,
272 const GlobalPosition& faceCenterPos)
273 {
274 // fill the forward axis data for this element in the stencil
275 axisData_.hasForwardNeighbor.set(distance, true);
276 axisData_.inAxisForwardDofs[distance] = gridView_.indexSet().subIndex(neighbor, localIntersectionIndex, codimIntersection);
277 const auto& forwardFacePos = getFacet_(localIntersectionIndex, neighbor).geometry().center();
278 axisData_.inAxisForwardDistances[distance] = (forwardFacePos - faceCenterPos).two_norm();
279
280 // recursion end if we added all axis data for the given stencil size (numParallelFaces)
281 const auto numForwardBackwardAxisDofs = axisData_.inAxisForwardDofs.size();
282 if (distance >= numForwardBackwardAxisDofs-1)
283 return;
284
285 for (const auto& intersection : intersections(gridView_, neighbor))
286 {
287 if (intersection.indexInInside() == localIntersectionIndex && intersection.neighbor())
288 {
289 addForwardNeighborAxisData_(intersection.outside(), distance+1, localIntersectionIndex, forwardFacePos);
290 break;
291 }
292 }
293 }
294
298 void addBackwardNeighborAxisData_(const Element& neighbor,
299 const std::size_t distance,
300 const unsigned int localIntersectionIndex,
301 const GlobalPosition& faceCenterPos)
302 {
303 // fill the forward axis data for this element in the stencil
304 axisData_.hasBackwardNeighbor.set(distance, true);
305 axisData_.inAxisBackwardDofs[distance] = gridView_.indexSet().subIndex(neighbor, localIntersectionIndex, codimIntersection);
306 const auto& backwardFacePos = getFacet_(localIntersectionIndex, neighbor).geometry().center();
307 axisData_.inAxisBackwardDistances[distance] = (backwardFacePos - faceCenterPos).two_norm();
308
309 // recursion end if we added all axis data for the given stencil size (numParallelFaces)
310 const auto numForwardBackwardAxisDofs = axisData_.inAxisForwardDofs.size();
311 if (distance >= numForwardBackwardAxisDofs-1)
312 return;
313
314 for (const auto& intersection : intersections(gridView_, neighbor))
315 {
316 if (intersection.indexInInside() == localIntersectionIndex && intersection.neighbor())
317 {
318 addBackwardNeighborAxisData_(intersection.outside(), distance+1, localIntersectionIndex, backwardFacePos);
319 break;
320 }
321 }
322 }
323
328 void fillPairData_()
329 {
330 // reset the pair data structs
331 pairData_ = {};
332
333 // set basic global positions
334 const auto& selfElementCenter = element_.geometry().center();
335 const auto& selfFacetCenter = intersection_.geometry().center();
336
337 // get the inner lateral Dof Index
338 SmallLocalIndexType numPairInnerLateralIdx = 0;
339 for (const auto& innerElementIntersection : intersections(gridView_, element_))
340 {
341 if (facetIsNormal_(innerElementIntersection.indexInInside(), intersection_.indexInInside()))
342 {
343 const auto innerElementIntersectionIdx = innerElementIntersection.indexInInside();
344 setLateralPairFirstInfo_(innerElementIntersectionIdx, element_, numPairInnerLateralIdx);
345
346 const auto distance = innerElementIntersection.geometry().center() - selfElementCenter;
347 pairData_[numPairInnerLateralIdx].lateralStaggeredFaceCenter = selfFacetCenter+distance;
348
349 numPairInnerLateralIdx++;
350 }
351 }
352
353 // get the outer lateral Dof Index
354 SmallLocalIndexType numPairOuterLateralIdx = 0;
355 if (intersection_.neighbor())
356 {
357 // the direct neighbor element and the respective intersection index
358 const auto& directNeighborElement = intersection_.outside();
359
360 for (const auto& directNeighborElementIntersection : intersections(gridView_, directNeighborElement))
361 {
362 // skip the directly neighboring face itself and its opposing one
363 if (facetIsNormal_(directNeighborElementIntersection.indexInInside(), intersection_.indexInOutside()))
364 {
365 const auto directNeighborElemIsIdx = directNeighborElementIntersection.indexInInside();
366 setLateralPairSecondInfo_(directNeighborElemIsIdx, directNeighborElement, numPairOuterLateralIdx);
367 numPairOuterLateralIdx++;
368 }
369 }
370 }
371 else // intersection is on boundary
372 {
373 for (const auto& intersection : intersections(gridView_, element_))
374 {
375 if (facetIsNormal_(intersection.indexInInside(), intersection_.indexInInside()))
376 {
377 assert(!pairData_[numPairOuterLateralIdx].hasOuterLateral);
378
379 const auto normalDistanceoffset = selfFacetCenter - selfElementCenter;
380 pairData_[numPairOuterLateralIdx].lateralDistance = normalDistanceoffset.two_norm();
381 numPairOuterLateralIdx++;
382 }
383 }
384 }
385
386 // fill the pair data with the parallel dofs and distances
387 const auto parallelLocalIdx = intersection_.indexInInside();
388 SmallLocalIndexType numPairParallelIdx = 0;
389 for (const auto& intersection : intersections(gridView_, element_))
390 {
391 if (facetIsNormal_(intersection.indexInInside(), parallelLocalIdx))
392 {
393 if (intersection.neighbor())
394 {
395 // recursively insert parallel neighbor faces into pair data
396 const auto parallelAxisIdx = directionIndex(intersection);
397 const auto localLateralIntersectionIndex = intersection.indexInInside();
398
399 /*
400 * ------------
401 * | |
402 * | |
403 * | |
404 * iiiiiiiiiii*bbbbbbbbbbb
405 * | o zzzzzzzz |
406 * | o zzzzzzzz |
407 * | o zzzzzzzz |
408 * -----------------------
409 *
410 * i:intersection,o:intersection_, b: outerIntersection, z: intersection_.outside()
411 */
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;
417
418 /* ------------
419 * | o
420 * | o
421 * | o
422 * iiiiiiiiiii------------
423 * | zzzzzzzz b |
424 * | zzzzzzzz b |
425 * | zzzzzzzz b |
426 * -----------------------
427 *
428 * i:intersection,o:intersection_, b: outerIntersection, z: intersection.outside()
429 */
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;
435
436
437 addParallelNeighborPairData_(intersection.outside(), 0, localLateralIntersectionIndex, parallelLocalIdx, parallelAxisIdx, numPairParallelIdx);
438 }
439
440 ++numPairParallelIdx;
441 }
442 }
443 }
444
445 // zero distance means direct neighbor to element_
446 void addParallelNeighborPairData_(const Element& neighbor,
447 const std::size_t distance,
448 const unsigned int localLateralIntersectionIndex,
449 const unsigned int parallelLocalIdx,
450 const unsigned int parallelAxisIdx,
451 const SmallLocalIndexType dataIdx)
452 {
453 // fill the pair data for this element in the stencil
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);
457
458 // recursion end if we added all parallel data for the given stencil size (numParallelFaces)
459 const auto numParallelFaces = pairData_[0].parallelCellWidths.size();
460 if (distance >= numParallelFaces-1)
461 return;
462
463 // continue with neighbor's neighbor element in the same direction, if we find one
464 for (const auto& lateralIntersection : intersections(gridView_, neighbor))
465 {
466 // we only expect to find one intersection matching this condition
467 if (lateralIntersection.indexInInside() == localLateralIntersectionIndex)
468 {
469 // if there is no neighbor recursion ends here
470 if (lateralIntersection.neighbor())
471 addParallelNeighborPairData_(lateralIntersection.outside(), distance+1,
472 localLateralIntersectionIndex, parallelLocalIdx, parallelAxisIdx, dataIdx);
473 break;
474 }
475 }
476 }
477
483 int localOppositeIdx_(const int idx) const
484 {
485 return (idx % 2) ? (idx - 1) : (idx + 1);
486 }
487
494 bool facetIsNormal_(const int selfIdx, const int otherIdx) const
495 {
496 return !(selfIdx == otherIdx || localOppositeIdx_(selfIdx) == otherIdx);
497 }
498
499 auto getFacet_(const int localFacetIdx, const Element& element) const
500 {
501 return element.template subEntity <1> (localFacetIdx);
502 }
503
505 void setLateralPairFirstInfo_(const int isIdx, const Element& element, const int numPairsIdx)
506 {
507 // store the inner lateral dofIdx
508 pairData_[numPairsIdx].lateralPair.first = gridView_.indexSet().subIndex(element, isIdx, codimIntersection);
509
510 // store the local lateral facet index
511 pairData_[numPairsIdx].localLateralFaceIdx = isIdx;
512 }
513
515 void setLateralPairSecondInfo_(const int isIdx, const Element& element, const int numPairsIdx)
516 {
517 // store the dofIdx
518 pairData_[numPairsIdx].hasOuterLateral = true;
519 pairData_[numPairsIdx].lateralPair.second = gridView_.indexSet().subIndex(element, isIdx, codimIntersection);
520
521 // set basic global positions
522 const auto& selfFacetCenter = intersection_.geometry().center();
523 const auto& selfElementCenter = element_.geometry().center();
524
525 const auto& neighborElement = intersection_.outside();
526 const auto& neighborElementCenter = neighborElement.geometry().center();
527 const auto& neighborFacetCenter = getFacet_(intersection_.indexInOutside(), neighborElement).geometry().center();
528
529 const Scalar insideLateralDistance = (selfFacetCenter - selfElementCenter).two_norm();
530 const Scalar outsideLateralDistance = (neighborFacetCenter - neighborElementCenter).two_norm();
531
532 pairData_[numPairsIdx].lateralDistance = insideLateralDistance + outsideLateralDistance;
533 }
534
536 Scalar setParallelPairCellWidths_(const Element& element, const int parallelAxisIdx) const
537 {
538 std::vector<GlobalPosition> faces;
539 faces.reserve(numFacets);
540 for (const auto& intersection : intersections(gridView_, element))
541 faces.push_back(intersection.geometry().center());
542
543 switch (parallelAxisIdx)
544 {
545 case 0:
546 return (faces[1] - faces[0]).two_norm();
547 case 1:
548 return (faces[3] - faces[2]).two_norm();
549 case 2:
550 {
551 assert(dim == 3);
552 return (faces[5] - faces[4]).two_norm();
553 }
554 default:
555 DUNE_THROW(Dune::InvalidStateException, "Something went terribly wrong");
556 }
557 }
558
559 Intersection intersection_;
560 const Element& element_;
561 const GridView gridView_;
562 AxisData axisData_;
563 std::array<PairData, numPairs> pairData_;
564};
565
566} // end namespace Dumux
567
568#endif
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
Defines the index types used for grid and local indices.
Define some often used mathematical functions.
Definition: adapt.hh:17
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