13#ifndef DUMUX_DISCRETIZATION_BOX_GEOMETRY_HELPER_HH
14#define DUMUX_DISCRETIZATION_BOX_GEOMETRY_HELPER_HH
18#include <dune/common/exceptions.hh>
20#include <dune/geometry/type.hh>
21#include <dune/geometry/typeindex.hh>
22#include <dune/geometry/referenceelements.hh>
23#include <dune/geometry/multilineargeometry.hh>
35 template<
int mydim,
int cdim >
38 using Type = std::array< Dune::FieldVector< ct, cdim >, (1<<(mydim)) >;
45 static const bool v =
true;
46 static const unsigned int topologyId = Dune::GeometryTypes::cube(mydim).id();
52template<Dune::GeometryType::Id gt>
58 using Key = std::pair<std::uint8_t, std::uint8_t>;
59 static constexpr std::array<std::array<Key, 2>, 2>
keys = {{
68 using Key = std::pair<std::uint8_t, std::uint8_t>;
69 static constexpr std::array<std::array<Key, 4>, 3>
keys = {{
79 using Key = std::pair<std::uint8_t, std::uint8_t>;
80 static constexpr std::array<std::array<Key, 4>, 4>
keys = {{
91 using Key = std::pair<std::uint8_t, std::uint8_t>;
92 static constexpr std::array<std::array<Key, 8>, 4>
keys = {{
93 {
Key{0, 3},
Key{0, 2},
Key{1, 2},
Key{0, 1},
Key{3, 2},
Key{1, 1},
Key{2, 1},
Key{0, 0} },
94 {
Key{1, 3},
Key{2, 2},
Key{0, 2},
Key{0, 1},
Key{4, 2},
Key{3, 1},
Key{1, 1},
Key{0, 0} },
95 {
Key{2, 3},
Key{1, 2},
Key{2, 2},
Key{0, 1},
Key{5, 2},
Key{2, 1},
Key{3, 1},
Key{0, 0} },
96 {
Key{3, 3},
Key{3, 2},
Key{5, 2},
Key{2, 1},
Key{4, 2},
Key{1, 1},
Key{3, 1},
Key{0, 0} }
103 using Key = std::pair<std::uint8_t, std::uint8_t>;
104 static constexpr std::array<std::array<Key, 8>, 6>
keys = {{
105 {
Key{0, 3},
Key{3, 2},
Key{4, 2},
Key{3, 1},
Key{0, 2},
Key{0, 1},
Key{1, 1},
Key{0, 0} },
106 {
Key{1, 3},
Key{5, 2},
Key{3, 2},
Key{3, 1},
Key{1, 2},
Key{2, 1},
Key{0, 1},
Key{0, 0} },
107 {
Key{2, 3},
Key{4, 2},
Key{5, 2},
Key{3, 1},
Key{2, 2},
Key{1, 1},
Key{2, 1},
Key{0, 0} },
108 {
Key{3, 3},
Key{7, 2},
Key{6, 2},
Key{4, 1},
Key{0, 2},
Key{1, 1},
Key{0, 1},
Key{0, 0} },
109 {
Key{4, 3},
Key{6, 2},
Key{8, 2},
Key{4, 1},
Key{1, 2},
Key{0, 1},
Key{2, 1},
Key{0, 0} },
110 {
Key{5, 3},
Key{8, 2},
Key{7, 2},
Key{4, 1},
Key{2, 2},
Key{2, 1},
Key{1, 1},
Key{0, 0} }
117 using Key = std::pair<std::uint8_t, std::uint8_t>;
118 static constexpr std::array<std::array<Key, 8>, 8>
keys = {{
119 {
Key{0, 3},
Key{6, 2},
Key{4, 2},
Key{4, 1},
Key{0, 2},
Key{2, 1},
Key{0, 1},
Key{0, 0} },
120 {
Key{1, 3},
Key{5, 2},
Key{6, 2},
Key{4, 1},
Key{1, 2},
Key{1, 1},
Key{2, 1},
Key{0, 0} },
121 {
Key{2, 3},
Key{4, 2},
Key{7, 2},
Key{4, 1},
Key{2, 2},
Key{0, 1},
Key{3, 1},
Key{0, 0} },
122 {
Key{3, 3},
Key{7, 2},
Key{5, 2},
Key{4, 1},
Key{3, 2},
Key{3, 1},
Key{1, 1},
Key{0, 0} },
123 {
Key{4, 3},
Key{8, 2},
Key{10, 2},
Key{5, 1},
Key{0, 2},
Key{0, 1},
Key{2, 1},
Key{0, 0} },
124 {
Key{5, 3},
Key{10, 2},
Key{9, 2},
Key{5, 1},
Key{1, 2},
Key{2, 1},
Key{1, 1},
Key{0, 0} },
125 {
Key{6, 3},
Key{11, 2},
Key{8, 2},
Key{5, 1},
Key{2, 2},
Key{3, 1},
Key{0, 1},
Key{0, 0} },
126 {
Key{7, 3},
Key{9, 2},
Key{11, 2},
Key{5, 1},
Key{3, 2},
Key{1, 1},
Key{3, 1},
Key{0, 0} }
130template<Dune::GeometryType::Id gt>
136 using Key = std::pair<std::uint8_t, std::uint8_t>;
137 static constexpr std::array<std::array<Key, 1>, 1>
keys = {{
145 using Key = std::pair<std::uint8_t, std::uint8_t>;
146 static constexpr std::array<std::array<Key, 2>, 3>
keys = {{
156 using Key = std::pair<std::uint8_t, std::uint8_t>;
157 static constexpr std::array<std::array<Key, 2>, 4>
keys = {{
168 using Key = std::pair<std::uint8_t, std::uint8_t>;
169 static constexpr std::array<std::array<Key, 4>, 6>
keys = {{
182 using Key = std::pair<std::uint8_t, std::uint8_t>;
183 static constexpr std::array<std::array<Key, 4>, 9>
keys = {{
199 using Key = std::pair<std::uint8_t, std::uint8_t>;
200 static constexpr std::array<std::array<Key, 4>, 12>
keys = {{
217template<
class S,
class Geo,
class KeyArray, std::size_t... I>
220 using Dune::referenceElement;
221 const auto ref = referenceElement(geo);
223 return { geo.global(ref.position(key[I].first, key[I].second))... };
227template<
class S,
class Geo,
class T, std::
size_t N,
class Indices = std::make_index_sequence<N>>
235template<
class S,
class Geo,
class KeyArray, std::size_t... I>
238 using Dune::referenceElement;
239 const auto ref = referenceElement(geo);
242 return { geo.global(ref.position(ref.subEntity(i, c, key[I].first, c+key[I].second), c+key[I].second))... };
247template<
class S,
class Geo,
class T, std::
size_t N,
class Indices = std::make_index_sequence<N>>
256template<
class Gr
idView,
int dim,
class ScvType,
class ScvfType>
260template <
class Gr
idView,
class ScvType,
class ScvfType>
264 using Scalar =
typename GridView::ctype;
265 using GlobalPosition =
typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
266 using ScvCornerStorage =
typename ScvType::Traits::CornerStorage;
267 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
268 using ScvGeometry =
typename ScvType::Traits::Geometry;
269 using ScvfGeometry =
typename ScvfType::Traits::Geometry;
271 using Element =
typename GridView::template Codim<0>::Entity;
272 using Intersection =
typename GridView::Intersection;
274 static constexpr int dim = 1;
290 return { Dune::GeometryTypes::line,
getScvCorners(localScvIdx) };
304 return ScvfCornerStorage{{ geo_.corner(localFacetIndex) }};
308 GlobalPosition
normal(
const ScvfCornerStorage& scvfCorners,
309 const std::vector<unsigned int>& scvIndices)
const
311 auto normal = geo_.corner(1) - geo_.corner(0);
319 return referenceElement(geo_).size(dim-1);
325 return referenceElement(geo_).size(dim);
333 const typename Element::Geometry& geo_;
337template <
class Gr
idView,
class ScvType,
class ScvfType>
338class BoxGeometryHelper<GridView, 2, ScvType, ScvfType>
340 using Scalar =
typename GridView::ctype;
341 using GlobalPosition =
typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
342 using ScvCornerStorage =
typename ScvType::Traits::CornerStorage;
343 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
345 using Element =
typename GridView::template Codim<0>::Entity;
346 using Intersection =
typename GridView::Intersection;
348 static constexpr auto dim = GridView::dimension;
349 static constexpr auto dimWorld = GridView::dimensionworld;
352 explicit BoxGeometryHelper(
const typename Element::Geometry& geometry)
357 ScvCornerStorage getScvCorners(
unsigned int localScvIdx)
const
360 const auto type = geo_.type();
361 if (type == Dune::GeometryTypes::triangle)
363 using Corners = Detail::Box::ScvCorners<Dune::GeometryTypes::triangle>;
366 else if (type == Dune::GeometryTypes::quadrilateral)
368 using Corners = Detail::Box::ScvCorners<Dune::GeometryTypes::quadrilateral>;
372 DUNE_THROW(Dune::NotImplemented,
"Box scv geometries for dim=" << dim
373 <<
" dimWorld=" << dimWorld
374 <<
" type=" << type);
378 ScvfCornerStorage getScvfCorners(
unsigned int localScvfIdx)
const
381 const auto type = geo_.type();
382 if (type == Dune::GeometryTypes::triangle)
384 using Corners = Detail::Box::ScvfCorners<Dune::GeometryTypes::triangle>;
387 else if (type == Dune::GeometryTypes::quadrilateral)
389 using Corners = Detail::Box::ScvfCorners<Dune::GeometryTypes::quadrilateral>;
393 DUNE_THROW(Dune::NotImplemented,
"Box scvf geometries for dim=" << dim
394 <<
" dimWorld=" << dimWorld
395 <<
" type=" << type);
399 ScvfCornerStorage getBoundaryScvfCorners(
unsigned int localFacetIndex,
400 unsigned int indexInFacet)
const
405 using Corners = Detail::Box::ScvCorners<Dune::GeometryTypes::line>;
406 constexpr int facetCodim = 1;
411 template <
int w = dimWorld>
412 typename std::enable_if<w == 3, GlobalPosition>::type
413 normal(
const ScvfCornerStorage& scvfCorners,
414 const std::vector<unsigned int>& scvIndices)
const
416 const auto v1 = geo_.corner(1) - geo_.corner(0);
417 const auto v2 = geo_.corner(2) - geo_.corner(0);
419 const auto t = scvfCorners[1] - scvfCorners[0];
424 const auto v = geo_.corner(scvIndices[1]) - geo_.corner(scvIndices[0]);
433 template <
int w = dimWorld>
434 typename std::enable_if<w == 2, GlobalPosition>::type
435 normal(
const ScvfCornerStorage& scvfCorners,
436 const std::vector<unsigned int>& scvIndices)
const
439 const auto t = scvfCorners[1] - scvfCorners[0];
440 GlobalPosition
normal({-t[1], t[0]});
444 const auto v = geo_.corner(scvIndices[1]) - geo_.corner(scvIndices[0]);
453 std::size_t numInteriorScvf()
const
455 return referenceElement(geo_).size(dim-1);
459 std::size_t numScv()
const
461 return referenceElement(geo_).size(dim);
465 const typename Element::Geometry& elementGeometry()
const
469 const typename Element::Geometry& geo_;
473template <
class Gr
idView,
class ScvType,
class ScvfType>
476 using Scalar =
typename GridView::ctype;
477 using GlobalPosition =
typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
478 using ScvCornerStorage =
typename ScvType::Traits::CornerStorage;
479 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
481 using Element =
typename GridView::template Codim<0>::Entity;
482 using Intersection =
typename GridView::Intersection;
484 static constexpr auto dim = GridView::dimension;
485 static constexpr auto dimWorld = GridView::dimensionworld;
488 explicit BoxGeometryHelper(
const typename Element::Geometry& geometry)
493 ScvCornerStorage getScvCorners(
unsigned int localScvIdx)
const
495 const auto type = geo_.type();
496 if (type == Dune::GeometryTypes::tetrahedron)
498 using Corners = Detail::Box::ScvCorners<Dune::GeometryTypes::tetrahedron>;
501 else if (type == Dune::GeometryTypes::prism)
503 using Corners = Detail::Box::ScvCorners<Dune::GeometryTypes::prism>;
506 else if (type == Dune::GeometryTypes::hexahedron)
508 using Corners = Detail::Box::ScvCorners<Dune::GeometryTypes::hexahedron>;
512 DUNE_THROW(Dune::NotImplemented,
"Box scv geometries for dim=" << dim
513 <<
" dimWorld=" << dimWorld
514 <<
" type=" << type);
518 ScvfCornerStorage getScvfCorners(
unsigned int localScvfIdx)
const
521 const auto type = geo_.type();
522 if (type == Dune::GeometryTypes::tetrahedron)
524 using Corners = Detail::Box::ScvfCorners<Dune::GeometryTypes::tetrahedron>;
527 else if (type == Dune::GeometryTypes::prism)
529 using Corners = Detail::Box::ScvfCorners<Dune::GeometryTypes::prism>;
532 else if (type == Dune::GeometryTypes::hexahedron)
534 using Corners = Detail::Box::ScvfCorners<Dune::GeometryTypes::hexahedron>;
538 DUNE_THROW(Dune::NotImplemented,
"Box scvf geometries for dim=" << dim
539 <<
" dimWorld=" << dimWorld
540 <<
" type=" << type);
544 ScvfCornerStorage getBoundaryScvfCorners(
unsigned localFacetIndex,
545 unsigned int indexInFacet)
const
547 constexpr int facetCodim = 1;
552 using Dune::referenceElement;
553 const auto type = referenceElement(geo_).type(localFacetIndex, facetCodim);
554 if (type == Dune::GeometryTypes::triangle)
556 using Corners = Detail::Box::ScvCorners<Dune::GeometryTypes::triangle>;
559 else if (type == Dune::GeometryTypes::quadrilateral)
561 using Corners = Detail::Box::ScvCorners<Dune::GeometryTypes::quadrilateral>;
565 DUNE_THROW(Dune::NotImplemented,
"Box boundary scvf geometries for dim=" << dim
566 <<
" dimWorld=" << dimWorld
567 <<
" type=" << type);
571 GlobalPosition
normal(
const ScvfCornerStorage& p,
572 const std::vector<unsigned int>& scvIndices)
const
577 const auto v = geo_.corner(scvIndices[1]) - geo_.corner(scvIndices[0]);
586 std::size_t numInteriorScvf()
const
588 return referenceElement(geo_).size(dim-1);
592 std::size_t numScv()
const
594 return referenceElement(geo_).size(dim);
598 const typename Element::Geometry& elementGeometry()
const
602 const typename Element::Geometry& geo_;
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition boxgeometryhelper.hh:329
std::size_t numInteriorScvf() const
number of sub control volume faces (number of edges)
Definition boxgeometryhelper.hh:317
BoxGeometryHelper(const typename Element::Geometry &geometry)
Definition boxgeometryhelper.hh:277
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition boxgeometryhelper.hh:282
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int) const
Create the sub control volume face geometries on the boundary.
Definition boxgeometryhelper.hh:301
GlobalPosition normal(const ScvfCornerStorage &scvfCorners, const std::vector< unsigned int > &scvIndices) const
get scvf normal vector
Definition boxgeometryhelper.hh:308
std::size_t numScv() const
number of sub control volumes (number of vertices)
Definition boxgeometryhelper.hh:323
ScvGeometry scvGeometry(unsigned int localScvIdx) const
Definition boxgeometryhelper.hh:288
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition boxgeometryhelper.hh:294
Create sub control volumes and sub control volume face geometries.
Definition boxgeometryhelper.hh:257
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:671
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition normal.hh:26
Define some often used mathematical functions.
Definition boxgeometryhelper.hh:50
S subEntityKeyToCornerStorage(const Geo &geo, unsigned int i, unsigned int c, const std::array< T, N > &key)
Definition boxgeometryhelper.hh:248
S subEntityKeyToCornerStorageImpl(const Geo &geo, unsigned int i, unsigned int c, const KeyArray &key, std::index_sequence< I... >)
Definition boxgeometryhelper.hh:236
S keyToCornerStorage(const Geo &geo, const std::array< T, N > &key)
Definition boxgeometryhelper.hh:228
S keyToCornerStorageImpl(const Geo &geo, const KeyArray &key, std::index_sequence< I... >)
Definition boxgeometryhelper.hh:218
Definition common/pdesolver.hh:24
Definition boxgeometryhelper.hh:37
std::array< Dune::FieldVector< ct, cdim >,(1<<(mydim)) > Type
Definition boxgeometryhelper.hh:38
Definition boxgeometryhelper.hh:44
static const bool v
Definition boxgeometryhelper.hh:45
static const unsigned int topologyId
Definition boxgeometryhelper.hh:46
Traits for an efficient corner storage for box method sub control volumes.
Definition boxgeometryhelper.hh:32
static constexpr std::array< std::array< Key, 8 >, 8 > keys
Definition boxgeometryhelper.hh:118
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:117
static constexpr std::array< std::array< Key, 2 >, 2 > keys
Definition boxgeometryhelper.hh:59
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:58
static constexpr std::array< std::array< Key, 8 >, 6 > keys
Definition boxgeometryhelper.hh:104
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:103
static constexpr std::array< std::array< Key, 4 >, 4 > keys
Definition boxgeometryhelper.hh:80
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:79
static constexpr std::array< std::array< Key, 8 >, 4 > keys
Definition boxgeometryhelper.hh:92
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:68
static constexpr std::array< std::array< Key, 4 >, 3 > keys
Definition boxgeometryhelper.hh:69
Definition boxgeometryhelper.hh:53
static constexpr std::array< std::array< Key, 4 >, 12 > keys
Definition boxgeometryhelper.hh:200
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:199
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:136
static constexpr std::array< std::array< Key, 1 >, 1 > keys
Definition boxgeometryhelper.hh:137
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:182
static constexpr std::array< std::array< Key, 4 >, 9 > keys
Definition boxgeometryhelper.hh:183
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:156
static constexpr std::array< std::array< Key, 2 >, 4 > keys
Definition boxgeometryhelper.hh:157
static constexpr std::array< std::array< Key, 4 >, 6 > keys
Definition boxgeometryhelper.hh:169
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:168
static constexpr std::array< std::array< Key, 2 >, 3 > keys
Definition boxgeometryhelper.hh:146
std::pair< std::uint8_t, std::uint8_t > Key
Definition boxgeometryhelper.hh:145
Definition boxgeometryhelper.hh:131