|
| | 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.
|
|
| 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.
|
|
|
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_ |
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
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
-
| intersection | Info on cell intersection |
| inside | Info describing inside face |
| outside | Info describing outside face |
| faceNormal | Face (unit) normal vector |
- Warning
- Not implemented!
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
-
| intersection | Info on cell intersection |
| inside | Info describing inside face |
| outside | Info describing outside face |
| faceNormal | Face (unit) normal |
- Warning
- LGR computations not implemented!