Class encapsulating basic saturation function behaviour, by which we mean constant, linear or quadratic relative permeability functions for a maximum of two phases, and zero capillary pressure.
More...
#include <SaturationPropsBasic.hpp>
|
enum | RelPermFunc { Constant,
Linear,
Quadratic
} |
|
|
| SaturationPropsBasic () |
| Default constructor.
|
|
void | init (const ParameterGroup ¶m) |
| Initialize from parameters. More...
|
|
void | init (const int num_phases, const RelPermFunc &relperm_func) |
| Initialize from arguments a basic Saturation property.
|
|
int | numPhases () const |
|
void | relperm (const int n, const double *s, double *kr, double *dkrds) const |
| Relative permeability. More...
|
|
void | capPress (const int n, const double *s, double *pc, double *dpcds) const |
| Capillary pressure. More...
|
|
void | satRange (const int n, double *smin, double *smax) const |
| Obtain the range of allowable saturation values. More...
|
|
Class encapsulating basic saturation function behaviour, by which we mean constant, linear or quadratic relative permeability functions for a maximum of two phases, and zero capillary pressure.
TODO: This class can easily be extended to three phases, by adding three-phase relperm behaviour.
◆ capPress()
void Opm::SaturationPropsBasic::capPress |
( |
const int |
n, |
|
|
const double * |
s, |
|
|
double * |
pc, |
|
|
double * |
dpcds |
|
) |
| const |
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} = {dpc_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m01 ...) |
◆ init()
Initialize from parameters.
The following parameters are accepted (defaults):
- num_phases (2) – Must be 1 or 2.
- relperm_func ("Linear") – Must be "Constant", "Linear" or "Quadratic".
◆ numPhases()
int Opm::SaturationPropsBasic::numPhases |
( |
| ) |
const |
- Returns
- P, the number of phases.
◆ relperm()
void Opm::SaturationPropsBasic::relperm |
( |
const int |
n, |
|
|
const double * |
s, |
|
|
double * |
kr, |
|
|
double * |
dkrds |
|
) |
| const |
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} = {dkr_i}{ds^j}, and is output in Fortran order (m_00 m_10 m_20 m01 ...) |
◆ satRange()
void Opm::SaturationPropsBasic::satRange |
( |
const int |
n, |
|
|
double * |
smin, |
|
|
double * |
smax |
|
) |
| const |
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. |
The documentation for this class was generated from the following files: