3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
boxgeometryhelper.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 *****************************************************************************/
25#ifndef DUMUX_DISCRETIZATION_BOX_GEOMETRY_HELPER_HH
26#define DUMUX_DISCRETIZATION_BOX_GEOMETRY_HELPER_HH
27
28#include <array>
29
30#include <dune/common/exceptions.hh>
31
32#include <dune/geometry/type.hh>
33#include <dune/geometry/typeindex.hh>
34#include <dune/geometry/referenceelements.hh>
35#include <dune/geometry/multilineargeometry.hh>
36
37#include <dumux/common/math.hh>
38
39namespace Dumux {
40
42template <class ct>
43struct BoxMLGeometryTraits : public Dune::MultiLinearGeometryTraits<ct>
44{
45 // we use static vectors to store the corners as we know
46 // the number of corners in advance (2^(mydim) corners (1<<(mydim))
47 template< int mydim, int cdim >
49 {
50 using Type = std::array< Dune::FieldVector< ct, cdim >, (1<<(mydim)) >;
51 };
52
53 // we know all scvfs will have the same geometry type
54 template< int mydim >
56 {
57 static const bool v = true;
58 static const unsigned int topologyId = Dune::GeometryTypes::cube(mydim).id();
59 };
60};
61
62namespace Detail::Box {
63
64template<Dune::GeometryType::Id gt>
66
67template<>
69{
70 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
71 static constexpr std::array<std::array<Key, 2>, 2> keys = {{
72 { Key{0, 1}, Key{0, 0} },
73 { Key{1, 1}, Key{0, 0} }
74 }};
75};
76
77template<>
78struct ScvCorners<Dune::GeometryTypes::triangle>
79{
80 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
81 static constexpr std::array<std::array<Key, 4>, 3> keys = {{
82 { Key{0, 2}, Key{0, 1}, Key{1, 1}, Key{0, 0} },
83 { Key{1, 2}, Key{2, 1}, Key{0, 1}, Key{0, 0} },
84 { Key{2, 2}, Key{1, 1}, Key{2, 1}, Key{0, 0} }
85 }};
86};
87
88template<>
89struct ScvCorners<Dune::GeometryTypes::quadrilateral>
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, 2}, Key{2, 1}, Key{0, 1}, Key{0, 0} },
94 { Key{1, 2}, Key{1, 1}, Key{2, 1}, Key{0, 0} },
95 { Key{2, 2}, Key{0, 1}, Key{3, 1}, Key{0, 0} },
96 { Key{3, 2}, Key{3, 1}, Key{1, 1}, Key{0, 0} }
97 }};
98};
99
100template<>
101struct ScvCorners<Dune::GeometryTypes::tetrahedron>
102{
103 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
104 static constexpr std::array<std::array<Key, 8>, 4> keys = {{
105 { Key{0, 3}, Key{0, 2}, Key{1, 2}, Key{0, 1}, Key{3, 2}, Key{1, 1}, Key{2, 1}, Key{0, 0} },
106 { Key{1, 3}, Key{2, 2}, Key{0, 2}, Key{0, 1}, Key{4, 2}, Key{3, 1}, Key{1, 1}, Key{0, 0} },
107 { Key{2, 3}, Key{1, 2}, Key{2, 2}, Key{0, 1}, Key{5, 2}, Key{2, 1}, Key{3, 1}, Key{0, 0} },
108 { Key{3, 3}, Key{3, 2}, Key{5, 2}, Key{2, 1}, Key{4, 2}, Key{1, 1}, Key{3, 1}, Key{0, 0} }
109 }};
110};
111
112template<>
113struct ScvCorners<Dune::GeometryTypes::prism>
114{
115 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
116 static constexpr std::array<std::array<Key, 8>, 6> keys = {{
117 { Key{0, 3}, Key{3, 2}, Key{4, 2}, Key{3, 1}, Key{0, 2}, Key{0, 1}, Key{1, 1}, Key{0, 0} },
118 { Key{1, 3}, Key{5, 2}, Key{3, 2}, Key{3, 1}, Key{1, 2}, Key{2, 1}, Key{0, 1}, Key{0, 0} },
119 { Key{2, 3}, Key{4, 2}, Key{5, 2}, Key{3, 1}, Key{2, 2}, Key{1, 1}, Key{2, 1}, Key{0, 0} },
120 { Key{3, 3}, Key{7, 2}, Key{6, 2}, Key{4, 1}, Key{0, 2}, Key{1, 1}, Key{0, 1}, Key{0, 0} },
121 { Key{4, 3}, Key{6, 2}, Key{8, 2}, Key{4, 1}, Key{1, 2}, Key{0, 1}, Key{2, 1}, Key{0, 0} },
122 { Key{5, 3}, Key{8, 2}, Key{7, 2}, Key{4, 1}, Key{2, 2}, Key{2, 1}, Key{1, 1}, Key{0, 0} }
123 }};
124};
125
126template<>
127struct ScvCorners<Dune::GeometryTypes::hexahedron>
128{
129 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
130 static constexpr std::array<std::array<Key, 8>, 8> keys = {{
131 { Key{0, 3}, Key{6, 2}, Key{4, 2}, Key{4, 1}, Key{0, 2}, Key{2, 1}, Key{0, 1}, Key{0, 0} },
132 { Key{1, 3}, Key{5, 2}, Key{6, 2}, Key{4, 1}, Key{1, 2}, Key{1, 1}, Key{2, 1}, Key{0, 0} },
133 { Key{2, 3}, Key{4, 2}, Key{7, 2}, Key{4, 1}, Key{2, 2}, Key{0, 1}, Key{3, 1}, Key{0, 0} },
134 { Key{3, 3}, Key{7, 2}, Key{5, 2}, Key{4, 1}, Key{3, 2}, Key{3, 1}, Key{1, 1}, Key{0, 0} },
135 { Key{4, 3}, Key{8, 2}, Key{10, 2}, Key{5, 1}, Key{0, 2}, Key{0, 1}, Key{2, 1}, Key{0, 0} },
136 { Key{5, 3}, Key{10, 2}, Key{9, 2}, Key{5, 1}, Key{1, 2}, Key{2, 1}, Key{1, 1}, Key{0, 0} },
137 { Key{6, 3}, Key{11, 2}, Key{8, 2}, Key{5, 1}, Key{2, 2}, Key{3, 1}, Key{0, 1}, Key{0, 0} },
138 { Key{7, 3}, Key{9, 2}, Key{11, 2}, Key{5, 1}, Key{3, 2}, Key{1, 1}, Key{3, 1}, Key{0, 0} }
139 }};
140};
141
142template<Dune::GeometryType::Id gt>
144
145template<>
147{
148 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
149 static constexpr std::array<std::array<Key, 1>, 1> keys = {{
150 { Key{0, 0} }
151 }};
152};
153
154template<>
155struct ScvfCorners<Dune::GeometryTypes::triangle>
156{
157 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
158 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
159 { Key{0, 0}, Key{0, 1} },
160 { Key{1, 1}, Key{0, 0} },
161 { Key{0, 0}, Key{2, 1} }
162 }};
163};
164
165template<>
166struct ScvfCorners<Dune::GeometryTypes::quadrilateral>
167{
168 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
169 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
170 { Key{0, 1}, Key{0, 0} },
171 { Key{0, 0}, Key{1, 1} },
172 { Key{0, 0}, Key{2, 1} },
173 { Key{3, 1}, Key{0, 0} }
174 }};
175};
176
177template<>
178struct ScvfCorners<Dune::GeometryTypes::tetrahedron>
179{
180 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
181 static constexpr std::array<std::array<Key, 4>, 6> keys = {{
182 { Key{0, 2}, Key{0, 1}, Key{1, 1}, Key{0, 0} },
183 { Key{0, 1}, Key{1, 2}, Key{0, 0}, Key{2, 1} },
184 { Key{2, 2}, Key{0, 1}, Key{3, 1}, Key{0, 0} },
185 { Key{2, 1}, Key{3, 2}, Key{0, 0}, Key{1, 1} },
186 { Key{3, 1}, Key{0, 0}, Key{4, 2}, Key{1, 1} },
187 { Key{5, 2}, Key{2, 1}, Key{3, 1}, Key{0, 0} }
188 }};
189};
190
191template<>
192struct ScvfCorners<Dune::GeometryTypes::prism>
193{
194 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
195 static constexpr std::array<std::array<Key, 4>, 9> keys = {{
196 { Key{0, 2}, Key{0, 1}, Key{1, 1}, Key{0, 0} },
197 { Key{1, 2}, Key{2, 1}, Key{0, 1}, Key{0, 0} },
198 { Key{2, 2}, Key{1, 1}, Key{2, 1}, Key{0, 0} },
199 { Key{3, 2}, Key{0, 1}, Key{3, 1}, Key{0, 0} },
200 { Key{4, 2}, Key{3, 1}, Key{1, 1}, Key{0, 0} },
201 { Key{5, 2}, Key{2, 1}, Key{3, 1}, Key{0, 0} },
202 { Key{6, 2}, Key{4, 1}, Key{0, 1}, Key{0, 0} },
203 { Key{7, 2}, Key{1, 1}, Key{4, 1}, Key{0, 0} },
204 { Key{8, 2}, Key{4, 1}, Key{2, 1}, Key{0, 0} }
205 }};
206};
207
208template<>
209struct ScvfCorners<Dune::GeometryTypes::hexahedron>
210{
211 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
212 static constexpr std::array<std::array<Key, 4>, 12> keys = {{
213 { Key{0, 1}, Key{0, 2}, Key{0, 0}, Key{2, 1} },
214 { Key{1, 1}, Key{0, 0}, Key{1, 2}, Key{2, 1} },
215 { Key{3, 1}, Key{2, 2}, Key{0, 0}, Key{0, 1} },
216 { Key{3, 2}, Key{3, 1}, Key{1, 1}, Key{0, 0} },
217 { Key{4, 1}, Key{4, 2}, Key{0, 0}, Key{0, 1} },
218 { Key{5, 2}, Key{4, 1}, Key{1, 1}, Key{0, 0} },
219 { Key{6, 2}, Key{4, 1}, Key{2, 1}, Key{0, 0} },
220 { Key{4, 1}, Key{7, 2}, Key{0, 0}, Key{3, 1} },
221 { Key{0, 0}, Key{0, 1}, Key{5, 1}, Key{8, 2} },
222 { Key{9, 2}, Key{1, 1}, Key{5, 1}, Key{0, 0} },
223 { Key{10, 2}, Key{2, 1}, Key{5, 1}, Key{0, 0} },
224 { Key{11, 2}, Key{5, 1}, Key{3, 1}, Key{0, 0} }
225 }};
226};
227
228// convert key array to global corner storage
229template<class S, class Geo, class KeyArray, std::size_t... I>
230S keyToCornerStorageImpl(const Geo& geo, const KeyArray& key, std::index_sequence<I...>)
231{
232 using Dune::referenceElement;
233 const auto ref = referenceElement(geo);
234 // key is a pair of a local sub-entity index (first) and the sub-entity's codim (second)
235 return { geo.global(ref.position(key[I].first, key[I].second))... };
236}
237
238// convert key array to global corner storage
239template<class S, class Geo, class T, std::size_t N, class Indices = std::make_index_sequence<N>>
240S keyToCornerStorage(const Geo& geo, const std::array<T, N>& key)
241{
242 return keyToCornerStorageImpl<S>(geo, key, Indices{});
243}
244
245// convert key array to global corner storage
246// for the i-th sub-entity of codim c (e.g. the i-th facet/codim-1-entity for boundaries)
247template<class S, class Geo, class KeyArray, std::size_t... I>
248S subEntityKeyToCornerStorageImpl(const Geo& geo, unsigned int i, unsigned int c, const KeyArray& key, std::index_sequence<I...>)
249{
250 using Dune::referenceElement;
251 const auto ref = referenceElement(geo);
252 // subEntity gives the subEntity number with respect to the codim-0 reference element
253 // key is a pair of a local sub-entity index (first) and the sub-entity's codim (second) but here w.r.t. the sub-entity i/c
254 return { geo.global(ref.position(ref.subEntity(i, c, key[I].first, c+key[I].second), c+key[I].second))... };
255}
256
257// convert key array to global corner storage
258// for the i-th sub-entity of codim c (e.g. the i-th facet/codim-1-entity for boundaries)
259template<class S, class Geo, class T, std::size_t N, class Indices = std::make_index_sequence<N>>
260S subEntityKeyToCornerStorage(const Geo& geo, unsigned int i, unsigned int c, const std::array<T, N>& key)
261{
262 return subEntityKeyToCornerStorageImpl<S>(geo, i, c, key, Indices{});
263}
264
265} // end namespace Detail::Box
266
268template<class GridView, int dim, class ScvType, class ScvfType>
270
272template <class GridView, class ScvType, class ScvfType>
273class BoxGeometryHelper<GridView, 1, ScvType, ScvfType>
274{
275private:
276 using Scalar = typename GridView::ctype;
277 using GlobalPosition = typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
278 using ScvCornerStorage = typename ScvType::Traits::CornerStorage;
279 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
280 using ScvGeometry = typename ScvType::Traits::Geometry;
281 using ScvfGeometry = typename ScvfType::Traits::Geometry;
282
283 using Element = typename GridView::template Codim<0>::Entity;
284 using Intersection = typename GridView::Intersection;
285
286 static constexpr int dim = 1;
287public:
288
289 BoxGeometryHelper(const typename Element::Geometry& geometry)
290 : geo_(geometry)
291 {}
292
294 ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
295 {
297 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localScvIdx]);
298 }
299
300 ScvGeometry scvGeometry(unsigned int localScvIdx) const
301 {
302 return { Dune::GeometryTypes::line, getScvCorners(localScvIdx) };
303 }
304
306 ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
307 {
309 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localScvfIdx]);
310 }
311
313 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex,
314 unsigned int) const
315 {
316 return ScvfCornerStorage{{ geo_.corner(localFacetIndex) }};
317 }
318
320 [[deprecated("Will be removed after release 3.6. Use other signature.")]]
321 ScvfCornerStorage getBoundaryScvfCorners(const Intersection& is,
322 const typename Intersection::Geometry& geometry,
323 unsigned int indexInIntersection) const
324 {
325 return getBoundaryScvfCorners(is.indexInInside(), indexInIntersection);
326 }
327
329 GlobalPosition normal(const ScvfCornerStorage& scvfCorners,
330 const std::vector<unsigned int>& scvIndices) const
331 {
332 auto normal = geo_.corner(1) - geo_.corner(0);
333 normal /= normal.two_norm();
334 return normal;
335 }
336
338 [[deprecated("Will be removed after 3.6. Use volume from geometry/volume.hh")]]
339 Scalar scvVolume(const ScvCornerStorage& scvCorners) const
340 {
341 return (scvCorners[1] - scvCorners[0]).two_norm();
342 }
343
345 [[deprecated("Will be removed after 3.6. Use volume from geometry/volume.hh")]]
346 Scalar scvfArea(const ScvfCornerStorage& scvfCorners) const
347 {
348 return 1.0;
349 }
350
352 std::size_t numInteriorScvf() const
353 {
354 return referenceElement(geo_).size(dim-1);
355 }
356
358 std::size_t numScv() const
359 {
360 return referenceElement(geo_).size(dim);
361 }
362
364 const typename Element::Geometry& elementGeometry() const
365 { return geo_; }
366
367private:
368 const typename Element::Geometry& geo_;
369};
370
372template <class GridView, class ScvType, class ScvfType>
373class BoxGeometryHelper<GridView, 2, ScvType, ScvfType>
374{
375 using Scalar = typename GridView::ctype;
376 using GlobalPosition = typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
377 using ScvCornerStorage = typename ScvType::Traits::CornerStorage;
378 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
379
380 using Element = typename GridView::template Codim<0>::Entity;
381 using Intersection = typename GridView::Intersection;
382
383 static constexpr auto dim = GridView::dimension;
384 static constexpr auto dimWorld = GridView::dimensionworld;
385public:
386
387 BoxGeometryHelper(const typename Element::Geometry& geometry)
388 : geo_(geometry)
389 {}
390
392 ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
393 {
394 // proceed according to number of corners of the element
395 const auto type = geo_.type();
396 if (type == Dune::GeometryTypes::triangle)
397 {
399 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localScvIdx]);
400 }
401 else if (type == Dune::GeometryTypes::quadrilateral)
402 {
404 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localScvIdx]);
405 }
406 else
407 DUNE_THROW(Dune::NotImplemented, "Box scv geometries for dim=" << dim
408 << " dimWorld=" << dimWorld
409 << " type=" << type);
410 }
411
413 ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
414 {
415 // proceed according to number of corners
416 const auto type = geo_.type();
417 if (type == Dune::GeometryTypes::triangle)
418 {
420 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localScvfIdx]);
421 }
422 else if (type == Dune::GeometryTypes::quadrilateral)
423 {
425 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localScvfIdx]);
426 }
427 else
428 DUNE_THROW(Dune::NotImplemented, "Box scvf geometries for dim=" << dim
429 << " dimWorld=" << dimWorld
430 << " type=" << type);
431 }
432
434 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex,
435 unsigned int indexInFacet) const
436 {
437 // we have to use the corresponding facet geometry as the intersection geometry
438 // might be rotated or flipped. This makes sure that the corners (dof location)
439 // and corresponding scvfs are sorted in the same way
441 constexpr int facetCodim = 1;
442 return Detail::Box::subEntityKeyToCornerStorage<ScvfCornerStorage>(geo_, localFacetIndex, facetCodim, Corners::keys[indexInFacet]);
443 }
444
446 [[deprecated("Will be removed after release 3.6. Use other signature.")]]
447 ScvfCornerStorage getBoundaryScvfCorners(const Intersection& is,
448 const typename Intersection::Geometry& isGeometry,
449 unsigned int indexInIntersection) const
450 {
451 return getBoundaryScvfCorners(is.indexInInside(), indexInIntersection);
452 }
453
455 template <int w = dimWorld>
456 typename std::enable_if<w == 3, GlobalPosition>::type
457 normal(const ScvfCornerStorage& scvfCorners,
458 const std::vector<unsigned int>& scvIndices) const
459 {
460 const auto v1 = geo_.corner(1) - geo_.corner(0);
461 const auto v2 = geo_.corner(2) - geo_.corner(0);
462 const auto v3 = Dumux::crossProduct(v1, v2);
463 const auto t = scvfCorners[1] - scvfCorners[0];
464 GlobalPosition normal = Dumux::crossProduct(v3, t);
465 normal /= normal.two_norm();
466
468 const auto v = geo_.corner(scvIndices[1]) - geo_.corner(scvIndices[0]);
469 const auto s = v*normal;
470 if (std::signbit(s))
471 normal *= -1;
472
473 return normal;
474 }
475
477 template <int w = dimWorld>
478 typename std::enable_if<w == 2, GlobalPosition>::type
479 normal(const ScvfCornerStorage& scvfCorners,
480 const std::vector<unsigned int>& scvIndices) const
481 {
483 const auto t = scvfCorners[1] - scvfCorners[0];
484 GlobalPosition normal({-t[1], t[0]});
485 normal /= normal.two_norm();
486
488 const auto v = geo_.corner(scvIndices[1]) - geo_.corner(scvIndices[0]);
489 const auto s = v*normal;
490 if (std::signbit(s))
491 normal *= -1;
492
493 return normal;
494 }
495
497 template <int w = dimWorld>
498 [[deprecated("Will be removed after 3.6. Use volume from geometry/volume.hh")]]
499 typename std::enable_if<w == 3, Scalar>::type
500 scvVolume(const ScvCornerStorage& p) const
501 {
502 return 0.5*Dumux::crossProduct(p[3]-p[0], p[2]-p[1]).two_norm();
503 }
504
506 template <int w = dimWorld>
507 [[deprecated("Will be removed after 3.6. Use volume from geometry/volume.hh")]]
508 typename std::enable_if<w == 2, Scalar>::type
509 scvVolume(const ScvCornerStorage& p) const
510 {
513 using std::abs;
514 return 0.5*abs(Dumux::crossProduct(p[3]-p[0], p[2]-p[1]));
515 }
516
518 [[deprecated("Will be removed after 3.6. Use volume from geometry/volume.hh")]]
519 Scalar scvfArea(const ScvfCornerStorage& p) const
520 {
521 return (p[1]-p[0]).two_norm();
522 }
523
525 std::size_t numInteriorScvf() const
526 {
527 return referenceElement(geo_).size(dim-1);
528 }
529
531 std::size_t numScv() const
532 {
533 return referenceElement(geo_).size(dim);
534 }
535
537 const typename Element::Geometry& elementGeometry() const
538 { return geo_; }
539
540private:
541 const typename Element::Geometry& geo_;
542};
543
545template <class GridView, class ScvType, class ScvfType>
546class BoxGeometryHelper<GridView, 3, ScvType, ScvfType>
547{
548 using Scalar = typename GridView::ctype;
549 using GlobalPosition = typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
550 using ScvCornerStorage = typename ScvType::Traits::CornerStorage;
551 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
552
553 using Element = typename GridView::template Codim<0>::Entity;
554 using Intersection = typename GridView::Intersection;
555
556 static constexpr auto dim = GridView::dimension;
557 static constexpr auto dimWorld = GridView::dimensionworld;
558
559public:
560 BoxGeometryHelper(const typename Element::Geometry& geometry)
561 : geo_(geometry)
562 {}
563
565 ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
566 {
567 const auto type = geo_.type();
568 if (type == Dune::GeometryTypes::tetrahedron)
569 {
571 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localScvIdx]);
572 }
573 else if (type == Dune::GeometryTypes::prism)
574 {
576 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localScvIdx]);
577 }
578 else if (type == Dune::GeometryTypes::hexahedron)
579 {
581 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localScvIdx]);
582 }
583 else
584 DUNE_THROW(Dune::NotImplemented, "Box scv geometries for dim=" << dim
585 << " dimWorld=" << dimWorld
586 << " type=" << type);
587 }
588
590 ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
591 {
592 // proceed according to number of corners
593 const auto type = geo_.type();
594 if (type == Dune::GeometryTypes::tetrahedron)
595 {
597 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localScvfIdx]);
598 }
599 else if (type == Dune::GeometryTypes::prism)
600 {
602 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localScvfIdx]);
603 }
604 else if (type == Dune::GeometryTypes::hexahedron)
605 {
607 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localScvfIdx]);
608 }
609 else
610 DUNE_THROW(Dune::NotImplemented, "Box scvf geometries for dim=" << dim
611 << " dimWorld=" << dimWorld
612 << " type=" << type);
613 }
614
616 ScvfCornerStorage getBoundaryScvfCorners(unsigned localFacetIndex,
617 unsigned int indexInFacet) const
618 {
619 constexpr int facetCodim = 1;
620
621 // we have to use the corresponding facet geometry as the intersection geometry
622 // might be rotated or flipped. This makes sure that the corners (dof location)
623 // and corresponding scvfs are sorted in the same way
624 using Dune::referenceElement;
625 const auto type = referenceElement(geo_).type(localFacetIndex, facetCodim);
626 if (type == Dune::GeometryTypes::triangle)
627 {
629 return Detail::Box::subEntityKeyToCornerStorage<ScvfCornerStorage>(geo_, localFacetIndex, facetCodim, Corners::keys[indexInFacet]);
630 }
631 else if (type == Dune::GeometryTypes::quadrilateral)
632 {
634 return Detail::Box::subEntityKeyToCornerStorage<ScvfCornerStorage>(geo_, localFacetIndex, facetCodim, Corners::keys[indexInFacet]);
635 }
636 else
637 DUNE_THROW(Dune::NotImplemented, "Box boundary scvf geometries for dim=" << dim
638 << " dimWorld=" << dimWorld
639 << " type=" << type);
640 }
641
643 [[deprecated("Will be removed after release 3.6. Use other signature.")]]
644 ScvfCornerStorage getBoundaryScvfCorners(const Intersection& is,
645 const typename Intersection::Geometry& isGeometry,
646 unsigned int indexInIntersection) const
647 {
648 return getBoundaryScvfCorners(is.indexInInside(), indexInIntersection);
649 }
650
652 GlobalPosition normal(const ScvfCornerStorage& p,
653 const std::vector<unsigned int>& scvIndices) const
654 {
655 auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]);
656 normal /= normal.two_norm();
657
658 const auto v = geo_.corner(scvIndices[1]) - geo_.corner(scvIndices[0]);
659 const auto s = v*normal;
660 if (std::signbit(s))
661 normal *= -1;
662
663 return normal;
664 }
665
667 [[deprecated("Will be removed after 3.6. Use volume from geometry/volume.hh")]]
668 Scalar scvVolume(const ScvCornerStorage& p) const
669 {
670 // after Grandy 1997, Efficient computation of volume of hexahedron
671 const auto v = p[7]-p[0];
672 return 1.0/6.0 * ( Dumux::tripleProduct(v, p[1]-p[0], p[3]-p[5])
673 + Dumux::tripleProduct(v, p[4]-p[0], p[5]-p[6])
674 + Dumux::tripleProduct(v, p[2]-p[0], p[6]-p[3]));
675 }
676
678 [[deprecated("Will be removed after 3.6. Use volume from geometry/volume.hh")]]
679 Scalar scvfArea(const ScvfCornerStorage& p) const
680 {
681 // after Wolfram alpha quadrilateral area
682 return 0.5*Dumux::crossProduct(p[3]-p[0], p[2]-p[1]).two_norm();
683 }
684
686 std::size_t numInteriorScvf() const
687 {
688 return referenceElement(geo_).size(dim-1);
689 }
690
692 std::size_t numScv() const
693 {
694 return referenceElement(geo_).size(dim);
695 }
696
698 const typename Element::Geometry& elementGeometry() const
699 { return geo_; }
700
701private:
702 const typename Element::Geometry& geo_;
703};
704
705} // end namespace Dumux
706
707#endif
Define some often used mathematical functions.
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:38
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
Scalar tripleProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2, const Dune::FieldVector< Scalar, 3 > &vec3)
Triple product of three vectors in three-dimensional Euclidean space retuning scalar.
Definition: math.hh:683
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: deprecated.hh:149
S subEntityKeyToCornerStorage(const Geo &geo, unsigned int i, unsigned int c, const std::array< T, N > &key)
Definition: boxgeometryhelper.hh:260
S subEntityKeyToCornerStorageImpl(const Geo &geo, unsigned int i, unsigned int c, const KeyArray &key, std::index_sequence< I... >)
Definition: boxgeometryhelper.hh:248
S keyToCornerStorage(const Geo &geo, const std::array< T, N > &key)
Definition: boxgeometryhelper.hh:240
S keyToCornerStorageImpl(const Geo &geo, const KeyArray &key, std::index_sequence< I... >)
Definition: boxgeometryhelper.hh:230
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
constexpr Line line
Definition: couplingmanager1d3d_line.hh:43
Traits for an efficient corner storage for box method sub control volumes.
Definition: boxgeometryhelper.hh:44
Definition: boxgeometryhelper.hh:49
std::array< Dune::FieldVector< ct, cdim >,(1<<(mydim)) > Type
Definition: boxgeometryhelper.hh:50
Definition: boxgeometryhelper.hh:56
static const bool v
Definition: boxgeometryhelper.hh:57
static const unsigned int topologyId
Definition: boxgeometryhelper.hh:58
Definition: boxgeometryhelper.hh:65
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:70
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:80
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:91
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:103
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:115
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:129
Definition: boxgeometryhelper.hh:143
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:148
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:157
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:168
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:180
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:194
std::pair< std::uint8_t, std::uint8_t > Key
Definition: boxgeometryhelper.hh:211
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:269
Scalar scvfArea(const ScvfCornerStorage &scvfCorners) const
get scvf area
Definition: boxgeometryhelper.hh:346
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition: boxgeometryhelper.hh:364
std::size_t numInteriorScvf() const
number of sub control volume faces (number of edges)
Definition: boxgeometryhelper.hh:352
BoxGeometryHelper(const typename Element::Geometry &geometry)
Definition: boxgeometryhelper.hh:289
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition: boxgeometryhelper.hh:294
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int) const
Create the sub control volume face geometries on the boundary.
Definition: boxgeometryhelper.hh:313
GlobalPosition normal(const ScvfCornerStorage &scvfCorners, const std::vector< unsigned int > &scvIndices) const
get scvf normal vector
Definition: boxgeometryhelper.hh:329
std::size_t numScv() const
number of sub control volumes (number of vertices)
Definition: boxgeometryhelper.hh:358
ScvGeometry scvGeometry(unsigned int localScvIdx) const
Definition: boxgeometryhelper.hh:300
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition: boxgeometryhelper.hh:306
ScvfCornerStorage getBoundaryScvfCorners(const Intersection &is, const typename Intersection::Geometry &geometry, unsigned int indexInIntersection) const
Create the sub control volume face geometries on the boundary.
Definition: boxgeometryhelper.hh:321
Scalar scvVolume(const ScvCornerStorage &scvCorners) const
get scv volume
Definition: boxgeometryhelper.hh:339
std::size_t numScv() const
number of sub control volumes (number of vertices)
Definition: boxgeometryhelper.hh:531
std::enable_if< w==2, GlobalPosition >::type normal(const ScvfCornerStorage &scvfCorners, const std::vector< unsigned int > &scvIndices) const
get scvf normal vector for dim == 2, dimworld == 2
Definition: boxgeometryhelper.hh:479
std::size_t numInteriorScvf() const
number of sub control volume faces (number of edges)
Definition: boxgeometryhelper.hh:525
Scalar scvfArea(const ScvfCornerStorage &p) const
get scvf area
Definition: boxgeometryhelper.hh:519
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition: boxgeometryhelper.hh:413
BoxGeometryHelper(const typename Element::Geometry &geometry)
Definition: boxgeometryhelper.hh:387
std::enable_if< w==3, Scalar >::type scvVolume(const ScvCornerStorage &p) const
get scv volume for dim == 2, dimworld == 3
Definition: boxgeometryhelper.hh:500
ScvfCornerStorage getBoundaryScvfCorners(const Intersection &is, const typename Intersection::Geometry &isGeometry, unsigned int indexInIntersection) const
Create the sub control volume face geometries on the boundary.
Definition: boxgeometryhelper.hh:447
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition: boxgeometryhelper.hh:434
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition: boxgeometryhelper.hh:537
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition: boxgeometryhelper.hh:392
std::enable_if< w==2, Scalar >::type scvVolume(const ScvCornerStorage &p) const
get scv volume for dim == 2, dimworld == 2
Definition: boxgeometryhelper.hh:509
std::enable_if< w==3, GlobalPosition >::type normal(const ScvfCornerStorage &scvfCorners, const std::vector< unsigned int > &scvIndices) const
get scvf normal vector for dim == 2, dimworld == 3
Definition: boxgeometryhelper.hh:457
Scalar scvfArea(const ScvfCornerStorage &p) const
get scvf area
Definition: boxgeometryhelper.hh:679
ScvfCornerStorage getBoundaryScvfCorners(unsigned localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition: boxgeometryhelper.hh:616
ScvfCornerStorage getBoundaryScvfCorners(const Intersection &is, const typename Intersection::Geometry &isGeometry, unsigned int indexInIntersection) const
Create the sub control volume face geometries on the boundary.
Definition: boxgeometryhelper.hh:644
GlobalPosition normal(const ScvfCornerStorage &p, const std::vector< unsigned int > &scvIndices) const
get scvf normal vector
Definition: boxgeometryhelper.hh:652
BoxGeometryHelper(const typename Element::Geometry &geometry)
Definition: boxgeometryhelper.hh:560
Scalar scvVolume(const ScvCornerStorage &p) const
get scv volume
Definition: boxgeometryhelper.hh:668
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition: boxgeometryhelper.hh:698
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition: boxgeometryhelper.hh:565
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the scvf corners.
Definition: boxgeometryhelper.hh:590
std::size_t numScv() const
number of sub control volumes (number of vertices)
Definition: boxgeometryhelper.hh:692
std::size_t numInteriorScvf() const
number of sub control volume faces (number of edges)
Definition: boxgeometryhelper.hh:686