Interface to saturation functions from deck.
More...
#include <SaturationPropsFromDeck.hpp>
|
| SaturationPropsFromDeck () |
| Default constructor. More...
|
|
void | init (const PhaseUsage &phaseUsage, std::shared_ptr< MaterialLawManager > materialLawManager) |
|
void | init (const Opm::Deck &deck, std::shared_ptr< MaterialLawManager > materialLawManager) |
|
int | numPhases () const |
|
void | relperm (const int n, const double *s, const int *cells, double *kr, double *dkrds) const |
|
void | capPress (const int n, const double *s, const int *cells, double *pc, double *dpcds) const |
|
void | satRange (const int n, const int *cells, double *smin, double *smax) const |
|
void | updateSatHyst (const int n, const int *cells, const double *s) |
|
void | setGasOilHystParams (const int n, const int *cells, const double *pcswmdc, const double *krnswdc) |
|
void | getGasOilHystParams (const int n, const int *cells, double *pcswmdc, double *krnswdc) const |
|
void | setOilWaterHystParams (const int n, const int *cells, const double *pcswmdc, const double *krnswdc) |
|
void | getOilWaterHystParams (const int n, const int *cells, double *pcswmdc, double *krnswdc) const |
|
void | swatInitScaling (const int cell, const double pcow, double &swat) |
|
const MaterialLawManager & | materialLawManager () const |
| Returns a reference to the MaterialLawManager. More...
|
|
Interface to saturation functions from deck.
◆ MaterialLawManager
◆ MaterialTraits
◆ PhaseIndex
Enumerator |
---|
Aqua | |
Liquid | |
Vapour | |
Solvent | |
Polymer | |
Energy | |
◆ SaturationPropsFromDeck()
Opm::SaturationPropsFromDeck::SaturationPropsFromDeck |
( |
| ) |
|
◆ capPress()
void Opm::SaturationPropsFromDeck::capPress |
( |
const int |
n, |
|
|
const double * |
s, |
|
|
const int * |
cells, |
|
|
double * |
pc, |
|
|
double * |
dpcds |
|
) |
| const |
|
virtual |
Capillary pressure. - Parameters
-
[in] | n | Number of data points. |
[in] | s | Array of nP saturation values. |
[out] | pc | Array of nP capillary pressure values, array must be valid before calling. |
[out] | dpcds | If non-null: array of nP^2 derivative values, array must be valid before calling. The P^2 derivative matrix is m_{ij} = \frac{dpc_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m01 ...) |
Implements Opm::SaturationPropsInterface.
◆ getGasOilHystParams()
void Opm::SaturationPropsFromDeck::getGasOilHystParams |
( |
const int |
n, |
|
|
const int * |
cells, |
|
|
double * |
pcswmdc, |
|
|
double * |
krnswdc |
|
) |
| const |
Get hysteresis parameters for gas-oil - Parameters
-
[in] | n | Number of data points. |
[in] | pcswmdc | Array of hysteresis parameters ( |
- See also
- EclHysteresisTwoPhaseLawParams::pcSwMdc(...))
- Parameters
-
[in] | krnswdc | Array of hysteresis parameters ( |
- See also
- EclHysteresisTwoPhaseLawParams::krnSwMdc(...))
◆ getOilWaterHystParams()
void Opm::SaturationPropsFromDeck::getOilWaterHystParams |
( |
const int |
n, |
|
|
const int * |
cells, |
|
|
double * |
pcswmdc, |
|
|
double * |
krnswdc |
|
) |
| const |
Get hysteresis parameters for oil-water - Parameters
-
[in] | n | Number of data points. |
[in] | pcswmdc | Array of hysteresis parameters ( |
- See also
- EclHysteresisTwoPhaseLawParams::pcSwMdc(...))
- Parameters
-
[in] | krnswdc | Array of hysteresis parameters ( |
- See also
- EclHysteresisTwoPhaseLawParams::krnSwMdc(...))
◆ init() [1/2]
void Opm::SaturationPropsFromDeck::init |
( |
const Opm::Deck & |
deck, |
|
|
std::shared_ptr< MaterialLawManager > |
materialLawManager |
|
) |
| |
|
inline |
◆ init() [2/2]
Initialize from a MaterialLawManager object. - Parameters
-
[in] | phaseUsage | Phase configuration |
[in] | materialLawManager | An initialized MaterialLawManager object |
Referenced by init().
◆ materialLawManager()
Returns a reference to the MaterialLawManager.
Referenced by init().
◆ numPhases()
int Opm::SaturationPropsFromDeck::numPhases |
( |
| ) |
const |
|
virtual |
◆ relperm()
void Opm::SaturationPropsFromDeck::relperm |
( |
const int |
n, |
|
|
const double * |
s, |
|
|
const int * |
cells, |
|
|
double * |
kr, |
|
|
double * |
dkrds |
|
) |
| const |
|
virtual |
Relative permeability. - Parameters
-
[in] | n | Number of data points. |
[in] | s | Array of nP saturation values. |
[out] | kr | Array of nP relperm values, array must be valid before calling. |
[out] | dkrds | If non-null: array of nP^2 relperm derivative values, array must be valid before calling. The P^2 derivative matrix is m_{ij} = \frac{dkr_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m01 ...) |
Implements Opm::SaturationPropsInterface.
◆ satRange()
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] | n | Number of data points. |
[out] | smin | Array of nP minimum s values, array must be valid before calling. |
[out] | smax | Array of nP maximum s values, array must be valid before calling. |
Implements Opm::SaturationPropsInterface.
◆ setGasOilHystParams()
void Opm::SaturationPropsFromDeck::setGasOilHystParams |
( |
const int |
n, |
|
|
const int * |
cells, |
|
|
const double * |
pcswmdc, |
|
|
const double * |
krnswdc |
|
) |
| |
Set hysteresis parameters for gas-oil - Parameters
-
[in] | n | Number of data points. |
[in] | pcswmdc | Array of hysteresis parameters ( |
- See also
- EclHysteresisTwoPhaseLawParams::pcSwMdc(...))
- Parameters
-
[in] | krnswdc | Array of hysteresis parameters ( |
- See also
- EclHysteresisTwoPhaseLawParams::krnSwMdc(...))
◆ setOilWaterHystParams()
void Opm::SaturationPropsFromDeck::setOilWaterHystParams |
( |
const int |
n, |
|
|
const int * |
cells, |
|
|
const double * |
pcswmdc, |
|
|
const double * |
krnswdc |
|
) |
| |
Set hysteresis parameters for oil-water - Parameters
-
[in] | n | Number of data points. |
[in] | pcswmdc | Array of hysteresis parameters ( |
- See also
- EclHysteresisTwoPhaseLawParams::pcSwMdc(...))
- Parameters
-
[in] | krnswdc | Array of hysteresis parameters ( |
- See also
- EclHysteresisTwoPhaseLawParams::krnSwMdc(...))
◆ swatInitScaling()
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] | cell | Cell index. |
[in] | pcow | P_oil - P_water. |
| [in/out] | swat Water saturation. / Possibly modified Water saturation.
|
Implements Opm::SaturationPropsInterface.
◆ updateSatHyst()
void Opm::SaturationPropsFromDeck::updateSatHyst |
( |
const int |
n, |
|
|
const int * |
cells, |
|
|
const double * |
s |
|
) |
| |
|
virtual |
Update saturation state for the hysteresis tracking - Parameters
-
[in] | n | Number of data points. |
[in] | s | Array of nP saturation values.
|
Implements Opm::SaturationPropsInterface.
◆ MaxNumPhases
const int Opm::BlackoilPhases::MaxNumPhases = 3 |
|
staticinherited |
◆ NumCryptoPhases
const int Opm::BlackoilPhases::NumCryptoPhases = 3 |
|
staticinherited |
The documentation for this class was generated from the following file:
|