version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MATERIAL_FLUIDMATRIX_BROOKS_COREY_HH
14#define DUMUX_MATERIAL_FLUIDMATRIX_BROOKS_COREY_HH
15
16#include <cmath>
17#include <algorithm>
18
22
24
38{
39public:
46 template<class Scalar>
47 struct Params
48 {
49 Params(Scalar pcEntry, Scalar lambda)
50 : pcEntry_(pcEntry), lambda_(lambda)
51 {}
52
53 Scalar pcEntry() const{ return pcEntry_; }
54 void setPcEntry(Scalar pcEntry){ pcEntry_ = pcEntry; }
55
56 Scalar lambda() const { return lambda_; }
57 void setLambda(Scalar lambda) { lambda_ = lambda; }
58
59 bool operator== (const Params& p) const
60 {
61 return Dune::FloatCmp::eq(pcEntry(), p.pcEntry(), 1e-6)
62 && Dune::FloatCmp::eq(lambda(), p.lambda(), 1e-6);
63 }
64
65 private:
66 Scalar pcEntry_, lambda_;
67 };
68
73 template<class Scalar = double>
74 static Params<Scalar> makeParams(const std::string& paramGroup)
75 {
76 const auto pcEntry = getParamFromGroup<Scalar>(paramGroup, "BrooksCoreyPcEntry");
77 const auto lambda = getParamFromGroup<Scalar>(paramGroup, "BrooksCoreyLambda");
78 return {pcEntry, lambda};
79 }
80
99 template<class Scalar>
100 static Scalar pc(Scalar swe, const Params<Scalar>& params)
101 {
102 using std::pow;
103 using std::clamp;
104
105 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
106
107 return params.pcEntry()*pow(swe, -1.0/params.lambda());
108 }
109
125 template<class Scalar>
126 static Scalar swe(Scalar pc, const Params<Scalar>& params)
127 {
128 using std::pow;
129 using std::max;
130
131 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
132
133 return pow(pc/params.pcEntry(), -params.lambda());
134 }
135
143 template<class Scalar>
144 static Scalar endPointPc(const Params<Scalar>& params)
145 { return params.pcEntry(); }
146
165 template<class Scalar>
166 static Scalar dpc_dswe(Scalar swe, const Params<Scalar>& params)
167 {
168 using std::pow;
169 using std::clamp;
170
171 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
172
173 return -params.pcEntry()/params.lambda() * pow(swe, -1.0/params.lambda() - 1.0);
174 }
175
189 template<class Scalar>
190 static Scalar dswe_dpc(Scalar pc, const Params<Scalar>& params)
191 {
192 using std::pow;
193 using std::max;
194
195 pc = max(pc, 0.0); // the equation below is undefined for negative pcs
196
197 return -params.lambda()/params.pcEntry() * pow(pc/params.pcEntry(), - params.lambda() - 1.0);
198 }
199
214 template<class Scalar>
215 static Scalar krw(Scalar swe, const Params<Scalar>& params)
216 {
217 using std::pow;
218 using std::clamp;
219
220 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
221
222 return pow(swe, 2.0/params.lambda() + 3.0);
223 }
224
240 template<class Scalar>
241 static Scalar dkrw_dswe(Scalar swe, const Params<Scalar>& params)
242 {
243 using std::pow;
244 using std::clamp;
245
246 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
247
248 return (2.0/params.lambda() + 3.0)*pow(swe, 2.0/params.lambda() + 2.0);
249 }
250
265 template<class Scalar>
266 static Scalar krn(Scalar swe, const Params<Scalar>& params)
267 {
268 using std::pow;
269 using std::clamp;
270
271 swe = clamp(swe, 0.0, 1.0); // the equation below is only defined for 0.0 <= sw <= 1.0
272
273 const Scalar exponent = 2.0/params.lambda() + 1.0;
274 const Scalar sne = 1.0 - swe;
275 return sne*sne*(1.0 - pow(swe, exponent));
276 }
277
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 lambdaInv = 1.0/params.lambda();
303 const auto swePow = pow(swe, 2*lambdaInv);
304 return 2.0*(swe - 1.0)*(1.0 + (0.5 + lambdaInv)*swePow - (1.5 + lambdaInv)*swePow*swe);
305 }
306};
307
315template <class Scalar>
317{
318public:
320 template<class S>
321 struct Params
322 {
323 Params(S pcLowSwe = 0.01) : pcLowSwe_(pcLowSwe) {}
324
325 S pcLowSwe() const { return pcLowSwe_; }
326 void setPcLowSwe(S pcLowSwe) { pcLowSwe_ = pcLowSwe; }
327
328 private:
329 S pcLowSwe_ = 0.01;
330 };
331
332 template<class MaterialLaw>
333 void init(const MaterialLaw* m, const std::string& paramGroup)
334 {
335 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup, "BrooksCoreyPcLowSweThreshold", 0.01);
336 entryPressure_ = getParamFromGroup<Scalar>(paramGroup, "BrooksCoreyPcEntry");
337
338 initPcParameters_(m, pcLowSwe_);
339 }
340
341 template<class MaterialLaw, class BaseParams, class EffToAbsParams>
342 void init(const MaterialLaw* m, const BaseParams& bp, const EffToAbsParams& etap, const Params<Scalar>& p)
343 {
344 pcLowSwe_ = p.pcLowSwe();
345 entryPressure_ = bp.pcEntry();
346
347 initPcParameters_(m, pcLowSwe_);
348 }
349
354 {
355 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
356 && Dune::FloatCmp::eq(entryPressure_, o.entryPressure_, 1e-6);
357 }
358
365 OptionalScalar<Scalar> pc(const Scalar swe) const
366 {
367 // make sure that the capillary pressure observes a derivative
368 // != 0 for 'illegal' saturations. This is favourable for the
369 // newton solver (if the derivative is calculated numerically)
370 // in order to get the saturation moving to the right
371 // direction if it temporarily is in an 'illegal' range.
372 if (swe <= pcLowSwe_)
373 return pcLowSwePcValue_ + pcDerivativeLowSw_*(swe - pcLowSwe_);
374
375 else if (swe >= 1.0)
376 return pcDerivativeHighSwEnd_*(swe - 1.0) + entryPressure_;
377
378 else
379 return {}; // no regularization
380 }
381
386 {
387 if (swe <= pcLowSwe_)
388 return pcDerivativeLowSw_;
389
390 else if (swe >= 1.0)
391 return pcDerivativeHighSwEnd_;
392
393 else
394 return {}; // no regularization
395 }
396
400 OptionalScalar<Scalar> swe(const Scalar pc) const
401 {
402 if (pc <= entryPressure_)
403 return 1.0 + (pc - entryPressure_)/pcDerivativeHighSwEnd_;
404
405 else if (pc >= pcLowSwePcValue_)
406 return (pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
407
408 else
409 return {}; // no regularization
410 }
411
416 {
417 if (pc <= entryPressure_)
418 return 1.0/pcDerivativeHighSwEnd_;
419
420 else if (pc >= pcLowSwePcValue_)
421 return 1.0/pcDerivativeLowSw_;
422
423 else
424 return {}; // no regularization
425 }
426
430 OptionalScalar<Scalar> krw(const Scalar swe) const
431 {
432 if (swe <= 0.0)
433 return 0.0;
434 else if (swe >= 1.0)
435 return 1.0;
436 else
437 return {}; // no regularization
438 }
439
444 {
445 if (swe <= 0.0)
446 return 0.0;
447 else if (swe >= 1.0)
448 return 0.0;
449 else
450 return {}; // no regularization
451 }
452
456 OptionalScalar<Scalar> krn(const Scalar swe) const
457 {
458 if (swe <= 0.0)
459 return 1.0;
460 else if (swe >= 1.0)
461 return 0.0;
462 else
463 return {}; // no regularization
464 }
465
470 {
471 if (swe <= 0.0)
472 return 0.0;
473 else if (swe >= 1.0)
474 return 0.0;
475 else
476 return {}; // no regularization
477 }
478
479private:
480 template<class MaterialLaw>
481 void initPcParameters_(const MaterialLaw* m, const Scalar lowSwe)
482 {
483 const auto lowSw = MaterialLaw::EffToAbs::sweToSw(lowSwe, m->effToAbsParams());
484 const auto highSw = MaterialLaw::EffToAbs::sweToSw(1.0, m->effToAbsParams());
485 const auto dsw_dswe = MaterialLaw::EffToAbs::dsw_dswe(m->effToAbsParams());
486 pcDerivativeLowSw_ = m->template dpc_dsw<false>(lowSw)*dsw_dswe;
487 pcDerivativeHighSwEnd_ = m->template dpc_dsw<false>(highSw)*dsw_dswe;
488 pcLowSwePcValue_ = m->template pc<false>(lowSw);
489 }
490
491 Scalar pcLowSwe_;
492 Scalar pcLowSwePcValue_;
493 Scalar entryPressure_;
494 Scalar pcDerivativeLowSw_;
495 Scalar pcDerivativeHighSwEnd_;
496};
497
502template<typename Scalar = double>
504
509template<typename Scalar = double>
511
512} // end namespace Dumux::FluidMatrix
513
514#endif
Implementation of the Brooks-Corey capillary pressure <-> saturation relation. This class bundles the...
Definition: brookscorey.hh:38
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:241
static Scalar pc(Scalar swe, const Params< Scalar > &params)
The capillary pressure-saturation curve according to Brooks & Corey.
Definition: brookscorey.hh:100
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:190
static Scalar swe(Scalar pc, const Params< Scalar > &params)
The saturation-capillary pressure curve according to Brooks & Corey.
Definition: brookscorey.hh:126
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:166
static Params< Scalar > makeParams(const std::string &paramGroup)
Construct from a subgroup from the global parameter tree.
Definition: brookscorey.hh:74
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:295
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:266
static Scalar endPointPc(const Params< Scalar > &params)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: brookscorey.hh:144
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:215
A regularization for the BrooksCorey material law.
Definition: brookscorey.hh:317
OptionalScalar< Scalar > pc(const Scalar swe) const
The regularized capillary pressure-saturation curve regularized part:
Definition: brookscorey.hh:365
OptionalScalar< Scalar > krn(const Scalar swe) const
The regularized relative permeability for the non-wetting phase.
Definition: brookscorey.hh:456
OptionalScalar< Scalar > dswe_dpc(const Scalar pc) const
The regularized partial derivative of the saturation to the capillary pressure.
Definition: brookscorey.hh:415
void init(const MaterialLaw *m, const std::string &paramGroup)
Definition: brookscorey.hh:333
OptionalScalar< Scalar > dkrn_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the non-wetting phase w....
Definition: brookscorey.hh:469
OptionalScalar< Scalar > dpc_dswe(const Scalar swe) const
The regularized partial derivative of the capillary pressure w.r.t. the saturation.
Definition: brookscorey.hh:385
OptionalScalar< Scalar > krw(const Scalar swe) const
The regularized relative permeability for the wetting phase.
Definition: brookscorey.hh:430
OptionalScalar< Scalar > dkrw_dswe(const Scalar swe) const
The regularized derivative of the relative permeability for the wetting phase w.r....
Definition: brookscorey.hh:443
OptionalScalar< Scalar > swe(const Scalar pc) const
The regularized saturation-capillary pressure curve.
Definition: brookscorey.hh:400
bool operator==(const BrooksCoreyRegularization &o) const
Equality comparison with another instance.
Definition: brookscorey.hh:353
void init(const MaterialLaw *m, const BaseParams &bp, const EffToAbsParams &etap, const Params< Scalar > &p)
Definition: brookscorey.hh:342
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
Definition: efftoabsdefaultpolicy.hh:34
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:45
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
Definition: brookscorey.hh:23
Scalar pcEntry(const Scalar surfaceTension, const Scalar contactAngle, const Scalar inscribedRadius, const Scalar shapeFactor) noexcept
The entry capillary pressure of a pore throat.
Definition: thresholdcapillarypressures.hh:82
A wrapper that can either contain a valid Scalar or NaN.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The parameter type.
Definition: brookscorey.hh:48
Scalar pcEntry() const
Definition: brookscorey.hh:53
void setPcEntry(Scalar pcEntry)
Definition: brookscorey.hh:54
bool operator==(const Params &p) const
Definition: brookscorey.hh:59
Scalar lambda() const
Definition: brookscorey.hh:56
Params(Scalar pcEntry, Scalar lambda)
Definition: brookscorey.hh:49
void setLambda(Scalar lambda)
Definition: brookscorey.hh:57
Regularization parameters.
Definition: brookscorey.hh:322
Params(S pcLowSwe=0.01)
Definition: brookscorey.hh:323
S pcLowSwe() const
Definition: brookscorey.hh:325
void setPcLowSwe(S pcLowSwe)
Definition: brookscorey.hh:326