Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar > Class Template Reference

#include <Transmissibility.hpp>

Public Types

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

Public Member Functions

 Transmissibility (const EclipseState &eclState, const GridView &gridView, const CartesianIndexMapper &cartMapper, const Grid &grid, std::function< std::array< double, dimWorld >(int)> centroids, bool enableEnergy, bool enableDiffusivity, bool enableDispersivity)
 
const DimMatrixpermeability (unsigned elemIdx) const
 Return the permeability for an element. More...
 
Scalar transmissibility (unsigned elemIdx1, unsigned elemIdx2) const
 Return the transmissibility for the intersection between two elements. More...
 
Scalar transmissibilityBoundary (unsigned elemIdx, unsigned boundaryFaceIdx) const
 Return the transmissibility for a given boundary segment. More...
 
Scalar thermalHalfTrans (unsigned insideElemIdx, unsigned outsideElemIdx) const
 Return the thermal "half transmissibility" for the intersection between two elements. More...
 
Scalar thermalHalfTransBoundary (unsigned insideElemIdx, unsigned boundaryFaceIdx) const
 
Scalar diffusivity (unsigned elemIdx1, unsigned elemIdx2) const
 Return the diffusivity for the intersection between two elements. More...
 
Scalar dispersivity (unsigned elemIdx1, unsigned elemIdx2) const
 Return the dispersivity for the intersection between two elements. More...
 
void finishInit (const std::function< unsigned int(unsigned int)> &map={})
 Actually compute the transmissibility over a face as a pre-compute step. More...
 
void update (bool global, const std::function< unsigned int(unsigned int)> &map={}, bool applyNncMultRegT=false)
 Compute all transmissibilities. More...
 

Protected Member Functions

void updateFromEclState_ (bool global)
 
void removeNonCartesianTransmissibilities_ (bool removeAll)
 
void applyAllZMultipliers_ (Scalar &trans, unsigned insideFaceIdx, unsigned outsideFaceIdx, unsigned insideCartElemIdx, unsigned outsideCartElemIdx, const TransMult &transMult, const std::array< int, dimWorld > &cartDims, bool pinchTop)
 Apply the Multipliers for the case PINCH(4)==TOPBOT. More...
 
std::array< std::vector< double >, 3 > createTransmissibilityArrays_ (const std::array< bool, 3 > &is_tran)
 Creates TRANS{XYZ} arrays for modification by FieldProps data. More...
 
void resetTransmissibilityFromArrays_ (const std::array< bool, 3 > &is_tran, const std::array< std::vector< double >, 3 > &trans)
 overwrites calculated transmissibilities More...
 
template<class Intersection >
void computeFaceProperties (const Intersection &intersection, const int, const int, const int, const int, DimVector &faceCenterInside, DimVector &faceCenterOutside, DimVector &faceAreaNormal, std::false_type) const
 
template<class Intersection >
void computeFaceProperties (const Intersection &intersection, const int insideElemIdx, const int insideFaceIdx, const int outsideElemIdx, const int outsideFaceIdx, DimVector &faceCenterInside, DimVector &faceCenterOutside, DimVector &faceAreaNormal, std::true_type) const
 
void applyNncToGridTrans_ (const std::unordered_map< std::size_t, int > &cartesianToCompressed)
 
void applyEditNncToGridTrans_ (const std::unordered_map< std::size_t, int > &globalToLocal)
 Multiplies the grid transmissibilities according to EDITNNC. More...
 
void applyEditNncrToGridTrans_ (const std::unordered_map< std::size_t, int > &globalToLocal)
 Resets the grid transmissibilities according to EDITNNCR. More...
 
void applyNncMultreg_ (const std::unordered_map< std::size_t, int > &globalToLocal)
 
void applyEditNncToGridTransHelper_ (const std::unordered_map< std::size_t, int > &globalToLocal, const std::string &keyword, const std::vector< NNCdata > &nncs, const std::function< KeywordLocation(const NNCdata &)> &getLocation, const std::function< void(double &, const double &)> &apply)
 
void extractPermeability_ ()
 
void extractPermeability_ (const std::function< unsigned int(unsigned int)> &map)
 
void extractPorosity_ ()
 
void extractDispersion_ ()
 
void computeHalfTrans_ (Scalar &halfTrans, const DimVector &areaNormal, int faceIdx, const DimVector &distance, const DimMatrix &perm) const
 
void computeHalfDiffusivity_ (Scalar &halfDiff, const DimVector &areaNormal, const DimVector &distance, const Scalar &poro) const
 
DimVector distanceVector_ (const DimVector &center, int faceIdx, unsigned elemIdx, const std::array< std::vector< DimVector >, dimWorld > &axisCentroids) const
 
void applyMultipliers_ (Scalar &trans, unsigned faceIdx, unsigned cartElemIdx, const TransMult &transMult) const
 
void applyNtg_ (Scalar &trans, unsigned faceIdx, unsigned elemIdx, const std::vector< double > &ntg) const
 

Protected Attributes

std::vector< DimMatrixpermeability_
 
std::vector< Scalar > porosity_
 
std::vector< Scalar > dispersion_
 
std::unordered_map< std::uint64_t, Scalar > trans_
 
const EclipseState & eclState_
 
const GridView & gridView_
 
const CartesianIndexMapper & cartMapper_
 
const Grid & grid_
 
std::function< std::array< double, dimWorld >(int)> centroids_
 
Scalar transmissibilityThreshold_
 
std::map< std::pair< unsigned, unsigned >, Scalar > transBoundary_
 
std::map< std::pair< unsigned, unsigned >, Scalar > thermalHalfTransBoundary_
 
bool enableEnergy_
 
bool enableDiffusivity_
 
bool enableDispersivity_
 
std::unordered_map< std::uint64_t, Scalar > thermalHalfTrans_
 
std::unordered_map< std::uint64_t, Scalar > diffusivity_
 
std::unordered_map< std::uint64_t, Scalar > dispersivity_
 
const LookUpData< Grid, GridView > lookUpData_
 
const LookUpCartesianData< Grid, GridView > lookUpCartesianData_
 

Member Typedef Documentation

◆ DimMatrix

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
using Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>

◆ DimVector

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
using Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::DimVector = Dune::FieldVector<Scalar, dimWorld>

Constructor & Destructor Documentation

◆ Transmissibility()

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

Member Function Documentation

◆ applyAllZMultipliers_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::applyAllZMultipliers_ ( Scalar &  trans,
unsigned  insideFaceIdx,
unsigned  outsideFaceIdx,
unsigned  insideCartElemIdx,
unsigned  outsideCartElemIdx,
const TransMult &  transMult,
const std::array< int, dimWorld > &  cartDims,
bool  pinchTop 
)
protected

Apply the Multipliers for the case PINCH(4)==TOPBOT.

Parameters
pinchTopWhether PINCH(5) is TOP, otherwise ALL is assumed.

◆ applyEditNncrToGridTrans_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::applyEditNncrToGridTrans_ ( const std::unordered_map< std::size_t, int > &  globalToLocal)
protected

Resets the grid transmissibilities according to EDITNNCR.

◆ applyEditNncToGridTrans_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::applyEditNncToGridTrans_ ( const std::unordered_map< std::size_t, int > &  globalToLocal)
protected

Multiplies the grid transmissibilities according to EDITNNC.

◆ applyEditNncToGridTransHelper_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::applyEditNncToGridTransHelper_ ( const std::unordered_map< std::size_t, int > &  globalToLocal,
const std::string &  keyword,
const std::vector< NNCdata > &  nncs,
const std::function< KeywordLocation(const NNCdata &)> &  getLocation,
const std::function< void(double &, const double &)> &  apply 
)
protected

References Opm::details::isId().

◆ applyMultipliers_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::applyMultipliers_ ( Scalar &  trans,
unsigned  faceIdx,
unsigned  cartElemIdx,
const TransMult &  transMult 
) const
protected

◆ applyNncMultreg_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::applyNncMultreg_ ( const std::unordered_map< std::size_t, int > &  globalToLocal)
protected

References Opm::details::isId().

◆ applyNncToGridTrans_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::applyNncToGridTrans_ ( const std::unordered_map< std::size_t, int > &  cartesianToCompressed)
protected

References Opm::details::isId().

◆ applyNtg_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::applyNtg_ ( Scalar &  trans,
unsigned  faceIdx,
unsigned  elemIdx,
const std::vector< double > &  ntg 
) const
protected

◆ computeFaceProperties() [1/2]

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
template<class Intersection >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::computeFaceProperties ( const Intersection &  intersection,
const int  insideElemIdx,
const int  insideFaceIdx,
const int  outsideElemIdx,
const int  outsideFaceIdx,
DimVector faceCenterInside,
DimVector faceCenterOutside,
DimVector faceAreaNormal,
std::true_type   
) const
protected

Alternatively, the actual area of the refined face can be computed as follows:

◆ computeFaceProperties() [2/2]

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
template<class Intersection >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::computeFaceProperties ( const Intersection &  intersection,
const int  ,
const int  ,
const int  ,
const int  ,
DimVector faceCenterInside,
DimVector faceCenterOutside,
DimVector faceAreaNormal,
std::false_type   
) const
protected

◆ computeHalfDiffusivity_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::computeHalfDiffusivity_ ( Scalar &  halfDiff,
const DimVector areaNormal,
const DimVector distance,
const Scalar &  poro 
) const
protected

◆ computeHalfTrans_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::computeHalfTrans_ ( Scalar &  halfTrans,
const DimVector areaNormal,
int  faceIdx,
const DimVector distance,
const DimMatrix perm 
) const
protected

◆ createTransmissibilityArrays_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::array< std::vector< double >, 3 > Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::createTransmissibilityArrays_ ( const std::array< bool, 3 > &  is_tran)
protected

Creates TRANS{XYZ} arrays for modification by FieldProps data.

Parameters
is_tranWhether TRAN{XYZ} will be modified. If entry is false the array will be empty
Returns
an array of vector (TRANX, TRANY, TRANZ}

References Opm::details::isId().

◆ diffusivity()

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

Return the diffusivity for the intersection between two elements.

References Opm::details::isId().

◆ dispersivity()

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

Return the dispersivity for the intersection between two elements.

References Opm::details::isId().

◆ distanceVector_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::DimVector Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::distanceVector_ ( const DimVector center,
int  faceIdx,
unsigned  elemIdx,
const std::array< std::vector< DimVector >, dimWorld > &  axisCentroids 
) const
protected

◆ extractDispersion_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::extractDispersion_
protected

◆ extractPermeability_() [1/2]

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::extractPermeability_
protected

◆ extractPermeability_() [2/2]

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::extractPermeability_ ( const std::function< unsigned int(unsigned int)> &  map)
protected

◆ extractPorosity_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::extractPorosity_
protected

◆ finishInit()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::finishInit ( const std::function< unsigned int(unsigned int)> &  map = {})
inline

Actually compute the transmissibility over a face as a pre-compute step.

This code actually uses the direction specific "centroids" of each element. These "centroids" are not the identical barycenter of the element, but the middle of the centers of the faces of the logical Cartesian cells, i.e., the centers of the faces of the reference elements. We do things this way because the barycenter of the element can be located outside of the element for sufficiently "ugly" (i.e., thin and "non-flat") elements which in turn leads to quite wrong permeabilities. This approach is probably not always correct either but at least it seems to be much better.

◆ permeability()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
const DimMatrix & Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::permeability ( unsigned  elemIdx) const
inline

◆ removeNonCartesianTransmissibilities_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::removeNonCartesianTransmissibilities_ ( bool  removeAll)
protected

◆ resetTransmissibilityFromArrays_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::resetTransmissibilityFromArrays_ ( const std::array< bool, 3 > &  is_tran,
const std::array< std::vector< double >, 3 > &  trans 
)
protected

overwrites calculated transmissibilities

Parameters
is_tranWhether TRAN{XYZ} have been modified.
transArrays with modified transmissibilities TRAN{XYZ}

c1 < c2 <=> gc1 < gc2 is no longer true when the grid is a CpGrid with LGRs. When cells c1 and c2 have the same parent cell on level zero, then gc1 == gc2.

References Opm::details::isId().

◆ thermalHalfTrans()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
Scalar Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::thermalHalfTrans ( unsigned  insideElemIdx,
unsigned  outsideElemIdx 
) const

Return the thermal "half transmissibility" for the intersection between two elements.

The "half transmissibility" features all sub-expressions of the "thermal transmissibility" which can be precomputed, i.e. they are not dependent on the current solution:

H_t = A * (n*d)/(d*d);

where A is the area of the intersection between the inside and outside elements, n is the outer unit normal, and d is the distance between the center of the inside cell and the center of the intersection.

References Opm::details::directionalIsId().

◆ thermalHalfTransBoundary()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
Scalar Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::thermalHalfTransBoundary ( unsigned  insideElemIdx,
unsigned  boundaryFaceIdx 
) const

◆ transmissibility()

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

Return the transmissibility for the intersection between two elements.

References Opm::details::isId().

◆ transmissibilityBoundary()

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

Return the transmissibility for a given boundary segment.

◆ update()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::update ( bool  global,
const std::function< unsigned int(unsigned int)> &  map = {},
bool  applyNncMultRegT = false 
)

Compute all transmissibilities.

Parameters
[in]globalWhether or not to call update() on all processes. Also, this updates the "thermal half transmissibilities" if energy is enabled.
[in]mapUndocumented.
[in]applyNncMultRegTWhether or not to apply regional multipliers to explicit NNCs. Explicit NNCs are those entered directly in the input data, e.g., through the NNC/EDITNNC/EDITNNCR keywords, or the result of generating connections to or within numerical aquifers. Default value: false, meaning do not apply regional multipliers to explicit NNCs.

References Opm::details::directionalIsId(), and Opm::details::isId().

◆ updateFromEclState_()

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
void Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::updateFromEclState_ ( bool  global)
protected

Member Data Documentation

◆ cartMapper_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
const CartesianIndexMapper& Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::cartMapper_
protected

◆ centroids_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::function<std::array<double,dimWorld>(int)> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::centroids_
protected

◆ diffusivity_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::unordered_map<std::uint64_t, Scalar> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::diffusivity_
protected

◆ dispersion_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::vector<Scalar> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::dispersion_
protected

◆ dispersivity_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::unordered_map<std::uint64_t, Scalar> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::dispersivity_
protected

◆ eclState_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
const EclipseState& Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::eclState_
protected

◆ enableDiffusivity_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
bool Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::enableDiffusivity_
protected

◆ enableDispersivity_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
bool Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::enableDispersivity_
protected

◆ enableEnergy_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
bool Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::enableEnergy_
protected

◆ grid_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
const Grid& Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::grid_
protected

◆ gridView_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
const GridView& Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::gridView_
protected

◆ lookUpCartesianData_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
const LookUpCartesianData<Grid,GridView> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::lookUpCartesianData_
protected

◆ lookUpData_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
const LookUpData<Grid,GridView> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::lookUpData_
protected

◆ permeability_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::vector<DimMatrix> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::permeability_
protected

◆ porosity_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::vector<Scalar> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::porosity_
protected

◆ thermalHalfTrans_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::unordered_map<std::uint64_t, Scalar> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::thermalHalfTrans_
protected

◆ thermalHalfTransBoundary_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::map<std::pair<unsigned, unsigned>, Scalar> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::thermalHalfTransBoundary_
protected

◆ trans_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::unordered_map<std::uint64_t, Scalar> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::trans_
protected

◆ transBoundary_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
std::map<std::pair<unsigned, unsigned>, Scalar> Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::transBoundary_
protected

◆ transmissibilityThreshold_

template<class Grid , class GridView , class ElementMapper , class CartesianIndexMapper , class Scalar >
Scalar Opm::Transmissibility< Grid, GridView, ElementMapper, CartesianIndexMapper, Scalar >::transmissibilityThreshold_
protected

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