13#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH
14#define DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH
18#include <dune/common/exceptions.hh>
20#include <dune/geometry/type.hh>
21#include <dune/geometry/referenceelements.hh>
22#include <dune/geometry/multilineargeometry.hh>
23#include <dune/common/reservedvector.hh>
38 template<
int mydim,
int cdim >
41 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)+1>;
47template<Dune::GeometryType::Id gt>
53 using Key = std::pair<std::uint8_t, std::uint8_t>;
54 static constexpr std::array<std::array<Key, 2>, 1> keys = {{
62 using Key = std::pair<std::uint8_t, std::uint8_t>;
63 static constexpr std::array<std::array<Key, 3>, 1> keys = {{
71 using Key = std::pair<std::uint8_t, std::uint8_t>;
72 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
80 using Key = std::pair<std::uint8_t, std::uint8_t>;
81 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
89 using Key = std::pair<std::uint8_t, std::uint8_t>;
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} }
96template<Dune::GeometryType::Id gt>
102 using Key = std::pair<std::uint8_t, std::uint8_t>;
103 static constexpr std::array<std::array<Key, 1>, 1> keys = {{
111 using Key = std::pair<std::uint8_t, std::uint8_t>;
112 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
122 using Key = std::pair<std::uint8_t, std::uint8_t>;
123 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
134 using Key = std::pair<std::uint8_t, std::uint8_t>;
135 static constexpr std::array<std::array<Key, 3>, 10> keys = {{
146 using Key = std::pair<std::uint8_t, std::uint8_t>;
147 static constexpr std::array<std::array<Key, 3>, 8> keys = {{
165template <
class Gr
idView,
class ScvType,
class ScvfType>
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;
174 using Element =
typename GridView::template Codim<0>::Entity;
175 using Intersection =
typename GridView::Intersection;
177 static constexpr auto dim = GridView::dimension;
178 static constexpr auto dimWorld = GridView::dimensionworld;
185 , boxHelper_(geometry)
191 return getScvCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvIdx);
195 template<
class Transformation>
196 static ScvCornerStorage
getScvCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvIdx)
199 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
200 const auto numBoxScv = ref.size(dim);
202 if (localScvIdx < numBoxScv)
203 return BoxHelper::getScvCorners(type, trans, localScvIdx);
205 const auto localOverlappingScvIdx = localScvIdx-numBoxScv;
206 if (type == Dune::GeometryTypes::triangle)
209 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localOverlappingScvIdx]);
211 else if (type == Dune::GeometryTypes::quadrilateral)
214 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localOverlappingScvIdx]);
216 else if (type == Dune::GeometryTypes::tetrahedron)
219 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localOverlappingScvIdx]);
221 else if (type == Dune::GeometryTypes::hexahedron)
224 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(ref, trans, Corners::keys[localOverlappingScvIdx]);
227 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv geometries for dim=" << dim
228 <<
" dimWorld=" << dimWorld
229 <<
" type=" << type);
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)
246 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv geometries for dim=" << dim
247 <<
" dimWorld=" << dimWorld
248 <<
" type=" << type);
254 return getScvfCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvfIdx);
258 template<
class Transformation>
259 static ScvfCornerStorage
getScvfCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvfIdx)
262 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
263 const auto numBoxScvf = ref.size(dim-1);
265 if (localScvfIdx < numBoxScvf)
266 return BoxHelper::getScvfCorners(type, trans, localScvfIdx);
268 const auto localOverlappingScvfIdx = localScvfIdx-numBoxScvf;
269 if (type == Dune::GeometryTypes::triangle)
272 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localOverlappingScvfIdx]);
274 else if (type == Dune::GeometryTypes::quadrilateral)
277 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localOverlappingScvfIdx]);
279 else if (type == Dune::GeometryTypes::tetrahedron)
282 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localOverlappingScvfIdx]);
284 else if (type == Dune::GeometryTypes::hexahedron)
287 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(ref, trans, Corners::keys[localOverlappingScvfIdx]);
290 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scvf geometries for dim=" << dim
291 <<
" dimWorld=" << dimWorld
292 <<
" type=" << type);
297 const auto numBoxScvf = boxHelper_.numInteriorScvf();
298 if (localScvfIdx < numBoxScvf)
299 return Dune::GeometryTypes::cube(dim-1);
301 return Dune::GeometryTypes::simplex(dim-1);
306 unsigned int indexInFacet)
const
308 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
313 return Dune::GeometryTypes::cube(dim-1);
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)
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)
335 const auto t = p[1] - p[0];
336 GlobalPosition
normal({-t[1], t[0]});
355 return BoxHelper::numInteriorScvf(type) + Dune::referenceElement<Scalar, dim>(type).size(dim);
361 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
367 return boxHelper_.numScv() + 1;
371 Scalar
scvVolume(
unsigned int localScvIdx,
const ScvCornerStorage& p)
const
374 if constexpr (dim == 3)
376 return octahedronVolume_(p);
378 return Dumux::convexPolytopeVolume<dim>(
380 [&](
unsigned int i){
return p[i]; }
387 return Dune::referenceElement<Scalar, dim>(type).size(dim) + 1;
399 return Dune::referenceElement<Scalar, dim>(type).size(iIdx, 1, dim);
405 return Dune::referenceElement<Scalar, dim>(type).subEntity(iIdx, 1, ilocalDofIdx, dim);
408 template<
class DofMapper>
409 static auto dofIndex(
const DofMapper& dofMapper,
const Element& element,
unsigned int localDofIdx)
412 return dofMapper.subIndex(element, localDofIdx, dim);
414 return dofMapper.index(element);
420 return geo_.corner(localDofIdx);
422 return geo_.center();
426 template<
class LocalKey>
427 static Element::Geometry::LocalCoordinate
localDofPosition(Dune::GeometryType type,
const LocalKey& localKey)
429 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
434 const auto numEdges = referenceElement(geo_).size(dim-1);
435 if (localScvfIndex < numEdges)
437 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
438 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
443 static_cast<LocalIndexType
>(localScvfIndex-numEdges)
449 const LocalIndexType insideScvIdx
450 =
static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
451 return { insideScvIdx, insideScvIdx };
456 if (localScvfIndex < boxHelper_.numInteriorScvf())
469 if (localScvIndex < boxHelper_.numScv())
476 static Element::Geometry::LocalCoordinate
localScvfCenter(Dune::GeometryType type,
unsigned int localScvfIdx)
482 static Element::Geometry::LocalCoordinate
localBoundaryScvfCenter(Dune::GeometryType type,
unsigned int localFacetIndex,
unsigned int indexInFace)
484 return Dumux::center(BoxHelper::getBoundaryScvfCorners(type, [&](
const auto& local){
return local; }, localFacetIndex, indexInFace));
488 Scalar octahedronVolume_(
const ScvCornerStorage& p)
const
497 const typename Element::Geometry& geo_;
498 BoxHelper boxHelper_;
501template <
class Gr
idView,
class ScvType,
class ScvfType>
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;
510 using Element =
typename GridView::template Codim<0>::Entity;
511 using Intersection =
typename GridView::Intersection;
513 static constexpr auto dim = GridView::dimension;
514 static constexpr auto dimWorld = GridView::dimensionworld;
521 , boxHelper_(geometry)
527 return getScvCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvIdx);
531 template<
class Transformation>
532 static ScvCornerStorage
getScvCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvIdx)
535 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
536 const auto numBoxScv = ref.size(dim);
538 if (localScvIdx < numBoxScv)
539 return BoxHelper::getScvCorners(type, trans, localScvIdx);
541 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv corners call for hybrid dofs");
547 const auto numBoxScv = boxHelper_.numScv();
549 if (localScvIdx < numBoxScv)
550 return Dune::GeometryTypes::cube(dim);
552 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv geometry call for hybrid dofs");
558 return getScvfCorners(geo_.type(), [&](
const auto& local){ return geo_.global(local); }, localScvfIdx);
562 template<
class Transformation>
563 static ScvfCornerStorage
getScvfCorners(Dune::GeometryType type, Transformation&& trans,
unsigned int localScvfIdx)
566 const auto& ref = Dune::referenceElement<Scalar, dim>(type);
567 const auto numBoxScvf = ref.size(dim-1);
569 if (localScvfIdx < numBoxScvf)
570 return BoxHelper::getScvfCorners(type, trans, localScvfIdx);
572 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scvf corners call for hybrid dofs");
577 const auto numBoxScvf = boxHelper_.numInteriorScvf();
578 if (localScvfIdx < numBoxScvf)
579 return Dune::GeometryTypes::cube(dim-1);
581 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble interior scvf geometry type call for hybrid dofs");
586 unsigned int indexInFacet)
const
588 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
593 return Dune::GeometryTypes::cube(dim-1);
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)
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)
615 const auto t = p[1] - p[0];
616 GlobalPosition
normal({-t[1], t[0]});
635 return BoxHelper::numInteriorScvf(type);
641 return Dune::referenceElement<Scalar, dim>(type).size(localFacetIndex, 1, dim);
647 return boxHelper_.numScv();
651 Scalar
scvVolume(
unsigned int localScvIdx,
const ScvCornerStorage& p)
const
655 return Dumux::convexPolytopeVolume<dim>(
657 [&](
unsigned int i){
return p[i]; }
664 return Dune::referenceElement<Scalar, dim>(type).size(dim) + 1;
676 return Dune::referenceElement<Scalar, dim>(type).size(iIdx, 1, dim);
682 return Dune::referenceElement<Scalar, dim>(type).subEntity(iIdx, 1, ilocalDofIdx, dim);
685 template<
class DofMapper>
686 static auto dofIndex(
const DofMapper& dofMapper,
const Element& element,
unsigned int localDofIdx)
689 return dofMapper.subIndex(element, localDofIdx, dim);
691 return dofMapper.index(element);
697 return geo_.corner(localDofIdx);
699 return geo_.center();
703 template<
class LocalKey>
704 static Element::Geometry::LocalCoordinate
localDofPosition(Dune::GeometryType type,
const LocalKey& localKey)
706 return Dune::referenceElement<Scalar, dim>(type).position(localKey.subEntity(), localKey.codim());
711 const auto numEdges = referenceElement(geo_).size(dim-1);
712 if (localScvfIndex < numEdges)
714 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
715 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
718 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv pair call for hybrid dofs");
723 const LocalIndexType insideScvIdx
724 =
static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
725 return { insideScvIdx, insideScvIdx };
738 static Element::Geometry::LocalCoordinate
localScvfCenter(Dune::GeometryType type,
unsigned int localScvfIdx)
744 static Element::Geometry::LocalCoordinate
localBoundaryScvfCenter(Dune::GeometryType type,
unsigned int localFacetIndex,
unsigned int indexInFace)
746 return Dumux::center(BoxHelper::getBoundaryScvfCorners(type, [&](
const auto& local){
return local; }, localFacetIndex, indexInFace));
750 const typename Element::Geometry& geo_;
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: 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.