3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Static Public Member Functions | List of all members
Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 > Class Template Reference

Dimension-specific mpfa helper class for dim == 3 & dimWorld == 3. More...

#include <dumux/discretization/cellcentered/mpfa/helper.hh>

Description

template<class GridGeometry>
class Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 >

Dimension-specific mpfa helper class for dim == 3 & dimWorld == 3.

Static Public Member Functions

template<class ScvBasis >
static ScvBasis calculateInnerNormals (const ScvBasis &scvBasis)
 Calculates the inner normal vectors to a given scv basis. More...
 
template<class ScvBasis >
static CoordScalar calculateDetX (const ScvBasis &scvBasis)
 Calculates the determinant of an scv basis. This is equal to the cross product for dim = dimWorld = 2. More...
 
static std::size_t getGlobalNumScvf (const GridView &gridView)
 Returns the global number of scvfs in the grid. Assumes the grid to be made up of only basic geometry types. Overload this function if you want to use different geometry types. More...
 
template<class ScvBasis >
static bool isRightHandSystem (const ScvBasis &scvBasis)
 Checks whether or not a given scv basis forms a right hand system. More...
 
template<class ElementGeometry , class ReferenceElement >
static ScvfPositionsOnIntersection computeScvfCornersOnIntersection (const ElementGeometry &eg, const ReferenceElement &refElement, unsigned int indexInElement, unsigned int numCorners)
 Returns a vector containing the positions on a given intersection that are relevant for scvf corner computation. Ordering -> 1: facet center, 2: facet corners, 3: edge centers. More...
 
static ScvfCornerVector getScvfCorners (const ScvfPositionsOnIntersection &p, unsigned int numIntersectionCorners, unsigned int cornerIdx)
 Returns the corners of the sub control volume face constructed in a corner (vertex) of an intersection. More...
 
static CoordScalar computeScvfArea (const ScvfCornerVector &scvfCorners)
 Calculates the area of an scvf. More...
 
static std::size_t getNumLocalScvfs (const Dune::GeometryType &gt)
 Calculates the number of scvfs in a given element geometry type. More...
 

Member Function Documentation

◆ calculateDetX()

template<class GridGeometry >
template<class ScvBasis >
static CoordScalar Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 >::calculateDetX ( const ScvBasis &  scvBasis)
inlinestatic

Calculates the determinant of an scv basis. This is equal to the cross product for dim = dimWorld = 2.

Parameters
scvBasisThe basis of an scv

◆ calculateInnerNormals()

template<class GridGeometry >
template<class ScvBasis >
static ScvBasis Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 >::calculateInnerNormals ( const ScvBasis &  scvBasis)
inlinestatic

Calculates the inner normal vectors to a given scv basis.

Parameters
scvBasisThe basis of an scv

◆ computeScvfArea()

template<class GridGeometry >
static CoordScalar Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 >::computeScvfArea ( const ScvfCornerVector &  scvfCorners)
inlinestatic

Calculates the area of an scvf.

Parameters
scvfCornersContainer with the corners of the scvf

◆ computeScvfCornersOnIntersection()

template<class GridGeometry >
template<class ElementGeometry , class ReferenceElement >
static ScvfPositionsOnIntersection Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 >::computeScvfCornersOnIntersection ( const ElementGeometry &  eg,
const ReferenceElement &  refElement,
unsigned int  indexInElement,
unsigned int  numCorners 
)
inlinestatic

Returns a vector containing the positions on a given intersection that are relevant for scvf corner computation. Ordering -> 1: facet center, 2: facet corners, 3: edge centers.

Parameters
egGeometry of the element the facet is embedded in
refElementReference element of the element the facet is embedded in
indexInElementThe local index of the facet in the element
numCornersThe number of corners on the facet

◆ getGlobalNumScvf()

template<class GridGeometry >
static std::size_t Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 >::getGlobalNumScvf ( const GridView &  gridView)
inlinestatic

Returns the global number of scvfs in the grid. Assumes the grid to be made up of only basic geometry types. Overload this function if you want to use different geometry types.

Parameters
gridViewThe grid view to be checked

◆ getNumLocalScvfs()

template<class GridGeometry >
static std::size_t Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 >::getNumLocalScvfs ( const Dune::GeometryType &  gt)
inlinestatic

Calculates the number of scvfs in a given element geometry type.

Parameters
gtThe element geometry type

◆ getScvfCorners()

template<class GridGeometry >
static ScvfCornerVector Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 >::getScvfCorners ( const ScvfPositionsOnIntersection &  p,
unsigned int  numIntersectionCorners,
unsigned int  cornerIdx 
)
inlinestatic

Returns the corners of the sub control volume face constructed in a corner (vertex) of an intersection.

Parameters
pContainer with all scvf corners of the intersection
numIntersectionCornersNumber of corners of the intersection
cornerIdxLocal vertex index on the intersection

◆ isRightHandSystem()

template<class GridGeometry >
template<class ScvBasis >
static bool Dumux::MpfaDimensionHelper< GridGeometry, 3, 3 >::isRightHandSystem ( const ScvBasis &  scvBasis)
inlinestatic

Checks whether or not a given scv basis forms a right hand system.

Parameters
scvBasisThe basis of an scv

The documentation for this class was generated from the following file: