version 3.8
discretization/pq1bubble/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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH
14#define DUMUX_DISCRETIZATION_PQ1BUBBLE_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>
28
29namespace Dumux {
30
32template <class ct>
33struct PQ1BubbleMLGeometryTraits : public Dune::MultiLinearGeometryTraits<ct>
34{
35 // we use static vectors to store the corners as we know
36 // the maximum number of corners in advance (2^dim)
37 template< int mydim, int cdim >
39 {
40 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)+1>;
41 };
42};
43
45
46template<Dune::GeometryType::Id gt>
48
49template<>
51{
52 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
53 static constexpr std::array<std::array<Key, 2>, 1> keys = {{
54 { Key{0, 1}, Key{1, 1} }
55 }};
56};
57
58template<>
59struct OverlappingScvCorners<Dune::GeometryTypes::triangle>
60{
61 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
62 static constexpr std::array<std::array<Key, 3>, 1> keys = {{
63 { Key{0, 1}, Key{1, 1}, Key{2, 1} }
64 }};
65};
66
67template<>
68struct OverlappingScvCorners<Dune::GeometryTypes::quadrilateral>
69{
70 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
71 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
72 { Key{2, 1}, Key{1, 1}, Key{0, 1}, Key{3, 1} }
73 }};
74};
75
76template<>
77struct OverlappingScvCorners<Dune::GeometryTypes::tetrahedron>
78{
79 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
80 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
81 { Key{0, 1}, Key{1, 1}, Key{2, 1}, Key{3, 1} }
82 }};
83};
84
85template<>
86struct OverlappingScvCorners<Dune::GeometryTypes::hexahedron>
87{
88 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
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} }
91 }};
92};
93
94
95template<Dune::GeometryType::Id gt>
97
98template<>
100{
101 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
102 static constexpr std::array<std::array<Key, 1>, 1> keys = {{
103 { Key{0, 0} }
104 }};
105};
106
107template<>
108struct OverlappingScvfCorners<Dune::GeometryTypes::triangle>
109{
110 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
111 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
112 { Key{0, 1}, Key{1, 1} },
113 { Key{0, 1}, Key{2, 1} },
114 { Key{1, 1}, Key{2, 1} }
115 }};
116};
117
118template<>
119struct OverlappingScvfCorners<Dune::GeometryTypes::quadrilateral>
120{
121 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
122 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
123 { Key{0, 1}, Key{2, 1} },
124 { Key{2, 1}, Key{1, 1} },
125 { Key{0, 1}, Key{3, 1} },
126 { Key{1, 1}, Key{3, 1} }
127 }};
128};
129
130template<>
131struct OverlappingScvfCorners<Dune::GeometryTypes::tetrahedron>
132{
133 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
134 static constexpr std::array<std::array<Key, 3>, 10> keys = {{
135 { Key{0, 1}, Key{1, 1}, Key{2, 1} },
136 { Key{0, 1}, Key{1, 1}, Key{3, 1} },
137 { Key{0, 1}, Key{2, 1}, Key{3, 1} },
138 { Key{1, 1}, Key{2, 1}, Key{3, 1} }
139 }};
140};
141
142template<>
143struct OverlappingScvfCorners<Dune::GeometryTypes::hexahedron>
144{
145 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
146 static constexpr std::array<std::array<Key, 3>, 8> keys = {{
147 { Key{4, 1}, Key{0, 1}, Key{2, 1} },
148 { Key{4, 1}, Key{2, 1}, Key{1, 1} },
149 { Key{4, 1}, Key{0, 1}, Key{3, 1} },
150 { Key{4, 1}, Key{1, 1}, Key{3, 1} },
151 { Key{5, 1}, Key{0, 1}, Key{2, 1} },
152 { Key{5, 1}, Key{2, 1}, Key{1, 1} },
153 { Key{5, 1}, Key{0, 1}, Key{3, 1} },
154 { Key{5, 1}, Key{1, 1}, Key{3, 1} }
155 }};
156};
157
158} // end namespace Detail::PQ1Bubble
159
164template <class GridView, class ScvType, class ScvfType>
166{
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;
172
173 using Element = typename GridView::template Codim<0>::Entity;
174 using Intersection = typename GridView::Intersection;
175
176 static constexpr auto dim = GridView::dimension;
177 static constexpr auto dimWorld = GridView::dimensionworld;
178public:
179
180 PQ1BubbleGeometryHelper(const typename Element::Geometry& geometry)
181 : geo_(geometry)
182 , boxHelper_(geometry)
183 {}
184
186 ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
187 {
188 // proceed according to number of corners of the element
189 const auto type = geo_.type();
190 const auto numBoxScv = boxHelper_.numScv();
191 // reuse box geometry helper for the corner scvs
192 if (localScvIdx < numBoxScv)
193 return boxHelper_.getScvCorners(localScvIdx);
194
195 const auto localOverlappingScvIdx = localScvIdx-numBoxScv;
196 if (type == Dune::GeometryTypes::triangle)
197 {
199 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
200 }
201 else if (type == Dune::GeometryTypes::quadrilateral)
202 {
204 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
205 }
206 else if (type == Dune::GeometryTypes::tetrahedron)
207 {
209 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
210 }
211 else if (type == Dune::GeometryTypes::hexahedron)
212 {
214 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
215 }
216 else
217 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv geometries for dim=" << dim
218 << " dimWorld=" << dimWorld
219 << " type=" << type);
220 }
221
222 Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
223 {
224 // proceed according to number of corners of the element
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); // octahedron
235 else
236 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv geometries for dim=" << dim
237 << " dimWorld=" << dimWorld
238 << " type=" << type);
239 }
240
242 ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
243 {
244 // proceed according to number of corners
245 const auto type = geo_.type();
246 const auto numBoxScvf = boxHelper_.numInteriorScvf();
247 // reuse box geometry helper for the corner scvs
248 if (localScvfIdx < numBoxScvf)
249 return boxHelper_.getScvfCorners(localScvfIdx);
250
251 const auto localOverlappingScvfIdx = localScvfIdx-numBoxScvf;
252 if (type == Dune::GeometryTypes::triangle)
253 {
255 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
256 }
257 else if (type == Dune::GeometryTypes::quadrilateral)
258 {
260 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
261 }
262 else if (type == Dune::GeometryTypes::tetrahedron)
263 {
265 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
266 }
267 else if (type == Dune::GeometryTypes::hexahedron)
268 {
270 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
271 }
272 else
273 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scvf geometries for dim=" << dim
274 << " dimWorld=" << dimWorld
275 << " type=" << type);
276 }
277
278 Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
279 {
280 const auto numBoxScvf = boxHelper_.numInteriorScvf();
281 if (localScvfIdx < numBoxScvf)
282 return Dune::GeometryTypes::cube(dim-1);
283 else
284 return Dune::GeometryTypes::simplex(dim-1);
285 }
286
288 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex,
289 unsigned int indexInFacet) const
290 {
291 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
292 }
293
294 template<int d = dimWorld, std::enable_if_t<(d==3), int> = 0>
295 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
296 {
297 auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]);
298 normal /= normal.two_norm();
299
300 GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]);
301
302 const auto s = v*normal;
303 if (std::signbit(s))
304 normal *= -1;
305
306 return normal;
307 }
308
309 template<int d = dimWorld, std::enable_if_t<(d==2), int> = 0>
310 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
311 {
313 const auto t = p[1] - p[0];
314 GlobalPosition normal({-t[1], t[0]});
315 normal /= normal.two_norm();
316
317 GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]);
318
319 const auto s = v*normal;
320 if (std::signbit(s))
321 normal *= -1;
322
323 return normal;
324 }
325
327 const typename Element::Geometry& elementGeometry() const
328 { return geo_; }
329
331 std::size_t numInteriorScvf() const
332 {
333 return boxHelper_.numInteriorScvf() + referenceElement(geo_).size(dim);
334 }
335
337 std::size_t numBoundaryScvf(unsigned int localFacetIndex) const
338 {
339 return referenceElement(geo_).size(localFacetIndex, 1, dim);
340 }
341
343 std::size_t numScv() const
344 {
345 return boxHelper_.numScv() + 1;
346 }
347
349 Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage& p) const
350 {
351 const auto scvType = getScvGeometryType(localScvIdx);
352 if constexpr (dim == 3)
353 if (scvType == Dune::GeometryTypes::none(dim))
354 return octahedronVolume_(p);
355
356 return Dumux::convexPolytopeVolume<dim>(
357 scvType,
358 [&](unsigned int i){ return p[i]; }
359 );
360 }
361
362 template<class DofMapper>
363 auto dofIndex(const DofMapper& dofMapper, const Element& element, unsigned int localScvIdx) const
364 {
365 if (localScvIdx < numScv()-1)
366 return dofMapper.subIndex(element, localScvIdx, dim);
367 else
368 return dofMapper.index(element);
369 }
370
371 GlobalPosition dofPosition(unsigned int localScvIdx) const
372 {
373 if (localScvIdx < numScv()-1)
374 return geo_.corner(localScvIdx);
375 else
376 return geo_.center();
377 }
378
379 std::array<LocalIndexType, 2> getScvPairForScvf(unsigned int localScvfIndex) const
380 {
381 const auto numEdges = referenceElement(geo_).size(dim-1);
382 if (localScvfIndex < numEdges)
383 return {
384 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
385 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
386 };
387 else
388 return {
389 static_cast<LocalIndexType>(numScv()-1),
390 static_cast<LocalIndexType>(localScvfIndex-numEdges)
391 };
392 }
393
394 std::array<LocalIndexType, 2> getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
395 {
396 const LocalIndexType insideScvIdx
397 = static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
398 return { insideScvIdx, insideScvIdx };
399 }
400
401 bool isOverlappingScvf(unsigned int localScvfIndex) const
402 {
403 if (localScvfIndex < boxHelper_.numInteriorScvf())
404 return false;
405 else
406 return true;
407 }
408
409 bool isOverlappingScv(unsigned int localScvIndex) const
410 {
411 if (localScvIndex < boxHelper_.numScv())
412 return false;
413 else
414 return true;
415 }
416
417private:
418 Scalar octahedronVolume_(const ScvCornerStorage& p) const
419 {
420 using std::abs;
421 return 1.0/6.0 * (
422 abs(Dumux::tripleProduct(p[4]-p[0], p[1]-p[0], p[2]-p[0]))
423 + abs(Dumux::tripleProduct(p[5]-p[0], p[1]-p[0], p[2]-p[0]))
424 );
425 }
426
427 const typename Element::Geometry& geo_;
429};
430
431} // end namespace Dumux
432
433#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
A class to create sub control volume and sub control volume face geometries per element.
Definition: discretization/pq1bubble/geometryhelper.hh:166
PQ1BubbleGeometryHelper(const typename Element::Geometry &geometry)
Definition: discretization/pq1bubble/geometryhelper.hh:180
bool isOverlappingScvf(unsigned int localScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:401
auto dofIndex(const DofMapper &dofMapper, const Element &element, unsigned int localScvIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:363
std::size_t numScv() const
number of sub control volumes (number of codim-1 entities)
Definition: discretization/pq1bubble/geometryhelper.hh:343
std::size_t numBoundaryScvf(unsigned int localFacetIndex) const
number of boundary sub control volume faces for face localFacetIndex
Definition: discretization/pq1bubble/geometryhelper.hh:337
std::size_t numInteriorScvf() const
number of interior sub control volume faces
Definition: discretization/pq1bubble/geometryhelper.hh:331
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:371
std::array< LocalIndexType, 2 > getScvPairForScvf(unsigned int localScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:379
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:278
bool isOverlappingScv(unsigned int localScvIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:409
Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage &p) const
get scv volume
Definition: discretization/pq1bubble/geometryhelper.hh:349
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition: discretization/pq1bubble/geometryhelper.hh:327
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:394
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition: discretization/pq1bubble/geometryhelper.hh:295
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:642
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:671
Define some often used mathematical functions.
constexpr None none
Definition: method.hh:153
CVFE< CVFEMethods::PQ1Bubble > PQ1Bubble
Definition: method.hh:108
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: adapt.hh:17
Definition: common/pdesolver.hh:24
Definition: discretization/pq1bubble/geometryhelper.hh:87
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
Definition: discretization/pq1bubble/geometryhelper.hh:69
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:70
Definition: discretization/pq1bubble/geometryhelper.hh:78
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:79
Definition: discretization/pq1bubble/geometryhelper.hh:60
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:61
Definition: discretization/pq1bubble/geometryhelper.hh:47
Definition: discretization/pq1bubble/geometryhelper.hh:144
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:145
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:101
Definition: discretization/pq1bubble/geometryhelper.hh:120
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:121
Definition: discretization/pq1bubble/geometryhelper.hh:132
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:133
Definition: discretization/pq1bubble/geometryhelper.hh:109
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.