My Project
Loading...
Searching...
No Matches
Opm::LinearMaterial< TraitsT, ParamsT > Class Template Reference

Implements a linear saturation-capillary pressure relation. More...

#include <LinearMaterial.hpp>

Inheritance diagram for Opm::LinearMaterial< TraitsT, ParamsT >:

Public Types

typedef TraitsT Traits
 
typedef ParamsT Params
 
typedef Traits::Scalar Scalar
 

Static Public Member Functions

template<class ContainerT , class FluidState >
static void capillaryPressures (ContainerT &values, const Params &params, const FluidState &state)
 The linear capillary pressure-saturation curve.
 
template<class ContainerT , class FluidState >
static void saturations (ContainerT &, const Params &, const FluidState &)
 The inverse of the capillary pressure.
 
template<class ContainerT , class FluidState >
static void relativePermeabilities (ContainerT &values, const Params &, const FluidState &state)
 The relative permeability of all phases.
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation pcnw (const Params &params, const FluidState &fs)
 The difference between the pressures of the non-wetting and wetting phase.
 
template<class Evaluation = Scalar>
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatPcnw (const Params &params, const Evaluation &Sw)
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Sw (const Params &, const FluidState &)
 Calculate wetting phase saturation given that the rest of the fluid state has been initialized.
 
template<class Evaluation = Scalar>
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatSw (const Params &, const Evaluation &)
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Sn (const Params &, const FluidState &)
 Calculate non-wetting liquid phase saturation given that the rest of the fluid state has been initialized.
 
template<class Evaluation = Scalar>
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatSn (const Params &, const Evaluation &)
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static std::enable_if< Traits::numPhases==3, Evaluation >::type Sg (const Params &, const FluidState &)
 Calculate gas phase saturation given that the rest of the fluid state has been initialized.
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation krw (const Params &, const FluidState &fs)
 The relative permability of the wetting phase.
 
template<class Evaluation = Scalar>
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatKrw (const Params &, const Evaluation &Sw)
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation krn (const Params &, const FluidState &fs)
 The relative permability of the liquid non-wetting phase.
 
template<class Evaluation = Scalar>
static std::enable_if< Traits::numPhases==2, Evaluation >::type twoPhaseSatKrn (const Params &, const Evaluation &Sw)
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static std::enable_if< Traits::numPhases==3, Evaluation >::type krg (const Params &, const FluidState &fs)
 The relative permability of the gas phase.
 
template<class FluidState , class Evaluation = Scalar>
static std::enable_if< Traits::numPhases==3, Evaluation >::type pcgn (const Params &params, const FluidState &fs)
 The difference between the pressures of the gas and the non-wetting phase.
 

Static Public Attributes

static const int numPhases = Traits::numPhases
 The number of fluid phases.
 
static const bool implementsTwoPhaseApi = (numPhases == 2)
 Specify whether this material law implements the two-phase convenience API.
 
static const bool implementsTwoPhaseSatApi = (numPhases == 2)
 Specify whether this material law implements the two-phase convenience API which only depends on the phase saturations.
 
static const bool isSaturationDependent = true
 Specify whether the quantities defined by this material law are saturation dependent.
 
static const bool isPressureDependent = false
 Specify whether the quantities defined by this material law are dependent on the absolute pressure.
 
static const bool isTemperatureDependent = false
 Specify whether the quantities defined by this material law are temperature dependent.
 
static const bool isCompositionDependent = false
 Specify whether the quantities defined by this material law are dependent on the phase composition.
 

Detailed Description

template<class TraitsT, class ParamsT = LinearMaterialParams<TraitsT>>
class Opm::LinearMaterial< TraitsT, ParamsT >

Implements a linear saturation-capillary pressure relation.

Implements a linear saturation-capillary pressure relation for M-phase fluid systems.

See also
LinearMaterialParams

Member Function Documentation

◆ capillaryPressures()

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class ContainerT , class FluidState >
static void Opm::LinearMaterial< TraitsT, ParamsT >::capillaryPressures ( ContainerT &  values,
const Params &  params,
const FluidState &  state 
)
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
valuesContainer for the return values
paramsParameters
stateThe fluid state

◆ krg()

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static std::enable_if< Traits::numPhases==3, Evaluation >::type Opm::LinearMaterial< TraitsT, ParamsT >::krg ( const Params &  ,
const FluidState &  fs 
)
inlinestatic

The relative permability of the gas phase.

This method is only available for at least three fluid phases

◆ pcgn()

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class FluidState , class Evaluation = Scalar>
static std::enable_if< Traits::numPhases==3, Evaluation >::type Opm::LinearMaterial< TraitsT, ParamsT >::pcgn ( const Params &  params,
const FluidState &  fs 
)
inlinestatic

The difference between the pressures of the gas and the non-wetting phase.

This method is only available for at least three fluid phases

◆ Sg()

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static std::enable_if< Traits::numPhases==3, Evaluation >::type Opm::LinearMaterial< TraitsT, ParamsT >::Sg ( const Params &  ,
const FluidState &   
)
inlinestatic

Calculate gas phase saturation given that the rest of the fluid state has been initialized.

This method is only available for at least three fluid phases


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