25#ifndef DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH
26#define DUMUX_DISCRETIZATION_PQ1BUBBLE_GEOMETRY_HELPER_HH
30#include <dune/common/exceptions.hh>
32#include <dune/geometry/type.hh>
33#include <dune/geometry/referenceelements.hh>
34#include <dune/geometry/multilineargeometry.hh>
35#include <dune/common/reservedvector.hh>
49 template<
int mydim,
int cdim >
52 using Type = Dune::ReservedVector< Dune::FieldVector< ct, cdim >, (1<<mydim)+1>;
58template<Dune::GeometryType::Id gt>
64 using Key = std::pair<std::uint8_t, std::uint8_t>;
65 static constexpr std::array<std::array<Key, 2>, 1> keys = {{
73 using Key = std::pair<std::uint8_t, std::uint8_t>;
74 static constexpr std::array<std::array<Key, 3>, 1> keys = {{
82 using Key = std::pair<std::uint8_t, std::uint8_t>;
83 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
91 using Key = std::pair<std::uint8_t, std::uint8_t>;
92 static constexpr std::array<std::array<Key, 4>, 1> keys = {{
100 using Key = std::pair<std::uint8_t, std::uint8_t>;
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} }
107template<Dune::GeometryType::Id gt>
113 using Key = std::pair<std::uint8_t, std::uint8_t>;
114 static constexpr std::array<std::array<Key, 1>, 1> keys = {{
122 using Key = std::pair<std::uint8_t, std::uint8_t>;
123 static constexpr std::array<std::array<Key, 2>, 3> keys = {{
133 using Key = std::pair<std::uint8_t, std::uint8_t>;
134 static constexpr std::array<std::array<Key, 2>, 4> keys = {{
145 using Key = std::pair<std::uint8_t, std::uint8_t>;
146 static constexpr std::array<std::array<Key, 3>, 10> keys = {{
157 using Key = std::pair<std::uint8_t, std::uint8_t>;
158 static constexpr std::array<std::array<Key, 3>, 8> keys = {{
176template <
class Gr
idView,
class ScvType,
class ScvfType>
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;
185 using Element =
typename GridView::template Codim<0>::Entity;
186 using Intersection =
typename GridView::Intersection;
188 static constexpr auto dim = GridView::dimension;
189 static constexpr auto dimWorld = GridView::dimensionworld;
194 , boxHelper_(geometry)
201 const auto type = geo_.type();
202 const auto numBoxScv = boxHelper_.numScv();
204 if (localScvIdx < numBoxScv)
205 return boxHelper_.getScvCorners(localScvIdx);
207 const auto localOverlappingScvIdx = localScvIdx-numBoxScv;
208 if (type == Dune::GeometryTypes::triangle)
211 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
213 else if (type == Dune::GeometryTypes::quadrilateral)
216 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
218 else if (type == Dune::GeometryTypes::tetrahedron)
221 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
223 else if (type == Dune::GeometryTypes::hexahedron)
226 return Detail::Box::keyToCornerStorage<ScvCornerStorage>(geo_, Corners::keys[localOverlappingScvIdx]);
229 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv geometries for dim=" << dim
230 <<
" dimWorld=" << dimWorld
231 <<
" type=" << type);
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)
248 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scv geometries for dim=" << dim
249 <<
" dimWorld=" << dimWorld
250 <<
" type=" << type);
257 const auto type = geo_.type();
258 const auto numBoxScvf = boxHelper_.numInteriorScvf();
260 if (localScvfIdx < numBoxScvf)
261 return boxHelper_.getScvfCorners(localScvfIdx);
263 const auto localOverlappingScvfIdx = localScvfIdx-numBoxScvf;
264 if (type == Dune::GeometryTypes::triangle)
267 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
269 else if (type == Dune::GeometryTypes::quadrilateral)
272 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
274 else if (type == Dune::GeometryTypes::tetrahedron)
277 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
279 else if (type == Dune::GeometryTypes::hexahedron)
282 return Detail::Box::keyToCornerStorage<ScvfCornerStorage>(geo_, Corners::keys[localOverlappingScvfIdx]);
285 DUNE_THROW(Dune::NotImplemented,
"PQ1Bubble scvf geometries for dim=" << dim
286 <<
" dimWorld=" << dimWorld
287 <<
" type=" << type);
292 const auto numBoxScvf = boxHelper_.numInteriorScvf();
293 if (localScvfIdx < numBoxScvf)
294 return Dune::GeometryTypes::cube(dim-1);
296 return Dune::GeometryTypes::simplex(dim-1);
301 unsigned int indexInFacet)
const
303 return boxHelper_.getBoundaryScvfCorners(localFacetIndex, indexInFacet);
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)
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)
325 const auto t = p[1] - p[0];
326 GlobalPosition
normal({-t[1], t[0]});
345 return boxHelper_.numInteriorScvf() + referenceElement(geo_).size(dim);
351 return referenceElement(geo_).size(localFacetIndex, 1, dim);
357 return boxHelper_.numScv() + 1;
361 Scalar
scvVolume(
unsigned int localScvIdx,
const ScvCornerStorage& p)
const
364 if constexpr (dim == 3)
366 return octahedronVolume_(p);
368 return Dumux::convexPolytopeVolume<dim>(
370 [&](
unsigned int i){
return p[i]; }
374 template<
class DofMapper>
375 auto dofIndex(
const DofMapper& dofMapper,
const Element& element,
unsigned int localScvIdx)
const
377 if (localScvIdx <
numScv()-1)
378 return dofMapper.subIndex(element, localScvIdx, dim);
380 return dofMapper.index(element);
385 if (localScvIdx <
numScv()-1)
386 return geo_.corner(localScvIdx);
388 return geo_.center();
393 const auto numEdges = referenceElement(geo_).size(dim-1);
394 if (localScvfIndex < numEdges)
396 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 0, dim)),
397 static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localScvfIndex, dim-1, 1, dim))
401 static_cast<LocalIndexType
>(
numScv()-1),
402 static_cast<LocalIndexType
>(localScvfIndex-numEdges)
408 const LocalIndexType insideScvIdx
409 =
static_cast<LocalIndexType
>(referenceElement(geo_).subEntity(localFacetIndex, 1, localIsScvfIndex, dim));
410 return { insideScvIdx, insideScvIdx };
415 if (localScvfIndex < boxHelper_.numInteriorScvf())
423 if (localScvIndex < boxHelper_.numScv())
430 Scalar octahedronVolume_(
const ScvCornerStorage& p)
const
439 const typename Element::Geometry& geo_;
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