25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_LINEAR_MATERIAL_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_LINEAR_MATERIAL_HH
47template <
class ScalarT,
class ParamsT = LinearMaterialParams<ScalarT> >
52 using Scalar =
typename Params::Scalar;
70 return (1 - swe)*(params.maxPc() - params.entryPc()) + params.entryPc();
89 return 1 - (pc - params.entryPc())/(params.maxPc() - params.entryPc());
100 {
return params.entryPc(); }
119 return - (params.maxPc() - params.entryPc());
134 return - 1/(params.maxPc() - params.entryPc());
149 return clamp(swe, 0.0, 1.0);
164 return clamp(1.0-swe, 0.0, 1.0);
194 template<
class Scalar>
204 Scalar
pcMax()
const {
return pcMax_; }
210 && Dune::FloatCmp::eq(
pcMax(), p.
pcMax(), 1e-6);
214 Scalar pcEntry_, pcMax_;
221 template<
class Scalar =
double>
224 const auto pcEntry = getParamFromGroup<Scalar>(paramGroup,
"LinearPcEntry");
225 const auto pcMax = getParamFromGroup<Scalar>(paramGroup,
"LinearPcMax");
226 return {pcEntry, pcMax};
232 template<
class Scalar>
241 template<
class Scalar>
250 template<
class Scalar>
257 template<
class Scalar>
266 template<
class Scalar>
275 template<
class Scalar>
279 return clamp(
swe, 0.0, 1.0);
285 template<
class Scalar>
294 template<
class Scalar>
298 return clamp(1.0-
swe, 0.0, 1.0);
304 template<
class Scalar>
311template<
typename Scalar =
double>
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Parameters for the linear capillary pressure and relative permeability <-> saturation relations.
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:286
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:49
static Scalar krw(const Params ¶ms, Scalar swe)
The relative permeability for the wetting phase.
Definition: linearmaterial.hh:146
static Scalar sw(const Params ¶ms, Scalar pc)
The saturation-capillary pressure curve.
Definition: linearmaterial.hh:87
static Scalar dpc_dswe(const Params ¶ms, Scalar swe)
Returns the partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: linearmaterial.hh:117
static Scalar endPointPc(const Params ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: linearmaterial.hh:99
static Scalar pc(const Params ¶ms, Scalar swe)
The linear capillary pressure-saturation curve.
Definition: linearmaterial.hh:68
typename Params::Scalar Scalar
Definition: linearmaterial.hh:52
static Scalar dsw_dpc(const Params ¶ms, Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: linearmaterial.hh:132
static Scalar krn(const Params ¶ms, Scalar swe)
The relative permeability for the nonwetting phase.
Definition: linearmaterial.hh:161
ParamsT Params
Definition: linearmaterial.hh:51
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:188
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:258
static Scalar endPointPc(const Params< Scalar > ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: linearmaterial.hh:251
static Scalar dkrw_dswe(Scalar swe, const Params< Scalar > ¶ms)
The derivative of the relative permeability.
Definition: linearmaterial.hh:286
static Scalar pc(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve.
Definition: linearmaterial.hh:233
static Scalar krn(Scalar swe, const Params< Scalar > ¶ms)
The relative permeability for the non-wetting phase.
Definition: linearmaterial.hh:295
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:267
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > ¶ms)
The derivative of the relative permeability.
Definition: linearmaterial.hh:305
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: linearmaterial.hh:222
static Scalar krw(Scalar swe, const Params< Scalar > ¶ms)
The relative permeability for the wetting phase.
Definition: linearmaterial.hh:276
static Scalar swe(Scalar pc, const Params< Scalar > ¶ms)
The inverse saturation-capillary pressure curve.
Definition: linearmaterial.hh:242
The parameter type.
Definition: linearmaterial.hh:196
Scalar pcMax() const
Definition: linearmaterial.hh:204
Params(Scalar pcEntry, Scalar pcMax)
Definition: linearmaterial.hh:197
Scalar pcEntry() const
Definition: linearmaterial.hh:201
void setPcMax(Scalar pcMax)
Definition: linearmaterial.hh:205
bool operator==(const Params &p) const
Definition: linearmaterial.hh:207
void setPcEntry(Scalar pcEntry)
Definition: linearmaterial.hh:202
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition: materiallaw.hh:57