version 3.10-dev
heatpipelaw.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
14#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
15
16#include <cmath>
17#include <algorithm>
18
23
24namespace Dumux::FluidMatrix {
25
31template<class ScalarType, class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
32class HeatPipeLaw : Adapter<HeatPipeLaw<ScalarType, EffToAbsPolicy>, PcKrSw>
33{
34
35public:
36 using Scalar = ScalarType;
37 using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>;
38 using EffToAbs = EffToAbsPolicy;
39
43 static constexpr int numFluidPhases()
44 { return 2; }
45
49 static constexpr bool isRegularized()
50 { return false; }
51
56 struct Params
57 {
59 : gamma_(gamma), p0_(p0)
60 {}
61
62 Scalar gamma() const { return gamma_; }
63 void setGamma(Scalar gamma) { gamma_ = gamma; }
64
65 Scalar p0() const { return p0_; }
66 void setP0(Scalar p0) { p0_ = p0; }
67
68 bool operator== (const Params& p) const
69 {
70 return Dune::FloatCmp::eq(gamma(), p.gamma(), 1e-6)
71 && Dune::FloatCmp::eq(p0(), p.p0(), 1e-6);
72 }
73
74 private:
75 Scalar gamma_, p0_;
76 };
77
82 HeatPipeLaw() = delete;
83
88 explicit HeatPipeLaw(const std::string& paramGroup)
89 {
90 params_ = makeParams(paramGroup);
91 effToAbsParams_ = EffToAbs::template makeParams<Scalar>(paramGroup);
92 }
93
98 HeatPipeLaw(const Params& params,
100 : params_(params)
101 , effToAbsParams_(effToAbsParams)
102 , spline_(eps_, 1.0, // x1, x2
103 eps_*eps_*eps_, 1.0, // y1, y2
104 3.0*eps_*eps_, 0.0) // m1, m2
105 {}
106
111 static Params makeParams(const std::string& paramGroup)
112 {
113 const auto gamma = getParamFromGroup<Scalar>(paramGroup, "HeatpipeLawGamma");
114 const auto p0 = getParamFromGroup<Scalar>(paramGroup, "HeatpipeLawP0");
115 return {gamma, p0};
116 }
117
123 Scalar pc(const Scalar sw) const
124 {
125 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
126 const Scalar sne = 1 - swe;
127 const Scalar p0Gamma = params_.p0()*params_.gamma();
128
129 // regularization
130 if (sne >= 1.0)
131 {
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;
135 }
136 else if (sne <= 0.0)
137 return sne * p0Gamma*1.417;
138 else
139 return p0Gamma*((1.263*sne - 2.120)*sne + 1.417) * sne;
140 }
141
146 { return 0.0; }
147
154 Scalar dpc_dsw(const Scalar sw) const
155 {
156 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
157 const Scalar sne = 1 - swe;
158 const Scalar p0Gamma = params_.p0()*params_.gamma();
159 if (sne > 1.0)
160 sne = 1.0;
161 else if (sne <= 0.0)
162 return -p0Gamma*1.417*EffToAbs::dswe_dsw(effToAbsParams_);
163 else
164 return - p0Gamma*((3*1.263*sne - 2*2.120)*sne + 1.417)*EffToAbs::dswe_dsw(effToAbsParams_);
165 }
166
173 Scalar krw(const Scalar sw) const
174 {
175 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
176 return kr_(swe);
177 }
178
185 Scalar krn(const Scalar sw) const
186 {
187 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
188 const Scalar sne = 1 - swe; // TODO does this make sense?
189 return kr_(sne);
190 }
191
196 {
197 return params_ == o.params_
198 && effToAbsParams_ == o.effToAbsParams_;
199 }
200
202 { return effToAbsParams_; }
203
204private:
205
206 Scalar kr_(Scalar s) const
207 {
208 if (s >= 1.0)
209 return 1;
210 else if (s <= 0.0)
211 return 0;
212 else if (s > eps_)
213 return spline_.eval(s); // regularize
214 else
215 return s*s*s;
216 }
217
218 Params params_;
219 EffToAbsParams effToAbsParams_;
220 static constexpr Scalar eps_ = 0.95;
221 Dumux::Spline<Scalar> spline_;
222};
223} // end namespace Dumux::FluidMatrix
224
225#endif
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 &paramGroup)
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 &paramGroup)
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 &params, 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