3.4
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#include <cmath>
29#include <algorithm>
30
35
36namespace Dumux::FluidMatrix {
37
47{
48
49public:
63 template<class Scalar>
64 struct Params
65 {
66 Params(Scalar alpha, Scalar n, Scalar l = 0.5)
67 : alpha_(alpha), n_(n), m_(1.0 - 1.0/n), l_(l)
68 {}
69
70 Scalar alpha() const { return alpha_; }
71 void setAlpha(Scalar alpha) { alpha_ = alpha; }
72
73 Scalar m() const { return m_; }
74 void setM(Scalar m) { m_ = m; n_ = 1.0/(1.0 - m); }
75
76 Scalar n() const{ return n_; }
77 void setN(Scalar n){ n_ = n; m_ = 1.0 - 1.0/n; }
78
79 Scalar l() const { return l_; }
80 void setL(Scalar l) { l_ = l; }
81
82 bool operator== (const Params& p) const
83 {
84 return Dune::FloatCmp::eq(alpha_, p.alpha_, 1e-6)
85 && Dune::FloatCmp::eq(n_, p.n_, 1e-6)
86 && Dune::FloatCmp::eq(m_, p.m_, 1e-6)
87 && Dune::FloatCmp::eq(l_, p.l_, 1e-6);
88 }
89
90 private:
91 Scalar alpha_, n_, m_, l_;
92 };
93
98 template<class Scalar = double>
99 static Params<Scalar> makeParams(const std::string& paramGroup)
100 {
101 const auto n = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenN");
102 const auto alpha = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenAlpha");
103 // l is usually chosen to be 0.5 (according to Mualem (1976), WRR)
104 const auto l = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenL", 0.5);
105 return Params<Scalar>(alpha, n, l);
106 }
107
121 template<class Scalar>
122 static Scalar pc(Scalar swe, const Params<Scalar>& params)
123 {
124 using std::pow;
125 using std::clamp;
126
127 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
128
129 const Scalar pc = pow(pow(swe, -1.0/params.m()) - 1, 1.0/params.n())/params.alpha();
130 return pc;
131 }
132
147 template<class Scalar>
148 static Scalar swe(Scalar pc, const Params<Scalar>& params)
149 {
150 using std::pow;
151 using std::max;
152
153 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
154
155 const Scalar sw = pow(pow(params.alpha()*pc, params.n()) + 1, -params.m());
156 return sw;
157 }
158
163 template<class Scalar>
164 static Scalar endPointPc(const Params<Scalar>& params)
165 { return 0.0; }
166
183 template<class Scalar>
184 static Scalar dpc_dswe(Scalar swe, const Params<Scalar>& params)
185 {
186 using std::pow;
187 using std::clamp;
188
189 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
190
191 const Scalar powSwe = pow(swe, -1/params.m());
192 return - 1.0/params.alpha() * pow(powSwe - 1, 1.0/params.n() - 1)/params.n()
193 * powSwe/swe/params.m();
194 }
195
205 template<class Scalar>
206 static Scalar dswe_dpc(Scalar pc, const Params<Scalar>& params)
207 {
208 using std::pow;
209 using std::max;
210
211 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
212
213 const Scalar powAlphaPc = pow(params.alpha()*pc, params.n());
214 return -pow(powAlphaPc + 1, -params.m()-1)*params.m()*powAlphaPc/pc*params.n();
215 }
216
227 template<class Scalar>
228 static Scalar krw(Scalar swe, const Params<Scalar>& params)
229 {
230 using std::pow;
231 using std::clamp;
232
233 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
234
235 const Scalar r = 1.0 - pow(1.0 - pow(swe, 1.0/params.m()), params.m());
236 return pow(swe, params.l())*r*r;
237 }
238
249 template<class Scalar>
250 static Scalar dkrw_dswe(Scalar swe, const Params<Scalar>& params)
251 {
252 using std::pow;
253 using std::clamp;
254
255 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
256
257 const Scalar x = 1.0 - pow(swe, 1.0/params.m());
258 const Scalar xToM = pow(x, params.m());
259 return (1.0 - xToM)*pow(swe, params.l()-1) * ( (1.0 - xToM)*params.l() + 2*xToM*(1.0-x)/x );
260 }
261
272 template<class Scalar>
273 static Scalar krn(Scalar swe, const Params<Scalar>& params)
274 {
275 using std::pow;
276 using std::clamp;
277
278 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
279
280 return pow(1 - swe, params.l()) * pow(1 - pow(swe, 1.0/params.m()), 2*params.m());
281 }
282
294 template<class Scalar>
295 static Scalar dkrn_dswe(Scalar swe, const Params<Scalar>& params)
296 {
297 using std::pow;
298 using std::clamp;
299
300 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
301
302 const auto sne = 1.0 - swe;
303 const auto x = 1.0 - pow(swe, 1.0/params.m());
304 return -pow(sne, params.l()-1.0) * pow(x, 2*params.m() - 1.0) * ( params.l()*x + 2.0*sne/swe*(1.0 - x) );
305 }
306};
307
315template <class Scalar>
317{
318public:
320 template<class S>
321 struct Params
322 {
330 { pcLowSwe_ = pcLowSwe; }
331
335 Scalar pcLowSwe() const
336 { return pcLowSwe_; }
337
342 { pcHighSwe_ = pcHighSwe; }
343
350 Scalar pcHighSwe() const
351 { return pcHighSwe_; }
352
358 { krnLowSwe_ = krnLowSwe; }
359
364 Scalar krnLowSwe() const
365 { return krnLowSwe_; }
366
372 { krwHighSwe_ = krwHighSwe; }
373
378 Scalar krwHighSwe() const
379 { return krwHighSwe_; }
380
381private:
382 S pcLowSwe_ = 0.01;
383 S pcHighSwe_ = 0.99;
384 S krnLowSwe_ = 0.1;
385 S krwHighSwe_ = 0.9;
386 };
387
389 template<class MaterialLaw>
390 void init(const MaterialLaw* m, const std::string& paramGroup)
391 {
392 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenPcLowSweThreshold", 0.01);
393 pcHighSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenPcHighSweThreshold", 0.99);
394 krwHighSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenKrwHighSweThreshold", 0.9);
395 krnLowSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenKrnLowSweThreshold", 0.1);
396
397 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
398 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
399 }
400
401 template<class MaterialLaw, class BaseParams, class EffToAbsParams>
402 void init(const MaterialLaw* m, const BaseParams& bp, const EffToAbsParams& etap, const Params<Scalar>& p)
403 {
404 pcLowSwe_ = p.pcLowSwe();
405 pcHighSwe_ = p.pcHighSwe();
406 krwHighSwe_ = p.krwHighSwe();
407 krnLowSwe_ = p.krnLowSwe();
408
409 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
410 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
411 }
412
417 {
418 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
419 && Dune::FloatCmp::eq(pcHighSwe_, o.pcHighSwe_, 1e-6)
420 && Dune::FloatCmp::eq(krwHighSwe_, o.krwHighSwe_, 1e-6)
421 && Dune::FloatCmp::eq(krnLowSwe_, o.krnLowSwe_, 1e-6);
422 }
423
431 OptionalScalar<Scalar> pc(const Scalar swe) const
432 {
433 // make sure that the capillary pressure observes a derivative
434 // != 0 for 'illegal' saturations. This is favourable for the
435 // newton solver (if the derivative is calculated numerically)
436 // in order to get the saturation moving to the right
437 // direction if it temporarily is in an 'illegal' range.
438 if (swe <= pcLowSwe_)
439 return pcLowSwePcValue_ + pcDerivativeLowSw_*(swe - pcLowSwe_);
440
441 else if (swe >= 1.0)
442 return pcDerivativeHighSweEnd_*(swe - 1.0);
443
444 else if (swe > pcHighSwe_)
445 return pcSpline_.eval(swe);
446
447 else
448 return {}; // no regularization
449 }
450
455 {
456 if (swe <= pcLowSwe_)
457 return pcDerivativeLowSw_;
458
459 else if (swe >= 1.0)
460 return pcDerivativeHighSweEnd_;
461
462 else if (swe > pcHighSwe_)
463 return pcSpline_.evalDerivative(swe);
464
465 else
466 return {}; // no regularization
467 }
468
472 OptionalScalar<Scalar> swe(const Scalar pc) const
473 {
474 if (pc <= 0.0)
475 {
476 if (pcHighSwe_ >= 1.0)
477 return 1.0;
478 else
479 return pc/pcDerivativeHighSweEnd_ + 1.0;
480 }
481
482 // invert spline
483 else if (pc <= pcHighSwePcValue_)
484 return pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0, pc);
485
486 else if (pc >= pcLowSwePcValue_)
487 return (pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
488
489 else
490 return {}; // no regularization
491 }
492
497 {
498 if (pc <= 0.0)
499 {
500 if (pcHighSwe_ >= 1.0)
501 return 0.0;
502 else
503 return 1.0/pcDerivativeHighSweEnd_;
504 }
505
506 // derivative of the inverse of the function is one over derivative of the function
507 else if (pc <= pcHighSwePcValue_)
508 return 1.0/pcSpline_.evalDerivative(pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0, pc));
509
510 else if (pc >= pcLowSwePcValue_)
511 return 1.0/pcDerivativeLowSw_;
512
513 else
514 return {}; // no regularization
515 }
516
520 OptionalScalar<Scalar> krw(const Scalar swe) const
521 {
522 if (swe <= 0.0)
523 return 0.0;
524 else if (swe >= 1.0)
525 return 1.0;
526 else if (swe >= krwHighSwe_)
527 return krwSpline_.eval(swe);
528 else
529 return {}; // no regularization
530 }
531
536 {
537 if (swe <= 0.0)
538 return 0.0;
539 else if (swe >= 1.0)
540 return 0.0;
541 else if (swe >= krwHighSwe_)
542 return krwSpline_.evalDerivative(swe);
543 else
544 return {}; // no regularization
545 }
546
550 OptionalScalar<Scalar> krn(const Scalar swe) const
551 {
552 if (swe <= 0.0)
553 return 1.0;
554 else if (swe >= 1.0)
555 return 0.0;
556 else if (swe <= krnLowSwe_)
557 return krnSpline_.eval(swe);
558 else
559 return {}; // no regularization
560 }
561
566 {
567 if (swe <= 0.0)
568 return 0.0;
569 else if (swe >= 1.0)
570 return 0.0;
571 else if (swe <= krnLowSwe_)
572 return krnSpline_.evalDerivative(swe);
573 else
574 return {}; // no regularization
575 }
576
577private:
578 template<class MaterialLaw>
579 void initPcParameters_(const MaterialLaw* m, const Scalar lowSwe, const Scalar highSwe)
580 {
581 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
582 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
583 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
584
585 pcDerivativeLowSw_ = m->template dpc_dsw<false>(lowSw)*dsw_dswe;
586
587 pcDerivativeHighSweThreshold_ = m->template dpc_dsw<false>(highSw)*dsw_dswe;
588 pcDerivativeHighSweEnd_ = 2.0*(0.0 - m->template pc<false>(highSw))/(1.0 - highSwe);
589
590 pcLowSwePcValue_ = m->template pc<false>(lowSw);
591 pcHighSwePcValue_ = m->template pc<false>(highSw);
592
593 // Only initialize regularization spline if given parameters are in
594 // the admissible range. When constructing with non-sensible parameters
595 // the spline construction might fail (e.g. highSwe == 1.0)
596 if (highSwe < 1.0)
597 pcSpline_ = Spline<Scalar>(highSwe, 1.0, // x0, x1
598 pcHighSwePcValue_, 0, // y0, y1
599 pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_); // m0, m1
600 }
601
602 template<class MaterialLaw>
603 void initKrParameters_(const MaterialLaw* m, const Scalar lowSwe, const Scalar highSwe)
604 {
605 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
606 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
607 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
608
609 const auto krwHighSw = m->template krw<false>(highSw);
610 const auto dkrwHighSw = m->template dkrw_dsw<false>(highSw)*dsw_dswe;
611
612 const auto krnLowSw = m->template krn<false>(lowSw);
613 const auto dkrnLowSw = m->template dkrn_dsw<false>(lowSw)*dsw_dswe;
614
615 if (highSwe < 1.0)
616 krwSpline_ = Spline<Scalar>(highSwe, 1.0, // x0, x1
617 krwHighSw, 1.0, // y0, y1
618 dkrwHighSw, 0.0); // m0, m1
619
620 if (lowSwe > 0.0)
621 krnSpline_ = Spline<Scalar>(0.0, lowSwe, // x0, x1
622 1.0, krnLowSw, // y0, y1
623 0.0, dkrnLowSw); // m0, m1
624 }
625
626 Scalar pcLowSwe_, pcHighSwe_;
627 Scalar pcLowSwePcValue_, pcHighSwePcValue_;
628 Scalar krwHighSwe_, krnLowSwe_;
629 Scalar pcDerivativeLowSw_;
630 Scalar pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_;
631
632 Spline<Scalar> pcSpline_;
633 Spline<Scalar> krwSpline_;
634 Spline<Scalar> krnSpline_;
635};
636
641template<typename Scalar = double>
643
648template<typename Scalar = double>
650
651} // end namespace Dumux::FluidMatrix
652
653#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...
Definition: brookscorey.hh:35
A 3rd order polynomial spline.
Definition: spline.hh:55
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
Definition: efftoabsdefaultpolicy.hh:46
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:57
Implementation of the van Genuchten capillary pressure <-> saturation relation, and relative permeabi...
Definition: vangenuchten.hh:47
static Scalar dswe_dpc(Scalar pc, const Params< Scalar > &params)
The partial derivative of the effective saturation to the capillary pressure according to van Genucht...
Definition: vangenuchten.hh:206
static Scalar endPointPc(const Params< Scalar > &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: vangenuchten.hh:164
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability for the non-wetting phase in regard to the wetting satura...
Definition: vangenuchten.hh:295
static Scalar swe(Scalar pc, const Params< Scalar > &params)
The saturation-capillary pressure curve according to van Genuchten.
Definition: vangenuchten.hh:148
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:228
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:250
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:184
static Scalar pc(Scalar swe, const Params< Scalar > &params)
The capillary pressure-saturation curve according to van Genuchten.
Definition: vangenuchten.hh:122
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:273
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: vangenuchten.hh:99
The parameter type.
Definition: vangenuchten.hh:65
Params(Scalar alpha, Scalar n, Scalar l=0.5)
Definition: vangenuchten.hh:66
Scalar alpha() const
Definition: vangenuchten.hh:70
bool operator==(const Params &p) const
Definition: vangenuchten.hh:82
Scalar l() const
Definition: vangenuchten.hh:79
void setAlpha(Scalar alpha)
Definition: vangenuchten.hh:71
void setL(Scalar l)
Definition: vangenuchten.hh:80
void setM(Scalar m)
Definition: vangenuchten.hh:74
Scalar n() const
Definition: vangenuchten.hh:76
Scalar m() const
Definition: vangenuchten.hh:73
void setN(Scalar n)
Definition: vangenuchten.hh:77
A regularization for the VanGenuchten material law.
Definition: vangenuchten.hh:317
bool operator==(const VanGenuchtenRegularization &o) const
Equality comparison with another instance.
Definition: vangenuchten.hh:416
void init(const MaterialLaw *m, const BaseParams &bp, const EffToAbsParams &etap, const Params< Scalar > &p)
Definition: vangenuchten.hh:402
OptionalScalar< Scalar > dpc_dswe(const Scalar swe) const
The regularized partial derivative of the capillary pressure w.r.t. the saturation.
Definition: vangenuchten.hh:454
void init(const MaterialLaw *m, const std::string &paramGroup)
Initialize the spline.
Definition: vangenuchten.hh:390
OptionalScalar< Scalar > krw(const Scalar swe) const
The regularized relative permeability for the wetting phase.
Definition: vangenuchten.hh:520
OptionalScalar< Scalar > krn(const Scalar swe) const
The regularized relative permeability for the non-wetting phase.
Definition: vangenuchten.hh:550
OptionalScalar< Scalar > swe(const Scalar pc) const
The regularized saturation-capillary pressure curve.
Definition: vangenuchten.hh:472
OptionalScalar< Scalar > dkrn_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the non-wetting phase w....
Definition: vangenuchten.hh:565
OptionalScalar< Scalar > pc(const Scalar swe) const
The regularized capillary pressure-saturation curve regularized part:
Definition: vangenuchten.hh:431
OptionalScalar< Scalar > dswe_dpc(const Scalar pc) const
The regularized partial derivative of the saturation to the capillary pressure.
Definition: vangenuchten.hh:496
OptionalScalar< Scalar > dkrw_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the wetting phase w.r....
Definition: vangenuchten.hh:535
Regularization parameters.
Definition: vangenuchten.hh:322
Scalar krnLowSwe() const
Threshold saturation below which the relative permeability of the non-wetting phase gets regularized.
Definition: vangenuchten.hh:364
void setKrwHighSwe(Scalar krwHighSwe)
Set the threshold saturation above which the relative permeability of the wetting phase gets regulari...
Definition: vangenuchten.hh:371
Scalar krwHighSwe() const
Threshold saturation above which the relative permeability of the wetting phase gets regularized.
Definition: vangenuchten.hh:378
void setKrnLowSwe(Scalar krnLowSwe)
Set the threshold saturation below which the relative permeability of the non-wetting phase gets regu...
Definition: vangenuchten.hh:357
void setPcHighSwe(Scalar pcHighSwe)
Set the threshold saturation above which the capillary pressure is regularized.
Definition: vangenuchten.hh:341
Scalar pcHighSwe() const
Threshold saturation above which the capillary pressure is regularized.
Definition: vangenuchten.hh:350
void setPcLowSwe(Scalar pcLowSwe)
Set the threshold saturation below which the capillary pressure is regularized.
Definition: vangenuchten.hh:329
Scalar pcLowSwe() const
Threshold saturation below which the capillary pressure is regularized.
Definition: vangenuchten.hh:335