Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup FaceCenteredStaggeredDiscretization | ||
10 | * \copydoc Dumux::FaceCenteredStaggeredLocalIntersectionIndexMapper | ||
11 | */ | ||
12 | #ifndef DUMUX_DISCRETIZATION_LOCAL_INTERSECTION_INDEX_MAPPER_HH | ||
13 | #define DUMUX_DISCRETIZATION_LOCAL_INTERSECTION_INDEX_MAPPER_HH | ||
14 | |||
15 | #include <array> | ||
16 | #include <numeric> | ||
17 | #include <dune/common/float_cmp.hh> | ||
18 | #include <dumux/common/indextraits.hh> | ||
19 | #include <dumux/common/math.hh> | ||
20 | #include <dumux/common/parameters.hh> | ||
21 | #include <dumux/discretization/facecentered/staggered/normalaxis.hh> | ||
22 | #include <dumux/discretization/facecentered/staggered/consistentlyorientedgrid.hh> | ||
23 | |||
24 | namespace Dumux { | ||
25 | |||
26 | namespace Detail { | ||
27 | |||
28 | template<class GridView, bool consistentlyOrientedGrid> | ||
29 | class FaceCenteredStaggeredLocalIntersectionIndexMapper; | ||
30 | |||
31 | /*! | ||
32 | * \ingroup FaceCenteredStaggeredDiscretization | ||
33 | * \brief Provides a mapping of local intersection indices (indexInInside) | ||
34 | * such that the local indices always follow the order of a reference element, | ||
35 | * regardless of how the element is oriented. | ||
36 | */ | ||
37 | template<class GridView> | ||
38 | class FaceCenteredStaggeredLocalIntersectionIndexMapper<GridView, false> | ||
39 | { | ||
40 | using SmallLocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex; | ||
41 | using Element = typename GridView::template Codim<0>::Entity; | ||
42 | static constexpr auto numElementFaces = GridView::Grid::dimension * 2; | ||
43 | public: | ||
44 | 247557 | void update(const GridView& gv, const Element& element) | |
45 | { | ||
46 |
4/6✓ Branch 0 taken 8 times.
✓ Branch 1 taken 247549 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 8 times.
✗ Branch 7 not taken.
|
247557 | static const bool makeConsistentlyOriented = getParam<bool>("Grid.MakeConsistentlyOriented", true); |
47 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 247557 times.
|
247557 | if (!makeConsistentlyOriented) |
48 | { | ||
49 | ✗ | std::iota(realToRefMap_.begin(), realToRefMap_.end(), 0); | |
50 | ✗ | refToRealMap_ = realToRefMap_; | |
51 | ✗ | return; | |
52 | } | ||
53 | |||
54 |
5/9✓ Branch 1 taken 109 times.
✓ Branch 2 taken 436 times.
✓ Branch 3 taken 436 times.
✓ Branch 4 taken 1237349 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 9 taken 989792 times.
✗ Branch 10 not taken.
|
2722909 | for (const auto& is : intersections(gv, element)) |
55 | { | ||
56 |
1/2✓ Branch 1 taken 989792 times.
✗ Branch 2 not taken.
|
990228 | const auto& otherOuterNormal = is.centerUnitOuterNormal(); |
57 | 990228 | const auto idx = normalAxis(otherOuterNormal); | |
58 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 990228 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 436 times.
✓ Branch 4 taken 989792 times.
✗ Branch 5 not taken.
|
1980456 | const int positveOrientation = !std::signbit(otherOuterNormal[idx]); |
59 | 990228 | const auto refIdx = idx * 2 + positveOrientation; | |
60 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 990228 times.
✗ Branch 2 not taken.
|
990228 | const auto realIdx = is.indexInInside(); |
61 |
1/3✗ Branch 0 not taken.
✓ Branch 1 taken 990228 times.
✗ Branch 2 not taken.
|
990228 | realToRefMap_[realIdx] = refIdx; |
62 |
3/6✗ Branch 0 not taken.
✓ Branch 1 taken 990228 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 436 times.
✓ Branch 4 taken 989792 times.
✗ Branch 5 not taken.
|
1980456 | refToRealMap_[refIdx] = realIdx; |
63 | } | ||
64 | } | ||
65 | |||
66 | //! Return the intersection's actual local indexInElement given a local reference index. | ||
67 | SmallLocalIndexType realToRefIdx(const SmallLocalIndexType localIsIdx) const | ||
68 |
14/20✓ Branch 1 taken 14144 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 14144 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 839584 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 839584 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 25984 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 25984 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 17040 times.
✓ Branch 19 taken 17040 times.
✓ Branch 20 taken 17040 times.
✓ Branch 21 taken 17040 times.
✓ Branch 22 taken 68160 times.
✓ Branch 23 taken 68160 times.
✓ Branch 24 taken 68160 times.
✓ Branch 25 taken 68160 times.
|
2101136 | { return realToRefMap_[localIsIdx]; } |
69 | |||
70 | //! Return the intersection's local reference indexInElement given an actual local index. | ||
71 | SmallLocalIndexType refToRealIdx(const SmallLocalIndexType localIsIdx) const | ||
72 |
10/24✗ Branch 0 not taken.
✓ Branch 1 taken 436 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 436 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 872 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 8 taken 872 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✓ Branch 13 taken 14144 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 14144 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 28288 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 28288 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 1678656 times.
✗ Branch 26 not taken.
✓ Branch 28 taken 1678656 times.
✗ Branch 29 not taken.
|
3444792 | { return refToRealMap_[localIsIdx]; } |
73 | |||
74 | private: | ||
75 | std::array<SmallLocalIndexType, numElementFaces> realToRefMap_ = {}; | ||
76 | std::array<SmallLocalIndexType, numElementFaces> refToRealMap_ = {}; | ||
77 | }; | ||
78 | |||
79 | /*! | ||
80 | * \ingroup FaceCenteredStaggeredDiscretization | ||
81 | * \brief Provides a mapping of local intersection indices (indexInInside) | ||
82 | * such that the local indices always follow the order of a reference element, | ||
83 | * regardless of how the element in oriented. | ||
84 | * \note This specialization is used for grids not supporting rotated elements. | ||
85 | * No mapping needs to be done here. | ||
86 | */ | ||
87 | template<class GridView> | ||
88 | class FaceCenteredStaggeredLocalIntersectionIndexMapper<GridView, true> | ||
89 | { | ||
90 | using SmallLocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex; | ||
91 | using Element = typename GridView::template Codim<0>::Entity; | ||
92 | public: | ||
93 | |||
94 | //! Update the map for getting the corresponding local face indices in another element. | ||
95 | ✗ | void update(const GridView&, const Element&) {} | |
96 | |||
97 | //! Return the intersection's actual local indexInElement given a local reference index. | ||
98 | SmallLocalIndexType realToRefIdx(const SmallLocalIndexType localIsIdx) const | ||
99 | { return localIsIdx; } | ||
100 | |||
101 | //! Return the intersection's local reference indexInElement given an actual local index. | ||
102 | SmallLocalIndexType refToRealIdx(const SmallLocalIndexType localIsIdx) const | ||
103 | { return localIsIdx; } | ||
104 | }; | ||
105 | |||
106 | } // end namespace Detail | ||
107 | |||
108 | /*! | ||
109 | * \ingroup FaceCenteredStaggeredDiscretization | ||
110 | * \brief Provides a mapping of local intersection indices (indexInInside) | ||
111 | * such that the local indices always follow the order of a reference element, | ||
112 | * regardless of how the element is oriented. Does not do anything for grids | ||
113 | * not supporting rotated elements (such as Dune::YaspGrid). | ||
114 | */ | ||
115 | template<class GridView> | ||
116 | using FaceCenteredStaggeredLocalIntersectionIndexMapper = | ||
117 | Detail::FaceCenteredStaggeredLocalIntersectionIndexMapper<GridView, ConsistentlyOrientedGrid<typename GridView::Grid>{}>; | ||
118 | |||
119 | } // end namespace Dumux | ||
120 | |||
121 | #endif | ||
122 |