version 3.11-dev
Loading...
Searching...
No Matches
discretization/pq3/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//
7// This file includes code adapted from dune-functions.
8// These parts are clearly marked with inline comments and were originally
9// licensed under LGPL-3.0-or-later. They are adapted and relicensed here
10// under the terms of the GNU GPL-3.0-or-later as permitted by LGPLv3 Sec. 3.
11//
12// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file AUTHORS.md
13//
20#ifndef DUMUX_DISCRETIZATION_PQ3_GEOMETRY_HELPER_HH
21#define DUMUX_DISCRETIZATION_PQ3_GEOMETRY_HELPER_HH
22
23#include <array>
24
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>
30
31#include <dumux/common/math.hh>
34// Reuse PQ2 corner storage traits (same structure, different order)
37
38namespace Dumux {
39
51template <class GridView, class ScvType, class ScvfType>
53{
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;
59
60 using Element = typename GridView::template Codim<0>::Entity;
61 using Intersection = typename GridView::Intersection;
62
63 static constexpr auto dim = GridView::dimension;
64 static constexpr auto dimWorld = GridView::dimensionworld;
65
67public:
68
69 HybridPQ3GeometryHelper(const typename Element::Geometry& geometry)
70 : geo_(geometry)
71 , boxHelper_(geometry)
72 {}
73
75 ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
76 {
77 return getScvCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localScvIdx);
78 }
79
81 template<class Transformation>
82 static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localScvIdx)
83 {
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);
88
89 DUNE_THROW(Dune::NotImplemented, "PQ3 scv corners call for hybrid dofs");
90 }
91
92 Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
93 {
94 const auto numBoxScv = boxHelper_.numScv();
95 if (localScvIdx < numBoxScv)
96 return Dune::GeometryTypes::cube(dim);
97
98 DUNE_THROW(Dune::NotImplemented, "PQ3 scv geometry call for hybrid dofs");
99 }
100
102 ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
103 {
104 return getScvfCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localScvfIdx);
105 }
106
108 template<class Transformation>
109 static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localScvfIdx)
110 {
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);
115
116 DUNE_THROW(Dune::NotImplemented, "PQ3 scvf corners call for hybrid dofs");
117 }
118
119 Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
120 {
121 const auto numBoxScvf = boxHelper_.numInteriorScvf();
122 if (localScvfIdx < numBoxScvf)
123 return Dune::GeometryTypes::cube(dim-1);
124
125 DUNE_THROW(Dune::NotImplemented, "PQ3 interior scvf geometry type call for hybrid dofs");
126 }
127
129 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
130 {
131 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
132 }
133
134 Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
135 {
136 return Dune::GeometryTypes::cube(dim-1);
137 }
138
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)
141 {
142 auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]);
143 normal /= normal.two_norm();
144
145 GlobalPosition v = geo_.corner(scvPair[1]) - geo_.corner(scvPair[0]);
146
147 const auto s = v*normal;
148 if (std::signbit(s))
149 normal *= -1;
150
151 return normal;
152 }
153
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)
156 {
157 const auto t = p[1] - p[0];
158 GlobalPosition normal({-t[1], t[0]});
159 normal /= normal.two_norm();
160
161 GlobalPosition v = geo_.corner(scvPair[1]) - geo_.corner(scvPair[0]);
162
163 const auto s = v*normal;
164 if (std::signbit(s))
165 normal *= -1;
166
167 return normal;
168 }
169
171 const typename Element::Geometry& elementGeometry() const
172 { return geo_; }
173
175 static auto numInteriorScvf(Dune::GeometryType type)
176 {
177 return BoxHelper::numInteriorScvf(type);
178 }
179
181 static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
182 {
183 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
184 }
185
187 std::size_t numScv() const
188 {
189 return boxHelper_.numScv();
190 }
191
193 Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage& p) const
194 {
195 const auto scvType = getScvGeometryType(localScvIdx);
197 scvType,
198 [&](unsigned int i){ return p[i]; }
199 );
200 }
201
203 template<class LocalKey>
204 static auto localDofOnIntersection(Dune::GeometryType type, unsigned int iIdx, const LocalKey& localKey)
205 {
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()))
210 return true;
211 return false;
212 }
213
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)
218 { return PQ3LagrangeDofHelper<GridView>::dofIndex(dofMapper, element, localKey, gidSet); }
219
220 template<class Geometry, class LocalKey>
221 static GlobalPosition dofPosition(const Geometry& geo, const LocalKey& localKey)
222 { return PQ3LagrangeDofHelper<GridView>::dofPosition(geo, localKey); }
223
224 template<class LocalKey>
225 GlobalPosition dofPosition(const LocalKey& localKey) const
226 { return dofPosition(geo_, localKey); }
227
229 template<class LocalKey>
230 static typename Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey& localKey)
231 { return PQ3LagrangeDofHelper<GridView>::localDofPos(type, localKey); }
232
233 std::array<LocalIndexType, 2> getScvPairForScvf(unsigned int localScvfIndex) const
234 {
235 const auto numBoxFaces = boxHelper_.numInteriorScvf();
236 if (localScvfIndex < numBoxFaces)
237 {
238 return {
239 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
240 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
241 };
242 }
243
244 DUNE_THROW(Dune::NotImplemented, "PQ3 scv pair call for hybrid dofs");
245 }
246
247 std::array<LocalIndexType, 2> getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
248 {
249 const auto numBoxScvf = referenceElement(geo_).size(localFacetIndex, 1, dim);
250 if (localIsScvfIndex < numBoxScvf)
251 {
252 const LocalIndexType insideScvIdx
253 = static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
254 return { insideScvIdx, insideScvIdx };
255 }
256
257 DUNE_THROW(Dune::NotImplemented, "PQ3 scv boundary pair call for hybrid dofs");
258 }
259
260 bool isOverlappingScvf(unsigned int localScvfIndex) const
261 { return false; }
262
263 bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
264 { return false; }
265
266 bool isOverlappingScv(unsigned int localScvIndex) const
267 { return false; }
268
269private:
270 const typename Element::Geometry& geo_;
271 BoxHelper boxHelper_;
272};
273
274} // end namespace Dumux
275
276#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/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.
Definition adapt.hh:17
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.