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 >:
Inheritance graph

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. More...
 
template<class ContainerT , class FluidState >
static void saturations (ContainerT &, const Params &, const FluidState &)
 The inverse of the capillary pressure. More...
 
template<class ContainerT , class FluidState >
static void relativePermeabilities (ContainerT &values, const Params &, const FluidState &state)
 The relative permeability of all phases. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation krw (const Params &, const FluidState &fs)
 The relative permability of the wetting phase. More...
 
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. More...
 
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. More...
 
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. More...
 

Static Public Attributes

static const int numPhases = Traits::numPhases
 The number of fluid phases. More...
 
static const bool implementsTwoPhaseApi = (numPhases == 2)
 
static const bool implementsTwoPhaseSatApi = (numPhases == 2)
 
static const bool isSaturationDependent = true
 
static const bool isPressureDependent = false
 
static const bool isTemperatureDependent = false
 
static const bool isCompositionDependent = false
 

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 Typedef Documentation

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
typedef ParamsT Opm::LinearMaterial< TraitsT, ParamsT >::Params
template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
typedef Traits::Scalar Opm::LinearMaterial< TraitsT, ParamsT >::Scalar
template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
typedef TraitsT Opm::LinearMaterial< TraitsT, ParamsT >::Traits

Member Function Documentation

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

References Valgrind::CheckDefined().

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

References Opm::LocalAd::max(), Opm::LocalAd::min(), and Opm::LinearMaterial< TraitsT, ParamsT >::Sg().

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Opm::LinearMaterial< TraitsT, ParamsT >::krn ( const Params ,
const FluidState &  fs 
)
inlinestatic

The relative permability of the liquid non-wetting phase.

References Opm::LocalAd::max(), Opm::LocalAd::min(), and Opm::LinearMaterial< TraitsT, ParamsT >::Sn().

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Opm::LinearMaterial< TraitsT, ParamsT >::krw ( const Params ,
const FluidState &  fs 
)
inlinestatic

The relative permability of the wetting phase.

References Opm::LocalAd::max(), Opm::LocalAd::min(), and Opm::LinearMaterial< TraitsT, ParamsT >::Sw().

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

References Valgrind::CheckDefined(), Opm::LinearMaterial< TraitsT, ParamsT >::Sg(), and Opm::LinearMaterial< TraitsT, ParamsT >::Sn().

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Opm::LinearMaterial< TraitsT, ParamsT >::pcnw ( const Params params,
const FluidState &  fs 
)
inlinestatic

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

References Valgrind::CheckDefined(), Opm::LinearMaterial< TraitsT, ParamsT >::Sn(), and Opm::LinearMaterial< TraitsT, ParamsT >::Sw().

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class ContainerT , class FluidState >
static void Opm::LinearMaterial< TraitsT, ParamsT >::relativePermeabilities ( ContainerT &  values,
const Params ,
const FluidState &  state 
)
inlinestatic

The relative permeability of all phases.

References Valgrind::CheckDefined(), Opm::LocalAd::max(), and Opm::LocalAd::min().

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class ContainerT , class FluidState >
static void Opm::LinearMaterial< TraitsT, ParamsT >::saturations ( ContainerT &  ,
const Params ,
const FluidState &   
)
inlinestatic

The inverse of the capillary pressure.

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

Referenced by Opm::LinearMaterial< TraitsT, ParamsT >::krg(), and Opm::LinearMaterial< TraitsT, ParamsT >::pcgn().

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Opm::LinearMaterial< TraitsT, ParamsT >::Sn ( const Params ,
const FluidState &   
)
inlinestatic

Calculate non-wetting liquid phase saturation given that the rest of the fluid state has been initialized.

Referenced by Opm::LinearMaterial< TraitsT, ParamsT >::krn(), Opm::LinearMaterial< TraitsT, ParamsT >::pcgn(), and Opm::LinearMaterial< TraitsT, ParamsT >::pcnw().

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Opm::LinearMaterial< TraitsT, ParamsT >::Sw ( const Params ,
const FluidState &   
)
inlinestatic

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

Referenced by Opm::LinearMaterial< TraitsT, ParamsT >::krw(), Opm::LinearMaterial< TraitsT, ParamsT >::pcnw(), and Opm::LinearMaterial< TraitsT, ParamsT >::twoPhaseSatPcnw().

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class Evaluation = Scalar>
static std::enable_if<Traits::numPhases == 2, Evaluation>::type Opm::LinearMaterial< TraitsT, ParamsT >::twoPhaseSatKrn ( const Params ,
const Evaluation &  Sw 
)
inlinestatic
template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class Evaluation = Scalar>
static std::enable_if<Traits::numPhases == 2, Evaluation>::type Opm::LinearMaterial< TraitsT, ParamsT >::twoPhaseSatKrw ( const Params ,
const Evaluation &  Sw 
)
inlinestatic
template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class Evaluation = Scalar>
static std::enable_if<Traits::numPhases == 2, Evaluation>::type Opm::LinearMaterial< TraitsT, ParamsT >::twoPhaseSatPcnw ( const Params params,
const Evaluation &  Sw 
)
inlinestatic
template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class Evaluation = Scalar>
static std::enable_if<Traits::numPhases == 2, Evaluation>::type Opm::LinearMaterial< TraitsT, ParamsT >::twoPhaseSatSn ( const Params ,
const Evaluation &   
)
inlinestatic
template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
template<class Evaluation = Scalar>
static std::enable_if<Traits::numPhases == 2, Evaluation>::type Opm::LinearMaterial< TraitsT, ParamsT >::twoPhaseSatSw ( const Params ,
const Evaluation &   
)
inlinestatic

Member Data Documentation

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
const bool Opm::LinearMaterial< TraitsT, ParamsT >::implementsTwoPhaseApi = (numPhases == 2)
static

Specify whether this material law implements the two-phase convenience API

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
const bool Opm::LinearMaterial< TraitsT, ParamsT >::implementsTwoPhaseSatApi = (numPhases == 2)
static

Specify whether this material law implements the two-phase convenience API which only depends on the phase saturations

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
const bool Opm::LinearMaterial< TraitsT, ParamsT >::isCompositionDependent = false
static

Specify whether the quantities defined by this material law are dependent on the phase composition

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
const bool Opm::LinearMaterial< TraitsT, ParamsT >::isPressureDependent = false
static

Specify whether the quantities defined by this material law are dependent on the absolute pressure

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
const bool Opm::LinearMaterial< TraitsT, ParamsT >::isSaturationDependent = true
static

Specify whether the quantities defined by this material law are saturation dependent

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
const bool Opm::LinearMaterial< TraitsT, ParamsT >::isTemperatureDependent = false
static

Specify whether the quantities defined by this material law are temperature dependent

template<class TraitsT , class ParamsT = LinearMaterialParams<TraitsT>>
const int Opm::LinearMaterial< TraitsT, ParamsT >::numPhases = Traits::numPhases
static

The number of fluid phases.


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