13#ifndef DUMUX_MATERIAL_FLUIDMATRIX_BROOKS_COREY_HH
14#define DUMUX_MATERIAL_FLUIDMATRIX_BROOKS_COREY_HH
46 template<
class Scalar>
53 Scalar
pcEntry()
const{
return pcEntry_; }
56 Scalar
lambda()
const {
return lambda_; }
66 Scalar pcEntry_, lambda_;
73 template<
class Scalar =
double>
76 const auto pcEntry = getParamFromGroup<Scalar>(paramGroup,
"BrooksCoreyPcEntry");
77 const auto lambda = getParamFromGroup<Scalar>(paramGroup,
"BrooksCoreyLambda");
99 template<
class Scalar>
105 swe = clamp(
swe, 0.0, 1.0);
125 template<
class Scalar>
143 template<
class Scalar>
165 template<
class Scalar>
171 swe = clamp(
swe, 0.0, 1.0);
189 template<
class Scalar>
214 template<
class Scalar>
220 swe = clamp(
swe, 0.0, 1.0);
222 return pow(
swe, 2.0/params.
lambda() + 3.0);
240 template<
class Scalar>
246 swe = clamp(
swe, 0.0, 1.0);
248 return (2.0/params.
lambda() + 3.0)*pow(
swe, 2.0/params.
lambda() + 2.0);
265 template<
class Scalar>
271 swe = clamp(
swe, 0.0, 1.0);
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));
294 template<
class Scalar>
300 swe = clamp(
swe, 0.0, 1.0);
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);
315template <
class Scalar>
332 template<
class MaterialLaw>
333 void init(
const MaterialLaw* m,
const std::string& paramGroup)
335 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"BrooksCoreyPcLowSweThreshold", 0.01);
336 entryPressure_ = getParamFromGroup<Scalar>(paramGroup,
"BrooksCoreyPcEntry");
338 initPcParameters_(m, pcLowSwe_);
341 template<
class MaterialLaw,
class BaseParams,
class EffToAbsParams>
342 void init(
const MaterialLaw* m,
const BaseParams& bp,
const EffToAbsParams& etap,
const Params<Scalar>& p)
345 entryPressure_ = bp.pcEntry();
347 initPcParameters_(m, pcLowSwe_);
355 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
356 && Dune::FloatCmp::eq(entryPressure_, o.entryPressure_, 1e-6);
372 if (
swe <= pcLowSwe_)
373 return pcLowSwePcValue_ + pcDerivativeLowSw_*(
swe - pcLowSwe_);
376 return pcDerivativeHighSwEnd_*(
swe - 1.0) + entryPressure_;
387 if (
swe <= pcLowSwe_)
388 return pcDerivativeLowSw_;
391 return pcDerivativeHighSwEnd_;
402 if (
pc <= entryPressure_)
403 return 1.0 + (
pc - entryPressure_)/pcDerivativeHighSwEnd_;
405 else if (
pc >= pcLowSwePcValue_)
406 return (
pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
417 if (
pc <= entryPressure_)
418 return 1.0/pcDerivativeHighSwEnd_;
420 else if (
pc >= pcLowSwePcValue_)
421 return 1.0/pcDerivativeLowSw_;
480 template<
class MaterialLaw>
481 void initPcParameters_(
const MaterialLaw* m,
const Scalar lowSwe)
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);
492 Scalar pcLowSwePcValue_;
493 Scalar entryPressure_;
494 Scalar pcDerivativeLowSw_;
495 Scalar pcDerivativeHighSwEnd_;
502template<
typename Scalar =
double>
509template<
typename Scalar =
double>
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 > ¶ms)
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 > ¶ms)
The capillary pressure-saturation curve according to Brooks & Corey.
Definition: brookscorey.hh:100
static Scalar dswe_dpc(Scalar pc, const Params< Scalar > ¶ms)
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 > ¶ms)
The saturation-capillary pressure curve according to Brooks & Corey.
Definition: brookscorey.hh:126
static Scalar dpc_dswe(Scalar swe, const Params< Scalar > ¶ms)
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 ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: brookscorey.hh:74
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > ¶ms)
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 > ¶ms)
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 > ¶ms)
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 > ¶ms)
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 ¶mGroup)
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