3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH
26#define DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH
27
28#include <array>
29
30#include <dune/common/exceptions.hh>
31
32#include <dune/geometry/type.hh>
33#include <dune/geometry/referenceelements.hh>
34#include <dune/geometry/multilineargeometry.hh>
35#include <dune/common/reservedvector.hh>
36
37#include <dumux/common/math.hh>
40
41namespace Dumux {
42
44template <class ct>
45struct PQ1BubbleMLGeometryTraits : public Dune::MultiLinearGeometryTraits<ct>
46{
47 // we use static vectors to store the corners as we know
48 // the maximum number of corners in advance (2^dim)
49 template< int mydim, int cdim >
51 {
52 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)+1>;
53 };
54};
55
56namespace Detail::PQ1Bubble {
57
58template<Dune::GeometryType::Id gt>
60
61template<>
63{
64 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
65 static constexpr std::array<std::array<Key, 2>, 1> keys = {{
66 { Key{0, 1}, Key{1, 1} }
67 }};
68};
69
70template<>
71struct OverlappingScvCorners<Dune::GeometryTypes::triangle>
72{
73 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
74 static constexpr std::array<std::array<Key, 3>, 1> keys = {{
75 { Key{0, 1}, Key{1, 1}, Key{2, 1} }
76 }};
77};
78
79template<>
80struct OverlappingScvCorners<Dune::GeometryTypes::quadrilateral>
81{
82 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
83 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
84 { Key{2, 1}, Key{1, 1}, Key{0, 1}, Key{3, 1} }
85 }};
86};
87
88template<>
89struct OverlappingScvCorners<Dune::GeometryTypes::tetrahedron>
90{
91 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
92 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
93 { Key{0, 1}, Key{1, 1}, Key{2, 1}, Key{3, 1} }
94 }};
95};
96
97template<>
98struct OverlappingScvCorners<Dune::GeometryTypes::hexahedron>
99{
100 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
101 static constexpr std::array<std::array<Key, 6>, 1> keys = {{
102 { Key{0, 1}, Key{2, 1}, Key{3, 1}, Key{1, 1}, Key{4, 1}, Key{5, 1} }
103 }};
104};
105
106
107template<Dune::GeometryType::Id gt>
109
110template<>
112{
113 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
114 static constexpr std::array<std::array<Key, 1>, 1> keys = {{
115 { Key{0, 0} }
116 }};
117};
118
119template<>
120struct OverlappingScvfCorners<Dune::GeometryTypes::triangle>
121{
122 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
123 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
124 { Key{0, 1}, Key{1, 1} },
125 { Key{0, 1}, Key{2, 1} },
126 { Key{1, 1}, Key{2, 1} }
127 }};
128};
129
130template<>
131struct OverlappingScvfCorners<Dune::GeometryTypes::quadrilateral>
132{
133 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
134 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
135 { Key{0, 1}, Key{2, 1} },
136 { Key{2, 1}, Key{1, 1} },
137 { Key{0, 1}, Key{3, 1} },
138 { Key{1, 1}, Key{3, 1} }
139 }};
140};
141
142template<>
143struct OverlappingScvfCorners<Dune::GeometryTypes::tetrahedron>
144{
145 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
146 static constexpr std::array<std::array<Key, 3>, 10> keys = {{
147 { Key{0, 1}, Key{1, 1}, Key{2, 1} },
148 { Key{0, 1}, Key{1, 1}, Key{3, 1} },
149 { Key{0, 1}, Key{2, 1}, Key{3, 1} },
150 { Key{1, 1}, Key{2, 1}, Key{3, 1} }
151 }};
152};
153
154template<>
155struct OverlappingScvfCorners<Dune::GeometryTypes::hexahedron>
156{
157 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
158 static constexpr std::array<std::array<Key, 3>, 8> keys = {{
159 { Key{4, 1}, Key{0, 1}, Key{2, 1} },
160 { Key{4, 1}, Key{2, 1}, Key{1, 1} },
161 { Key{4, 1}, Key{0, 1}, Key{3, 1} },
162 { Key{4, 1}, Key{1, 1}, Key{3, 1} },
163 { Key{5, 1}, Key{0, 1}, Key{2, 1} },
164 { Key{5, 1}, Key{2, 1}, Key{1, 1} },
165 { Key{5, 1}, Key{0, 1}, Key{3, 1} },
166 { Key{5, 1}, Key{1, 1}, Key{3, 1} }
167 }};
168};
169
170} // end namespace Detail::PQ1Bubble
171
176template <class GridView, class ScvType, class ScvfType>
178{
179 using Scalar = typename GridView::ctype;
180 using GlobalPosition = typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
181 using ScvCornerStorage = typename ScvType::Traits::CornerStorage;
182 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
183 using LocalIndexType = typename ScvType::Traits::LocalIndexType;
184
185 using Element = typename GridView::template Codim<0>::Entity;
186 using Intersection = typename GridView::Intersection;
187
188 static constexpr auto dim = GridView::dimension;
189 static constexpr auto dimWorld = GridView::dimensionworld;
190public:
191
192 PQ1BubbleGeometryHelper(const typename Element::Geometry& geometry)
193 : geo_(geometry)
194 , boxHelper_(geometry)
195 {}
196
198 ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
199 {
200 // proceed according to number of corners of the element
201 const auto type = geo_.type();
202 const auto numBoxScv = boxHelper_.numScv();
203 // reuse box geometry helper for the corner scvs
204 if (localScvIdx < numBoxScv)
205 return boxHelper_.getScvCorners(localScvIdx);
206
207 const auto localOverlappingScvIdx = localScvIdx-numBoxScv;
208 if (type == Dune::GeometryTypes::triangle)
209 {
211 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
212 }
213 else if (type == Dune::GeometryTypes::quadrilateral)
214 {
216 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
217 }
218 else if (type == Dune::GeometryTypes::tetrahedron)
219 {
221 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
222 }
223 else if (type == Dune::GeometryTypes::hexahedron)
224 {
226 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
227 }
228 else
229 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv geometries for dim=" << dim
230 << " dimWorld=" << dimWorld
231 << " type=" << type);
232 }
233
234 Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
235 {
236 // proceed according to number of corners of the element
237 const auto type = geo_.type();
238 const auto numBoxScv = boxHelper_.numScv();
239 if (localScvIdx < numBoxScv)
240 return Dune::GeometryTypes::cube(dim);
241 else if (type == Dune::GeometryTypes::simplex(dim))
242 return Dune::GeometryTypes::simplex(dim);
243 else if (type == Dune::GeometryTypes::quadrilateral)
244 return Dune::GeometryTypes::quadrilateral;
245 else if (type == Dune::GeometryTypes::hexahedron)
246 return Dune::GeometryTypes::none(dim); // octahedron
247 else
248 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv geometries for dim=" << dim
249 << " dimWorld=" << dimWorld
250 << " type=" << type);
251 }
252
254 ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
255 {
256 // proceed according to number of corners
257 const auto type = geo_.type();
258 const auto numBoxScvf = boxHelper_.numInteriorScvf();
259 // reuse box geometry helper for the corner scvs
260 if (localScvfIdx < numBoxScvf)
261 return boxHelper_.getScvfCorners(localScvfIdx);
262
263 const auto localOverlappingScvfIdx = localScvfIdx-numBoxScvf;
264 if (type == Dune::GeometryTypes::triangle)
265 {
267 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
268 }
269 else if (type == Dune::GeometryTypes::quadrilateral)
270 {
272 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
273 }
274 else if (type == Dune::GeometryTypes::tetrahedron)
275 {
277 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
278 }
279 else if (type == Dune::GeometryTypes::hexahedron)
280 {
282 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
283 }
284 else
285 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scvf geometries for dim=" << dim
286 << " dimWorld=" << dimWorld
287 << " type=" << type);
288 }
289
290 Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
291 {
292 const auto numBoxScvf = boxHelper_.numInteriorScvf();
293 if (localScvfIdx < numBoxScvf)
294 return Dune::GeometryTypes::cube(dim-1);
295 else
296 return Dune::GeometryTypes::simplex(dim-1);
297 }
298
300 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex,
301 unsigned int indexInFacet) const
302 {
303 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
304 }
305
306 template<int d = dimWorld, std::enable_if_t<(d==3), int> = 0>
307 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
308 {
309 auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]);
310 normal /= normal.two_norm();
311
312 GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]);
313
314 const auto s = v*normal;
315 if (std::signbit(s))
316 normal *= -1;
317
318 return normal;
319 }
320
321 template<int d = dimWorld, std::enable_if_t<(d==2), int> = 0>
322 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
323 {
325 const auto t = p[1] - p[0];
326 GlobalPosition normal({-t[1], t[0]});
327 normal /= normal.two_norm();
328
329 GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]);
330
331 const auto s = v*normal;
332 if (std::signbit(s))
333 normal *= -1;
334
335 return normal;
336 }
337
339 const typename Element::Geometry& elementGeometry() const
340 { return geo_; }
341
343 std::size_t numInteriorScvf() const
344 {
345 return boxHelper_.numInteriorScvf() + referenceElement(geo_).size(dim);
346 }
347
349 std::size_t numBoundaryScvf(unsigned int localFacetIndex) const
350 {
351 return referenceElement(geo_).size(localFacetIndex, 1, dim);
352 }
353
355 std::size_t numScv() const
356 {
357 return boxHelper_.numScv() + 1;
358 }
359
361 Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage& p) const
362 {
363 const auto scvType = getScvGeometryType(localScvIdx);
364 if constexpr (dim == 3)
365 if (scvType == Dune::GeometryTypes::none(dim))
366 return octahedronVolume_(p);
367
368 return Dumux::convexPolytopeVolume<dim>(
369 scvType,
370 [&](unsigned int i){ return p[i]; }
371 );
372 }
373
374 template<class DofMapper>
375 auto dofIndex(const DofMapper& dofMapper, const Element& element, unsigned int localScvIdx) const
376 {
377 if (localScvIdx < numScv()-1)
378 return dofMapper.subIndex(element, localScvIdx, dim);
379 else
380 return dofMapper.index(element);
381 }
382
383 GlobalPosition dofPosition(unsigned int localScvIdx) const
384 {
385 if (localScvIdx < numScv()-1)
386 return geo_.corner(localScvIdx);
387 else
388 return geo_.center();
389 }
390
391 std::array<LocalIndexType, 2> getScvPairForScvf(unsigned int localScvfIndex) const
392 {
393 const auto numEdges = referenceElement(geo_).size(dim-1);
394 if (localScvfIndex < numEdges)
395 return {
396 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
397 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
398 };
399 else
400 return {
401 static_cast<LocalIndexType>(numScv()-1),
402 static_cast<LocalIndexType>(localScvfIndex-numEdges)
403 };
404 }
405
406 std::array<LocalIndexType, 2> getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
407 {
408 const LocalIndexType insideScvIdx
409 = static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
410 return { insideScvIdx, insideScvIdx };
411 }
412
413 bool isOverlappingScvf(unsigned int localScvfIndex) const
414 {
415 if (localScvfIndex < boxHelper_.numInteriorScvf())
416 return false;
417 else
418 return true;
419 }
420
421 bool isOverlappingScv(unsigned int localScvIndex) const
422 {
423 if (localScvIndex < boxHelper_.numScv())
424 return false;
425 else
426 return true;
427 }
428
429private:
430 Scalar octahedronVolume_(const ScvCornerStorage& p) const
431 {
432 using std::abs;
433 return 1.0/6.0 * (
434 abs(Dumux::tripleProduct(p[4]-p[0], p[1]-p[0], p[2]-p[0]))
435 + abs(Dumux::tripleProduct(p[5]-p[0], p[1]-p[0], p[2]-p[0]))
436 );
437 }
438
439 const typename Element::Geometry& geo_;
441};
442
443} // end namespace Dumux
444
445#endif
Define some often used mathematical functions.
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Compute the volume of several common geometry types.
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:654
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:683
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Definition: deprecated.hh:149
constexpr None none
Definition: method.hh:142
CVFE< CVFEMethods::PQ1Bubble > PQ1Bubble
Definition: method.hh:97
constexpr Line line
Definition: couplingmanager1d3d_line.hh:43
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:269
Traits for an efficient corner storage for the PQ1Bubble method.
Definition: discretization/pq1bubble/geometryhelper.hh:46
Definition: discretization/pq1bubble/geometryhelper.hh:51
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<< mydim)+1 > Type
Definition: discretization/pq1bubble/geometryhelper.hh:52
Definition: discretization/pq1bubble/geometryhelper.hh:59
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:64
Definition: discretization/pq1bubble/geometryhelper.hh:72
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:73
Definition: discretization/pq1bubble/geometryhelper.hh:81
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:82
Definition: discretization/pq1bubble/geometryhelper.hh:90
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:91
Definition: discretization/pq1bubble/geometryhelper.hh:99
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:100
Definition: discretization/pq1bubble/geometryhelper.hh:108
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:113
Definition: discretization/pq1bubble/geometryhelper.hh:121
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:122
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:144
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:145
Definition: discretization/pq1bubble/geometryhelper.hh:156
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:157
A class to create sub control volume and sub control volume face geometries per element.
Definition: discretization/pq1bubble/geometryhelper.hh:178
PQ1BubbleGeometryHelper(const typename Element::Geometry &geometry)
Definition: discretization/pq1bubble/geometryhelper.hh:192
bool isOverlappingScvf(unsigned int localScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:413
auto dofIndex(const DofMapper &dofMapper, const Element &element, unsigned int localScvIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:375
std::size_t numScv() const
number of sub control volumes (number of codim-1 entities)
Definition: discretization/pq1bubble/geometryhelper.hh:355
std::size_t numBoundaryScvf(unsigned int localFacetIndex) const
number of boundary sub control volume faces for face localFacetIndex
Definition: discretization/pq1bubble/geometryhelper.hh:349
std::size_t numInteriorScvf() const
number of interior sub control volume faces
Definition: discretization/pq1bubble/geometryhelper.hh:343
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition: discretization/pq1bubble/geometryhelper.hh:300
GlobalPosition dofPosition(unsigned int localScvIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:383
std::array< LocalIndexType, 2 > getScvPairForScvf(unsigned int localScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:391
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:290
bool isOverlappingScv(unsigned int localScvIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:421
Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage &p) const
get scv volume
Definition: discretization/pq1bubble/geometryhelper.hh:361
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition: discretization/pq1bubble/geometryhelper.hh:339
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition: discretization/pq1bubble/geometryhelper.hh:254
std::array< LocalIndexType, 2 > getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:406
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition: discretization/pq1bubble/geometryhelper.hh:307
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition: discretization/pq1bubble/geometryhelper.hh:198
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:234