13#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
14#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
31template<
class ScalarType,
class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
70 return Dune::FloatCmp::eq(
gamma(), p.
gamma(), 1e-6)
71 && Dune::FloatCmp::eq(
p0(), p.
p0(), 1e-6);
91 effToAbsParams_ = EffToAbs::template makeParams<Scalar>(paramGroup);
113 const auto gamma = getParamFromGroup<Scalar>(paramGroup,
"HeatpipeLawGamma");
114 const auto p0 = getParamFromGroup<Scalar>(paramGroup,
"HeatpipeLawP0");
125 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
126 const Scalar sne = 1 - swe;
132 const Scalar y = p0Gamma*( (1.263*1.0 - 2.120)*1.0 + 1.417)*1.0;
133 const Scalar m = p0Gamma*((3*1.263*1.0 - 2*2.120)*1.0 + 1.417);
134 return (sne - 1)*m + y;
137 return sne * p0Gamma*1.417;
139 return p0Gamma*((1.263*sne - 2.120)*sne + 1.417) * sne;
156 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
157 const Scalar sne = 1 - swe;
162 return -p0Gamma*1.417*EffToAbs::dswe_dsw(effToAbsParams_);
164 return - p0Gamma*((3*1.263*sne - 2*2.120)*sne + 1.417)*EffToAbs::dswe_dsw(effToAbsParams_);
175 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
187 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
188 const Scalar sne = 1 - swe;
197 return params_ == o.params_
198 && effToAbsParams_ == o.effToAbsParams_;
202 {
return effToAbsParams_; }
213 return spline_.eval(s);
220 static constexpr Scalar eps_ = 0.95;
Implementation of the capillary pressure <-> saturation relation for the heatpipe problem.
Definition: heatpipelaw.hh:33
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: heatpipelaw.hh:123
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase of the medium.
Definition: heatpipelaw.hh:185
static Params makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: heatpipelaw.hh:111
ScalarType Scalar
Definition: heatpipelaw.hh:36
EffToAbsPolicy EffToAbs
Definition: heatpipelaw.hh:38
HeatPipeLaw(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: heatpipelaw.hh:88
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase of the medium.
Definition: heatpipelaw.hh:173
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: heatpipelaw.hh:154
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: heatpipelaw.hh:37
HeatPipeLaw()=delete
Deleted default constructor (so we are never in an undefined state)
const EffToAbsParams & effToAbsParams() const
Definition: heatpipelaw.hh:201
bool operator==(const HeatPipeLaw< Scalar, EffToAbs > &o) const
Equality comparison with another instance.
Definition: heatpipelaw.hh:195
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: heatpipelaw.hh:43
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: heatpipelaw.hh:49
Scalar endPointPc() const
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: heatpipelaw.hh:145
HeatPipeLaw(const Params ¶ms, const EffToAbsParams &effToAbsParams={})
Construct from parameter structs.
Definition: heatpipelaw.hh:98
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
void setP0(Scalar p0)
Definition: heatpipelaw.hh:66
Scalar gamma() const
Definition: heatpipelaw.hh:62
void setGamma(Scalar gamma)
Definition: heatpipelaw.hh:63
bool operator==(const Params &p) const
Definition: heatpipelaw.hh:68
Scalar p0() const
Definition: heatpipelaw.hh:65
Params(Scalar gamma, Scalar p0)
Definition: heatpipelaw.hh:58