26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_SCV_GEOMETRY_HELPER_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_O_SCV_GEOMETRY_HELPER_HH
31#include <dune/common/exceptions.hh>
32#include <dune/geometry/type.hh>
33#include <dune/geometry/multilineargeometry.hh>
42template<
class LocalScvType>
45 using ctype =
typename LocalScvType::ctype;
46 using LocalIndexType =
typename LocalScvType::LocalIndexType;
48 static constexpr int dim = LocalScvType::myDimension;
49 static constexpr int dimWorld = LocalScvType::worldDimension;
51 struct MLGTraits :
public Dune::MultiLinearGeometryTraits<ctype>
54 template<
int mydim,
int cdim >
56 {
using Type = std::array<
typename LocalScvType::GlobalCoordinate, (1<<dim) >; };
62 static const bool v =
true;
63 static const unsigned int topologyId = Dune::GeometryTypes::cube(d).id();
69 using ScvGeometry = Dune::MultiLinearGeometry<ctype, dim, dimWorld, MLGTraits>;
72 template<
class InteractionVolume,
class FVElementGeometry>
74 const InteractionVolume& iv,
75 const FVElementGeometry& fvGeometry)
77 const auto& scv = iv.localScv(ivLocalScvIdx);
81 const auto& firstGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(0)).gridScvfIndex());
82 const auto& secondGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(1)).gridScvfIndex());
84 typename MLGTraits::template CornerStorage<dim, dimWorld>::Type corners;
85 corners[0] = fvGeometry.scv( scv.gridScvIndex() ).center();
86 corners[1] = firstGridScvf.facetCorner();
87 corners[2] = secondGridScvf.facetCorner();
88 corners[3] = secondGridScvf.vertexCorner();
91 typename LocalScvType::LocalBasis basis;
92 basis[0] = corners[1] - corners[0];
93 basis[1] = corners[2] - corners[0];
94 if ( !fvGeometry.gridGeometry().mpfaHelper().isRightHandSystem(basis) )
95 swap(corners[1], corners[2]);
97 return ScvGeometry(Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners);
100 DUNE_THROW(Dune::NotImplemented,
"Mpfa-o local scv geometry computation in 3d");
102 DUNE_THROW(Dune::InvalidStateException,
"Mpfa only works in 2d or 3d");
Helper class providing functionality to compute the geometry of the interaction-volume local sub-cont...
Definition: scvgeometryhelper.hh:44
static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const InteractionVolume &iv, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: scvgeometryhelper.hh:73
Dune::MultiLinearGeometry< ctype, dim, dimWorld, MLGTraits > ScvGeometry
export the geometry type of the local scvs
Definition: scvgeometryhelper.hh:69
Definition: scvgeometryhelper.hh:56
std::array< typename LocalScvType::GlobalCoordinate,(1<< dim) > Type
Definition: scvgeometryhelper.hh:56
Definition: scvgeometryhelper.hh:61
static const bool v
Definition: scvgeometryhelper.hh:62
static const unsigned int topologyId
Definition: scvgeometryhelper.hh:63