25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH
51template <
class ScalarT,
class ParamsT = VanGenuchtenParams<ScalarT> >
56 using Scalar =
typename Params::Scalar;
78 swe = clamp(swe, 0.0, 1.0);
80 const Scalar pc = pow(pow(swe, -1.0/params.vgm()) - 1, 1.0/params.vgn())/params.vgAlpha();
107 const Scalar sw = pow(pow(params.vgAlpha()*pc, params.vgn()) + 1, -params.vgm());
147 swe = clamp(swe, 0.0, 1.0);
149 const Scalar powSwe = pow(swe, -1/params.vgm());
150 return - 1.0/params.vgAlpha() * pow(powSwe - 1, 1.0/params.vgn() - 1)/params.vgn()
151 * powSwe/swe/params.vgm();
174 const Scalar powAlphaPc = pow(params.vgAlpha()*pc, params.vgn());
175 return -pow(powAlphaPc + 1, -params.vgm()-1)*params.vgm()*powAlphaPc/pc*params.vgn();
197 swe = clamp(swe, 0.0, 1.0);
199 const Scalar r = 1.0 - pow(1.0 - pow(swe, 1.0/params.vgm()), params.vgm());
200 return pow(swe, params.vgl())*r*r;
222 swe = clamp(swe, 0.0, 1.0);
224 const Scalar x = 1.0 - pow(swe, 1.0/params.vgm());
225 const Scalar xToM = pow(x, params.vgm());
226 return (1.0 - xToM)*pow(swe, params.vgl()-1) * ( (1.0 - xToM)*params.vgl() + 2*xToM*(1.0-x)/x );
249 swe = clamp(swe, 0.0, 1.0);
251 return pow(1 - swe, params.vgl()) * pow(1 - pow(swe, 1.0/params.vgm()), 2*params.vgm());
274 swe = clamp(swe, 0.0, 1.0);
276 const auto sne = 1.0 - swe;
277 const auto x = 1.0 - pow(swe, 1.0/params.vgm());
278 return -pow(sne, params.vgl()-1.0) * pow(x, 2*params.vgm() - 1.0) * ( params.vgl()*x + 2.0*sne/swe*(1.0 - x) );
321 template<
class Scalar>
325 : alpha_(
alpha), n_(
n), m_(1.0 - 1.0/
n), l_(
l)
328 Scalar
alpha()
const {
return alpha_; }
331 Scalar
m()
const {
return m_; }
332 void setM(Scalar
m) { m_ =
m; n_ = 1.0/(1.0 -
m); }
334 Scalar
n()
const{
return n_; }
335 void setN(Scalar
n){ n_ =
n; m_ = 1.0 - 1.0/
n; }
337 Scalar
l()
const {
return l_; }
342 return Dune::FloatCmp::eq(alpha_, p.alpha_, 1e-6)
343 && Dune::FloatCmp::eq(n_, p.n_, 1e-6)
344 && Dune::FloatCmp::eq(m_, p.m_, 1e-6)
345 && Dune::FloatCmp::eq(l_, p.l_, 1e-6);
349 Scalar alpha_, n_, m_, l_;
356 template<
class Scalar =
double>
359 const auto n = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenN");
360 const auto alpha = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenAlpha");
362 const auto l = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenL", 0.5);
379 template<
class Scalar>
385 swe = clamp(
swe, 0.0, 1.0);
387 const Scalar
pc = pow(pow(
swe, -1.0/params.
m()) - 1, 1.0/params.
n())/params.
alpha();
405 template<
class Scalar>
413 const Scalar sw = pow(pow(params.
alpha()*
pc, params.
n()) + 1, -params.
m());
421 template<
class Scalar>
441 template<
class Scalar>
447 swe = clamp(
swe, 0.0, 1.0);
449 const Scalar powSwe = pow(
swe, -1/params.
m());
450 return - 1.0/params.
alpha() * pow(powSwe - 1, 1.0/params.
n() - 1)/params.
n()
451 * powSwe/
swe/params.
m();
463 template<
class Scalar>
471 const Scalar powAlphaPc = pow(params.
alpha()*
pc, params.
n());
472 return -pow(powAlphaPc + 1, -params.
m()-1)*params.
m()*powAlphaPc/
pc*params.
n();
485 template<
class Scalar>
491 swe = clamp(
swe, 0.0, 1.0);
493 const Scalar r = 1.0 - pow(1.0 - pow(
swe, 1.0/params.
m()), params.
m());
494 return pow(
swe, params.
l())*r*r;
507 template<
class Scalar>
513 swe = clamp(
swe, 0.0, 1.0);
515 const Scalar x = 1.0 - pow(
swe, 1.0/params.
m());
516 const Scalar xToM = pow(x, params.
m());
517 return (1.0 - xToM)*pow(
swe, params.
l()-1) * ( (1.0 - xToM)*params.
l() + 2*xToM*(1.0-x)/x );
530 template<
class Scalar>
536 swe = clamp(
swe, 0.0, 1.0);
538 return pow(1 -
swe, params.
l()) * pow(1 - pow(
swe, 1.0/params.
m()), 2*params.
m());
552 template<
class Scalar>
558 swe = clamp(
swe, 0.0, 1.0);
560 const auto sne = 1.0 -
swe;
561 const auto x = 1.0 - pow(
swe, 1.0/params.
m());
562 return -pow(sne, params.
l()-1.0) * pow(x, 2*params.
m() - 1.0) * ( params.
l()*x + 2.0*sne/
swe*(1.0 - x) );
573template <
class Scalar>
594 {
return pcLowSwe_; }
609 {
return pcHighSwe_; }
623 {
return krnLowSwe_; }
637 {
return krwHighSwe_; }
647 template<
class MaterialLaw>
648 void init(
const MaterialLaw* m,
const std::string& paramGroup)
650 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenPcLowSweThreshold", 0.01);
651 pcHighSwe_ = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenPcHighSweThreshold", 0.99);
652 krwHighSwe_ = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenKrwHighSweThreshold", 0.9);
653 krnLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenKrnLowSweThreshold", 0.1);
655 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
656 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
659 template<
class MaterialLaw,
class BaseParams,
class EffToAbsParams>
660 void init(
const MaterialLaw* m,
const BaseParams& bp,
const EffToAbsParams& etap,
const Params<Scalar>& p)
667 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
668 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
676 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
677 && Dune::FloatCmp::eq(pcHighSwe_, o.pcHighSwe_, 1e-6)
678 && Dune::FloatCmp::eq(krwHighSwe_, o.krwHighSwe_, 1e-6)
679 && Dune::FloatCmp::eq(krnLowSwe_, o.krnLowSwe_, 1e-6);
697 return pcLowSwePcValue_ + pcDerivativeLowSw_*(
swe - pcLowSwe_);
700 return pcDerivativeHighSweEnd_*(
swe - 1.0);
702 else if (
swe > pcHighSwe_)
703 return pcSpline_.eval(
swe);
715 return pcDerivativeLowSw_;
718 return pcDerivativeHighSweEnd_;
720 else if (
swe > pcHighSwe_)
721 return pcSpline_.evalDerivative(
swe);
734 if (pcHighSwe_ > 1.0 - std::numeric_limits<Scalar>::epsilon())
737 return pc/pcDerivativeHighSweEnd_ + 1.0;
741 else if (
pc < pcHighSwePcValue_)
742 return pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0,
pc);
744 else if (
pc >= pcLowSwePcValue_)
745 return (
pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
758 if (pcHighSwe_ > 1.0 - std::numeric_limits<Scalar>::epsilon())
761 return 1.0/pcDerivativeHighSweEnd_;
765 else if (
pc < pcHighSwePcValue_)
766 return 1.0/pcSpline_.evalDerivative(pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0,
pc));
768 else if (
pc >= pcLowSwePcValue_)
769 return 1.0/pcDerivativeLowSw_;
782 else if (
swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
784 else if (
swe > krwHighSwe_)
785 return krwSpline_.eval(
swe);
797 else if (
swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
799 else if (
swe > krwHighSwe_)
800 return krwSpline_.evalDerivative(
swe);
812 else if (
swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
814 else if (
swe < krnLowSwe_)
815 return krnSpline_.eval(
swe);
827 else if (
swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
829 else if (
swe < krnLowSwe_)
830 return krnSpline_.evalDerivative(
swe);
836 template<
class MaterialLaw>
837 void initPcParameters_(
const MaterialLaw* m,
const Scalar lowSwe,
const Scalar highSwe)
839 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
840 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
841 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
843 pcDerivativeLowSw_ = m->template dpc_dsw<false>(lowSw)*dsw_dswe;
845 pcDerivativeHighSweThreshold_ = m->template dpc_dsw<false>(highSw)*dsw_dswe;
846 pcDerivativeHighSweEnd_ = 2.0*(0.0 - m->template pc<false>(highSw))/(1.0 - highSwe);
848 pcLowSwePcValue_ = m->template pc<false>(lowSw);
849 pcHighSwePcValue_ = m->template pc<false>(highSw);
852 pcHighSwePcValue_, 0,
853 pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_);
858 template<
class MaterialLaw>
859 void initKrParameters_(
const MaterialLaw* m,
const Scalar lowSwe,
const Scalar highSwe)
861 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
862 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
863 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
865 const auto krwHighSw = m->template krw<false>(highSw);
866 const auto dkrwHighSw = m->template dkrw_dsw<false>(highSw)*dsw_dswe;
868 const auto krnLowSw = m->template krn<false>(lowSw);
869 const auto dkrnLowSw = m->template dkrn_dsw<false>(lowSw)*dsw_dswe;
881 Scalar pcLowSwe_, pcHighSwe_;
882 Scalar pcLowSwePcValue_, pcHighSwePcValue_;
883 Scalar krwHighSwe_, krnLowSwe_;
884 Scalar pcDerivativeLowSw_;
885 Scalar pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_;
887 Spline<Scalar> pcSpline_;
888 Spline<Scalar> krwSpline_;
889 Spline<Scalar> krnSpline_;
896template<
typename Scalar =
double>
903template<
typename Scalar =
double>
Provides 3rd order polynomial splines.
A wrapper that can either contain a valid Scalar or NaN.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
Specification of the material parameters for the van Genuchten-Mualem constitutive relations.
Definition: brookscorey.hh:286
A 3rd order polynomial spline.
Definition: spline.hh:55
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
Definition: efftoabsdefaultpolicy.hh:46
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:57
Implementation of the van Genuchten capillary pressure <-> saturation relation. This class bundles th...
Definition: vangenuchten.hh:53
static Scalar dkrw_dswe(const Params ¶ms, Scalar swe)
The derivative of the relative permeability for the wetting phase in regard to the wetting saturation...
Definition: vangenuchten.hh:217
static Scalar sw(const Params ¶ms, Scalar pc)
The saturation-capillary pressure curve according to van Genuchten.
Definition: vangenuchten.hh:100
typename Params::Scalar Scalar
Definition: vangenuchten.hh:56
static Scalar dpc_dswe(const Params ¶ms, Scalar swe)
The partial derivative of the capillary pressure w.r.t. the effective saturation according to van Gen...
Definition: vangenuchten.hh:142
static Scalar krw(const Params ¶ms, Scalar swe)
The relative permeability for the wetting phase of the medium implied by van Genuchten / Mualem param...
Definition: vangenuchten.hh:192
static Scalar dswe_dpc(const Params ¶ms, Scalar pc)
The partial derivative of the effective saturation to the capillary pressure according to van Genucht...
Definition: vangenuchten.hh:167
static Scalar endPointPc(const Params ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: vangenuchten.hh:119
static Scalar krn(const Params ¶ms, Scalar swe)
The relative permeability for the nonwetting phase of the medium implied by van Genuchten's parameter...
Definition: vangenuchten.hh:244
ParamsT Params
Definition: vangenuchten.hh:55
static Scalar pc(const Params ¶ms, Scalar swe)
The capillary pressure-saturation curve according to van Genuchten.
Definition: vangenuchten.hh:73
static Scalar dkrn_dswe(const Params ¶ms, Scalar swe)
The derivative of the relative permeability for the nonwetting phase in regard to the wetting saturat...
Definition: vangenuchten.hh:269
Implementation of the van Genuchten capillary pressure <-> saturation relation, and relative permeabi...
Definition: vangenuchten.hh:305
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:464
static Scalar endPointPc(const Params< Scalar > ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: vangenuchten.hh:422
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:553
static Scalar swe(Scalar pc, const Params< Scalar > ¶ms)
The saturation-capillary pressure curve according to van Genuchten.
Definition: vangenuchten.hh:406
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:486
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:508
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:442
static Scalar pc(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve according to van Genuchten.
Definition: vangenuchten.hh:380
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:531
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: vangenuchten.hh:357
The parameter type.
Definition: vangenuchten.hh:323
Params(Scalar alpha, Scalar n, Scalar l=0.5)
Definition: vangenuchten.hh:324
Scalar alpha() const
Definition: vangenuchten.hh:328
bool operator==(const Params &p) const
Definition: vangenuchten.hh:340
Scalar l() const
Definition: vangenuchten.hh:337
void setAlpha(Scalar alpha)
Definition: vangenuchten.hh:329
void setL(Scalar l)
Definition: vangenuchten.hh:338
void setM(Scalar m)
Definition: vangenuchten.hh:332
Scalar n() const
Definition: vangenuchten.hh:334
Scalar m() const
Definition: vangenuchten.hh:331
void setN(Scalar n)
Definition: vangenuchten.hh:335
A regularization for the VanGenuchten material law.
Definition: vangenuchten.hh:575
bool operator==(const VanGenuchtenRegularization &o) const
Equality comparison with another instance.
Definition: vangenuchten.hh:674
void init(const MaterialLaw *m, const BaseParams &bp, const EffToAbsParams &etap, const Params< Scalar > &p)
Definition: vangenuchten.hh:660
OptionalScalar< Scalar > dpc_dswe(const Scalar swe) const
The regularized partial derivative of the capillary pressure w.r.t. the saturation.
Definition: vangenuchten.hh:712
void init(const MaterialLaw *m, const std::string ¶mGroup)
Initialize the spline.
Definition: vangenuchten.hh:648
OptionalScalar< Scalar > krw(const Scalar swe) const
The regularized relative permeability for the wetting phase.
Definition: vangenuchten.hh:778
OptionalScalar< Scalar > krn(const Scalar swe) const
The regularized relative permeability for the non-wetting phase.
Definition: vangenuchten.hh:808
OptionalScalar< Scalar > swe(const Scalar pc) const
The regularized saturation-capillary pressure curve.
Definition: vangenuchten.hh:730
OptionalScalar< Scalar > dkrn_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the non-wetting phase w....
Definition: vangenuchten.hh:823
OptionalScalar< Scalar > pc(const Scalar swe) const
The regularized capillary pressure-saturation curve regularized part:
Definition: vangenuchten.hh:689
OptionalScalar< Scalar > dswe_dpc(const Scalar pc) const
The regularized partial derivative of the saturation to the capillary pressure.
Definition: vangenuchten.hh:754
OptionalScalar< Scalar > dkrw_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the wetting phase w.r....
Definition: vangenuchten.hh:793
Regularization parameters.
Definition: vangenuchten.hh:580
Scalar krnLowSwe() const
Threshold saturation below which the relative permeability of the non-wetting phase gets regularized.
Definition: vangenuchten.hh:622
void setKrwHighSwe(Scalar krwHighSwe)
Set the threshold saturation above which the relative permeability of the wetting phase gets regulari...
Definition: vangenuchten.hh:629
Scalar krwHighSwe() const
Threshold saturation above which the relative permeability of the wetting phase gets regularized.
Definition: vangenuchten.hh:636
void setKrnLowSwe(Scalar krnLowSwe)
Set the threshold saturation below which the relative permeability of the non-wetting phase gets regu...
Definition: vangenuchten.hh:615
void setPcHighSwe(Scalar pcHighSwe)
Set the threshold saturation above which the capillary pressure is regularized.
Definition: vangenuchten.hh:599
Scalar pcHighSwe() const
Threshold saturation above which the capillary pressure is regularized.
Definition: vangenuchten.hh:608
void setPcLowSwe(Scalar pcLowSwe)
Set the threshold saturation below which the capillary pressure is regularized.
Definition: vangenuchten.hh:587
Scalar pcLowSwe() const
Threshold saturation below which the capillary pressure is regularized.
Definition: vangenuchten.hh:593