13#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH
14#define DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH
18#include <dune/common/exceptions.hh>
20#include <dune/geometry/type.hh>
21#include <dune/geometry/referenceelements.hh>
22#include <dune/geometry/multilineargeometry.hh>
23#include <dune/common/reservedvector.hh>
39 template<
int mydim,
int cdim >
42 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)>;
49template<Dune::GeometryType::Id gt>
55 using Key = std::pair<std::uint8_t, std::uint8_t>;
56 static constexpr std::array<std::array<Key, 2>, 1>
keys = {{
64 using Key = std::pair<std::uint8_t, std::uint8_t>;
65 static constexpr std::array<std::array<Key, 3>, 1>
keys = {{
73 using Key = std::pair<std::uint8_t, std::uint8_t>;
74 static constexpr std::array<std::array<Key, 4>, 1>
keys = {{
82 using Key = std::pair<std::uint8_t, std::uint8_t>;
83 static constexpr std::array<std::array<Key, 4>, 1>
keys = {{
91 using Key = std::pair<std::uint8_t, std::uint8_t>;
92 static constexpr std::array<std::array<Key, 6>, 1>
keys = {{
93 {
Key{0, 1},
Key{2, 1},
Key{3, 1},
Key{1, 1},
Key{4, 1},
Key{5, 1} }
98template<Dune::GeometryType::Id gt>
104 using Key = std::pair<std::uint8_t, std::uint8_t>;
105 static constexpr std::array<std::array<Key, 1>, 1>
keys = {{
113 using Key = std::pair<std::uint8_t, std::uint8_t>;
114 static constexpr std::array<std::array<Key, 2>, 3>
keys = {{
124 using Key = std::pair<std::uint8_t, std::uint8_t>;
125 static constexpr std::array<std::array<Key, 2>, 4>
keys = {{
136 using Key = std::pair<std::uint8_t, std::uint8_t>;
137 static constexpr std::array<std::array<Key, 3>, 4>
keys = {{
148 using Key = std::pair<std::uint8_t, std::uint8_t>;
149 static constexpr std::array<std::array<Key, 3>, 8>
keys = {{
167template <
class Gr
idView,
class ScvType,
class ScvfType>
170 using Scalar =
typename GridView::ctype;
171 using GlobalPosition =
typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
172 using ScvCornerStorage =
typename ScvType::Traits::CornerStorage;
173 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
174 using LocalIndexType =
typename ScvType::Traits::LocalIndexType;
176 using Element =
typename GridView::template Codim<0>::Entity;
177 using Intersection =
typename GridView::Intersection;
179 static constexpr auto dim = GridView::dimension;
180 static constexpr auto dimWorld = GridView::dimensionworld;
187 , boxHelper_(geometry)
193 return getScvCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvIdx);
197 template<
class Transformation>
198 static ScvCornerStorage
getScvCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvIdx)
201 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
202 const auto numBoxScv = ref.size(dim);
204 if (localScvIdx < numBoxScv)
205 return BoxHelper::getScvCorners(type, trans, localScvIdx);
207 const auto localOverlappingScvIdx = localScvIdx-numBoxScv;
208 if (type == Dune::GeometryTypes::triangle)
213 else if (type == Dune::GeometryTypes::quadrilateral)
218 else if (type == Dune::GeometryTypes::tetrahedron)
223 else if (type == Dune::GeometryTypes::hexahedron)
229 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv geometries for dim=" << dim
230 <<
" dimWorld=" << dimWorld
231 <<
" type=" << type);
237 const auto type = geo_.type();
238 const auto numBoxScv = boxHelper_.numScv();
239 if (localScvIdx < numBoxScv)
240 return Dune::GeometryTypes::cube(dim);
241 else if (type == Dune::GeometryTypes::simplex(dim))
242 return Dune::GeometryTypes::simplex(dim);
243 else if (type == Dune::GeometryTypes::quadrilateral)
244 return Dune::GeometryTypes::quadrilateral;
245 else if (type == Dune::GeometryTypes::hexahedron)
246 return Dune::GeometryTypes::none(dim);
248 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv geometries for dim=" << dim
249 <<
" dimWorld=" << dimWorld
250 <<
" type=" << type);
256 return getScvfCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvfIdx);
260 template<
class Transformation>
261 static ScvfCornerStorage
getScvfCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvfIdx)
264 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
265 const auto numBoxScvf = ref.size(dim-1);
267 if (localScvfIdx < numBoxScvf)
268 return BoxHelper::getScvfCorners(type, trans, localScvfIdx);
270 const auto localOverlappingScvfIdx = localScvfIdx-numBoxScvf;
271 if (type == Dune::GeometryTypes::triangle)
276 else if (type == Dune::GeometryTypes::quadrilateral)
281 else if (type == Dune::GeometryTypes::tetrahedron)
286 else if (type == Dune::GeometryTypes::hexahedron)
292 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scvf geometries for dim=" << dim
293 <<
" dimWorld=" << dimWorld
294 <<
" type=" << type);
299 const auto numBoxScvf = boxHelper_.numInteriorScvf();
300 if (localScvfIdx < numBoxScvf)
301 return Dune::GeometryTypes::cube(dim-1);
303 return Dune::GeometryTypes::simplex(dim-1);
308 unsigned int indexInFacet)
const
310 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
315 return Dune::GeometryTypes::cube(dim-1);
318 template<
int d = dimWorld, std::enable_if_t<(d==3),
int> = 0>
319 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
333 template<
int d = dimWorld, std::enable_if_t<(d==2),
int> = 0>
334 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
337 const auto t = p[1] - p[0];
338 GlobalPosition
normal({-t[1], t[0]});
357 return BoxHelper::numInteriorScvf(type) + Dune::referenceElement<Scalar, dim>(type).size(dim);
363 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
369 return boxHelper_.numScv() + 1;
373 Scalar
scvVolume(
unsigned int localScvIdx,
const ScvCornerStorage& p)
const
376 if constexpr (dim == 3)
377 if (scvType == Dune::GeometryTypes::none(dim))
378 return octahedronVolume_(p);
382 [&](
unsigned int i){
return p[i]; }
389 return Dune::referenceElement<Scalar, dim>(type).size(dim) + 1;
401 return Dune::referenceElement<Scalar, dim>(type).size(iIdx, 1, dim);
407 return Dune::referenceElement<Scalar, dim>(type).subEntity(iIdx, 1, ilocalDofIdx, dim);
410 template<
class DofMapper,
class LocalKey>
411 static auto dofIndex(
const DofMapper& dofMapper,
const Element& element,
const LocalKey& localKey)
413 return dofMapper.subIndex(element, localKey.subEntity(), localKey.codim()) + localKey.index();
419 return geo_.corner(localDofIdx);
421 return geo_.center();
425 template<
class LocalKey>
426 static Element::Geometry::LocalCoordinate
localDofPosition(Dune::GeometryType type,
const LocalKey& localKey)
428 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
433 const auto numEdges = referenceElement(geo_).size(dim-1);
434 if (localScvfIndex < numEdges)
436 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
437 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
442 static_cast<LocalIndexType
>(localScvfIndex-numEdges)
448 const LocalIndexType insideScvIdx
449 =
static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
450 return { insideScvIdx, insideScvIdx };
455 if (localScvfIndex < boxHelper_.numInteriorScvf())
468 if (localScvIndex < boxHelper_.numScv())
475 static Element::Geometry::LocalCoordinate
localScvfCenter(Dune::GeometryType type,
unsigned int localScvfIdx)
481 static Element::Geometry::LocalCoordinate
localBoundaryScvfCenter(Dune::GeometryType type,
unsigned int localFacetIndex,
unsigned int indexInFace)
483 return Dumux::center(BoxHelper::getBoundaryScvfCorners(type, [&](
const auto& local){
return local; }, localFacetIndex, indexInFace));
487 Scalar octahedronVolume_(
const ScvCornerStorage& p)
const
496 const typename Element::Geometry& geo_;
497 BoxHelper boxHelper_;
500template <
class Gr
idView,
class ScvType,
class ScvfType, std::
size_t numCubeBubbleDofs>
503 using Scalar =
typename GridView::ctype;
504 using GlobalPosition =
typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
505 using ScvCornerStorage =
typename ScvType::Traits::CornerStorage;
506 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
507 using LocalIndexType =
typename ScvType::Traits::LocalIndexType;
509 using Element =
typename GridView::template Codim<0>::Entity;
510 using Intersection =
typename GridView::Intersection;
512 static constexpr auto dim = GridView::dimension;
513 static constexpr auto dimWorld = GridView::dimensionworld;
520 , boxHelper_(geometry)
526 return getScvCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvIdx);
530 template<
class Transformation>
531 static ScvCornerStorage
getScvCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvIdx)
534 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
535 const auto numBoxScv = ref.size(dim);
537 if (localScvIdx < numBoxScv)
538 return BoxHelper::getScvCorners(type, trans, localScvIdx);
540 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv corners call for hybrid dofs");
546 const auto numBoxScv = boxHelper_.numScv();
548 if (localScvIdx < numBoxScv)
549 return Dune::GeometryTypes::cube(dim);
551 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv geometry call for hybrid dofs");
557 return getScvfCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvfIdx);
561 template<
class Transformation>
562 static ScvfCornerStorage
getScvfCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvfIdx)
565 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
566 const auto numBoxScvf = ref.size(dim-1);
568 if (localScvfIdx < numBoxScvf)
569 return BoxHelper::getScvfCorners(type, trans, localScvfIdx);
571 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scvf corners call for hybrid dofs");
576 const auto numBoxScvf = boxHelper_.numInteriorScvf();
577 if (localScvfIdx < numBoxScvf)
578 return Dune::GeometryTypes::cube(dim-1);
580 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble interior scvf geometry type call for hybrid dofs");
585 unsigned int indexInFacet)
const
587 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
592 return Dune::GeometryTypes::cube(dim-1);
595 template<
int d = dimWorld, std::enable_if_t<(d==3),
int> = 0>
596 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
610 template<
int d = dimWorld, std::enable_if_t<(d==2),
int> = 0>
611 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
614 const auto t = p[1] - p[0];
615 GlobalPosition
normal({-t[1], t[0]});
634 return BoxHelper::numInteriorScvf(type);
640 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
646 return boxHelper_.numScv();
650 Scalar
scvVolume(
unsigned int localScvIdx,
const ScvCornerStorage& p)
const
656 [&](
unsigned int i){
return p[i]; }
663 const auto numVertexDofs = Dune::referenceElement<Scalar, dim>(type).size(dim);
664 return numVertexDofs + (type.isCube() ? numCubeBubbleDofs : 1);
670 return type.isCube() ? numCubeBubbleDofs : 1;
676 return Dune::referenceElement<Scalar, dim>(type).size(iIdx, 1, dim);
682 return Dune::referenceElement<Scalar, dim>(type).subEntity(iIdx, 1, ilocalDofIdx, dim);
685 template<
class DofMapper,
class LocalKey>
686 static auto dofIndex(
const DofMapper& dofMapper,
const Element& element,
const LocalKey& localKey)
689 return dofMapper.subIndex(element, localKey.subEntity(), localKey.codim()) + localKey.index();
694 const auto numVertexDofs = Dune::referenceElement<Scalar, dim>(geo_.type()).size(dim);
695 if (localDofIdx < numVertexDofs)
696 return geo_.corner(localDofIdx);
698 return geo_.center();
702 template<
class LocalKey>
703 static Element::Geometry::LocalCoordinate
localDofPosition(Dune::GeometryType type,
const LocalKey& localKey)
705 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
710 const auto numEdges = referenceElement(geo_).size(dim-1);
711 if (localScvfIndex < numEdges)
713 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
714 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
717 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv pair call for hybrid dofs");
722 const LocalIndexType insideScvIdx
723 =
static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
724 return { insideScvIdx, insideScvIdx };
737 static Element::Geometry::LocalCoordinate
localScvfCenter(Dune::GeometryType type,
unsigned int localScvfIdx)
743 static Element::Geometry::LocalCoordinate
localBoundaryScvfCenter(Dune::GeometryType type,
unsigned int localFacetIndex,
unsigned int indexInFace)
745 return Dumux::center(BoxHelper::getBoundaryScvfCorners(type, [&](
const auto& local){
return local; }, localFacetIndex, indexInFace));
749 const typename Element::Geometry& geo_;
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Compute the center point of a convex polytope geometry or a random-access container of corner points.
Create sub control volumes and sub control volume face geometries.
Definition boxgeometryhelper.hh:257
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition discretization/pq1bubble/geometryhelper.hh:596
static auto numInteriorScvf(Dune::GeometryType type)
number of interior sub control volume faces
Definition discretization/pq1bubble/geometryhelper.hh:632
static auto numLocalDofsIntersection(Dune::GeometryType type, unsigned int iIdx)
Number of local dofs related to an intersection with index iIdx.
Definition discretization/pq1bubble/geometryhelper.hh:674
std::array< LocalIndexType, 2 > getScvPairForScvf(unsigned int localScvfIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:708
static Element::Geometry::LocalCoordinate localScvfCenter(Dune::GeometryType type, unsigned int localScvfIdx)
local scvf center
Definition discretization/pq1bubble/geometryhelper.hh:737
bool isOverlappingScv(unsigned int localScvIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:733
HybridPQ1BubbleGeometryHelper(const typename Element::Geometry &geometry)
Definition discretization/pq1bubble/geometryhelper.hh:518
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition discretization/pq1bubble/geometryhelper.hh:628
static Element::Geometry::LocalCoordinate localBoundaryScvfCenter(Dune::GeometryType type, unsigned int localFacetIndex, unsigned int indexInFace)
local boundary scvf center
Definition discretization/pq1bubble/geometryhelper.hh:743
static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvfIdx)
Create a vector with the corners of sub control volume faces.
Definition discretization/pq1bubble/geometryhelper.hh:562
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition discretization/pq1bubble/geometryhelper.hh:584
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:543
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:574
static auto localDofIndexIntersection(Dune::GeometryType type, unsigned int iIdx, unsigned int ilocalDofIdx)
Local dof index related to a localDof, with index ilocalDofIdx, on an intersection with index iIdx.
Definition discretization/pq1bubble/geometryhelper.hh:680
static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvIdx)
Create a vector with the scv corners.
Definition discretization/pq1bubble/geometryhelper.hh:531
bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:730
static std::size_t numNonCVLocalDofs(Dune::GeometryType type)
number of hybrid dofs
Definition discretization/pq1bubble/geometryhelper.hh:668
static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey &localKey)
local dof position
Definition discretization/pq1bubble/geometryhelper.hh:703
static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
number of boundary sub control volume faces for face localFacetIndex
Definition discretization/pq1bubble/geometryhelper.hh:638
static auto dofIndex(const DofMapper &dofMapper, const Element &element, const LocalKey &localKey)
Definition discretization/pq1bubble/geometryhelper.hh:686
Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage &p) const
get scv volume
Definition discretization/pq1bubble/geometryhelper.hh:650
std::array< LocalIndexType, 2 > getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:720
GlobalPosition dofPosition(unsigned int localDofIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:692
bool isOverlappingScvf(unsigned int localScvfIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:727
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition discretization/pq1bubble/geometryhelper.hh:555
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition discretization/pq1bubble/geometryhelper.hh:524
std::size_t numScv() const
number of sub control volumes (number of codim-1 entities)
Definition discretization/pq1bubble/geometryhelper.hh:644
static std::size_t numElementDofs(Dune::GeometryType type)
number of element dofs
Definition discretization/pq1bubble/geometryhelper.hh:661
Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:590
PQ1BubbleGeometryHelper(const typename Element::Geometry &geometry)
Definition discretization/pq1bubble/geometryhelper.hh:185
static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
number of boundary sub control volume faces for face localFacetIndex
Definition discretization/pq1bubble/geometryhelper.hh:361
bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:461
bool isOverlappingScvf(unsigned int localScvfIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:453
static std::size_t numNonCVLocalDofs(Dune::GeometryType type)
number of hybrid dofs
Definition discretization/pq1bubble/geometryhelper.hh:393
std::size_t numScv() const
number of sub control volumes (number of codim-1 entities)
Definition discretization/pq1bubble/geometryhelper.hh:367
static auto numInteriorScvf(Dune::GeometryType type)
number of interior sub control volume faces
Definition discretization/pq1bubble/geometryhelper.hh:355
static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvfIdx)
Create a vector with the corners of sub control volume faces.
Definition discretization/pq1bubble/geometryhelper.hh:261
static std::size_t numElementDofs(Dune::GeometryType type)
number of element dofs
Definition discretization/pq1bubble/geometryhelper.hh:387
static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvIdx)
Create a vector with the scv corners.
Definition discretization/pq1bubble/geometryhelper.hh:198
Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:313
static Element::Geometry::LocalCoordinate localBoundaryScvfCenter(Dune::GeometryType type, unsigned int localFacetIndex, unsigned int indexInFace)
local boundary scvf center
Definition discretization/pq1bubble/geometryhelper.hh:481
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition discretization/pq1bubble/geometryhelper.hh:307
std::array< LocalIndexType, 2 > getScvPairForScvf(unsigned int localScvfIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:431
GlobalPosition dofPosition(unsigned int localDofIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:416
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:297
static auto localDofIndexIntersection(Dune::GeometryType type, unsigned int iIdx, unsigned int ilocalDofIdx)
Local dof index related to a localDof, with index ilocalDofIdx, on an intersection with index iIdx.
Definition discretization/pq1bubble/geometryhelper.hh:405
bool isOverlappingScv(unsigned int localScvIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:466
static auto numLocalDofsIntersection(Dune::GeometryType type, unsigned int iIdx)
Number of local dofs related to an intersection with index iIdx.
Definition discretization/pq1bubble/geometryhelper.hh:399
Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage &p) const
get scv volume
Definition discretization/pq1bubble/geometryhelper.hh:373
static auto dofIndex(const DofMapper &dofMapper, const Element &element, const LocalKey &localKey)
Definition discretization/pq1bubble/geometryhelper.hh:411
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition discretization/pq1bubble/geometryhelper.hh:351
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition discretization/pq1bubble/geometryhelper.hh:254
std::array< LocalIndexType, 2 > getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:446
static Element::Geometry::LocalCoordinate localScvfCenter(Dune::GeometryType type, unsigned int localScvfIdx)
local scvf center
Definition discretization/pq1bubble/geometryhelper.hh:475
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition discretization/pq1bubble/geometryhelper.hh:319
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition discretization/pq1bubble/geometryhelper.hh:191
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:234
static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey &localKey)
local dof position
Definition discretization/pq1bubble/geometryhelper.hh:426
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
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:700
auto convexPolytopeVolume(Dune::GeometryType type, const CornerF &c)
Compute the volume of several common geometry types.
Definition volume.hh:41
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 keyToCornerStorage(const ReferenceElement &ref, Transformation &&trans, const std::array< T, N > &key)
Definition boxgeometryhelper.hh:228
Definition discretization/pq1bubble/geometryhelper.hh:47
Definition common/pdesolver.hh:24
static constexpr std::array< std::array< Key, 6 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:92
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:91
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:55
static constexpr std::array< std::array< Key, 2 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:56
static constexpr std::array< std::array< Key, 4 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:74
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:73
static constexpr std::array< std::array< Key, 4 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:83
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:82
static constexpr std::array< std::array< Key, 3 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:65
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:64
Definition discretization/pq1bubble/geometryhelper.hh:50
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:148
static constexpr std::array< std::array< Key, 3 >, 8 > keys
Definition discretization/pq1bubble/geometryhelper.hh:149
static constexpr std::array< std::array< Key, 1 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:105
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:104
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:124
static constexpr std::array< std::array< Key, 2 >, 4 > keys
Definition discretization/pq1bubble/geometryhelper.hh:125
static constexpr std::array< std::array< Key, 3 >, 4 > keys
Definition discretization/pq1bubble/geometryhelper.hh:137
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:136
static constexpr std::array< std::array< Key, 2 >, 3 > keys
Definition discretization/pq1bubble/geometryhelper.hh:114
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:113
Definition discretization/pq1bubble/geometryhelper.hh:99
Definition discretization/pq1bubble/geometryhelper.hh:41
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<< mydim)> Type
Definition discretization/pq1bubble/geometryhelper.hh:42
Traits for an efficient corner storage for the PQ1Bubble method.
Definition discretization/pq1bubble/geometryhelper.hh:35
Compute the volume of several common geometry types.