3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
brookscorey.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_BROOKS_COREY_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_BROOKS_COREY_HH
27
28#include <cmath>
29#include <algorithm>
30
34
36
50{
51public:
58 template<class Scalar>
59 struct Params
60 {
61 Params(Scalar pcEntry, Scalar lambda)
62 : pcEntry_(pcEntry), lambda_(lambda)
63 {}
64
65 Scalar pcEntry() const{ return pcEntry_; }
66 void setPcEntry(Scalar pcEntry){ pcEntry_ = pcEntry; }
67
68 Scalar lambda() const { return lambda_; }
69 void setLambda(Scalar lambda) { lambda_ = lambda; }
70
71 bool operator== (const Params& p) const
72 {
73 return Dune::FloatCmp::eq(pcEntry(), p.pcEntry(), 1e-6)
74 && Dune::FloatCmp::eq(lambda(), p.lambda(), 1e-6);
75 }
76
77 private:
78 Scalar pcEntry_, lambda_;
79 };
80
85 template<class Scalar = double>
86 static Params<Scalar> makeParams(const std::string& paramGroup)
87 {
88 const auto pcEntry = getParamFromGroup<Scalar>(paramGroup, "BrooksCoreyPcEntry");
89 const auto lambda = getParamFromGroup<Scalar>(paramGroup, "BrooksCoreyLambda");
90 return {pcEntry, lambda};
91 }
92
111 template<class Scalar>
112 static Scalar pc(Scalar swe, const Params<Scalar>& params)
113 {
114 using std::pow;
115 using std::clamp;
116
117 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
118
119 return params.pcEntry()*pow(swe, -1.0/params.lambda());
120 }
121
137 template<class Scalar>
138 static Scalar swe(Scalar pc, const Params<Scalar>& params)
139 {
140 using std::pow;
141 using std::max;
142
143 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
144
145 return pow(pc/params.pcEntry(), -params.lambda());
146 }
147
155 template<class Scalar>
156 static Scalar endPointPc(const Params<Scalar>& params)
157 { return params.pcEntry(); }
158
177 template<class Scalar>
178 static Scalar dpc_dswe(Scalar swe, const Params<Scalar>& params)
179 {
180 using std::pow;
181 using std::clamp;
182
183 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
184
185 return -params.pcEntry()/params.lambda() * pow(swe, -1.0/params.lambda() - 1.0);
186 }
187
201 template<class Scalar>
202 static Scalar dswe_dpc(Scalar pc, const Params<Scalar>& params)
203 {
204 using std::pow;
205 using std::max;
206
207 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
208
209 return -params.lambda()/params.pcEntry() * pow(pc/params.pcEntry(), - params.lambda() - 1.0);
210 }
211
226 template<class Scalar>
227 static Scalar krw(Scalar swe, const Params<Scalar>& params)
228 {
229 using std::pow;
230 using std::clamp;
231
232 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
233
234 return pow(swe, 2.0/params.lambda() + 3.0);
235 }
236
252 template<class Scalar>
253 static Scalar dkrw_dswe(Scalar swe, const Params<Scalar>& params)
254 {
255 using std::pow;
256 using std::clamp;
257
258 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
259
260 return (2.0/params.lambda() + 3.0)*pow(swe, 2.0/params.lambda() + 2.0);
261 }
262
277 template<class Scalar>
278 static Scalar krn(Scalar swe, const Params<Scalar>& params)
279 {
280 using std::pow;
281 using std::clamp;
282
283 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
284
285 const Scalar exponent = 2.0/params.lambda() + 1.0;
286 const Scalar sne = 1.0 - swe;
287 return sne*sne*(1.0 - pow(swe, exponent));
288 }
289
306 template<class Scalar>
307 static Scalar dkrn_dswe(Scalar swe, const Params<Scalar>& params)
308 {
309 using std::pow;
310 using std::clamp;
311
312 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
313
314 const auto lambdaInv = 1.0/params.lambda();
315 const auto swePow = pow(swe, 2*lambdaInv);
316 return 2.0*(swe - 1.0)*(1.0 + (0.5 + lambdaInv)*swePow - (1.5 + lambdaInv)*swePow*swe);
317 }
318};
319
327template <class Scalar>
329{
330public:
332 template<class S>
333 struct Params
334 {
335 Params(S pcLowSwe = 0.01) : pcLowSwe_(pcLowSwe) {}
336
337 S pcLowSwe() const { return pcLowSwe_; }
338 void setPcLowSwe(S pcLowSwe) { pcLowSwe_ = pcLowSwe; }
339
340 private:
341 S pcLowSwe_ = 0.01;
342 };
343
344 template<class MaterialLaw>
345 void init(const MaterialLaw* m, const std::string& paramGroup)
346 {
347 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup, "BrooksCoreyPcLowSweThreshold", 0.01);
348 entryPressure_ = getParamFromGroup<Scalar>(paramGroup, "BrooksCoreyPcEntry");
349
350 initPcParameters_(m, pcLowSwe_);
351 }
352
353 template<class MaterialLaw, class BaseParams, class EffToAbsParams>
354 void init(const MaterialLaw* m, const BaseParams& bp, const EffToAbsParams& etap, const Params<Scalar>& p)
355 {
356 pcLowSwe_ = p.pcLowSwe();
357 entryPressure_ = bp.pcEntry();
358
359 initPcParameters_(m, pcLowSwe_);
360 }
361
366 {
367 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
368 && Dune::FloatCmp::eq(entryPressure_, o.entryPressure_, 1e-6);
369 }
370
377 OptionalScalar<Scalar> pc(const Scalar swe) const
378 {
379 // make sure that the capillary pressure observes a derivative
380 // != 0 for 'illegal' saturations. This is favourable for the
381 // newton solver (if the derivative is calculated numerically)
382 // in order to get the saturation moving to the right
383 // direction if it temporarily is in an 'illegal' range.
384 if (swe <= pcLowSwe_)
385 return pcLowSwePcValue_ + pcDerivativeLowSw_*(swe - pcLowSwe_);
386
387 else if (swe >= 1.0)
388 return pcDerivativeHighSwEnd_*(swe - 1.0) + entryPressure_;
389
390 else
391 return {}; // no regularization
392 }
393
398 {
399 if (swe <= pcLowSwe_)
400 return pcDerivativeLowSw_;
401
402 else if (swe >= 1.0)
403 return pcDerivativeHighSwEnd_;
404
405 else
406 return {}; // no regularization
407 }
408
412 OptionalScalar<Scalar> swe(const Scalar pc) const
413 {
414 if (pc <= entryPressure_)
415 return 1.0 + (pc - entryPressure_)/pcDerivativeHighSwEnd_;
416
417 else if (pc >= pcLowSwePcValue_)
418 return (pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
419
420 else
421 return {}; // no regularization
422 }
423
428 {
429 if (pc <= entryPressure_)
430 return 1.0/pcDerivativeHighSwEnd_;
431
432 else if (pc >= pcLowSwePcValue_)
433 return 1.0/pcDerivativeLowSw_;
434
435 else
436 return {}; // no regularization
437 }
438
442 OptionalScalar<Scalar> krw(const Scalar swe) const
443 {
444 if (swe <= 0.0)
445 return 0.0;
446 else if (swe >= 1.0)
447 return 1.0;
448 else
449 return {}; // no regularization
450 }
451
456 {
457 if (swe <= 0.0)
458 return 0.0;
459 else if (swe >= 1.0)
460 return 0.0;
461 else
462 return {}; // no regularization
463 }
464
468 OptionalScalar<Scalar> krn(const Scalar swe) const
469 {
470 if (swe <= 0.0)
471 return 1.0;
472 else if (swe >= 1.0)
473 return 0.0;
474 else
475 return {}; // no regularization
476 }
477
482 {
483 if (swe <= 0.0)
484 return 0.0;
485 else if (swe >= 1.0)
486 return 0.0;
487 else
488 return {}; // no regularization
489 }
490
491private:
492 template<class MaterialLaw>
493 void initPcParameters_(const MaterialLaw* m, const Scalar lowSwe)
494 {
495 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
496 const auto highSw = MaterialLaw::EffToAbs::sweToSw(1.0, m->effToAbsParams());
497 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
498 pcDerivativeLowSw_ = m->template dpc_dsw<false>(lowSw)*dsw_dswe;
499 pcDerivativeHighSwEnd_ = m->template dpc_dsw<false>(highSw)*dsw_dswe;
500 pcLowSwePcValue_ = m->template pc<false>(lowSw);
501 }
502
503 Scalar pcLowSwe_;
504 Scalar pcLowSwePcValue_;
505 Scalar entryPressure_;
506 Scalar pcDerivativeLowSw_;
507 Scalar pcDerivativeHighSwEnd_;
508};
509
514template<typename Scalar = double>
516
521template<typename Scalar = double>
523
524} // end namespace Dumux::FluidMatrix
525
526#endif
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
Implementation of the Brooks-Corey capillary pressure <-> saturation relation. This class bundles the...
Definition: brookscorey.hh:50
static Scalar dkrw_dswe(Scalar swe, const Params< Scalar > &params)
The derivative of the relative permeability for the wetting phase with regard to the wetting saturati...
Definition: brookscorey.hh:253
static Scalar pc(Scalar swe, const Params< Scalar > &params)
The capillary pressure-saturation curve according to Brooks & Corey.
Definition: brookscorey.hh:112
static Scalar dswe_dpc(Scalar pc, const Params< Scalar > &params)
The partial derivative of the effective saturation w.r.t. the capillary pressure according to Brooks ...
Definition: brookscorey.hh:202
static Scalar swe(Scalar pc, const Params< Scalar > &params)
The saturation-capillary pressure curve according to Brooks & Corey.
Definition: brookscorey.hh:138
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 Brooks ...
Definition: brookscorey.hh:178
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: brookscorey.hh:86
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: brookscorey.hh:307
static Scalar krn(Scalar swe, const Params< Scalar > &params)
The relative permeability for the non-wetting phase of the medium as implied by the Brooks-Corey para...
Definition: brookscorey.hh:278
static Scalar endPointPc(const Params< Scalar > &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: brookscorey.hh:156
static Scalar krw(Scalar swe, const Params< Scalar > &params)
The relative permeability for the wetting phase of the medium implied by the Brooks-Corey parameteriz...
Definition: brookscorey.hh:227
The parameter type.
Definition: brookscorey.hh:60
Scalar pcEntry() const
Definition: brookscorey.hh:65
void setPcEntry(Scalar pcEntry)
Definition: brookscorey.hh:66
bool operator==(const Params &p) const
Definition: brookscorey.hh:71
Scalar lambda() const
Definition: brookscorey.hh:68
Params(Scalar pcEntry, Scalar lambda)
Definition: brookscorey.hh:61
void setLambda(Scalar lambda)
Definition: brookscorey.hh:69
A regularization for the BrooksCorey material law.
Definition: brookscorey.hh:329
OptionalScalar< Scalar > pc(const Scalar swe) const
The regularized capillary pressure-saturation curve regularized part:
Definition: brookscorey.hh:377
OptionalScalar< Scalar > krn(const Scalar swe) const
The regularized relative permeability for the non-wetting phase.
Definition: brookscorey.hh:468
OptionalScalar< Scalar > dswe_dpc(const Scalar pc) const
The regularized partial derivative of the saturation to the capillary pressure.
Definition: brookscorey.hh:427
void init(const MaterialLaw *m, const std::string &paramGroup)
Definition: brookscorey.hh:345
OptionalScalar< Scalar > dkrn_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the non-wetting phase w....
Definition: brookscorey.hh:481
OptionalScalar< Scalar > dpc_dswe(const Scalar swe) const
The regularized partial derivative of the capillary pressure w.r.t. the saturation.
Definition: brookscorey.hh:397
OptionalScalar< Scalar > krw(const Scalar swe) const
The regularized relative permeability for the wetting phase.
Definition: brookscorey.hh:442
OptionalScalar< Scalar > dkrw_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the wetting phase w.r....
Definition: brookscorey.hh:455
OptionalScalar< Scalar > swe(const Scalar pc) const
The regularized saturation-capillary pressure curve.
Definition: brookscorey.hh:412
bool operator==(const BrooksCoreyRegularization &o) const
Equality comparison with another instance.
Definition: brookscorey.hh:365
void init(const MaterialLaw *m, const BaseParams &bp, const EffToAbsParams &etap, const Params< Scalar > &p)
Definition: brookscorey.hh:354
Regularization parameters.
Definition: brookscorey.hh:334
Params(S pcLowSwe=0.01)
Definition: brookscorey.hh:335
S pcLowSwe() const
Definition: brookscorey.hh:337
void setPcLowSwe(S pcLowSwe)
Definition: brookscorey.hh:338
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