25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
48template <
class ScalarT,
class ParamsT = HeatPipeLawParams<ScalarT> >
53 using Scalar =
typename Params::Scalar;
64 Scalar p0Gamma = params.p0()*params.gamma();
68 Scalar y = p0Gamma*( (1.263*1.0 - 2.120)*1.0 + 1.417)*1.0;
69 Scalar m = p0Gamma*((3*1.263*1.0 - 2*2.120)*1.0 + 1.417);
70 return (Sn - 1)*m + y;
78 return p0Gamma*((1.263*Sn - 2.120)*Sn + 1.417) * Sn;
90 DUNE_THROW(Dune::NotImplemented,
"HeatPipeLaw::Sw");
102 Scalar p0Gamma = params.p0()*params.gamma();
105 else if (Sn <= 0.0) {
106 Scalar m = -p0Gamma*1.417;
110 Scalar m = - p0Gamma*((3*1.263*Sn - 2*2.120)*Sn + 1.417);
122 DUNE_THROW(Dune::NotImplemented,
"HeatPipeLaw::dSw_dpC");
149 static Scalar kr_(Scalar S)
151 const Scalar eps = 0.95;
187template<
class ScalarType,
class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
226 return Dune::FloatCmp::eq(
gamma(), p.
gamma(), 1e-6)
227 && Dune::FloatCmp::eq(
p0(), p.
p0(), 1e-6);
247 effToAbsParams_ = EffToAbs::template makeParams<Scalar>(paramGroup);
269 const auto gamma = getParamFromGroup<Scalar>(paramGroup,
"HeatpipeLawGamma");
270 const auto p0 = getParamFromGroup<Scalar>(paramGroup,
"HeatpipeLawP0");
279 template<
class Scalar>
282 const Scalar sne = 1 - swe;
288 const Scalar y = p0Gamma*( (1.263*1.0 - 2.120)*1.0 + 1.417)*1.0;
289 const Scalar m = p0Gamma*((3*1.263*1.0 - 2*2.120)*1.0 + 1.417);
290 return (sne - 1)*m + y;
293 return sne * p0Gamma*1.417;
295 return p0Gamma*((1.263*sne - 2.120)*sne + 1.417) * sne;
301 template<
class Scalar>
311 template<
class Scalar>
314 const Scalar sne = 1 - swe;
319 return -p0Gamma*1.417;
321 return - p0Gamma*((3*1.263*sne - 2*2.120)*sne + 1.417);
330 template<
class Scalar>
342 template<
class Scalar>
345 const Scalar sne = 1 - swe;
354 return params_ == o.params_
355 && effToAbsParams_ == o.effToAbsParams_;
359 {
return effToAbsParams_; }
370 return spline_.eval(s);
377 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.
Specification of the material params for the heat pipe's capillary pressure model.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
Definition: brookscorey.hh:286
A 3rd order polynomial spline.
Definition: spline.hh:55
Implementation of the capillary pressure <-> saturation relation for the heatpipe problem.
Definition: heatpipelaw.hh:50
static Scalar dpC_dSw(const Params ¶ms, Scalar Sw)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: heatpipelaw.hh:99
static Scalar krn(const Params ¶ms, Scalar Sw)
The relative permeability for the nonwetting phase.
Definition: heatpipelaw.hh:142
typename Params::Scalar Scalar
Definition: heatpipelaw.hh:53
static Scalar krw(const Params ¶ms, Scalar Sw)
The relative permeability for the wetting phase.
Definition: heatpipelaw.hh:131
ParamsT Params
Definition: heatpipelaw.hh:52
static Scalar pc(const Params ¶ms, Scalar Sw)
The capillary pressure-saturation curve.
Definition: heatpipelaw.hh:61
static Scalar Sw(const Params ¶ms, Scalar pC)
The saturation-capillary pressure curve.
Definition: heatpipelaw.hh:88
static Scalar dSw_dpC(const Params ¶ms, Scalar pC)
Returns the partial derivative of the effective saturation to the capillary pressure.
Definition: heatpipelaw.hh:120
Implementation of the capillary pressure <-> saturation relation for the heatpipe problem.
Definition: heatpipelaw.hh:189
static Params makeParams(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: heatpipelaw.hh:267
ScalarType Scalar
Definition: heatpipelaw.hh:192
Scalar pc(Scalar swe) const
The capillary pressure-saturation curve.
Definition: heatpipelaw.hh:280
EffToAbsPolicy EffToAbs
Definition: heatpipelaw.hh:194
HeatPipeLaw(const std::string ¶mGroup)
Construct from a subgroup from the global parameter tree.
Definition: heatpipelaw.hh:244
Scalar krw(Scalar swe) const
The relative permeability for the wetting phase of the medium.
Definition: heatpipelaw.hh:331
Scalar dpc_dswe(Scalar swe) const
The partial derivative of the capillary pressure w.r.t. the effective saturation.
Definition: heatpipelaw.hh:312
Scalar endPointPc() const
The capillary pressure at Swe = 1.0 also called end point capillary pressure.
Definition: heatpipelaw.hh:302
typename EffToAbsPolicy::template Params< Scalar > EffToAbsParams
Definition: heatpipelaw.hh:193
HeatPipeLaw()=delete
Deleted default constructor (so we are never in an undefined state)
const EffToAbsParams & effToAbsParams() const
Definition: heatpipelaw.hh:358
Scalar krn(Scalar swe) const
The relative permeability for the non-wetting phase of the medium.
Definition: heatpipelaw.hh:343
bool operator==(const HeatPipeLaw< Scalar, EffToAbs > &o) const
Equality comparison with another instance.
Definition: heatpipelaw.hh:352
static constexpr int numFluidPhases()
Return the number of fluid phases.
Definition: heatpipelaw.hh:199
static constexpr bool isRegularized()
Return whether this law is regularized.
Definition: heatpipelaw.hh:205
HeatPipeLaw(const Params ¶ms, const EffToAbsParams &effToAbsParams={})
Construct from parameter structs.
Definition: heatpipelaw.hh:254
void setP0(Scalar p0)
Definition: heatpipelaw.hh:222
Scalar gamma() const
Definition: heatpipelaw.hh:218
void setGamma(Scalar gamma)
Definition: heatpipelaw.hh:219
bool operator==(const Params &p) const
Definition: heatpipelaw.hh:224
Scalar p0() const
Definition: heatpipelaw.hh:221
Params(Scalar gamma, Scalar p0)
Definition: heatpipelaw.hh:214
Adapter to inherit from, allowing the inheriting class to be wrapped by the makeFluidMatrixInteractio...
Definition: fluidmatrixinteraction.hh:67