3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
vangenuchten.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH
27
28// remove from here after release 3.3 /////////////
29#include "vangenuchtenparams.hh"
30
31#include <algorithm>
32#include <cmath>
33#include <cassert>
34
35namespace Dumux {
36
51template <class ScalarT, class ParamsT = VanGenuchtenParams<ScalarT> >
52class [[deprecated("Use new material laws and FluidMatrix::VanGenuchten instead!")]] VanGenuchten
53{
54public:
55 using Params = ParamsT;
56 using Scalar = typename Params::Scalar;
57
73 static Scalar pc(const Params &params, Scalar swe)
74 {
75 using std::pow;
76 using std::clamp;
77
78 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
79
80 const Scalar pc = pow(pow(swe, -1.0/params.vgm()) - 1, 1.0/params.vgn())/params.vgAlpha();
81 return pc;
82 }
83
100 static Scalar sw(const Params &params, Scalar pc)
101 {
102 using std::pow;
103 using std::max;
104
105 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
106
107 const Scalar sw = pow(pow(params.vgAlpha()*pc, params.vgn()) + 1, -params.vgm());
108 return sw;
109 }
110
119 static Scalar endPointPc(const Params &params)
120 { return 0.0; }
121
142 static Scalar dpc_dswe(const Params &params, Scalar swe)
143 {
144 using std::pow;
145 using std::clamp;
146
147 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
148
149 const Scalar powSwe = pow(swe, -1/params.vgm());
150 return - 1.0/params.vgAlpha() * pow(powSwe - 1, 1.0/params.vgn() - 1)/params.vgn()
151 * powSwe/swe/params.vgm();
152 }
153
167 static Scalar dswe_dpc(const Params &params, Scalar pc)
168 {
169 using std::pow;
170 using std::max;
171
172 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
173
174 const Scalar powAlphaPc = pow(params.vgAlpha()*pc, params.vgn());
175 return -pow(powAlphaPc + 1, -params.vgm()-1)*params.vgm()*powAlphaPc/pc*params.vgn();
176 }
177
192 static Scalar krw(const Params &params, Scalar swe)
193 {
194 using std::pow;
195 using std::clamp;
196
197 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
198
199 const Scalar r = 1.0 - pow(1.0 - pow(swe, 1.0/params.vgm()), params.vgm());
200 return pow(swe, params.vgl())*r*r;
201 }
202
217 static Scalar dkrw_dswe(const Params &params, Scalar swe)
218 {
219 using std::pow;
220 using std::clamp;
221
222 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
223
224 const Scalar x = 1.0 - pow(swe, 1.0/params.vgm());
225 const Scalar xToM = pow(x, params.vgm());
226 return (1.0 - xToM)*pow(swe, params.vgl()-1) * ( (1.0 - xToM)*params.vgl() + 2*xToM*(1.0-x)/x );
227 }
228
244 static Scalar krn(const Params &params, Scalar swe)
245 {
246 using std::pow;
247 using std::clamp;
248
249 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
250
251 return pow(1 - swe, params.vgl()) * pow(1 - pow(swe, 1.0/params.vgm()), 2*params.vgm());
252 }
253
269 static Scalar dkrn_dswe(const Params &params, Scalar swe)
270 {
271 using std::pow;
272 using std::clamp;
273
274 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
275
276 const auto sne = 1.0 - swe;
277 const auto x = 1.0 - pow(swe, 1.0/params.vgm());
278 return -pow(sne, params.vgl()-1.0) * pow(x, 2*params.vgm() - 1.0) * ( params.vgl()*x + 2.0*sne/swe*(1.0 - x) );
279 }
280
281};
282
283} // end namespace Dumux
284// remove until here after release 3.3 /////////////
285
286#include <cmath>
287#include <algorithm>
288
290#include <dumux/common/spline.hh>
293
294namespace Dumux::FluidMatrix {
295
305{
306
307public:
321 template<class Scalar>
322 struct Params
323 {
324 Params(Scalar alpha, Scalar n, Scalar l = 0.5)
325 : alpha_(alpha), n_(n), m_(1.0 - 1.0/n), l_(l)
326 {}
327
328 Scalar alpha() const { return alpha_; }
329 void setAlpha(Scalar alpha) { alpha_ = alpha; }
330
331 Scalar m() const { return m_; }
332 void setM(Scalar m) { m_ = m; n_ = 1.0/(1.0 - m); }
333
334 Scalar n() const{ return n_; }
335 void setN(Scalar n){ n_ = n; m_ = 1.0 - 1.0/n; }
336
337 Scalar l() const { return l_; }
338 void setL(Scalar l) { l_ = l; }
339
340 bool operator== (const Params& p) const
341 {
342 return Dune::FloatCmp::eq(alpha_, p.alpha_, 1e-6)
343 && Dune::FloatCmp::eq(n_, p.n_, 1e-6)
344 && Dune::FloatCmp::eq(m_, p.m_, 1e-6)
345 && Dune::FloatCmp::eq(l_, p.l_, 1e-6);
346 }
347
348 private:
349 Scalar alpha_, n_, m_, l_;
350 };
351
356 template<class Scalar = double>
357 static Params<Scalar> makeParams(const std::string& paramGroup)
358 {
359 const auto n = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenN");
360 const auto alpha = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenAlpha");
361 // l is usually chosen to be 0.5 (according to Mualem (1976), WRR)
362 const auto l = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenL", 0.5);
363 return Params<Scalar>(alpha, n, l);
364 }
365
379 template<class Scalar>
380 static Scalar pc(Scalar swe, const Params<Scalar>& params)
381 {
382 using std::pow;
383 using std::clamp;
384
385 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
386
387 const Scalar pc = pow(pow(swe, -1.0/params.m()) - 1, 1.0/params.n())/params.alpha();
388 return pc;
389 }
390
405 template<class Scalar>
406 static Scalar swe(Scalar pc, const Params<Scalar>& params)
407 {
408 using std::pow;
409 using std::max;
410
411 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
412
413 const Scalar sw = pow(pow(params.alpha()*pc, params.n()) + 1, -params.m());
414 return sw;
415 }
416
421 template<class Scalar>
422 static Scalar endPointPc(const Params<Scalar>& params)
423 { return 0.0; }
424
441 template<class Scalar>
442 static Scalar dpc_dswe(Scalar swe, const Params<Scalar>& params)
443 {
444 using std::pow;
445 using std::clamp;
446
447 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
448
449 const Scalar powSwe = pow(swe, -1/params.m());
450 return - 1.0/params.alpha() * pow(powSwe - 1, 1.0/params.n() - 1)/params.n()
451 * powSwe/swe/params.m();
452 }
453
463 template<class Scalar>
464 static Scalar dswe_dpc(Scalar pc, const Params<Scalar>& params)
465 {
466 using std::pow;
467 using std::max;
468
469 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
470
471 const Scalar powAlphaPc = pow(params.alpha()*pc, params.n());
472 return -pow(powAlphaPc + 1, -params.m()-1)*params.m()*powAlphaPc/pc*params.n();
473 }
474
485 template<class Scalar>
486 static Scalar krw(Scalar swe, const Params<Scalar>& params)
487 {
488 using std::pow;
489 using std::clamp;
490
491 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
492
493 const Scalar r = 1.0 - pow(1.0 - pow(swe, 1.0/params.m()), params.m());
494 return pow(swe, params.l())*r*r;
495 }
496
507 template<class Scalar>
508 static Scalar dkrw_dswe(Scalar swe, const Params<Scalar>& params)
509 {
510 using std::pow;
511 using std::clamp;
512
513 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
514
515 const Scalar x = 1.0 - pow(swe, 1.0/params.m());
516 const Scalar xToM = pow(x, params.m());
517 return (1.0 - xToM)*pow(swe, params.l()-1) * ( (1.0 - xToM)*params.l() + 2*xToM*(1.0-x)/x );
518 }
519
530 template<class Scalar>
531 static Scalar krn(Scalar swe, const Params<Scalar>& params)
532 {
533 using std::pow;
534 using std::clamp;
535
536 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
537
538 return pow(1 - swe, params.l()) * pow(1 - pow(swe, 1.0/params.m()), 2*params.m());
539 }
540
552 template<class Scalar>
553 static Scalar dkrn_dswe(Scalar swe, const Params<Scalar>& params)
554 {
555 using std::pow;
556 using std::clamp;
557
558 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
559
560 const auto sne = 1.0 - swe;
561 const auto x = 1.0 - pow(swe, 1.0/params.m());
562 return -pow(sne, params.l()-1.0) * pow(x, 2*params.m() - 1.0) * ( params.l()*x + 2.0*sne/swe*(1.0 - x) );
563 }
564};
565
573template <class Scalar>
575{
576public:
578 template<class S>
579 struct Params
580 {
588 { pcLowSwe_ = pcLowSwe; }
589
593 Scalar pcLowSwe() const
594 { return pcLowSwe_; }
595
600 { pcHighSwe_ = pcHighSwe; }
601
608 Scalar pcHighSwe() const
609 { return pcHighSwe_; }
610
616 { krnLowSwe_ = krnLowSwe; }
617
622 Scalar krnLowSwe() const
623 { return krnLowSwe_; }
624
630 { krwHighSwe_ = krwHighSwe; }
631
636 Scalar krwHighSwe() const
637 { return krwHighSwe_; }
638
639private:
640 S pcLowSwe_ = 0.01;
641 S pcHighSwe_ = 0.99;
642 S krnLowSwe_ = 0.1;
643 S krwHighSwe_ = 0.9;
644 };
645
647 template<class MaterialLaw>
648 void init(const MaterialLaw* m, const std::string& paramGroup)
649 {
650 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenPcLowSweThreshold", 0.01);
651 pcHighSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenPcHighSweThreshold", 0.99);
652 krwHighSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenKrwHighSweThreshold", 0.9);
653 krnLowSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenKrnLowSweThreshold", 0.1);
654
655 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
656 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
657 }
658
659 template<class MaterialLaw, class BaseParams, class EffToAbsParams>
660 void init(const MaterialLaw* m, const BaseParams& bp, const EffToAbsParams& etap, const Params<Scalar>& p)
661 {
662 pcLowSwe_ = p.pcLowSwe();
663 pcHighSwe_ = p.pcHighSwe();
664 krwHighSwe_ = p.krwHighSwe();
665 krnLowSwe_ = p.krnLowSwe();
666
667 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
668 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
669 }
670
675 {
676 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
677 && Dune::FloatCmp::eq(pcHighSwe_, o.pcHighSwe_, 1e-6)
678 && Dune::FloatCmp::eq(krwHighSwe_, o.krwHighSwe_, 1e-6)
679 && Dune::FloatCmp::eq(krnLowSwe_, o.krnLowSwe_, 1e-6);
680 }
681
689 OptionalScalar<Scalar> pc(const Scalar swe) const
690 {
691 // make sure that the capillary pressure observes a derivative
692 // != 0 for 'illegal' saturations. This is favourable for the
693 // newton solver (if the derivative is calculated numerically)
694 // in order to get the saturation moving to the right
695 // direction if it temporarily is in an 'illegal' range.
696 if (swe < pcLowSwe_)
697 return pcLowSwePcValue_ + pcDerivativeLowSw_*(swe - pcLowSwe_);
698
699 else if (swe > 1.0)
700 return pcDerivativeHighSweEnd_*(swe - 1.0);
701
702 else if (swe > pcHighSwe_)
703 return pcSpline_.eval(swe);
704
705 else
706 return {}; // no regularization
707 }
708
713 {
714 if (swe < pcLowSwe_)
715 return pcDerivativeLowSw_;
716
717 else if (swe > 1.0)
718 return pcDerivativeHighSweEnd_;
719
720 else if (swe > pcHighSwe_)
721 return pcSpline_.evalDerivative(swe);
722
723 else
724 return {}; // no regularization
725 }
726
730 OptionalScalar<Scalar> swe(const Scalar pc) const
731 {
732 if (pc <= 0.0)
733 {
734 if (pcHighSwe_ > 1.0 - std::numeric_limits<Scalar>::epsilon())
735 return 1.0;
736 else
737 return pc/pcDerivativeHighSweEnd_ + 1.0;
738 }
739
740 // invert spline
741 else if (pc < pcHighSwePcValue_)
742 return pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0, pc);
743
744 else if (pc >= pcLowSwePcValue_)
745 return (pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
746
747 else
748 return {}; // no regularization
749 }
750
755 {
756 if (pc <= 0.0)
757 {
758 if (pcHighSwe_ > 1.0 - std::numeric_limits<Scalar>::epsilon())
759 return 0.0;
760 else
761 return 1.0/pcDerivativeHighSweEnd_;
762 }
763
764 // derivative of the inverse of the function is one over derivative of the function
765 else if (pc < pcHighSwePcValue_)
766 return 1.0/pcSpline_.evalDerivative(pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0, pc));
767
768 else if (pc >= pcLowSwePcValue_)
769 return 1.0/pcDerivativeLowSw_;
770
771 else
772 return {}; // no regularization
773 }
774
778 OptionalScalar<Scalar> krw(const Scalar swe) const
779 {
780 if (swe < 0.0)
781 return 0.0;
782 else if (swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
783 return 1.0;
784 else if (swe > krwHighSwe_)
785 return krwSpline_.eval(swe);
786 else
787 return {}; // no regularization
788 }
789
794 {
795 if (swe < 0.0)
796 return 0.0;
797 else if (swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
798 return 0.0;
799 else if (swe > krwHighSwe_)
800 return krwSpline_.evalDerivative(swe);
801 else
802 return {}; // no regularization
803 }
804
808 OptionalScalar<Scalar> krn(const Scalar swe) const
809 {
810 if (swe < 0.0)
811 return 1.0;
812 else if (swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
813 return 0.0;
814 else if (swe < krnLowSwe_)
815 return krnSpline_.eval(swe);
816 else
817 return {}; // no regularization
818 }
819
824 {
825 if (swe < 0.0)
826 return 0.0;
827 else if (swe > 1.0 - std::numeric_limits<Scalar>::epsilon())
828 return 0.0;
829 else if (swe < krnLowSwe_)
830 return krnSpline_.evalDerivative(swe);
831 else
832 return {}; // no regularization
833 }
834
835private:
836 template<class MaterialLaw>
837 void initPcParameters_(const MaterialLaw* m, const Scalar lowSwe, const Scalar highSwe)
838 {
839 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
840 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
841 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
842
843 pcDerivativeLowSw_ = m->template dpc_dsw<false>(lowSw)*dsw_dswe;
844
845 pcDerivativeHighSweThreshold_ = m->template dpc_dsw<false>(highSw)*dsw_dswe;
846 pcDerivativeHighSweEnd_ = 2.0*(0.0 - m->template pc<false>(highSw))/(1.0 - highSwe);
847
848 pcLowSwePcValue_ = m->template pc<false>(lowSw);
849 pcHighSwePcValue_ = m->template pc<false>(highSw);
850
851 pcSpline_ = Spline<Scalar>(highSwe, 1.0, // x0, x1
852 pcHighSwePcValue_, 0, // y0, y1
853 pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_); // m0, m1
854
855
856 }
857
858 template<class MaterialLaw>
859 void initKrParameters_(const MaterialLaw* m, const Scalar lowSwe, const Scalar highSwe)
860 {
861 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
862 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
863 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
864
865 const auto krwHighSw = m->template krw<false>(highSw);
866 const auto dkrwHighSw = m->template dkrw_dsw<false>(highSw)*dsw_dswe;
867
868 const auto krnLowSw = m->template krn<false>(lowSw);
869 const auto dkrnLowSw = m->template dkrn_dsw<false>(lowSw)*dsw_dswe;
870
871 krwSpline_ = Spline<Scalar>(highSwe, 1.0, // x0, x1
872 krwHighSw, 1.0, // y0, y1
873 dkrwHighSw, 0.0); // m0, m1
874
875 krnSpline_ = Spline<Scalar>(0.0, lowSwe, // x0, x1
876 1.0, krnLowSw, // y0, y1
877 0.0, dkrnLowSw); // m0, m1
878
879 }
880
881 Scalar pcLowSwe_, pcHighSwe_;
882 Scalar pcLowSwePcValue_, pcHighSwePcValue_;
883 Scalar krwHighSwe_, krnLowSwe_;
884 Scalar pcDerivativeLowSw_;
885 Scalar pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_;
886
887 Spline<Scalar> pcSpline_;
888 Spline<Scalar> krwSpline_;
889 Spline<Scalar> krnSpline_;
890};
891
896template<typename Scalar = double>
898
903template<typename Scalar = double>
905
906} // end namespace Dumux::FluidMatrix
907
908#endif
Provides 3rd order polynomial splines.
A wrapper that can either contain a valid Scalar or NaN.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
Specification of the material parameters for the van Genuchten-Mualem constitutive relations.
Definition: adapt.hh:29
Definition: brookscorey.hh:286
A 3rd order polynomial spline.
Definition: spline.hh:55
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
Definition: efftoabsdefaultpolicy.hh:46
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:57
Implementation of the van Genuchten capillary pressure <-> saturation relation. This class bundles th...
Definition: vangenuchten.hh:53
static Scalar dkrw_dswe(const Params &params, Scalar swe)
The derivative of the relative permeability for the wetting phase in regard to the wetting saturation...
Definition: vangenuchten.hh:217
static Scalar sw(const Params &params, Scalar pc)
The saturation-capillary pressure curve according to van Genuchten.
Definition: vangenuchten.hh:100
typename Params::Scalar Scalar
Definition: vangenuchten.hh:56
static Scalar dpc_dswe(const Params &params, Scalar swe)
The partial derivative of the capillary pressure w.r.t. the effective saturation according to van Gen...
Definition: vangenuchten.hh:142
static Scalar krw(const Params &params, Scalar swe)
The relative permeability for the wetting phase of the medium implied by van Genuchten / Mualem param...
Definition: vangenuchten.hh:192
static Scalar dswe_dpc(const Params &params, Scalar pc)
The partial derivative of the effective saturation to the capillary pressure according to van Genucht...
Definition: vangenuchten.hh:167
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: vangenuchten.hh:119
static Scalar krn(const Params &params, Scalar swe)
The relative permeability for the nonwetting phase of the medium implied by van Genuchten's parameter...
Definition: vangenuchten.hh:244
ParamsT Params
Definition: vangenuchten.hh:55
static Scalar pc(const Params &params, Scalar swe)
The capillary pressure-saturation curve according to van Genuchten.
Definition: vangenuchten.hh:73
static Scalar dkrn_dswe(const Params &params, Scalar swe)
The derivative of the relative permeability for the nonwetting phase in regard to the wetting saturat...
Definition: vangenuchten.hh:269
Implementation of the van Genuchten capillary pressure <-> saturation relation, and relative permeabi...
Definition: vangenuchten.hh:305
static Scalar dswe_dpc(Scalar pc, const Params< Scalar > &params)
The partial derivative of the effective saturation to the capillary pressure according to van Genucht...
Definition: vangenuchten.hh:464
static Scalar endPointPc(const Params< Scalar > &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: vangenuchten.hh:422
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability for the non-wetting phase in regard to the wetting satura...
Definition: vangenuchten.hh:553
static Scalar swe(Scalar pc, const Params< Scalar > &params)
The saturation-capillary pressure curve according to van Genuchten.
Definition: vangenuchten.hh:406
static Scalar krw(Scalar swe, const Params< Scalar > &params)
The relative permeability for the wetting phase of the medium implied by van Genuchten / Mualem param...
Definition: vangenuchten.hh:486
static Scalar dkrw_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability for the wetting phase in regard to the wetting saturation...
Definition: vangenuchten.hh:508
static Scalar dpc_dswe(Scalar swe, const Params< Scalar > &params)
The partial derivative of the capillary pressure w.r.t. the effective saturation according to van Gen...
Definition: vangenuchten.hh:442
static Scalar pc(Scalar swe, const Params< Scalar > &params)
The capillary pressure-saturation curve according to van Genuchten.
Definition: vangenuchten.hh:380
static Scalar krn(Scalar swe, const Params< Scalar > &params)
The relative permeability for the non-wetting phase of the medium implied by van Genuchten's paramete...
Definition: vangenuchten.hh:531
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: vangenuchten.hh:357
The parameter type.
Definition: vangenuchten.hh:323
Params(Scalar alpha, Scalar n, Scalar l=0.5)
Definition: vangenuchten.hh:324
Scalar alpha() const
Definition: vangenuchten.hh:328
bool operator==(const Params &p) const
Definition: vangenuchten.hh:340
Scalar l() const
Definition: vangenuchten.hh:337
void setAlpha(Scalar alpha)
Definition: vangenuchten.hh:329
void setL(Scalar l)
Definition: vangenuchten.hh:338
void setM(Scalar m)
Definition: vangenuchten.hh:332
Scalar n() const
Definition: vangenuchten.hh:334
Scalar m() const
Definition: vangenuchten.hh:331
void setN(Scalar n)
Definition: vangenuchten.hh:335
A regularization for the VanGenuchten material law.
Definition: vangenuchten.hh:575
bool operator==(const VanGenuchtenRegularization &o) const
Equality comparison with another instance.
Definition: vangenuchten.hh:674
void init(const MaterialLaw *m, const BaseParams &bp, const EffToAbsParams &etap, const Params< Scalar > &p)
Definition: vangenuchten.hh:660
OptionalScalar< Scalar > dpc_dswe(const Scalar swe) const
The regularized partial derivative of the capillary pressure w.r.t. the saturation.
Definition: vangenuchten.hh:712
void init(const MaterialLaw *m, const std::string &paramGroup)
Initialize the spline.
Definition: vangenuchten.hh:648
OptionalScalar< Scalar > krw(const Scalar swe) const
The regularized relative permeability for the wetting phase.
Definition: vangenuchten.hh:778
OptionalScalar< Scalar > krn(const Scalar swe) const
The regularized relative permeability for the non-wetting phase.
Definition: vangenuchten.hh:808
OptionalScalar< Scalar > swe(const Scalar pc) const
The regularized saturation-capillary pressure curve.
Definition: vangenuchten.hh:730
OptionalScalar< Scalar > dkrn_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the non-wetting phase w....
Definition: vangenuchten.hh:823
OptionalScalar< Scalar > pc(const Scalar swe) const
The regularized capillary pressure-saturation curve regularized part:
Definition: vangenuchten.hh:689
OptionalScalar< Scalar > dswe_dpc(const Scalar pc) const
The regularized partial derivative of the saturation to the capillary pressure.
Definition: vangenuchten.hh:754
OptionalScalar< Scalar > dkrw_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the wetting phase w.r....
Definition: vangenuchten.hh:793
Regularization parameters.
Definition: vangenuchten.hh:580
Scalar krnLowSwe() const
Threshold saturation below which the relative permeability of the non-wetting phase gets regularized.
Definition: vangenuchten.hh:622
void setKrwHighSwe(Scalar krwHighSwe)
Set the threshold saturation above which the relative permeability of the wetting phase gets regulari...
Definition: vangenuchten.hh:629
Scalar krwHighSwe() const
Threshold saturation above which the relative permeability of the wetting phase gets regularized.
Definition: vangenuchten.hh:636
void setKrnLowSwe(Scalar krnLowSwe)
Set the threshold saturation below which the relative permeability of the non-wetting phase gets regu...
Definition: vangenuchten.hh:615
void setPcHighSwe(Scalar pcHighSwe)
Set the threshold saturation above which the capillary pressure is regularized.
Definition: vangenuchten.hh:599
Scalar pcHighSwe() const
Threshold saturation above which the capillary pressure is regularized.
Definition: vangenuchten.hh:608
void setPcLowSwe(Scalar pcLowSwe)
Set the threshold saturation below which the capillary pressure is regularized.
Definition: vangenuchten.hh:587
Scalar pcLowSwe() const
Threshold saturation below which the capillary pressure is regularized.
Definition: vangenuchten.hh:593