25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_BROOKS_COREY_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_BROOKS_COREY_HH
58 template<
class Scalar>
65 Scalar
pcEntry()
const{
return pcEntry_; }
68 Scalar
lambda()
const {
return lambda_; }
78 Scalar pcEntry_, lambda_;
85 template<
class Scalar =
double>
88 const auto pcEntry = getParamFromGroup<Scalar>(paramGroup,
"BrooksCoreyPcEntry");
89 const auto lambda = getParamFromGroup<Scalar>(paramGroup,
"BrooksCoreyLambda");
111 template<
class Scalar>
117 swe = clamp(
swe, 0.0, 1.0);
137 template<
class Scalar>
155 template<
class Scalar>
177 template<
class Scalar>
183 swe = clamp(
swe, 0.0, 1.0);
201 template<
class Scalar>
226 template<
class Scalar>
232 swe = clamp(
swe, 0.0, 1.0);
234 return pow(
swe, 2.0/params.
lambda() + 3.0);
252 template<
class Scalar>
258 swe = clamp(
swe, 0.0, 1.0);
260 return (2.0/params.
lambda() + 3.0)*pow(
swe, 2.0/params.
lambda() + 2.0);
277 template<
class Scalar>
283 swe = clamp(
swe, 0.0, 1.0);
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));
306 template<
class Scalar>
312 swe = clamp(
swe, 0.0, 1.0);
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);
327template <
class Scalar>
344 template<
class MaterialLaw>
345 void init(
const MaterialLaw* m,
const std::string& paramGroup)
347 pcLowSwe_ = getParamFromGroup<Scalar>(paramGroup,
"BrooksCoreyPcLowSweThreshold", 0.01);
348 entryPressure_ = getParamFromGroup<Scalar>(paramGroup,
"BrooksCoreyPcEntry");
350 initPcParameters_(m, pcLowSwe_);
353 template<
class MaterialLaw,
class BaseParams,
class EffToAbsParams>
354 void init(
const MaterialLaw* m,
const BaseParams& bp,
const EffToAbsParams& etap,
const Params<Scalar>& p)
357 entryPressure_ = bp.pcEntry();
359 initPcParameters_(m, pcLowSwe_);
367 return Dune::FloatCmp::eq(pcLowSwe_, o.pcLowSwe_, 1e-6)
368 && Dune::FloatCmp::eq(entryPressure_, o.entryPressure_, 1e-6);
384 if (
swe <= pcLowSwe_)
385 return pcLowSwePcValue_ + pcDerivativeLowSw_*(
swe - pcLowSwe_);
388 return pcDerivativeHighSwEnd_*(
swe - 1.0) + entryPressure_;
399 if (
swe <= pcLowSwe_)
400 return pcDerivativeLowSw_;
403 return pcDerivativeHighSwEnd_;
414 if (
pc <= entryPressure_)
415 return 1.0 + (
pc - entryPressure_)/pcDerivativeHighSwEnd_;
417 else if (
pc >= pcLowSwePcValue_)
418 return (
pc - pcLowSwePcValue_)/pcDerivativeLowSw_ + pcLowSwe_;
429 if (
pc <= entryPressure_)
430 return 1.0/pcDerivativeHighSwEnd_;
432 else if (
pc >= pcLowSwePcValue_)
433 return 1.0/pcDerivativeLowSw_;
492 template<
class MaterialLaw>
493 void initPcParameters_(
const MaterialLaw* m,
const Scalar lowSwe)
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);
504 Scalar pcLowSwePcValue_;
505 Scalar entryPressure_;
506 Scalar pcDerivativeLowSw_;
507 Scalar pcDerivativeHighSwEnd_;
514template<
typename Scalar =
double>
521template<
typename Scalar =
double>
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
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:94
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 > ¶ms)
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 > ¶ms)
The capillary pressure-saturation curve according to Brooks & Corey.
Definition: brookscorey.hh:112
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:202
static Scalar swe(Scalar pc, const Params< Scalar > ¶ms)
The saturation-capillary pressure curve according to Brooks & Corey.
Definition: brookscorey.hh:138
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:178
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: brookscorey.hh:86
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:307
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:278
static Scalar endPointPc(const Params< Scalar > ¶ms)
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 > ¶ms)
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 ¶mGroup)
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