3.5-git
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#include <cmath>
29#include <algorithm>
30
35
36namespace Dumux::FluidMatrix {
37
43template<class ScalarType, class EffToAbsPolicy = TwoPEffToAbsDefaultPolicy>
44class HeatPipeLaw : Adapter<HeatPipeLaw<ScalarType, EffToAbsPolicy>, PcKrSw>
45{
46
47public:
48 using Scalar = ScalarType;
49 using EffToAbsParams = typename EffToAbsPolicy::template Params<Scalar>;
50 using EffToAbs = EffToAbsPolicy;
51
55 static constexpr int numFluidPhases()
56 { return 2; }
57
61 static constexpr bool isRegularized()
62 { return false; }
63
68 struct Params
69 {
71 : gamma_(gamma), p0_(p0)
72 {}
73
74 Scalar gamma() const { return gamma_; }
75 void setGamma(Scalar gamma) { gamma_ = gamma; }
76
77 Scalar p0() const { return p0_; }
78 void setP0(Scalar p0) { p0_ = p0; }
79
80 bool operator== (const Params& p) const
81 {
82 return Dune::FloatCmp::eq(gamma(), p.gamma(), 1e-6)
83 && Dune::FloatCmp::eq(p0(), p.p0(), 1e-6);
84 }
85
86 private:
87 Scalar gamma_, p0_;
88 };
89
94 HeatPipeLaw() = delete;
95
100 explicit HeatPipeLaw(const std::string& paramGroup)
101 {
102 params_ = makeParams(paramGroup);
103 effToAbsParams_ = EffToAbs::template makeParams<Scalar>(paramGroup);
104 }
105
110 HeatPipeLaw(const Params& params,
111 const EffToAbsParams& effToAbsParams = {})
112 : params_(params)
113 , effToAbsParams_(effToAbsParams)
114 , spline_(eps_, 1.0, // x1, x2
115 eps_*eps_*eps_, 1.0, // y1, y2
116 3.0*eps_*eps_, 0.0) // m1, m2
117 {}
118
123 static Params makeParams(const std::string& paramGroup)
124 {
125 const auto gamma = getParamFromGroup<Scalar>(paramGroup, "HeatpipeLawGamma");
126 const auto p0 = getParamFromGroup<Scalar>(paramGroup, "HeatpipeLawP0");
127 return {gamma, p0};
128 }
129
135 Scalar pc(const Scalar sw) const
136 {
137 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
138 const Scalar sne = 1 - swe;
139 const Scalar p0Gamma = params_.p0()*params_.gamma();
140
141 // regularization
142 if (sne >= 1.0)
143 {
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;
147 }
148 else if (sne <= 0.0)
149 return sne * p0Gamma*1.417;
150 else
151 return p0Gamma*((1.263*sne - 2.120)*sne + 1.417) * sne;
152 }
153
158 { return 0.0; }
159
166 Scalar dpc_dsw(const Scalar sw) const
167 {
168 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
169 const Scalar sne = 1 - swe;
170 const Scalar p0Gamma = params_.p0()*params_.gamma();
171 if (sne > 1.0)
172 sne = 1.0;
173 else if (sne <= 0.0)
174 return -p0Gamma*1.417*EffToAbs::dswe_dsw(effToAbsParams_);
175 else
176 return - p0Gamma*((3*1.263*sne - 2*2.120)*sne + 1.417)*EffToAbs::dswe_dsw(effToAbsParams_);
177 }
178
185 Scalar krw(const Scalar sw) const
186 {
187 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
188 return kr_(swe);
189 }
190
197 Scalar krn(const Scalar sw) const
198 {
199 const Scalar swe = EffToAbs::swToSwe(sw, effToAbsParams_);
200 const Scalar sne = 1 - swe; // TODO does this make sense?
201 return kr_(sne);
202 }
203
208 {
209 return params_ == o.params_
210 && effToAbsParams_ == o.effToAbsParams_;
211 }
212
214 { return effToAbsParams_; }
215
216private:
217
218 Scalar kr_(Scalar s) const
219 {
220 if (s >= 1.0)
221 return 1;
222 else if (s <= 0.0)
223 return 0;
224 else if (s > eps_)
225 return spline_.eval(s); // regularize
226 else
227 return s*s*s;
228 }
229
230 Params params_;
231 EffToAbsParams effToAbsParams_;
232 static constexpr Scalar eps_ = 0.95;
233 Dumux::Spline<Scalar> spline_;
234};
235} // end namespace Dumux::FluidMatrix
236
237#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.
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 &paramGroup)
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 &paramGroup)
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 &params, 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