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::ParkerVanGen3P< ScalarT, ParamsT > Class Template Reference

Implementation of van Genuchten's capillary pressure <-> saturation relation. This class bundles the "raw" curves as static members and doesn't concern itself converting absolute to effective saturations and vince versa. More...

#include <dumux/material/fluidmatrixinteractions/3p/parkervangen3p.hh>

Description

template<class ScalarT, class ParamsT = ParkerVanGen3PParams<ScalarT>>
class Dumux::ParkerVanGen3P< ScalarT, ParamsT >

Implementation of van Genuchten's capillary pressure <-> saturation relation. This class bundles the "raw" curves as static members and doesn't concern itself converting absolute to effective saturations and vince versa.

See also
VanGenuchten, VanGenuchtenThreephase

Public Types

using Params = ParamsT
 
using Scalar = typename Params::Scalar
 

Static Public Member Functions

static Scalar pc (const Params &params, const Scalar sw)
 The capillary pressure-saturation curve. More...
 
static Scalar pcgw (const Params &params, const Scalar swe)
 The capillary pressure-saturation curve for the gas and wetting phase. More...
 
static Scalar pcnw (const Params &params, const Scalar swe)
 The capillary pressure-saturation curve for the non-wettigng and wetting phase. More...
 
static Scalar pcgn (const Params &params, const Scalar ste)
 The capillary pressure-saturation curve for the gas and non-wetting phase. More...
 
static Scalar pcAlpha (const Params &params, Scalar sne)
 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_dswe (const Params &params, const Scalar swe)
 Returns the partial derivative of the capillary pressure to the effective saturation. More...
 
static Scalar dpcgw_dswe (const Params &params, const Scalar seRegu)
 Returns the partial derivative of the capillary pressure to the effective saturation. More...
 
static Scalar dpcnw_dswe (const Params &params, const Scalar seRegu)
 Returns the partial derivative of the capillary pressure to the effective saturation. More...
 
static Scalar dpcgn_dste (const Params &params, const Scalar seRegu)
 Returns the partial derivative of the capillary pressure to the effective saturation. More...
 
static Scalar dswe_dpc (const Params &params, const Scalar pc)
 Returns the partial derivative of the effective saturation to the capillary pressure. More...
 
static Scalar krw (const Params &params, const Scalar swe)
 The relative permeability for the wetting phase of the medium implied by van Genuchten's parameterization. More...
 
static Scalar krn (const Params &params, const Scalar swe, const Scalar sn, const Scalar ste)
 The relative permeability for the non-wetting phase after the Model of Parker et al. (1987). More...
 
static Scalar krg (const Params &params, const Scalar ste)
 The relative permeability for the non-wetting phase of the medium implied by van Genuchten's parameterization. More...
 
static Scalar dkrg_dste (const Params &params, Scalar ste)
 The derivative of the relative permeability for the gas phase in regard to the total liquid saturation of the medium as implied by the van Genuchten parameterization. More...
 
static Scalar kr (const Params &params, const int phaseIdx, const Scalar swe, const Scalar sn, const Scalar ste)
 The relative permeability for a phase. More...
 
static Scalar bulkDensTimesAdsorpCoeff (const Params &params)
 the basis for calculating adsorbed NAPL in storage term More...
 

Member Typedef Documentation

◆ Params

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
using Dumux::ParkerVanGen3P< ScalarT, ParamsT >::Params = ParamsT

◆ Scalar

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
using Dumux::ParkerVanGen3P< ScalarT, ParamsT >::Scalar = typename Params::Scalar

Member Function Documentation

◆ bulkDensTimesAdsorpCoeff()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::bulkDensTimesAdsorpCoeff ( const Params params)
inlinestatic

the basis for calculating adsorbed NAPL in storage term

Parameters
paramsArray of parameters

◆ dkrg_dste()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::dkrg_dste ( const Params params,
Scalar  ste 
)
inlinestatic

The derivative of the relative permeability for the gas phase in regard to the total liquid saturation of the medium as implied by the van Genuchten parameterization.

Parameters
steThe mobile total liquid saturation.
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.

◆ dpc_dswe()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::dpc_dswe ( const Params params,
const Scalar  swe 
)
inlinestatic

Returns the partial derivative of the capillary pressure to the effective saturation.

Parameters
paramsArray of parameters
sweEffective wetting liquid saturation

◆ dpcgn_dste()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::dpcgn_dste ( const Params params,
const Scalar  seRegu 
)
inlinestatic

Returns the partial derivative of the capillary pressure to the effective saturation.

Parameters
paramsArray of parameters
seReguEffective wetting phase saturation for regularization

◆ dpcgw_dswe()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::dpcgw_dswe ( const Params params,
const Scalar  seRegu 
)
inlinestatic

Returns the partial derivative of the capillary pressure to the effective saturation.

Parameters
paramsArray of parameters
seReguEffective wetting phase saturation for regularization

◆ dpcnw_dswe()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::dpcnw_dswe ( const Params params,
const Scalar  seRegu 
)
inlinestatic

Returns the partial derivative of the capillary pressure to the effective saturation.

Parameters
paramsArray of parameters
seReguEffective wetting phase saturation for regularization

◆ dswe_dpc()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::dswe_dpc ( const Params params,
const Scalar  pc 
)
inlinestatic

Returns the partial derivative of the effective saturation to the capillary pressure.

Parameters
paramsArray of parameters
pcCapillary pressure in \(\mathrm{[Pa]}\)

◆ kr()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::kr ( const Params params,
const int  phaseIdx,
const Scalar  swe,
const Scalar  sn,
const Scalar  ste 
)
inlinestatic

The relative permeability for a phase.

Parameters
paramsArray of parameters.
phaseIdxIndicator, The saturation of all phases.
sweEffective wetting phase saturation
snAbsolute non-wetting liquid saturation
steEffective total liquid (wetting + non-wetting) saturation

◆ krg()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::krg ( const Params params,
const Scalar  ste 
)
inlinestatic

The relative permeability for the non-wetting phase of the medium implied by van Genuchten's parameterization.

The permeability of gas in a 3p system equals the standard 2p description. (see p61. in "Comparison of the Three-Phase Oil Relative Permeability Models" MOJDEH DELSHAD and GARY A. POPE, Transport in Porous Media 4 (1989), 59-83.) [17]

Parameters
paramsArray of parameters.
steEffective total liquid (wetting + non-wetting) saturation

◆ krn()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::krn ( const Params params,
const Scalar  swe,
const Scalar  sn,
const Scalar  ste 
)
inlinestatic

The relative permeability for the non-wetting phase after the Model of Parker et al. (1987).

See model 7 in "Comparison of the Three-Phase Oil Relative Permeability Models" MOJDEH DELSHAD and GARY A. POPE, Transport in Porous Media 4 (1989), 59-83 [17]
or more comprehensive in "Estimation of primary drainage three-phase relative permeability for organic liquid transport in the vadose zone", Leonardo I. Oliveira, Avery H. Demond, Journal of Contaminant Hydrology 66 (2003), 261-285 [46]

Parameters
paramsArray of parameters.
sweEffective wetting phase saturation
snAbsolute non-wetting liquid saturation
steEffective total liquid (wetting + non-wetting) saturation

◆ krw()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::krw ( const Params params,
const Scalar  swe 
)
inlinestatic

The relative permeability for the wetting phase of the medium implied by van Genuchten's parameterization.

The permeability of water in a 3p system equals the standard 2p description. (see p61. in "Comparison of the Three-Phase Oil Relative Permeability Models" MOJDEH DELSHAD and GARY A. POPE, Transport in Porous Media 4 (1989), 59-83.) [17]

Parameters
paramsArray of parameters.
sweEffective wetting phase saturation

◆ pc()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::pc ( const Params params,
const Scalar  sw 
)
inlinestatic

The capillary pressure-saturation curve.

Parameters
paramsArray of parameters
swwetting phase saturation

◆ pcAlpha()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::pcAlpha ( const Params params,
Scalar  sne 
)
inlinestatic

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

Parameters
paramsArray of parameters
sneNon-wetting liquid saturation

◆ pcgn()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::pcgn ( const Params params,
const Scalar  ste 
)
inlinestatic

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

Parameters
paramsArray of parameters
steEffective total liquid (wetting + non-wetting) saturation

◆ pcgw()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::pcgw ( const Params params,
const Scalar  swe 
)
inlinestatic

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

Parameters
paramsArray of parameters
sweEffective wetting phase saturation

◆ pcnw()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::pcnw ( const Params params,
const Scalar  swe 
)
inlinestatic

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

Parameters
paramsArray of parameters
sweEffective wetting phase saturation

◆ sw()

template<class ScalarT , class ParamsT = ParkerVanGen3PParams<ScalarT>>
static Scalar Dumux::ParkerVanGen3P< ScalarT, ParamsT >::sw ( const Params params,
const Scalar  pc 
)
inlinestatic

The saturation-capillary pressure curve.

Parameters
paramsArray of parameters
pcCapillary pressure in \(\mathrm{[Pa]}\)

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