version 3.8
discretization/facecentered/staggered/geometryhelper.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_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>
20
21namespace Dumux {
26template<class GridView>
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
34public:
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 FaceCenteredStaggeredGeometryHelper(const GridView& gridView) : gridView_(gridView) {}
40
42 static constexpr SmallLocalIndexType localOppositeIdx(const SmallLocalIndexType ownLocalFaceIndex)
43 {
44 return isOdd_(ownLocalFaceIndex) ? (ownLocalFaceIndex - 1) : (ownLocalFaceIndex + 1);
45 }
46
48 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 switch (ownLocalScvfIndex)
59 {
60 case 1: return 7;
61 case 7: return 1;
62 case 2: return 10;
63 case 10: return 2;
64 case 4: return 8;
65 case 8: return 4;
66 case 5: return 11;
67 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 switch (ownLocalScvfIndex)
78 {
79 case 1: return 11;
80 case 11: return 1;
81 case 2: return 16;
82 case 16: return 2;
83 case 3: return 21;
84 case 21: return 3;
85 case 4: return 26;
86 case 26: return 4;
87 case 6: return 12;
88 case 12: return 6;
89 case 7: return 17;
90 case 17: return 7;
91 case 8: return 22;
92 case 22: return 8;
93 case 9: return 27;
94 case 27: return 9;
95 case 13: return 23;
96 case 23: return 13;
97 case 14: return 28;
98 case 28: return 14;
99 case 18: return 24;
100 case 24: return 18;
101 case 19: return 29;
102 case 29: return 19;
103
104 default:
105 {
106 assert(false && "No lateral orthogonal scvf found");
107 return -1;
108 }
109 }
110 }
111 }
112
114 static constexpr auto localLaterFaceIndices(const SmallLocalIndexType ownLocalFaceIndex)
115 {
116 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 return table[ownLocalFaceIndex];
128 }
129
131 static Facet facet(const SmallLocalIndexType localFacetIdx, const Element& element)
132 {
133 return element.template subEntity <1> (localFacetIdx);
134 }
135
137 auto intersection(const SmallLocalIndexType localFacetIdx, const Element& element) const
138 {
139 for (const auto& intersection : intersections(gridView(), element))
140 {
141 if (intersection.indexInInside() == localFacetIdx)
142 return intersection;
143 }
144 DUNE_THROW(Dune::InvalidStateException, "localFacetIdx " << localFacetIdx << " out of range");
145 }
146
148 template<class FVElementGeometry, class SubControlVolumeFace>
149 static bool scvfIntegrationPointInConcaveCorner(const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
150 {
151 assert (scvf.isLateral());
152
153 using Grid = typename FVElementGeometry::GridGeometry::Grid;
155 return false;
156 else
157 {
158 if (scvf.boundary())
159 return false;
160
161 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
162 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
163 int onBoundaryCounter = 0;
164 onBoundaryCounter += static_cast<int>(insideScv.boundary());
165 onBoundaryCounter += static_cast<int>(outsideScv.boundary());
166 return onBoundaryCounter == 1;
167 }
168 }
169
170 template<class SubControlVolumeFace, class SubControlVolume>
171 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 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 return isOdd_(scv.indexInElement()) ? scvf.localIndex() - offset : scvf.localIndex() + offset;
180 }
181
182 const GridView& gridView() const
183 { return gridView_; }
184
185private:
186
187 static constexpr bool isOdd_(int number)
188 { return number % 2; }
189
190 GridView gridView_;
191};
192
193} // end namespace Dumux
194
195#endif
Face centered staggered geometry helper.
Definition: discretization/facecentered/staggered/geometryhelper.hh:28
static Facet facet(const SmallLocalIndexType localFacetIdx, const Element &element)
Returns an element's facet based on the local facet index.
Definition: discretization/facecentered/staggered/geometryhelper.hh:131
static constexpr auto localLaterFaceIndices(const SmallLocalIndexType ownLocalFaceIndex)
Returns the local indices of the faces lateral to the own one.
Definition: discretization/facecentered/staggered/geometryhelper.hh:114
static bool scvfIntegrationPointInConcaveCorner(const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Returns true if the IP of an scvf lies on a concave corner.
Definition: discretization/facecentered/staggered/geometryhelper.hh:149
static constexpr SmallLocalIndexType localOppositeIdx(const SmallLocalIndexType ownLocalFaceIndex)
Returns the local index of the opposing face.
Definition: discretization/facecentered/staggered/geometryhelper.hh:42
static constexpr auto numLateralFacesPerScv
Definition: discretization/facecentered/staggered/geometryhelper.hh:37
const GridView & gridView() const
Definition: discretization/facecentered/staggered/geometryhelper.hh:182
static constexpr auto numElementFaces
Definition: discretization/facecentered/staggered/geometryhelper.hh:36
FaceCenteredStaggeredGeometryHelper(const GridView &gridView)
Definition: discretization/facecentered/staggered/geometryhelper.hh:39
static SmallLocalIndexType localIndexOutsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/geometryhelper.hh:171
auto intersection(const SmallLocalIndexType localFacetIdx, const Element &element) const
Returns an element's intersection based on the local facet index.
Definition: discretization/facecentered/staggered/geometryhelper.hh:137
static constexpr int lateralOrthogonalScvfLocalIndex(const SmallLocalIndexType ownLocalScvfIndex)
Return the local index of a lateral orthogonal scvf.
Definition: discretization/facecentered/staggered/geometryhelper.hh:48
static constexpr auto dim
Definition: discretization/facecentered/staggered/geometryhelper.hh:35
Defines the index types used for grid and local indices.
Definition: adapt.hh:17
Type trait to determine if a grid supports concave corners (e.g. by cutting out a hole from the domai...
Definition: gridsupportsconcavecorners.hh:31
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:29