12#ifndef PARKER_VANGENUCHTEN_3P_HH
13#define PARKER_VANGENUCHTEN_3P_HH
32 template<
class Scalar>
77 return Dune::FloatCmp::eq(
swr(), p.
swr(), 1e-6)
78 && Dune::FloatCmp::eq(
snr(), p.
snr(), 1e-6)
79 && Dune::FloatCmp::eq(
sgr(), p.
sgr(), 1e-6);
91 template<
class Scalar>
95 params.
setSwr(getParamFromGroup<Scalar>(paramGroup,
"Swr", 0.0));
96 params.
setSnr(getParamFromGroup<Scalar>(paramGroup,
"Snr", 0.0));
97 params.
setSgr(getParamFromGroup<Scalar>(paramGroup,
"Sgr", 0.0));
110 template<
class Scalar>
113 return (sw - params.
swr())/(1.0 - params.
swr());
125 template<
class Scalar>
128 return swe*(1.0 - params.
swr()) + params.
swr();
139 template<
class Scalar>
142 return 1.0/(1.0 - params.
swr());
153 template<
class Scalar>
156 return 1.0 - params.
swr();
168 template<
class Scalar>
183 template<
class Scalar>
186 return (st-params.
swr()) / (1-params.
swr());
198 template<
class Scalar>
201 return ste*(1.0 - params.
swr()) + params.
swr();
212 template<
class Scalar>
215 return 1.0/(1.0 - params.
swr() );
226 template<
class Scalar>
229 return 1.0 - params.
swr();
258 template<
class Scalar>
262 Scalar
betaNw = 1.0, Scalar
betaGn = 1.0, Scalar
betaGw = 1.0,
bool regardSnr =
false)
263 : alpha_(
alpha), n_(
n), m_(1.0 - 1.0/
n), swr_(
swr), snr_(
snr)
267 Scalar
alpha()
const {
return alpha_; }
270 Scalar
m()
const {
return m_; }
271 void setM(Scalar
m) { m_ =
m; n_ = 1.0/(1.0 -
m); }
273 Scalar
n()
const{
return n_; }
274 void setN(Scalar
n){ n_ =
n; m_ = 1.0 - 1.0/
n; }
276 Scalar
swr()
const {
return swr_; }
279 Scalar
snr()
const {
return snr_; }
282 Scalar
betaNw()
const {
return betaNw_; }
285 Scalar
betaGn()
const {
return betaGn_; }
288 Scalar
betaGw()
const {
return betaGw_; }
296 return Dune::FloatCmp::eq(alpha_, p.alpha_, 1e-6)
297 && Dune::FloatCmp::eq(n_, p.n_, 1e-6)
298 && Dune::FloatCmp::eq(m_, p.m_, 1e-6)
299 && Dune::FloatCmp::eq(swr_, p.swr_, 1e-6)
300 && Dune::FloatCmp::eq(snr_, p.snr_, 1e-6)
301 && Dune::FloatCmp::eq(betaNw_, p.betaNw_, 1e-6)
302 && Dune::FloatCmp::eq(betaGn_, p.betaGn_, 1e-6)
303 && Dune::FloatCmp::eq(betaGw_, p.betaGw_, 1e-6)
304 && regardSnr_ == p.regardSnr_;
308 Scalar alpha_, n_, m_, swr_, snr_, betaNw_, betaGn_, betaGw_;
316 template<
class Scalar =
double>
319 const auto n = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenN");
320 const auto alpha = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenAlpha");
321 const auto swr = getParamFromGroup<Scalar>(paramGroup,
"Swr", 0.0);
322 const auto snr = getParamFromGroup<Scalar>(paramGroup,
"Snr", 0.0);
323 const auto betaNw = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenBetaNw", 1.0);
324 const auto betaGn = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenBetaGn", 1.0);
325 const auto betaGw = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenBetaGw", 1.0);
326 const auto regardSnr = getParamFromGroup<bool>(paramGroup,
"ParkerVanGenuchtenRegardSnrForKrn",
false);
328 betaNw, betaGn, betaGw, regardSnr );
336 template<
class Scalar>
339 assert(0 <= swe && swe <= 1);
340 return pc_(swe, params);
348 template<
class Scalar>
351 assert(0 <= swe && swe <= 1);
352 return pc_(swe, params)/params.
betaNw();
360 template<
class Scalar>
363 assert(0 <= ste && ste <= 1);
364 return pc_(ste, params)/params.
betaGn();
372 template<
class Scalar>
381 if (sne > params.
snr())
385 if (params.
snr() >= 0.001)
386 return sne/params.
snr();
398 template<
class Scalar>
402 const Scalar powSeRegu = pow(swe, -1/params.
m());
403 return - 1.0/params.
alpha() * pow(powSeRegu - 1, 1.0/params.
n() - 1)/params.
n()
404 * powSeRegu/swe/params.
m()/params.
betaGw();
413 template<
class Scalar>
417 const Scalar powSeRegu = pow(swe, -1/params.
m());
418 return - 1.0/params.
alpha() * pow(powSeRegu - 1, 1.0/params.
n() - 1)/params.
n()
419 * powSeRegu/swe/params.
m()/params.
betaNw();
428 template<
class Scalar>
432 const Scalar powSeRegu = pow(ste, -1/params.
m());
433 return - 1.0/params.
alpha() * pow(powSeRegu - 1, 1.0/params.
n() - 1)/params.
n()
434 * powSeRegu/ste/params.
m()/params.
betaGn();
449 template<
class Scalar>
454 const Scalar r = 1.0 - pow(1 - pow(swe, 1/params.
m()), params.
m());
455 return sqrt(swe)*r*r;
475 template<
class Scalar>
476 static Scalar
krn(
const Scalar swe,
const Scalar sn,
const Scalar ste,
const Params<Scalar>& params)
480 krn = pow(1 - pow(swe, 1/params.
m()), params.
m());
481 krn -= pow(1 - pow(ste, 1/params.
m()), params.
m());
489 const Scalar resIncluded = clamp(sn - params.
snr()/ (1-params.
swr()), 0.0, 1.0);
490 krn *= sqrt(resIncluded);
493 krn *= sqrt(sn / (1 - params.
swr()));
510 template<
class Scalar>
513 assert(0 <= ste && ste <= 1);
516 return cbrt(1 - ste) * pow(1 - pow(ste, 1/params.
m()), 2*params.
m());
528 template<
class Scalar>
531 assert(0 < ste && ste <= 1);
534 const Scalar x = pow(ste, 1.0/params.
m());
536 -pow(1.0 - x, 2*params.
m())
537 *pow(1.0 - ste, -2.0/3)
548 template<
class Scalar>
549 static Scalar
kr(
const int phaseIdx,
const Scalar swe,
const Scalar sne,
const Params<Scalar>& params)
554 return krw(params, swe, sne);
556 return krn(params, swe, sne);
558 return krg(params, swe, sne);
560 DUNE_THROW(Dune::InvalidStateException,
561 "Invalid phase index ");
572 template<
class Scalar>
573 const static Scalar pc_(
const Scalar se,
const Params<Scalar>& params)
576 return pow(pow(se, -1/params.m()) - 1, 1/params.n())/params.alpha();
588template <
class Scalar>
611 {
return pcLowSwe_; }
626 {
return pcHighSwe_; }
640 {
return krnLowSwe_; }
654 {
return krgLowSte_; }
668 {
return krwHighSwe_; }
677 { constRegularization_ = input; }
684 {
return constRegularization_; }
692 bool constRegularization_ =
false;
696 template<
class MaterialLaw>
697 void init(
const MaterialLaw* m,
const std::string& paramGroup)
699 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenPcLowSweThreshold", 0.01);
700 pcHighSwe_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenPcHighSweThreshold", 0.99);
701 krwHighSwe_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenKrwHighSweThreshold", 0.9);
702 krnLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenKrnLowSweThreshold", 0.1);
703 krgLowSte_ = getParamFromGroup<Scalar>(paramGroup,
"ParkerVanGenuchtenKrgLowSteThreshold", 1e-3);
704 constRegularization_ = getParamFromGroup<bool>(paramGroup,
"VanGenuchtenConstantRegularization",
false);
706 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
707 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
710 template<
class MaterialLaw,
class BaseParams,
class EffToAbsParams>
711 void init(
const MaterialLaw* m,
const BaseParams& bp,
const EffToAbsParams& etap,
const Params<Scalar>& p)
720 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
721 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
729 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
730 && Dune::FloatCmp::eq(pcHighSwe_, o.pcHighSwe_, 1e-6)
731 && Dune::FloatCmp::eq(krwHighSwe_, o.krwHighSwe_, 1e-6)
732 && Dune::FloatCmp::eq(krnLowSwe_, o.krnLowSwe_, 1e-6)
733 && constRegularization_ == o.constRegularization_;
748 if (constRegularization_)
749 swe = clamp(swe, 0.0, 1.0);
757 return pcgwLowSwePcgwValue_ + pcgwDerivativeLowSw_*(swe - pcLowSwe_);
760 return pcgwDerivativeHighSweEnd_*(swe - 1.0);
762 else if (swe > pcHighSwe_)
763 return pcgwSpline_.eval(swe);
781 if (constRegularization_)
782 swe = clamp(swe, 0.0, 1.0);
790 return pcnwLowSwePcnwValue_ + pcnwDerivativeLowSw_*(swe - pcLowSwe_);
793 return pcnwDerivativeHighSweEnd_*(swe - 1.0);
795 else if (swe > pcHighSwe_)
796 return pcnwSpline_.eval(swe);
814 if (constRegularization_)
815 ste = clamp(ste, 0.0, 1.0);
823 const Scalar pcLowSte = pcLowSwe_;
824 const Scalar pcHighSte = pcHighSwe_;
826 return pcgnLowStePcgnValue_ + pcgnDerivativeLowSt_*(ste - pcLowSte);
829 return pcgnDerivativeHighSteEnd_*(ste - 1.0);
831 else if (ste > pcHighSte)
832 return pcgnSpline_.eval(ste);
856 else if (swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
872 swe = clamp(swe, 0.0, 1.0);
873 ste = clamp(ste, 0.0, 1.0);
875 if (ste - swe <= 0.0)
888 if (ste > 1.0 - std::numeric_limits<Scalar>::epsilon())
894 if (ste <= krgLowSte_)
895 return krgLowStkrgValue_ + krgDerivativeLowSt_*(ste - krgLowSte_);
905 const Scalar st = ste*(1 - swr_) + swr_;
906 const Scalar sg = 1.0 - st;
914 const Scalar scalFact = max(0.0, (sg - sgr_)/(0.1 - sgr_));
931 return krw(swe, sne);
933 return krn(swe, sne);
935 return krg(swe, sne);
937 DUNE_THROW(Dune::InvalidStateException,
938 "Invalid phase index ");
942 template<
class MaterialLaw>
943 void initPcParameters_(
const MaterialLaw* m,
const Scalar lowSwe,
const Scalar highSwe)
945 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
946 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
947 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
948 const auto dst_dste = MaterialLaw::EffToAbs::dst_dste(m->effToAbsParams());
950 baseLawParamsPtr_ = &m->basicParams();
953 pcgwLowSwePcgwValue_ = m->template pcgw<false>(lowSw, 0.0);
954 pcgwDerivativeLowSw_ = m->template dpcgw_dsw<false>(lowSw, 0.0)*dsw_dswe;
955 pcgwHighSwePcgwValue_ = m->template pcgw<false>(highSw, 0.0);
956 pcgwDerivativeHighSweThreshold_ = m->template dpcgw_dsw<false>(highSw, 0.0)*dsw_dswe;
957 pcgwDerivativeHighSweEnd_ = 2.0*(0.0 - m->template pcgw<false>(highSw, 0.0))/(1.0 - highSwe);
959 pcgwHighSwePcgwValue_, 0,
960 pcgwDerivativeHighSweThreshold_, pcgwDerivativeHighSweEnd_);
963 pcnwLowSwePcnwValue_ = m->template pcnw<false>(lowSw, 0.0);
964 pcnwDerivativeLowSw_ = m->template dpcnw_dsw<false>(lowSw, 0.0)*dsw_dswe;
965 pcnwHighSwePcnwValue_ = m->template pcnw<false>(highSw, 0.0);
966 pcnwDerivativeHighSweThreshold_ = m->template dpcnw_dsw<false>(highSw, 0.0);
967 pcnwDerivativeHighSweEnd_ = 2.0*(0.0 - m->template pcnw<false>(highSw, 0.0))/(1.0 - highSwe);
969 pcnwHighSwePcnwValue_, 0,
970 pcnwDerivativeHighSweThreshold_, pcnwDerivativeHighSweEnd_);
973 pcgnLowStePcgnValue_ = m->template pcgn<false>(lowSw, 0.0);
974 pcgnDerivativeLowSt_ = m->template dpcgn_dst<false>(lowSw, 0.0)*dst_dste;
975 pcgnHighSwePcgnValue_ = m->template pcgn<false>(highSw, 0.0);
976 pcgnDerivativeHighSteThreshold_ = m->template dpcgn_dst<false>(highSw, 0.0);
977 pcgnDerivativeHighSteEnd_ = 2.0*(0.0 - m->template pcgn<false>(highSw, 0.0))/(1.0 - highSwe);
979 pcgnHighSwePcgnValue_, 0,
980 pcgnDerivativeHighSteThreshold_, pcgnDerivativeHighSteEnd_);
984 template<
class MaterialLaw>
985 void initKrParameters_(
const MaterialLaw* m,
const Scalar lowSwe,
const Scalar highSwe)
990 swr_ = m->effToAbsParams().swr();
991 sgr_ = m->effToAbsParams().sgr();
994 Scalar krgLowStkrgValue_;
995 Scalar krgDerivativeLowSt_;
997 Scalar pcLowSwe_, pcHighSwe_;
998 Scalar krwHighSwe_, krnLowSwe_, krgLowSte_;
1001 Scalar pcgwLowSwePcgwValue_;
1002 Scalar pcgwHighSwePcgwValue_;
1003 Scalar pcgwDerivativeLowSw_;
1004 Scalar pcgwDerivativeHighSweThreshold_;
1005 Scalar pcgwDerivativeHighSweEnd_;
1008 Scalar pcgnLowStePcgnValue_;
1009 Scalar pcgnHighSwePcgnValue_;
1010 Scalar pcgnDerivativeLowSt_;
1011 Scalar pcgnDerivativeHighSteThreshold_;
1012 Scalar pcgnDerivativeHighSteEnd_;
1015 Scalar pcnwLowSwePcnwValue_;
1016 Scalar pcnwHighSwePcnwValue_;
1017 Scalar pcnwDerivativeLowSw_;
1018 Scalar pcnwDerivativeHighSweThreshold_;
1019 Scalar pcnwDerivativeHighSweEnd_;
1021 Spline<Scalar> pcgwSpline_;
1022 Spline<Scalar> pcnwSpline_;
1023 Spline<Scalar> pcgnSpline_;
1024 Spline<Scalar> krwSpline_;
1025 Spline<Scalar> krnSpline_;
1029 bool constRegularization_;
1031 const BaseLawParams* baseLawParamsPtr_;
1038template<
class ScalarType,
1040 class Regularization = NoRegularization,
1041 class EffToAbsPolicy = ParkerVanGenuchten3PEffToAbsPolicy>
1058 {
return !std::is_same<Regularization, NoRegularization>::value; }
1075 regularization_.init(
this, paramGroup);
1085 : basicParams_(baseParams)
1089 regularization_.init(
this, baseParams,
effToAbsParams, regParams);
1095 template<
bool enableRegularization = isRegularized()>
1098 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1099 if constexpr (enableRegularization)
1101 const auto regularized = regularization_.pcgw(swe);
1103 return regularized.value();
1106 return BaseLaw::pcgw(swe, basicParams_);
1112 template<
bool enableRegularization = isRegularized()>
1115 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1116 if constexpr (enableRegularization)
1118 const auto regularized = regularization_.pcnw(swe);
1120 return regularized.value();
1123 return BaseLaw::pcnw(swe, basicParams_);
1131 template<
bool enableRegularization = isRegularized()>
1134 const auto swe = EffToAbs::swToSwe(sw + sn, effToAbsParams_);
1135 if constexpr (enableRegularization)
1137 const auto regularized = regularization_.pcgn(swe);
1139 return regularized.value();
1142 return BaseLaw::pcgn(swe, basicParams_);
1148 template<
bool enableRegularization = isRegularized()>
1151 const auto sne = EffToAbs::snToSne(sn, effToAbsParams_);
1152 if constexpr (enableRegularization)
1154 const auto regularized = regularization_.pcAlpha(sne);
1156 return regularized.value();
1158 return BaseLaw::pcAlpha(sne, basicParams_);
1164 template<
bool enableRegularization = isRegularized()>
1167 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1168 if constexpr (enableRegularization)
1170 const auto regularized = regularization_.dpcgw_dswe(swe);
1172 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
1175 return BaseLaw::dpcgw_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
1181 template<
bool enableRegularization = isRegularized()>
1184 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1185 if constexpr (enableRegularization)
1187 const auto regularized = regularization_.dpcnw_dswe(swe);
1189 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
1192 return BaseLaw::dpcnw_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
1198 template<
bool enableRegularization = isRegularized()>
1201 const auto ste = EffToAbs::stToSte(st, effToAbsParams_);
1202 if constexpr (enableRegularization)
1204 const auto regularized = regularization_.dpcgn_dste(ste);
1206 return regularized.value()*EffToAbs::dswte_dst(effToAbsParams_);
1209 return BaseLaw::dpcgn_dste(ste, basicParams_)*EffToAbs::dste_dst(effToAbsParams_);
1217 template<
bool enableRegularization = isRegularized()>
1220 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1221 if constexpr (enableRegularization)
1223 const auto regularized = regularization_.krw(swe);
1225 return regularized.value();
1228 return BaseLaw::krw(swe, basicParams_);
1236 template<
bool enableRegularization = isRegularized()>
1239 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1240 const auto ste = EffToAbs::stToSte(sw + sn, effToAbsParams_);
1241 if constexpr (enableRegularization)
1243 const auto regularized = regularization_.krn(swe, sn, ste);
1245 return regularized.value();
1248 return BaseLaw::krn(swe, sn, ste, basicParams_);
1256 template<
bool enableRegularization = isRegularized()>
1259 const auto ste = EffToAbs::stToSte(sw + sn, effToAbsParams_);
1260 if constexpr (enableRegularization)
1262 const auto regularized = regularization_.krg(ste);
1264 return regularized.value();
1267 return BaseLaw::krg(ste, basicParams_);
1276 template<
bool enableRegularization = isRegularized()>
1288 DUNE_THROW(Dune::InvalidStateException,
1289 "Invalid phase index ");
1296 template<
bool enableRegularization = isRegularized()>
1299 const auto ste = EffToAbs::stToSte(st, effToAbsParams_);
1300 if constexpr (enableRegularization)
1302 const auto regularized = regularization_.dkrg_dste(ste);
1304 return regularized.value()*EffToAbs::dste_dst(effToAbsParams_);
1307 return BaseLaw::dkrg_dste(ste, basicParams_)*EffToAbs::dste_dst(effToAbsParams_);
1315 return basicParams_ == o.basicParams_
1316 && effToAbsParams_ == o.effToAbsParams_
1317 && regularization_ == o.regularization_;
1326 return BaseLaw::template makeParams<Scalar>(paramGroup);
1333 {
return basicParams_; }
1341 return EffToAbs::template makeParams<Scalar>(paramGroup);
1348 {
return effToAbsParams_; }
1353 Regularization regularization_;
1360template<
class Scalar>
1367template<
class Scalar>
Implementation of Parker/vanGenuchten's capillary pressure <-> saturation relation for three phases....
Definition: parkervangenuchten.hh:241
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: parkervangenuchten.hh:317
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:476
static Scalar pcnw(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve for the non-wettigng and wetting phase.
Definition: parkervangenuchten.hh:349
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:511
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:414
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:549
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:529
static Scalar pcgw(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: parkervangenuchten.hh:337
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:429
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:450
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:373
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:399
static Scalar pcgn(const Scalar ste, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve for the gas and nonwetting phase.
Definition: parkervangenuchten.hh:361
A regularization for the ParkerVanGenuchten3PRegularization material law.
Definition: parkervangenuchten.hh:590
void init(const MaterialLaw *m, const BaseParams &bp, const EffToAbsParams &etap, const Params< Scalar > &p)
Definition: parkervangenuchten.hh:711
OptionalScalar< Scalar > pcgn(Scalar ste) const
The regularized capillary pressure-saturation curve for the gas and nonwetting phase regularized part...
Definition: parkervangenuchten.hh:810
OptionalScalar< Scalar > pcnw(Scalar swe) const
The regularized capillary pressure-saturation curve for the nonwetting and wetting phase regularized ...
Definition: parkervangenuchten.hh:777
OptionalScalar< Scalar > pcgw(Scalar swe) const
The regularized capillary pressure-saturation curve for the gas and wetting phase regularized part:
Definition: parkervangenuchten.hh:744
bool operator==(const ParkerVanGenuchten3PRegularization &o) const
Equality comparison with another instance.
Definition: parkervangenuchten.hh:727
OptionalScalar< Scalar > kr(const int phaseIdx, const Scalar swe, const Scalar sne) const
The relative permeability for a phase.
Definition: parkervangenuchten.hh:926
void init(const MaterialLaw *m, const std::string ¶mGroup)
Initialize the spline.
Definition: parkervangenuchten.hh:697
OptionalScalar< Scalar > krn(Scalar swe, const Scalar sn, Scalar ste) const
The regularized relative permeability for the nonwetting phase.
Definition: parkervangenuchten.hh:869
OptionalScalar< Scalar > krg(const Scalar ste) const
The regularized relative permeability for the gas phase.
Definition: parkervangenuchten.hh:885
OptionalScalar< Scalar > pcAlpha(Scalar sne) const
This function ensures a continuous transition from 2 to 3 phases and vice versa.
Definition: parkervangenuchten.hh:842
OptionalScalar< Scalar > krw(const Scalar swe) const
The regularized relative permeability for the wetting phase.
Definition: parkervangenuchten.hh:852
Parker van Genuchten material law.
Definition: parkervangenuchten.hh:1043
typename Regularization::template Params< Scalar > RegularizationParams
Definition: parkervangenuchten.hh:1050
Scalar kr(const int phaseIdx, const Scalar sw, const Scalar sn) const
The relative permeability for the nonwetting phase.
Definition: parkervangenuchten.hh:1277
ParkerVanGenuchtenMaterialLaw(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: parkervangenuchten.hh:1070
Scalar pcnw(const Scalar sw, const Scalar) const
The capillary pressure-saturation curve for the nonwetting and wetting phase.
Definition: parkervangenuchten.hh:1113
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: parkervangenuchten.hh:1049
const EffToAbsParams & effToAbsParams() const
Return the parameters of the EffToAbs policy.
Definition: parkervangenuchten.hh:1347
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:1182
ScalarType Scalar
Definition: parkervangenuchten.hh:1046
Scalar krn(const Scalar sw, const Scalar sn) const
The relative permeability for the nonwetting phase.
Definition: parkervangenuchten.hh:1237
static BasicParams makeBasicParams(const std::string ¶mGroup)
Create the base law's parameters using input file parameters.
Definition: parkervangenuchten.hh:1324
Scalar pcgn(const Scalar sw, const Scalar sn) const
The capillary pressure-saturation curve for the gas and nonwetting phase.
Definition: parkervangenuchten.hh:1132
Scalar krw(const Scalar sw, const Scalar sn) const
The relative permeability for the wetting phase.
Definition: parkervangenuchten.hh:1218
Scalar dpcgw_dsw(const Scalar sw, const Scalar) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: parkervangenuchten.hh:1165
Scalar krg(const Scalar sw, const Scalar sn) const
The relative permeability for the nonwetting phase.
Definition: parkervangenuchten.hh:1257
Scalar pcgw(const Scalar sw, const Scalar) const
The capillary pressure-saturation curve for the gas and wetting phase.
Definition: parkervangenuchten.hh:1096
static EffToAbsParams makeEffToAbsParams(const std::string ¶mGroup)
Create the parameters of the EffToAbs policy using input file parameters.
Definition: parkervangenuchten.hh:1339
const BasicParams & basicParams() const
Return the base law's parameters.
Definition: parkervangenuchten.hh:1332
Scalar dkrg_dst(const Scalar st) const
The derivative of the relative permeability for the nonwetting phase w.r.t. saturation.
Definition: parkervangenuchten.hh:1297
ParkerVanGenuchtenMaterialLaw(const BasicParams &baseParams, const EffToAbsParams &effToAbsParams={}, const RegularizationParams ®Params={})
Construct from parameter structs.
Definition: parkervangenuchten.hh:1082
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: parkervangenuchten.hh:1057
bool operator==(const ParkerVanGenuchtenMaterialLaw &o) const
Equality comparison with another instance.
Definition: parkervangenuchten.hh:1313
typename BaseLaw::template Params< Scalar > BasicParams
Definition: parkervangenuchten.hh:1048
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:1149
Scalar dpcgn_dst(const Scalar st, const Scalar) const
The partial derivative of the capillary pressure w.r.t. the saturation.
Definition: parkervangenuchten.hh:1199
EffToAbsPolicy EffToAbs
Definition: parkervangenuchten.hh:1052
A 3rd order polynomial spline.
Definition: spline.hh:43
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: brookscorey.hh:23
A tag to turn off regularization and it's overhead.
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.
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55
The parameter type.
Definition: parkervangenuchten.hh:260
Scalar betaGn() const
Definition: parkervangenuchten.hh:285
Scalar snr() const
Definition: parkervangenuchten.hh:279
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:261
Scalar swr() const
Definition: parkervangenuchten.hh:276
void setAlpha(Scalar alpha)
Definition: parkervangenuchten.hh:268
void setSwr(Scalar swr)
Definition: parkervangenuchten.hh:277
void setM(Scalar m)
Definition: parkervangenuchten.hh:271
void setBetaGw(Scalar betaGw)
Definition: parkervangenuchten.hh:289
void setBetaNw(Scalar betaNw)
Definition: parkervangenuchten.hh:283
void setBetaGn(Scalar betaGn)
Definition: parkervangenuchten.hh:286
Scalar alpha() const
Definition: parkervangenuchten.hh:267
void setN(Scalar n)
Definition: parkervangenuchten.hh:274
bool regardSnrForKrn() const
Definition: parkervangenuchten.hh:291
Scalar betaNw() const
Definition: parkervangenuchten.hh:282
void setRegardSnrForKrn(bool v)
Definition: parkervangenuchten.hh:292
Scalar m() const
Definition: parkervangenuchten.hh:270
Scalar n() const
Definition: parkervangenuchten.hh:273
void setSnr(Scalar swr)
Definition: parkervangenuchten.hh:280
Scalar betaGw() const
Definition: parkervangenuchten.hh:288
bool operator==(const Params &p) const
Definition: parkervangenuchten.hh:294
The parameter type.
Definition: parkervangenuchten.hh:34
Scalar swr() const
Return the residual wetting saturation.
Definition: parkervangenuchten.hh:42
Params(const Scalar swr=0.0, const Scalar snr=0.0, const Scalar sgr=0.0)
Definition: parkervangenuchten.hh:35
void setSwr(Scalar v)
Set the residual wetting saturation.
Definition: parkervangenuchten.hh:48
Scalar snr() const
Return the residual nonwetting saturation.
Definition: parkervangenuchten.hh:54
void setSgr(Scalar v)
Set the residual gas phase saturation.
Definition: parkervangenuchten.hh:72
bool operator==(const Params &p) const
Definition: parkervangenuchten.hh:75
void setSnr(Scalar v)
Set the residual nonwetting saturation.
Definition: parkervangenuchten.hh:60
Scalar sgr() const
Return the residual gas phase saturation.
Definition: parkervangenuchten.hh:66
Definition: parkervangenuchten.hh:25
static Scalar snToSne(const Scalar sn, const Params< Scalar > ¶ms)
Convert an absolute nonwetting saturation to an effective one.
Definition: parkervangenuchten.hh:169
static Scalar dswe_dsw(const Params< Scalar > ¶ms)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: parkervangenuchten.hh:140
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: parkervangenuchten.hh:92
static Scalar steToSt(const Scalar ste, const Params< Scalar > ¶ms)
Convert an effective wetting saturation to an absolute one.
Definition: parkervangenuchten.hh:199
static Scalar dste_dst(const Params< Scalar > ¶ms)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: parkervangenuchten.hh:213
static Scalar stToSte(const Scalar st, const Params< Scalar > ¶ms)
Convert an absolute total liquid saturation to an effective one.
Definition: parkervangenuchten.hh:184
static Scalar swToSwe(const Scalar sw, const Params< Scalar > ¶ms)
Convert an absolute wetting saturation to an effective one.
Definition: parkervangenuchten.hh:111
static Scalar dsw_dswe(const Params< Scalar > ¶ms)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: parkervangenuchten.hh:154
static Scalar sweToSw(const Scalar swe, const Params< Scalar > ¶ms)
Convert an effective wetting saturation to an absolute one.
Definition: parkervangenuchten.hh:126
static Scalar dst_dste(const Params< Scalar > ¶ms)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: parkervangenuchten.hh:227
Regularization parameters.
Definition: parkervangenuchten.hh:597
Scalar krnLowSwe() const
Threshold saturation below which the relative permeability of the nonwetting phase gets regularized.
Definition: parkervangenuchten.hh:639
void setKrnLowSwe(Scalar krnLowSwe)
Set the threshold saturation below which the relative permeability of the nonwetting phase gets regul...
Definition: parkervangenuchten.hh:632
Scalar pcLowSwe() const
Threshold saturation below which the capillary pressure is regularized.
Definition: parkervangenuchten.hh:610
void setPcLowSwe(Scalar pcLowSwe)
Set the threshold saturation below which the capillary pressure is regularized.
Definition: parkervangenuchten.hh:604
Scalar pcHighSwe() const
Threshold saturation above which the capillary pressure is regularized.
Definition: parkervangenuchten.hh:625
void setConstRegularization(const bool input)
Choose whether to use a constant value for regularization of the pc-S curves or not.
Definition: parkervangenuchten.hh:676
void setKrwHighSwe(Scalar krwHighSwe)
Set the threshold saturation above which the relative permeability of the wetting phase gets regulari...
Definition: parkervangenuchten.hh:660
bool constRegularization() const
Returns whether to use a constant value for regularization of the pc-S curves or not.
Definition: parkervangenuchten.hh:683
Scalar krwHighSwe() const
Threshold saturation above which the relative permeability of the wetting phase gets regularized.
Definition: parkervangenuchten.hh:667
void setPcHighSwe(Scalar pcHighSwe)
Set the threshold saturation above which the capillary pressure is regularized.
Definition: parkervangenuchten.hh:616
void setKrgLowSte(Scalar krgLowSte)
Set the threshold saturation below which the relative permeability of the nonwetting phase gets regul...
Definition: parkervangenuchten.hh:646
Scalar krgLowSte() const
Threshold saturation below which the relative permeability of the nonwetting phase gets regularized.
Definition: parkervangenuchten.hh:653