Opm::SaturationPropsBasic Class Reference

#include <SaturationPropsBasic.hpp>

Public Types

enum  RelPermFunc { Constant, Linear, Quadratic }
 

Public Member Functions

 SaturationPropsBasic ()
 Default constructor. More...
 
void init (const parameter::ParameterGroup &param)
 
void init (const int num_phases, const RelPermFunc &relperm_func)
 Initialize from arguments a basic Saturation property. More...
 
int numPhases () const
 
void relperm (const int n, const double *s, double *kr, double *dkrds) const
 
void capPress (const int n, const double *s, double *pc, double *dpcds) const
 
void satRange (const int n, double *smin, double *smax) const
 

Detailed Description

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.

Member Enumeration Documentation

Enumerator
Constant 
Linear 
Quadratic 

Constructor & Destructor Documentation

Opm::SaturationPropsBasic::SaturationPropsBasic ( )

Default constructor.

Member Function Documentation

void Opm::SaturationPropsBasic::capPress ( const int  n,
const double *  s,
double *  pc,
double *  dpcds 
) const

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 ...)
void Opm::SaturationPropsBasic::init ( const parameter::ParameterGroup param)

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".
void Opm::SaturationPropsBasic::init ( const int  num_phases,
const RelPermFunc relperm_func 
)
inline

Initialize from arguments a basic Saturation property.

int Opm::SaturationPropsBasic::numPhases ( ) const
Returns
P, the number of phases.
void Opm::SaturationPropsBasic::relperm ( const int  n,
const double *  s,
double *  kr,
double *  dkrds 
) const

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 ...)
void Opm::SaturationPropsBasic::satRange ( const int  n,
double *  smin,
double *  smax 
) const

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.

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