25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH
63 template<
class Scalar>
67 : alpha_(
alpha), n_(
n), m_(1.0 - 1.0/
n), l_(
l)
70 Scalar
alpha()
const {
return alpha_; }
73 Scalar
m()
const {
return m_; }
74 void setM(Scalar
m) { m_ =
m; n_ = 1.0/(1.0 -
m); }
76 Scalar
n()
const{
return n_; }
77 void setN(Scalar
n){ n_ =
n; m_ = 1.0 - 1.0/
n; }
79 Scalar
l()
const {
return l_; }
84 return Dune::FloatCmp::eq(alpha_, p.alpha_, 1e-6)
85 && Dune::FloatCmp::eq(n_, p.n_, 1e-6)
86 && Dune::FloatCmp::eq(m_, p.m_, 1e-6)
87 && Dune::FloatCmp::eq(l_, p.l_, 1e-6);
91 Scalar alpha_, n_, m_, l_;
98 template<
class Scalar =
double>
101 const auto n = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenN");
102 const auto alpha = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenAlpha");
104 const auto l = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenL", 0.5);
121 template<
class Scalar>
127 swe = clamp(
swe, 0.0, 1.0);
129 const Scalar
pc = pow(pow(
swe, -1.0/params.
m()) - 1, 1.0/params.
n())/params.
alpha();
147 template<
class Scalar>
155 const Scalar sw = pow(pow(params.
alpha()*
pc, params.
n()) + 1, -params.
m());
163 template<
class Scalar>
183 template<
class Scalar>
189 swe = clamp(
swe, 0.0, 1.0);
191 const Scalar powSwe = pow(
swe, -1/params.
m());
192 return - 1.0/params.
alpha() * pow(powSwe - 1, 1.0/params.
n() - 1)/params.
n()
193 * powSwe/
swe/params.
m();
205 template<
class Scalar>
213 const Scalar powAlphaPc = pow(params.
alpha()*
pc, params.
n());
214 return -pow(powAlphaPc + 1, -params.
m()-1)*params.
m()*powAlphaPc/
pc*params.
n();
227 template<
class Scalar>
233 swe = clamp(
swe, 0.0, 1.0);
235 const Scalar r = 1.0 - pow(1.0 - pow(
swe, 1.0/params.
m()), params.
m());
236 return pow(
swe, params.
l())*r*r;
249 template<
class Scalar>
255 swe = clamp(
swe, 0.0, 1.0);
257 const Scalar x = 1.0 - pow(
swe, 1.0/params.
m());
258 const Scalar xToM = pow(x, params.
m());
259 return (1.0 - xToM)*pow(
swe, params.
l()-1) * ( (1.0 - xToM)*params.
l() + 2*xToM*(1.0-x)/x );
272 template<
class Scalar>
278 swe = clamp(
swe, 0.0, 1.0);
280 return pow(1 -
swe, params.
l()) * pow(1 - pow(
swe, 1.0/params.
m()), 2*params.
m());
294 template<
class Scalar>
300 swe = clamp(
swe, 0.0, 1.0);
302 const auto sne = 1.0 -
swe;
303 const auto x = 1.0 - pow(
swe, 1.0/params.
m());
304 return -pow(sne, params.
l()-1.0) * pow(x, 2*params.
m() - 1.0) * ( params.
l()*x + 2.0*sne/
swe*(1.0 - x) );
315template <
class Scalar>
336 {
return pcLowSwe_; }
351 {
return pcHighSwe_; }
365 {
return krnLowSwe_; }
379 {
return krwHighSwe_; }
389 template<
class MaterialLaw>
390 void init(
const MaterialLaw* m,
const std::string& paramGroup)
392 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenPcLowSweThreshold", 0.01);
393 pcHighSwe_ = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenPcHighSweThreshold", 0.99);
394 krwHighSwe_ = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenKrwHighSweThreshold", 0.9);
395 krnLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"VanGenuchtenKrnLowSweThreshold", 0.1);
397 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
398 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
401 template<
class MaterialLaw,
class BaseParams,
class EffToAbsParams>
402 void init(
const MaterialLaw* m,
const BaseParams& bp,
const EffToAbsParams& etap,
const Params<Scalar>& p)
409 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
410 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
418 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
419 && Dune::FloatCmp::eq(pcHighSwe_, o.pcHighSwe_, 1e-6)
420 && Dune::FloatCmp::eq(krwHighSwe_, o.krwHighSwe_, 1e-6)
421 && Dune::FloatCmp::eq(krnLowSwe_, o.krnLowSwe_, 1e-6);
438 if (
swe <= pcLowSwe_)
439 return pcLowSwePcValue_ + pcDerivativeLowSw_*(
swe - pcLowSwe_);
442 return pcDerivativeHighSweEnd_*(
swe - 1.0);
444 else if (
swe > pcHighSwe_)
445 return pcSpline_.eval(
swe);
456 if (
swe <= pcLowSwe_)
457 return pcDerivativeLowSw_;
460 return pcDerivativeHighSweEnd_;
462 else if (
swe > pcHighSwe_)
463 return pcSpline_.evalDerivative(
swe);
476 if (pcHighSwe_ >= 1.0)
479 return pc/pcDerivativeHighSweEnd_ + 1.0;
483 else if (
pc <= pcHighSwePcValue_)
484 return pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0,
pc);
486 else if (
pc >= pcLowSwePcValue_)
487 return (
pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
500 if (pcHighSwe_ >= 1.0)
503 return 1.0/pcDerivativeHighSweEnd_;
507 else if (
pc <= pcHighSwePcValue_)
508 return 1.0/pcSpline_.evalDerivative(pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0,
pc));
510 else if (
pc >= pcLowSwePcValue_)
511 return 1.0/pcDerivativeLowSw_;
526 else if (
swe >= krwHighSwe_)
527 return krwSpline_.eval(
swe);
541 else if (
swe >= krwHighSwe_)
542 return krwSpline_.evalDerivative(
swe);
556 else if (
swe <= krnLowSwe_)
557 return krnSpline_.eval(
swe);
571 else if (
swe <= krnLowSwe_)
572 return krnSpline_.evalDerivative(
swe);
578 template<
class MaterialLaw>
579 void initPcParameters_(
const MaterialLaw* m,
const Scalar lowSwe,
const Scalar highSwe)
581 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
582 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
583 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
585 pcDerivativeLowSw_ = m->template dpc_dsw<false>(lowSw)*dsw_dswe;
587 pcDerivativeHighSweThreshold_ = m->template dpc_dsw<false>(highSw)*dsw_dswe;
588 pcDerivativeHighSweEnd_ = 2.0*(0.0 - m->template pc<false>(highSw))/(1.0 - highSwe);
590 pcLowSwePcValue_ = m->template pc<false>(lowSw);
591 pcHighSwePcValue_ = m->template pc<false>(highSw);
598 pcHighSwePcValue_, 0,
599 pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_);
602 template<
class MaterialLaw>
603 void initKrParameters_(
const MaterialLaw* m,
const Scalar lowSwe,
const Scalar highSwe)
605 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
606 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
607 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
609 const auto krwHighSw = m->template krw<false>(highSw);
610 const auto dkrwHighSw = m->template dkrw_dsw<false>(highSw)*dsw_dswe;
612 const auto krnLowSw = m->template krn<false>(lowSw);
613 const auto dkrnLowSw = m->template dkrn_dsw<false>(lowSw)*dsw_dswe;
626 Scalar pcLowSwe_, pcHighSwe_;
627 Scalar pcLowSwePcValue_, pcHighSwePcValue_;
628 Scalar krwHighSwe_, krnLowSwe_;
629 Scalar pcDerivativeLowSw_;
630 Scalar pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_;
632 Spline<Scalar> pcSpline_;
633 Spline<Scalar> krwSpline_;
634 Spline<Scalar> krnSpline_;
641template<
typename Scalar =
double>
648template<
typename Scalar =
double>
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.
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
Definition: brookscorey.hh:35
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, and relative permeabi...
Definition: vangenuchten.hh:47
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:206
static Scalar endPointPc(const Params< Scalar > ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: vangenuchten.hh:164
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:295
static Scalar swe(Scalar pc, const Params< Scalar > ¶ms)
The saturation-capillary pressure curve according to van Genuchten.
Definition: vangenuchten.hh:148
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:228
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:250
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:184
static Scalar pc(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve according to van Genuchten.
Definition: vangenuchten.hh:122
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:273
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: vangenuchten.hh:99
The parameter type.
Definition: vangenuchten.hh:65
Params(Scalar alpha, Scalar n, Scalar l=0.5)
Definition: vangenuchten.hh:66
Scalar alpha() const
Definition: vangenuchten.hh:70
bool operator==(const Params &p) const
Definition: vangenuchten.hh:82
Scalar l() const
Definition: vangenuchten.hh:79
void setAlpha(Scalar alpha)
Definition: vangenuchten.hh:71
void setL(Scalar l)
Definition: vangenuchten.hh:80
void setM(Scalar m)
Definition: vangenuchten.hh:74
Scalar n() const
Definition: vangenuchten.hh:76
Scalar m() const
Definition: vangenuchten.hh:73
void setN(Scalar n)
Definition: vangenuchten.hh:77
A regularization for the VanGenuchten material law.
Definition: vangenuchten.hh:317
bool operator==(const VanGenuchtenRegularization &o) const
Equality comparison with another instance.
Definition: vangenuchten.hh:416
void init(const MaterialLaw *m, const BaseParams &bp, const EffToAbsParams &etap, const Params< Scalar > &p)
Definition: vangenuchten.hh:402
OptionalScalar< Scalar > dpc_dswe(const Scalar swe) const
The regularized partial derivative of the capillary pressure w.r.t. the saturation.
Definition: vangenuchten.hh:454
void init(const MaterialLaw *m, const std::string ¶mGroup)
Initialize the spline.
Definition: vangenuchten.hh:390
OptionalScalar< Scalar > krw(const Scalar swe) const
The regularized relative permeability for the wetting phase.
Definition: vangenuchten.hh:520
OptionalScalar< Scalar > krn(const Scalar swe) const
The regularized relative permeability for the non-wetting phase.
Definition: vangenuchten.hh:550
OptionalScalar< Scalar > swe(const Scalar pc) const
The regularized saturation-capillary pressure curve.
Definition: vangenuchten.hh:472
OptionalScalar< Scalar > dkrn_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the non-wetting phase w....
Definition: vangenuchten.hh:565
OptionalScalar< Scalar > pc(const Scalar swe) const
The regularized capillary pressure-saturation curve regularized part:
Definition: vangenuchten.hh:431
OptionalScalar< Scalar > dswe_dpc(const Scalar pc) const
The regularized partial derivative of the saturation to the capillary pressure.
Definition: vangenuchten.hh:496
OptionalScalar< Scalar > dkrw_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the wetting phase w.r....
Definition: vangenuchten.hh:535
Regularization parameters.
Definition: vangenuchten.hh:322
Scalar krnLowSwe() const
Threshold saturation below which the relative permeability of the non-wetting phase gets regularized.
Definition: vangenuchten.hh:364
void setKrwHighSwe(Scalar krwHighSwe)
Set the threshold saturation above which the relative permeability of the wetting phase gets regulari...
Definition: vangenuchten.hh:371
Scalar krwHighSwe() const
Threshold saturation above which the relative permeability of the wetting phase gets regularized.
Definition: vangenuchten.hh:378
void setKrnLowSwe(Scalar krnLowSwe)
Set the threshold saturation below which the relative permeability of the non-wetting phase gets regu...
Definition: vangenuchten.hh:357
void setPcHighSwe(Scalar pcHighSwe)
Set the threshold saturation above which the capillary pressure is regularized.
Definition: vangenuchten.hh:341
Scalar pcHighSwe() const
Threshold saturation above which the capillary pressure is regularized.
Definition: vangenuchten.hh:350
void setPcLowSwe(Scalar pcLowSwe)
Set the threshold saturation below which the capillary pressure is regularized.
Definition: vangenuchten.hh:329
Scalar pcLowSwe() const
Threshold saturation below which the capillary pressure is regularized.
Definition: vangenuchten.hh:335