25#ifndef DUMUX_REGULARIZED_LINEAR_MATERIAL_HH
26#define DUMUX_REGULARIZED_LINEAR_MATERIAL_HH
56template <
class ScalarT,
class ParamsT = RegularizedLinearMaterialParams<ScalarT> >
63 using Scalar =
typename Params::Scalar;
110 {
return params.entryPc(); }
142 return LinearMaterial::dswe_dpc(params,
pc);
155 return relperm_(params, swe);
169 return relperm_(params, sne);
175 const Scalar lowS = params.krLowS();
176 const Scalar highS = params.krHighS();
179 const Scalar m = (1 - ((1 - highS) + lowS)/2 ) / (1 - (1 - highS) - lowS);
194 else if (S > highS) {
197 1 - (1 - highS)/2, 1,
203 return lowS/2 + m*(S - lowS);
Provides 3rd order polynomial splines.
Parameters that are necessary for the regularization of the linear constitutive relations.
Linear capillary pressure and relative permeability <-> saturation relations.
A 3rd order polynomial spline.
Definition: spline.hh:55
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:48
static Scalar sw(const Params ¶ms, Scalar pc)
The saturation-capillary pressure curve.
Definition: linearmaterial.hh:86
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:116
static Scalar pc(const Params ¶ms, Scalar swe)
The linear capillary pressure-saturation curve.
Definition: linearmaterial.hh:67
Implements a linear saturation-capillary pressure relation.
Definition: regularizedlinearmaterial.hh:58
static Scalar endPointPc(const Params ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: regularizedlinearmaterial.hh:109
ParamsT Params
Definition: regularizedlinearmaterial.hh:62
static Scalar dswe_dpc(const Params ¶ms, Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: regularizedlinearmaterial.hh:140
static Scalar sw(const Params ¶ms, Scalar pc)
The saturation-capillary pressure curve.
Definition: regularizedlinearmaterial.hh:97
static Scalar dpc_dswe(const Params ¶ms, Scalar swe)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: regularizedlinearmaterial.hh:126
static Scalar pc(const Params ¶ms, Scalar swe)
The linear capillary pressure-saturation curve.
Definition: regularizedlinearmaterial.hh:78
static Scalar krn(const Params ¶ms, Scalar swe)
The relative permeability for the non-wetting phase.
Definition: regularizedlinearmaterial.hh:166
typename Params::Scalar Scalar
Definition: regularizedlinearmaterial.hh:63
static Scalar krw(const Params ¶ms, Scalar swe)
The relative permeability for the wetting phase.
Definition: regularizedlinearmaterial.hh:153