25#ifndef LINEAR_MATERIAL_HH
26#define LINEAR_MATERIAL_HH
46template <
class ScalarT,
class ParamsT = LinearMaterialParams<ScalarT> >
51 using Scalar =
typename Params::Scalar;
69 return (1 - swe)*(params.maxPc() - params.entryPc()) + params.entryPc();
88 return 1 - (
pc - params.entryPc())/(params.maxPc() - params.entryPc());
99 {
return params.entryPc(); }
118 return - (params.maxPc() - params.entryPc());
133 return - 1/(params.maxPc() - params.entryPc());
149 return max(min(swe,1.0),0.0);
166 return max(min(sne,1.0),0.0);
Parameters for the linear capillary pressure and relative permeability <-> saturation relations.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Linear capillary pressure and relative permeability <-> saturation relations.
Definition: linearmaterial.hh:48
static Scalar krw(const Params ¶ms, Scalar swe)
The relative permeability for the wetting phase.
Definition: linearmaterial.hh:145
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 endPointPc(const Params ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: linearmaterial.hh:98
static Scalar pc(const Params ¶ms, Scalar swe)
The linear capillary pressure-saturation curve.
Definition: linearmaterial.hh:67
typename Params::Scalar Scalar
Definition: linearmaterial.hh:51
static Scalar dsw_dpc(const Params ¶ms, Scalar pc)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: linearmaterial.hh:131
static Scalar krn(const Params ¶ms, Scalar swe)
The relative permeability for the non-wetting phase.
Definition: linearmaterial.hh:161
ParamsT Params
Definition: linearmaterial.hh:50