Opm::FluidMatrixInteractionBlackoil< ScalarT, ParamsT > Class Template Reference

Capillary pressures and relative permeabilities for a black oil system. More...

#include <FluidMatrixInteractionBlackoil.hpp>

Inheritance diagram for Opm::FluidMatrixInteractionBlackoil< ScalarT, ParamsT >:
Inheritance graph

Public Types

typedef ParamsT Params
 
typedef Params::Scalar Scalar
 
enum  { numComponents = 3 }
 
enum  { numPhases = 3 }
 
enum  ComponentIndex { Water = 0 , Oil = 1 , Gas = 2 }
 
enum  PhaseIndex { Aqua = 0 , Liquid = 1 , Vapour = 2 }
 
typedef Dune::FieldVector< Scalar, numComponentsCompVec
 
typedef Dune::FieldVector< Scalar, numPhasesPhaseVec
 
typedef Dune::FieldMatrix< Scalar, numComponents, numPhasesPhaseToCompMatrix
 
typedef Dune::FieldMatrix< Scalar, numPhases, numPhasesPhaseJacobian
 

Static Public Member Functions

template<class pcContainerT , class SatContainerT >
static void pC (pcContainerT &pc, const Params &params, const SatContainerT &saturations, Scalar)
 The linear capillary pressure-saturation curve. More...
 
template<class SatContainerT , class pcContainerT >
static void S (SatContainerT &saturations, const Params &params, const pcContainerT &pc, Scalar)
 The saturation-capillary pressure curve. More...
 
template<class krContainerT , class SatContainerT >
static void kr (krContainerT &kr, const Params &params, const SatContainerT &saturations, Scalar)
 The relative permeability of all phases. More...
 
template<class krContainerT , class SatContainerT >
static void dkr (krContainerT &dkr, const Params &params, const SatContainerT &saturations, Scalar)
 The saturation derivatives of relative permeability of all phases. More...
 

Detailed Description

template<class ScalarT, class ParamsT = FluidMatrixInteractionBlackoilParams<ScalarT>>
class Opm::FluidMatrixInteractionBlackoil< ScalarT, ParamsT >

Capillary pressures and relative permeabilities for a black oil system.

Member Typedef Documentation

◆ CompVec

typedef Dune::FieldVector<Scalar, numComponents> Opm::BlackoilDefs::CompVec
inherited

◆ Params

template<class ScalarT , class ParamsT = FluidMatrixInteractionBlackoilParams<ScalarT>>
typedef ParamsT Opm::FluidMatrixInteractionBlackoil< ScalarT, ParamsT >::Params

◆ PhaseJacobian

typedef Dune::FieldMatrix<Scalar, numPhases, numPhases> Opm::BlackoilDefs::PhaseJacobian
inherited

◆ PhaseToCompMatrix

typedef Dune::FieldMatrix<Scalar, numComponents, numPhases> Opm::BlackoilDefs::PhaseToCompMatrix
inherited

◆ PhaseVec

typedef Dune::FieldVector<Scalar, numPhases> Opm::BlackoilDefs::PhaseVec
inherited

◆ Scalar

template<class ScalarT , class ParamsT = FluidMatrixInteractionBlackoilParams<ScalarT>>
typedef Params::Scalar Opm::FluidMatrixInteractionBlackoil< ScalarT, ParamsT >::Scalar

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
inherited
Enumerator
numComponents 

◆ anonymous enum

anonymous enum
inherited
Enumerator
numPhases 

◆ ComponentIndex

Enumerator
Water 
Oil 
Gas 

◆ PhaseIndex

Enumerator
Aqua 
Liquid 
Vapour 

Member Function Documentation

◆ dkr()

template<class ScalarT , class ParamsT = FluidMatrixInteractionBlackoilParams<ScalarT>>
template<class krContainerT , class SatContainerT >
static void Opm::FluidMatrixInteractionBlackoil< ScalarT, ParamsT >::dkr ( krContainerT &  dkr,
const Params params,
const SatContainerT &  saturations,
Scalar   
)
inlinestatic

◆ kr()

template<class ScalarT , class ParamsT = FluidMatrixInteractionBlackoilParams<ScalarT>>
template<class krContainerT , class SatContainerT >
static void Opm::FluidMatrixInteractionBlackoil< ScalarT, ParamsT >::kr ( krContainerT &  kr,
const Params params,
const SatContainerT &  saturations,
Scalar   
)
inlinestatic

◆ pC()

template<class ScalarT , class ParamsT = FluidMatrixInteractionBlackoilParams<ScalarT>>
template<class pcContainerT , class SatContainerT >
static void Opm::FluidMatrixInteractionBlackoil< ScalarT, ParamsT >::pC ( pcContainerT &  pc,
const Params params,
const SatContainerT &  saturations,
Scalar   
)
inlinestatic

The linear capillary pressure-saturation curve.

This material law is linear:

\[
p_C = (1 - \overline{S}_w) (p_{C,max} - p_{C,entry}) + p_{C,entry}
\]

Parameters
SweEffective saturation of of the wetting phase $\overline{S}_w$

References Opm::BlackoilDefs::Aqua, Opm::BlackoilDefs::Liquid, and Opm::BlackoilDefs::Vapour.

◆ S()

template<class ScalarT , class ParamsT = FluidMatrixInteractionBlackoilParams<ScalarT>>
template<class SatContainerT , class pcContainerT >
static void Opm::FluidMatrixInteractionBlackoil< ScalarT, ParamsT >::S ( SatContainerT &  saturations,
const Params params,
const pcContainerT &  pc,
Scalar   
)
inlinestatic

The saturation-capillary pressure curve.

This is the inverse of the capillary pressure-saturation curve:

\[
S_w = 1 - \frac{p_C - p_{C,entry}}{p_{C,max} - p_{C,entry}}
\]

Parameters
pCCapillary pressure $\p_C$
Returns
The effective saturaion of the wetting phase $\overline{S}_w$

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