version 3.11-dev
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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_GEOMETRY_HELPER_HH
13#define DUMUX_DISCRETIZATION_FACECENTERED_DIAMOND_GEOMETRY_HELPER_HH
14
15#include <array>
16
17#include <dune/common/reservedvector.hh>
18#include <dune/common/fvector.hh>
19#include <dune/geometry/multilineargeometry.hh>
20#include <dune/geometry/type.hh>
21
22#include <dumux/common/math.hh>
24
25namespace Dumux {
26
28template <class ct>
29struct FCDiamondMLGeometryTraits : public Dune::MultiLinearGeometryTraits<ct>
30{
31 // we use static vectors to store the corners as we know
32 // the maximum number of corners in advance (2^dim)
33 template< int mydim, int cdim >
35 {
36 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)>;
37 };
38};
39
41
42template<Dune::GeometryType::Id gt>
44
45template<>
46struct ScvCorners<Dune::GeometryTypes::triangle>
47{
48 static constexpr Dune::GeometryType type()
49 { return Dune::GeometryTypes::triangle; }
50
51 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
52 static constexpr std::array<std::array<Key, 3>, 3> keys = {{
53 { Key{0, 2}, Key{1, 2}, Key{0, 0} },
54 { Key{2, 2}, Key{0, 2}, Key{0, 0} },
55 { Key{1, 2}, Key{2, 2}, Key{0, 0} }
56 }};
57};
58
59template<>
60struct ScvCorners<Dune::GeometryTypes::quadrilateral>
61{
62 static constexpr Dune::GeometryType type()
63 { return Dune::GeometryTypes::triangle; }
64
65 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
66 static constexpr std::array<std::array<Key, 3>, 4> keys = {{
67 { Key{2, 2}, Key{0, 2}, Key{0, 0} },
68 { Key{1, 2}, Key{3, 2}, Key{0, 0} },
69 { Key{0, 2}, Key{1, 2}, Key{0, 0} },
70 { Key{3, 2}, Key{2, 2}, Key{0, 0} }
71 }};
72};
73
74template<>
75struct ScvCorners<Dune::GeometryTypes::tetrahedron>
76{
77 static constexpr Dune::GeometryType type()
78 { return Dune::GeometryTypes::tetrahedron; }
79
80 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
81 static constexpr std::array<std::array<Key, 4>, 4> keys = {{
82 { Key{0, 3}, Key{1, 3}, Key{2, 3}, Key{0, 0} },
83 { Key{1, 3}, Key{0, 3}, Key{3, 3}, Key{0, 0} },
84 { Key{0, 3}, Key{2, 3}, Key{3, 3}, Key{0, 0} },
85 { Key{2, 3}, Key{1, 3}, Key{3, 3}, Key{0, 0} }
86 }};
87};
88
89template<>
90struct ScvCorners<Dune::GeometryTypes::hexahedron>
91{
92 static constexpr Dune::GeometryType type()
93 { return Dune::GeometryTypes::pyramid; }
94
95 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
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} },
102 { Key{6, 3}, Key{7, 3}, Key{4, 3}, Key{5, 3}, Key{0, 0} }
103 }};
104};
105
106template<Dune::GeometryType::Id gt>
108
109template<>
110struct ScvfCorners<Dune::GeometryTypes::triangle>
111{
112 static constexpr Dune::GeometryType type()
113 { return Dune::GeometryTypes::line; }
114
115 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
116 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
117 { Key{0, 2}, Key{0, 0} },
118 { Key{1, 2}, Key{0, 0} },
119 { Key{2, 2}, Key{0, 0} }
120 }};
121};
122
123template<>
124struct ScvfCorners<Dune::GeometryTypes::quadrilateral>
125{
126 static constexpr Dune::GeometryType type()
127 { return Dune::GeometryTypes::line; }
128
129 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
130 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
131 { Key{0, 2}, Key{0, 0} },
132 { Key{1, 2}, Key{0, 0} },
133 { Key{2, 2}, Key{0, 0} },
134 { Key{3, 2}, Key{0, 0} }
135 }};
136};
137
138template<>
139struct ScvfCorners<Dune::GeometryTypes::tetrahedron>
140{
141 static constexpr Dune::GeometryType type()
142 { return Dune::GeometryTypes::triangle; }
143
144 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
145 static constexpr std::array<std::array<Key, 3>, 6> keys = {{
146 { Key{0, 3}, Key{1, 3}, Key{0, 0} },
147 { Key{2, 3}, Key{0, 3}, Key{0, 0} },
148 { Key{1, 3}, Key{2, 3}, Key{0, 0} },
149 { Key{0, 3}, Key{3, 3}, Key{0, 0} },
150 { Key{1, 3}, Key{3, 3}, Key{0, 0} },
151 { Key{2, 3}, Key{3, 3}, Key{0, 0} }
152 }};
153};
154
155template<>
156struct ScvfCorners<Dune::GeometryTypes::hexahedron>
157{
158 static constexpr Dune::GeometryType type()
159 { return Dune::GeometryTypes::triangle; }
160
161 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
162 static constexpr std::array<std::array<Key, 3>, 12> keys = {{
163 { Key{0, 3}, Key{4, 3}, Key{0, 0} },
164 { Key{1, 3}, Key{5, 3}, Key{0, 0} },
165 { Key{2, 3}, Key{6, 3}, Key{0, 0} },
166 { Key{3, 3}, Key{7, 3}, Key{0, 0} },
167 { Key{0, 3}, Key{2, 3}, Key{0, 0} },
168 { Key{1, 3}, Key{3, 3}, Key{0, 0} },
169 { Key{0, 3}, Key{1, 3}, Key{0, 0} },
170 { Key{2, 3}, Key{3, 3}, Key{0, 0} },
171 { Key{4, 3}, Key{6, 3}, Key{0, 0} },
172 { Key{5, 3}, Key{7, 3}, Key{0, 0} },
173 { Key{4, 3}, Key{5, 3}, Key{0, 0} },
174 { Key{6, 3}, Key{7, 3}, Key{0, 0} }
175 }};
176};
177
178// convert key array to corner storage
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...>)
181{
182 // key is a pair of a local sub-entity index (first) and the sub-entity's codim (second)
183 return { trans(ref.position(key[I].first, key[I].second))... };
184}
185
186// convert key array to corner storage
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)
189{
190 return keyToCornerStorageImpl<S>(ref, trans, key, Indices{});
191}
192
193// boundary corners for the i-th facet
194template<class S, int dim, class ReferenceElement, class Transformation, std::size_t... ii>
195S boundaryCornerStorageImpl(const ReferenceElement& ref, Transformation&& trans, unsigned int i, std::index_sequence<ii...>)
196{
197 // simply the vertices of the facet i
198 return { trans(ref.position(ref.subEntity(i, 1, ii, dim), dim))... };
199}
200
201// boundary corners for the i-th facet
202template<class S, std::size_t numCorners, int dim, class ReferenceElement, class Transformation>
203S boundaryCornerStorage(const ReferenceElement& ref, Transformation&& trans, unsigned int i)
204{
205 return boundaryCornerStorageImpl<S, dim>(ref, trans, i, std::make_index_sequence<numCorners>{});
206}
207
208template<class IndexType, Dune::GeometryType::Id gt>
210
211template<class IndexType>
212struct InsideOutsideScv<IndexType, Dune::GeometryTypes::triangle>
213{
214 static constexpr std::array<std::array<IndexType, 2>, 3> pairs = {{
215 {0, 1}, {0, 2}, {1, 2}
216 }};
217};
218
219template<class IndexType>
220struct InsideOutsideScv<IndexType, Dune::GeometryTypes::quadrilateral>
221{
222 static constexpr std::array<std::array<IndexType, 2>, 4> pairs = {{
223 {0, 2}, {1, 2}, {0, 3}, {1, 3}
224 }};
225};
226
227template<class IndexType>
228struct InsideOutsideScv<IndexType, Dune::GeometryTypes::tetrahedron>
229{
230 static constexpr std::array<std::array<IndexType, 2>, 6> pairs = {{
231 {0, 1}, {0, 2}, {0, 3}, {1, 2}, {1, 3}, {2, 3},
232 }};
233};
234
235template<class IndexType>
236struct InsideOutsideScv<IndexType, Dune::GeometryTypes::hexahedron>
237{
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}
241 }};
242};
243
244} // end namespace Detail::FCDiamond
245
250template <class GridView, class ScvType, class ScvfType>
252{
253 using Element = typename GridView::template Codim<0>::Entity;
254
255 static constexpr auto dim = GridView::dimension;
256 static constexpr auto dimWorld = GridView::dimensionworld;
257
258 using ScvCornerStorage = typename ScvType::Traits::CornerStorage;
259 using LocalIndexType = typename ScvType::Traits::LocalIndexType;
260 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
261
262 // for the normal
263 using Scalar = typename GridView::ctype;
264 using GlobalPosition = typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
265
266public:
267 explicit DiamondGeometryHelper(const typename Element::Geometry& geo)
268 : geo_(geo)
269 {}
270
272 ScvCornerStorage getScvCorners(unsigned int localFacetIndex) const
273 {
274 return getScvCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localFacetIndex);
275 }
276
278 template<class Transformation>
279 static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localFacetIndex)
280 {
281 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
282 if (type == Dune::GeometryTypes::triangle)
283 {
285 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localFacetIndex]);
286 }
287 else if (type == Dune::GeometryTypes::quadrilateral)
288 {
290 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localFacetIndex]);
291 }
292 else if (type == Dune::GeometryTypes::tetrahedron)
293 {
295 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localFacetIndex]);
296 }
297 else if (type == Dune::GeometryTypes::hexahedron)
298 {
300 return Detail::FCDiamond::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localFacetIndex]);
301 }
302 else
303 DUNE_THROW(Dune::NotImplemented, "Scv geometries for type " << type);
304 }
305
307 ScvfCornerStorage getScvfCorners(unsigned int localEdgeIndex) const
308 {
309 return getScvfCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localEdgeIndex);
310 }
311
313 template<class Transformation>
314 static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localEdgeIndex)
315 {
316 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
317 if (type == Dune::GeometryTypes::triangle)
318 {
320 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localEdgeIndex]);
321 }
322 else if (type == Dune::GeometryTypes::quadrilateral)
323 {
325 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localEdgeIndex]);
326 }
327 else if (type == Dune::GeometryTypes::tetrahedron)
328 {
330 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localEdgeIndex]);
331 }
332 else if (type == Dune::GeometryTypes::hexahedron)
333 {
335 return Detail::FCDiamond::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localEdgeIndex]);
336 }
337 else
338 DUNE_THROW(Dune::NotImplemented, "Scvf geometries for type " << type);
339 }
340
342 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex) const
343 {
344 return getBoundaryScvfCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localFacetIndex);
345 }
346
348 template<class Transformation>
349 static ScvfCornerStorage getBoundaryScvfCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localFacetIndex)
350 {
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);
358 else
359 DUNE_THROW(Dune::NotImplemented, "Boundary scvf geometries for type " << type);
360 }
361
362 GlobalPosition facetCenter(unsigned int localFacetIndex) const
363 {
364 return geo_.global(referenceElement(geo_).position(localFacetIndex, 1));
365 }
366
367 std::array<LocalIndexType, 2> getInsideOutsideScvForScvf(unsigned int localEdgeIndex)
368 {
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)
378 else
379 DUNE_THROW(Dune::NotImplemented, "Inside outside scv pairs for type " << type);
380 }
381
383 std::size_t numInteriorScvf()
384 {
385 return referenceElement(geo_).size(2);
386 }
387
389 std::size_t numScv()
390 {
391 return referenceElement(geo_).size(1);
392 }
393
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)
396 {
397 auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]);
398 normal /= normal.two_norm();
399
400 const auto ref = referenceElement(geo_);
401 const auto v = facetCenter_(scvPair[1], ref) - facetCenter_(scvPair[0], ref);
402
403 const auto s = v*normal;
404 if (std::signbit(s))
405 normal *= -1;
406
407 return normal;
408 }
409
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)
412 {
414 const auto t = p[1] - p[0];
415 GlobalPosition normal({-t[1], t[0]});
416 normal /= normal.two_norm();
417
418 const auto ref = referenceElement(geo_);
419 const auto v = facetCenter_(scvPair[1], ref) - facetCenter_(scvPair[0], ref);
420
421 const auto s = v*normal;
422 if (std::signbit(s))
423 normal *= -1;
424
425 return normal;
426 }
427
428 const typename Element::Geometry& elementGeometry() const
429 { return geo_; }
430
432 template<class LocalKey>
433 static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey& localKey)
434 {
435 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
436 }
437
439 static Element::Geometry::LocalCoordinate localScvfCenter(Dune::GeometryType type, unsigned int localScvfIdx)
440 {
441 return Dumux::center(getScvfCorners(type, [&](const auto& local){ return local; }, localScvfIdx));
442 }
443
445 static Element::Geometry::LocalCoordinate localBoundaryScvfCenter(Dune::GeometryType type, unsigned int localFacetIndex)
446 {
447 return Dumux::center(getBoundaryScvfCorners(type, [&](const auto& local){ return local; }, localFacetIndex));
448 }
449
450private:
451 template<class RefElement>
452 GlobalPosition facetCenter_(unsigned int localFacetIndex, const RefElement& ref) const
453 {
454 return geo_.global(ref.position(localFacetIndex, 1));
455 }
456
457 const typename Element::Geometry& geo_;
458};
459
460} // end namespace Dumux
461
462#endif
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: adapt.hh:17
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