13#ifndef DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH 
   14#define DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH 
   51    template<
class Scalar>
 
   55        : alpha_(
alpha), n_(
n), m_(1.0 - 1.0/
n), l_(
l)
 
 
   58        Scalar 
alpha()
 const { 
return alpha_; }
 
   61        Scalar 
m()
 const { 
return m_; }
 
   62        void setM(Scalar 
m) { m_ = 
m; n_ = 1.0/(1.0 - 
m); }
 
   64        Scalar 
n()
 const{ 
return n_; }
 
   65        void setN(Scalar 
n){ n_ = 
n; m_ = 1.0 - 1.0/
n; }
 
   67        Scalar 
l()
 const { 
return l_; }
 
   72            return Dune::FloatCmp::eq(alpha_, p.alpha_, 1e-6)
 
   73                   && Dune::FloatCmp::eq(n_, p.n_, 1e-6)
 
   74                   && Dune::FloatCmp::eq(m_, p.m_, 1e-6)
 
   75                   && Dune::FloatCmp::eq(l_, p.l_, 1e-6);
 
 
   79        Scalar alpha_, n_, m_, l_;
 
 
   86    template<
class Scalar = 
double>
 
  109    template<
class Scalar>
 
  115        swe = clamp(
swe, 0.0, 1.0); 
 
  117        const Scalar 
pc = pow(pow(
swe, -1.0/params.
m()) - 1, 1.0/params.
n())/params.
alpha();
 
 
  135    template<
class Scalar>
 
  143        const Scalar sw = pow(pow(params.
alpha()*
pc, params.
n()) + 1, -params.
m());
 
 
  151    template<
class Scalar>
 
  171    template<
class Scalar>
 
  177        swe = clamp(
swe, 0.0, 1.0); 
 
  179        const Scalar powSwe = pow(
swe, -1/params.
m());
 
  180        return - 1.0/params.
alpha() * pow(powSwe - 1, 1.0/params.
n() - 1)/params.
n()
 
  181                                  * powSwe/
swe/params.
m();
 
 
  193    template<
class Scalar>
 
  201        const Scalar powAlphaPc = pow(params.
alpha()*
pc, params.
n());
 
  202        return -pow(powAlphaPc + 1, -params.
m()-1)*params.
m()*powAlphaPc/
pc*params.
n();
 
 
  215    template<
class Scalar>
 
  221        swe = clamp(
swe, 0.0, 1.0); 
 
  223        const Scalar r = 1.0 - pow(1.0 - pow(
swe, 1.0/params.
m()), params.
m());
 
  224        return pow(
swe, params.
l())*r*r;
 
 
  237    template<
class Scalar>
 
  243        swe = clamp(
swe, 0.0, 1.0); 
 
  245        const Scalar x = 1.0 - pow(
swe, 1.0/params.
m());
 
  246        const Scalar xToM = pow(x, params.
m());
 
  247        return (1.0 - xToM)*pow(
swe, params.
l()-1) * ( (1.0 - xToM)*params.
l() + 2*xToM*(1.0-x)/x );
 
 
  260    template<
class Scalar>
 
  266        swe = clamp(
swe, 0.0, 1.0); 
 
  268        return pow(1 - 
swe, params.
l()) * pow(1 - pow(
swe, 1.0/params.
m()), 2*params.
m());
 
 
  282    template<
class Scalar>
 
  288        swe = clamp(
swe, 0.0, 1.0); 
 
  290        const auto sne = 1.0 - 
swe;
 
  291        const auto x = 1.0 - pow(
swe, 1.0/params.
m());
 
  292        return -pow(sne, params.
l()-1.0) * pow(x, 2*params.
m() - 1.0) * ( params.
l()*x + 2.0*sne/
swe*(1.0 - x) );
 
 
 
  303template <
class Scalar>
 
  324        { 
return pcLowSwe_; }
 
 
  339        { 
return pcHighSwe_; }
 
 
  353        { 
return krnLowSwe_; }
 
 
  367        { 
return krwHighSwe_; }
 
 
 
  377    template<
class MaterialLaw>
 
  378    void init(
const MaterialLaw* m, 
const std::string& paramGroup)
 
  385        initPcParameters_(m, pcLowSwe_, pcHighSwe_);
 
  386        initKrParameters_(m, krnLowSwe_, krwHighSwe_);
 
 
  389    template<
class MaterialLaw, 
class BaseParams, 
class EffToAbsParams>
 
  390    void init(
const MaterialLaw* m, 
const BaseParams& bp, 
const EffToAbsParams& etap, 
const Params<Scalar>& p)
 
  397        initPcParameters_(m, pcLowSwe_, pcHighSwe_);
 
  398        initKrParameters_(m, krnLowSwe_, krwHighSwe_);
 
 
  406        return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
 
  407               && Dune::FloatCmp::eq(pcHighSwe_, o.pcHighSwe_, 1e-6)
 
  408               && Dune::FloatCmp::eq(krwHighSwe_, o.krwHighSwe_, 1e-6)
 
  409               && Dune::FloatCmp::eq(krnLowSwe_, o.krnLowSwe_, 1e-6);
 
 
  426        if (
swe <= pcLowSwe_)
 
  427            return pcLowSwePcValue_ + pcDerivativeLowSw_*(
swe - pcLowSwe_);
 
  430            return pcDerivativeHighSweEnd_*(
swe - 1.0);
 
  432        else if (
swe > pcHighSwe_)
 
  433            return pcSpline_.eval(
swe);
 
 
  444        if (
swe <= pcLowSwe_)
 
  445            return pcDerivativeLowSw_;
 
  448            return pcDerivativeHighSweEnd_;
 
  450        else if (
swe > pcHighSwe_)
 
  451            return pcSpline_.evalDerivative(
swe);
 
 
  464            if (pcHighSwe_ >= 1.0)
 
  467                return pc/pcDerivativeHighSweEnd_ + 1.0;
 
  471        else if (
pc <= pcHighSwePcValue_)
 
  472            return pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0, 
pc);
 
  474        else if (
pc >= pcLowSwePcValue_)
 
  475            return (
pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
 
 
  488            if (pcHighSwe_ >= 1.0)
 
  491                return 1.0/pcDerivativeHighSweEnd_;
 
  495        else if (
pc <= pcHighSwePcValue_)
 
  496            return 1.0/pcSpline_.evalDerivative(pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0, 
pc));
 
  498        else if (
pc >= pcLowSwePcValue_)
 
  499            return 1.0/pcDerivativeLowSw_;
 
 
  514        else if (
swe >= krwHighSwe_)
 
  515            return krwSpline_.eval(
swe);
 
 
  529        else if (
swe >= krwHighSwe_)
 
  530            return krwSpline_.evalDerivative(
swe);
 
 
  544        else if (
swe <= krnLowSwe_)
 
  545            return krnSpline_.eval(
swe);
 
 
  559        else if (
swe <= krnLowSwe_)
 
  560            return krnSpline_.evalDerivative(
swe);
 
 
  566    template<
class MaterialLaw>
 
  567    void initPcParameters_(
const MaterialLaw* m, 
const Scalar lowSwe, 
const Scalar highSwe)
 
  569        const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
 
  570        const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
 
  571        const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
 
  573        pcDerivativeLowSw_ = m->template dpc_dsw<false>(lowSw)*dsw_dswe;
 
  575        pcDerivativeHighSweThreshold_ = m->template dpc_dsw<false>(highSw)*dsw_dswe;
 
  576        pcDerivativeHighSweEnd_ = 2.0*(0.0 - m->template 
pc<false>(highSw))/(1.0 - highSwe);
 
  578        pcLowSwePcValue_ = m->template 
pc<false>(lowSw);
 
  579        pcHighSwePcValue_ = m->template 
pc<false>(highSw);
 
  586                                       pcHighSwePcValue_, 0, 
 
  587                                       pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_); 
 
  590    template<
class MaterialLaw>
 
  591    void initKrParameters_(
const MaterialLaw* m, 
const Scalar lowSwe, 
const Scalar highSwe)
 
  593        const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
 
  594        const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
 
  595        const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
 
  597        const auto krwHighSw = m->template 
krw<false>(highSw);
 
  598        const auto dkrwHighSw = m->template dkrw_dsw<false>(highSw)*dsw_dswe;
 
  600        const auto krnLowSw = m->template 
krn<false>(lowSw);
 
  601        const auto dkrnLowSw = m->template dkrn_dsw<false>(lowSw)*dsw_dswe;
 
  614    Scalar pcLowSwe_, pcHighSwe_;
 
  615    Scalar pcLowSwePcValue_, pcHighSwePcValue_;
 
  616    Scalar krwHighSwe_, krnLowSwe_;
 
  617    Scalar pcDerivativeLowSw_;
 
  618    Scalar pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_;
 
  620    Spline<Scalar> pcSpline_;
 
  621    Spline<Scalar> krwSpline_;
 
  622    Spline<Scalar> krnSpline_;
 
 
  629template<
typename Scalar = 
double>
 
  636template<
typename Scalar = 
double>
 
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
Definition efftoabsdefaultpolicy.hh:34
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition materiallaw.hh:45
Implementation of the van Genuchten capillary pressure <-> saturation relation, and relative permeabi...
Definition vangenuchten.hh:35
static Scalar dswe_dpc(Scalar pc, const Params< Scalar > ¶ms)
The partial derivative of the effective saturation to the capillary pressure according to van Genucht...
Definition vangenuchten.hh:194
static Scalar endPointPc(const Params< Scalar > ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition vangenuchten.hh:152
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > ¶ms)
The derivative of the relative permeability for the non-wetting phase in regard to the wetting satura...
Definition vangenuchten.hh:283
static Scalar swe(Scalar pc, const Params< Scalar > ¶ms)
The saturation-capillary pressure curve according to van Genuchten.
Definition vangenuchten.hh:136
static Scalar krw(Scalar swe, const Params< Scalar > ¶ms)
The relative permeability for the wetting phase of the medium implied by van Genuchten / Mualem param...
Definition vangenuchten.hh:216
static Scalar dkrw_dswe(Scalar swe, const Params< Scalar > ¶ms)
The derivative of the relative permeability for the wetting phase in regard to the wetting saturation...
Definition vangenuchten.hh:238
static Scalar dpc_dswe(Scalar swe, const Params< Scalar > ¶ms)
The partial derivative of the capillary pressure w.r.t. the effective saturation according to van Gen...
Definition vangenuchten.hh:172
static Scalar pc(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve according to van Genuchten.
Definition vangenuchten.hh:110
static Scalar krn(Scalar swe, const Params< Scalar > ¶ms)
The relative permeability for the non-wetting phase of the medium implied by van Genuchten's paramete...
Definition vangenuchten.hh:261
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition vangenuchten.hh:87
A regularization for the VanGenuchten material law.
Definition vangenuchten.hh:305
bool operator==(const VanGenuchtenRegularization &o) const
Equality comparison with another instance.
Definition vangenuchten.hh:404
void init(const MaterialLaw *m, const BaseParams &bp, const EffToAbsParams &etap, const Params< Scalar > &p)
Definition vangenuchten.hh:390
OptionalScalar< Scalar > dpc_dswe(const Scalar swe) const
The regularized partial derivative of the capillary pressure w.r.t. the saturation.
Definition vangenuchten.hh:442
void init(const MaterialLaw *m, const std::string ¶mGroup)
Initialize the spline.
Definition vangenuchten.hh:378
OptionalScalar< Scalar > krw(const Scalar swe) const
The regularized relative permeability for the wetting phase.
Definition vangenuchten.hh:508
OptionalScalar< Scalar > krn(const Scalar swe) const
The regularized relative permeability for the non-wetting phase.
Definition vangenuchten.hh:538
OptionalScalar< Scalar > swe(const Scalar pc) const
The regularized saturation-capillary pressure curve.
Definition vangenuchten.hh:460
OptionalScalar< Scalar > dkrn_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the non-wetting phase w....
Definition vangenuchten.hh:553
OptionalScalar< Scalar > pc(const Scalar swe) const
The regularized capillary pressure-saturation curve regularized part:
Definition vangenuchten.hh:419
OptionalScalar< Scalar > dswe_dpc(const Scalar pc) const
The regularized partial derivative of the saturation to the capillary pressure.
Definition vangenuchten.hh:484
OptionalScalar< Scalar > dkrw_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the wetting phase w.r....
Definition vangenuchten.hh:523
A 3rd order polynomial spline.
Definition spline.hh:43
TwoPMaterialLaw< Scalar, VanGenuchten, NoRegularization, TwoPEffToAbsDefaultPolicy > VanGenuchtenNoReg
A default configuration without regularization for using the VanGenuchten material law.
Definition vangenuchten.hh:637
TwoPMaterialLaw< Scalar, VanGenuchten, VanGenuchtenRegularization< Scalar >, TwoPEffToAbsDefaultPolicy > VanGenuchtenDefault
A default configuration for using the VanGenuchten material law.
Definition vangenuchten.hh:630
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
Definition brookscorey.hh:23
A wrapper that can either contain a valid Scalar or NaN.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Provides 3rd order polynomial splines.
The parameter type.
Definition vangenuchten.hh:53
Params(Scalar alpha, Scalar n, Scalar l=0.5)
Definition vangenuchten.hh:54
Scalar alpha() const
Definition vangenuchten.hh:58
bool operator==(const Params &p) const
Definition vangenuchten.hh:70
Scalar l() const
Definition vangenuchten.hh:67
void setAlpha(Scalar alpha)
Definition vangenuchten.hh:59
void setL(Scalar l)
Definition vangenuchten.hh:68
void setM(Scalar m)
Definition vangenuchten.hh:62
Scalar n() const
Definition vangenuchten.hh:64
Scalar m() const
Definition vangenuchten.hh:61
void setN(Scalar n)
Definition vangenuchten.hh:65
Regularization parameters.
Definition vangenuchten.hh:310
Scalar krnLowSwe() const
Threshold saturation below which the relative permeability of the non-wetting phase gets regularized.
Definition vangenuchten.hh:352
void setKrwHighSwe(Scalar krwHighSwe)
Set the threshold saturation above which the relative permeability of the wetting phase gets regulari...
Definition vangenuchten.hh:359
Scalar krwHighSwe() const
Threshold saturation above which the relative permeability of the wetting phase gets regularized.
Definition vangenuchten.hh:366
void setKrnLowSwe(Scalar krnLowSwe)
Set the threshold saturation below which the relative permeability of the non-wetting phase gets regu...
Definition vangenuchten.hh:345
void setPcHighSwe(Scalar pcHighSwe)
Set the threshold saturation above which the capillary pressure is regularized.
Definition vangenuchten.hh:329
Scalar pcHighSwe() const
Threshold saturation above which the capillary pressure is regularized.
Definition vangenuchten.hh:338
void setPcLowSwe(Scalar pcLowSwe)
Set the threshold saturation below which the capillary pressure is regularized.
Definition vangenuchten.hh:317
Scalar pcLowSwe() const
Threshold saturation below which the capillary pressure is regularized.
Definition vangenuchten.hh:323
A type for an optional scalar (contains either a valid number or NaN)
Definition optionalscalar.hh:27