version 3.8
parkervangenuchten.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef PARKER_VANGENUCHTEN_3P_HH
13#define PARKER_VANGENUCHTEN_3P_HH
14
15#include <algorithm>
21
22namespace Dumux::FluidMatrix {
23
25{
32 template<class Scalar>
33 struct Params
34 {
35 Params(const Scalar swr = 0.0, const Scalar snr = 0.0, const Scalar sgr = 0.0)
36 : swr_(swr), snr_(snr), sgr_(sgr)
37 {}
38
42 Scalar swr() const
43 { return swr_; }
44
48 void setSwr(Scalar v)
49 { swr_ = v; }
50
54 Scalar snr() const
55 { return snr_; }
56
60 void setSnr(Scalar v)
61 { snr_ = v; }
66 Scalar sgr() const
67 { return sgr_; }
68
72 void setSgr(Scalar v)
73 { sgr_ = v; }
74
75 bool operator== (const Params& p) const
76 {
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);
80 }
81 private:
82 Scalar swr_;
83 Scalar snr_;
84 Scalar sgr_;
85 };
86
91 template<class Scalar>
92 static Params<Scalar> makeParams(const std::string& paramGroup)
93 {
94 Params<Scalar> params;
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));
98 return params;
99 }
100
110 template<class Scalar>
111 static Scalar swToSwe(const Scalar sw, const Params<Scalar>& params)
112 {
113 return (sw - params.swr())/(1.0 - params.swr()); // TODO other residual saturations?
114 }
115
125 template<class Scalar>
126 static Scalar sweToSw(const Scalar swe, const Params<Scalar>& params)
127 {
128 return swe*(1.0 - params.swr()) + params.swr(); // TODO other residual saturations?
129 }
130
139 template<class Scalar>
140 static Scalar dswe_dsw(const Params<Scalar>& params)
141 {
142 return 1.0/(1.0 - params.swr()); // TODO other residual saturations?
143 }
144
153 template<class Scalar>
154 static Scalar dsw_dswe(const Params<Scalar>& params)
155 {
156 return 1.0 - params.swr(); // TODO other residual saturations?
157 }
158
168 template<class Scalar>
169 static Scalar snToSne(const Scalar sn, const Params<Scalar>& params)
170 {
171 return sn; // sne equals sn // TODO other residual saturations?
172 }
173
183 template<class Scalar>
184 static Scalar stToSte(const Scalar st, const Params<Scalar>& params)
185 {
186 return (st-params.swr()) / (1-params.swr()); // TODO other residual saturations?
187 }
188
198 template<class Scalar>
199 static Scalar steToSt(const Scalar ste, const Params<Scalar>& params)
200 {
201 return ste*(1.0 - params.swr()) + params.swr(); // TODO other residual saturations?
202 }
203
212 template<class Scalar>
213 static Scalar dste_dst(const Params<Scalar>& params)
214 {
215 return 1.0/(1.0 - params.swr() /*- params.snr() - params.sgr()*/); // TODO other residual saturations?
216 }
217
226 template<class Scalar>
227 static Scalar dst_dste(const Params<Scalar>& params)
228 {
229 return 1.0 - params.swr(); // TODO other residual saturations?
230 }
231};
232
241{
242
243public:
258 template<class Scalar>
259 struct Params
260 {
261 Params(Scalar alpha, Scalar n, Scalar swr = 0.0, Scalar snr = 0.0,
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)
264 , betaNw_(betaNw), betaGn_(betaGn), betaGw_(betaGw), regardSnr_(regardSnr)
265 {}
266
267 Scalar alpha() const { return alpha_; }
268 void setAlpha(Scalar alpha) { alpha_ = alpha; }
269
270 Scalar m() const { return m_; }
271 void setM(Scalar m) { m_ = m; n_ = 1.0/(1.0 - m); }
272
273 Scalar n() const{ return n_; }
274 void setN(Scalar n){ n_ = n; m_ = 1.0 - 1.0/n; }
275
276 Scalar swr() const { return swr_; }
277 void setSwr(Scalar swr) { swr_ = swr; }
278
279 Scalar snr() const { return snr_; }
280 void setSnr(Scalar swr) { snr_ = swr; }
281
282 Scalar betaNw() const { return betaNw_; }
283 void setBetaNw(Scalar betaNw) { betaNw_ = betaNw; }
284
285 Scalar betaGn() const { return betaGn_; }
286 void setBetaGn(Scalar betaGn) { betaGn_ = betaGn; }
287
288 Scalar betaGw() const { return betaGw_; }
289 void setBetaGw(Scalar betaGw) { betaGw_ = betaGw; }
290
291 bool regardSnrForKrn() const { return regardSnr_; }
292 void setRegardSnrForKrn(bool v) {regardSnr_ = v; }
293
294 bool operator== (const Params& p) const
295 {
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_;
305 }
306
307 private:
308 Scalar alpha_, n_, m_, swr_, snr_, betaNw_, betaGn_, betaGw_;
309 bool regardSnr_;
310 };
311
316 template<class Scalar = double>
317 static Params<Scalar> makeParams(const std::string& paramGroup)
318 {
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);
327 return Params<Scalar>(alpha, n, swr, snr,
328 betaNw, betaGn, betaGw, regardSnr );
329 }
330
336 template<class Scalar>
337 static Scalar pcgw(Scalar swe, const Params<Scalar>& params)
338 {
339 assert(0 <= swe && swe <= 1);
340 return pc_(swe, params);
341 }
342
348 template<class Scalar>
349 static Scalar pcnw(Scalar swe, const Params<Scalar>& params)
350 {
351 assert(0 <= swe && swe <= 1);
352 return pc_(swe, params)/params.betaNw();
353 }
354
360 template<class Scalar>
361 static Scalar pcgn(const Scalar ste, const Params<Scalar>& params)
362 {
363 assert(0 <= ste && ste <= 1);
364 return pc_(ste, params)/params.betaGn();
365 }
366
372 template<class Scalar>
373 static Scalar pcAlpha(Scalar sne, const Params<Scalar>& params)
374 {
375 /* regularization */
376 if (sne <= 0.001)
377 sne = 0.0;
378 if (sne >= 1.0)
379 sne = 1.0;
380
381 if (sne > params.snr())
382 return 1.0;
383 else
384 {
385 if (params.snr() >= 0.001)
386 return sne/params.snr();
387 else
388 return 0.0;
389 };
390 }
391
398 template<class Scalar>
399 static Scalar dpcgw_dswe(const Scalar swe, const Params<Scalar>& params)
400 {
401 using std::pow;
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();
405 }
406
413 template<class Scalar>
414 static Scalar dpcnw_dswe(const Scalar swe, const Params<Scalar>& params)
415 {
416 using std::pow;
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();
420 }
421
428 template<class Scalar>
429 static Scalar dpcgn_dste(const Scalar ste, const Params<Scalar>& params)
430 {
431 using std::pow;
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();
435 }
436
449 template<class Scalar>
450 static Scalar krw(const Scalar swe, const Params<Scalar>& params)
451 {
452 using std::pow;
453 using std::sqrt;
454 const Scalar r = 1.0 - pow(1 - pow(swe, 1/params.m()), params.m());
455 return sqrt(swe)*r*r;
456 }
457
475 template<class Scalar>
476 static Scalar krn(const Scalar swe, const Scalar sn, const Scalar ste, const Params<Scalar>& params)
477 {
478 Scalar krn;
479 using std::pow;
480 krn = pow(1 - pow(swe, 1/params.m()), params.m());
481 krn -= pow(1 - pow(ste, 1/params.m()), params.m());
482 krn *= krn;
483
484 using std::clamp;
485 using std::sqrt;
486 if (params.regardSnrForKrn())
487 {
488 // regard Snr in the permeability of the n-phase, see Helmig1997
489 const Scalar resIncluded = clamp(sn - params.snr()/ (1-params.swr()), 0.0, 1.0);
490 krn *= sqrt(resIncluded);
491 }
492 else
493 krn *= sqrt(sn / (1 - params.swr())); // Hint: (ste - swe) = sn / (1-Swr)
494
495 return krn;
496 }
497
510 template<class Scalar>
511 static Scalar krg(const Scalar ste, const Params<Scalar>& params)
512 {
513 assert(0 <= ste && ste <= 1);
514 using std::cbrt;
515 using std::pow;
516 return cbrt(1 - ste) * pow(1 - pow(ste, 1/params.m()), 2*params.m());
517 }
518
528 template<class Scalar>
529 static Scalar dkrg_dste(const Scalar ste, const Params<Scalar>& params)
530 {
531 assert(0 < ste && ste <= 1);
532
533 using std::pow;
534 const Scalar x = pow(ste, 1.0/params.m());
535 return
536 -pow(1.0 - x, 2*params.m())
537 *pow(1.0 - ste, -2.0/3)
538 *(1.0/3 + 2*x/ste);
539 }
540
548 template<class Scalar>
549 static Scalar kr(const int phaseIdx, const Scalar swe, const Scalar sne, const Params<Scalar>& params)
550 {
551 switch (phaseIdx)
552 {
553 case 0:
554 return krw(params, swe, sne);
555 case 1:
556 return krn(params, swe, sne);
557 case 2:
558 return krg(params, swe, sne);
559 }
560 DUNE_THROW(Dune::InvalidStateException,
561 "Invalid phase index ");
562 }
563
564private:
565
572 template<class Scalar>
573 const static Scalar pc_(const Scalar se, const Params<Scalar>& params)
574 {
575 using std::pow;
576 return pow(pow(se, -1/params.m()) - 1, 1/params.n())/params.alpha();
577 }
578
579};
580
588template <class Scalar>
590{
591 using BaseLawParams = typename ParkerVanGenuchten3P::Params<Scalar>;
592
593public:
595 template<class S>
596 struct Params
597 {
605 { pcLowSwe_ = pcLowSwe; }
606
610 Scalar pcLowSwe() const
611 { return pcLowSwe_; }
612
617 { pcHighSwe_ = pcHighSwe; }
618
625 Scalar pcHighSwe() const
626 { return pcHighSwe_; }
627
633 { krnLowSwe_ = krnLowSwe; }
634
639 Scalar krnLowSwe() const
640 { return krnLowSwe_; }
641
647 { krgLowSte_ = krgLowSte; }
648
653 Scalar krgLowSte() const
654 { return krgLowSte_; }
655
661 { krwHighSwe_ = krwHighSwe; }
662
667 Scalar krwHighSwe() const
668 { return krwHighSwe_; }
669
670
676 void setConstRegularization(const bool input)
677 { constRegularization_ = input; }
678
684 { return constRegularization_; }
685
686private:
687 S pcLowSwe_ = 0.01;
688 S pcHighSwe_ = 0.99;
689 S krnLowSwe_ = 0.1;
690 S krwHighSwe_ = 0.9;
691 S krgLowSte_ = 1e-3;
692 bool constRegularization_ = false;
693 };
694
696 template<class MaterialLaw>
697 void init(const MaterialLaw* m, const std::string& paramGroup)
698 {
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);
705
706 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
707 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
708 }
709
710 template<class MaterialLaw, class BaseParams, class EffToAbsParams>
711 void init(const MaterialLaw* m, const BaseParams& bp, const EffToAbsParams& etap, const Params<Scalar>& p)
712 {
713 pcLowSwe_ = p.pcLowSwe();
714 pcHighSwe_ = p.pcHighSwe();
715 krwHighSwe_ = p.krwHighSwe();
716 krnLowSwe_ = p.krnLowSwe();
717 krgLowSte_ = p.krgLowSte();
718 constRegularization_ = p.constRegularization();
719
720 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
721 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
722 }
723
728 {
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_;
734 }
735
744 OptionalScalar<Scalar> pcgw(Scalar swe) const
745 {
746 // if specified, a constant value is used for regularization
747 using std::clamp;
748 if (constRegularization_)
749 swe = clamp(swe, 0.0, 1.0);
750
751 // make sure that the capillary pressure observes a derivative
752 // != 0 for 'illegal' saturations. This is favourable for the
753 // newton solver (if the derivative is calculated numerically)
754 // in order to get the saturation moving to the right
755 // direction if it temporarily is in an 'illegal' range.
756 if (swe < pcLowSwe_)
757 return pcgwLowSwePcgwValue_ + pcgwDerivativeLowSw_*(swe - pcLowSwe_);
758
759 else if (swe > 1.0)
760 return pcgwDerivativeHighSweEnd_*(swe - 1.0);
761
762 else if (swe > pcHighSwe_)
763 return pcgwSpline_.eval(swe);
764
765 else
766 return {}; // no regularization
767 }
768
777 OptionalScalar<Scalar> pcnw(Scalar swe) const
778 {
779 // if specified, a constant value is used for regularization
780 using std::clamp;
781 if (constRegularization_)
782 swe = clamp(swe, 0.0, 1.0);
783
784 // make sure that the capillary pressure observes a derivative
785 // != 0 for 'illegal' saturations. This is favourable for the
786 // newton solver (if the derivative is calculated numerically)
787 // in order to get the saturation moving to the right
788 // direction if it temporarily is in an 'illegal' range.
789 if (swe < pcLowSwe_)
790 return pcnwLowSwePcnwValue_ + pcnwDerivativeLowSw_*(swe - pcLowSwe_);
791
792 else if (swe > 1.0)
793 return pcnwDerivativeHighSweEnd_*(swe - 1.0);
794
795 else if (swe > pcHighSwe_)
796 return pcnwSpline_.eval(swe);
797
798 else
799 return {}; // no regularization
800 }
801
810 OptionalScalar<Scalar> pcgn(Scalar ste) const
811 {
812 // if specified, a constant value is used for regularization
813 using std::clamp;
814 if (constRegularization_)
815 ste = clamp(ste, 0.0, 1.0);
816
817
818 // make sure that the capillary pressure observes a derivative
819 // != 0 for 'illegal' saturations. This is favourable for the
820 // newton solver (if the derivative is calculated numerically)
821 // in order to get the saturation moving to the right
822 // direction if it temporarily is in an 'illegal' range.
823 const Scalar pcLowSte = pcLowSwe_;
824 const Scalar pcHighSte = pcHighSwe_;
825 if (ste < pcLowSte)
826 return pcgnLowStePcgnValue_ + pcgnDerivativeLowSt_*(ste - pcLowSte);
827
828 else if (ste > 1.0)
829 return pcgnDerivativeHighSteEnd_*(ste - 1.0);
830
831 else if (ste > pcHighSte)
832 return pcgnSpline_.eval(ste);
833
834 else
835 return {}; // no regularization
836 }
837
843 {
844 // no regularization
845 return {};
846 }
847
852 OptionalScalar<Scalar> krw(const Scalar swe) const
853 {
854 if (swe < 0.0)
855 return 0.0;
856 else if (swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
857 return 1.0;
858 else
859 return {}; // no regularization
860 }
861
862
869 OptionalScalar<Scalar> krn(Scalar swe, const Scalar sn, Scalar ste) const
870 {
871 using std::clamp;
872 swe = clamp(swe, 0.0, 1.0);
873 ste = clamp(ste, 0.0, 1.0);
874
875 if (ste - swe <= 0.0)
876 return 0.0;
877 else
878 return ParkerVanGenuchten3P::krn(swe, sn, ste, *baseLawParamsPtr_);
879 }
880
885 OptionalScalar<Scalar> krg(const Scalar ste) const
886 {
887 //return 0 if there is no gas
888 if (ste > 1.0 - std::numeric_limits<Scalar>::epsilon())
889 return 0.0;
890
891 // use linear regularization for very high gas saturations
892 // to avoid a kink in the curve and to maintain a slope for
893 // the Newton solver
894 if (ste <= krgLowSte_)
895 return krgLowStkrgValue_ + krgDerivativeLowSt_*(ste - krgLowSte_);
896 else
897 {
898 // For very low gas saturations:
899 // We use a scaling factor that decreases the gas phase permeability quite fast a very low gas phase
900 // saturations, thus making that phase virtually immobile.
901 // This prevents numerical issues related to the degeneration of the gas phase mass balance for the 3p3c model
902 // at very low gas phase saturations.
903
904 // get the absolute gas phase saturation
905 const Scalar st = ste*(1 - swr_) + swr_;
906 const Scalar sg = 1.0 - st;
907
908 // do not regularize
909 if (sg > 0.1)
910 return {};
911
912 // return original curve scaled by factor
913 using std::max;
914 const Scalar scalFact = max(0.0, (sg - sgr_)/(0.1 - sgr_));
915
916 return ParkerVanGenuchten3P::krg(ste, *baseLawParamsPtr_) * scalFact;
917 }
918 }
919
926 OptionalScalar<Scalar> kr(const int phaseIdx, const Scalar swe, const Scalar sne) const
927 {
928 switch (phaseIdx)
929 {
930 case 0:
931 return krw(swe, sne);
932 case 1:
933 return krn(swe, sne);
934 case 2:
935 return krg(swe, sne);
936 }
937 DUNE_THROW(Dune::InvalidStateException,
938 "Invalid phase index ");
939 }
940
941private:
942 template<class MaterialLaw>
943 void initPcParameters_(const MaterialLaw* m, const Scalar lowSwe, const Scalar highSwe)
944 {
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());
949
950 baseLawParamsPtr_ = &m->basicParams();
951
952 // pcgw
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);
958 pcgwSpline_ = Spline<Scalar>(highSwe, 1.0, // x0, x1
959 pcgwHighSwePcgwValue_, 0, // y0, y1
960 pcgwDerivativeHighSweThreshold_, pcgwDerivativeHighSweEnd_); // m0, m1
961
962 // pcnw
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);
968 pcnwSpline_ = Spline<Scalar>(highSwe, 1.0, // x0, x1
969 pcnwHighSwePcnwValue_, 0, // y0, y1
970 pcnwDerivativeHighSweThreshold_, pcnwDerivativeHighSweEnd_); // m0, m1
971
972 // pcgn
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);
978 pcgnSpline_ = Spline<Scalar>(highSwe, 1.0, // x0, x1
979 pcgnHighSwePcgnValue_, 0, // y0, y1
980 pcgnDerivativeHighSteThreshold_, pcgnDerivativeHighSteEnd_); // m0, m1
981
982 }
983
984 template<class MaterialLaw>
985 void initKrParameters_(const MaterialLaw* m, const Scalar lowSwe, const Scalar highSwe)
986 {
987 krgLowStkrgValue_ = ParkerVanGenuchten3P::krg(krgLowSte_, *baseLawParamsPtr_);
988 krgDerivativeLowSt_ = ParkerVanGenuchten3P::dkrg_dste(krgLowSte_, *baseLawParamsPtr_);
989
990 swr_ = m->effToAbsParams().swr();
991 sgr_ = m->effToAbsParams().sgr();
992 }
993
994 Scalar krgLowStkrgValue_;
995 Scalar krgDerivativeLowSt_;
996
997 Scalar pcLowSwe_, pcHighSwe_;
998 Scalar krwHighSwe_, krnLowSwe_, krgLowSte_;
999
1000 // pcgw
1001 Scalar pcgwLowSwePcgwValue_;
1002 Scalar pcgwHighSwePcgwValue_;
1003 Scalar pcgwDerivativeLowSw_;
1004 Scalar pcgwDerivativeHighSweThreshold_;
1005 Scalar pcgwDerivativeHighSweEnd_;
1006
1007 // pcgn
1008 Scalar pcgnLowStePcgnValue_;
1009 Scalar pcgnHighSwePcgnValue_;
1010 Scalar pcgnDerivativeLowSt_;
1011 Scalar pcgnDerivativeHighSteThreshold_;
1012 Scalar pcgnDerivativeHighSteEnd_;
1013
1014 // pcnw
1015 Scalar pcnwLowSwePcnwValue_;
1016 Scalar pcnwHighSwePcnwValue_;
1017 Scalar pcnwDerivativeLowSw_;
1018 Scalar pcnwDerivativeHighSweThreshold_;
1019 Scalar pcnwDerivativeHighSweEnd_;
1020
1021 Spline<Scalar> pcgwSpline_;
1022 Spline<Scalar> pcnwSpline_;
1023 Spline<Scalar> pcgnSpline_;
1024 Spline<Scalar> krwSpline_;
1025 Spline<Scalar> krnSpline_;
1026
1027 Scalar swr_, sgr_;
1028
1029 bool constRegularization_;
1030
1031 const BaseLawParams* baseLawParamsPtr_;
1032};
1033
1038template<class ScalarType,
1039 class BaseLaw,
1040 class Regularization = NoRegularization,
1041 class EffToAbsPolicy = ParkerVanGenuchten3PEffToAbsPolicy>
1042class ParkerVanGenuchtenMaterialLaw : public Adapter<ParkerVanGenuchtenMaterialLaw<ScalarType, BaseLaw, Regularization, EffToAbsPolicy>, ThreePhasePcKrSw>
1043{
1044public:
1045
1046 using Scalar = ScalarType;
1047
1048 using BasicParams = typename BaseLaw::template Params<Scalar>;
1049 using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>;
1050 using RegularizationParams = typename Regularization::template Params<Scalar>;
1051
1052 using EffToAbs = EffToAbsPolicy;
1053
1057 static constexpr bool isRegularized()
1058 { return !std::is_same<Regularization, NoRegularization>::value; }
1059
1065
1070 explicit ParkerVanGenuchtenMaterialLaw(const std::string& paramGroup)
1071 : basicParams_(makeBasicParams(paramGroup))
1072 , effToAbsParams_(makeEffToAbsParams(paramGroup))
1073 {
1074 if constexpr (isRegularized())
1075 regularization_.init(this, paramGroup);
1076 }
1077
1083 const EffToAbsParams& effToAbsParams = {},
1084 const RegularizationParams& regParams = {})
1085 : basicParams_(baseParams)
1086 , effToAbsParams_(effToAbsParams)
1087 {
1088 if constexpr (isRegularized())
1089 regularization_.init(this, baseParams, effToAbsParams, regParams);
1090 }
1091
1095 template<bool enableRegularization = isRegularized()>
1096 Scalar pcgw(const Scalar sw, const Scalar /*dummySn*/) const
1097 {
1098 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1099 if constexpr (enableRegularization)
1100 {
1101 const auto regularized = regularization_.pcgw(swe);
1102 if (regularized)
1103 return regularized.value();
1104 }
1105
1106 return BaseLaw::pcgw(swe, basicParams_);
1107 }
1108
1112 template<bool enableRegularization = isRegularized()>
1113 Scalar pcnw(const Scalar sw, const Scalar /*dummySn*/) const
1114 {
1115 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1116 if constexpr (enableRegularization)
1117 {
1118 const auto regularized = regularization_.pcnw(swe);
1119 if (regularized)
1120 return regularized.value();
1121 }
1122
1123 return BaseLaw::pcnw(swe, basicParams_);
1124 }
1125
1131 template<bool enableRegularization = isRegularized()>
1132 Scalar pcgn(const Scalar sw, const Scalar sn) const
1133 {
1134 const auto swe = EffToAbs::swToSwe(sw + sn, effToAbsParams_);
1135 if constexpr (enableRegularization)
1136 {
1137 const auto regularized = regularization_.pcgn(swe);
1138 if (regularized)
1139 return regularized.value();
1140 }
1141
1142 return BaseLaw::pcgn(swe, basicParams_);
1143 }
1144
1148 template<bool enableRegularization = isRegularized()>
1149 Scalar pcAlpha(const Scalar /*dummySw*/, const Scalar sn) const
1150 {
1151 const auto sne = EffToAbs::snToSne(sn, effToAbsParams_);
1152 if constexpr (enableRegularization)
1153 {
1154 const auto regularized = regularization_.pcAlpha(sne);
1155 if (regularized)
1156 return regularized.value();
1157 }
1158 return BaseLaw::pcAlpha(sne, basicParams_);
1159 }
1160
1164 template<bool enableRegularization = isRegularized()>
1165 Scalar dpcgw_dsw(const Scalar sw, const Scalar /*dummyS*/) const
1166 {
1167 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1168 if constexpr (enableRegularization)
1169 {
1170 const auto regularized = regularization_.dpcgw_dswe(swe);
1171 if (regularized)
1172 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
1173 }
1174
1175 return BaseLaw::dpcgw_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
1176 }
1177
1181 template<bool enableRegularization = isRegularized()>
1182 Scalar dpcnw_dsw(const Scalar sw, const Scalar /*dummySw*/) const
1183 {
1184 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1185 if constexpr (enableRegularization)
1186 {
1187 const auto regularized = regularization_.dpcnw_dswe(swe);
1188 if (regularized)
1189 return regularized.value()*EffToAbs::dswe_dsw(effToAbsParams_);
1190 }
1191
1192 return BaseLaw::dpcnw_dswe(swe, basicParams_)*EffToAbs::dswe_dsw(effToAbsParams_);
1193 }
1194
1198 template<bool enableRegularization = isRegularized()>
1199 Scalar dpcgn_dst(const Scalar st, const Scalar /*dummySw*/) const
1200 {
1201 const auto ste = EffToAbs::stToSte(st, effToAbsParams_);
1202 if constexpr (enableRegularization)
1203 {
1204 const auto regularized = regularization_.dpcgn_dste(ste);
1205 if (regularized)
1206 return regularized.value()*EffToAbs::dswte_dst(effToAbsParams_);
1207 }
1208
1209 return BaseLaw::dpcgn_dste(ste, basicParams_)*EffToAbs::dste_dst(effToAbsParams_);
1210 }
1211
1217 template<bool enableRegularization = isRegularized()>
1218 Scalar krw(const Scalar sw, const Scalar sn) const
1219 {
1220 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1221 if constexpr (enableRegularization)
1222 {
1223 const auto regularized = regularization_.krw(swe);
1224 if (regularized)
1225 return regularized.value();
1226 }
1227
1228 return BaseLaw::krw(swe, basicParams_);
1229 }
1230
1236 template<bool enableRegularization = isRegularized()>
1237 Scalar krn(const Scalar sw, const Scalar sn) const
1238 {
1239 const auto swe = EffToAbs::swToSwe(sw, effToAbsParams_);
1240 const auto ste = EffToAbs::stToSte(sw + sn, effToAbsParams_);
1241 if constexpr (enableRegularization)
1242 {
1243 const auto regularized = regularization_.krn(swe, sn, ste);
1244 if (regularized)
1245 return regularized.value();
1246 }
1247
1248 return BaseLaw::krn(swe, sn, ste, basicParams_);
1249 }
1250
1256 template<bool enableRegularization = isRegularized()>
1257 Scalar krg(const Scalar sw, const Scalar sn) const
1258 {
1259 const auto ste = EffToAbs::stToSte(sw + sn, effToAbsParams_);
1260 if constexpr (enableRegularization)
1261 {
1262 const auto regularized = regularization_.krg(ste);
1263 if (regularized)
1264 return regularized.value();
1265 }
1266
1267 return BaseLaw::krg(ste, basicParams_);
1268 }
1269
1276 template<bool enableRegularization = isRegularized()>
1277 Scalar kr(const int phaseIdx, const Scalar sw, const Scalar sn) const
1278 {
1279 switch (phaseIdx)
1280 {
1281 case 0:
1282 return krw(sw, sn);
1283 case 1:
1284 return krn(sw, sn);
1285 case 2:
1286 return krg(sw, sn);
1287 }
1288 DUNE_THROW(Dune::InvalidStateException,
1289 "Invalid phase index ");
1290 }
1291
1296 template<bool enableRegularization = isRegularized()>
1297 Scalar dkrg_dst(const Scalar st) const
1298 {
1299 const auto ste = EffToAbs::stToSte(st, effToAbsParams_);
1300 if constexpr (enableRegularization)
1301 {
1302 const auto regularized = regularization_.dkrg_dste(ste);
1303 if (regularized)
1304 return regularized.value()*EffToAbs::dste_dst(effToAbsParams_);
1305 }
1306
1307 return BaseLaw::dkrg_dste(ste, basicParams_)*EffToAbs::dste_dst(effToAbsParams_);
1308 }
1309
1314 {
1315 return basicParams_ == o.basicParams_
1316 && effToAbsParams_ == o.effToAbsParams_
1317 && regularization_ == o.regularization_;
1318 }
1319
1324 static BasicParams makeBasicParams(const std::string& paramGroup)
1325 {
1326 return BaseLaw::template makeParams<Scalar>(paramGroup);
1327 }
1328
1333 { return basicParams_; }
1334
1339 static EffToAbsParams makeEffToAbsParams(const std::string& paramGroup)
1340 {
1341 return EffToAbs::template makeParams<Scalar>(paramGroup);
1342 }
1343
1348 { return effToAbsParams_; }
1349
1350private:
1351 BasicParams basicParams_;
1352 EffToAbsParams effToAbsParams_;
1353 Regularization regularization_;
1354};
1355
1360template<class Scalar>
1362
1367template<class Scalar>
1369
1370} // end namespace Dumux::FluidMatrix
1371
1372#endif
Implementation of Parker/vanGenuchten's capillary pressure <-> saturation relation for three phases....
Definition: parkervangenuchten.hh:241
static Params< Scalar > makeParams(const std::string &paramGroup)
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 > &params)
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 > &params)
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 > &params)
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 > &params)
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 > &params)
The relative permeability for a phase.
Definition: parkervangenuchten.hh:549
static Scalar dkrg_dste(const Scalar ste, const Params< Scalar > &params)
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 > &params)
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 > &params)
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 > &params)
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 > &params)
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 > &params)
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 > &params)
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 &paramGroup)
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 &paramGroup)
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 &paramGroup)
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 &paramGroup)
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 &regParams={})
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 > &params)
Convert an absolute nonwetting saturation to an effective one.
Definition: parkervangenuchten.hh:169
static Scalar dswe_dsw(const Params< Scalar > &params)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: parkervangenuchten.hh:140
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: parkervangenuchten.hh:92
static Scalar steToSt(const Scalar ste, const Params< Scalar > &params)
Convert an effective wetting saturation to an absolute one.
Definition: parkervangenuchten.hh:199
static Scalar dste_dst(const Params< Scalar > &params)
Derivative of the effective saturation w.r.t. the absolute saturation.
Definition: parkervangenuchten.hh:213
static Scalar stToSte(const Scalar st, const Params< Scalar > &params)
Convert an absolute total liquid saturation to an effective one.
Definition: parkervangenuchten.hh:184
static Scalar swToSwe(const Scalar sw, const Params< Scalar > &params)
Convert an absolute wetting saturation to an effective one.
Definition: parkervangenuchten.hh:111
static Scalar dsw_dswe(const Params< Scalar > &params)
Derivative of the absolute saturation w.r.t. the effective saturation.
Definition: parkervangenuchten.hh:154
static Scalar sweToSw(const Scalar swe, const Params< Scalar > &params)
Convert an effective wetting saturation to an absolute one.
Definition: parkervangenuchten.hh:126
static Scalar dst_dste(const Params< Scalar > &params)
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