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> 
   37    template< 
int mydim, 
int cdim >
 
   40        using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)+1>;
 
 
 
   46template<Dune::GeometryType::Id gt>
 
   52    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
   53    static constexpr std::array<std::array<Key, 2>, 1> 
keys = {{
 
 
 
   61    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
   62    static constexpr std::array<std::array<Key, 3>, 1> 
keys = {{
 
 
 
   70    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
   71    static constexpr std::array<std::array<Key, 4>, 1> 
keys = {{
 
 
 
   79    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
   80    static constexpr std::array<std::array<Key, 4>, 1> 
keys = {{
 
 
 
   88    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
   89    static constexpr std::array<std::array<Key, 6>, 1> 
keys = {{
 
   90        { 
Key{0, 1}, 
Key{2, 1}, 
Key{3, 1}, 
Key{1, 1}, 
Key{4, 1}, 
Key{5, 1} }
 
 
 
   95template<Dune::GeometryType::Id gt>
 
  101    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
  102    static constexpr std::array<std::array<Key, 1>, 1> 
keys = {{
 
 
 
  110    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
  111    static constexpr std::array<std::array<Key, 2>, 3> 
keys = {{
 
 
 
  121    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
  122    static constexpr std::array<std::array<Key, 2>, 4> 
keys = {{
 
 
 
  133    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
  134    static constexpr std::array<std::array<Key, 3>, 10> 
keys = {{
 
 
 
  145    using Key = std::pair<std::uint8_t, std::uint8_t>; 
 
  146    static constexpr std::array<std::array<Key, 3>, 8> 
keys = {{
 
 
 
 
  164template <
class Gr
idView, 
class ScvType, 
class ScvfType>
 
  167    using Scalar = 
typename GridView::ctype;
 
  168    using GlobalPosition = 
typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
 
  169    using ScvCornerStorage = 
typename ScvType::Traits::CornerStorage;
 
  170    using ScvfCornerStorage = 
typename ScvfType::Traits::CornerStorage;
 
  171    using LocalIndexType = 
typename ScvType::Traits::LocalIndexType;
 
  173    using Element = 
typename GridView::template Codim<0>::Entity;
 
  174    using Intersection = 
typename GridView::Intersection;
 
  176    static constexpr auto dim = GridView::dimension;
 
  177    static constexpr auto dimWorld = GridView::dimensionworld;
 
  182    , boxHelper_(geometry)
 
 
  189        const auto type = geo_.type();
 
  190        const auto numBoxScv = boxHelper_.numScv();
 
  192        if (localScvIdx < numBoxScv)
 
  193            return boxHelper_.getScvCorners(localScvIdx);
 
  195        const auto localOverlappingScvIdx = localScvIdx-numBoxScv;
 
  196        if (type == Dune::GeometryTypes::triangle)
 
  201        else if (type == Dune::GeometryTypes::quadrilateral)
 
  206        else if (type == Dune::GeometryTypes::tetrahedron)
 
  211        else if (type == Dune::GeometryTypes::hexahedron)
 
  217            DUNE_THROW(Dune::NotImplemented, 
"PQ1Bubble scv geometries for dim=" << dim
 
  218                                                            << 
" dimWorld=" << dimWorld
 
  219                                                            << 
" type=" << type);
 
 
  225        const auto type = geo_.type();
 
  226        const auto numBoxScv = boxHelper_.numScv();
 
  227        if (localScvIdx < numBoxScv)
 
  228            return Dune::GeometryTypes::cube(dim);
 
  229        else if (type == Dune::GeometryTypes::simplex(dim))
 
  230            return Dune::GeometryTypes::simplex(dim);
 
  231        else if (type == Dune::GeometryTypes::quadrilateral)
 
  232            return Dune::GeometryTypes::quadrilateral;
 
  233        else if (type == Dune::GeometryTypes::hexahedron)
 
  234            return Dune::GeometryTypes::none(dim); 
 
  236            DUNE_THROW(Dune::NotImplemented, 
"PQ1Bubble scv geometries for dim=" << dim
 
  237                                                            << 
" dimWorld=" << dimWorld
 
  238                                                            << 
" type=" << type);
 
 
  245        const auto type = geo_.type();
 
  246        const auto numBoxScvf = boxHelper_.numInteriorScvf();
 
  248        if (localScvfIdx < numBoxScvf)
 
  249            return boxHelper_.getScvfCorners(localScvfIdx);
 
  251        const auto localOverlappingScvfIdx = localScvfIdx-numBoxScvf;
 
  252        if (type == Dune::GeometryTypes::triangle)
 
  257        else if (type == Dune::GeometryTypes::quadrilateral)
 
  262        else if (type == Dune::GeometryTypes::tetrahedron)
 
  267        else if (type == Dune::GeometryTypes::hexahedron)
 
  273            DUNE_THROW(Dune::NotImplemented, 
"PQ1Bubble scvf geometries for dim=" << dim
 
  274                                                            << 
" dimWorld=" << dimWorld
 
  275                                                            << 
" type=" << type);
 
 
  280        const auto numBoxScvf = boxHelper_.numInteriorScvf();
 
  281        if (localScvfIdx < numBoxScvf)
 
  282            return Dune::GeometryTypes::cube(dim-1);
 
  284            return Dune::GeometryTypes::simplex(dim-1);
 
 
  289                                             unsigned int indexInFacet)
 const 
  291        return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
 
 
  296        return Dune::GeometryTypes::cube(dim-1);
 
 
  299    template<
int d = dimWorld, std::enable_if_t<(d==3), 
int> = 0>
 
  300    GlobalPosition 
normal(
const ScvfCornerStorage& p, 
const std::array<LocalIndexType, 2>& scvPair)
 
 
  314    template<
int d = dimWorld, std::enable_if_t<(d==2), 
int> = 0>
 
  315    GlobalPosition 
normal(
const ScvfCornerStorage& p, 
const std::array<LocalIndexType, 2>& scvPair)
 
  318        const auto t = p[1] - p[0];
 
  319        GlobalPosition 
normal({-t[1], t[0]});
 
 
  338        return boxHelper_.numInteriorScvf() + referenceElement(geo_).size(dim);
 
 
  344        return referenceElement(geo_).size(localFacetIndex, 1, dim);
 
 
  350        return boxHelper_.numScv() + 1;
 
 
  354    Scalar 
scvVolume(
unsigned int localScvIdx, 
const ScvCornerStorage& p)
 const 
  357        if constexpr (dim == 3)
 
  358            if (scvType == Dune::GeometryTypes::none(dim))
 
  359                return octahedronVolume_(p);
 
  363            [&](
unsigned int i){ 
return p[i]; }
 
 
  367    template<
class DofMapper>
 
  368    auto dofIndex(
const DofMapper& dofMapper, 
const Element& element, 
unsigned int localScvIdx)
 const 
  370        if (localScvIdx < 
numScv()-1)
 
  371            return dofMapper.subIndex(element, localScvIdx, dim);
 
  373            return dofMapper.index(element);
 
 
  378        if (localScvIdx < 
numScv()-1)
 
  379            return geo_.corner(localScvIdx);
 
  381            return geo_.center();
 
 
  386        const auto numEdges = referenceElement(geo_).size(dim-1);
 
  387        if (localScvfIndex < numEdges)
 
  389                static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
 
  390                static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
 
  394                static_cast<LocalIndexType
>(
numScv()-1),
 
  395                static_cast<LocalIndexType
>(localScvfIndex-numEdges)
 
 
  401        const LocalIndexType insideScvIdx
 
  402            = 
static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
 
  403        return { insideScvIdx, insideScvIdx };
 
 
  408        if (localScvfIndex < boxHelper_.numInteriorScvf())
 
 
  421        if (localScvIndex < boxHelper_.numScv())
 
 
  428    Scalar octahedronVolume_(
const ScvCornerStorage& p)
 const 
  437    const typename Element::Geometry& geo_; 
 
 
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Create sub control volumes and sub control volume face geometries.
Definition boxgeometryhelper.hh:257
PQ1BubbleGeometryHelper(const typename Element::Geometry &geometry)
Definition discretization/pq1bubble/geometryhelper.hh:180
bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:414
bool isOverlappingScvf(unsigned int localScvfIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:406
auto dofIndex(const DofMapper &dofMapper, const Element &element, unsigned int localScvIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:368
std::size_t numScv() const
number of sub control volumes (number of codim-1 entities)
Definition discretization/pq1bubble/geometryhelper.hh:348
std::size_t numBoundaryScvf(unsigned int localFacetIndex) const
number of boundary sub control volume faces for face localFacetIndex
Definition discretization/pq1bubble/geometryhelper.hh:342
Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:294
std::size_t numInteriorScvf() const
number of interior sub control volume faces
Definition discretization/pq1bubble/geometryhelper.hh:336
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition discretization/pq1bubble/geometryhelper.hh:288
GlobalPosition dofPosition(unsigned int localScvIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:376
std::array< LocalIndexType, 2 > getScvPairForScvf(unsigned int localScvfIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:384
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:278
bool isOverlappingScv(unsigned int localScvIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:419
Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage &p) const
get scv volume
Definition discretization/pq1bubble/geometryhelper.hh:354
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition discretization/pq1bubble/geometryhelper.hh:332
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition discretization/pq1bubble/geometryhelper.hh:242
std::array< LocalIndexType, 2 > getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
Definition discretization/pq1bubble/geometryhelper.hh:399
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition discretization/pq1bubble/geometryhelper.hh:300
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition discretization/pq1bubble/geometryhelper.hh:186
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition discretization/pq1bubble/geometryhelper.hh:222
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
Define some often used mathematical functions.
S keyToCornerStorage(const Geo &geo, const std::array< T, N > &key)
Definition boxgeometryhelper.hh:228
Definition discretization/pq1bubble/geometryhelper.hh:44
Definition common/pdesolver.hh:24
static constexpr std::array< std::array< Key, 6 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:89
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:88
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:52
static constexpr std::array< std::array< Key, 2 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:53
static constexpr std::array< std::array< Key, 4 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:71
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:70
static constexpr std::array< std::array< Key, 4 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:80
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:79
static constexpr std::array< std::array< Key, 3 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:62
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:61
Definition discretization/pq1bubble/geometryhelper.hh:47
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:145
static constexpr std::array< std::array< Key, 3 >, 8 > keys
Definition discretization/pq1bubble/geometryhelper.hh:146
static constexpr std::array< std::array< Key, 1 >, 1 > keys
Definition discretization/pq1bubble/geometryhelper.hh:102
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:101
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:121
static constexpr std::array< std::array< Key, 2 >, 4 > keys
Definition discretization/pq1bubble/geometryhelper.hh:122
static constexpr std::array< std::array< Key, 3 >, 10 > keys
Definition discretization/pq1bubble/geometryhelper.hh:134
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:133
static constexpr std::array< std::array< Key, 2 >, 3 > keys
Definition discretization/pq1bubble/geometryhelper.hh:111
std::pair< std::uint8_t, std::uint8_t > Key
Definition discretization/pq1bubble/geometryhelper.hh:110
Definition discretization/pq1bubble/geometryhelper.hh:96
Definition discretization/pq1bubble/geometryhelper.hh:39
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<< mydim)+1 > Type
Definition discretization/pq1bubble/geometryhelper.hh:40
Traits for an efficient corner storage for the PQ1Bubble method.
Definition discretization/pq1bubble/geometryhelper.hh:34
Compute the volume of several common geometry types.