version 3.11-dev
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-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_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>
29
30namespace Dumux {
31
33template <class ct>
34struct PQ1BubbleMLGeometryTraits : 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)
38 template< int mydim, int cdim >
40 {
41 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)+1>;
42 };
43};
44
46
47template<Dune::GeometryType::Id gt>
49
50template<>
52{
53 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
54 static constexpr std::array<std::array<Key, 2>, 1> keys = {{
55 { Key{0, 1}, Key{1, 1} }
56 }};
57};
58
59template<>
60struct OverlappingScvCorners<Dune::GeometryTypes::triangle>
61{
62 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
63 static constexpr std::array<std::array<Key, 3>, 1> keys = {{
64 { Key{0, 1}, Key{1, 1}, Key{2, 1} }
65 }};
66};
67
68template<>
69struct OverlappingScvCorners<Dune::GeometryTypes::quadrilateral>
70{
71 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
72 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
73 { Key{2, 1}, Key{1, 1}, Key{0, 1}, Key{3, 1} }
74 }};
75};
76
77template<>
78struct OverlappingScvCorners<Dune::GeometryTypes::tetrahedron>
79{
80 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
81 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
82 { Key{0, 1}, Key{1, 1}, Key{2, 1}, Key{3, 1} }
83 }};
84};
85
86template<>
87struct OverlappingScvCorners<Dune::GeometryTypes::hexahedron>
88{
89 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
90 static constexpr std::array<std::array<Key, 6>, 1> keys = {{
91 { Key{0, 1}, Key{2, 1}, Key{3, 1}, Key{1, 1}, Key{4, 1}, Key{5, 1} }
92 }};
93};
94
95
96template<Dune::GeometryType::Id gt>
98
99template<>
101{
102 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
103 static constexpr std::array<std::array<Key, 1>, 1> keys = {{
104 { Key{0, 0} }
105 }};
106};
107
108template<>
109struct OverlappingScvfCorners<Dune::GeometryTypes::triangle>
110{
111 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
112 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
113 { Key{0, 1}, Key{1, 1} },
114 { Key{0, 1}, Key{2, 1} },
115 { Key{1, 1}, Key{2, 1} }
116 }};
117};
118
119template<>
120struct OverlappingScvfCorners<Dune::GeometryTypes::quadrilateral>
121{
122 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
123 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
124 { Key{0, 1}, Key{2, 1} },
125 { Key{2, 1}, Key{1, 1} },
126 { Key{0, 1}, Key{3, 1} },
127 { Key{1, 1}, Key{3, 1} }
128 }};
129};
130
131template<>
132struct OverlappingScvfCorners<Dune::GeometryTypes::tetrahedron>
133{
134 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
135 static constexpr std::array<std::array<Key, 3>, 10> keys = {{
136 { Key{0, 1}, Key{1, 1}, Key{2, 1} },
137 { Key{0, 1}, Key{1, 1}, Key{3, 1} },
138 { Key{0, 1}, Key{2, 1}, Key{3, 1} },
139 { Key{1, 1}, Key{2, 1}, Key{3, 1} }
140 }};
141};
142
143template<>
144struct OverlappingScvfCorners<Dune::GeometryTypes::hexahedron>
145{
146 using Key = std::pair<std::uint8_t, std::uint8_t>; // (i, codim)
147 static constexpr std::array<std::array<Key, 3>, 8> keys = {{
148 { Key{4, 1}, Key{0, 1}, Key{2, 1} },
149 { Key{4, 1}, Key{2, 1}, Key{1, 1} },
150 { Key{4, 1}, Key{0, 1}, Key{3, 1} },
151 { Key{4, 1}, Key{1, 1}, Key{3, 1} },
152 { Key{5, 1}, Key{0, 1}, Key{2, 1} },
153 { Key{5, 1}, Key{2, 1}, Key{1, 1} },
154 { Key{5, 1}, Key{0, 1}, Key{3, 1} },
155 { Key{5, 1}, Key{1, 1}, Key{3, 1} }
156 }};
157};
158
159} // end namespace Detail::PQ1Bubble
160
165template <class GridView, class ScvType, class ScvfType>
167{
168 using Scalar = typename GridView::ctype;
169 using GlobalPosition = typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
170 using ScvCornerStorage = typename ScvType::Traits::CornerStorage;
171 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
172 using LocalIndexType = typename ScvType::Traits::LocalIndexType;
173
174 using Element = typename GridView::template Codim<0>::Entity;
175 using Intersection = typename GridView::Intersection;
176
177 static constexpr auto dim = GridView::dimension;
178 static constexpr auto dimWorld = GridView::dimensionworld;
179
181public:
182
183 PQ1BubbleGeometryHelper(const typename Element::Geometry& geometry)
184 : geo_(geometry)
185 , boxHelper_(geometry)
186 {}
187
189 ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
190 {
191 return getScvCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localScvIdx);
192 }
193
195 template<class Transformation>
196 static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localScvIdx)
197 {
198 // proceed according to number of corners of the element
199 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
200 const auto numBoxScv = ref.size(dim);
201 // reuse box geometry helper for the corner scvs
202 if (localScvIdx < numBoxScv)
203 return BoxHelper::getScvCorners(type, trans, localScvIdx);
204
205 const auto localOverlappingScvIdx = localScvIdx-numBoxScv;
206 if (type == Dune::GeometryTypes::triangle)
207 {
209 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localOverlappingScvIdx]);
210 }
211 else if (type == Dune::GeometryTypes::quadrilateral)
212 {
214 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localOverlappingScvIdx]);
215 }
216 else if (type == Dune::GeometryTypes::tetrahedron)
217 {
219 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localOverlappingScvIdx]);
220 }
221 else if (type == Dune::GeometryTypes::hexahedron)
222 {
224 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localOverlappingScvIdx]);
225 }
226 else
227 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv geometries for dim=" << dim
228 << " dimWorld=" << dimWorld
229 << " type=" << type);
230 }
231
232 Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
233 {
234 // proceed according to number of corners of the element
235 const auto type = geo_.type();
236 const auto numBoxScv = boxHelper_.numScv();
237 if (localScvIdx < numBoxScv)
238 return Dune::GeometryTypes::cube(dim);
239 else if (type == Dune::GeometryTypes::simplex(dim))
240 return Dune::GeometryTypes::simplex(dim);
241 else if (type == Dune::GeometryTypes::quadrilateral)
242 return Dune::GeometryTypes::quadrilateral;
243 else if (type == Dune::GeometryTypes::hexahedron)
244 return Dune::GeometryTypes::none(dim); // octahedron
245 else
246 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv geometries for dim=" << dim
247 << " dimWorld=" << dimWorld
248 << " type=" << type);
249 }
250
252 ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
253 {
254 return getScvfCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localScvfIdx);
255 }
256
258 template<class Transformation>
259 static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localScvfIdx)
260 {
261 // proceed according to number of corners
262 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
263 const auto numBoxScvf = ref.size(dim-1);
264 // reuse box geometry helper for the corner scvs
265 if (localScvfIdx < numBoxScvf)
266 return BoxHelper::getScvfCorners(type, trans, localScvfIdx);
267
268 const auto localOverlappingScvfIdx = localScvfIdx-numBoxScvf;
269 if (type == Dune::GeometryTypes::triangle)
270 {
272 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localOverlappingScvfIdx]);
273 }
274 else if (type == Dune::GeometryTypes::quadrilateral)
275 {
277 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localOverlappingScvfIdx]);
278 }
279 else if (type == Dune::GeometryTypes::tetrahedron)
280 {
282 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localOverlappingScvfIdx]);
283 }
284 else if (type == Dune::GeometryTypes::hexahedron)
285 {
287 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localOverlappingScvfIdx]);
288 }
289 else
290 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scvf geometries for dim=" << dim
291 << " dimWorld=" << dimWorld
292 << " type=" << type);
293 }
294
295 Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
296 {
297 const auto numBoxScvf = boxHelper_.numInteriorScvf();
298 if (localScvfIdx < numBoxScvf)
299 return Dune::GeometryTypes::cube(dim-1);
300 else
301 return Dune::GeometryTypes::simplex(dim-1);
302 }
303
305 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex,
306 unsigned int indexInFacet) const
307 {
308 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
309 }
310
311 Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
312 {
313 return Dune::GeometryTypes::cube(dim-1);
314 }
315
316 template<int d = dimWorld, std::enable_if_t<(d==3), int> = 0>
317 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
318 {
319 auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]);
320 normal /= normal.two_norm();
321
322 GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]);
323
324 const auto s = v*normal;
325 if (std::signbit(s))
326 normal *= -1;
327
328 return normal;
329 }
330
331 template<int d = dimWorld, std::enable_if_t<(d==2), int> = 0>
332 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
333 {
335 const auto t = p[1] - p[0];
336 GlobalPosition normal({-t[1], t[0]});
337 normal /= normal.two_norm();
338
339 GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]);
340
341 const auto s = v*normal;
342 if (std::signbit(s))
343 normal *= -1;
344
345 return normal;
346 }
347
349 const typename Element::Geometry& elementGeometry() const
350 { return geo_; }
351
353 static auto numInteriorScvf(Dune::GeometryType type)
354 {
355 return BoxHelper::numInteriorScvf(type) + Dune::referenceElement<Scalar, dim>(type).size(dim);
356 }
357
359 static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
360 {
361 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
362 }
363
365 std::size_t numScv() const
366 {
367 return boxHelper_.numScv() + 1;
368 }
369
371 Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage& p) const
372 {
373 const auto scvType = getScvGeometryType(localScvIdx);
374 if constexpr (dim == 3)
375 if (scvType == Dune::GeometryTypes::none(dim))
376 return octahedronVolume_(p);
377
378 return Dumux::convexPolytopeVolume<dim>(
379 scvType,
380 [&](unsigned int i){ return p[i]; }
381 );
382 }
383
385 static std::size_t numElementDofs(Dune::GeometryType type)
386 {
387 return Dune::referenceElement<Scalar, dim>(type).size(dim) + 1;
388 }
389
391 static std::size_t numNonCVLocalDofs(Dune::GeometryType type)
392 {
393 return 0;
394 }
395
397 static auto numLocalDofsIntersection(Dune::GeometryType type, unsigned int iIdx)
398 {
399 return Dune::referenceElement<Scalar, dim>(type).size(iIdx, 1, dim);
400 }
401
403 static auto localDofIndexIntersection(Dune::GeometryType type, unsigned int iIdx, unsigned int ilocalDofIdx)
404 {
405 return Dune::referenceElement<Scalar, dim>(type).subEntity(iIdx, 1, ilocalDofIdx, dim);
406 }
407
408 template<class DofMapper>
409 static auto dofIndex(const DofMapper& dofMapper, const Element& element, unsigned int localDofIdx)
410 {
411 if (localDofIdx < numElementDofs(element.type())-1)
412 return dofMapper.subIndex(element, localDofIdx, dim);
413 else
414 return dofMapper.index(element);
415 }
416
417 GlobalPosition dofPosition(unsigned int localDofIdx) const
418 {
419 if (localDofIdx < numElementDofs(geo_.type())-1)
420 return geo_.corner(localDofIdx);
421 else
422 return geo_.center();
423 }
424
426 template<class LocalKey>
427 static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey& localKey)
428 {
429 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
430 }
431
432 std::array<LocalIndexType, 2> getScvPairForScvf(unsigned int localScvfIndex) const
433 {
434 const auto numEdges = referenceElement(geo_).size(dim-1);
435 if (localScvfIndex < numEdges)
436 return {
437 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
438 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
439 };
440 else
441 return {
442 static_cast<LocalIndexType>(numElementDofs(geo_.type())-1),
443 static_cast<LocalIndexType>(localScvfIndex-numEdges)
444 };
445 }
446
447 std::array<LocalIndexType, 2> getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
448 {
449 const LocalIndexType insideScvIdx
450 = static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
451 return { insideScvIdx, insideScvIdx };
452 }
453
454 bool isOverlappingScvf(unsigned int localScvfIndex) const
455 {
456 if (localScvfIndex < boxHelper_.numInteriorScvf())
457 return false;
458 else
459 return true;
460 }
461
462 bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
463 {
464 return false;
465 }
466
467 bool isOverlappingScv(unsigned int localScvIndex) const
468 {
469 if (localScvIndex < boxHelper_.numScv())
470 return false;
471 else
472 return true;
473 }
474
476 static Element::Geometry::LocalCoordinate localScvfCenter(Dune::GeometryType type, unsigned int localScvfIdx)
477 {
478 return Dumux::center(getScvfCorners(type, [&](const auto& local){ return local; }, localScvfIdx));
479 }
480
482 static Element::Geometry::LocalCoordinate localBoundaryScvfCenter(Dune::GeometryType type, unsigned int localFacetIndex, unsigned int indexInFace)
483 {
484 return Dumux::center(BoxHelper::getBoundaryScvfCorners(type, [&](const auto& local){ return local; }, localFacetIndex, indexInFace));
485 }
486
487private:
488 Scalar octahedronVolume_(const ScvCornerStorage& p) const
489 {
490 using std::abs;
491 return 1.0/6.0 * (
492 abs(Dumux::tripleProduct(p[4]-p[0], p[1]-p[0], p[2]-p[0]))
493 + abs(Dumux::tripleProduct(p[5]-p[0], p[1]-p[0], p[2]-p[0]))
494 );
495 }
496
497 const typename Element::Geometry& geo_;
498 BoxHelper boxHelper_;
499};
500
501template <class GridView, class ScvType, class ScvfType>
503{
504 using Scalar = typename GridView::ctype;
505 using GlobalPosition = typename Dune::FieldVector<Scalar, GridView::dimensionworld>;
506 using ScvCornerStorage = typename ScvType::Traits::CornerStorage;
507 using ScvfCornerStorage = typename ScvfType::Traits::CornerStorage;
508 using LocalIndexType = typename ScvType::Traits::LocalIndexType;
509
510 using Element = typename GridView::template Codim<0>::Entity;
511 using Intersection = typename GridView::Intersection;
512
513 static constexpr auto dim = GridView::dimension;
514 static constexpr auto dimWorld = GridView::dimensionworld;
515
517public:
518
519 HybridPQ1BubbleGeometryHelper(const typename Element::Geometry& geometry)
520 : geo_(geometry)
521 , boxHelper_(geometry)
522 {}
523
525 ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
526 {
527 return getScvCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localScvIdx);
528 }
529
531 template<class Transformation>
532 static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localScvIdx)
533 {
534 // proceed according to number of corners of the element
535 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
536 const auto numBoxScv = ref.size(dim);
537 // reuse box geometry helper for the corner scvs
538 if (localScvIdx < numBoxScv)
539 return BoxHelper::getScvCorners(type, trans, localScvIdx);
540
541 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv corners call for hybrid dofs");
542 }
543
544 Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
545 {
546 // proceed according to number of corners of the element
547 const auto numBoxScv = boxHelper_.numScv();
548
549 if (localScvIdx < numBoxScv)
550 return Dune::GeometryTypes::cube(dim);
551
552 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv geometry call for hybrid dofs");
553 }
554
556 ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
557 {
558 return getScvfCorners(geo_.type(), [&](const auto& local){ return geo_.global(local); }, localScvfIdx);
559 }
560
562 template<class Transformation>
563 static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation&& trans, unsigned int localScvfIdx)
564 {
565 // proceed according to number of corners
566 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
567 const auto numBoxScvf = ref.size(dim-1);
568 // reuse box geometry helper for scvfs
569 if (localScvfIdx < numBoxScvf)
570 return BoxHelper::getScvfCorners(type, trans, localScvfIdx);
571
572 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scvf corners call for hybrid dofs");
573 }
574
575 Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
576 {
577 const auto numBoxScvf = boxHelper_.numInteriorScvf();
578 if (localScvfIdx < numBoxScvf)
579 return Dune::GeometryTypes::cube(dim-1);
580
581 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble interior scvf geometry type call for hybrid dofs");
582 }
583
585 ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex,
586 unsigned int indexInFacet) const
587 {
588 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
589 }
590
591 Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
592 {
593 return Dune::GeometryTypes::cube(dim-1);
594 }
595
596 template<int d = dimWorld, std::enable_if_t<(d==3), int> = 0>
597 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
598 {
599 auto normal = Dumux::crossProduct(p[1]-p[0], p[2]-p[0]);
600 normal /= normal.two_norm();
601
602 GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]);
603
604 const auto s = v*normal;
605 if (std::signbit(s))
606 normal *= -1;
607
608 return normal;
609 }
610
611 template<int d = dimWorld, std::enable_if_t<(d==2), int> = 0>
612 GlobalPosition normal(const ScvfCornerStorage& p, const std::array<LocalIndexType, 2>& scvPair)
613 {
615 const auto t = p[1] - p[0];
616 GlobalPosition normal({-t[1], t[0]});
617 normal /= normal.two_norm();
618
619 GlobalPosition v = dofPosition(scvPair[1]) - dofPosition(scvPair[0]);
620
621 const auto s = v*normal;
622 if (std::signbit(s))
623 normal *= -1;
624
625 return normal;
626 }
627
629 const typename Element::Geometry& elementGeometry() const
630 { return geo_; }
631
633 static auto numInteriorScvf(Dune::GeometryType type)
634 {
635 return BoxHelper::numInteriorScvf(type);
636 }
637
639 static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
640 {
641 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
642 }
643
645 std::size_t numScv() const
646 {
647 return boxHelper_.numScv();
648 }
649
651 Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage& p) const
652 {
653 const auto scvType = getScvGeometryType(localScvIdx);
654
655 return Dumux::convexPolytopeVolume<dim>(
656 scvType,
657 [&](unsigned int i){ return p[i]; }
658 );
659 }
660
662 static std::size_t numElementDofs(Dune::GeometryType type)
663 {
664 return Dune::referenceElement<Scalar, dim>(type).size(dim) + 1;
665 }
666
668 static std::size_t numNonCVLocalDofs(Dune::GeometryType type)
669 {
670 return 1;
671 }
672
674 static auto numLocalDofsIntersection(Dune::GeometryType type, unsigned int iIdx)
675 {
676 return Dune::referenceElement<Scalar, dim>(type).size(iIdx, 1, dim);
677 }
678
680 static auto localDofIndexIntersection(Dune::GeometryType type, unsigned int iIdx, unsigned int ilocalDofIdx)
681 {
682 return Dune::referenceElement<Scalar, dim>(type).subEntity(iIdx, 1, ilocalDofIdx, dim);
683 }
684
685 template<class DofMapper>
686 static auto dofIndex(const DofMapper& dofMapper, const Element& element, unsigned int localDofIdx)
687 {
688 if (localDofIdx < numElementDofs(element.type())-1)
689 return dofMapper.subIndex(element, localDofIdx, dim);
690 else
691 return dofMapper.index(element);
692 }
693
694 GlobalPosition dofPosition(unsigned int localDofIdx) const
695 {
696 if (localDofIdx < numElementDofs(geo_.type())-1)
697 return geo_.corner(localDofIdx);
698 else
699 return geo_.center();
700 }
701
703 template<class LocalKey>
704 static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey& localKey)
705 {
706 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
707 }
708
709 std::array<LocalIndexType, 2> getScvPairForScvf(unsigned int localScvfIndex) const
710 {
711 const auto numEdges = referenceElement(geo_).size(dim-1);
712 if (localScvfIndex < numEdges)
713 return {
714 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
715 static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
716 };
717
718 DUNE_THROW(Dune::NotImplemented, "PQ1Bubble scv pair call for hybrid dofs");
719 }
720
721 std::array<LocalIndexType, 2> getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
722 {
723 const LocalIndexType insideScvIdx
724 = static_cast<LocalIndexType>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
725 return { insideScvIdx, insideScvIdx };
726 }
727
728 bool isOverlappingScvf(unsigned int localScvfIndex) const
729 { return false; }
730
731 bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
732 { return false; }
733
734 bool isOverlappingScv(unsigned int localScvIndex) const
735 { return false; }
736
738 static Element::Geometry::LocalCoordinate localScvfCenter(Dune::GeometryType type, unsigned int localScvfIdx)
739 {
740 return Dumux::center(getScvfCorners(type, [&](const auto& local){ return local; }, localScvfIdx));
741 }
742
744 static Element::Geometry::LocalCoordinate localBoundaryScvfCenter(Dune::GeometryType type, unsigned int localFacetIndex, unsigned int indexInFace)
745 {
746 return Dumux::center(BoxHelper::getBoundaryScvfCorners(type, [&](const auto& local){ return local; }, localFacetIndex, indexInFace));
747 }
748
749private:
750 const typename Element::Geometry& geo_;
752};
753
754} // end namespace Dumux
755
756#endif
Helper class constructing the dual grid finite volume geometries for the box discretizazion method.
Compute the center point of a convex polytope geometry or a random-access container of corner points.
Create sub control volumes and sub control volume face geometries.
Definition: boxgeometryhelper.hh:257
Definition: discretization/pq1bubble/geometryhelper.hh:503
bool isOverlappingScvf(unsigned int localScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:728
std::size_t numScv() const
number of sub control volumes (number of codim-1 entities)
Definition: discretization/pq1bubble/geometryhelper.hh:645
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition: discretization/pq1bubble/geometryhelper.hh:597
bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:731
bool isOverlappingScv(unsigned int localScvIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:734
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition: discretization/pq1bubble/geometryhelper.hh:629
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:544
static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey &localKey)
local dof position
Definition: discretization/pq1bubble/geometryhelper.hh:704
static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvfIdx)
Create a vector with the corners of sub control volume faces.
Definition: discretization/pq1bubble/geometryhelper.hh:563
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition: discretization/pq1bubble/geometryhelper.hh:585
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:575
static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvIdx)
Create a vector with the scv corners.
Definition: discretization/pq1bubble/geometryhelper.hh:532
static Element::Geometry::LocalCoordinate localBoundaryScvfCenter(Dune::GeometryType type, unsigned int localFacetIndex, unsigned int indexInFace)
local boundary scvf center
Definition: discretization/pq1bubble/geometryhelper.hh:744
static std::size_t numNonCVLocalDofs(Dune::GeometryType type)
number of hybrid dofs
Definition: discretization/pq1bubble/geometryhelper.hh:668
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition: discretization/pq1bubble/geometryhelper.hh:525
static auto localDofIndexIntersection(Dune::GeometryType type, unsigned int iIdx, unsigned int ilocalDofIdx)
Local dof index related to a localDof, with index ilocalDofIdx, on an intersection with index iIdx.
Definition: discretization/pq1bubble/geometryhelper.hh:680
std::array< LocalIndexType, 2 > getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:721
static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
number of boundary sub control volume faces for face localFacetIndex
Definition: discretization/pq1bubble/geometryhelper.hh:639
Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage &p) const
get scv volume
Definition: discretization/pq1bubble/geometryhelper.hh:651
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition: discretization/pq1bubble/geometryhelper.hh:556
GlobalPosition dofPosition(unsigned int localDofIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:694
std::array< LocalIndexType, 2 > getScvPairForScvf(unsigned int localScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:709
static std::size_t numElementDofs(Dune::GeometryType type)
number of element dofs
Definition: discretization/pq1bubble/geometryhelper.hh:662
static Element::Geometry::LocalCoordinate localScvfCenter(Dune::GeometryType type, unsigned int localScvfIdx)
local scvf center
Definition: discretization/pq1bubble/geometryhelper.hh:738
HybridPQ1BubbleGeometryHelper(const typename Element::Geometry &geometry)
Definition: discretization/pq1bubble/geometryhelper.hh:519
static auto dofIndex(const DofMapper &dofMapper, const Element &element, unsigned int localDofIdx)
Definition: discretization/pq1bubble/geometryhelper.hh:686
static auto numInteriorScvf(Dune::GeometryType type)
number of interior sub control volume faces
Definition: discretization/pq1bubble/geometryhelper.hh:633
static auto numLocalDofsIntersection(Dune::GeometryType type, unsigned int iIdx)
Number of local dofs related to an intersection with index iIdx.
Definition: discretization/pq1bubble/geometryhelper.hh:674
Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:591
A class to create sub control volume and sub control volume face geometries per element.
Definition: discretization/pq1bubble/geometryhelper.hh:167
PQ1BubbleGeometryHelper(const typename Element::Geometry &geometry)
Definition: discretization/pq1bubble/geometryhelper.hh:183
static auto numBoundaryScvf(Dune::GeometryType type, unsigned int localFacetIndex)
number of boundary sub control volume faces for face localFacetIndex
Definition: discretization/pq1bubble/geometryhelper.hh:359
bool isOverlappingBoundaryScvf(unsigned int localFacetIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:462
bool isOverlappingScvf(unsigned int localScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:454
static std::size_t numNonCVLocalDofs(Dune::GeometryType type)
number of hybrid dofs
Definition: discretization/pq1bubble/geometryhelper.hh:391
std::size_t numScv() const
number of sub control volumes (number of codim-1 entities)
Definition: discretization/pq1bubble/geometryhelper.hh:365
static auto numInteriorScvf(Dune::GeometryType type)
number of interior sub control volume faces
Definition: discretization/pq1bubble/geometryhelper.hh:353
static ScvfCornerStorage getScvfCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvfIdx)
Create a vector with the corners of sub control volume faces.
Definition: discretization/pq1bubble/geometryhelper.hh:259
static std::size_t numElementDofs(Dune::GeometryType type)
number of element dofs
Definition: discretization/pq1bubble/geometryhelper.hh:385
static ScvCornerStorage getScvCorners(Dune::GeometryType type, Transformation &&trans, unsigned int localScvIdx)
Create a vector with the scv corners.
Definition: discretization/pq1bubble/geometryhelper.hh:196
Dune::GeometryType getBoundaryScvfGeometryType(unsigned int localScvfIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:311
static Element::Geometry::LocalCoordinate localBoundaryScvfCenter(Dune::GeometryType type, unsigned int localFacetIndex, unsigned int indexInFace)
local boundary scvf center
Definition: discretization/pq1bubble/geometryhelper.hh:482
ScvfCornerStorage getBoundaryScvfCorners(unsigned int localFacetIndex, unsigned int indexInFacet) const
Create the sub control volume face geometries on the boundary.
Definition: discretization/pq1bubble/geometryhelper.hh:305
std::array< LocalIndexType, 2 > getScvPairForScvf(unsigned int localScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:432
static auto dofIndex(const DofMapper &dofMapper, const Element &element, unsigned int localDofIdx)
Definition: discretization/pq1bubble/geometryhelper.hh:409
GlobalPosition dofPosition(unsigned int localDofIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:417
Dune::GeometryType getInteriorScvfGeometryType(unsigned int localScvfIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:295
static auto localDofIndexIntersection(Dune::GeometryType type, unsigned int iIdx, unsigned int ilocalDofIdx)
Local dof index related to a localDof, with index ilocalDofIdx, on an intersection with index iIdx.
Definition: discretization/pq1bubble/geometryhelper.hh:403
bool isOverlappingScv(unsigned int localScvIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:467
static auto numLocalDofsIntersection(Dune::GeometryType type, unsigned int iIdx)
Number of local dofs related to an intersection with index iIdx.
Definition: discretization/pq1bubble/geometryhelper.hh:397
Scalar scvVolume(unsigned int localScvIdx, const ScvCornerStorage &p) const
get scv volume
Definition: discretization/pq1bubble/geometryhelper.hh:371
const Element::Geometry & elementGeometry() const
the wrapped element geometry
Definition: discretization/pq1bubble/geometryhelper.hh:349
ScvfCornerStorage getScvfCorners(unsigned int localScvfIdx) const
Create a vector with the corners of sub control volume faces.
Definition: discretization/pq1bubble/geometryhelper.hh:252
std::array< LocalIndexType, 2 > getScvPairForBoundaryScvf(unsigned int localFacetIndex, unsigned int localIsScvfIndex) const
Definition: discretization/pq1bubble/geometryhelper.hh:447
static Element::Geometry::LocalCoordinate localScvfCenter(Dune::GeometryType type, unsigned int localScvfIdx)
local scvf center
Definition: discretization/pq1bubble/geometryhelper.hh:476
GlobalPosition normal(const ScvfCornerStorage &p, const std::array< LocalIndexType, 2 > &scvPair)
Definition: discretization/pq1bubble/geometryhelper.hh:317
ScvCornerStorage getScvCorners(unsigned int localScvIdx) const
Create a vector with the scv corners.
Definition: discretization/pq1bubble/geometryhelper.hh:189
Dune::GeometryType getScvGeometryType(unsigned int localScvIdx) const
Definition: discretization/pq1bubble/geometryhelper.hh:232
static Element::Geometry::LocalCoordinate localDofPosition(Dune::GeometryType type, const LocalKey &localKey)
local dof position
Definition: discretization/pq1bubble/geometryhelper.hh:427
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
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:700
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
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:88
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:89
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:53
Definition: discretization/pq1bubble/geometryhelper.hh:70
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:71
Definition: discretization/pq1bubble/geometryhelper.hh:79
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:80
Definition: discretization/pq1bubble/geometryhelper.hh:61
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:62
Definition: discretization/pq1bubble/geometryhelper.hh:48
Definition: discretization/pq1bubble/geometryhelper.hh:145
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:146
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:102
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:133
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:134
Definition: discretization/pq1bubble/geometryhelper.hh:110
std::pair< std::uint8_t, std::uint8_t > Key
Definition: discretization/pq1bubble/geometryhelper.hh:111
Definition: discretization/pq1bubble/geometryhelper.hh:97
Definition: discretization/pq1bubble/geometryhelper.hh:40
Dune::ReservedVector< Dune::FieldVector< ct, cdim >,(1<< mydim)+1 > Type
Definition: discretization/pq1bubble/geometryhelper.hh:41
Traits for an efficient corner storage for the PQ1Bubble method.
Definition: discretization/pq1bubble/geometryhelper.hh:35
Compute the volume of several common geometry types.