All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Pages
Opm::SaturationPropsFromDeck Class Reference

Interface to saturation functions from deck. More...

#include <SaturationPropsFromDeck.hpp>

Inheritance diagram for Opm::SaturationPropsFromDeck:
Opm::SaturationPropsInterface Opm::BlackoilPhases

Public Types

typedef
Opm::ThreePhaseMaterialTraits
< double, BlackoilPhases::Aqua,
BlackoilPhases::Liquid,
BlackoilPhases::Vapour > 
MaterialTraits
 
typedef
Opm::EclMaterialLawManager
< MaterialTraits
MaterialLawManager
 
- Public Types inherited from Opm::BlackoilPhases
enum  PhaseIndex { Aqua = 0, Liquid = 1, Vapour = 2 }
 

Public Member Functions

 SaturationPropsFromDeck ()
 Default constructor.
 
void init (const PhaseUsage &phaseUsage, std::shared_ptr< MaterialLawManager > materialLawManager)
 Initialize from a MaterialLawManager object. More...
 
void init (const Opm::Deck &deck, std::shared_ptr< MaterialLawManager > materialLawManager)
 Initialize from deck and MaterialLawManager. More...
 
int numPhases () const
 
void relperm (const int n, const double *s, const int *cells, double *kr, double *dkrds) const
 Relative permeability. More...
 
void capPress (const int n, const double *s, const int *cells, double *pc, double *dpcds) const
 Capillary pressure. More...
 
void satRange (const int n, const int *cells, double *smin, double *smax) const
 Obtain the range of allowable saturation values. More...
 
void updateSatHyst (const int n, const int *cells, const double *s)
 Update saturation state for the hysteresis tracking. More...
 
void setGasOilHystParams (const int n, const int *cells, const double *pcswmdc, const double *krnswdc)
 Set hysteresis parameters for gas-oil. More...
 
void getGasOilHystParams (const int n, const int *cells, double *pcswmdc, double *krnswdc) const
 Get hysteresis parameters for gas-oil. More...
 
void setOilWaterHystParams (const int n, const int *cells, const double *pcswmdc, const double *krnswdc)
 Set hysteresis parameters for oil-water. More...
 
void getOilWaterHystParams (const int n, const int *cells, double *pcswmdc, double *krnswdc) const
 Get hysteresis parameters for oil-water. More...
 
void swatInitScaling (const int cell, const double pcow, double &swat)
 Update capillary pressure scaling according to pressure diff. More...
 
const MaterialLawManagermaterialLawManager () const
 Returns a reference to the MaterialLawManager.
 
- Public Member Functions inherited from Opm::SaturationPropsInterface
virtual ~SaturationPropsInterface ()
 Virtual destructor.
 

Additional Inherited Members

- Static Public Attributes inherited from Opm::BlackoilPhases
static const int MaxNumPhases = 3
 

Detailed Description

Interface to saturation functions from deck.

Member Function Documentation

void Opm::SaturationPropsFromDeck::capPress ( const int  n,
const double *  s,
const int *  cells,
double *  pc,
double *  dpcds 
) const
virtual

Capillary pressure.

Parameters
[in]nNumber of data points.
[in]sArray of nP saturation values.
[out]pcArray of nP capillary pressure values, array must be valid before calling.
[out]dpcdsIf non-null: array of nP^2 derivative values, array must be valid before calling. The P^2 derivative matrix is m_{ij} = {dpc_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m01 ...)
[in]nNumber of data points.
[in]sArray of nP saturation values.
[in]cellsArray of n cell indices to be associated with the s values.
[out]pcArray of nP capillary pressure values, array must be valid before calling.
[out]dpcdsIf non-null: array of nP^2 derivative values, array must be valid before calling. The P^2 derivative matrix is m_{ij} = {dpc_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m01 ...)

Implements Opm::SaturationPropsInterface.

void Opm::SaturationPropsFromDeck::getGasOilHystParams ( const int  n,
const int *  cells,
double *  pcswmdc,
double *  krnswdc 
) const

Get hysteresis parameters for gas-oil.

Parameters
[in]nNumber of data points.
[in]pcswmdcArray of hysteresis parameters (
See Also
EclHysteresisTwoPhaseLawParams::pcSwMdc(...))
Parameters
[in]krnswdcArray of hysteresis parameters (
See Also
EclHysteresisTwoPhaseLawParams::krnSwMdc(...))
void Opm::SaturationPropsFromDeck::getOilWaterHystParams ( const int  n,
const int *  cells,
double *  pcswmdc,
double *  krnswdc 
) const

Get hysteresis parameters for oil-water.

Parameters
[in]nNumber of data points.
[in]pcswmdcArray of hysteresis parameters (
See Also
EclHysteresisTwoPhaseLawParams::pcSwMdc(...))
Parameters
[in]krnswdcArray of hysteresis parameters (
See Also
EclHysteresisTwoPhaseLawParams::krnSwMdc(...))
void Opm::SaturationPropsFromDeck::init ( const PhaseUsage phaseUsage,
std::shared_ptr< MaterialLawManager materialLawManager 
)

Initialize from a MaterialLawManager object.

Initialize from deck.

Parameters
[in]phaseUsagePhase configuration
[in]materialLawManagerAn initialized MaterialLawManager object
void Opm::SaturationPropsFromDeck::init ( const Opm::Deck &  deck,
std::shared_ptr< MaterialLawManager materialLawManager 
)
inline

Initialize from deck and MaterialLawManager.

Parameters
[in]deckInput deck
[in]materialLawManagerAn initialized MaterialLawManager object
int Opm::SaturationPropsFromDeck::numPhases ( ) const
virtual
Returns
P, the number of phases.

Implements Opm::SaturationPropsInterface.

void Opm::SaturationPropsFromDeck::relperm ( const int  n,
const double *  s,
const int *  cells,
double *  kr,
double *  dkrds 
) const
virtual

Relative permeability.

Parameters
[in]nNumber of data points.
[in]sArray of nP saturation values.
[out]krArray of nP relperm values, array must be valid before calling.
[out]dkrdsIf non-null: array of nP^2 relperm derivative values, array must be valid before calling. The P^2 derivative matrix is m_{ij} = {dkr_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m01 ...)
[in]nNumber of data points.
[in]sArray of nP saturation values.
[in]cellsArray of n cell indices to be associated with the s values.
[out]krArray of nP relperm values, array must be valid before calling.
[out]dkrdsIf non-null: array of nP^2 relperm derivative values, array must be valid before calling. The P^2 derivative matrix is m_{ij} = {dkr_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m01 ...)

Implements Opm::SaturationPropsInterface.

void Opm::SaturationPropsFromDeck::satRange ( const int  n,
const int *  cells,
double *  smin,
double *  smax 
) const
virtual

Obtain the range of allowable saturation values.

Parameters
[in]nNumber of data points.
[out]sminArray of nP minimum s values, array must be valid before calling.
[out]smaxArray of nP maximum s values, array must be valid before calling.
[in]nNumber of data points.
[in]cellsArray of n cell indices.
[out]sminArray of nP minimum s values, array must be valid before calling.
[out]smaxArray of nP maximum s values, array must be valid before calling.

Implements Opm::SaturationPropsInterface.

void Opm::SaturationPropsFromDeck::setGasOilHystParams ( const int  n,
const int *  cells,
const double *  pcswmdc,
const double *  krnswdc 
)

Set hysteresis parameters for gas-oil.

Parameters
[in]nNumber of data points.
[in]pcswmdcArray of hysteresis parameters (
See Also
EclHysteresisTwoPhaseLawParams::pcSwMdc(...))
Parameters
[in]krnswdcArray of hysteresis parameters (
See Also
EclHysteresisTwoPhaseLawParams::krnSwMdc(...))
void Opm::SaturationPropsFromDeck::setOilWaterHystParams ( const int  n,
const int *  cells,
const double *  pcswmdc,
const double *  krnswdc 
)

Set hysteresis parameters for oil-water.

Parameters
[in]nNumber of data points.
[in]pcswmdcArray of hysteresis parameters (
See Also
EclHysteresisTwoPhaseLawParams::pcSwMdc(...))
Parameters
[in]krnswdcArray of hysteresis parameters (
See Also
EclHysteresisTwoPhaseLawParams::krnSwMdc(...))
void Opm::SaturationPropsFromDeck::swatInitScaling ( const int  cell,
const double  pcow,
double &  swat 
)
virtual

Update capillary pressure scaling according to pressure diff.

and initial water saturation.

Parameters
[in]cellCell index.
[in]pcowP_oil - P_water.
in/out]swat Water saturation. / Possibly modified Water saturation.

Implements Opm::SaturationPropsInterface.

void Opm::SaturationPropsFromDeck::updateSatHyst ( const int  n,
const int *  cells,
const double *  s 
)
virtual

Update saturation state for the hysteresis tracking.

Parameters
[in]nNumber of data points.
[in]sArray of nP saturation values.

Implements Opm::SaturationPropsInterface.


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