12#ifndef DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_GEOMETRY_HELPER_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_GEOMETRY_HELPER_HH
17#include <dune/common/reservedvector.hh>
18#include <dune/common/fvector.hh>
19#include <dune/geometry/multilineargeometry.hh>
20#include <dune/geometry/type.hh>
33 template<
int mydim,
int cdim >
36 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)>;
42template<Dune::GeometryType::Id gt>
48 static constexpr Dune::GeometryType
type()
49 {
return Dune::GeometryTypes::triangle; }
51 using Key = std::pair<std::uint8_t, std::uint8_t>;
52 static constexpr std::array<std::array<Key, 3>, 3> keys = {{
62 static constexpr Dune::GeometryType
type()
63 {
return Dune::GeometryTypes::triangle; }
65 using Key = std::pair<std::uint8_t, std::uint8_t>;
66 static constexpr std::array<std::array<Key, 3>, 4> keys = {{
77 static constexpr Dune::GeometryType
type()
78 {
return Dune::GeometryTypes::tetrahedron; }
80 using Key = std::pair<std::uint8_t, std::uint8_t>;
81 static constexpr std::array<std::array<Key, 4>, 4> keys = {{
92 static constexpr Dune::GeometryType
type()
93 {
return Dune::GeometryTypes::pyramid; }
95 using Key = std::pair<std::uint8_t, std::uint8_t>;
96 static constexpr std::array<std::array<Key, 5>, 6> keys = {{
97 {
Key{4, 3},
Key{0, 3},
Key{6, 3},
Key{2, 3},
Key{0, 0} },
98 {
Key{1, 3},
Key{5, 3},
Key{3, 3},
Key{7, 3},
Key{0, 0} },
99 {
Key{4, 3},
Key{5, 3},
Key{0, 3},
Key{1, 3},
Key{0, 0} },
100 {
Key{2, 3},
Key{3, 3},
Key{6, 3},
Key{7, 3},
Key{0, 0} },
101 {
Key{0, 3},
Key{1, 3},
Key{2, 3},
Key{3, 3},
Key{0, 0} },
106template<Dune::GeometryType::Id gt>
112 static constexpr Dune::GeometryType
type()
115 using Key = std::pair<std::uint8_t, std::uint8_t>;
116 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
126 static constexpr Dune::GeometryType
type()
129 using Key = std::pair<std::uint8_t, std::uint8_t>;
130 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
141 static constexpr Dune::GeometryType
type()
142 {
return Dune::GeometryTypes::triangle; }
144 using Key = std::pair<std::uint8_t, std::uint8_t>;
145 static constexpr std::array<std::array<Key, 3>, 6> keys = {{
158 static constexpr Dune::GeometryType
type()
159 {
return Dune::GeometryTypes::triangle; }
161 using Key = std::pair<std::uint8_t, std::uint8_t>;
162 static constexpr std::array<std::array<Key, 3>, 12> keys = {{
179template<
class S,
class ReferenceElement,
class Transformation,
class KeyArray, std::size_t... I>
180S
keyToCornerStorageImpl(
const ReferenceElement& ref, Transformation&& trans,
const KeyArray& key, std::index_sequence<I...>)
183 return { trans(ref.position(key[I].first, key[I].second))... };
187template<
class S,
class ReferenceElement,
class Transformation,
class T, std::
size_t N,
class Indices = std::make_index_sequence<N>>
188S
keyToCornerStorage(
const ReferenceElement& ref, Transformation&& trans,
const std::array<T, N>& key)
190 return keyToCornerStorageImpl<S>(ref, trans, key, Indices{});
194template<
class S,
int dim,
class ReferenceElement,
class Transformation, std::size_t... ii>
198 return { trans(ref.position(ref.subEntity(i, 1, ii, dim), dim))... };
202template<
class S, std::
size_t numCorners,
int dim,
class ReferenceElement,
class Transformation>
205 return boundaryCornerStorageImpl<S, dim>(ref, trans, i, std::make_index_sequence<numCorners>{});
208template<
class IndexType, Dune::GeometryType::Id gt>
211template<
class IndexType>
214 static constexpr std::array<std::array<IndexType, 2>, 3> pairs = {{
215 {0, 1}, {0, 2}, {1, 2}
219template<
class IndexType>
222 static constexpr std::array<std::array<IndexType, 2>, 4> pairs = {{
223 {0, 2}, {1, 2}, {0, 3}, {1, 3}
227template<
class IndexType>
230 static constexpr std::array<std::array<IndexType, 2>, 6> pairs = {{
231 {0, 1}, {0, 2}, {0, 3}, {1, 2}, {1, 3}, {2, 3},
235template<
class IndexType>
238 static constexpr std::array<std::array<IndexType, 2>, 12> pairs = {{
239 {0, 2}, {1, 2}, {0, 3}, {1, 3}, {0, 4}, {1, 4},
240 {2, 4}, {3, 4}, {0, 5}, {1, 5}, {2, 5}, {3, 5}
250template <
class Gr
idView,
class ScvType,
class ScvfType>
253 using Element =
typename GridView::template Codim<0>::Entity;
255 static constexpr auto dim = GridView::dimension;
256 static constexpr auto dimWorld = GridView::dimensionworld;
258 using ScvCornerStorage =
typename ScvType::Traits::CornerStorage;
259 using LocalIndexType =
typename ScvType::Traits::LocalIndexType;
260 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
263 using Scalar =
typename GridView::ctype;
264 using GlobalPosition =
typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
274 return getScvCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localFacetIndex);
278 template<
class Transformation>
279 static ScvCornerStorage
getScvCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localFacetIndex)
281 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
282 if (type == Dune::GeometryTypes::triangle)
285 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localFacetIndex]);
287 else if (type == Dune::GeometryTypes::quadrilateral)
290 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localFacetIndex]);
292 else if (type == Dune::GeometryTypes::tetrahedron)
295 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localFacetIndex]);
297 else if (type == Dune::GeometryTypes::hexahedron)
300 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localFacetIndex]);
303 DUNE_THROW(Dune::NotImplemented,
"Scv geometries for type " << type);
309 return getScvfCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localEdgeIndex);
313 template<
class Transformation>
314 static ScvfCornerStorage
getScvfCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localEdgeIndex)
316 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
317 if (type == Dune::GeometryTypes::triangle)
320 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localEdgeIndex]);
322 else if (type == Dune::GeometryTypes::quadrilateral)
325 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localEdgeIndex]);
327 else if (type == Dune::GeometryTypes::tetrahedron)
330 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localEdgeIndex]);
332 else if (type == Dune::GeometryTypes::hexahedron)
335 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localEdgeIndex]);
338 DUNE_THROW(Dune::NotImplemented,
"Scvf geometries for type " << type);
344 return getBoundaryScvfCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localFacetIndex);
348 template<
class Transformation>
349 static ScvfCornerStorage
getBoundaryScvfCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localFacetIndex)
351 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
352 if (type == Dune::GeometryTypes::triangle || type == Dune::GeometryTypes::quadrilateral)
353 return Detail::FCDiamond::boundaryCornerStorage<ScvfCornerStorage, 2, dim>(ref, trans, localFacetIndex);
354 else if (type == Dune::GeometryTypes::tetrahedron)
355 return Detail::FCDiamond::boundaryCornerStorage<ScvfCornerStorage, 3, dim>(ref, trans, localFacetIndex);
356 else if (type == Dune::GeometryTypes::hexahedron)
357 return Detail::FCDiamond::boundaryCornerStorage<ScvfCornerStorage, 4, dim>(ref, trans, localFacetIndex);
359 DUNE_THROW(Dune::NotImplemented,
"Boundary scvf geometries for type " << type);
364 return geo_.global(referenceElement(geo_).position(localFacetIndex, 1));
369 const auto type = geo_.type();
370 if (type == Dune::GeometryTypes::triangle)
372 else if (type == Dune::GeometryTypes::quadrilateral)
374 else if (type == Dune::GeometryTypes::tetrahedron)
376 else if (type == Dune::GeometryTypes::hexahedron)
379 DUNE_THROW(Dune::NotImplemented,
"Inside outside scv pairs for type " << type);
385 return referenceElement(geo_).size(2);
391 return referenceElement(geo_).size(1);
394 template<
int d = dimWorld, std::enable_if_t<(d==3),
int> = 0>
395 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
400 const auto ref = referenceElement(geo_);
401 const auto v = facetCenter_(scvPair[1], ref) - facetCenter_(scvPair[0], ref);
410 template<
int d = dimWorld, std::enable_if_t<(d==2),
int> = 0>
411 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
414 const auto t = p[1] - p[0];
415 GlobalPosition
normal({-t[1], t[0]});
418 const auto ref = referenceElement(geo_);
419 const auto v = facetCenter_(scvPair[1], ref) - facetCenter_(scvPair[0], ref);
432 template<
class LocalKey>
433 static Element::Geometry::LocalCoordinate
localDofPosition(Dune::GeometryType type,
const LocalKey& localKey)
435 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
439 static Element::Geometry::LocalCoordinate
localScvfCenter(Dune::GeometryType type,
unsigned int localScvfIdx)
451 template<
class RefElement>
452 GlobalPosition facetCenter_(
unsigned int localFacetIndex,
const RefElement& ref)
const
454 return geo_.global(ref.position(localFacetIndex, 1));
457 const typename Element::Geometry& geo_;
Compute the center point of a convex polytope geometry or a random-access container of corner points.
Helper class to construct SCVs and SCVFs for the diamond scheme.
Definition: discretization/facecentered/diamond/geometryhelper.hh:252
GlobalPosition facetCenter(unsigned int localFacetIndex) const
Definition: discretization/facecentered/diamond/geometryhelper.hh:362
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:307
static ScvfCornerStorage getBoundaryScvfCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localFacetIndex)
Create the sub control volume face geometries on the boundary for a given face index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:349
std::size_t numInteriorScvf()
number of interior sub control volume faces (number of codim-2 entities)
Definition: discretization/facecentered/diamond/geometryhelper.hh:383
static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey &localKey)
local dof position
Definition: discretization/facecentered/diamond/geometryhelper.hh:433
std::array< LocalIndexType, 2 > getInsideOutsideScvForScvf(unsigned int localEdgeIndex)
Definition: discretization/facecentered/diamond/geometryhelper.hh:367
static Element::Geometry::LocalCoordinate localScvfCenter(Dune::GeometryType type, unsigned int localScvfIdx)
local scvf center
Definition: discretization/facecentered/diamond/geometryhelper.hh:439
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:272
DiamondGeometryHelper(const typename Element::Geometry &geo)
Definition: discretization/facecentered/diamond/geometryhelper.hh:267
static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localFacetIndex)
Create a corner storage with the scv corners for a given face (codim-1) index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:279
static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localEdgeIndex)
Create a corner storage with the scvf corners for a given edge (codim-2) index.
Definition: discretization/facecentered/diamond/geometryhelper.hh:314
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition: discretization/facecentered/diamond/geometryhelper.hh:395
static Element::Geometry::LocalCoordinate localBoundaryScvfCenter(Dune::GeometryType type, unsigned int localFacetIndex)
local boundary scvf center
Definition: discretization/facecentered/diamond/geometryhelper.hh:445
const Element::Geometry & elementGeometry() const
Definition: discretization/facecentered/diamond/geometryhelper.hh:428
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:342
std::size_t numScv()
number of sub control volumes (number of codim-1 entities)
Definition: discretization/facecentered/diamond/geometryhelper.hh:389
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
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
Define some often used mathematical functions.
S boundaryCornerStorageImpl(const ReferenceElement &ref, Transformation &&trans, unsigned int i, std::index_sequence< ii... >)
Definition: discretization/facecentered/diamond/geometryhelper.hh:195
S keyToCornerStorage(const ReferenceElement &ref, Transformation &&trans, const std::array< T, N > &key)
Definition: discretization/facecentered/diamond/geometryhelper.hh:188
S boundaryCornerStorage(const ReferenceElement &ref, Transformation &&trans, unsigned int i)
Definition: discretization/facecentered/diamond/geometryhelper.hh:203
S keyToCornerStorageImpl(const ReferenceElement &ref, Transformation &&trans, const KeyArray &key, std::index_sequence< I... >)
Definition: discretization/facecentered/diamond/geometryhelper.hh:180
CVFE< CVFEMethods::CR_RT > FCDiamond
Definition: method.hh:101
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: common/pdesolver.hh:24
Definition: discretization/facecentered/diamond/geometryhelper.hh:209
Definition: discretization/facecentered/diamond/geometryhelper.hh:91
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:92
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:95
Definition: discretization/facecentered/diamond/geometryhelper.hh:61
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:62
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:65
Definition: discretization/facecentered/diamond/geometryhelper.hh:76
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:80
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:77
Definition: discretization/facecentered/diamond/geometryhelper.hh:47
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:48
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:51
Definition: discretization/facecentered/diamond/geometryhelper.hh:43
Definition: discretization/facecentered/diamond/geometryhelper.hh:157
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:158
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:161
Definition: discretization/facecentered/diamond/geometryhelper.hh:125
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:126
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:129
Definition: discretization/facecentered/diamond/geometryhelper.hh:140
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:144
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:141
Definition: discretization/facecentered/diamond/geometryhelper.hh:111
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/facecentered/diamond/geometryhelper.hh:115
static constexpr Dune::GeometryType type()
Definition: discretization/facecentered/diamond/geometryhelper.hh:112
Definition: discretization/facecentered/diamond/geometryhelper.hh:107
Definition: discretization/facecentered/diamond/geometryhelper.hh:35
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<< mydim)> Type
Definition: discretization/facecentered/diamond/geometryhelper.hh:36
Traits for an efficient corner storage for fc diamond method.
Definition: discretization/facecentered/diamond/geometryhelper.hh:30