24#ifndef DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_GEOMETRY_HELPER_HH
25#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_GEOMETRY_HELPER_HH
29#include <dune/common/reservedvector.hh>
30#include <dune/common/fvector.hh>
31#include <dune/geometry/multilineargeometry.hh>
32#include <dune/geometry/type.hh>
44 template<
int mydim,
int cdim >
47 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)>;
53template<Dune::GeometryType::Id gt>
59 static constexpr Dune::GeometryType
type()
60 {
return Dune::GeometryTypes::triangle; }
62 using Key = std::pair<std::uint8_t, std::uint8_t>;
63 static constexpr std::array<std::array<Key, 3>, 3> keys = {{
73 static constexpr Dune::GeometryType
type()
74 {
return Dune::GeometryTypes::triangle; }
76 using Key = std::pair<std::uint8_t, std::uint8_t>;
77 static constexpr std::array<std::array<Key, 3>, 4> keys = {{
88 static constexpr Dune::GeometryType
type()
89 {
return Dune::GeometryTypes::tetrahedron; }
91 using Key = std::pair<std::uint8_t, std::uint8_t>;
92 static constexpr std::array<std::array<Key, 4>, 4> keys = {{
103 static constexpr Dune::GeometryType
type()
104 {
return Dune::GeometryTypes::pyramid; }
106 using Key = std::pair<std::uint8_t, std::uint8_t>;
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} },
117template<Dune::GeometryType::Id gt>
123 static constexpr Dune::GeometryType
type()
126 using Key = std::pair<std::uint8_t, std::uint8_t>;
127 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
137 static constexpr Dune::GeometryType
type()
140 using Key = std::pair<std::uint8_t, std::uint8_t>;
141 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
152 static constexpr Dune::GeometryType
type()
153 {
return Dune::GeometryTypes::triangle; }
155 using Key = std::pair<std::uint8_t, std::uint8_t>;
156 static constexpr std::array<std::array<Key, 3>, 6> keys = {{
169 static constexpr Dune::GeometryType
type()
170 {
return Dune::GeometryTypes::triangle; }
172 using Key = std::pair<std::uint8_t, std::uint8_t>;
173 static constexpr std::array<std::array<Key, 3>, 12> keys = {{
190template<
class S,
class Geo,
class KeyArray, std::size_t... I>
193 using Dune::referenceElement;
194 const auto ref = referenceElement(geo);
196 return { geo.global(ref.position(key[I].first, key[I].second))... };
200template<
class S,
class Geo,
class T, std::
size_t N,
class Indices = std::make_index_sequence<N>>
203 return keyToCornerStorageImpl<S>(geo, key, Indices{});
207template<
class S,
class Geo, std::size_t... ii>
210 using Dune::referenceElement;
211 const auto ref = referenceElement(geo);
213 return { geo.global(ref.position(ref.subEntity(i, 1, ii, Geo::mydimension), Geo::mydimension))... };
217template<
class S, std::
size_t numCorners,
class Geo>
220 return boundaryCornerStorageImpl<S>(geo, i, std::make_index_sequence<numCorners>{});
223template<
class IndexType, Dune::GeometryType::Id gt>
226template<
class IndexType>
229 static constexpr std::array<std::array<IndexType, 2>, 3> pairs = {{
230 {0, 1}, {0, 2}, {1, 2}
234template<
class IndexType>
237 static constexpr std::array<std::array<IndexType, 2>, 4> pairs = {{
238 {0, 2}, {1, 2}, {0, 3}, {1, 3}
242template<
class IndexType>
245 static constexpr std::array<std::array<IndexType, 2>, 6> pairs = {{
246 {0, 1}, {0, 2}, {0, 3}, {1, 2}, {1, 3}, {2, 3},
250template<
class IndexType>
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}
265template <
class Gr
idView,
class ScvType,
class ScvfType>
268 using Element =
typename GridView::template Codim<0>::Entity;
270 static constexpr auto dim = GridView::dimension;
271 static constexpr auto dimWorld = GridView::dimensionworld;
273 using ScvCornerStorage =
typename ScvType::Traits::CornerStorage;
274 using LocalIndexType =
typename ScvType::Traits::LocalIndexType;
275 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
278 using ctype =
typename GridView::ctype;
279 using GlobalPosition =
typename Dune::FieldVector<ctype, GridView::dimensionworld>;
289 const auto type = geo_.type();
290 if (type == Dune::GeometryTypes::triangle)
293 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localFacetIndex]);
295 else if (type == Dune::GeometryTypes::quadrilateral)
298 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localFacetIndex]);
300 else if (type == Dune::GeometryTypes::tetrahedron)
303 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localFacetIndex]);
305 else if (type == Dune::GeometryTypes::hexahedron)
308 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localFacetIndex]);
311 DUNE_THROW(Dune::NotImplemented,
"Scv geometries for type " << type);
317 const auto type = geo_.type();
318 if (type == Dune::GeometryTypes::triangle)
321 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localEdgeIndex]);
323 else if (type == Dune::GeometryTypes::quadrilateral)
326 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localEdgeIndex]);
328 else if (type == Dune::GeometryTypes::tetrahedron)
331 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localEdgeIndex]);
333 else if (type == Dune::GeometryTypes::hexahedron)
336 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localEdgeIndex]);
339 DUNE_THROW(Dune::NotImplemented,
"Scvf geometries for type " << type);
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);
353 DUNE_THROW(Dune::NotImplemented,
"Boundary scvf geometries for type " << type);
358 return geo_.global(referenceElement(geo_).position(localFacetIndex, 1));
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)
373 DUNE_THROW(Dune::NotImplemented,
"Inside outside scv pairs for type " << type);
379 return referenceElement(geo_).size(2);
385 return referenceElement(geo_).size(1);
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)
394 const auto ref = referenceElement(geo_);
395 const auto v = facetCenter_(scvPair[1], ref) - facetCenter_(scvPair[0], ref);
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)
408 const auto t = p[1] - p[0];
409 GlobalPosition
normal({-t[1], t[0]});
412 const auto ref = referenceElement(geo_);
413 const auto v = facetCenter_(scvPair[1], ref) - facetCenter_(scvPair[0], ref);
426 template<
class RefElement>
427 GlobalPosition facetCenter_(
unsigned int localFacetIndex,
const RefElement& ref)
const
429 return geo_.global(ref.position(localFacetIndex, 1));
432 const typename Element::Geometry& geo_;
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