25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_LINEAR_MATERIAL_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_LINEAR_MATERIAL_HH
52 template<
class Scalar>
59 Scalar
pcEntry()
const {
return pcEntry_; }
62 Scalar
pcMax()
const {
return pcMax_; }
68 && Dune::FloatCmp::eq(
pcMax(), p.
pcMax(), 1e-6);
72 Scalar pcEntry_, pcMax_;
79 template<
class Scalar =
double>
82 const auto pcEntry = getParamFromGroup<Scalar>(paramGroup,
"LinearPcEntry");
83 const auto pcMax = getParamFromGroup<Scalar>(paramGroup,
"LinearPcMax");
90 template<
class Scalar>
99 template<
class Scalar>
108 template<
class Scalar>
115 template<
class Scalar>
124 template<
class Scalar>
133 template<
class Scalar>
137 return clamp(
swe, 0.0, 1.0);
143 template<
class Scalar>
152 template<
class Scalar>
156 return clamp(1.0-
swe, 0.0, 1.0);
162 template<
class Scalar>
169template<
typename Scalar =
double>
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
A tag to turn off regularization and it's overhead.
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
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:46
static Scalar dpc_dswe(Scalar swe, const Params< Scalar > ¶ms)
The partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: linearmaterial.hh:116
static Scalar endPointPc(const Params< Scalar > ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: linearmaterial.hh:109
static Scalar dkrw_dswe(Scalar swe, const Params< Scalar > ¶ms)
The derivative of the relative permeability.
Definition: linearmaterial.hh:144
static Scalar pc(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve.
Definition: linearmaterial.hh:91
static Scalar krn(Scalar swe, const Params< Scalar > ¶ms)
The relative permeability for the non-wetting phase.
Definition: linearmaterial.hh:153
static Scalar dswe_dpc(Scalar pc, const Params< Scalar > ¶ms)
The partial derivative of the effective saturation w.r.t. the capillary pressure.
Definition: linearmaterial.hh:125
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > ¶ms)
The derivative of the relative permeability.
Definition: linearmaterial.hh:163
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: linearmaterial.hh:80
static Scalar krw(Scalar swe, const Params< Scalar > ¶ms)
The relative permeability for the wetting phase.
Definition: linearmaterial.hh:134
static Scalar swe(Scalar pc, const Params< Scalar > ¶ms)
The inverse saturation-capillary pressure curve.
Definition: linearmaterial.hh:100
The parameter type.
Definition: linearmaterial.hh:54
Scalar pcMax() const
Definition: linearmaterial.hh:62
Params(Scalar pcEntry, Scalar pcMax)
Definition: linearmaterial.hh:55
Scalar pcEntry() const
Definition: linearmaterial.hh:59
void setPcMax(Scalar pcMax)
Definition: linearmaterial.hh:63
bool operator==(const Params &p) const
Definition: linearmaterial.hh:65
void setPcEntry(Scalar pcEntry)
Definition: linearmaterial.hh:60
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:57