3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Public Types | Static Public Member Functions | List of all members
Dumux::EffToAbsLaw< EffLawT, AbsParamsT > Class Template Reference

This material law takes a material law defined for effective saturations and converts it to a material law defined on absolute saturations. More...

#include <dumux/material/fluidmatrixinteractions/2p/efftoabslaw.hh>

Description

template<class EffLawT, class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
class Dumux::EffToAbsLaw< EffLawT, AbsParamsT >

This material law takes a material law defined for effective saturations and converts it to a material law defined on absolute saturations.

The idea: "material laws" (like VanGenuchten or BrooksCorey) are defined for effective saturations. The numeric calculations however are performed with absolute saturations. The EffToAbsLaw class gets the "material laws" actually used as well as the corresponding parameter container as template arguments.

Subsequently, the desired function (pc, sw... ) of the actually used "material laws" are called but with the saturations already converted from absolute to effective.

This approach makes sure that in the "material laws" only effective saturations are considered, which makes sense, as these laws only deal with effective saturations. This also allows for changing the calculation of the effective saturations easily, as this is subject of discussion / may be problem specific.

Additionally, handing over effective saturations to the "material laws" in stead of them calculating effective saturations prevents accidently "converting twice".

This boils down to:

Public Types

using Params = AbsParamsT
 
using Scalar = typename EffLaw::Scalar
 
using Params = AbsParamsT
 
using Scalar = typename EffLaw::Scalar
 

Static Public Member Functions

static Scalar pc (const Params &params, Scalar sw)
 The capillary pressure-saturation curve. More...
 
static Scalar sw (const Params &params, Scalar pc)
 The saturation-capillary pressure curve. More...
 
static Scalar endPointPc (const Params &params)
 The capillary pressure at Swe = 1.0 also called end point capillary pressure. More...
 
static Scalar dpc_dsw (const Params &params, Scalar sw)
 Returns the partial derivative of the capillary pressure w.r.t the absolute saturation. More...
 
static Scalar dsw_dpc (const Params &params, Scalar pc)
 Returns the partial derivative of the absolute saturation w.r.t. the capillary pressure. More...
 
static Scalar krw (const Params &params, Scalar sw)
 The relative permeability for the wetting phase. More...
 
static Scalar dkrw_dsw (const Params &params, Scalar sw)
 Returns the partial derivative of the relative permeability of the wetting phase with respect to the wetting saturation. More...
 
static Scalar krn (const Params &params, Scalar sw)
 The relative permeability for the non-wetting phase. More...
 
static Scalar dkrn_dsw (const Params &params, Scalar sw)
 Returns the partial derivative of the relative permeability of the non-wetting phase with respect to the wetting saturation. More...
 
static Scalar swToSwe (const Params &params, Scalar sw)
 Convert an absolute wetting saturation to an effective one. More...
 
static Scalar snToSne (const Params &params, Scalar sn)
 Convert an absolute non-wetting saturation to an effective one. More...
 
static Scalar sweToSw (const Params &params, Scalar swe)
 Convert an effective wetting saturation to an absolute one. More...
 
static Scalar sweToSw_ (const Params &params, Scalar swe)
 
static Scalar sneToSn (const Params &params, Scalar sne)
 Convert an effective non-wetting saturation to an absolute one. More...
 
static Scalar dswe_dsw (const Params &params)
 Derivative of the effective saturation w.r.t. the absolute saturation. More...
 
static Scalar dswe_dsw_ (const Params &params)
 
static Scalar dsw_dswe (const Params &params)
 Derivative of the absolute saturation w.r.t. the effective saturation. More...
 
static Scalar dsw_dswe_ (const Params &params)
 
static Scalar pc (const Params &params, const Scalar sw)
 The capillary pressure-saturation curve. More...
 
static Scalar pcgw (const Params &params, const Scalar sw)
 The capillary pressure-saturation curve for the gas and wetting phase. More...
 
static Scalar pcnw (const Params &params, const Scalar sw)
 The capillary pressure-saturation curve the non-wetting and wetting phase. More...
 
static Scalar pcgn (const Params &params, const Scalar st)
 The capillary pressure-saturation curve for the gas and non-wetting phase. More...
 
static Scalar pcAlpha (const Params &params, const Scalar sn)
 This function ensures a continuous transition from 2 to 3 phases and vice versa. More...
 
static Scalar sw (const Params &params, const Scalar pc)
 The saturation-capillary pressure curve. More...
 
static Scalar dpc_dsw (const Params &params, const Scalar sw)
 Returns the partial derivative of the capillary pressure w.r.t the absolute saturation. In this case the chain rule needs to be applied: \(\mathrm{ p_c = p_c( \overline{S}_w (S_w)) \rightarrow p_c ^\prime = \frac{\partial p_c}{\partial \overline S_w} \frac{\partial \overline{S}_w}{\partial S_w} }\). More...
 
static Scalar dsw_dpc (const Params &params, const Scalar pc)
 Returns the partial derivative of the absolute saturation w.r.t. the capillary pressure. In this case the chain rule needs to be applied: \(\mathrm{ S_w = S_w(\overline{S}_w (p_c) ) \rightarrow S_w^\prime = \frac{\partial S_w}{\partial \overline{S}_w} \frac{\partial \overline{S}_w}{\partial p_c} }\). More...
 
static Scalar krw (const Params &params, Scalar sw, const Scalar sn)
 The relative permeability for the wetting phase. More...
 
static Scalar krn (const Params &params, const Scalar sw, const Scalar sn)
 The relative permeability for the non-wetting phase. More...
 
static Scalar krg (const Params &params, const Scalar sw, const Scalar sn)
 The relative permeability for the gas phase. More...
 
static Scalar kr (const Params &params, const int phaseIdx, const Scalar sw, const Scalar sn, const Scalar sg)
 The relative permeability for a phase. More...
 
static Scalar bulkDensTimesAdsorpCoeff (const Params &params)
 the basis for calculating adsorbed NAPL in storage term More...
 
static Scalar swToSwe (const Params &params, const Scalar sw)
 Convert an absolute wetting saturation to an effective one. More...
 
static Scalar snToSne (const Params &params, const Scalar sn)
 Convert an absolute non-wetting saturation to an effective one. More...
 
static Scalar stToSte (const Params &params, const Scalar st)
 Convert an absolute total liquid saturation to an effective one. More...
 
static Scalar sgToSge (const Params &params, Scalar sg)
 Convert an absolute gas saturation to an effective one. More...
 
static Scalar sweToSw_ (const Params &params, Scalar swe)
 Convert an effective wetting saturation to an absolute one. More...
 
static Scalar sneToSn_ (const Params &params, Scalar sne)
 
static Scalar sgeToSg_ (const Params &params, Scalar sge)
 
static Scalar dswe_dsw_ (const Params &params)
 Derivative of the effective saturation w.r.t. the absolute saturation. More...
 
static Scalar dsw_dswe_ (const Params &params)
 Derivative of the absolute saturation w.r.t. the effective saturation. More...
 

Member Typedef Documentation

◆ Params [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
using Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::Params = AbsParamsT

◆ Params [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
using Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::Params = AbsParamsT

◆ Scalar [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
using Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::Scalar = typename EffLaw::Scalar

◆ Scalar [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
using Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::Scalar = typename EffLaw::Scalar

Member Function Documentation

◆ bulkDensTimesAdsorpCoeff()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::bulkDensTimesAdsorpCoeff ( const Params params)
inlinestatic

the basis for calculating adsorbed NAPL in storage term

Parameters
paramsArray of parameters

◆ dkrn_dsw()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dkrn_dsw ( const Params params,
Scalar  sw 
)
inlinestatic

Returns the partial derivative of the relative permeability of the non-wetting phase with respect to the wetting saturation.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law.

◆ dkrw_dsw()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dkrw_dsw ( const Params params,
Scalar  sw 
)
inlinestatic

Returns the partial derivative of the relative permeability of the wetting phase with respect to the wetting saturation.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law.

◆ dpc_dsw() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dpc_dsw ( const Params params,
const Scalar  sw 
)
inlinestatic

Returns the partial derivative of the capillary pressure w.r.t the absolute saturation. In this case the chain rule needs to be applied: \(\mathrm{ p_c = p_c( \overline{S}_w (S_w)) \rightarrow p_c ^\prime = \frac{\partial p_c}{\partial \overline S_w} \frac{\partial \overline{S}_w}{\partial S_w} }\).

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Partial derivative of \(\mathrm{[p_c]}\) w.r.t. effective saturation according to EffLaw e.g. Brooks & Corey, van Genuchten, linear... .

◆ dpc_dsw() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dpc_dsw ( const Params params,
Scalar  sw 
)
inlinestatic

Returns the partial derivative of the capillary pressure w.r.t the absolute saturation.

In this case the chain rule needs to be applied: \(\mathrm{ p_c = p_c( \overline{S}_w (S_w)) \rightarrow p_c ^\prime = \frac{\partial p_c}{\partial \overline S_w} \frac{\partial \overline{S}_w}{\partial S_w} }\)

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Partial derivative of \(\mathrm{[p_c]}\) w.r.t. effective saturation according to EffLaw e.g. Brooks & Corey, van Genuchten, linear... .

◆ dsw_dpc() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dsw_dpc ( const Params params,
const Scalar  pc 
)
inlinestatic

Returns the partial derivative of the absolute saturation w.r.t. the capillary pressure. In this case the chain rule needs to be applied: \(\mathrm{ S_w = S_w(\overline{S}_w (p_c) ) \rightarrow S_w^\prime = \frac{\partial S_w}{\partial \overline{S}_w} \frac{\partial \overline{S}_w}{\partial p_c} }\).

Parameters
pcCapillary pressure \(\mathrm{[p_c]}\) in \(\mathrm{[Pa]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Partial derivative of effective saturation w.r.t. \(\mathrm{[p_c]}\) according to EffLaw e.g. Brooks & Corey, van Genuchten, linear... .

◆ dsw_dpc() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dsw_dpc ( const Params params,
Scalar  pc 
)
inlinestatic

Returns the partial derivative of the absolute saturation w.r.t. the capillary pressure.

In this case the chain rule needs to be applied: \(\mathrm{ S_w = S_w(\overline{S}_w (p_c) ) \rightarrow S_w^\prime = \frac{\partial S_w}{\partial \overline{S}_w} \frac{\partial \overline{S}_w}{\partial p_c} }\)

Parameters
pcCapillary pressure \(\mathrm{[p_c]}\) in \(\mathrm{[Pa]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Partial derivative of effective saturation w.r.t. \(\mathrm{[p_c]}\) according to EffLaw e.g. Brooks & Corey, van Genuchten, linear... .

◆ dsw_dswe()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dsw_dswe ( const Params params)
inlinestatic

Derivative of the absolute saturation w.r.t. the effective saturation.

Parameters
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Derivative of the absolute saturation w.r.t. the effective saturation.

◆ dsw_dswe_() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dsw_dswe_ ( const Params params)
inlinestatic

◆ dsw_dswe_() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dsw_dswe_ ( const Params params)
inlinestatic

Derivative of the absolute saturation w.r.t. the effective saturation.

Parameters
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Derivative of the absolute saturation w.r.t. the effective saturation.

◆ dswe_dsw()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dswe_dsw ( const Params params)
inlinestatic

Derivative of the effective saturation w.r.t. the absolute saturation.

Parameters
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Derivative of the effective saturation w.r.t. the absolute saturation.

◆ dswe_dsw_() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dswe_dsw_ ( const Params params)
inlinestatic

◆ dswe_dsw_() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::dswe_dsw_ ( const Params params)
inlinestatic

Derivative of the effective saturation w.r.t. the absolute saturation.

Parameters
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Derivative of the effective saturation w.r.t. the absolute saturation.

◆ endPointPc()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::endPointPc ( const Params params)
inlinestatic

The capillary pressure at Swe = 1.0 also called end point capillary pressure.

Parameters
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.

◆ kr()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::kr ( const Params params,
const int  phaseIdx,
const Scalar  sw,
const Scalar  sn,
const Scalar  sg 
)
inlinestatic

The relative permeability for a phase.

Parameters
swWetting liquid saturation
sgGas saturation
snNon-wetting liquid saturation
paramsArray of parameters.
phaseIdxindicator, The saturation of all phases.

◆ krg()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::krg ( const Params params,
const Scalar  sw,
const Scalar  sn 
)
inlinestatic

The relative permeability for the gas phase.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[{S}_w]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
snAbsolute saturation of the non-wetting phase \(\mathrm{[{S}_n]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Relative permeability of the non-wetting phase calculated as implied by EffLaw e.g. Brooks & Corey, van Genuchten, linear... .

◆ krn() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::krn ( const Params params,
const Scalar  sw,
const Scalar  sn 
)
inlinestatic

The relative permeability for the non-wetting phase.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[{S}_w]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
snAbsolute saturation of the non-wetting phase \(\mathrm{[{S}_n]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Relative permeability of the non-wetting phase calculated as implied by EffLaw e.g. Brooks & Corey, van Genuchten, linear... .

◆ krn() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::krn ( const Params params,
Scalar  sw 
)
inlinestatic

The relative permeability for the non-wetting phase.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[{S}_w]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Relative permeability of the non-wetting phase calculated as implied by EffLaw e.g. Brooks & Corey, van Genuchten, linear... .

◆ krw() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::krw ( const Params params,
Scalar  sw 
)
inlinestatic

The relative permeability for the wetting phase.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Relative permeability of the wetting phase calculated as implied by EffLaw e.g. Brooks & Corey, van Genuchten, linear... .

◆ krw() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::krw ( const Params params,
Scalar  sw,
const Scalar  sn 
)
inlinestatic

The relative permeability for the wetting phase.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
snAbsolute saturation of the non-wetting phase \(\mathrm{[{S}_n]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Relative permeability of the wetting phase calculated as implied by EffLaw e.g. Brooks & Corey, van Genuchten, linear... .

◆ pc() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::pc ( const Params params,
const Scalar  sw 
)
inlinestatic

The capillary pressure-saturation curve.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Capillary pressure calculated by specific constitutive relation (EffLaw e.g. Brooks & Corey, van Genuchten, linear...)

◆ pc() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::pc ( const Params params,
Scalar  sw 
)
inlinestatic

The capillary pressure-saturation curve.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[\overline{S}_w]}\). It is converted to effective saturation and then handed over to the material law actually used for calculation.
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Capillary pressure calculated by specific constitutive relation (EffLaw e.g. Brooks & Corey, van Genuchten, linear...)

◆ pcAlpha()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::pcAlpha ( const Params params,
const Scalar  sn 
)
inlinestatic

This function ensures a continuous transition from 2 to 3 phases and vice versa.

Parameters
paramsArray of parameters
snNon-wetting liquid saturation

◆ pcgn()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::pcgn ( const Params params,
const Scalar  st 
)
inlinestatic

The capillary pressure-saturation curve for the gas and non-wetting phase.

Parameters
paramsArray of parameters
stsum of wetting (liquid) phase saturations

◆ pcgw()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::pcgw ( const Params params,
const Scalar  sw 
)
inlinestatic

The capillary pressure-saturation curve for the gas and wetting phase.

Parameters
paramsArray of parameters
swwetting phase saturation or sum of wetting phase saturations

◆ pcnw()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::pcnw ( const Params params,
const Scalar  sw 
)
inlinestatic

The capillary pressure-saturation curve the non-wetting and wetting phase.

Parameters
paramsArray of parameters
swwetting phase saturation or sum of wetting phase saturations

◆ sgeToSg_()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::sgeToSg_ ( const Params params,
Scalar  sge 
)
inlinestatic

◆ sgToSge()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::sgToSge ( const Params params,
Scalar  sg 
)
inlinestatic

Convert an absolute gas saturation to an effective one.

Parameters
sgAbsolute saturation of the gas phase \(\mathrm{[{S}_n]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Effective saturation of the non-wetting phase.

◆ sneToSn()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::sneToSn ( const Params params,
Scalar  sne 
)
inlinestatic

Convert an effective non-wetting saturation to an absolute one.

Parameters
sneEffective saturation of the non-wetting phase \(\mathrm{[{S}_n]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Absolute saturation of the non-wetting phase.

◆ sneToSn_()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::sneToSn_ ( const Params params,
Scalar  sne 
)
inlinestatic

◆ snToSne() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::snToSne ( const Params params,
const Scalar  sn 
)
inlinestatic

Convert an absolute non-wetting saturation to an effective one.

Parameters
snAbsolute saturation of the non-wetting phase \(\mathrm{[{S}_n]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Effective saturation of the non-wetting phase.

◆ snToSne() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::snToSne ( const Params params,
Scalar  sn 
)
inlinestatic

Convert an absolute non-wetting saturation to an effective one.

Parameters
snAbsolute saturation of the non-wetting phase \(\mathrm{[{S}_n]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Effective saturation of the non-wetting phase.

◆ stToSte()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::stToSte ( const Params params,
const Scalar  st 
)
inlinestatic

Convert an absolute total liquid saturation to an effective one.

Parameters
stAbsolute saturation of the total liquid phase (sw+sn) \(\mathrm{[{S}_n]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Effective saturation of the non-wetting phase.

◆ sw() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::sw ( const Params params,
const Scalar  pc 
)
inlinestatic

The saturation-capillary pressure curve.

Parameters
pcCapillary pressure \(\mathrm{[p_c]}\) in \(\mathrm{[Pa]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Absolute wetting phase saturation calculated as inverse of (EffLaw e.g. Brooks & Corey, van Genuchten, linear...) constitutive relation.
The absolute saturation of the wetting phase \(\mathrm{[S_w]}\)

◆ sw() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::sw ( const Params params,
Scalar  pc 
)
inlinestatic

The saturation-capillary pressure curve.

Parameters
pcCapillary pressure \(\mathrm{[p_c]}\) in \(\mathrm{[Pa]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Absolute wetting phase saturation \(\mathrm{[S_w]}\) calculated as inverse of (EffLaw e.g. Brooks & Corey, van Genuchten, linear...) constitutive relation.

◆ sweToSw()

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::sweToSw ( const Params params,
Scalar  swe 
)
inlinestatic

Convert an effective wetting saturation to an absolute one.

Parameters
sweEffective saturation of the non-wetting phase \(\mathrm{[\overline{S}_n]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Absolute saturation of the non-wetting phase.

◆ sweToSw_() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::sweToSw_ ( const Params params,
Scalar  swe 
)
inlinestatic

◆ sweToSw_() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::sweToSw_ ( const Params params,
Scalar  swe 
)
inlinestatic

Convert an effective wetting saturation to an absolute one.

Parameters
sweEffective saturation of the non-wetting phase \(\mathrm{[\overline{S}_n]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Absolute saturation of the non-wetting phase.

◆ swToSwe() [1/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::swToSwe ( const Params params,
const Scalar  sw 
)
inlinestatic

Convert an absolute wetting saturation to an effective one.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[{S}_w]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Effective saturation of the wetting phase.

◆ swToSwe() [2/2]

template<class EffLawT , class AbsParamsT = EffToAbsLawParams<typename EffLawT::Params>>
static Scalar Dumux::EffToAbsLaw< EffLawT, AbsParamsT >::swToSwe ( const Params params,
Scalar  sw 
)
inlinestatic

Convert an absolute wetting saturation to an effective one.

Parameters
swAbsolute saturation of the wetting phase \(\mathrm{[{S}_w]}\).
paramsA container object that is populated with the appropriate coefficients for the respective law. Therefore, in the (problem specific) spatialParameters first, the material law is chosen, and then the params container is constructed accordingly. Afterwards the values are set there, too.
Returns
Effective saturation of the wetting phase.

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