Region1.hpp
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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_IAPWS_REGION1_HPP
28#define OPM_IAPWS_REGION1_HPP
29
31
32#include <cmath>
33
34namespace Opm {
35namespace IAPWS {
49template <class Scalar>
51{
52public:
60 template <class Evaluation>
61 static bool isValid(const Evaluation& temperature, const Evaluation& pressure)
62 {
63 return temperature <= 623.15 && pressure <= 100e6;
64
65 // actually this is:
66 /*
67 return
68 (
69 273.15 <= temperature &&
70 temperature <= 623.15 &&
71 pressure >= vaporPressure(temperature) &&
72 pressure <= 100e6
73 );
74 */
75 }
76
82 template <class Evaluation>
83 static Evaluation tau(const Evaluation& temperature)
84 { return 1386.0 / temperature; }
85
92 template <class Evaluation>
93 static Evaluation dtau_dT(const Evaluation& temperature)
94 { return - 1386.0 / (temperature*temperature); }
95
101 template <class Evaluation>
102 static Evaluation pi(const Evaluation& pressure)
103 { return pressure / 16.53e6; }
104
111 template <class Evaluation>
112 static Scalar dpi_dp(const Evaluation& /*pressure*/)
113 { return 1.0 / 16.53e6; }
114
121 template <class Evaluation>
122 static Scalar dp_dpi(const Evaluation& /*pressure*/)
123 { return 16.53e6; }
124
135 template <class Evaluation>
136 static Evaluation gamma(const Evaluation& temperature, const Evaluation& pressure)
137 {
138 const Evaluation tau_ = tau(temperature); /* reduced temperature */
139 const Evaluation pi_ = pi(pressure); /* reduced pressure */
140
141 Evaluation result = 0;
142 for (int i = 0; i < 34; ++i) {
143 result += n(i)*pow(7.1 - pi_, I(i))*pow(tau_ - 1.222, J(i));
144 }
145
146 return result;
147 }
148
149
161 template <class Evaluation>
162 static Evaluation dgamma_dtau(const Evaluation& temperature, const Evaluation& pressure)
163 {
164 const Evaluation tau_ = tau(temperature); /* reduced temperature */
165 const Evaluation pi_ = pi(pressure); /* reduced pressure */
166
167 Evaluation result = 0.0;
168 for (int i = 0; i < 34; i++) {
169 result +=
170 n(i) *
171 pow(7.1 - pi_, static_cast<Scalar>(I(i))) *
172 pow(tau_ - 1.222, static_cast<Scalar>(J(i)-1)) *
173 J(i);
174 }
175
176 return result;
177 }
178
190 template <class Evaluation>
191 static Evaluation dgamma_dpi(const Evaluation& temperature, const Evaluation& pressure)
192 {
193 const Evaluation tau_ = tau(temperature); /* reduced temperature */
194 const Evaluation pi_ = pi(pressure); /* reduced pressure */
195
196 Evaluation result = 0.0;
197 for (int i = 0; i < 34; i++) {
198 result +=
199 -n(i) *
200 I(i) *
201 pow(7.1 - pi_, static_cast<Scalar>(I(i) - 1)) *
202 pow(tau_ - 1.222, static_cast<Scalar>(J(i)));
203 }
204
205 return result;
206 }
207
220 template <class Evaluation>
221 static Evaluation ddgamma_dtaudpi(const Evaluation& temperature, const Evaluation& pressure)
222 {
223 const Evaluation tau_ = tau(temperature); /* reduced temperature */
224 const Evaluation pi_ = pi(pressure); /* reduced pressure */
225
226 Evaluation result = 0.0;
227 for (int i = 0; i < 34; i++) {
228 result +=
229 -n(i) *
230 I(i) *
231 J(i) *
232 pow(7.1 - pi_, static_cast<Scalar>(I(i) - 1)) *
233 pow(tau_ - 1.222, static_cast<Scalar>(J(i) - 1));
234 }
235
236 return result;
237 }
238
251 template <class Evaluation>
252 static Evaluation ddgamma_ddpi(const Evaluation& temperature, const Evaluation& pressure)
253 {
254 const Evaluation tau_ = tau(temperature); /* reduced temperature */
255 const Evaluation pi_ = pi(pressure); /* reduced pressure */
256
257 Evaluation result = 0.0;
258 for (int i = 0; i < 34; i++) {
259 result +=
260 n(i) *
261 I(i) *
262 (I(i) - 1) *
263 pow(7.1 - pi_, I(i) - 2) *
264 pow(tau_ - 1.222, J(i));
265 }
266
267 return result;
268 }
269
281 template <class Evaluation>
282 static Evaluation ddgamma_ddtau(const Evaluation& temperature, const Evaluation& pressure)
283 {
284 const Evaluation tau_ = tau(temperature); /* reduced temperature */
285 const Evaluation pi_ = pi(pressure); /* reduced pressure */
286
287 Evaluation result = 0.0;
288 for (int i = 0; i < 34; i++) {
289 result +=
290 n(i) *
291 pow(7.1 - pi_, I(i)) *
292 J(i) *
293 (J(i) - 1) *
294 pow(tau_ - 1.222, J(i) - 2);
295 }
296
297 return result;
298 }
299
300private:
301 static Scalar n(int i)
302 {
303 static const Scalar n[34] = {
304 0.14632971213167, -0.84548187169114, -0.37563603672040e1,
305 0.33855169168385e1, -0.95791963387872, 0.15772038513228,
306 -0.16616417199501e-1, 0.81214629983568e-3, 0.28319080123804e-3,
307 -0.60706301565874e-3, -0.18990068218419e-1, -0.32529748770505e-1,
308 -0.21841717175414e-1, -0.52838357969930e-4, -0.47184321073267e-3,
309 -0.30001780793026e-3, 0.47661393906987e-4, -0.44141845330846e-5,
310 -0.72694996297594e-15,-0.31679644845054e-4, -0.28270797985312e-5,
311 -0.85205128120103e-9, -0.22425281908000e-5, -0.65171222895601e-6,
312 -0.14341729937924e-12,-0.40516996860117e-6, -0.12734301741641e-8,
313 -0.17424871230634e-9, -0.68762131295531e-18, 0.14478307828521e-19,
314 0.26335781662795e-22,-0.11947622640071e-22, 0.18228094581404e-23,
315 -0.93537087292458e-25
316 };
317 return n[i];
318 }
319
320 static Scalar I(int i)
321 {
322 static const short int I[34] = {
323 0, 0, 0,
324 0, 0, 0,
325 0, 0, 1,
326 1, 1, 1,
327 1, 1, 2,
328 2, 2, 2,
329 2, 3, 3,
330 3, 4, 4,
331 4, 5, 8,
332 8, 21, 23,
333 29, 30, 31,
334 32
335 };
336 return I[i];
337 }
338
339 static Scalar J(int i)
340 {
341 static const short int J[34] = {
342 -2, -1, 0,
343 1, 2, 3,
344 4, 5, -9,
345 -7, -1, 0,
346 1, 3, -3,
347 0, 1, 3,
348 17, -4, 0,
349 6, -5, -2,
350 10, -8, -11,
351 -6, -29, -31,
352 -38, -39, -40,
353 -41
354 };
355 return J[i];
356 }
357
358};
359
360} // namespace IAPWS
361} // namespace Opm
362
363#endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Implements the equations for region 1 of the IAPWS '97 formulation.
Definition: Region1.hpp:51
static Evaluation ddgamma_ddpi(const Evaluation &temperature, const Evaluation &pressure)
The second partial derivative of the Gibbs free energy to the normalized pressure for IAPWS region 1 ...
Definition: Region1.hpp:252
static Evaluation dtau_dT(const Evaluation &temperature)
Returns the derivative of the reduced temperature to the temperature for IAPWS region 1 in .
Definition: Region1.hpp:93
static Evaluation tau(const Evaluation &temperature)
Returns the reduced temperature for IAPWS region 1.
Definition: Region1.hpp:83
static Evaluation ddgamma_ddtau(const Evaluation &temperature, const Evaluation &pressure)
The second partial derivative of the Gibbs free energy to the normalized temperature for IAPWS region...
Definition: Region1.hpp:282
static Scalar dp_dpi(const Evaluation &)
Returns the derivative of the pressure to the reduced pressure for IAPWS region 1 in .
Definition: Region1.hpp:122
static Evaluation pi(const Evaluation &pressure)
Returns the reduced pressure for IAPWS region 1.
Definition: Region1.hpp:102
static Evaluation gamma(const Evaluation &temperature, const Evaluation &pressure)
The Gibbs free energy (dimensionless) for IAPWS region 1 (i.e. liquid)
Definition: Region1.hpp:136
static Evaluation dgamma_dtau(const Evaluation &temperature, const Evaluation &pressure)
The partial derivative of the Gibbs free energy to the normalized temperature for IAPWS region 1 (i....
Definition: Region1.hpp:162
static bool isValid(const Evaluation &temperature, const Evaluation &pressure)
Returns true if IAPWS region 1 applies for a (temperature in , pressure in ) pair.
Definition: Region1.hpp:61
static Evaluation ddgamma_dtaudpi(const Evaluation &temperature, const Evaluation &pressure)
The partial derivative of the Gibbs free energy to the normalized pressure and to the normalized temp...
Definition: Region1.hpp:221
static Scalar dpi_dp(const Evaluation &)
Returns the derivative of the reduced pressure to the pressure for IAPWS region 1 in .
Definition: Region1.hpp:112
static Evaluation dgamma_dpi(const Evaluation &temperature, const Evaluation &pressure)
The partial derivative of the Gibbs free energy to the normalized pressure for IAPWS region 1 (i....
Definition: Region1.hpp:191
Definition: Air_Mesitylene.hpp:34
ReturnEval_< Evaluation1, Evaluation2 >::type pow(const Evaluation1 &base, const Evaluation2 &exp)
Definition: MathToolbox.hpp:416