25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
43template<
class ScalarType,
class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
82 return Dune::FloatCmp::eq(
gamma(), p.
gamma(), 1e-6)
83 && Dune::FloatCmp::eq(
p0(), p.
p0(), 1e-6);
103 effToAbsParams_ = EffToAbs::template makeParams<Scalar>(paramGroup);
125 const auto gamma = getParamFromGroup<Scalar>(paramGroup,
"HeatpipeLawGamma");
126 const auto p0 = getParamFromGroup<Scalar>(paramGroup,
"HeatpipeLawP0");
137 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
138 const Scalar sne = 1 - swe;
144 const Scalar y = p0Gamma*( (1.263*1.0 - 2.120)*1.0 + 1.417)*1.0;
145 const Scalar m = p0Gamma*((3*1.263*1.0 - 2*2.120)*1.0 + 1.417);
146 return (sne - 1)*m + y;
149 return sne * p0Gamma*1.417;
151 return p0Gamma*((1.263*sne - 2.120)*sne + 1.417) * sne;
168 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
169 const Scalar sne = 1 - swe;
174 return -p0Gamma*1.417*EffToAbs::dswe_dsw(effToAbsParams_);
176 return - p0Gamma*((3*1.263*sne - 2*2.120)*sne + 1.417)*EffToAbs::dswe_dsw(effToAbsParams_);
187 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
199 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
200 const Scalar sne = 1 - swe;
209 return params_ == o.params_
210 && effToAbsParams_ == o.effToAbsParams_;
214 {
return effToAbsParams_; }
225 return spline_.eval(s);
232 static constexpr Scalar eps_ = 0.95;
Provides 3rd order polynomial splines.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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:35
A 3rd order polynomial spline.
Definition: spline.hh:55
Implementation of the capillary pressure <-> saturation relation for the heatpipe problem.
Definition: heatpipelaw.hh:45
Scalar pc(const Scalar sw) const
The capillary pressure-saturation curve.
Definition: heatpipelaw.hh:135
Scalar krn(const Scalar sw) const
The relative permeability for the non-wetting phase of the medium.
Definition: heatpipelaw.hh:197
static Params makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: heatpipelaw.hh:123
ScalarType Scalar
Definition: heatpipelaw.hh:48
EffToAbsPolicy EffToAbs
Definition: heatpipelaw.hh:50
HeatPipeLaw(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: heatpipelaw.hh:100
Scalar krw(const Scalar sw) const
The relative permeability for the wetting phase of the medium.
Definition: heatpipelaw.hh:185
Scalar dpc_dsw(const Scalar sw) const
The partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: heatpipelaw.hh:166
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: heatpipelaw.hh:49
HeatPipeLaw()=delete
Deleted default constructor (so we are never in an undefined state)
const EffToAbsParams & effToAbsParams() const
Definition: heatpipelaw.hh:213
bool operator==(const HeatPipeLaw< Scalar, EffToAbs > &o) const
Equality comparison with another instance.
Definition: heatpipelaw.hh:207
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: heatpipelaw.hh:55
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: heatpipelaw.hh:61
Scalar endPointPc() const
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: heatpipelaw.hh:157
HeatPipeLaw(const Params ¶ms, const EffToAbsParams &effToAbsParams={})
Construct from parameter structs.
Definition: heatpipelaw.hh:110
void setP0(Scalar p0)
Definition: heatpipelaw.hh:78
Scalar gamma() const
Definition: heatpipelaw.hh:74
void setGamma(Scalar gamma)
Definition: heatpipelaw.hh:75
bool operator==(const Params &p) const
Definition: heatpipelaw.hh:80
Scalar p0() const
Definition: heatpipelaw.hh:77
Params(Scalar gamma, Scalar p0)
Definition: heatpipelaw.hh:70
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67