3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
regularizedvangenuchten.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 REGULARIZED_VAN_GENUCHTEN_HH
26#define REGULARIZED_VAN_GENUCHTEN_HH
27
28#warning "This header is deprecated. Removal after 3.3. Use new material laws."
29
30#include "vangenuchten.hh"
32
33#include <algorithm>
34
36
37namespace Dumux {
38
71template <class ScalarT, class ParamsT = RegularizedVanGenuchtenParams<ScalarT> >
73{
75
76public:
77 using Params = ParamsT;
78 using Scalar = typename Params::Scalar;
79
93 static Scalar pc(const Params &params, Scalar swe)
94 {
95 // retrieve the low and the high threshold saturations for the
96 // unregularized capillary pressure curve from the parameters
97 const Scalar swThLow = params.pcLowSw();
98 const Scalar swThHigh = params.pcHighSw();
99
100 // make sure that the capillary pressure observes a derivative
101 // != 0 for 'illegal' saturations. This is favourable for the
102 // newton solver (if the derivative is calculated numerically)
103 // in order to get the saturation moving to the right
104 // direction if it temporarily is in an 'illegal' range.
105 if (swe < swThLow) {
106 return VanGenuchten::pc(params, swThLow) + mLow_(params)*(swe - swThLow);
107 }
108 else if (swe > swThHigh)
109 {
110 Scalar yTh = VanGenuchten::pc(params, swThHigh);
111 Scalar m1 = (0.0 - yTh)/(1.0 - swThHigh)*2;
112
113 if (swe < 1.0) {
114 // use spline between threshold swe and 1.0
115 Scalar mTh = VanGenuchten::dpc_dswe(params, swThHigh);
116 Spline<Scalar> sp(swThHigh, 1.0, // x0, x1
117 yTh, 0, // y0, y1
118 mTh, m1); // m0, m1
119 return sp.eval(swe);
120 }
121 else {
122 // straight line for swe > 1.0
123 return m1*(swe - 1.0) + 0.0;
124 }
125 }
126
127 // if the effective saturation is in an 'reasonable'
128 // range, we use the real van genuchten law...
129 return VanGenuchten::pc(params, swe);
130 }
131
147 static Scalar sw(const Params &params, Scalar pc)
148 {
149 // retrieve the low and the high threshold saturations for the
150 // unregularized capillary pressure curve from the parameters
151 const Scalar swThLow = params.pcLowSw();
152 const Scalar swThHigh = params.pcHighSw();
153
154 // calculate the saturation which corrosponds to the
155 // saturation in the non-regularized verision of van
156 // Genuchten's law
157 Scalar sw;
158 if (pc <= 0) {
159 // for swThHigh = 1.0 the slope would get infinity
160 // swThHigh > 1.0 are not sensible threshold values
161 // setting swThHigh = 1.0 is a way to disable regularization
162 if (swThHigh > 1.0 - std::numeric_limits<Scalar>::epsilon())
163 return 1.0;
164
165 // invert straight line for swe > 1.0
166 Scalar yTh = VanGenuchten::pc(params, swThHigh);
167 Scalar m1 = (0.0 - yTh)/(1.0 - swThHigh)*2;
168 return pc/m1 + 1.0;
169 }
170 else
171 sw = VanGenuchten::sw(params, pc);
172
173 // invert the regularization if necessary
174 if (sw <= swThLow) {
175 // invert the low saturation regularization of pc()
176 Scalar pcswLow = VanGenuchten::pc(params, swThLow);
177 return (pc - pcswLow)/mLow_(params) + swThLow;
178 }
179 else if (sw > swThHigh)
180 {
181 Scalar yTh = VanGenuchten::pc(params, swThHigh);
182 Scalar m1 = (0.0 - yTh)/(1.0 - swThHigh)*2;
183
184 // invert spline between threshold swe and 1.0
185 Scalar mTh = VanGenuchten::dpc_dswe(params, swThHigh);
186 Spline<Scalar> sp(swThHigh, 1.0, // x0, x1
187 yTh, 0, // m0, m1
188 mTh, m1); // m0, m1
189 return sp.intersectInterval(swThHigh, 1.0,
190 0, 0, 0, pc);
191 }
192
193 return sw;
194 }
195
203 static Scalar endPointPc(const Params &params)
204 { return 0.0; }
205
220 static Scalar dpc_dswe(const Params &params, Scalar swe)
221 {
222 // retrieve the low and the high threshold saturations for the
223 // unregularized capillary pressure curve from the parameters
224 const Scalar swThLow = params.pcLowSw();
225 const Scalar swThHigh = params.pcHighSw();
226
227 // derivative of the regularization
228 if (swe < swThLow) {
229 // the slope of the straight line used in pc()
230 return mLow_(params);
231 }
232 else if (swe > swThHigh)
233 {
234 Scalar yTh = VanGenuchten::pc(params, swThHigh);
235 Scalar m1 = (0.0 - yTh)/(1.0 - swThHigh)*2;
236
237 if (swe < 1.0) {
238 // use spline between threshold swe and 1.0
239 Scalar mTh = VanGenuchten::dpc_dswe(params, swThHigh);
240 Spline<Scalar> sp(swThHigh, 1.0, // x0, x1
241 yTh, 0, // y0, y1
242 mTh, m1); // m0, m1
243 return sp.evalDerivative(swe);
244 }
245 else {
246 // straight line for swe > 1.0
247 return m1;
248 }
249 }
250
251 return VanGenuchten::dpc_dswe(params, swe);
252 }
253
267 static Scalar dswe_dpc(const Params &params, Scalar pc)
268 {
269 const Scalar swThLow = params.pcLowSw();
270 const Scalar swThHigh = params.pcHighSw();
271
272 if (pc <= 0)
273 {
274 // for swThHigh = 1.0 the slope gets infinity
275 // swThHigh > 1.0 are not sensible threshold values
276 // setting swThHigh = 1.0 is a way to disable regularization
277 // return the maximum representable number with the given precision
278 if (swThHigh > 1.0 - std::numeric_limits<Scalar>::epsilon())
279 return std::numeric_limits<Scalar>::max();
280
281 Scalar yTh = VanGenuchten::pc(params, swThHigh);
282 Scalar m1 = (0.0 - yTh)/(1.0 - swThHigh)*2;
283 return 1.0/m1;
284 }
285
286 const auto sw = VanGenuchten::sw(params, pc);
287
288 // derivative of the regularization
289 if (sw <= swThLow)
290 return 1/mLow_(params);
291
292 if (sw > swThHigh)
293 {
294 const Scalar yTh = VanGenuchten::pc(params, swThHigh);
295 const Scalar m1 = (0.0 - yTh)/(1.0 - swThHigh)*2;
296
297 // invert spline between threshold swe and 1.0
298 const Scalar mTh = VanGenuchten::dpc_dswe(params, swThHigh);
299 const Spline<Scalar> sp(swThHigh, 1.0, // x0, x1
300 yTh, 0, // m0, m1
301 mTh, m1); // m0, m1
302 const auto swReg = sp.intersectInterval(swThHigh, 1.0,
303 0, 0, 0, pc);
304 // derivative of the inverse of the function is one over derivative of the function
305 return 1.0/sp.evalDerivative(swReg);
306 }
307
308 return VanGenuchten::dswe_dpc(params, pc);
309 }
310
325 static Scalar krw(const Params &params, Scalar swe)
326 {
327 // retrieve the high threshold saturation for the
328 // unregularized relative permeability curve of the wetting
329 // phase from the parameters
330 const Scalar swThHigh = params.krwHighSw();
331
332 if (swe < 0)
333 return 0;
334 else if (swe > 1 - std::numeric_limits<Scalar>::epsilon())
335 return 1;
336 else if (swe > swThHigh) {
338 Spline sp(swThHigh, 1.0, // x1, x2
339 VanGenuchten::krw(params, swThHigh), 1.0, // y1, y2
340 VanGenuchten::dkrw_dswe(params, swThHigh), 0); // m1, m2
341 return sp.eval(swe);
342 }
343
344 return VanGenuchten::krw(params, swe);
345 }
346
354 static Scalar dkrw_dswe(const Params &params, Scalar swe)
355 {
356 // retrieve the high threshold saturation for the
357 // unregularized relative permeability curve of the wetting
358 // phase from the parameters
359 const Scalar swThHigh = params.krwHighSw();
360
361 // derivative of the regularization
362 // the slope is zero below sw=0.0 and above sw=1.0
363 if (swe < 0)
364 return 0.0;
365 else if (swe > 1 - std::numeric_limits<Scalar>::epsilon())
366 return 0.0;
367
368 // for high sw we need the slope of the interpolation spline
369 else if (swe > swThHigh) {
371 Spline sp(swThHigh, 1.0, // x1, x2
372 VanGenuchten::krw(params, swThHigh), 1.0, // y1, y2
373 VanGenuchten::dkrw_dswe(params, swThHigh), 0); // m1, m2
374 return sp.evalDerivative(swe);
375 }
376
377 return VanGenuchten::dkrw_dswe(params, swe);
378 }
379
394 static Scalar krn(const Params &params, Scalar swe)
395 {
396 // retrieve the low threshold saturation for the unregularized
397 // relative permeability curve of the nonwetting phase from
398 // the parameters
399 const Scalar swThLow = params.krnLowSw();
400
401 if (swe <= 0)
402 return 1;
403 else if (swe >= 1)
404 return 0;
405 else if (swe < swThLow) {
407 Spline sp(0.0, swThLow, // x1, x2
408 1.0, VanGenuchten::krn(params, swThLow), // y1, y2
409 0.0, VanGenuchten::dkrn_dswe(params, swThLow)); // m1, m2
410 return sp.eval(swe);
411 }
412
413 return VanGenuchten::krn(params, swe);
414 }
415
423 static Scalar dkrn_dswe(const Params &params, Scalar swe)
424 {
425 // retrieve the low threshold saturation for the unregularized
426 // relative permeability curve of the nonwetting phase from
427 // the parameters
428 const Scalar swThLow = params.krnLowSw();
429
430 if (swe <= 0)
431 return 0.0;
432 else if (swe >= 1)
433 return 0.0;
434 else if (swe < swThLow) {
436 Spline sp(0.0, swThLow, // x1, x2
437 1.0, VanGenuchten::krn(params, swThLow), // y1, y2
438 0.0, VanGenuchten::dkrn_dswe(params, swThLow)); // m1, m2
439 return sp.evalDerivative(swe);
440 }
441
442 return VanGenuchten::dkrn_dswe(params, swe);
443 }
444
445private:
446 // the slope of the straight line used to regularize saturations
447 // below the minimum saturation
448
457 static Scalar mLow_(const Params &params)
458 {
459 const Scalar swThLow = params.pcLowSw();
460
461 return VanGenuchten::dpc_dswe(params, swThLow);
462 }
463
472 static Scalar mHigh_(const Params &params)
473 {
474 const Scalar swThHigh = params.pcHighSw();
475
476 // for swThHigh = 1.0 the slope would get infinity
477 // swThHigh > 1.0 are not sensible threshold values
478 // setting swThHigh = 1.0 is a way to disable regularization
479 if (swThHigh > 1.0 - std::numeric_limits<Scalar>::epsilon())
480 return 0.0;
481
482 Scalar pcswHigh = VanGenuchten::pc(params, swThHigh);
483 return (0 - pcswHigh)/(1.0 - swThHigh);
484 }
485};
486
487} // end namespace Dumux
488
489#endif
Provides 3rd order polynomial splines.
Parameters that are necessary for the regularization of VanGenuchten "material law".
Implementation of the capillary pressure and relative permeability <-> saturation relations according...
Definition: adapt.hh:29
A 3rd order polynomial spline.
Definition: spline.hh:55
Scalar evalDerivative(Scalar x, bool extrapolate=false) const
Evaluate the spline's derivative at a given position.
Definition: splinecommon_.hh:175
Scalar intersectInterval(Scalar x0, Scalar x1, Scalar a, Scalar b, Scalar c, Scalar d) const
Find the intersections of the spline with a cubic polynomial in a sub-intervall of the spline,...
Definition: splinecommon_.hh:205
Scalar eval(Scalar x, bool extrapolate=false) const
Evaluate the spline at a given position.
Definition: splinecommon_.hh:141
Implementation of the regularized van Genuchten's capillary pressure / relative permeability <-> satu...
Definition: regularizedvangenuchten.hh:73
static Scalar dpc_dswe(const Params &params, Scalar swe)
A regularized version of the partial derivative of the w.r.t. effective saturation according to van ...
Definition: regularizedvangenuchten.hh:220
static Scalar dkrw_dswe(const Params &params, Scalar swe)
A regularized version of the derivative of the relative permeability for the wetting phase in regard ...
Definition: regularizedvangenuchten.hh:354
static Scalar endPointPc(const Params &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: regularizedvangenuchten.hh:203
typename Params::Scalar Scalar
Definition: regularizedvangenuchten.hh:78
static Scalar sw(const Params &params, Scalar pc)
A regularized van Genuchten saturation-capillary pressure curve.
Definition: regularizedvangenuchten.hh:147
static Scalar pc(const Params &params, Scalar swe)
A regularized van Genuchten capillary pressure-saturation curve.
Definition: regularizedvangenuchten.hh:93
static Scalar krn(const Params &params, Scalar swe)
Regularized version of the relative permeability for the nonwetting phase of the medium implied by th...
Definition: regularizedvangenuchten.hh:394
static Scalar dswe_dpc(const Params &params, Scalar pc)
A regularized version of the partial derivative of the w.r.t. cap.pressure according to van Genuchte...
Definition: regularizedvangenuchten.hh:267
ParamsT Params
Definition: regularizedvangenuchten.hh:77
static Scalar krw(const Params &params, Scalar swe)
Regularized version of the relative permeability for the wetting phase of the medium implied by the v...
Definition: regularizedvangenuchten.hh:325
static Scalar dkrn_dswe(const Params &params, Scalar swe)
A regularized version of the derivative of the relative permeability for the nonwetting phase in rega...
Definition: regularizedvangenuchten.hh:423
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
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 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
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