24#ifndef PARKER_VANGENUCHTEN_3P_HH
25#define PARKER_VANGENUCHTEN_3P_HH
44 template<
class Scalar>
89 return Dune::FloatCmp::eq(
swr(), p.
swr(), 1e-6)
90 && Dune::FloatCmp::eq(
snr(), p.
snr(), 1e-6)
91 && Dune::FloatCmp::eq(
sgr(), p.
sgr(), 1e-6);
103 template<
class Scalar>
107 params.
setSwr(getParamFromGroup<Scalar>(paramGroup,
"Swr", 0.0));
108 params.
setSnr(getParamFromGroup<Scalar>(paramGroup,
"Snr", 0.0));
109 params.
setSgr(getParamFromGroup<Scalar>(paramGroup,
"Sgr", 0.0));
122 template<
class Scalar>
125 return (sw - params.
swr())/(1.0 - params.
swr());
137 template<
class Scalar>
140 return swe*(1.0 - params.
swr()) + params.
swr();
151 template<
class Scalar>
154 return 1.0/(1.0 - params.
swr());
165 template<
class Scalar>
168 return 1.0 - params.
swr();
180 template<
class Scalar>
195 template<
class Scalar>
198 return (st-params.
swr()) / (1-params.
swr());
210 template<
class Scalar>
213 return ste*(1.0 - params.
swr()) + params.
swr();
224 template<
class Scalar>
227 return 1.0/(1.0 - params.
swr() );
238 template<
class Scalar>
241 return 1.0 - params.
swr();
270 template<
class Scalar>
274 Scalar
betaNw = 1.0, Scalar
betaGn = 1.0, Scalar
betaGw = 1.0,
bool regardSnr =
false)
275 : alpha_(
alpha), n_(
n), m_(1.0 - 1.0/
n), swr_(
swr), snr_(
snr)
279 Scalar
alpha()
const {
return alpha_; }
282 Scalar
m()
const {
return m_; }
283 void setM(Scalar
m) { m_ =
m; n_ = 1.0/(1.0 -
m); }
285 Scalar
n()
const{
return n_; }
286 void setN(Scalar
n){ n_ =
n; m_ = 1.0 - 1.0/
n; }
288 Scalar
swr()
const {
return swr_; }
291 Scalar
snr()
const {
return snr_; }
294 Scalar
betaNw()
const {
return betaNw_; }
297 Scalar
betaGn()
const {
return betaGn_; }
300 Scalar
betaGw()
const {
return betaGw_; }
308 return Dune::FloatCmp::eq(alpha_, p.alpha_, 1e-6)
309 && Dune::FloatCmp::eq(n_, p.n_, 1e-6)
310 && Dune::FloatCmp::eq(m_, p.m_, 1e-6)
311 && Dune::FloatCmp::eq(swr_, p.swr_, 1e-6)
312 && Dune::FloatCmp::eq(snr_, p.snr_, 1e-6)
313 && Dune::FloatCmp::eq(betaNw_, p.betaNw_, 1e-6)
314 && Dune::FloatCmp::eq(betaGn_, p.betaGn_, 1e-6)
315 && Dune::FloatCmp::eq(betaGw_, p.betaGw_, 1e-6)
316 && regardSnr_ == p.regardSnr_;
320 Scalar alpha_, n_, m_, swr_, snr_, betaNw_, betaGn_, betaGw_;
328 template<
class Scalar =
double>
331 const auto n = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenN");
332 const auto alpha = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenAlpha");
333 const auto swr = getParamFromGroup<Scalar>(paramGroup,
"Swr", 0.0);
334 const auto snr = getParamFromGroup<Scalar>(paramGroup,
"Snr", 0.0);
335 const auto betaNw = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenBetaNw", 1.0);
336 const auto betaGn = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenBetaGn", 1.0);
337 const auto betaGw = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenBetaGw", 1.0);
338 const auto regardSnr = getParamFromGroup<bool>(paramGroup,
"ParkerVanGenuchtenRegardSnrForKrn",
false);
340 betaNw, betaGn, betaGw, regardSnr );
348 template<
class Scalar>
351 assert(0 <= swe && swe <= 1);
352 return pc_(swe, params);
360 template<
class Scalar>
363 assert(0 <= swe && swe <= 1);
364 return pc_(swe, params)/params.
betaNw();
372 template<
class Scalar>
375 assert(0 <= ste && ste <= 1);
376 return pc_(ste, params)/params.
betaGn();
384 template<
class Scalar>
393 if (sne > params.
snr())
397 if (params.
snr() >= 0.001)
398 return sne/params.
snr();
410 template<
class Scalar>
414 const Scalar powSeRegu = pow(swe, -1/params.
m());
415 return - 1.0/params.
alpha() * pow(powSeRegu - 1, 1.0/params.
n() - 1)/params.
n()
416 * powSeRegu/swe/params.
m()/params.
betaGw();
425 template<
class Scalar>
429 const Scalar powSeRegu = pow(swe, -1/params.
m());
430 return - 1.0/params.
alpha() * pow(powSeRegu - 1, 1.0/params.
n() - 1)/params.
n()
431 * powSeRegu/swe/params.
m()/params.
betaNw();
440 template<
class Scalar>
444 const Scalar powSeRegu = pow(ste, -1/params.
m());
445 return - 1.0/params.
alpha() * pow(powSeRegu - 1, 1.0/params.
n() - 1)/params.
n()
446 * powSeRegu/ste/params.
m()/params.
betaGn();
461 template<
class Scalar>
466 const Scalar r = 1.0 - pow(1 - pow(swe, 1/params.
m()), params.
m());
467 return sqrt(swe)*r*r;
487 template<
class Scalar>
488 static Scalar
krn(
const Scalar swe,
const Scalar sn,
const Scalar ste,
const Params<Scalar>& params)
492 krn = pow(1 - pow(swe, 1/params.
m()), params.
m());
493 krn -= pow(1 - pow(ste, 1/params.
m()), params.
m());
501 const Scalar resIncluded = clamp(sn - params.
snr()/ (1-params.
swr()), 0.0, 1.0);
502 krn *= sqrt(resIncluded);
505 krn *= sqrt(sn / (1 - params.
swr()));
522 template<
class Scalar>
525 assert(0 <= ste && ste <= 1);
528 return cbrt(1 - ste) * pow(1 - pow(ste, 1/params.
m()), 2*params.
m());
540 template<
class Scalar>
543 assert(0 < ste && ste <= 1);
546 const Scalar x = pow(ste, 1.0/params.
m());
548 -pow(1.0 - x, 2*params.
m())
549 *pow(1.0 - ste, -2.0/3)
560 template<
class Scalar>
561 static Scalar
kr(
const int phaseIdx,
const Scalar swe,
const Scalar sne,
const Params<Scalar>& params)
566 return krw(params, swe, sne);
568 return krn(params, swe, sne);
570 return krg(params, swe, sne);
572 DUNE_THROW(Dune::InvalidStateException,
573 "Invalid phase index ");
584 template<
class Scalar>
585 const static Scalar pc_(
const Scalar se,
const Params<Scalar>& params)
588 return pow(pow(se, -1/params.m()) - 1, 1/params.n())/params.alpha();
600template <
class Scalar>
623 {
return pcLowSwe_; }
638 {
return pcHighSwe_; }
652 {
return krnLowSwe_; }
666 {
return krgLowSte_; }
680 {
return krwHighSwe_; }
689 { constRegularization_ = input; }
696 {
return constRegularization_; }
704 bool constRegularization_ =
false;
708 template<
class MaterialLaw>
709 void init(
const MaterialLaw* m,
const std::string& paramGroup)
711 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenPcLowSweThreshold", 0.01);
712 pcHighSwe_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenPcHighSweThreshold", 0.99);
713 krwHighSwe_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenKrwHighSweThreshold", 0.9);
714 krnLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenKrnLowSweThreshold", 0.1);
715 krgLowSte_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenKrgLowSteThreshold", 1e-3);
716 constRegularization_ = getParamFromGroup<bool>(paramGroup,
"VanGenuchtenConstantRegularization",
false);
718 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
719 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
722 template<
class MaterialLaw,
class BaseParams,
class EffToAbsParams>
723 void init(
const MaterialLaw* m,
const BaseParams& bp,
const EffToAbsParams& etap,
const Params<Scalar>& p)
732 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
733 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
741 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
742 && Dune::FloatCmp::eq(pcHighSwe_, o.pcHighSwe_, 1e-6)
743 && Dune::FloatCmp::eq(krwHighSwe_, o.krwHighSwe_, 1e-6)
744 && Dune::FloatCmp::eq(krnLowSwe_, o.krnLowSwe_, 1e-6)
745 && constRegularization_ == o.constRegularization_;
760 if (constRegularization_)
761 swe = clamp(swe, 0.0, 1.0);
769 return pcgwLowSwePcgwValue_ + pcgwDerivativeLowSw_*(swe - pcLowSwe_);
772 return pcgwDerivativeHighSweEnd_*(swe - 1.0);
774 else if (swe > pcHighSwe_)
775 return pcgwSpline_.eval(swe);
793 if (constRegularization_)
794 swe = clamp(swe, 0.0, 1.0);
802 return pcnwLowSwePcnwValue_ + pcnwDerivativeLowSw_*(swe - pcLowSwe_);
805 return pcnwDerivativeHighSweEnd_*(swe - 1.0);
807 else if (swe > pcHighSwe_)
808 return pcnwSpline_.eval(swe);
826 if (constRegularization_)
827 ste = clamp(ste, 0.0, 1.0);
835 const Scalar pcLowSte = pcLowSwe_;
836 const Scalar pcHighSte = pcHighSwe_;
838 return pcgnLowStePcgnValue_ + pcgnDerivativeLowSt_*(ste - pcLowSte);
841 return pcgnDerivativeHighSteEnd_*(ste - 1.0);
843 else if (ste > pcHighSte)
844 return pcgnSpline_.eval(ste);
868 else if (swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
884 swe = clamp(swe, 0.0, 1.0);
885 ste = clamp(ste, 0.0, 1.0);
887 if (ste - swe <= 0.0)
900 if (ste > 1.0 - std::numeric_limits<Scalar>::epsilon())
906 if (ste <= krgLowSte_)
907 return krgLowStkrgValue_ + krgDerivativeLowSt_*(ste - krgLowSte_);
917 const Scalar st = ste*(1 - swr_) + swr_;
918 const Scalar sg = 1.0 - st;
926 const Scalar scalFact = max(0.0, (sg - sgr_)/(0.1 - sgr_));
943 return krw(swe, sne);
945 return krn(swe, sne);
947 return krg(swe, sne);
949 DUNE_THROW(Dune::InvalidStateException,
950 "Invalid phase index ");
954 template<
class MaterialLaw>
955 void initPcParameters_(
const MaterialLaw* m,
const Scalar lowSwe,
const Scalar highSwe)
957 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
958 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
959 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
960 const auto dst_dste = MaterialLaw::EffToAbs::dst_dste(m->effToAbsParams());
962 baseLawParamsPtr_ = &m->basicParams();
965 pcgwLowSwePcgwValue_ = m->template pcgw<false>(lowSw, 0.0);
966 pcgwDerivativeLowSw_ = m->template dpcgw_dsw<false>(lowSw, 0.0)*dsw_dswe;
967 pcgwHighSwePcgwValue_ = m->template pcgw<false>(highSw, 0.0);
968 pcgwDerivativeHighSweThreshold_ = m->template dpcgw_dsw<false>(highSw, 0.0)*dsw_dswe;
969 pcgwDerivativeHighSweEnd_ = 2.0*(0.0 - m->template pcgw<false>(highSw, 0.0))/(1.0 - highSwe);
971 pcgwHighSwePcgwValue_, 0,
972 pcgwDerivativeHighSweThreshold_, pcgwDerivativeHighSweEnd_);
975 pcnwLowSwePcnwValue_ = m->template pcnw<false>(lowSw, 0.0);
976 pcnwDerivativeLowSw_ = m->template dpcnw_dsw<false>(lowSw, 0.0)*dsw_dswe;
977 pcnwHighSwePcnwValue_ = m->template pcnw<false>(highSw, 0.0);
978 pcnwDerivativeHighSweThreshold_ = m->template dpcnw_dsw<false>(highSw, 0.0);
979 pcnwDerivativeHighSweEnd_ = 2.0*(0.0 - m->template pcnw<false>(highSw, 0.0))/(1.0 - highSwe);
981 pcnwHighSwePcnwValue_, 0,
982 pcnwDerivativeHighSweThreshold_, pcnwDerivativeHighSweEnd_);
985 pcgnLowStePcgnValue_ = m->template pcgn<false>(lowSw, 0.0);
986 pcgnDerivativeLowSt_ = m->template dpcgn_dst<false>(lowSw, 0.0)*dst_dste;
987 pcgnHighSwePcgnValue_ = m->template pcgn<false>(highSw, 0.0);
988 pcgnDerivativeHighSteThreshold_ = m->template dpcgn_dst<false>(highSw, 0.0);
989 pcgnDerivativeHighSteEnd_ = 2.0*(0.0 - m->template pcgn<false>(highSw, 0.0))/(1.0 - highSwe);
991 pcgnHighSwePcgnValue_, 0,
992 pcgnDerivativeHighSteThreshold_, pcgnDerivativeHighSteEnd_);
996 template<
class MaterialLaw>
997 void initKrParameters_(
const MaterialLaw* m,
const Scalar lowSwe,
const Scalar highSwe)
1002 swr_ = m->effToAbsParams().swr();
1003 sgr_ = m->effToAbsParams().sgr();
1006 Scalar krgLowStkrgValue_;
1007 Scalar krgDerivativeLowSt_;
1009 Scalar pcLowSwe_, pcHighSwe_;
1010 Scalar krwHighSwe_, krnLowSwe_, krgLowSte_;
1013 Scalar pcgwLowSwePcgwValue_;
1014 Scalar pcgwHighSwePcgwValue_;
1015 Scalar pcgwDerivativeLowSw_;
1016 Scalar pcgwDerivativeHighSweThreshold_;
1017 Scalar pcgwDerivativeHighSweEnd_;
1020 Scalar pcgnLowStePcgnValue_;
1021 Scalar pcgnHighSwePcgnValue_;
1022 Scalar pcgnDerivativeLowSt_;
1023 Scalar pcgnDerivativeHighSteThreshold_;
1024 Scalar pcgnDerivativeHighSteEnd_;
1027 Scalar pcnwLowSwePcnwValue_;
1028 Scalar pcnwHighSwePcnwValue_;
1029 Scalar pcnwDerivativeLowSw_;
1030 Scalar pcnwDerivativeHighSweThreshold_;
1031 Scalar pcnwDerivativeHighSweEnd_;
1033 Spline<Scalar> pcgwSpline_;
1034 Spline<Scalar> pcnwSpline_;
1035 Spline<Scalar> pcgnSpline_;
1036 Spline<Scalar> krwSpline_;
1037 Spline<Scalar> krnSpline_;
1041 bool constRegularization_;
1043 const BaseLawParams* baseLawParamsPtr_;
1046template<
class ScalarType,
1048 class Regularization = NoRegularization,
1049 class EffToAbsPolicy = ParkerVanGenuchten3PEffToAbsPolicy>
1066 {
return !std::is_same<Regularization, NoRegularization>::value; }
1083 regularization_.init(
this, paramGroup);
1093 : basicParams_(baseParams)
1097 regularization_.init(
this, baseParams,
effToAbsParams, regParams);
1103 template<
bool enableRegularization = isRegularized()>
1106 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1107 if constexpr (enableRegularization)
1109 const auto regularized = regularization_.pcgw(swe);
1111 return regularized.value();
1114 return BaseLaw::pcgw(swe, basicParams_);
1120 template<
bool enableRegularization = isRegularized()>
1123 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1124 if constexpr (enableRegularization)
1126 const auto regularized = regularization_.pcnw(swe);
1128 return regularized.value();
1131 return BaseLaw::pcnw(swe, basicParams_);
1139 template<
bool enableRegularization = isRegularized()>
1142 const auto swe = EffToAbs::swToSwe(sw + sn, effToAbsParams_);
1143 if constexpr (enableRegularization)
1145 const auto regularized = regularization_.pcgn(swe);
1147 return regularized.value();
1150 return BaseLaw::pcgn(swe, basicParams_);
1156 template<
bool enableRegularization = isRegularized()>
1159 const auto sne = EffToAbs::snToSne(sn, effToAbsParams_);
1160 if constexpr (enableRegularization)
1162 const auto regularized = regularization_.pcAlpha(sne);
1164 return regularized.value();
1166 return BaseLaw::pcAlpha(sne, basicParams_);
1172 template<
bool enableRegularization = isRegularized()>
1175 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1176 if constexpr (enableRegularization)
1178 const auto regularized = regularization_.dpcgw_dswe(swe);
1180 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
1183 return BaseLaw::dpcgw_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
1189 template<
bool enableRegularization = isRegularized()>
1192 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1193 if constexpr (enableRegularization)
1195 const auto regularized = regularization_.dpcnw_dswe(swe);
1197 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
1200 return BaseLaw::dpcnw_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
1206 template<
bool enableRegularization = isRegularized()>
1209 const auto ste = EffToAbs::stToSte(st, effToAbsParams_);
1210 if constexpr (enableRegularization)
1212 const auto regularized = regularization_.dpcgn_dste(ste);
1214 return regularized.value()*EffToAbs::dswte_dst(effToAbsParams_);
1217 return BaseLaw::dpcgn_dste(ste, basicParams_)*EffToAbs::dste_dst(effToAbsParams_);
1225 template<
bool enableRegularization = isRegularized()>
1228 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1229 if constexpr (enableRegularization)
1231 const auto regularized = regularization_.krw(swe);
1233 return regularized.value();
1236 return BaseLaw::krw(swe, basicParams_);
1244 template<
bool enableRegularization = isRegularized()>
1247 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1248 const auto ste = EffToAbs::stToSte(sw + sn, effToAbsParams_);
1249 if constexpr (enableRegularization)
1251 const auto regularized = regularization_.krn(swe, sn, ste);
1253 return regularized.value();
1256 return BaseLaw::krn(swe, sn, ste, basicParams_);
1264 template<
bool enableRegularization = isRegularized()>
1267 const auto ste = EffToAbs::stToSte(sw + sn, effToAbsParams_);
1268 if constexpr (enableRegularization)
1270 const auto regularized = regularization_.krg(ste);
1272 return regularized.value();
1275 return BaseLaw::krg(ste, basicParams_);
1284 template<
bool enableRegularization = isRegularized()>
1296 DUNE_THROW(Dune::InvalidStateException,
1297 "Invalid phase index ");
1304 template<
bool enableRegularization = isRegularized()>
1307 const auto ste = EffToAbs::stToSte(st, effToAbsParams_);
1308 if constexpr (enableRegularization)
1310 const auto regularized = regularization_.dkrg_dste(ste);
1312 return regularized.value()*EffToAbs::dste_dst(effToAbsParams_);
1315 return BaseLaw::dkrg_dste(ste, basicParams_)*EffToAbs::dste_dst(effToAbsParams_);
1323 return basicParams_ == o.basicParams_
1324 && effToAbsParams_ == o.effToAbsParams_
1325 && regularization_ == o.regularization_;
1334 return BaseLaw::template makeParams<Scalar>(paramGroup);
1341 {
return basicParams_; }
1349 return EffToAbs::template makeParams<Scalar>(paramGroup);
1356 {
return effToAbsParams_; }
1361 Regularization regularization_;
1368template<
class Scalar>
1375template<
class Scalar>
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.
A tag to turn off regularization and it's overhead.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: brookscorey.hh:286
A 3rd order polynomial spline.
Definition: spline.hh:55
Definition: parkervangenuchten.hh:37
static Scalar snToSne(const Scalar sn, const Params< Scalar > ¶ms)
Convert an absolute nonwetting saturation to an effective one.
Definition: parkervangenuchten.hh:181
static Scalar dswe_dsw(const Params< Scalar > ¶ms)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: parkervangenuchten.hh:152
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: parkervangenuchten.hh:104
static Scalar steToSt(const Scalar ste, const Params< Scalar > ¶ms)
Convert an effective wetting saturation to an absolute one.
Definition: parkervangenuchten.hh:211
static Scalar dste_dst(const Params< Scalar > ¶ms)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: parkervangenuchten.hh:225
static Scalar stToSte(const Scalar st, const Params< Scalar > ¶ms)
Convert an absolute total liquid saturation to an effective one.
Definition: parkervangenuchten.hh:196
static Scalar swToSwe(const Scalar sw, const Params< Scalar > ¶ms)
Convert an absolute wetting saturation to an effective one.
Definition: parkervangenuchten.hh:123
static Scalar dsw_dswe(const Params< Scalar > ¶ms)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: parkervangenuchten.hh:166
static Scalar sweToSw(const Scalar swe, const Params< Scalar > ¶ms)
Convert an effective wetting saturation to an absolute one.
Definition: parkervangenuchten.hh:138
static Scalar dst_dste(const Params< Scalar > ¶ms)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: parkervangenuchten.hh:239
The parameter type.
Definition: parkervangenuchten.hh:46
Scalar swr() const
Return the residual wetting saturation.
Definition: parkervangenuchten.hh:54
Params(const Scalar swr=0.0, const Scalar snr=0.0, const Scalar sgr=0.0)
Definition: parkervangenuchten.hh:47
void setSwr(Scalar v)
Set the residual wetting saturation.
Definition: parkervangenuchten.hh:60
Scalar snr() const
Return the residual nonwetting saturation.
Definition: parkervangenuchten.hh:66
void setSgr(Scalar v)
Set the residual gas phase saturation.
Definition: parkervangenuchten.hh:84
bool operator==(const Params &p) const
Definition: parkervangenuchten.hh:87
void setSnr(Scalar v)
Set the residual nonwetting saturation.
Definition: parkervangenuchten.hh:72
Scalar sgr() const
Return the residual gas phase saturation.
Definition: parkervangenuchten.hh:78
Implementation of Parker/vanGenuchten's capillary pressure <-> saturation relation for three phases....
Definition: parkervangenuchten.hh:253
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: parkervangenuchten.hh:329
static Scalar krn(const Scalar swe, const Scalar sn, const Scalar ste, const Params< Scalar > ¶ms)
The relative permeability for the nonwetting phase after the Model of Parker et al....
Definition: parkervangenuchten.hh:488
static Scalar pcnw(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve for the non-wettigng and wetting phase.
Definition: parkervangenuchten.hh:361
static Scalar krg(const Scalar ste, const Params< Scalar > ¶ms)
The relative permeability for the nonwetting phase of the medium implied by van Genuchten's parameter...
Definition: parkervangenuchten.hh:523
static Scalar dpcnw_dswe(const Scalar swe, const Params< Scalar > ¶ms)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: parkervangenuchten.hh:426
static Scalar kr(const int phaseIdx, const Scalar swe, const Scalar sne, const Params< Scalar > ¶ms)
The relative permeability for a phase.
Definition: parkervangenuchten.hh:561
static Scalar dkrg_dste(const Scalar ste, const Params< Scalar > ¶ms)
The derivative of the relative permeability for the gas phase in regard to the total liquid saturatio...
Definition: parkervangenuchten.hh:541
static Scalar pcgw(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: parkervangenuchten.hh:349
static Scalar dpcgn_dste(const Scalar ste, const Params< Scalar > ¶ms)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: parkervangenuchten.hh:441
static Scalar krw(const Scalar swe, const Params< Scalar > ¶ms)
The relative permeability for the wetting phase of the medium implied by van Genuchten's parameteriza...
Definition: parkervangenuchten.hh:462
static Scalar pcAlpha(Scalar sne, const Params< Scalar > ¶ms)
This function ensures a continuous transition from 2 to 3 phases and vice versa.
Definition: parkervangenuchten.hh:385
static Scalar dpcgw_dswe(const Scalar swe, const Params< Scalar > ¶ms)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: parkervangenuchten.hh:411
static Scalar pcgn(const Scalar ste, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve for the gas and nonwetting phase.
Definition: parkervangenuchten.hh:373
The parameter type.
Definition: parkervangenuchten.hh:272
Scalar betaGn() const
Definition: parkervangenuchten.hh:297
Scalar snr() const
Definition: parkervangenuchten.hh:291
Params(Scalar alpha, Scalar n, Scalar swr=0.0, Scalar snr=0.0, Scalar betaNw=1.0, Scalar betaGn=1.0, Scalar betaGw=1.0, bool regardSnr=false)
Definition: parkervangenuchten.hh:273
Scalar swr() const
Definition: parkervangenuchten.hh:288
void setAlpha(Scalar alpha)
Definition: parkervangenuchten.hh:280
void setSwr(Scalar swr)
Definition: parkervangenuchten.hh:289
void setM(Scalar m)
Definition: parkervangenuchten.hh:283
void setBetaGw(Scalar betaGw)
Definition: parkervangenuchten.hh:301
void setBetaNw(Scalar betaNw)
Definition: parkervangenuchten.hh:295
void setBetaGn(Scalar betaGn)
Definition: parkervangenuchten.hh:298
Scalar alpha() const
Definition: parkervangenuchten.hh:279
void setN(Scalar n)
Definition: parkervangenuchten.hh:286
bool regardSnrForKrn() const
Definition: parkervangenuchten.hh:303
Scalar betaNw() const
Definition: parkervangenuchten.hh:294
void setRegardSnrForKrn(bool v)
Definition: parkervangenuchten.hh:304
Scalar m() const
Definition: parkervangenuchten.hh:282
Scalar n() const
Definition: parkervangenuchten.hh:285
void setSnr(Scalar swr)
Definition: parkervangenuchten.hh:292
Scalar betaGw() const
Definition: parkervangenuchten.hh:300
bool operator==(const Params &p) const
Definition: parkervangenuchten.hh:306
A regularization for the ParkerVanGenuchten3PRegularization material law.
Definition: parkervangenuchten.hh:602
void init(const MaterialLaw *m, const BaseParams &bp, const EffToAbsParams &etap, const Params< Scalar > &p)
Definition: parkervangenuchten.hh:723
OptionalScalar< Scalar > pcgn(Scalar ste) const
The regularized capillary pressure-saturation curve for the gas and nonwetting phase regularized part...
Definition: parkervangenuchten.hh:822
OptionalScalar< Scalar > pcnw(Scalar swe) const
The regularized capillary pressure-saturation curve for the nonwetting and wetting phase regularized ...
Definition: parkervangenuchten.hh:789
OptionalScalar< Scalar > pcgw(Scalar swe) const
The regularized capillary pressure-saturation curve for the gas and wetting phase regularized part:
Definition: parkervangenuchten.hh:756
bool operator==(const ParkerVanGenuchten3PRegularization &o) const
Equality comparison with another instance.
Definition: parkervangenuchten.hh:739
OptionalScalar< Scalar > kr(const int phaseIdx, const Scalar swe, const Scalar sne) const
The relative permeability for a phase.
Definition: parkervangenuchten.hh:938
void init(const MaterialLaw *m, const std::string ¶mGroup)
Initialize the spline.
Definition: parkervangenuchten.hh:709
OptionalScalar< Scalar > krn(Scalar swe, const Scalar sn, Scalar ste) const
The regularized relative permeability for the nonwetting phase.
Definition: parkervangenuchten.hh:881
OptionalScalar< Scalar > krg(const Scalar ste) const
The regularized relative permeability for the gas phase.
Definition: parkervangenuchten.hh:897
OptionalScalar< Scalar > pcAlpha(Scalar sne) const
This function ensures a continuous transition from 2 to 3 phases and vice versa.
Definition: parkervangenuchten.hh:854
OptionalScalar< Scalar > krw(const Scalar swe) const
The regularized relative permeability for the wetting phase.
Definition: parkervangenuchten.hh:864
Regularization parameters.
Definition: parkervangenuchten.hh:609
Scalar krnLowSwe() const
Threshold saturation below which the relative permeability of the nonwetting phase gets regularized.
Definition: parkervangenuchten.hh:651
void setKrnLowSwe(Scalar krnLowSwe)
Set the threshold saturation below which the relative permeability of the nonwetting phase gets regul...
Definition: parkervangenuchten.hh:644
Scalar pcLowSwe() const
Threshold saturation below which the capillary pressure is regularized.
Definition: parkervangenuchten.hh:622
void setPcLowSwe(Scalar pcLowSwe)
Set the threshold saturation below which the capillary pressure is regularized.
Definition: parkervangenuchten.hh:616
Scalar pcHighSwe() const
Threshold saturation above which the capillary pressure is regularized.
Definition: parkervangenuchten.hh:637
void setConstRegularization(const bool input)
Choose whether to use a constant value for regularization of the pc-S curves or not.
Definition: parkervangenuchten.hh:688
void setKrwHighSwe(Scalar krwHighSwe)
Set the threshold saturation above which the relative permeability of the wetting phase gets regulari...
Definition: parkervangenuchten.hh:672
bool constRegularization() const
Returns whether to use a constant value for regularization of the pc-S curves or not.
Definition: parkervangenuchten.hh:695
Scalar krwHighSwe() const
Threshold saturation above which the relative permeability of the wetting phase gets regularized.
Definition: parkervangenuchten.hh:679
void setPcHighSwe(Scalar pcHighSwe)
Set the threshold saturation above which the capillary pressure is regularized.
Definition: parkervangenuchten.hh:628
void setKrgLowSte(Scalar krgLowSte)
Set the threshold saturation below which the relative permeability of the nonwetting phase gets regul...
Definition: parkervangenuchten.hh:658
Scalar krgLowSte() const
Threshold saturation below which the relative permeability of the nonwetting phase gets regularized.
Definition: parkervangenuchten.hh:665
Definition: parkervangenuchten.hh:1051
typename Regularization::template Params< Scalar > RegularizationParams
Definition: parkervangenuchten.hh:1058
Scalar kr(const int phaseIdx, const Scalar sw, const Scalar sn) const
The relative permeability for the nonwetting phase.
Definition: parkervangenuchten.hh:1285
ParkerVanGenuchtenMaterialLaw(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: parkervangenuchten.hh:1078
Scalar pcnw(const Scalar sw, const Scalar) const
The capillary pressure-saturation curve for the nonwetting and wetting phase.
Definition: parkervangenuchten.hh:1121
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: parkervangenuchten.hh:1057
const EffToAbsParams & effToAbsParams() const
Return the parameters of the EffToAbs policy.
Definition: parkervangenuchten.hh:1355
ParkerVanGenuchtenMaterialLaw()=delete
Deleted default constructor (so we are never in an undefined state)
Scalar dpcnw_dsw(const Scalar sw, const Scalar) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: parkervangenuchten.hh:1190
ScalarType Scalar
Definition: parkervangenuchten.hh:1054
Scalar krn(const Scalar sw, const Scalar sn) const
The relative permeability for the nonwetting phase.
Definition: parkervangenuchten.hh:1245
static BasicParams makeBasicParams(const std::string ¶mGroup)
Create the base law's parameters using input file parameters.
Definition: parkervangenuchten.hh:1332
Scalar pcgn(const Scalar sw, const Scalar sn) const
The capillary pressure-saturation curve for the gas and nonwetting phase.
Definition: parkervangenuchten.hh:1140
Scalar krw(const Scalar sw, const Scalar sn) const
The relative permeability for the wetting phase.
Definition: parkervangenuchten.hh:1226
Scalar dpcgw_dsw(const Scalar sw, const Scalar) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: parkervangenuchten.hh:1173
Scalar krg(const Scalar sw, const Scalar sn) const
The relative permeability for the nonwetting phase.
Definition: parkervangenuchten.hh:1265
Scalar pcgw(const Scalar sw, const Scalar) const
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: parkervangenuchten.hh:1104
static EffToAbsParams makeEffToAbsParams(const std::string ¶mGroup)
Create the parameters of the EffToAbs policy using input file parameters.
Definition: parkervangenuchten.hh:1347
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: parkervangenuchten.hh:1340
Scalar dkrg_dst(const Scalar st) const
The derivative of the relative permeability for the nonwetting phase w.r.t. saturation.
Definition: parkervangenuchten.hh:1305
ParkerVanGenuchtenMaterialLaw(const BasicParams &baseParams, const EffToAbsParams &effToAbsParams={}, const RegularizationParams ®Params={})
Construct from parameter structs.
Definition: parkervangenuchten.hh:1090
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: parkervangenuchten.hh:1065
bool operator==(const ParkerVanGenuchtenMaterialLaw &o) const
Equality comparison with another instance.
Definition: parkervangenuchten.hh:1321
typename BaseLaw::template Params< Scalar > BasicParams
Definition: parkervangenuchten.hh:1056
Scalar pcAlpha(const Scalar, const Scalar sn) const
This function ensures a continuous transition from 2 to 3 phases and vice versa.
Definition: parkervangenuchten.hh:1157
Scalar dpcgn_dst(const Scalar st, const Scalar) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: parkervangenuchten.hh:1207
EffToAbsPolicy EffToAbs
Definition: parkervangenuchten.hh:1060
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67