GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/material/fluidmatrixinteractions/2p/vangenuchten.hh
Date: 2024-05-04 19:09:25
Exec Total Coverage
Lines: 153 183 83.6%
Functions: 22 36 61.1%
Branches: 84 106 79.2%

Line Branch Exec Source
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 //
7 /*!
8 * \file
9 * \ingroup Fluidmatrixinteractions
10 * \brief Implementation of the capillary pressure and
11 * relative permeability <-> saturation relations according to van Genuchten.
12 */
13 #ifndef DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH
14 #define DUMUX_MATERIAL_FLUIDMATRIX_VAN_GENUCHTEN_HH
15
16 #include <cmath>
17 #include <algorithm>
18
19 #include <dumux/common/parameters.hh>
20 #include <dumux/common/spline.hh>
21 #include <dumux/common/optionalscalar.hh>
22 #include <dumux/material/fluidmatrixinteractions/2p/materiallaw.hh>
23
24 namespace Dumux::FluidMatrix {
25
26 /*!
27 * \ingroup Fluidmatrixinteractions
28 * \brief Implementation of the van Genuchten capillary pressure <->
29 * saturation relation, and relative permeability.
30 *
31 * \note Capillary pressure model from van Genuchten (1980),
32 * relative permeability model from Mualem (1976)
33 */
34 class VanGenuchten
35 {
36
37 public:
38 /*!
39 * \brief The parameter type
40 * \tparam Scalar The scalar type
41 * \note The van Genuchten laws are parameterized with four parameters: \f$\mathrm{n, m, \alpha, l}\f$.
42 *
43 * - \f$\mathrm{\alpha}\f$ shape parameter \f$\mathrm{[1/Pa]}\f$
44 * - \f$\mathrm{m}\f$ shape parameter \f$\mathrm{[-]}\f$
45 * - \f$\mathrm{n}\f$ shape parameter \f$\mathrm{[-]}\f$
46 * - \f$\mathrm{l}\f$ pore-connectivity parameter \f$\mathrm{[-]}\f$ of Mualem's relative permeability curve
47 *
48 * \note In the original Mualem (1976) paper the pore-connectivity parameter is called "n". It's referred to as "l" in
49 * several later publication of van Genuchten, e.g. van Genuchten (1991), Shaap & van Genuchten (2006).
50 */
51 template<class Scalar>
52 struct Params
53 {
54 135 Params(Scalar alpha, Scalar n, Scalar l = 0.5)
55 135 : alpha_(alpha), n_(n), m_(1.0 - 1.0/n), l_(l)
56 {}
57
58 Scalar alpha() const { return alpha_; }
59 void setAlpha(Scalar alpha) { alpha_ = alpha; }
60
61 Scalar m() const { return m_; }
62 void setM(Scalar m) { m_ = m; n_ = 1.0/(1.0 - m); }
63
64 Scalar n() const{ return n_; }
65 void setN(Scalar n){ n_ = n; m_ = 1.0 - 1.0/n; }
66
67 Scalar l() const { return l_; }
68 void setL(Scalar l) { l_ = l; }
69
70 4527 bool operator== (const Params& p) const
71 {
72
2/2
✓ Branch 0 taken 63 times.
✓ Branch 1 taken 4464 times.
4527 return Dune::FloatCmp::eq(alpha_, p.alpha_, 1e-6)
73
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 4401 times.
✓ Branch 2 taken 4401 times.
✗ Branch 3 not taken.
4401 && Dune::FloatCmp::eq(n_, p.n_, 1e-6)
74
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 4401 times.
✓ Branch 2 taken 4401 times.
✗ Branch 3 not taken.
4401 && Dune::FloatCmp::eq(m_, p.m_, 1e-6)
75
4/6
✓ Branch 0 taken 4401 times.
✓ Branch 1 taken 126 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 4401 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 4401 times.
8928 && Dune::FloatCmp::eq(l_, p.l_, 1e-6);
76 }
77
78 private:
79 Scalar alpha_, n_, m_, l_;
80 };
81
82 /*!
83 * \brief Construct from a subgroup from the global parameter tree
84 * \note This will give you nice error messages if a mandatory parameter is missing
85 */
86 template<class Scalar = double>
87 134 static Params<Scalar> makeParams(const std::string& paramGroup)
88 {
89 134 const auto n = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenN");
90 134 const auto alpha = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenAlpha");
91 // l is usually chosen to be 0.5 (according to Mualem (1976), WRR)
92 134 const auto l = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenL", 0.5);
93 268 return Params<Scalar>(alpha, n, l);
94 }
95
96 /*!
97 * \brief The capillary pressure-saturation curve according to van Genuchten.
98 *
99 * Van Genuchten's empirical capillary pressure <-> saturation
100 * function is given by
101 * \f$\mathrm{
102 p_c = (\overline{S}_w^{-1/m} - 1)^{1/n}/\alpha
103 }\f$
104 * \param swe Effective saturation of the wetting phase \f$\mathrm{\overline{S}_w}\f$
105 * \param params A container object that is populated with the appropriate coefficients for the respective law.
106 * \note Instead of undefined behaviour if swe is not in the valid range, we return a valid number,
107 * by clamping the input. Note that pc(swe = 0.0) = inf, have a look at RegularizedVanGenuchten if this is a problem.
108 */
109 template<class Scalar>
110 13014759 static Scalar pc(Scalar swe, const Params<Scalar>& params)
111 {
112 using std::pow;
113 using std::clamp;
114
115
1/2
✓ Branch 0 taken 13014759 times.
✗ Branch 1 not taken.
13014759 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
116
117 13014759 const Scalar pc = pow(pow(swe, -1.0/params.m()) - 1, 1.0/params.n())/params.alpha();
118 13014759 return pc;
119 }
120
121 /*!
122 * \brief The saturation-capillary pressure curve according to van Genuchten.
123 *
124 * This is the inverse of the capillary pressure-saturation curve:
125 * \f$\mathrm{
126 \overline{S}_w = {p_c}^{-1} = ((\alpha p_c)^n + 1)^{-m}
127 }\f$
128 *
129 * \param pc Capillary pressure \f$\mathrm{p_c}\f$ in \f$\mathrm{[Pa]}\f$
130 * \param params A container object that is populated with the appropriate coefficients for the respective law.
131 * \return The effective saturation of the wetting phase \f$\mathrm{\overline{S}_w}\f$
132 * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
133 * i.e. sw(pc < 0.0) = 0.0, by clamping the input to the physical bounds.
134 */
135 template<class Scalar>
136 793449694 static Scalar swe(Scalar pc, const Params<Scalar>& params)
137 {
138 using std::pow;
139 using std::max;
140
141
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 793449694 times.
793449694 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
142
143 793449694 const Scalar sw = pow(pow(params.alpha()*pc, params.n()) + 1, -params.m());
144 793449694 return sw;
145 }
146
147 /*!
148 * \brief The capillary pressure at Swe = 1.0 also called end point capillary pressure
149 * \param params A container object that is populated with the appropriate coefficients for the respective law.
150 */
151 template<class Scalar>
152 static Scalar endPointPc(const Params<Scalar>& params)
153 { return 0.0; }
154
155 /*!
156 * \brief The partial derivative of the capillary
157 * pressure w.r.t. the effective saturation according to van Genuchten.
158 *
159 * This is equivalent to
160 * \f$\mathrm{
161 \frac{\partial p_c}{\partial \overline{S}_w} =
162 -\frac{1}{\alpha} (\overline{S}_w^{-1/m} - 1)^{1/n - }
163 \overline{S}_w^{-1/m} / \overline{S}_w / m
164 }\f$
165 *
166 * \param swe Effective saturation of the wetting phase \f$\mathrm{\overline{S}_w}\f$
167 * \param params A container object that is populated with the appropriate coefficients for the respective law.
168 * \note Instead of undefined behaviour if swe is not in the valid range, we return a valid number,
169 * by clamping the input.
170 */
171 template<class Scalar>
172 11658390 static Scalar dpc_dswe(Scalar swe, const Params<Scalar>& params)
173 {
174 using std::pow;
175 using std::clamp;
176
177
1/2
✓ Branch 0 taken 11658390 times.
✗ Branch 1 not taken.
11658390 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
178
179 11658390 const Scalar powSwe = pow(swe, -1/params.m());
180 11658390 return - 1.0/params.alpha() * pow(powSwe - 1, 1.0/params.n() - 1)/params.n()
181 11658390 * powSwe/swe/params.m();
182 }
183
184 /*!
185 * \brief The partial derivative of the effective
186 * saturation to the capillary pressure according to van Genuchten.
187 *
188 * \param pc Capillary pressure \f$\mathrm{p_C}\f$ in \f$\mathrm{[Pa]}\f$
189 * \param params A container object that is populated with the appropriate coefficients for the respective law.
190 * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
191 * by clamping the input.
192 */
193 template<class Scalar>
194 38360413 static Scalar dswe_dpc(Scalar pc, const Params<Scalar>& params)
195 {
196 using std::pow;
197 using std::max;
198
199
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 38360413 times.
38360413 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
200
201 38360413 const Scalar powAlphaPc = pow(params.alpha()*pc, params.n());
202 38360413 return -pow(powAlphaPc + 1, -params.m()-1)*params.m()*powAlphaPc/pc*params.n();
203 }
204
205 /*!
206 * \brief The relative permeability for the wetting phase of
207 * the medium implied by van Genuchten / Mualem
208 * parameterization.
209 *
210 * \param swe The mobile saturation of the wetting phase.
211 * \param params A container object that is populated with the appropriate coefficients for the respective law.
212 * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
213 * by clamping the input.
214 */
215 template<class Scalar>
216 820907417 static Scalar krw(Scalar swe, const Params<Scalar>& params)
217 {
218 using std::pow;
219 using std::clamp;
220
221
1/2
✓ Branch 0 taken 820907417 times.
✗ Branch 1 not taken.
820907417 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
222
223 820907417 const Scalar r = 1.0 - pow(1.0 - pow(swe, 1.0/params.m()), params.m());
224 820907417 return pow(swe, params.l())*r*r;
225 }
226
227 /*!
228 * \brief The derivative of the relative permeability for the
229 * wetting phase in regard to the wetting saturation of the
230 * medium implied by the van Genuchten parameterization.
231 *
232 * \param swe The mobile saturation of the wetting phase.
233 * \param params A container object that is populated with the appropriate coefficients for the respective law.
234 * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
235 * by clamping the input.
236 */
237 template<class Scalar>
238 29567654 static Scalar dkrw_dswe(Scalar swe, const Params<Scalar>& params)
239 {
240 using std::pow;
241 using std::clamp;
242
243
1/2
✓ Branch 0 taken 29567654 times.
✗ Branch 1 not taken.
29567654 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
244
245 29567654 const Scalar x = 1.0 - pow(swe, 1.0/params.m());
246 29567654 const Scalar xToM = pow(x, params.m());
247 29567654 return (1.0 - xToM)*pow(swe, params.l()-1) * ( (1.0 - xToM)*params.l() + 2*xToM*(1.0-x)/x );
248 }
249
250 /*!
251 * \brief The relative permeability for the non-wetting phase
252 * of the medium implied by van Genuchten's
253 * parameterization.
254 *
255 * \param swe The mobile saturation of the wetting phase.
256 * \param params A container object that is populated with the appropriate coefficients for the respective law.
257 * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
258 * by clamping the input.
259 */
260 template<class Scalar>
261 17103047 static Scalar krn(Scalar swe, const Params<Scalar>& params)
262 {
263 using std::pow;
264 using std::clamp;
265
266
1/2
✓ Branch 0 taken 17103047 times.
✗ Branch 1 not taken.
17103047 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
267
268 17103047 return pow(1 - swe, params.l()) * pow(1 - pow(swe, 1.0/params.m()), 2*params.m());
269 }
270
271 /*!
272 * \brief The derivative of the relative permeability for the
273 * non-wetting phase in regard to the wetting saturation of
274 * the medium as implied by the van Genuchten
275 * parameterization.
276 *
277 * \param swe The mobile saturation of the wetting phase.
278 * \param params A container object that is populated with the appropriate coefficients for the respective law.
279 * \note Instead of undefined behaviour if pc is not in the valid range, we return a valid number,
280 * by clamping the input.
281 */
282 template<class Scalar>
283 1123664 static Scalar dkrn_dswe(Scalar swe, const Params<Scalar>& params)
284 {
285 using std::pow;
286 using std::clamp;
287
288
1/2
✓ Branch 0 taken 1123664 times.
✗ Branch 1 not taken.
1123664 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
289
290 1123664 const auto sne = 1.0 - swe;
291 1123664 const auto x = 1.0 - pow(swe, 1.0/params.m());
292 1123664 return -pow(sne, params.l()-1.0) * pow(x, 2*params.m() - 1.0) * ( params.l()*x + 2.0*sne/swe*(1.0 - x) );
293 }
294 };
295
296 /*!
297 * \ingroup Fluidmatrixinteractions
298 * \brief A regularization for the VanGenuchten material law
299 * \note Regularization can be turned of by setting the threshold parameters
300 * out of range (runtime) or by replacing
301 * this class by NoRegularization (compile time).
302 */
303 template <class Scalar>
304 480 class VanGenuchtenRegularization
305 {
306 public:
307 //! Regularization parameters
308 template<class S>
309 struct Params
310 {
311 /*!
312 * \brief Set the threshold saturation below which the capillary pressure is regularized.
313 *
314 * Most problems are very sensitive to this value (e.g. making it smaller might
315 * result in very high capillary pressures)
316 */
317 void setPcLowSwe(Scalar pcLowSwe)
318 { pcLowSwe_ = pcLowSwe; }
319
320 /*!
321 * \brief Threshold saturation below which the capillary pressure is regularized.
322 */
323 Scalar pcLowSwe() const
324 { return pcLowSwe_; }
325
326 /*!
327 * \brief Set the threshold saturation above which the capillary pressure is regularized.
328 */
329 void setPcHighSwe(Scalar pcHighSwe)
330 { pcHighSwe_ = pcHighSwe; }
331
332 /*!
333 * \brief Threshold saturation above which the capillary pressure is regularized.
334 *
335 * Most problems are very sensitive to this value (e.g. making it smaller might
336 * result in negative capillary pressures).
337 */
338 Scalar pcHighSwe() const
339 { return pcHighSwe_; }
340
341 /*!
342 * \brief Set the threshold saturation below which the relative
343 * permeability of the non-wetting phase gets regularized.
344 */
345 void setKrnLowSwe(Scalar krnLowSwe)
346 { krnLowSwe_ = krnLowSwe; }
347
348 /*!
349 * \brief Threshold saturation below which the relative
350 * permeability of the non-wetting phase gets regularized.
351 */
352 Scalar krnLowSwe() const
353 { return krnLowSwe_; }
354
355 /*!
356 * \brief Set the threshold saturation above which the relative
357 * permeability of the wetting phase gets regularized.
358 */
359 void setKrwHighSwe(Scalar krwHighSwe)
360 { krwHighSwe_ = krwHighSwe; }
361
362 /*!
363 * \brief Threshold saturation above which the relative
364 * permeability of the wetting phase gets regularized.
365 */
366 Scalar krwHighSwe() const
367 { return krwHighSwe_; }
368
369 private:
370 S pcLowSwe_ = 0.01;
371 S pcHighSwe_ = 0.99;
372 S krnLowSwe_ = 0.1;
373 S krwHighSwe_ = 0.9;
374 };
375
376 //! Initialize the spline
377 template<class MaterialLaw>
378 119 void init(const MaterialLaw* m, const std::string& paramGroup)
379 {
380 119 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenPcLowSweThreshold", 0.01);
381 119 pcHighSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenPcHighSweThreshold", 0.99);
382 119 krwHighSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenKrwHighSweThreshold", 0.9);
383 119 krnLowSwe_ = getParamFromGroup<Scalar>(paramGroup, "VanGenuchtenKrnLowSweThreshold", 0.1);
384
385 119 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
386 119 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
387 119 }
388
389 template<class MaterialLaw, class BaseParams, class EffToAbsParams>
390 void init(const MaterialLaw* m, const BaseParams& bp, const EffToAbsParams& etap, const Params<Scalar>& p)
391 {
392 pcLowSwe_ = p.pcLowSwe();
393 pcHighSwe_ = p.pcHighSwe();
394 krwHighSwe_ = p.krwHighSwe();
395 krnLowSwe_ = p.krnLowSwe();
396
397 initPcParameters_(m, pcLowSwe_, pcHighSwe_);
398 initKrParameters_(m, krnLowSwe_, krwHighSwe_);
399 }
400
401 /*!
402 * \brief Equality comparison with another instance
403 */
404 4401 bool operator== (const VanGenuchtenRegularization& o) const
405 {
406
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 4401 times.
4401 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
407
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 4401 times.
✓ Branch 2 taken 4401 times.
✗ Branch 3 not taken.
4401 && Dune::FloatCmp::eq(pcHighSwe_, o.pcHighSwe_, 1e-6)
408
2/4
✗ Branch 0 not taken.
✓ Branch 1 taken 4401 times.
✓ Branch 2 taken 4401 times.
✗ Branch 3 not taken.
4401 && Dune::FloatCmp::eq(krwHighSwe_, o.krwHighSwe_, 1e-6)
409
3/6
✓ Branch 0 taken 4401 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 4401 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 4401 times.
8802 && Dune::FloatCmp::eq(krnLowSwe_, o.krnLowSwe_, 1e-6);
410 }
411
412 /*!
413 * \brief The regularized capillary pressure-saturation curve
414 * regularized part:
415 * - low saturation: extend the \f$\mathrm{p_c(S_w)}\f$ curve with the slope at the regularization point (i.e. no kink).
416 * - high saturation: connect the high regularization point with \f$\mathrm{\overline{S}_w =1}\f$
417 * with a spline and continue linearly for \f$\mathrm{\overline{S}_w > 1}\f$
418 */
419 35990053 OptionalScalar<Scalar> pc(const Scalar swe) const
420 {
421 // make sure that the capillary pressure observes a derivative
422 // != 0 for 'illegal' saturations. This is favourable for the
423 // newton solver (if the derivative is calculated numerically)
424 // in order to get the saturation moving to the right
425 // direction if it temporarily is in an 'illegal' range.
426
2/2
✓ Branch 0 taken 1001636 times.
✓ Branch 1 taken 34988417 times.
35990053 if (swe <= pcLowSwe_)
427 1001636 return pcLowSwePcValue_ + pcDerivativeLowSw_*(swe - pcLowSwe_);
428
429
2/2
✓ Branch 0 taken 18421512 times.
✓ Branch 1 taken 16566905 times.
34988417 else if (swe >= 1.0)
430 18421512 return pcDerivativeHighSweEnd_*(swe - 1.0);
431
432
2/2
✓ Branch 0 taken 6922229 times.
✓ Branch 1 taken 9644676 times.
16566905 else if (swe > pcHighSwe_)
433 6922229 return pcSpline_.eval(swe);
434
435 else
436 9644676 return {}; // no regularization
437 }
438
439 /*!
440 * \brief The regularized partial derivative of the capillary pressure w.r.t. the saturation
441 */
442 1743816 OptionalScalar<Scalar> dpc_dswe(const Scalar swe) const
443 {
444
2/2
✓ Branch 0 taken 22 times.
✓ Branch 1 taken 1743794 times.
1743816 if (swe <= pcLowSwe_)
445 22 return pcDerivativeLowSw_;
446
447
2/2
✓ Branch 0 taken 1566190 times.
✓ Branch 1 taken 177604 times.
1743794 else if (swe >= 1.0)
448 1566190 return pcDerivativeHighSweEnd_;
449
450
2/2
✓ Branch 0 taken 85488 times.
✓ Branch 1 taken 92116 times.
177604 else if (swe > pcHighSwe_)
451 85488 return pcSpline_.evalDerivative(swe);
452
453 else
454 92116 return {}; // no regularization
455 }
456
457 /*!
458 * \brief The regularized saturation-capillary pressure curve
459 */
460 810745599 OptionalScalar<Scalar> swe(const Scalar pc) const
461 {
462
2/2
✓ Branch 0 taken 19760766 times.
✓ Branch 1 taken 790984833 times.
810745599 if (pc <= 0.0)
463 {
464
2/2
✓ Branch 0 taken 19120175 times.
✓ Branch 1 taken 640591 times.
19760766 if (pcHighSwe_ >= 1.0)
465 19120175 return 1.0;
466 else
467 640591 return pc/pcDerivativeHighSweEnd_ + 1.0;
468 }
469
470 // invert spline
471
2/2
✓ Branch 0 taken 55773 times.
✓ Branch 1 taken 790929060 times.
790984833 else if (pc <= pcHighSwePcValue_)
472 55773 return pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0, pc);
473
474
2/2
✓ Branch 0 taken 3441979 times.
✓ Branch 1 taken 787487081 times.
790929060 else if (pc >= pcLowSwePcValue_)
475 3441979 return (pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
476
477 else
478 787487081 return {}; // no regularization
479 }
480
481 /*!
482 * \brief The regularized partial derivative of the saturation to the capillary pressure
483 */
484 41607302 OptionalScalar<Scalar> dswe_dpc(const Scalar pc) const
485 {
486
2/2
✓ Branch 0 taken 3752283 times.
✓ Branch 1 taken 37855019 times.
41607302 if (pc <= 0.0)
487 {
488
2/2
✓ Branch 0 taken 3752263 times.
✓ Branch 1 taken 20 times.
3752283 if (pcHighSwe_ >= 1.0)
489 3752263 return 0.0;
490 else
491 20 return 1.0/pcDerivativeHighSweEnd_;
492 }
493
494 // derivative of the inverse of the function is one over derivative of the function
495
2/2
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 37855017 times.
37855019 else if (pc <= pcHighSwePcValue_)
496 2 return 1.0/pcSpline_.evalDerivative(pcSpline_.intersectInterval(pcHighSwe_, 1.0, 0.0, 0.0, 0.0, pc));
497
498
2/2
✓ Branch 0 taken 494606 times.
✓ Branch 1 taken 37360411 times.
37855017 else if (pc >= pcLowSwePcValue_)
499 494606 return 1.0/pcDerivativeLowSw_;
500
501 else
502 37360411 return {}; // no regularization
503 }
504
505 /*!
506 * \brief The regularized relative permeability for the wetting phase
507 */
508 843935702 OptionalScalar<Scalar> krw(const Scalar swe) const
509 {
510
2/2
✓ Branch 0 taken 2917436 times.
✓ Branch 1 taken 841018266 times.
843935702 if (swe <= 0.0)
511 2917436 return 0.0;
512
2/2
✓ Branch 0 taken 36957808 times.
✓ Branch 1 taken 804060458 times.
841018266 else if (swe >= 1.0)
513 36957808 return 1.0;
514
2/2
✓ Branch 0 taken 10613293 times.
✓ Branch 1 taken 793447165 times.
804060458 else if (swe >= krwHighSwe_)
515 10613293 return krwSpline_.eval(swe);
516 else
517 793447165 return {}; // no regularization
518 }
519
520 /*!
521 * \brief The regularized derivative of the relative permeability for the wetting phase w.r.t. saturation
522 */
523 34307830 OptionalScalar<Scalar> dkrw_dswe(const Scalar swe) const
524 {
525
2/2
✓ Branch 0 taken 375168 times.
✓ Branch 1 taken 33932662 times.
34307830 if (swe <= 0.0)
526 375168 return 0.0;
527
2/2
✓ Branch 0 taken 4057722 times.
✓ Branch 1 taken 29874940 times.
33932662 else if (swe >= 1.0)
528 4057722 return 0.0;
529
2/2
✓ Branch 0 taken 1307408 times.
✓ Branch 1 taken 28567532 times.
29874940 else if (swe >= krwHighSwe_)
530 1307408 return krwSpline_.evalDerivative(swe);
531 else
532 28567532 return {}; // no regularization
533 }
534
535 /*!
536 * \brief The regularized relative permeability for the non-wetting phase
537 */
538 33824794 OptionalScalar<Scalar> krn(const Scalar swe) const
539 {
540
2/2
✓ Branch 0 taken 485115 times.
✓ Branch 1 taken 33339679 times.
33824794 if (swe <= 0.0)
541 485115 return 1.0;
542
2/2
✓ Branch 0 taken 17298757 times.
✓ Branch 1 taken 16040922 times.
33339679 else if (swe >= 1.0)
543 17298757 return 0.0;
544
2/2
✓ Branch 0 taken 307609 times.
✓ Branch 1 taken 15733313 times.
16040922 else if (swe <= krnLowSwe_)
545 307609 return krnSpline_.eval(swe);
546 else
547 15733313 return {}; // no regularization
548 }
549
550 /*!
551 * \brief The regularized derivative of the relative permeability for the non-wetting phase w.r.t. saturation
552 */
553 1190756 OptionalScalar<Scalar> dkrn_dswe(const Scalar swe) const
554 {
555
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 1190746 times.
1190756 if (swe <= 0.0)
556 10 return 0.0;
557
2/2
✓ Branch 0 taken 1067196 times.
✓ Branch 1 taken 123550 times.
1190746 else if (swe >= 1.0)
558 1067196 return 0.0;
559
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 123542 times.
123550 else if (swe <= krnLowSwe_)
560 8 return krnSpline_.evalDerivative(swe);
561 else
562 123542 return {}; // no regularization
563 }
564
565 private:
566 template<class MaterialLaw>
567 120 void initPcParameters_(const MaterialLaw* m, const Scalar lowSwe, const Scalar highSwe)
568 {
569 240 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
570 240 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
571 240 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
572
573 120 pcDerivativeLowSw_ = m->template dpc_dsw<false>(lowSw)*dsw_dswe;
574
575 120 pcDerivativeHighSweThreshold_ = m->template dpc_dsw<false>(highSw)*dsw_dswe;
576 120 pcDerivativeHighSweEnd_ = 2.0*(0.0 - m->template pc<false>(highSw))/(1.0 - highSwe);
577
578 120 pcLowSwePcValue_ = m->template pc<false>(lowSw);
579 120 pcHighSwePcValue_ = m->template pc<false>(highSw);
580
581 // Only initialize regularization spline if given parameters are in
582 // the admissible range. When constructing with non-sensible parameters
583 // the spline construction might fail (e.g. highSwe == 1.0)
584
2/2
✓ Branch 0 taken 107 times.
✓ Branch 1 taken 13 times.
120 if (highSwe < 1.0)
585 107 pcSpline_ = Spline<Scalar>(highSwe, 1.0, // x0, x1
586 pcHighSwePcValue_, 0, // y0, y1
587 pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_); // m0, m1
588 120 }
589
590 template<class MaterialLaw>
591 120 void initKrParameters_(const MaterialLaw* m, const Scalar lowSwe, const Scalar highSwe)
592 {
593 240 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
594 240 const auto highSw = MaterialLaw::EffToAbs::sweToSw(highSwe, m->effToAbsParams());
595 240 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
596
597 120 const auto krwHighSw = m->template krw<false>(highSw);
598 120 const auto dkrwHighSw = m->template dkrw_dsw<false>(highSw)*dsw_dswe;
599
600 120 const auto krnLowSw = m->template krn<false>(lowSw);
601 120 const auto dkrnLowSw = m->template dkrn_dsw<false>(lowSw)*dsw_dswe;
602
603
2/2
✓ Branch 0 taken 112 times.
✓ Branch 1 taken 8 times.
120 if (highSwe < 1.0)
604 112 krwSpline_ = Spline<Scalar>(highSwe, 1.0, // x0, x1
605 krwHighSw, 1.0, // y0, y1
606 dkrwHighSw, 0.0); // m0, m1
607
608
2/2
✓ Branch 0 taken 118 times.
✓ Branch 1 taken 2 times.
120 if (lowSwe > 0.0)
609 118 krnSpline_ = Spline<Scalar>(0.0, lowSwe, // x0, x1
610 1.0, krnLowSw, // y0, y1
611 0.0, dkrnLowSw); // m0, m1
612 120 }
613
614 Scalar pcLowSwe_, pcHighSwe_;
615 Scalar pcLowSwePcValue_, pcHighSwePcValue_;
616 Scalar krwHighSwe_, krnLowSwe_;
617 Scalar pcDerivativeLowSw_;
618 Scalar pcDerivativeHighSweThreshold_, pcDerivativeHighSweEnd_;
619
620 Spline<Scalar> pcSpline_;
621 Spline<Scalar> krwSpline_;
622 Spline<Scalar> krnSpline_;
623 };
624
625 /*!
626 * \ingroup Fluidmatrixinteractions
627 * \brief A default configuration for using the VanGenuchten material law
628 */
629 template<typename Scalar = double>
630 using VanGenuchtenDefault = TwoPMaterialLaw<Scalar, VanGenuchten, VanGenuchtenRegularization<Scalar>, TwoPEffToAbsDefaultPolicy>;
631
632 /*!
633 * \ingroup Fluidmatrixinteractions
634 * \brief A default configuration without regularization for using the VanGenuchten material law
635 */
636 template<typename Scalar = double>
637 using VanGenuchtenNoReg = TwoPMaterialLaw<Scalar, VanGenuchten, NoRegularization, TwoPEffToAbsDefaultPolicy>;
638
639 } // end namespace Dumux::FluidMatrix
640
641 #endif
642