3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
26#define DUMUX_MATERIAL_FLUIDMATRIX_TWOP_HEATPIPELAW_HH
27
28// remove from here after release 3.3 /////////////
29#include "heatpipelawparams.hh"
30
32
33#include <algorithm>
34
35#include <math.h>
36#include <assert.h>
37
38namespace Dumux {
39
48template <class ScalarT, class ParamsT = HeatPipeLawParams<ScalarT> >
49class [[deprecated("Use new material laws and FluidMatrix::HeatPipeLaw instead!")]] HeatPipeLaw
50{
51public:
52 using Params = ParamsT;
53 using Scalar = typename Params::Scalar;
54
61 static Scalar pc(const Params &params, Scalar Sw)
62 {
63 Scalar Sn = 1 - Sw;
64 Scalar p0Gamma = params.p0()*params.gamma();
65
66 // regularization
67 if (Sn >= 1.0) {
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;
71 }
72 else if (Sn <= 0.0) {
73 Scalar y = 0.0;
74 Scalar m = p0Gamma*1.417;
75 return Sn*m + y;
76 }
77
78 return p0Gamma*((1.263*Sn - 2.120)*Sn + 1.417) * Sn;
79 }
80
88 static Scalar Sw(const Params &params, Scalar pC)
89 {
90 DUNE_THROW(Dune::NotImplemented, "HeatPipeLaw::Sw");
91 }
92
99 static Scalar dpC_dSw(const Params &params, Scalar Sw)
100 {
101 Scalar Sn = 1 - Sw;
102 Scalar p0Gamma = params.p0()*params.gamma();
103 if (Sn > 1.0)
104 Sn = 1.0;
105 else if (Sn <= 0.0) {
106 Scalar m = -p0Gamma*1.417;
107 return m;
108 }
109
110 Scalar m = - p0Gamma*((3*1.263*Sn - 2*2.120)*Sn + 1.417);
111 return m;
112 }
113
120 static Scalar dSw_dpC(const Params &params, Scalar pC)
121 {
122 DUNE_THROW(Dune::NotImplemented, "HeatPipeLaw::dSw_dpC");
123 }
124
131 static Scalar krw(const Params &params, Scalar Sw)
132 {
133 return kr_(Sw);
134 }
135
142 static Scalar krn(const Params &params, Scalar Sw)
143 {
144 Scalar Sn = 1 - Sw;
145 return kr_(Sn);
146 }
147
148private:
149 static Scalar kr_(Scalar S)
150 {
151 const Scalar eps = 0.95;
152 if (S >= 1)
153 return 1;
154 else if (S <= 0)
155 return 0;
156 else if (S > eps) {
157 // regularize
159 Spline sp(eps, 1.0, // x1, x2
160 eps*eps*eps, 1, // y1, y2
161 3*eps*eps, 0); // m1, m2
162 return sp.eval(S);
163 }
164
165 return S*S*S;
166 }
167};
168
169} // end namespace Dumux
170// remove until here after release 3.3 /////////////
171
172#include <cmath>
173#include <algorithm>
174
176#include <dumux/common/spline.hh>
179
180namespace Dumux::FluidMatrix {
181
187template<class ScalarType, class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
188class HeatPipeLaw : Adapter<HeatPipeLaw<ScalarType, EffToAbsPolicy>, PcKrSw>
189{
190
191public:
192 using Scalar = ScalarType;
193 using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>;
194 using EffToAbs = EffToAbsPolicy;
195
199 static constexpr int numFluidPhases()
200 { return 2; }
201
205 static constexpr bool isRegularized()
206 { return false; }
207
212 struct Params
213 {
215 : gamma_(gamma), p0_(p0)
216 {}
217
218 Scalar gamma() const { return gamma_; }
219 void setGamma(Scalar gamma) { gamma_ = gamma; }
220
221 Scalar p0() const { return p0_; }
222 void setP0(Scalar p0) { p0_ = p0; }
223
224 bool operator== (const Params& p) const
225 {
226 return Dune::FloatCmp::eq(gamma(), p.gamma(), 1e-6)
227 && Dune::FloatCmp::eq(p0(), p.p0(), 1e-6);
228 }
229
230 private:
231 Scalar gamma_, p0_;
232 };
233
238 HeatPipeLaw() = delete;
239
244 explicit HeatPipeLaw(const std::string& paramGroup)
245 {
246 params_ = makeParams(paramGroup);
247 effToAbsParams_ = EffToAbs::template makeParams<Scalar>(paramGroup);
248 }
249
254 HeatPipeLaw(const Params& params,
255 const EffToAbsParams& effToAbsParams = {})
256 : params_(params)
257 , effToAbsParams_(effToAbsParams)
258 , spline_(eps_, 1.0, // x1, x2
259 eps_*eps_*eps_, 1.0, // y1, y2
260 3.0*eps_*eps_, 0.0) // m1, m2
261 {}
262
267 static Params makeParams(const std::string& paramGroup)
268 {
269 const auto gamma = getParamFromGroup<Scalar>(paramGroup, "HeatpipeLawGamma");
270 const auto p0 = getParamFromGroup<Scalar>(paramGroup, "HeatpipeLawP0");
271 return {gamma, p0};
272 }
273
279 template<class Scalar>
280 Scalar pc(Scalar swe) const
281 {
282 const Scalar sne = 1 - swe;
283 const Scalar p0Gamma = params_.p0()*params_.gamma();
284
285 // regularization
286 if (sne >= 1.0)
287 {
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;
291 }
292 else if (sne <= 0.0)
293 return sne * p0Gamma*1.417;
294 else
295 return p0Gamma*((1.263*sne - 2.120)*sne + 1.417) * sne;
296 }
297
301 template<class Scalar>
303 { return 0.0; }
304
311 template<class Scalar>
313 {
314 const Scalar sne = 1 - swe;
315 const Scalar p0Gamma = params_.p0()*params_.gamma();
316 if (sne > 1.0)
317 sne = 1.0;
318 else if (sne <= 0.0)
319 return -p0Gamma*1.417;
320 else
321 return - p0Gamma*((3*1.263*sne - 2*2.120)*sne + 1.417);
322 }
323
330 template<class Scalar>
331 Scalar krw(Scalar swe) const
332 {
333 return kr_(swe);
334 }
335
342 template<class Scalar>
343 Scalar krn(Scalar swe) const
344 {
345 const Scalar sne = 1 - swe; // TODO does this make sense?
346 return kr_(sne);
347 }
348
353 {
354 return params_ == o.params_
355 && effToAbsParams_ == o.effToAbsParams_;
356 }
357
359 { return effToAbsParams_; }
360
361private:
362
363 Scalar kr_(Scalar s) const
364 {
365 if (s >= 1.0)
366 return 1;
367 else if (s <= 0.0)
368 return 0;
369 else if (s > eps_)
370 return spline_.eval(s); // regularize
371 else
372 return s*s*s;
373 }
374
375 Params params_;
376 EffToAbsParams effToAbsParams_;
377 static constexpr Scalar eps_ = 0.95;
378 Dumux::Spline<Scalar> spline_;
379};
380} // end namespace Dumux::FluidMatrix
381
382#endif
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: adapt.hh:29
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 &params, Scalar Sw)
Returns the partial derivative of the capillary pressure to the effective saturation.
Definition: heatpipelaw.hh:99
static Scalar krn(const Params &params, 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 &params, Scalar Sw)
The relative permeability for the wetting phase.
Definition: heatpipelaw.hh:131
ParamsT Params
Definition: heatpipelaw.hh:52
static Scalar pc(const Params &params, Scalar Sw)
The capillary pressure-saturation curve.
Definition: heatpipelaw.hh:61
static Scalar Sw(const Params &params, Scalar pC)
The saturation-capillary pressure curve.
Definition: heatpipelaw.hh:88
static Scalar dSw_dpC(const Params &params, 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 &paramGroup)
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 &paramGroup)
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 &params, 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