opm-simulators
Loading...
Searching...
No Matches
Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar > Class Template Reference

Cell face properties needed in TPSA equation calculations. More...

#include <FacePropertiesTPSA.hpp>

Classes

struct  FaceInfo

Public Types

using DimVector = Dune::FieldVector<Scalar, dimWorld>

Public Member Functions

 FacePropertiesTPSA (const EclipseState &eclState, const GridView &gridView, const CartesianIndexMapper &cartMapper, const Grid &grid, std::function< std::array< double, dimWorld >(int)> centroids)
 Constructor.
void finishInit ()
 Compute TPSA face properties.
void update ()
 Compute TPSA face properties.
Scalar weightAverage (unsigned elemIdx1, unsigned elemIdx2) const
 Average (half-)weight at interface between two elements.
Scalar weightAverageBoundary (unsigned elemIdx1, unsigned boundaryFaceIdx) const
 Average (half-)weight at boundary interface.
Scalar weightProduct (unsigned elemIdx1, unsigned elemIdx2) const
 Product of weights at interface between two elements.
Scalar weightProductBoundary (unsigned elemIdx1, unsigned boundaryFaceIdx) const
 Product of weights at boundary interface.
Scalar normalDistance (unsigned elemIdx1, unsigned elemIdx2) const
 Distance between two elements.
Scalar normalDistanceBoundary (unsigned elemIdx1, unsigned boundaryFaceIdx) const
 Distance to boundary interface.
DimVector cellFaceNormal (unsigned elemIdx1, unsigned elemIdx2)
 Cell face normal at interface between two elements.
const DimVector & cellFaceNormalBoundary (unsigned elemIdx1, unsigned boundaryFaceIdx) const
 Cell face normal of boundary interface.
const Scalar shearModulus (unsigned elemIdx) const
 Return shear modulus of an element.

Protected Member Functions

template<class Intersection>
void computeCellProperties (const Intersection &intersection, FaceInfo &inside, FaceInfo &outside, DimVector &faceNormal, std::false_type) const
 Compute face properties from general DUNE grid.
template<class Intersection>
void computeCellProperties (const Intersection &intersection, FaceInfo &inside, FaceInfo &outside, DimVector &faceNormal, std::true_type) const
 Compute face properties from DUNE CpGrid.
Scalar computeDistance_ (const DimVector &distVec, const DimVector &faceNormal)
 Compute normal distance from cell center to face center.
Scalar computeWeight_ (const Scalar distance, const Scalar smod)
 Compute weight ratio between distance and shear modulus.
DimVector distanceVector_ (const DimVector &faceCenter, const unsigned &cellIdx) const
 Distance vector from cell center to face center.
void extractSModulus_ ()
 Extract shear modulus from eclState.

Protected Attributes

std::vector< Scalar > sModulus_
std::unordered_map< std::uint64_t, Scalar > weightsAvg_
std::unordered_map< std::uint64_t, Scalar > weightsProd_
std::unordered_map< std::uint64_t, Scalar > distance_
std::unordered_map< std::uint64_t, DimVector > faceNormal_
std::map< std::pair< unsigned, unsigned >, Scalar > weightsAvgBoundary_
std::map< std::pair< unsigned, unsigned >, Scalar > weightsProdBoundary_
std::map< std::pair< unsigned, unsigned >, Scalar > distanceBoundary_
std::map< std::pair< unsigned, unsigned >, DimVector > faceNormalBoundary_
const EclipseState & eclState_
const GridView & gridView_
const CartesianIndexMapper & cartMapper_
const Grid & grid_
std::function< std::array< double, dimWorld >(int)> centroids_
std::vector< std::array< double, dimWorld > > centroids_cache_
const LookUpData< Grid, GridView > lookUpData_
const LookUpCartesianData< Grid, GridView > lookUpCartesianData_

Detailed Description

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
class Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >

Cell face properties needed in TPSA equation calculations.

Similar calculations as done in Transmissibility class for TPFA

Constructor & Destructor Documentation

◆ FacePropertiesTPSA()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::FacePropertiesTPSA ( const EclipseState & eclState,
const GridView & gridView,
const CartesianIndexMapper & cartMapper,
const Grid & grid,
std::function< std::array< double, dimWorld >(int)> centroids )

Constructor.

Parameters
eclStateEclipse state made from a deck
gridViewThe view on the DUNE grid which ought to be used (normally the leaf grid view)
cartMapperThe cartesian index mapper for lookup of cartesian indices
gridThe grid to lookup cell properties
centroidsFunction to lookup cell centroids

Member Function Documentation

◆ cellFaceNormal()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::DimVector Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::cellFaceNormal ( unsigned elemIdx1,
unsigned elemIdx2 )

Cell face normal at interface between two elements.

Parameters
elemIdx1Cell index 1
elemIdx2Cell index 2
Returns
Face normals
Note
It is assumed that normals point from lower index to higher index!

◆ cellFaceNormalBoundary()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
const FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::DimVector & Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::cellFaceNormalBoundary ( unsigned elemIdx,
unsigned boundaryFaceIdx ) const

Cell face normal of boundary interface.

Parameters
elemIdxCell index
boundaryFaceIdxFace index
Returns
Face normals
Note
Boundary normals points outwards

◆ computeCellProperties() [1/2]

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
template<class Intersection>
void Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::computeCellProperties ( const Intersection & intersection,
FaceInfo & inside,
FaceInfo & outside,
DimVector & faceNormal,
std::false_type  ) const
protected

Compute face properties from general DUNE grid.

Parameters
intersectionInfo on cell intersection
insideInfo describing inside face
outsideInfo describing outside face
faceNormalFace (unit) normal vector
Warning
Not implemented!

◆ computeCellProperties() [2/2]

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
template<class Intersection>
void Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::computeCellProperties ( const Intersection & intersection,
FaceInfo & inside,
FaceInfo & outside,
DimVector & faceNormal,
std::true_type  ) const
protected

Compute face properties from DUNE CpGrid.

Parameters
intersectionInfo on cell intersection
insideInfo describing inside face
outsideInfo describing outside face
faceNormalFace (unit) normal
Warning
LGR computations not implemented!

◆ computeDistance_()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
Scalar Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::computeDistance_ ( const DimVector & distVec,
const DimVector & faceNormal )
protected

Compute normal distance from cell center to face center.

Parameters
distVecDistance vector from cell center to face center
faceNormalFace (unit) normal vector
Returns
Distance cell center -> face center

◆ computeWeight_()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
Scalar Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::computeWeight_ ( const Scalar distance,
const Scalar smod )
protected

Compute weight ratio between distance and shear modulus.

Parameters
distanceNormal distance from cell to face center
smodShear modulus
Returns
Weight ratio

◆ distanceVector_()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::DimVector Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::distanceVector_ ( const DimVector & faceCenter,
const unsigned & cellIdx ) const
protected

Distance vector from cell center to face center.

Parameters
faceCenterFace center coordinates
cellIdxCell index
Returns
Distance vector cell center -> face center

◆ extractSModulus_()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
void Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::extractSModulus_ ( )
protected

Extract shear modulus from eclState.

Note
(from Transmissibility::extractPorosity()): All arrays provided by eclState are one-per-cell of "uncompressed" grid, whereas the simulation grid might remove a few elements (e.g. because it is distributed over several processes).

◆ normalDistance()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
Scalar Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::normalDistance ( unsigned elemIdx1,
unsigned elemIdx2 ) const

Distance between two elements.

Parameters
elemIdx1Cell index 1
elemIdx2Cell index 2
Returns
Distance

◆ normalDistanceBoundary()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
Scalar Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::normalDistanceBoundary ( unsigned elemIdx,
unsigned boundaryFaceIdx ) const

Distance to boundary interface.

Parameters
elemIdxCell index
boundaryFaceIdxFace index
Returns
Distance to boundary

◆ weightAverage()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
Scalar Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::weightAverage ( unsigned elemIdx1,
unsigned elemIdx2 ) const

Average (half-)weight at interface between two elements.

Parameters
elemIdx1Cell index 1
elemIdx2Cell index 2
Returns
Average weight

◆ weightAverageBoundary()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
Scalar Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::weightAverageBoundary ( unsigned elemIdx,
unsigned boundaryFaceIdx ) const

Average (half-)weight at boundary interface.

Parameters
elemIdxCell index
boundaryFaceIdxFace index
Returns
Average weight at boundary

◆ weightProduct()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
Scalar Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::weightProduct ( unsigned elemIdx1,
unsigned elemIdx2 ) const

Product of weights at interface between two elements.

Parameters
elemIdx1Cell index 1
elemIdx2Cell index 2
Returns
Weight product

◆ weightProductBoundary()

template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
Scalar Opm::FacePropertiesTPSA< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::weightProductBoundary ( unsigned elemIdx,
unsigned boundaryFaceIdx ) const

Product of weights at boundary interface.

Parameters
elemIdxCell index
boundaryFaceIdxFace index
Returns
Weight product at boundary

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