13#ifndef DUMUX_MATERIAL_FLUIDMATRIX_LINEAR_MATERIAL_HH 
   14#define DUMUX_MATERIAL_FLUIDMATRIX_LINEAR_MATERIAL_HH 
   40    template<
class Scalar>
 
   47        Scalar 
pcEntry()
 const { 
return pcEntry_; }
 
   50        Scalar 
pcMax()
 const { 
return pcMax_; }
 
   56                   && Dune::FloatCmp::eq(
pcMax(), p.
pcMax(), 1e-6);
 
 
   60        Scalar pcEntry_, pcMax_;
 
 
   67    template<
class Scalar = 
double>
 
   72        return {pcEntry, pcMax};
 
 
   78    template<
class Scalar>
 
   87    template<
class Scalar>
 
   96    template<
class Scalar>
 
  103    template<
class Scalar>
 
  112    template<
class Scalar>
 
  121    template<
class Scalar>
 
  125        return clamp(
swe, 0.0, 1.0);
 
 
  131    template<
class Scalar>
 
  140    template<
class Scalar>
 
  144        return clamp(1.0-
swe, 0.0, 1.0);
 
 
  150    template<
class Scalar>
 
 
  157template<
typename Scalar = 
double>
 
Linear capillary pressure and relative permeability <-> saturation relations.
Definition linearmaterial.hh:34
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:104
static Scalar endPointPc(const Params< Scalar > ¶ms)
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition linearmaterial.hh:97
static Scalar dkrw_dswe(Scalar swe, const Params< Scalar > ¶ms)
The derivative of the relative permeability.
Definition linearmaterial.hh:132
static Scalar pc(Scalar swe, const Params< Scalar > ¶ms)
The capillary pressure-saturation curve.
Definition linearmaterial.hh:79
static Scalar krn(Scalar swe, const Params< Scalar > ¶ms)
The relative permeability for the non-wetting phase.
Definition linearmaterial.hh:141
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:113
static Scalar dkrn_dswe(Scalar swe, const Params< Scalar > ¶ms)
The derivative of the relative permeability.
Definition linearmaterial.hh:151
static Params< Scalar > makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition linearmaterial.hh:68
static Scalar krw(Scalar swe, const Params< Scalar > ¶ms)
The relative permeability for the wetting phase.
Definition linearmaterial.hh:122
static Scalar swe(Scalar pc, const Params< Scalar > ¶ms)
The inverse saturation-capillary pressure curve.
Definition linearmaterial.hh:88
Wrapper class to implement regularized material laws (pc-sw, kr-sw) with a conversion policy between ...
Definition materiallaw.hh:45
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
Implementation helper for capillary pressure and relative permeability <-> saturation relations for t...
Definition brookscorey.hh:23
TwoPMaterialLaw< Scalar, LinearMaterial, NoRegularization, TwoPEffToAbsDefaultPolicy > LinearMaterialDefault
Definition linearmaterial.hh:158
A tag to turn off regularization and it's overhead.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The parameter type.
Definition linearmaterial.hh:42
Scalar pcMax() const
Definition linearmaterial.hh:50
Params(Scalar pcEntry, Scalar pcMax)
Definition linearmaterial.hh:43
Scalar pcEntry() const
Definition linearmaterial.hh:47
void setPcMax(Scalar pcMax)
Definition linearmaterial.hh:51
bool operator==(const Params &p) const
Definition linearmaterial.hh:53
void setPcEntry(Scalar pcEntry)
Definition linearmaterial.hh:48