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-FileCopyrightText: 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 | * \brief Geometry helper for face-centered staggered scheme. | ||
11 | */ | ||
12 | #ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_GEOMETRY_HELPER_HH | ||
13 | #define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_GEOMETRY_HELPER_HH | ||
14 | |||
15 | #include <array> | ||
16 | #include <cassert> | ||
17 | #include <dune/common/exceptions.hh> | ||
18 | #include <dumux/common/indextraits.hh> | ||
19 | #include <dumux/discretization/facecentered/staggered/gridsupportsconcavecorners.hh> | ||
20 | |||
21 | namespace Dumux { | ||
22 | /*! | ||
23 | * \ingroup FaceCenteredStaggeredDiscretization | ||
24 | * \brief Face centered staggered geometry helper | ||
25 | */ | ||
26 | template<class GridView> | ||
27 |
1/4✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
4 | class FaceCenteredStaggeredGeometryHelper |
28 | { | ||
29 | using GridIndexType = typename IndexTraits<GridView>::GridIndex; | ||
30 | using SmallLocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex; | ||
31 | using Element = typename GridView::template Codim<0>::Entity; | ||
32 | using Facet = typename GridView::template Codim<1>::Entity; | ||
33 | |||
34 | public: | ||
35 | static constexpr auto dim = GridView::Grid::dimension; | ||
36 | static constexpr auto numElementFaces = dim * 2; | ||
37 | static constexpr auto numLateralFacesPerScv = 2 * (dim - 1); | ||
38 | |||
39 | 74 | FaceCenteredStaggeredGeometryHelper(const GridView& gridView) : gridView_(gridView) {} | |
40 | |||
41 | //! Returns the local index of the opposing face. | ||
42 | 9969044 | static constexpr SmallLocalIndexType localOppositeIdx(const SmallLocalIndexType ownLocalFaceIndex) | |
43 | { | ||
44 | 9969044 | return isOdd_(ownLocalFaceIndex) ? (ownLocalFaceIndex - 1) : (ownLocalFaceIndex + 1); | |
45 | } | ||
46 | |||
47 | //! Return the local index of a lateral orthogonal scvf | ||
48 | 348555706 | static constexpr int lateralOrthogonalScvfLocalIndex(const SmallLocalIndexType ownLocalScvfIndex) | |
49 | { | ||
50 | if constexpr(GridView::Grid::dimension == 1) | ||
51 | { | ||
52 | ✗ | assert(false && "There are no lateral scvfs in 1D"); | |
53 | return -1; | ||
54 | } | ||
55 | |||
56 | if constexpr (GridView::Grid::dimension == 2) | ||
57 | { | ||
58 |
8/9✓ Branch 0 taken 42975978 times.
✓ Branch 1 taken 43122721 times.
✓ Branch 2 taken 42931680 times.
✓ Branch 3 taken 43325253 times.
✓ Branch 4 taken 43026963 times.
✓ Branch 5 taken 43176930 times.
✓ Branch 6 taken 42991097 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 43279404 times.
|
344830026 | switch (ownLocalScvfIndex) |
59 | { | ||
60 | case 1: return 7; | ||
61 | 42975978 | case 7: return 1; | |
62 | 43122721 | case 2: return 10; | |
63 | 42931680 | case 10: return 2; | |
64 | 43325253 | case 4: return 8; | |
65 | 43026963 | case 8: return 4; | |
66 | 43176930 | case 5: return 11; | |
67 | 42991097 | case 11: return 5; | |
68 | ✗ | default: | |
69 | { | ||
70 | ✗ | assert(false && "No lateral orthogonal scvf found"); | |
71 | return -1; | ||
72 | } | ||
73 | } | ||
74 | } | ||
75 | else | ||
76 | { | ||
77 |
24/25✓ Branch 0 taken 164613 times.
✓ Branch 1 taken 146005 times.
✓ Branch 2 taken 165301 times.
✓ Branch 3 taken 146005 times.
✓ Branch 4 taken 157473 times.
✓ Branch 5 taken 146005 times.
✓ Branch 6 taken 156801 times.
✓ Branch 7 taken 142605 times.
✓ Branch 8 taken 164223 times.
✓ Branch 9 taken 145013 times.
✓ Branch 10 taken 167255 times.
✓ Branch 11 taken 144893 times.
✓ Branch 12 taken 159315 times.
✓ Branch 13 taken 142605 times.
✓ Branch 14 taken 156411 times.
✓ Branch 15 taken 162183 times.
✓ Branch 16 taken 155043 times.
✓ Branch 17 taken 162183 times.
✓ Branch 18 taken 154371 times.
✓ Branch 19 taken 165975 times.
✓ Branch 20 taken 158155 times.
✓ Branch 21 taken 162871 times.
✓ Branch 22 taken 154371 times.
✗ Branch 23 not taken.
✓ Branch 24 taken 146005 times.
|
3725680 | switch (ownLocalScvfIndex) |
78 | { | ||
79 | case 1: return 11; | ||
80 | 164613 | case 11: return 1; | |
81 | 146005 | case 2: return 16; | |
82 | 165301 | case 16: return 2; | |
83 | 146005 | case 3: return 21; | |
84 | 157473 | case 21: return 3; | |
85 | 146005 | case 4: return 26; | |
86 | 156801 | case 26: return 4; | |
87 | 142605 | case 6: return 12; | |
88 | 164223 | case 12: return 6; | |
89 | 145013 | case 7: return 17; | |
90 | 167255 | case 17: return 7; | |
91 | 144893 | case 8: return 22; | |
92 | 159315 | case 22: return 8; | |
93 | 142605 | case 9: return 27; | |
94 | 156411 | case 27: return 9; | |
95 | 162183 | case 13: return 23; | |
96 | 155043 | case 23: return 13; | |
97 | 162183 | case 14: return 28; | |
98 | 154371 | case 28: return 14; | |
99 | 165975 | case 18: return 24; | |
100 | 158155 | case 24: return 18; | |
101 | 162871 | case 19: return 29; | |
102 | 154371 | case 29: return 19; | |
103 | |||
104 | ✗ | default: | |
105 | { | ||
106 | ✗ | assert(false && "No lateral orthogonal scvf found"); | |
107 | return -1; | ||
108 | } | ||
109 | } | ||
110 | } | ||
111 | } | ||
112 | |||
113 | //! Returns the local indices of the faces lateral to the own one. | ||
114 | 9825714 | static constexpr auto localLaterFaceIndices(const SmallLocalIndexType ownLocalFaceIndex) | |
115 | { | ||
116 | 9825650 | constexpr auto table = [] | |
117 | { | ||
118 | using Table = std::array<std::array<SmallLocalIndexType, numLateralFacesPerScv>, numElementFaces>; | ||
119 | if constexpr (dim == 1) | ||
120 | return Table{}; | ||
121 | else if constexpr (dim == 2) | ||
122 | return Table {{ {2,3}, {2,3}, {0,1}, {0,1} }}; | ||
123 | else | ||
124 | return Table {{ {2,3,4,5}, {2,3,4,5}, {0,1,4,5}, {0,1,4,5}, {0,1,2,3}, {0,1,2,3} }}; | ||
125 | }(); | ||
126 | |||
127 | 9825650 | return table[ownLocalFaceIndex]; | |
128 | } | ||
129 | |||
130 | //! Returns an element's facet based on the local facet index. | ||
131 | 19649656 | static Facet facet(const SmallLocalIndexType localFacetIdx, const Element& element) | |
132 | { | ||
133 |
1/2✓ Branch 1 taken 1692800 times.
✗ Branch 2 not taken.
|
21342456 | return element.template subEntity <1> (localFacetIdx); |
134 | } | ||
135 | |||
136 | //! Returns an element's intersection based on the local facet index. | ||
137 | 19673800 | auto intersection(const SmallLocalIndexType localFacetIdx, const Element& element) const | |
138 | { | ||
139 |
10/15✓ Branch 1 taken 43387720 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 26441640 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1019032 times.
✓ Branch 7 taken 1573548 times.
✓ Branch 11 taken 1573548 times.
✗ Branch 12 not taken.
✓ Branch 13 taken 2592580 times.
✗ Branch 14 not taken.
✓ Branch 0 taken 872 times.
✓ Branch 3 taken 21215620 times.
✓ Branch 8 taken 2560416 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 4267360 times.
|
62575484 | for (const auto& intersection : intersections(gridView(), element)) |
140 | { | ||
141 |
2/2✓ Branch 0 taken 19673800 times.
✓ Branch 1 taken 29555700 times.
|
49229500 | if (intersection.indexInInside() == localFacetIdx) |
142 | 22399776 | return intersection; | |
143 | } | ||
144 | ✗ | DUNE_THROW(Dune::InvalidStateException, "localFacetIdx " << localFacetIdx << " out of range"); | |
145 | } | ||
146 | |||
147 | //! Returns true if the IP of an scvf lies on a concave corner | ||
148 | template<class FVElementGeometry, class SubControlVolumeFace> | ||
149 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 132555420 times.
|
132555420 | static bool scvfIntegrationPointInConcaveCorner(const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf) |
150 | { | ||
151 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 132555420 times.
|
132555420 | assert (scvf.isLateral()); |
152 | |||
153 | using Grid = typename FVElementGeometry::GridGeometry::Grid; | ||
154 | if constexpr (!GridSupportsConcaveCorners<Grid>::value) | ||
155 | return false; | ||
156 | else | ||
157 | { | ||
158 |
1/2✓ Branch 0 taken 4614292 times.
✗ Branch 1 not taken.
|
4614292 | if (scvf.boundary()) |
159 | return false; | ||
160 | |||
161 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4614292 times.
|
4614292 | const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx()); |
162 | 4614292 | const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx()); | |
163 | 4614292 | int onBoundaryCounter = 0; | |
164 | 4614292 | onBoundaryCounter += static_cast<int>(insideScv.boundary()); | |
165 | 4614292 | onBoundaryCounter += static_cast<int>(outsideScv.boundary()); | |
166 | 4614292 | return onBoundaryCounter == 1; | |
167 | } | ||
168 | } | ||
169 | |||
170 | template<class SubControlVolumeFace, class SubControlVolume> | ||
171 | 94928 | static SmallLocalIndexType localIndexOutsideScvfWithSameIntegrationPoint(const SubControlVolumeFace& scvf, | |
172 | const SubControlVolume& scv) | ||
173 | { | ||
174 | // In 2D there are 3 non-boundary faces per scv. In 3D, there are 5. | ||
175 | // This number of scvfs per scv is used as an offset to find the indexes in the outside half-scv. | ||
176 |
2/2✓ Branch 0 taken 49784 times.
✓ Branch 1 taken 45144 times.
|
94928 | const SmallLocalIndexType offset = (dim == 2) ? 3 : 5; |
177 | // For half-scvs with odd indexes, the outside half-scv has scvf local indexes with + offset. | ||
178 | // For half-scvs with even indexes, the outside half-scv has scvf local indexes have a - offset. | ||
179 | 94928 | return isOdd_(scv.indexInElement()) ? scvf.localIndex() - offset : scvf.localIndex() + offset; | |
180 | } | ||
181 | |||
182 | 19673800 | const GridView& gridView() const | |
183 | 19673800 | { return gridView_; } | |
184 | |||
185 | private: | ||
186 | |||
187 | 10063972 | static constexpr bool isOdd_(int number) | |
188 |
6/6✓ Branch 0 taken 784397 times.
✓ Branch 1 taken 779757 times.
✓ Branch 2 taken 4232839 times.
✓ Branch 3 taken 4232839 times.
✓ Branch 4 taken 17070 times.
✓ Branch 5 taken 17070 times.
|
10063972 | { return number % 2; } |
189 | |||
190 | GridView gridView_; | ||
191 | }; | ||
192 | |||
193 | } // end namespace Dumux | ||
194 | |||
195 | #endif | ||
196 |