3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_FACECENTERED_STAGGERED_GEOMETRY_HELPER_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_GEOMETRY_HELPER_HH
26
27#include <array>
28#include <cassert>
29#include <dune/common/exceptions.hh>
32
33namespace Dumux {
38template<class GridView>
40{
41 using GridIndexType = typename IndexTraits<GridView>::GridIndex;
42 using SmallLocalIndexType = typename IndexTraits<GridView>::SmallLocalIndex;
43 using Element = typename GridView::template Codim<0>::Entity;
44 using Facet = typename GridView::template Codim<1>::Entity;
45
46public:
47 static constexpr auto dim = GridView::Grid::dimension;
48 static constexpr auto numElementFaces = dim * 2;
49 static constexpr auto numLateralFacesPerScv = 2 * (dim - 1);
50
51 FaceCenteredStaggeredGeometryHelper(const GridView& gridView) : gridView_(gridView) {}
52
54 static constexpr SmallLocalIndexType localOppositeIdx(const SmallLocalIndexType ownLocalFaceIndex)
55 {
56 return isOdd_(ownLocalFaceIndex) ? (ownLocalFaceIndex - 1) : (ownLocalFaceIndex + 1);
57 }
58
60 static constexpr int lateralOrthogonalScvfLocalIndex(const SmallLocalIndexType ownLocalScvfIndex)
61 {
62 if constexpr(GridView::Grid::dimension == 1)
63 {
64 assert(false && "There are no lateral scvfs in 1D");
65 return -1;
66 }
67
68 if constexpr (GridView::Grid::dimension == 2)
69 {
70 switch (ownLocalScvfIndex)
71 {
72 case 1: return 7;
73 case 7: return 1;
74 case 2: return 10;
75 case 10: return 2;
76 case 4: return 8;
77 case 8: return 4;
78 case 5: return 11;
79 case 11: return 5;
80 default:
81 {
82 assert(false && "No lateral orthogonal scvf found");
83 return -1;
84 }
85 }
86 }
87 else
88 {
89 switch (ownLocalScvfIndex)
90 {
91 case 1: return 11;
92 case 11: return 1;
93 case 2: return 16;
94 case 16: return 2;
95 case 3: return 21;
96 case 21: return 3;
97 case 4: return 26;
98 case 26: return 4;
99 case 6: return 12;
100 case 12: return 6;
101 case 7: return 17;
102 case 17: return 7;
103 case 8: return 22;
104 case 22: return 8;
105 case 9: return 27;
106 case 27: return 9;
107 case 13: return 23;
108 case 23: return 13;
109 case 14: return 28;
110 case 28: return 14;
111 case 18: return 24;
112 case 24: return 18;
113 case 19: return 29;
114 case 29: return 19;
115
116 default:
117 {
118 assert(false && "No lateral orthogonal scvf found");
119 return -1;
120 }
121 }
122 }
123 }
124
126 static constexpr auto localLaterFaceIndices(const SmallLocalIndexType ownLocalFaceIndex)
127 {
128 constexpr auto table = []
129 {
130 using Table = std::array<std::array<SmallLocalIndexType, numLateralFacesPerScv>, numElementFaces>;
131 if constexpr (dim == 1)
132 return Table{};
133 else if constexpr (dim == 2)
134 return Table {{ {2,3}, {2,3}, {0,1}, {0,1} }};
135 else
136 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} }};
137 }();
138
139 return table[ownLocalFaceIndex];
140 }
141
143 static Facet facet(const SmallLocalIndexType localFacetIdx, const Element& element)
144 {
145 return element.template subEntity <1> (localFacetIdx);
146 }
147
149 auto intersection(const SmallLocalIndexType localFacetIdx, const Element& element) const
150 {
151 for (const auto& intersection : intersections(gridView(), element))
152 {
153 if (intersection.indexInInside() == localFacetIdx)
154 return intersection;
155 }
156 DUNE_THROW(Dune::InvalidStateException, "localFacetIdx " << localFacetIdx << " out of range");
157 }
158
160 template<class FVElementGeometry, class SubControlVolumeFace>
161 static bool scvfIntegrationPointInConcaveCorner(const FVElementGeometry& fvGeometry, const SubControlVolumeFace& scvf)
162 {
163 assert (scvf.isLateral());
164
165 using Grid = typename FVElementGeometry::GridGeometry::Grid;
167 return false;
168 else
169 {
170 if (scvf.boundary())
171 return false;
172
173 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
174 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
175 int onBoundaryCounter = 0;
176 onBoundaryCounter += static_cast<int>(insideScv.boundary());
177 onBoundaryCounter += static_cast<int>(outsideScv.boundary());
178 return onBoundaryCounter == 1;
179 }
180 }
181
182 template<class SubControlVolumeFace, class SubControlVolume>
183 static SmallLocalIndexType localIndexOutsideScvfWithSameIntegrationPoint(const SubControlVolumeFace& scvf,
184 const SubControlVolume& scv)
185 {
186 // In 2D there are 3 non-boundary faces per scv. In 3D, there are 5.
187 // This number of scvfs per scv is used as an offset to find the indexes in the outside half-scv.
188 const SmallLocalIndexType offset = (dim == 2) ? 3 : 5;
189 // For half-scvs with odd indexes, the outside half-scv has scvf local indexes with + offset.
190 // For half-scvs with even indexes, the outside half-scv has scvf local indexes have a - offset.
191 return isOdd_(scv.indexInElement()) ? scvf.localIndex() - offset : scvf.localIndex() + offset;
192 }
193
194 const GridView& gridView() const
195 { return gridView_; }
196
197private:
198
199 static constexpr bool isOdd_(int number)
200 { return number % 2; }
201
202 GridView gridView_;
203};
204
205} // end namespace Dumux
206
207#endif
Defines the index types used for grid and local indices.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
Face centered staggered geometry helper.
Definition: discretization/facecentered/staggered/geometryhelper.hh:40
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:143
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:126
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:161
static constexpr SmallLocalIndexType localOppositeIdx(const SmallLocalIndexType ownLocalFaceIndex)
Returns the local index of the opposing face.
Definition: discretization/facecentered/staggered/geometryhelper.hh:54
static constexpr auto numLateralFacesPerScv
Definition: discretization/facecentered/staggered/geometryhelper.hh:49
const GridView & gridView() const
Definition: discretization/facecentered/staggered/geometryhelper.hh:194
static constexpr auto numElementFaces
Definition: discretization/facecentered/staggered/geometryhelper.hh:48
FaceCenteredStaggeredGeometryHelper(const GridView &gridView)
Definition: discretization/facecentered/staggered/geometryhelper.hh:51
static SmallLocalIndexType localIndexOutsideScvfWithSameIntegrationPoint(const SubControlVolumeFace &scvf, const SubControlVolume &scv)
Definition: discretization/facecentered/staggered/geometryhelper.hh:183
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:149
static constexpr int lateralOrthogonalScvfLocalIndex(const SmallLocalIndexType ownLocalScvfIndex)
Return the local index of a lateral orthogonal scvf.
Definition: discretization/facecentered/staggered/geometryhelper.hh:60
static constexpr auto dim
Definition: discretization/facecentered/staggered/geometryhelper.hh:47
Type trait to determine if a grid supports concave corners (e.g. by cutting out a hole from the domai...
Definition: gridsupportsconcavecorners.hh:43