13#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_SMOOTHED_LINEAR_LAW_HH
14#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_SMOOTHED_LINEAR_LAW_HH
42template<
class ScalarType,
class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
87 return Dune::FloatCmp::eq(
pe(), p.
pe(), 1e-6)
88 && Dune::FloatCmp::eq(
pcMax(), p.
pcMax(), 1e-6)
94 Scalar pe_, pcMax_, krLowS_, krHighS_;
119 , splineM_((1.0 - ((1.0 - params_.krHighS()) + params_.krLowS())/2.0 )
120 / (1.0 - (1.0 - params_.krHighS()) - params_.krLowS()))
121 , splineLowS_(0.0, params_.krLowS(),
122 0.0, params_.krLowS()/2.0,
124 , splineHighS_(params_.krHighS(), 1.0,
125 1.0 - (1.0 - params_.krHighS())/2.0, 1.0,
135 const auto pe = getParamFromGroup<Scalar>(paramGroup,
"SmoothedLinearLawPe");
136 const auto pcMax = getParamFromGroup<Scalar>(paramGroup,
"SmoothedLinearLawPcMax");
137 const auto krLowS = getParamFromGroup<Scalar>(paramGroup,
"SmoothedLinearLawKrLowS");
138 const auto krHighS = getParamFromGroup<Scalar>(paramGroup,
"SmoothedLinearLawKrHighS");
139 return {pe, pcMax, krLowS, krHighS};
147 return (1.0 -
swe)*(params_.
pcMax() - params_.
pe()) + params_.
pe();
155 return 1.0 - (
pc - params_.
pe())/(params_.
pcMax() - params_.
pe());
162 {
return params_.
pe(); }
169 return params_.
pe() - params_.
pcMax();
177 return 1.0/(params_.
pe() - params_.
pcMax());
185 return relperm_(
swe);
194 return relperm_(sne);
202 return params_ == o.params_
203 && effToAbsParams_ == o.effToAbsParams_;
207 {
return effToAbsParams_; }
223 return splineLowS_.eval(S);
225 return splineHighS_.eval(S);
228 return lowS/2.0 + splineM_*(S - lowS);
Implements a linear saturation-capillary pressure relation.
Definition: smoothedlinearlaw.hh:44
Scalar krw(Scalar swe) const
The relative permeability for the wetting phase.
Definition: smoothedlinearlaw.hh:183
EffToAbsPolicy EffToAbs
Definition: smoothedlinearlaw.hh:49
Scalar endPointPc() const
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: smoothedlinearlaw.hh:161
SmoothedLinearLaw(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: smoothedlinearlaw.hh:107
SmoothedLinearLaw()=delete
Deleted default constructor (so we are never in an undefined state)
static Params makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: smoothedlinearlaw.hh:133
SmoothedLinearLaw(const Params ¶ms, const EffToAbsParams &effToAbsParams={})
Construct from parameter structs.
Definition: smoothedlinearlaw.hh:115
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: smoothedlinearlaw.hh:60
Scalar swe(Scalar pc) const
The inverse saturation-capillary pressure curve.
Definition: smoothedlinearlaw.hh:153
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: smoothedlinearlaw.hh:48
const EffToAbsParams & effToAbsParams() const
Definition: smoothedlinearlaw.hh:206
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: smoothedlinearlaw.hh:54
bool operator==(const SmoothedLinearLaw< Scalar, EffToAbs > &o) const
Equality comparison with another instance.
Definition: smoothedlinearlaw.hh:200
Scalar pc(Scalar swe) const
The capillary pressure-saturation curve.
Definition: smoothedlinearlaw.hh:145
ScalarType Scalar
Definition: smoothedlinearlaw.hh:47
Scalar dswe_dpc(Scalar pc) const
The partial derivative of the effective saturation w.r.t. the capillary pressure.
Definition: smoothedlinearlaw.hh:175
Scalar dpc_dswe(Scalar swe) const
The partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: smoothedlinearlaw.hh:167
Scalar krn(Scalar swe) const
The relative permeability for the non-wetting phase.
Definition: smoothedlinearlaw.hh:191
A 3rd order polynomial spline.
Definition: spline.hh:43
This is a policy for 2p material laws how to convert absolute to relative saturations and vice versa.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: brookscorey.hh:23
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Provides 3rd order polynomial splines.
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:55
The parameter type.
Definition: smoothedlinearlaw.hh:68
void setKrLowS(Scalar krLowS)
Definition: smoothedlinearlaw.hh:80
void setPe(Scalar pe)
Definition: smoothedlinearlaw.hh:74
Scalar krHighS() const
Definition: smoothedlinearlaw.hh:82
Scalar pcMax() const
Definition: smoothedlinearlaw.hh:76
Scalar krLowS() const
Definition: smoothedlinearlaw.hh:79
Scalar pe() const
Definition: smoothedlinearlaw.hh:73
Params(Scalar pe, Scalar pcMax, Scalar krLowS, Scalar krHighS)
Definition: smoothedlinearlaw.hh:69
void setPcMax(Scalar pcMax)
Definition: smoothedlinearlaw.hh:77
void setKrHighS(Scalar krHighS)
Definition: smoothedlinearlaw.hh:83
bool operator==(const Params &p) const
Definition: smoothedlinearlaw.hh:85