version 3.11-dev
Loading...
Searching...
No Matches
discretization/pq2/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//
13#ifndef DUMUX_DISCRETIZATION_PQ2_GEOMETRY_HELPER_HH
14#define DUMUX_DISCRETIZATION_PQ2_GEOMETRY_HELPER_HH
15
16#include <array>
17
18#include <dune/common/exceptions.hh>
19
20#include <dune/geometry/type.hh>
21#include <dune/geometry/referenceelements.hh>
22#include <dune/geometry/multilineargeometry.hh>
23#include <dune/common/reservedvector.hh>
24
25#include <dumux/common/math.hh>
29
30namespace Dumux {
31
33template <class ct>
34struct PQ2MLGeometryTraits : public Dune::MultiLinearGeometryTraits<ct>
35{
36 // we use static vectors to store the corners as we know
37 // the maximum number of corners in advance (2^dim for box-type sub-cells)
38 template< int mydim, int cdim >
40 {
41 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)>;
42 };
43};
44
45
50template <class GridView, class ScvType, class ScvfType>
52{
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;
58
59 using Element = typename GridView::template Codim<0>::Entity;
60 using Intersection = typename GridView::Intersection;
61
62 static constexpr auto dim = GridView::dimension;
63 static constexpr auto dimWorld = GridView::dimensionworld;
64
66public:
67
68 HybridPQ2GeometryHelper(const typename Element::Geometry& geometry)
69 : geo_(geometry)
70 , boxHelper_(geometry)
71 {}
72
74 ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
75 {
76 return getScvCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localScvIdx);
77 }
78
80 template<class Transformation>
81 static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localScvIdx)
82 {
83 // proceed according to number of corners of the element
84 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
85 const auto numBoxScv = ref.size(dim);
86 // reuse box geometry helper for the corner scvs
87 if (localScvIdx < numBoxScv)
88 return BoxHelper::getScvCorners(type, trans, localScvIdx);
89
90 DUNE_THROW(Dune::NotImplemented, "PQ2 scv corners call for hybrid dofs");
91 }
92
93 Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
94 {
95 // proceed according to number of corners of the element
96 const auto numBoxScv = boxHelper_.numScv();
97
98 if (localScvIdx < numBoxScv)
99 return Dune::GeometryTypes::cube(dim);
100
101 DUNE_THROW(Dune::NotImplemented, "PQ2 scv geometry call for hybrid dofs");
102 }
103
105 ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
106 {
107 return getScvfCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localScvfIdx);
108 }
109
111 template<class Transformation>
112 static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localScvfIdx)
113 {
114 // proceed according to number of corners
115 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
116 const auto numBoxScvf = ref.size(dim-1);
117 // reuse box geometry helper for scvfs
118 if (localScvfIdx < numBoxScvf)
119 return BoxHelper::getScvfCorners(type, trans, localScvfIdx);
120
121 DUNE_THROW(Dune::NotImplemented, "PQ2 scvf corners call for hybrid dofs");
122 }
123
124 Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
125 {
126 const auto numBoxScvf = boxHelper_.numInteriorScvf();
127 if (localScvfIdx < numBoxScvf)
128 return Dune::GeometryTypes::cube(dim-1);
129
130 DUNE_THROW(Dune::NotImplemented, "PQ2 interior scvf geometry type call for hybrid dofs");
131 }
132
134 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex,
135 unsigned int indexInFacet) const
136 {
137 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
138 }
139
140 Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
141 {
142 return Dune::GeometryTypes::cube(dim-1);
143 }
144
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)
147 {
148 auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]);
149 normal /= normal.two_norm();
150
151 GlobalPosition v = geo_.corner(scvPair[1]) - geo_.corner(scvPair[0]);
152
153 const auto s = v*normal;
154 if (std::signbit(s))
155 normal *= -1;
156
157 return normal;
158 }
159
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)
162 {
164 const auto t = p[1] - p[0];
165 GlobalPosition normal({-t[1], t[0]});
166 normal /= normal.two_norm();
167
168 GlobalPosition v = geo_.corner(scvPair[1]) - geo_.corner(scvPair[0]);
169
170 const auto s = v*normal;
171 if (std::signbit(s))
172 normal *= -1;
173
174 return normal;
175 }
176
178 const typename Element::Geometry& elementGeometry() const
179 { return geo_; }
180
182 static auto numInteriorScvf(Dune::GeometryType type)
183 {
184 return BoxHelper::numInteriorScvf(type);
185 }
186
188 static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
189 {
190 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
191 }
192
194 std::size_t numScv() const
195 {
196 return boxHelper_.numScv();
197 }
198
200 Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage& p) const
201 {
202 const auto scvType = getScvGeometryType(localScvIdx);
203
205 scvType,
206 [&](unsigned int i){ return p[i]; }
207 );
208 }
209
211 template<class LocalKey>
212 static auto localDofOnIntersection(Dune::GeometryType type, unsigned int iIdx, const LocalKey& localKey)
213 {
214 const auto& refElement = Dune::referenceElement<Scalar, dim>(type);
215
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()))
219 return true;
220
221 return false;
222 }
223
224 template<class DofMapper, class LocalKey>
225 static auto dofIndex(const DofMapper& dofMapper, const Element& element, const LocalKey& localKey)
226 { return PQ2LagrangeDofHelper<GridView>::dofIndex(dofMapper, element, localKey); }
227
228 template<class Geometry, class LocalKey>
229 static GlobalPosition dofPosition(const Geometry& geo, const LocalKey& localKey)
230 { return PQ2LagrangeDofHelper<GridView>::dofPosition(geo, localKey); }
231
232 template<class LocalKey>
233 GlobalPosition dofPosition(const LocalKey& localKey) const
234 { return dofPosition(geo_, localKey); }
235
237 template<class LocalKey>
238 static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey& localKey)
239 {
240 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
241 }
242
243 std::array<LocalIndexType, 2> getScvPairForScvf(unsigned int localScvfIndex) const
244 {
245 const auto numBoxFaces = boxHelper_.numInteriorScvf();
246 if (localScvfIndex < numBoxFaces)
247 {
248 return {
249 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
250 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
251 };
252 }
253
254 DUNE_THROW(Dune::NotImplemented, "PQ2 scv pair call for hybrid dofs");
255 }
256
257 std::array<LocalIndexType, 2> getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
258 {
259 const auto numBoxScvf = referenceElement(geo_).size(localFacetIndex, 1, dim);
260 if (localIsScvfIndex < numBoxScvf)
261 {
262 const LocalIndexType insideScvIdx
263 = static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
264 return { insideScvIdx, insideScvIdx };
265 }
266
267 DUNE_THROW(Dune::NotImplemented, "PQ2 scv boundary pair call for hybrid dofs");
268 }
269
270 // For hybrid dofs we don't construct scvfs
271 bool isOverlappingScvf(unsigned int localScvfIndex) const
272 { return false; }
273
274 bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
275 { return false; }
276
277 // For hybrid dofs we don't construct scvs
278 bool isOverlappingScv(unsigned int localScvIndex) const
279 { return false; }
280
281private:
282 const typename Element::Geometry& geo_;
283 BoxHelper boxHelper_;
284};
285
286} // end namespace Dumux
287
288#endif
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.
Definition adapt.hh:17
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.