3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_STAGGERED_GEOMETRY_HELPER_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_GEOMETRY_HELPER_HH
26
27#include <dune/geometry/multilineargeometry.hh>
28
29#include <dumux/common/math.hh>
31#include <type_traits>
32#include <algorithm>
33#include <array>
34#include <bitset>
35#include <stack>
36
37namespace Dumux
38{
39
40namespace Detail {
41
46template<class GridView, int upwindSchemeOrder>
48{
49 using Scalar = typename GridView::ctype;
50 using GlobalPosition = typename GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate;
53
54 std::bitset<upwindSchemeOrder> hasParallelNeighbor;
55 std::array<GridIndexType, upwindSchemeOrder> parallelDofs;
56 std::array<Scalar, upwindSchemeOrder> parallelCellWidths;
57 bool hasOuterLateral = false;
58 std::pair<GridIndexType, GridIndexType> lateralPair;
62};
63
64
69template<class GridView, int upwindSchemeOrder>
70struct AxisData;
75template<class GridView, int upwindSchemeOrder>
77{
78 using Scalar = typename GridView::ctype;
80
83 std::bitset<upwindSchemeOrder-1> hasForwardNeighbor;
84 std::bitset<upwindSchemeOrder-1> hasBackwardNeighbor;
85 std::array<GridIndexType, upwindSchemeOrder-1> inAxisForwardDofs;
86 std::array<GridIndexType, upwindSchemeOrder-1> inAxisBackwardDofs;
88 std::array<Scalar, upwindSchemeOrder-1> inAxisForwardDistances;
89 std::array<Scalar, upwindSchemeOrder-1> inAxisBackwardDistances;
90};
91
96template<class GridView>
97struct AxisData<GridView, 1>
98{
99 using Scalar = typename GridView::ctype;
101
105};
106
107} // namespace Detail
108
113template<class Vector>
114inline static unsigned int directionIndex(Vector&& vector)
115{
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();
118}
119
125template<class GridView, int upwindSchemeOrder>
127{
128 using Scalar = typename GridView::ctype;
129 static constexpr auto dim = GridView::dimension;
130
131 using Element = typename GridView::template Codim<0>::Entity;
132 using Intersection = typename GridView::Intersection;
133 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
134
135 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
136 using SmallLocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
137
138 //TODO include assert that checks for quad geometry
139 static constexpr auto codimIntersection = 1;
140 static constexpr auto numfacets = dim * 2;
141 static constexpr auto numPairs = 2 * (dim - 1);
142
143 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
144
145public:
148
149 FreeFlowStaggeredGeometryHelper(const Element& element, const GridView& gridView) : element_(element), elementGeometry_(element.geometry()), gridView_(gridView)
150 { }
151
153 template<class IntersectionMapper>
154 void updateLocalFace(const IntersectionMapper& intersectionMapper_, const Intersection& intersection)
155 {
156 intersection_ = intersection;
157 fillAxisData_();
158 fillPairData_();
159 }
160
164 SmallLocalIndexType localFaceIndex() const
165 {
166 return intersection_.indexInInside();
167 }
168
173 {
174 return axisData_;
175 }
176
180 std::array<PairData, numPairs> pairData() const
181 {
182 return pairData_;
183 }
184
188 unsigned int directionIndex() const
189 {
190 return Dumux::directionIndex(intersection_.centerUnitOuterNormal());
191 }
192
196 unsigned int directionIndex(const Intersection& intersection) const
197 {
198 return Dumux::directionIndex(std::move(intersection.centerUnitOuterNormal()));
199 }
200
201private:
202
207 void fillAxisData_()
208 {
209 fillOuterAxisData_(std::integral_constant<bool, useHigherOrder>{});
210
211 const auto inIdx = intersection_.indexInInside();
212 const auto oppIdx = localOppositeIdx_(inIdx);
213
214 // Set the self Dof
215 axisData_.selfDof = gridView_.indexSet().subIndex(intersection_.inside(), inIdx, codimIntersection);
216
217 // Set the opposite Dof
218 axisData_.oppositeDof = gridView_.indexSet().subIndex(intersection_.inside(), oppIdx, codimIntersection);
219
220 // Set the Self to Opposite Distance
221 const auto self = getFacet_(inIdx, element_);
222 const auto opposite = getFacet_(oppIdx, element_);
223 axisData_.selfToOppositeDistance = (self.geometry().center() - opposite.geometry().center()).two_norm();
224
225 }
226
231 void fillOuterAxisData_(std::false_type) {}
232
237 void fillOuterAxisData_(std::true_type)
238 {
239 // reset the axis data struct
240 axisData_ = {};
241
242 const auto numForwardBackwardAxisDofs = axisData_.inAxisForwardDofs.size();
243
244 // Set the Forward Dofs
245 std::stack<Element> inAxisForwardElementStack;
246 const auto inIdx = intersection_.indexInInside();
247 auto selfFace = getFacet_(inIdx, element_);
248 if(intersection_.neighbor())
249 {
250 inAxisForwardElementStack.push(intersection_.outside());
251 bool keepStackingForward = (inAxisForwardElementStack.size() < numForwardBackwardAxisDofs);
252 while(keepStackingForward)
253 {
254 auto e = inAxisForwardElementStack.top();
255 for(const auto& intersection : intersections(gridView_,e))
256 {
257 if( (intersection.indexInInside() == inIdx ) )
258 {
259 if( intersection.neighbor())
260 {
261 inAxisForwardElementStack.push(intersection.outside());
262 keepStackingForward = (inAxisForwardElementStack.size() < numForwardBackwardAxisDofs);
263 }
264 else
265 {
266 keepStackingForward = false;
267 }
268 }
269 }
270 }
271 }
272
273 std::vector<GlobalPosition> forwardFaceCoordinates(inAxisForwardElementStack.size(), selfFace.geometry().center());
274 while(!inAxisForwardElementStack.empty())
275 {
276 int forwardIdx = inAxisForwardElementStack.size()-1;
277 axisData_.hasForwardNeighbor.set(forwardIdx, true);
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();
282 }
283
284 const auto self = getFacet_(inIdx, element_);
285 for (int i = 0; i< forwardFaceCoordinates.size(); i++)
286 {
287 if (i == 0)
288 axisData_.inAxisForwardDistances[i] = (forwardFaceCoordinates[i] - self.geometry().center()).two_norm();
289 else
290 axisData_.inAxisForwardDistances[i] = (forwardFaceCoordinates[i]- forwardFaceCoordinates[i-1]).two_norm();
291 }
292
293 // Set the Backward Dofs
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_))
298 {
299 if(intersection.indexInInside() == oppIdx)
300 {
301 if(intersection.neighbor())
302 {
303 inAxisBackwardElementStack.push(intersection.outside());
304 bool keepStackingBackward = (inAxisBackwardElementStack.size() < numForwardBackwardAxisDofs);
305 while(keepStackingBackward)
306 {
307 auto e = inAxisBackwardElementStack.top();
308 for(const auto& intersectionOut : intersections(gridView_,e))
309 {
310
311 if( (intersectionOut.indexInInside() == oppIdx ) )
312 {
313 if( intersectionOut.neighbor())
314 {
315 inAxisBackwardElementStack.push(intersectionOut.outside());
316 keepStackingBackward = (inAxisBackwardElementStack.size() < numForwardBackwardAxisDofs);
317 }
318 else
319 {
320 keepStackingBackward = false;
321 }
322 }
323 }
324 }
325 }
326 }
327 }
328
329 std::vector<GlobalPosition> backwardFaceCoordinates(inAxisBackwardElementStack.size(), oppFace.geometry().center());
330 while(!inAxisBackwardElementStack.empty())
331 {
332 int backwardIdx = inAxisBackwardElementStack.size()-1;
333 axisData_.hasBackwardNeighbor.set(backwardIdx, true);
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();
338 }
339
340 const auto opposite = getFacet_(oppIdx, element_);
341 for (int i = 0; i< backwardFaceCoordinates.size(); i++)
342 {
343 if (i == 0)
344 axisData_.inAxisBackwardDistances[i] = (backwardFaceCoordinates[i] - opposite.geometry().center()).two_norm();
345 else
346 axisData_.inAxisBackwardDistances[i] = (backwardFaceCoordinates[i] - backwardFaceCoordinates[i-1]).two_norm();
347 }
348 }
349
354 void fillPairData_()
355 {
356 // reset the pair data structs
357 pairData_ = {};
358
359 // set basic global positions
360 const auto& selfElementCenter = element_.geometry().center();
361 const auto& selfFacetCenter = intersection_.geometry().center();
362
363 // get the inner lateral Dof Index
364 SmallLocalIndexType numPairInnerLateralIdx = 0;
365 for (const auto& innerElementIntersection : intersections(gridView_, element_))
366 {
367 if (facetIsNormal_(innerElementIntersection.indexInInside(), intersection_.indexInInside()))
368 {
369 const auto innerElementIntersectionIdx = innerElementIntersection.indexInInside();
370 setLateralPairFirstInfo_(innerElementIntersectionIdx, element_, numPairInnerLateralIdx);
371
372 const auto distance = innerElementIntersection.geometry().center() - selfElementCenter;
373 pairData_[numPairInnerLateralIdx].lateralStaggeredFaceCenter = std::move(selfFacetCenter + distance);
374
375 numPairInnerLateralIdx++;
376 }
377 }
378
379 // get the outer lateral Dof Index
380 SmallLocalIndexType numPairOuterLateralIdx = 0;
381 if (intersection_.neighbor())
382 {
383 // the direct neighbor element and the respective intersection index
384 const auto& directNeighborElement = intersection_.outside();
385
386 for (const auto& directNeighborElementIntersection : intersections(gridView_, directNeighborElement))
387 {
388 // skip the directly neighboring face itself and its opposing one
389 if (facetIsNormal_(directNeighborElementIntersection.indexInInside(), intersection_.indexInOutside()))
390 {
391 const auto directNeighborElemIsIdx = directNeighborElementIntersection.indexInInside();
392 setLateralPairSecondInfo_(directNeighborElemIsIdx, directNeighborElement, numPairOuterLateralIdx);
393 numPairOuterLateralIdx++;
394 }
395 }
396 }
397 else // intersection is on boundary
398 {
399 for (const auto& intersection : intersections(gridView_, element_))
400 {
401 if (facetIsNormal_(intersection.indexInInside(), intersection_.indexInInside()))
402 {
403 assert(!pairData_[numPairOuterLateralIdx].hasOuterLateral);
404
405 const auto normalDistanceoffset = selfFacetCenter - selfElementCenter;
406 pairData_[numPairOuterLateralIdx].lateralDistance = std::move(normalDistanceoffset.two_norm());
407 numPairOuterLateralIdx++;
408 }
409 }
410 }
411
412 fillParallelPairData_(std::integral_constant<bool, useHigherOrder>{});
413 }
414
419 void fillParallelPairData_(std::false_type)
420 {
421 // set basic global positions and stencil size definitions
422 // get the parallel Dofs
423 const auto parallelLocalIdx = intersection_.indexInInside();
424 SmallLocalIndexType numPairParallelIdx = 0;
425
426 for (const auto& intersection : intersections(gridView_, element_))
427 {
428 if (facetIsNormal_(intersection.indexInInside(), parallelLocalIdx))
429 {
430 auto parallelAxisIdx = directionIndex(intersection);
431 if (intersection.neighbor())
432 {
433 // If the lateral intersection has a neighboring cell, go in and store the parallel information.
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);
438 }
439
440 numPairParallelIdx++;
441 }
442 }
443 }
444
449 void fillParallelPairData_(std::true_type)
450 {
451 // set basic global positions and stencil size definitions
452 const auto numParallelFaces = pairData_[0].parallelCellWidths.size();
453
454 // get the parallel Dofs
455 const auto parallelLocalIdx = intersection_.indexInInside();
456 SmallLocalIndexType numPairParallelIdx = 0;
457 std::stack<Element> parallelElementStack;
458 for(const auto& intersection : intersections(gridView_, element_))
459 {
460 if( facetIsNormal_(intersection.indexInInside(), parallelLocalIdx) )
461 {
462 if( intersection.neighbor() )
463 {
464 auto parallelAxisIdx = directionIndex(intersection);
465 auto localLateralIntersectionIndex = intersection.indexInInside();
466 auto e = element_;
467
468 bool keepStacking = (parallelElementStack.size() < numParallelFaces);
469 while(keepStacking)
470 {
471 for(const auto& lateralIntersection : intersections(gridView_, e))
472 {
473 if( lateralIntersection.indexInInside() == localLateralIntersectionIndex )
474 {
475 if( lateralIntersection.neighbor() )
476 {
477 parallelElementStack.push(lateralIntersection.outside());
478 keepStacking = (parallelElementStack.size() < numParallelFaces);
479 }
480 else
481 {
482 keepStacking = false;
483 }
484 }
485 }
486 e = parallelElementStack.top();
487 }
488
489 while(!parallelElementStack.empty())
490 {
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();
495 }
496
497 }
498
499 numPairParallelIdx++;
500 }
501 }
502 }
503
509 int localOppositeIdx_(const int idx) const
510 {
511 return (idx % 2) ? (idx - 1) : (idx + 1);
512 }
513
520 bool facetIsNormal_(const int selfIdx, const int otherIdx) const
521 {
522 return !(selfIdx == otherIdx || localOppositeIdx_(selfIdx) == otherIdx);
523 };
524
525 auto getFacet_(const int localFacetIdx, const Element& element) const
526 {
527 return element.template subEntity <1> (localFacetIdx);
528 };
529
531 void setLateralPairFirstInfo_(const int isIdx, const Element& element, const int numPairsIdx)
532 {
533 // store the inner lateral dofIdx
534 pairData_[numPairsIdx].lateralPair.first = gridView_.indexSet().subIndex(element, isIdx, codimIntersection);
535
536 // store the local lateral facet index
537 pairData_[numPairsIdx].localLateralFaceIdx = isIdx;
538 }
539
541 void setLateralPairSecondInfo_(const int isIdx, const Element& element, const int numPairsIdx)
542 {
543 // store the dofIdx
544 pairData_[numPairsIdx].hasOuterLateral = true;
545 pairData_[numPairsIdx].lateralPair.second = gridView_.indexSet().subIndex(element, isIdx, codimIntersection);
546
547 // set basic global positions
548 const auto& selfFacetCenter = intersection_.geometry().center();
549 const auto& selfElementCenter = element_.geometry().center();
550
551 const auto& neighborElement = intersection_.outside();
552 const auto& neighborElementCenter = neighborElement.geometry().center();
553 const auto& neighborFacetCenter = getFacet_(intersection_.indexInOutside(), neighborElement).geometry().center();
554
555 const Scalar insideLateralDistance = (selfFacetCenter - selfElementCenter).two_norm();
556 const Scalar outsideLateralDistance = (neighborFacetCenter - neighborElementCenter).two_norm();
557
558 pairData_[numPairsIdx].lateralDistance = insideLateralDistance + outsideLateralDistance;
559 }
560
562 Scalar setParallelPairCellWidths_(const Element& element, const int parallelAxisIdx) const
563 {
564 std::vector<GlobalPosition> faces;
565 faces.reserve(numfacets);
566 for (const auto& intElement : intersections(gridView_, element))
567 {
568 faces.push_back(intElement.geometry().center());
569 }
570 switch(parallelAxisIdx)
571 {
572 case 0:
573 return (faces[1] - faces[0]).two_norm();
574 case 1:
575 return (faces[3] - faces[2]).two_norm();
576 case 2:
577 {
578 assert(dim == 3);
579 return (faces[5] - faces[4]).two_norm();
580 }
581 default:
582 DUNE_THROW(Dune::InvalidStateException, "Something went terribly wrong");
583 }
584 }
585
586 Intersection intersection_;
587 const Element element_;
588 const typename Element::Geometry elementGeometry_;
589 const GridView gridView_;
590 AxisData axisData_;
591 std::array<PairData, numPairs> pairData_;
592};
593
594} // end namespace Dumux
595
596#endif
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
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44
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