13#ifndef DUMUX_DISCRETIZATION_PQ2_GEOMETRY_HELPER_HH
14#define DUMUX_DISCRETIZATION_PQ2_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>
38 template<
int mydim,
int cdim >
41 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)>;
50template <
class Gr
idView,
class ScvType,
class ScvfType>
53 using Scalar =
typename GridView::ctype;
54 using GlobalPosition =
typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
55 using ScvCornerStorage =
typename ScvType::Traits::CornerStorage;
56 using ScvfCornerStorage =
typename ScvfType::Traits::CornerStorage;
57 using LocalIndexType =
typename ScvType::Traits::LocalIndexType;
59 using Element =
typename GridView::template Codim<0>::Entity;
60 using Intersection =
typename GridView::Intersection;
62 static constexpr auto dim = GridView::dimension;
63 static constexpr auto dimWorld = GridView::dimensionworld;
70 , boxHelper_(geometry)
76 return getScvCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvIdx);
80 template<
class Transformation>
81 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);
87 if (localScvIdx < numBoxScv)
88 return BoxHelper::getScvCorners(type, trans, localScvIdx);
90 DUNE_THROW(Dune::NotImplemented,
"PQ2 scv corners call for hybrid dofs");
96 const auto numBoxScv = boxHelper_.numScv();
98 if (localScvIdx < numBoxScv)
99 return Dune::GeometryTypes::cube(dim);
101 DUNE_THROW(Dune::NotImplemented,
"PQ2 scv geometry call for hybrid dofs");
107 return getScvfCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvfIdx);
111 template<
class Transformation>
112 static ScvfCornerStorage
getScvfCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvfIdx)
115 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
116 const auto numBoxScvf = ref.size(dim-1);
118 if (localScvfIdx < numBoxScvf)
119 return BoxHelper::getScvfCorners(type, trans, localScvfIdx);
121 DUNE_THROW(Dune::NotImplemented,
"PQ2 scvf corners call for hybrid dofs");
126 const auto numBoxScvf = boxHelper_.numInteriorScvf();
127 if (localScvfIdx < numBoxScvf)
128 return Dune::GeometryTypes::cube(dim-1);
130 DUNE_THROW(Dune::NotImplemented,
"PQ2 interior scvf geometry type call for hybrid dofs");
135 unsigned int indexInFacet)
const
137 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
142 return Dune::GeometryTypes::cube(dim-1);
145 template<
int d = dimWorld, std::enable_if_t<(d==3),
int> = 0>
146 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
151 GlobalPosition v = geo_.corner(scvPair[1]) - geo_.corner(scvPair[0]);
160 template<
int d = dimWorld, std::enable_if_t<(d==2),
int> = 0>
161 GlobalPosition
normal(
const ScvfCornerStorage& p,
const std::array<LocalIndexType, 2>& scvPair)
164 const auto t = p[1] - p[0];
165 GlobalPosition
normal({-t[1], t[0]});
168 GlobalPosition v = geo_.corner(scvPair[1]) - geo_.corner(scvPair[0]);
184 return BoxHelper::numInteriorScvf(type);
190 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
196 return boxHelper_.numScv();
200 Scalar
scvVolume(
unsigned int localScvIdx,
const ScvCornerStorage& p)
const
206 [&](
unsigned int i){
return p[i]; }
211 template<
class LocalKey>
214 const auto& refElement = Dune::referenceElement<Scalar, dim>(type);
216 const auto numEntitiesIntersection = refElement.size(iIdx, 1, localKey.codim());
217 for(std::size_t idx=0; idx < numEntitiesIntersection; idx++)
218 if(localKey.subEntity() == refElement.subEntity(iIdx, 1, idx, localKey.codim()))
224 template<
class DofMapper,
class LocalKey>
225 static auto dofIndex(
const DofMapper& dofMapper,
const Element& element,
const LocalKey& localKey)
228 template<
class Geometry,
class LocalKey>
229 static GlobalPosition
dofPosition(
const Geometry& geo,
const LocalKey& localKey)
232 template<
class LocalKey>
237 template<
class LocalKey>
238 static Element::Geometry::LocalCoordinate
localDofPosition(Dune::GeometryType type,
const LocalKey& localKey)
240 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
245 const auto numBoxFaces = boxHelper_.numInteriorScvf();
246 if (localScvfIndex < numBoxFaces)
249 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
250 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
254 DUNE_THROW(Dune::NotImplemented,
"PQ2 scv pair call for hybrid dofs");
259 const auto numBoxScvf = referenceElement(geo_).size(localFacetIndex, 1, dim);
260 if (localIsScvfIndex < numBoxScvf)
262 const LocalIndexType insideScvIdx
263 =
static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
264 return { insideScvIdx, insideScvIdx };
267 DUNE_THROW(Dune::NotImplemented,
"PQ2 scv boundary pair call for hybrid dofs");
282 const typename Element::Geometry& geo_;
283 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/pq2/geometryhelper.hh:140
std::size_t numScv() const
number of sub control volumes (number of codim-1 entities)
Definition discretization/pq2/geometryhelper.hh:194
bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
Definition discretization/pq2/geometryhelper.hh:274
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition discretization/pq2/geometryhelper.hh:105
GlobalPosition dofPosition(const LocalKey &localKey) const
Definition discretization/pq2/geometryhelper.hh:233
bool isOverlappingScvf(unsigned int localScvfIndex) const
Definition discretization/pq2/geometryhelper.hh:271
static auto dofIndex(const DofMapper &dofMapper, const Element &element, const LocalKey &localKey)
Definition discretization/pq2/geometryhelper.hh:225
static auto localDofOnIntersection(Dune::GeometryType type, unsigned int iIdx, const LocalKey &localKey)
Local dof index related to a localDof, with index ilocalDofIdx, on an intersection with index iIdx.
Definition discretization/pq2/geometryhelper.hh:212
std::array< LocalIndexType, 2 > getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
Definition discretization/pq2/geometryhelper.hh:257
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition discretization/pq2/geometryhelper.hh:74
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition discretization/pq2/geometryhelper.hh:178
std::array< LocalIndexType, 2 > getScvPairForScvf(unsigned int localScvfIndex) const
Definition discretization/pq2/geometryhelper.hh:243
static auto numInteriorScvf(Dune::GeometryType type)
number of interior sub control volume faces
Definition discretization/pq2/geometryhelper.hh:182
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition discretization/pq2/geometryhelper.hh:134
Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage &p) const
get scv volume
Definition discretization/pq2/geometryhelper.hh:200
static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
number of boundary sub control volume faces for face localFacetIndex
Definition discretization/pq2/geometryhelper.hh:188
static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey &localKey)
local dof position
Definition discretization/pq2/geometryhelper.hh:238
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition discretization/pq2/geometryhelper.hh:124
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition discretization/pq2/geometryhelper.hh:93
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition discretization/pq2/geometryhelper.hh:146
bool isOverlappingScv(unsigned int localScvIndex) const
Definition discretization/pq2/geometryhelper.hh:278
HybridPQ2GeometryHelper(const typename Element::Geometry &geometry)
Definition discretization/pq2/geometryhelper.hh:68
static GlobalPosition dofPosition(const Geometry &geo, const LocalKey &localKey)
Definition discretization/pq2/geometryhelper.hh:229
static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvfIdx)
Create a vector with the corners of sub control volume faces.
Definition discretization/pq2/geometryhelper.hh:112
static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvIdx)
Create a vector with the scv corners.
Definition discretization/pq2/geometryhelper.hh:81
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-2 Lagrange elements.
static std::size_t dofIndex(const DofMapper &m, const Element &e, const LocalKey &lk, const IdSet &)
Definition pq2/dofhelper.hh:44
static GlobalPosition dofPosition(const Geometry &geo, const LocalKey &lk)
Physical position of the DOF in global coordinates.
Definition pq2/dofhelper.hh:55
Definition discretization/pq2/geometryhelper.hh:40
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<< mydim)> Type
Definition discretization/pq2/geometryhelper.hh:41
Traits for an efficient corner storage for the PQ2 method.
Definition discretization/pq2/geometryhelper.hh:35
Compute the volume of several common geometry types.