Opm::spu_2p::ModelParameterStorage Class Reference

#include <SinglePointUpwindTwoPhase.hpp>

Public Member Functions

 ModelParameterStorage (int nc, int totconn)
 
double & drho ()
 
double drho () const
 
double * mob (int c)
 
const double * mob (int c) const
 
double * dmob (int c)
 
const double * dmob (int c) const
 
double * porevol ()
 
double porevol (int c) const
 
double & dg (int i)
 
double dg (int i) const
 
double & ds (int c)
 
double ds (int c) const
 
double & pc (int c)
 
double pc (int c) const
 
double & dpc (int c)
 
double dpc (int c) const
 
double & trans (int f)
 
double trans (int f) const
 

Detailed Description

Internal class to manage the direct and derived quantities needed to formulate the fluid transport system.

Note: This class elides off-diagonal elements of the phase mobility Jacobian, $(\partial_{s_\beta} \lambda_\alpha)_{\alpha\beta}$. These elements are assumed to be strictly equal to zero. In other words, the relative permeability of phase $\alpha$ is assumed to depend only on the saturation of phase $\alpha$. This convention allows storing only the two diagonals of the mobility Jacobian per grid cell.

The static gravity term is the scalar value

\[ \Delta G_i = \mathsf{T}_f\, \vec{g}\cdot(\Bar{x}_f - \Bar{x}_c) \]

in which i is the half face index corresponding to the cell-face pair (f,c) and $\mathsf{T}_f$ is the absolute (bacground) two-point transmissibility of face f.

The fluid transport problem is formulated in terms of saturation changes, $\Delta s$, per cell. These changes are the primary degrees of freedom in this model.

Capillary pressures are defined by the fluid model, but usually correspond to $p_w - p_n$ (e.g., $p_\mathit{oil} - p_\mathit{water}$).

Constructor & Destructor Documentation

Opm::spu_2p::ModelParameterStorage::ModelParameterStorage ( int  nc,
int  totconn 
)
inline

Constructor.

Parameters
[in]ncTotal number of grid cells.
[in]totconnTotal number of connections, accumulated per cell (``half faces'').

Member Function Documentation

double& Opm::spu_2p::ModelParameterStorage::dg ( int  i)
inline

Static gravity term associated to single half face.

Parameters
[in]iHalf face index corresponding to particular cell-face pair.
Returns
Read-write reference to static gravity term of single half face.
double Opm::spu_2p::ModelParameterStorage::dg ( int  i) const
inline

Static gravity term associated to single half face.

Parameters
[in]iHalf face index corresponding to particular cell-face pair.
Returns
Read-only reference to static gravity term of single half face.
double* Opm::spu_2p::ModelParameterStorage::dmob ( int  c)
inline

Diagonal elements of phase mobility derivative (Jacobian).

Parameters
[in]cCell.
Returns
Read-write reference to diagonal elements of phase mobility Jacobian in cell c.

Referenced by Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::initIteration(), and Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::sourceTerms().

const double* Opm::spu_2p::ModelParameterStorage::dmob ( int  c) const
inline

Diagonal elements of phase mobility derivative (Jacobian).

Parameters
[in]cCell.
Returns
Read-only reference to two consecutive diagonal elements of phase mobility Jacobian in cell c.
double& Opm::spu_2p::ModelParameterStorage::dpc ( int  c)
inline

Derivative of capillary pressure with respect to saturation.

Parameters
[in]cCell
Returns
Read-write reference to capillary pressure derivative with respect to primary saturation in cell c.

Referenced by Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::initIteration().

double Opm::spu_2p::ModelParameterStorage::dpc ( int  c) const
inline

Derivative of capillary pressure with respect to saturation.

Parameters
[in]cCell
Returns
Read-only reference to capillary pressure derivative with respect to primary saturation in cell c.
double& Opm::spu_2p::ModelParameterStorage::drho ( )
inline

Modifiable density difference.

Returns
Reference to modifiable internal representation of fluid phase density difference.

Referenced by Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::SinglePointUpwindTwoPhase().

double Opm::spu_2p::ModelParameterStorage::drho ( ) const
inline

Read-only density difference.

Returns
Read-only value of current fluid phase difference value.
double& Opm::spu_2p::ModelParameterStorage::ds ( int  c)
inline

Saturation change in particular cell.

Parameters
[in]c
Returns
Read-write reference to saturation change (scalar) in cell c.

Referenced by Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::accumulation(), and Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::initIteration().

double Opm::spu_2p::ModelParameterStorage::ds ( int  c) const
inline

Saturation change in particular cell.

Parameters
[in]c
Returns
Read-only reference to saturation change (scalar) in cell c.
double* Opm::spu_2p::ModelParameterStorage::mob ( int  c)
inline

Phase mobility in cell.

Parameters
[in]cCell.
Returns
Read-write reference to two consecutive phase mobilities in cell c.

Referenced by Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::initIteration(), and Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::sourceTerms().

const double* Opm::spu_2p::ModelParameterStorage::mob ( int  c) const
inline

Phase mobility in cell.

Parameters
[in]cCell.
Returns
Read-only reference to two consecutive phase mobilities in cell c.
double& Opm::spu_2p::ModelParameterStorage::pc ( int  c)
inline

Capillary pressure in particular cell.

Parameters
[in]cCell.
Returns
Read-write reference to capillary pressure in cell c.

Referenced by Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::initIteration().

double Opm::spu_2p::ModelParameterStorage::pc ( int  c) const
inline

Capillary pressure in particular cell.

Parameters
[in]cCell
Returns
Read-only reference to capillary pressure in cell c.
double* Opm::spu_2p::ModelParameterStorage::porevol ( )
inline
double Opm::spu_2p::ModelParameterStorage::porevol ( int  c) const
inline

Pore volume of single cell.

Parameters
[in]cCell.
Returns
Pore volume of cell c.
double& Opm::spu_2p::ModelParameterStorage::trans ( int  f)
inline

Background (absolute) face transmissibility of particular face.

Parameters
[in]fFace
Returns
Read-write reference to background face transmissibility of face f.

Referenced by Opm::SinglePointUpwindTwoPhase< Opm::Opm::SimpleFluid2pWrappingProps >::initGravityTrans().

double Opm::spu_2p::ModelParameterStorage::trans ( int  f) const
inline

Background (absolute) face transmissibility of particular face.

Parameters
[in]fFace
Returns
Read-only reference to bacground face transmissibility of face f.

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