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