Opm::IncompPropertiesBasic Class Reference

#include <IncompPropertiesBasic.hpp>

Inheritance diagram for Opm::IncompPropertiesBasic:
Inheritance graph

Public Member Functions

 IncompPropertiesBasic (const parameter::ParameterGroup &param, const int dim, const int num_cells)
 
 IncompPropertiesBasic (const int num_phases, const SaturationPropsBasic::RelPermFunc &relpermfunc, const std::vector< double > &rho, const std::vector< double > &mu, const double porosity, const double permeability, const int dim, const int num_cells)
 
virtual ~IncompPropertiesBasic ()
 Destructor. More...
 
virtual int numDimensions () const
 
virtual int numCells () const
 
virtual const double * porosity () const
 
virtual const double * permeability () const
 
virtual int numPhases () const
 
virtual const double * viscosity () const
 
virtual const double * density () const
 
virtual const double * surfaceDensity () const
 
virtual void relperm (const int n, const double *s, const int *cells, double *kr, double *dkrds) const
 
virtual void capPress (const int n, const double *s, const int *cells, double *pc, double *dpcds) const
 
virtual void satRange (const int n, const int *cells, double *smin, double *smax) const
 

Detailed Description

Concrete class implementing the incompressible property interface, reading all necessary input from parameters.

Supports variable number of spatial dimensions, called D. Supports variable number of phases, called P. In general, when arguments call for n values of some vector or matrix property, such as saturation, they shall always be ordered cellwise:

\[ [s^1_0, s^2_0, s^3_0, s^1_1, s^2_2, \ldots ] \]

in which $ s^i_j $ denotes saturation of phase i in cell j.

Constructor & Destructor Documentation

Opm::IncompPropertiesBasic::IncompPropertiesBasic ( const parameter::ParameterGroup param,
const int  dim,
const int  num_cells 
)

Construct from parameters. Note that all values passed through param should be in convenient units, as documented below. The following parameters are accepted (defaults):

  • num_phases (2) – Must be 1 or 2.
  • relperm_func ("Linear") – Must be "Constant", "Linear" or "Quadratic".
  • rho1 rho2, rho3 (1.0e3) – Density in kg/m^3.
  • mu1 mu2, mu3 (1.0) – Viscosity in cP.
  • porosity (1.0) – Porosity.
  • permeability (100.0) – Permeability in mD.
Opm::IncompPropertiesBasic::IncompPropertiesBasic ( const int  num_phases,
const SaturationPropsBasic::RelPermFunc relpermfunc,
const std::vector< double > &  rho,
const std::vector< double > &  mu,
const double  porosity,
const double  permeability,
const int  dim,
const int  num_cells 
)

Construct properties from arguments. Note that all values should be given in SI units for this constructor.

Parameters
[in]num_phasesMust be 1 or 2.
[in]rhoPhase densities in kg/m^3.
[in]muPhase viscosities in Pa*s.
[in]porosityMust be in [0,1].
[in]permeabilityPermeability in m^2.
[in]dimMust be 2 or 3.
[in]num_cellsThe number of grid cells.
virtual Opm::IncompPropertiesBasic::~IncompPropertiesBasic ( )
virtual

Destructor.

Member Function Documentation

virtual void Opm::IncompPropertiesBasic::capPress ( const int  n,
const double *  s,
const int *  cells,
double *  pc,
double *  dpcds 
) const
virtual
Parameters
[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 m_01 ...)

Implements Opm::IncompPropertiesInterface.

virtual const double* Opm::IncompPropertiesBasic::density ( ) const
virtual

Densities of fluid phases at reservoir conditions.

Returns
Array of P density values.

Implements Opm::IncompPropertiesInterface.

virtual int Opm::IncompPropertiesBasic::numCells ( ) const
virtual
Returns
N, the number of cells.

Implements Opm::IncompPropertiesInterface.

virtual int Opm::IncompPropertiesBasic::numDimensions ( ) const
virtual
Returns
D, the number of spatial dimensions.

Implements Opm::IncompPropertiesInterface.

virtual int Opm::IncompPropertiesBasic::numPhases ( ) const
virtual
Returns
P, the number of phases (also the number of components).

Implements Opm::IncompPropertiesInterface.

virtual const double* Opm::IncompPropertiesBasic::permeability ( ) const
virtual
Returns
Array of ND^2 permeability values. The D^2 permeability values for a cell are organized as a matrix, which is symmetric (so ordering does not matter).

Implements Opm::IncompPropertiesInterface.

virtual const double* Opm::IncompPropertiesBasic::porosity ( ) const
virtual
Returns
Array of N porosity values.

Implements Opm::IncompPropertiesInterface.

virtual void Opm::IncompPropertiesBasic::relperm ( const int  n,
const double *  s,
const int *  cells,
double *  kr,
double *  dkrds 
) const
virtual
Parameters
[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 m_01 ...)

Implements Opm::IncompPropertiesInterface.

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

Obtain the range of allowable saturation values. In cell cells[i], saturation of phase p is allowed to be in the interval [smin[i*P + p], smax[i*P + p]].

Parameters
[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::IncompPropertiesInterface.

virtual const double* Opm::IncompPropertiesBasic::surfaceDensity ( ) const
virtual

Densities of fluid phases at surface conditions.

Returns
Array of P density values.

Implements Opm::IncompPropertiesInterface.

virtual const double* Opm::IncompPropertiesBasic::viscosity ( ) const
virtual
Returns
Array of P viscosity values.

Implements Opm::IncompPropertiesInterface.


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