20#ifndef DUMUX_DISCRETIZATION_PQ3_GEOMETRY_HELPER_HH
21#define DUMUX_DISCRETIZATION_PQ3_GEOMETRY_HELPER_HH
25#include <dune/common/exceptions.hh>
26#include <dune/geometry/type.hh>
27#include <dune/geometry/referenceelements.hh>
28#include <dune/geometry/multilineargeometry.hh>
29#include <dune/common/reservedvector.hh>
51template <
class Gr
idView,
class ScvType,
class ScvfType>
54 using Scalar =
typename GridView::ctype;
55 using GlobalPosition =
typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
56 using ScvCornerStorage =
typename ScvType::Traits::CornerStorage;
57 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
58 using LocalIndexType =
typename ScvType::Traits::LocalIndexType;
60 using Element =
typename GridView::template Codim<0>::Entity;
61 using Intersection =
typename GridView::Intersection;
63 static constexpr auto dim = GridView::dimension;
64 static constexpr auto dimWorld = GridView::dimensionworld;
71 , boxHelper_(geometry)
77 return getScvCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvIdx);
81 template<
class Transformation>
82 static ScvCornerStorage
getScvCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvIdx)
84 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
85 const auto numBoxScv = ref.size(dim);
86 if (localScvIdx < numBoxScv)
87 return BoxHelper::getScvCorners(type, trans, localScvIdx);
89 DUNE_THROW(Dune::NotImplemented,
"PQ3 scv corners call for hybrid dofs");
94 const auto numBoxScv = boxHelper_.numScv();
95 if (localScvIdx < numBoxScv)
96 return Dune::GeometryTypes::cube(dim);
98 DUNE_THROW(Dune::NotImplemented,
"PQ3 scv geometry call for hybrid dofs");
104 return getScvfCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvfIdx);
108 template<
class Transformation>
109 static ScvfCornerStorage
getScvfCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvfIdx)
111 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
112 const auto numBoxScvf = ref.size(dim-1);
113 if (localScvfIdx < numBoxScvf)
114 return BoxHelper::getScvfCorners(type, trans, localScvfIdx);
116 DUNE_THROW(Dune::NotImplemented,
"PQ3 scvf corners call for hybrid dofs");
121 const auto numBoxScvf = boxHelper_.numInteriorScvf();
122 if (localScvfIdx < numBoxScvf)
123 return Dune::GeometryTypes::cube(dim-1);
125 DUNE_THROW(Dune::NotImplemented,
"PQ3 interior scvf geometry type call for hybrid dofs");
131 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
136 return Dune::GeometryTypes::cube(dim-1);
139 template<
int d = dimWorld, std::enable_if_t<(d==3),
int> = 0>
140 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
145 GlobalPosition v = geo_.corner(scvPair[1]) - geo_.corner(scvPair[0]);
154 template<
int d = dimWorld, std::enable_if_t<(d==2),
int> = 0>
155 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
157 const auto t = p[1] - p[0];
158 GlobalPosition
normal({-t[1], t[0]});
161 GlobalPosition v = geo_.corner(scvPair[1]) - geo_.corner(scvPair[0]);
177 return BoxHelper::numInteriorScvf(type);
183 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
189 return boxHelper_.numScv();
193 Scalar
scvVolume(
unsigned int localScvIdx,
const ScvCornerStorage& p)
const
198 [&](
unsigned int i){
return p[i]; }
203 template<
class LocalKey>
206 const auto& refElement = Dune::referenceElement<Scalar, dim>(type);
207 const auto numEntitiesIntersection = refElement.size(iIdx, 1, localKey.codim());
208 for (std::size_t idx = 0; idx < numEntitiesIntersection; idx++)
209 if (localKey.subEntity() == refElement.subEntity(iIdx, 1, idx, localKey.codim()))
215 template<
class DofMapper,
class LocalKey,
class GlobalIdSet>
216 static auto dofIndex(
const DofMapper& dofMapper,
const Element& element,
217 const LocalKey& localKey,
const GlobalIdSet& gidSet)
220 template<
class Geometry,
class LocalKey>
221 static GlobalPosition
dofPosition(
const Geometry& geo,
const LocalKey& localKey)
224 template<
class LocalKey>
229 template<
class LocalKey>
230 static typename Element::Geometry::LocalCoordinate
localDofPosition(Dune::GeometryType type,
const LocalKey& localKey)
235 const auto numBoxFaces = boxHelper_.numInteriorScvf();
236 if (localScvfIndex < numBoxFaces)
239 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
240 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
244 DUNE_THROW(Dune::NotImplemented,
"PQ3 scv pair call for hybrid dofs");
249 const auto numBoxScvf = referenceElement(geo_).size(localFacetIndex, 1, dim);
250 if (localIsScvfIndex < numBoxScvf)
252 const LocalIndexType insideScvIdx
253 =
static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
254 return { insideScvIdx, insideScvIdx };
257 DUNE_THROW(Dune::NotImplemented,
"PQ3 scv boundary pair call for hybrid dofs");
270 const typename Element::Geometry& geo_;
271 BoxHelper boxHelper_;
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
Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
Definition discretization/pq3/geometryhelper.hh:134
std::array< LocalIndexType, 2 > getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
Definition discretization/pq3/geometryhelper.hh:247
bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
Definition discretization/pq3/geometryhelper.hh:263
static auto localDofOnIntersection(Dune::GeometryType type, unsigned int iIdx, const LocalKey &localKey)
Local dof index related to a localDof, with index localDofIdx, on an intersection with index iIdx.
Definition discretization/pq3/geometryhelper.hh:204
HybridPQ3GeometryHelper(const typename Element::Geometry &geometry)
Definition discretization/pq3/geometryhelper.hh:69
bool isOverlappingScvf(unsigned int localScvfIndex) const
Definition discretization/pq3/geometryhelper.hh:260
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition discretization/pq3/geometryhelper.hh:102
bool isOverlappingScv(unsigned int localScvIndex) const
Definition discretization/pq3/geometryhelper.hh:266
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition discretization/pq3/geometryhelper.hh:75
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition discretization/pq3/geometryhelper.hh:140
std::size_t numScv() const
number of sub control volumes (one per vertex)
Definition discretization/pq3/geometryhelper.hh:187
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition discretization/pq3/geometryhelper.hh:119
static GlobalPosition dofPosition(const Geometry &geo, const LocalKey &localKey)
Definition discretization/pq3/geometryhelper.hh:221
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition discretization/pq3/geometryhelper.hh:92
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition discretization/pq3/geometryhelper.hh:129
Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage &p) const
get scv volume
Definition discretization/pq3/geometryhelper.hh:193
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition discretization/pq3/geometryhelper.hh:171
GlobalPosition dofPosition(const LocalKey &localKey) const
Definition discretization/pq3/geometryhelper.hh:225
std::array< LocalIndexType, 2 > getScvPairForScvf(unsigned int localScvfIndex) const
Definition discretization/pq3/geometryhelper.hh:233
static auto dofIndex(const DofMapper &dofMapper, const Element &element, const LocalKey &localKey, const GlobalIdSet &gidSet)
Orientation-consistent global DOF index — delegates to PQ3LagrangeDofHelper.
Definition discretization/pq3/geometryhelper.hh:216
static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvfIdx)
Create a vector with the corners of sub control volume faces.
Definition discretization/pq3/geometryhelper.hh:109
static auto numInteriorScvf(Dune::GeometryType type)
number of interior sub control volume faces
Definition discretization/pq3/geometryhelper.hh:175
static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey &localKey)
Local coordinate of a DOF for order-3 Lagrange basis.
Definition discretization/pq3/geometryhelper.hh:230
static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
number of boundary sub control volume faces for face localFacetIndex
Definition discretization/pq3/geometryhelper.hh:181
static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvIdx)
Create a vector with the scv corners.
Definition discretization/pq3/geometryhelper.hh:82
Helper class constructing the dual grid finite volume geometries for the cvfe discretizazion method.
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
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.
Lightweight DOF helper for order-3 Lagrange elements.
static GlobalPosition dofPosition(const Geometry &geo, const LocalKey &lk)
Physical position of a DOF in global coordinates.
Definition pq3/dofhelper.hh:100
static GridView::template Codim< 0 >::Entity::Geometry::LocalCoordinate localDofPos(Dune::GeometryType gt, const LocalKey &lk)
Reference-element position of a DOF for order-3 Lagrange basis.
Definition pq3/dofhelper.hh:106
static std::size_t dofIndex(const DofMapper &m, const Element &e, const LocalKey &lk, const IdSet &idSet)
Orientation-consistent global DOF index.
Definition pq3/dofhelper.hh:64
Compute the volume of several common geometry types.