14#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_SCV_GEOMETRY_HELPER_HH
15#define DUMUX_DISCRETIZATION_CC_MPFA_O_SCV_GEOMETRY_HELPER_HH
19#include <dune/common/exceptions.hh>
20#include <dune/geometry/type.hh>
21#include <dune/geometry/multilineargeometry.hh>
30template<
class LocalScvType>
33 using ctype =
typename LocalScvType::ctype;
34 using LocalIndexType =
typename LocalScvType::LocalIndexType;
36 static constexpr int dim = LocalScvType::myDimension;
37 static constexpr int dimWorld = LocalScvType::worldDimension;
39 struct MLGTraits :
public Dune::MultiLinearGeometryTraits<ctype>
42 template<
int mydim,
int cdim >
44 {
using Type = std::array<
typename LocalScvType::GlobalCoordinate, (1<<dim) >; };
50 static const bool v =
true;
51 static const unsigned int topologyId = Dune::GeometryTypes::cube(d).id();
57 using ScvGeometry = Dune::MultiLinearGeometry<ctype, dim, dimWorld, MLGTraits>;
60 template<
class InteractionVolume,
class FVElementGeometry>
62 const InteractionVolume& iv,
63 const FVElementGeometry& fvGeometry)
65 const auto& scv = iv.localScv(ivLocalScvIdx);
69 const auto& firstGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(0)).gridScvfIndex());
70 const auto& secondGridScvf = fvGeometry.scvf(iv.localScvf(scv.localScvfIndex(1)).gridScvfIndex());
72 typename MLGTraits::template CornerStorage<dim, dimWorld>::Type corners;
73 corners[0] = fvGeometry.scv( scv.gridScvIndex() ).center();
74 corners[1] = fvGeometry.facetCorner(firstGridScvf);
75 corners[2] = fvGeometry.facetCorner(secondGridScvf);
76 corners[3] = fvGeometry.vertexCorner(secondGridScvf);
79 typename LocalScvType::LocalBasis basis;
80 basis[0] = corners[1] - corners[0];
81 basis[1] = corners[2] - corners[0];
82 if ( !fvGeometry.gridGeometry().mpfaHelper().isRightHandSystem(basis) )
83 swap(corners[1], corners[2]);
85 return ScvGeometry(Dune::GeometryTypes::cube(ScvGeometry::mydimension), corners);
88 DUNE_THROW(Dune::NotImplemented,
"Mpfa-o local scv geometry computation in 3d");
90 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:32
static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const InteractionVolume &iv, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: scvgeometryhelper.hh:61
Dune::MultiLinearGeometry< ctype, dim, dimWorld, MLGTraits > ScvGeometry
export the geometry type of the local scvs
Definition: scvgeometryhelper.hh:57
Definition: scvgeometryhelper.hh:44
std::array< typename LocalScvType::GlobalCoordinate,(1<< dim) > Type
Definition: scvgeometryhelper.hh:44
Definition: scvgeometryhelper.hh:49
static const bool v
Definition: scvgeometryhelper.hh:50
static const unsigned int topologyId
Definition: scvgeometryhelper.hh:51