version 3.8
discretization/facecentered/diamond/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_DIAMOND_GEOMETRY_HELPER_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_GEOMETRY_HELPER_HH
14
15#include <array>
16
17#include <dune/common/reservedvector.hh>
18#include <dune/common/fvector.hh>
19#include <dune/geometry/multilineargeometry.hh>
20#include <dune/geometry/type.hh>
21
22#include <dumux/common/math.hh>
23
24namespace Dumux {
25
27template <class ct>
28struct FCDiamondMLGeometryTraits : public Dune::MultiLinearGeometryTraits<ct>
29{
30 // we use static vectors to store the corners as we know
31 // the maximum number of corners in advance (2^dim)
32 template< int mydim, int cdim >
34 {
35 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)>;
36 };
37};
38
40
41template<Dune::GeometryType::Id gt>
43
44template<>
45struct ScvCorners<Dune::GeometryTypes::triangle>
46{
47 static constexpr Dune::GeometryType type()
48 { return Dune::GeometryTypes::triangle; }
49
50 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
51 static constexpr std::array<std::array<Key, 3>, 3> keys = {{
52 { Key{0, 2}, Key{1, 2}, Key{0, 0} },
53 { Key{2, 2}, Key{0, 2}, Key{0, 0} },
54 { Key{1, 2}, Key{2, 2}, Key{0, 0} }
55 }};
56};
57
58template<>
59struct ScvCorners<Dune::GeometryTypes::quadrilateral>
60{
61 static constexpr Dune::GeometryType type()
62 { return Dune::GeometryTypes::triangle; }
63
64 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
65 static constexpr std::array<std::array<Key, 3>, 4> keys = {{
66 { Key{2, 2}, Key{0, 2}, Key{0, 0} },
67 { Key{1, 2}, Key{3, 2}, Key{0, 0} },
68 { Key{0, 2}, Key{1, 2}, Key{0, 0} },
69 { Key{3, 2}, Key{2, 2}, Key{0, 0} }
70 }};
71};
72
73template<>
74struct ScvCorners<Dune::GeometryTypes::tetrahedron>
75{
76 static constexpr Dune::GeometryType type()
77 { return Dune::GeometryTypes::tetrahedron; }
78
79 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
80 static constexpr std::array<std::array<Key, 4>, 4> keys = {{
81 { Key{0, 3}, Key{1, 3}, Key{2, 3}, Key{0, 0} },
82 { Key{1, 3}, Key{0, 3}, Key{3, 3}, Key{0, 0} },
83 { Key{0, 3}, Key{2, 3}, Key{3, 3}, Key{0, 0} },
84 { Key{2, 3}, Key{1, 3}, Key{3, 3}, Key{0, 0} }
85 }};
86};
87
88template<>
89struct ScvCorners<Dune::GeometryTypes::hexahedron>
90{
91 static constexpr Dune::GeometryType type()
92 { return Dune::GeometryTypes::pyramid; }
93
94 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
95 static constexpr std::array<std::array<Key, 5>, 6> keys = {{
96 { Key{4, 3}, Key{0, 3}, Key{6, 3}, Key{2, 3}, Key{0, 0} },
97 { Key{1, 3}, Key{5, 3}, Key{3, 3}, Key{7, 3}, Key{0, 0} },
98 { Key{4, 3}, Key{5, 3}, Key{0, 3}, Key{1, 3}, Key{0, 0} },
99 { Key{2, 3}, Key{3, 3}, Key{6, 3}, Key{7, 3}, Key{0, 0} },
100 { Key{0, 3}, Key{1, 3}, Key{2, 3}, Key{3, 3}, Key{0, 0} },
101 { Key{6, 3}, Key{7, 3}, Key{4, 3}, Key{5, 3}, Key{0, 0} }
102 }};
103};
104
105template<Dune::GeometryType::Id gt>
107
108template<>
109struct ScvfCorners<Dune::GeometryTypes::triangle>
110{
111 static constexpr Dune::GeometryType type()
112 { return Dune::GeometryTypes::line; }
113
114 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
115 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
116 { Key{0, 2}, Key{0, 0} },
117 { Key{1, 2}, Key{0, 0} },
118 { Key{2, 2}, Key{0, 0} }
119 }};
120};
121
122template<>
123struct ScvfCorners<Dune::GeometryTypes::quadrilateral>
124{
125 static constexpr Dune::GeometryType type()
126 { return Dune::GeometryTypes::line; }
127
128 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
129 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
130 { Key{0, 2}, Key{0, 0} },
131 { Key{1, 2}, Key{0, 0} },
132 { Key{2, 2}, Key{0, 0} },
133 { Key{3, 2}, Key{0, 0} }
134 }};
135};
136
137template<>
138struct ScvfCorners<Dune::GeometryTypes::tetrahedron>
139{
140 static constexpr Dune::GeometryType type()
141 { return Dune::GeometryTypes::triangle; }
142
143 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
144 static constexpr std::array<std::array<Key, 3>, 6> keys = {{
145 { Key{0, 3}, Key{1, 3}, Key{0, 0} },
146 { Key{2, 3}, Key{0, 3}, Key{0, 0} },
147 { Key{1, 3}, Key{2, 3}, Key{0, 0} },
148 { Key{0, 3}, Key{3, 3}, Key{0, 0} },
149 { Key{1, 3}, Key{3, 3}, Key{0, 0} },
150 { Key{2, 3}, Key{3, 3}, Key{0, 0} }
151 }};
152};
153
154template<>
155struct ScvfCorners<Dune::GeometryTypes::hexahedron>
156{
157 static constexpr Dune::GeometryType type()
158 { return Dune::GeometryTypes::triangle; }
159
160 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
161 static constexpr std::array<std::array<Key, 3>, 12> keys = {{
162 { Key{0, 3}, Key{4, 3}, Key{0, 0} },
163 { Key{1, 3}, Key{5, 3}, Key{0, 0} },
164 { Key{2, 3}, Key{6, 3}, Key{0, 0} },
165 { Key{3, 3}, Key{7, 3}, Key{0, 0} },
166 { Key{0, 3}, Key{2, 3}, Key{0, 0} },
167 { Key{1, 3}, Key{3, 3}, Key{0, 0} },
168 { Key{0, 3}, Key{1, 3}, Key{0, 0} },
169 { Key{2, 3}, Key{3, 3}, Key{0, 0} },
170 { Key{4, 3}, Key{6, 3}, Key{0, 0} },
171 { Key{5, 3}, Key{7, 3}, Key{0, 0} },
172 { Key{4, 3}, Key{5, 3}, Key{0, 0} },
173 { Key{6, 3}, Key{7, 3}, Key{0, 0} }
174 }};
175};
176
177// convert key array to global corner storage
178template<class S, class Geo, class KeyArray, std::size_t... I>
179S keyToCornerStorageImpl(const Geo& geo, const KeyArray& key, std::index_sequence<I...>)
180{
181 using Dune::referenceElement;
182 const auto ref = referenceElement(geo);
183 // key is a pair of a local sub-entity index (first) and the sub-entity's codim (second)
184 return { geo.global(ref.position(key[I].first, key[I].second))... };
185}
186
187// convert key array to global corner storage
188template<class S, class Geo, class T, std::size_t N, class Indices = std::make_index_sequence<N>>
189S keyToCornerStorage(const Geo& geo, const std::array<T, N>& key)
190{
191 return keyToCornerStorageImpl<S>(geo, key, Indices{});
192}
193
194// boundary corners for the i-th facet
195template<class S, class Geo, std::size_t... ii>
196S boundaryCornerStorageImpl(const Geo& geo, unsigned int i, std::index_sequence<ii...>)
197{
198 using Dune::referenceElement;
199 const auto ref = referenceElement(geo);
200 // simply the vertices of the facet i
201 return { geo.global(ref.position(ref.subEntity(i, 1, ii, Geo::mydimension), Geo::mydimension))... };
202}
203
204// boundary corners for the i-th facet
205template<class S, std::size_t numCorners, class Geo>
206S boundaryCornerStorage(const Geo& geo, unsigned int i)
207{
208 return boundaryCornerStorageImpl<S>(geo, i, std::make_index_sequence<numCorners>{});
209}
210
211template<class IndexType, Dune::GeometryType::Id gt>
213
214template<class IndexType>
215struct InsideOutsideScv<IndexType, Dune::GeometryTypes::triangle>
216{
217 static constexpr std::array<std::array<IndexType, 2>, 3> pairs = {{
218 {0, 1}, {0, 2}, {1, 2}
219 }};
220};
221
222template<class IndexType>
223struct InsideOutsideScv<IndexType, Dune::GeometryTypes::quadrilateral>
224{
225 static constexpr std::array<std::array<IndexType, 2>, 4> pairs = {{
226 {0, 2}, {1, 2}, {0, 3}, {1, 3}
227 }};
228};
229
230template<class IndexType>
231struct InsideOutsideScv<IndexType, Dune::GeometryTypes::tetrahedron>
232{
233 static constexpr std::array<std::array<IndexType, 2>, 6> pairs = {{
234 {0, 1}, {0, 2}, {0, 3}, {1, 2}, {1, 3}, {2, 3},
235 }};
236};
237
238template<class IndexType>
239struct InsideOutsideScv<IndexType, Dune::GeometryTypes::hexahedron>
240{
241 static constexpr std::array<std::array<IndexType, 2>, 12> pairs = {{
242 {0, 2}, {1, 2}, {0, 3}, {1, 3}, {0, 4}, {1, 4},
243 {2, 4}, {3, 4}, {0, 5}, {1, 5}, {2, 5}, {3, 5}
244 }};
245};
246
247} // end namespace Detail::FCDiamond
248
253template <class GridView, class ScvType, class ScvfType>
255{
256 using Element = typename GridView::template Codim<0>::Entity;
257
258 static constexpr auto dim = GridView::dimension;
259 static constexpr auto dimWorld = GridView::dimensionworld;
260
261 using ScvCornerStorage = typename ScvType::Traits::CornerStorage;
262 using LocalIndexType = typename ScvType::Traits::LocalIndexType;
263 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
264
265 // for the normal
266 using ctype = typename GridView::ctype;
267 using GlobalPosition = typename Dune::FieldVector<ctype, GridView::dimensionworld>;
268
269public:
270 DiamondGeometryHelper(const typename Element::Geometry& geo)
271 : geo_(geo)
272 {}
273
275 ScvCornerStorage getScvCorners(unsigned int localFacetIndex) const
276 {
277 const auto type = geo_.type();
278 if (type == Dune::GeometryTypes::triangle)
279 {
281 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localFacetIndex]);
282 }
283 else if (type == Dune::GeometryTypes::quadrilateral)
284 {
286 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localFacetIndex]);
287 }
288 else if (type == Dune::GeometryTypes::tetrahedron)
289 {
291 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localFacetIndex]);
292 }
293 else if (type == Dune::GeometryTypes::hexahedron)
294 {
296 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localFacetIndex]);
297 }
298 else
299 DUNE_THROW(Dune::NotImplemented, "Scv geometries for type " << type);
300 }
301
303 ScvfCornerStorage getScvfCorners(unsigned int localEdgeIndex) const
304 {
305 const auto type = geo_.type();
306 if (type == Dune::GeometryTypes::triangle)
307 {
309 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localEdgeIndex]);
310 }
311 else if (type == Dune::GeometryTypes::quadrilateral)
312 {
314 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localEdgeIndex]);
315 }
316 else if (type == Dune::GeometryTypes::tetrahedron)
317 {
319 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localEdgeIndex]);
320 }
321 else if (type == Dune::GeometryTypes::hexahedron)
322 {
324 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localEdgeIndex]);
325 }
326 else
327 DUNE_THROW(Dune::NotImplemented, "Scvf geometries for type " << type);
328 }
329
331 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex) const
332 {
333 const auto type = geo_.type();
334 if (type == Dune::GeometryTypes::triangle || type == Dune::GeometryTypes::quadrilateral)
335 return Detail::FCDiamond::boundaryCornerStorage<ScvfCornerStorage, 2>(geo_, localFacetIndex);
336 else if (type == Dune::GeometryTypes::tetrahedron)
337 return Detail::FCDiamond::boundaryCornerStorage<ScvfCornerStorage, 3>(geo_, localFacetIndex);
338 else if (type == Dune::GeometryTypes::hexahedron)
339 return Detail::FCDiamond::boundaryCornerStorage<ScvfCornerStorage, 4>(geo_, localFacetIndex);
340 else
341 DUNE_THROW(Dune::NotImplemented, "Boundary scvf geometries for type " << type);
342 }
343
344 GlobalPosition facetCenter(unsigned int localFacetIndex) const
345 {
346 return geo_.global(referenceElement(geo_).position(localFacetIndex, 1));
347 }
348
349 std::array<LocalIndexType, 2> getInsideOutsideScvForScvf(unsigned int localEdgeIndex)
350 {
351 const auto type = geo_.type();
352 if (type == Dune::GeometryTypes::triangle)
354 else if (type == Dune::GeometryTypes::quadrilateral)
356 else if (type == Dune::GeometryTypes::tetrahedron)
358 else if (type == Dune::GeometryTypes::hexahedron)
360 else
361 DUNE_THROW(Dune::NotImplemented, "Inside outside scv pairs for type " << type);
362 }
363
365 std::size_t numInteriorScvf()
366 {
367 return referenceElement(geo_).size(2);
368 }
369
371 std::size_t numScv()
372 {
373 return referenceElement(geo_).size(1);
374 }
375
376 template<int d = dimWorld, std::enable_if_t<(d==3), int> = 0>
377 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
378 {
379 auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]);
380 normal /= normal.two_norm();
381
382 const auto ref = referenceElement(geo_);
383 const auto v = facetCenter_(scvPair[1], ref) - facetCenter_(scvPair[0], ref);
384
385 const auto s = v*normal;
386 if (std::signbit(s))
387 normal *= -1;
388
389 return normal;
390 }
391
392 template<int d = dimWorld, std::enable_if_t<(d==2), int> = 0>
393 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
394 {
396 const auto t = p[1] - p[0];
397 GlobalPosition normal({-t[1], t[0]});
398 normal /= normal.two_norm();
399
400 const auto ref = referenceElement(geo_);
401 const auto v = facetCenter_(scvPair[1], ref) - facetCenter_(scvPair[0], ref);
402
403 const auto s = v*normal;
404 if (std::signbit(s))
405 normal *= -1;
406
407 return normal;
408 }
409
410 const typename Element::Geometry& elementGeometry() const
411 { return geo_; }
412
413private:
414 template<class RefElement>
415 GlobalPosition facetCenter_(unsigned int localFacetIndex, const RefElement& ref) const
416 {
417 return geo_.global(ref.position(localFacetIndex, 1));
418 }
419
420 const typename Element::Geometry& geo_;
421};
422
423} // end namespace Dumux
424
425#endif
Helper class to construct SCVs and SCVFs for the diamond scheme.
Definition: discretization/facecentered/diamond/geometryhelper.hh:255
GlobalPosition facetCenter(unsigned int localFacetIndex) const
Definition: discretization/facecentered/diamond/geometryhelper.hh:344
ScvfCornerStorage getScvfCorners(unsigned int localEdgeIndex) const
Create a corner storage with the scvf corners for a given edge (codim-2) index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:303
std::size_t numInteriorScvf()
number of interior sub control volume faces (number of codim-2 entities)
Definition: discretization/facecentered/diamond/geometryhelper.hh:365
std::array< LocalIndexType, 2 > getInsideOutsideScvForScvf(unsigned int localEdgeIndex)
Definition: discretization/facecentered/diamond/geometryhelper.hh:349
ScvCornerStorage getScvCorners(unsigned int localFacetIndex) const
Create a corner storage with the scv corners for a given face (codim-1) index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:275
DiamondGeometryHelper(const typename Element::Geometry &geo)
Definition: discretization/facecentered/diamond/geometryhelper.hh:270
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition: discretization/facecentered/diamond/geometryhelper.hh:377
const Element::Geometry & elementGeometry() const
Definition: discretization/facecentered/diamond/geometryhelper.hh:410
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex) const
Create the sub control volume face geometries on the boundary for a given face index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:331
std::size_t numScv()
number of sub control volumes (number of codim-1 entities)
Definition: discretization/facecentered/diamond/geometryhelper.hh:371
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:642
Define some often used mathematical functions.
S keyToCornerStorage(const Geo &geo, const std::array< T, N > &key)
Definition: discretization/facecentered/diamond/geometryhelper.hh:189
S keyToCornerStorageImpl(const Geo &geo, const KeyArray &key, std::index_sequence< I... >)
Definition: discretization/facecentered/diamond/geometryhelper.hh:179
S boundaryCornerStorageImpl(const Geo &geo, unsigned int i, std::index_sequence< ii... >)
Definition: discretization/facecentered/diamond/geometryhelper.hh:196
S boundaryCornerStorage(const Geo &geo, unsigned int i)
Definition: discretization/facecentered/diamond/geometryhelper.hh:206
CVFE< CVFEMethods::CR_RT > FCDiamond
Definition: method.hh:101
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
Definition: discretization/facecentered/diamond/geometryhelper.hh:212
Definition: discretization/facecentered/diamond/geometryhelper.hh:90
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:91
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:94
Definition: discretization/facecentered/diamond/geometryhelper.hh:60
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:61
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:64
Definition: discretization/facecentered/diamond/geometryhelper.hh:75
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:79
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:76
Definition: discretization/facecentered/diamond/geometryhelper.hh:46
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:47
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:50
Definition: discretization/facecentered/diamond/geometryhelper.hh:42
Definition: discretization/facecentered/diamond/geometryhelper.hh:156
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:157
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:160
Definition: discretization/facecentered/diamond/geometryhelper.hh:124
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:125
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:128
Definition: discretization/facecentered/diamond/geometryhelper.hh:139
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:143
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:140
Definition: discretization/facecentered/diamond/geometryhelper.hh:110
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:114
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:111
Definition: discretization/facecentered/diamond/geometryhelper.hh:106
Definition: discretization/facecentered/diamond/geometryhelper.hh:34
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<< mydim)> Type
Definition: discretization/facecentered/diamond/geometryhelper.hh:35
Traits for an efficient corner storage for fc diamond method.
Definition: discretization/facecentered/diamond/geometryhelper.hh:29