Opm::RegularizedVanGenuchten< TraitsT, ParamsT > Class Template Reference

Implementation of the regularized van Genuchten's capillary pressure / relative permeability <-> saturation relation. More...

#include <RegularizedVanGenuchten.hpp>

Inheritance diagram for Opm::RegularizedVanGenuchten< TraitsT, ParamsT >:
Inheritance graph

Public Types

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

Static Public Member Functions

template<class Container , class FluidState >
static void capillaryPressures (Container &values, const Params &params, const FluidState &fs)
 Calculate the pressure difference of the phases in the most generic way. More...
 
template<class Container , class FluidState >
static void saturations (Container &values, const Params &params, const FluidState &fs)
 Calculate the saturations of the phases starting from their pressure differences. More...
 
template<class Container , class FluidState >
static void relativePermeabilities (Container &values, const Params &params, const FluidState &fs)
 Returns the relative permeabilities of the phases dependening on the phase saturations. More...
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation pcnw (const Params &params, const FluidState &fs)
 A regularized van Genuchten capillary pressure-saturation curve. More...
 
template<class Evaluation >
static Evaluation twoPhaseSatPcnw (const Params &params, const Evaluation &Sw)
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Sw (const Params &params, const FluidState &fs)
 A regularized van Genuchten saturation-capillary pressure curve. More...
 
template<class Evaluation >
static Evaluation twoPhaseSatSw (const Params &params, const Evaluation &pC)
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Sn (const Params &params, const FluidState &fs)
 Calculate the non-wetting phase saturations depending on the phase pressures. More...
 
template<class Evaluation >
static Evaluation twoPhaseSatSn (const Params &params, const Evaluation &pc)
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation krw (const Params &params, const FluidState &fs)
 Regularized version of the relative permeability for the wetting phase of the medium implied by the van Genuchten parameterization. More...
 
template<class Evaluation >
static Evaluation twoPhaseSatKrw (const Params &params, const Evaluation &Sw)
 
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation krn (const Params &params, const FluidState &fs)
 Regularized version of the relative permeability for the non-wetting phase of the medium implied by the van Genuchten parameterization. More...
 
template<class Evaluation >
static Evaluation twoPhaseSatKrn (const Params &params, const Evaluation &Sw)
 

Static Public Attributes

static const int numPhases = Traits::numPhases
 The number of fluid phases. More...
 
static const bool implementsTwoPhaseApi = true
 
static const bool implementsTwoPhaseSatApi = true
 
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 = RegularizedVanGenuchtenParams<TraitsT>>
class Opm::RegularizedVanGenuchten< TraitsT, ParamsT >

Implementation of the regularized van Genuchten's capillary pressure / relative permeability <-> saturation relation.

This class bundles the "raw" curves as static members and doesn't concern itself converting absolute to effective saturations and vice versa.

In order to avoid very steep gradients the marginal values are "regularized". This means that in stead of following the curve of the material law in these regions, some linear approximation is used. Doing this is not worse than following the material law. E.g. for very low wetting phase values the material laws predict infinite values for $p_c$ which is completely unphysical. In case of very high wetting phase saturations the difference between regularized and "pure" material law is not big.

Regularizing has the additional benefit of being numerically friendly: Newton's method does not like infinite gradients.

The implementation is accomplished as follows:

  • check whether we are in the range of regularization
    • yes: use the regularization
    • no: forward to the standard material law.

An example of the regularization of the capillary pressure curve is shown below:

See also
VanGenuchten

Member Typedef Documentation

◆ Params

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
typedef ParamsT Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::Params

◆ Scalar

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
typedef Traits::Scalar Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::Scalar

◆ Traits

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
typedef TraitsT Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::Traits

Member Function Documentation

◆ capillaryPressures()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class Container , class FluidState >
static void Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::capillaryPressures ( Container &  values,
const Params params,
const FluidState &  fs 
)
inlinestatic

Calculate the pressure difference of the phases in the most generic way.

◆ krn()

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

Regularized version of the relative permeability for the non-wetting phase of the medium implied by the van Genuchten parameterization.

regularized part:

  • below $ \overline S_w =0$: set relative permeability to zero
  • above $ \overline S_w =1$: set relative permeability to one
  • for $ 0 \leq \overline S_w \leq 0.05 $: use a spline as interpolation
Parameters
paramsThe parameter object expressing the coefficients required by the van Genuchten law.
fsThe fluid state for which the derivative ought to be calculated

References Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::Sw(), and Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatKrn().

◆ krw()

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

Regularized version of the relative permeability for the wetting phase of the medium implied by the van Genuchten parameterization.

regularized part:

  • below $ \overline S_w =0$: set relative permeability to zero
  • above $ \overline S_w =1$: set relative permeability to one
  • between $ 0.95 \leq \overline S_w \leq 1$: use a spline as interpolation

For not-regularized part:

Parameters
paramsThe parameter object expressing the coefficients required by the van Genuchten law.
fsThe fluid state for which the relative permeability ought to be calculated

References Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::Sw(), and Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatKrw().

◆ pcnw()

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

A regularized van Genuchten capillary pressure-saturation curve.

regularized part:

  • low saturation: extend the $p_c(S_w)$ curve with the slope at the regularization point (i.e. no kink).
  • high saturation: connect the high regularization point with $ \overline S_w =1$ by a straight line (yes, there is a kink :-( ).

For not-regularized part:

References Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::Sw(), and Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatPcnw().

◆ relativePermeabilities()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class Container , class FluidState >
static void Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::relativePermeabilities ( Container &  values,
const Params params,
const FluidState &  fs 
)
inlinestatic

Returns the relative permeabilities of the phases dependening on the phase saturations.

◆ saturations()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class Container , class FluidState >
static void Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::saturations ( Container &  values,
const Params params,
const FluidState &  fs 
)
inlinestatic

Calculate the saturations of the phases starting from their pressure differences.

◆ Sn()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::Sn ( const Params params,
const FluidState &  fs 
)
inlinestatic

Calculate the non-wetting phase saturations depending on the phase pressures.

◆ Sw()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class FluidState , class Evaluation = typename FluidState::Scalar>
static Evaluation Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::Sw ( const Params params,
const FluidState &  fs 
)
inlinestatic

A regularized van Genuchten saturation-capillary pressure curve.

regularized part:

  • low saturation: extend the $p_c(S_w)$ curve with the slope at the regularization point (i.e. no kink).
  • high saturation: connect the high regularization point with $ \overline S_w =1$ by a straight line (yes, there is a kink :-( ).

The according quantities are obtained by exploiting theorem of intersecting lines.

For not-regularized part:

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

\[
S_w = {p_C}^{-1} = ((\alpha p_C)^n + 1)^{-m}
\]

Parameters
paramsThe parameter object expressing the coefficients required by the van Genuchten law.
fsThe fluid state containing valid phase pressures

References Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatSw().

Referenced by Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::krn(), Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::krw(), Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::pcnw(), Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatKrn(), Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatKrw(), Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatPcnw(), and Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatSw().

◆ twoPhaseSatKrn()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class Evaluation >
static Evaluation Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatKrn ( const Params params,
const Evaluation &  Sw 
)
inlinestatic

◆ twoPhaseSatKrw()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class Evaluation >
static Evaluation Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatKrw ( const Params params,
const Evaluation &  Sw 
)
inlinestatic

◆ twoPhaseSatPcnw()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class Evaluation >
static Evaluation Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatPcnw ( const Params params,
const Evaluation &  Sw 
)
inlinestatic

◆ twoPhaseSatSn()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class Evaluation >
static Evaluation Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatSn ( const Params params,
const Evaluation &  pc 
)
inlinestatic

◆ twoPhaseSatSw()

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
template<class Evaluation >
static Evaluation Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::twoPhaseSatSw ( const Params params,
const Evaluation &  pC 
)
inlinestatic

Member Data Documentation

◆ implementsTwoPhaseApi

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
const bool Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::implementsTwoPhaseApi = true
static

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

◆ implementsTwoPhaseSatApi

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
const bool Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::implementsTwoPhaseSatApi = true
static

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

◆ isCompositionDependent

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
const bool Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::isCompositionDependent = false
static

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

◆ isPressureDependent

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
const bool Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::isPressureDependent = false
static

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

◆ isSaturationDependent

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
const bool Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::isSaturationDependent = true
static

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

◆ isTemperatureDependent

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
const bool Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::isTemperatureDependent = false
static

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

◆ numPhases

template<class TraitsT , class ParamsT = RegularizedVanGenuchtenParams<TraitsT>>
const int Opm::RegularizedVanGenuchten< TraitsT, ParamsT >::numPhases = Traits::numPhases
static

The number of fluid phases.


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