Region2.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_REGION2_HPP
28#define OPM_IAPWS_REGION2_HPP
29
31
32#include <cmath>
33
34namespace Opm {
35namespace IAPWS {
50template <class Scalar>
52{
53public:
61 template <class Evaluation>
62 static bool isValid(const Evaluation& temperature, const Evaluation& pressure)
63 {
64 return
65 temperature <= 623.15 && pressure <= 100e6;
66
67 // actually this is:
68 /*
69 return
70 (273.15 <= temperature && temperature <= 623.15 && pressure <= vaporPressure(temperature)) ||
71 (623.15 < temperature && temperature <= 863.15 && pressure <= auxPressure(temperature)) ||
72 (863.15 < temperature && temperature <= 1073.15 && pressure < 100e6);
73 */
74 }
75
81 template <class Evaluation>
82 static Evaluation tau(const Evaluation& temperature)
83 { return 540.0 / temperature; }
84
91 template <class Evaluation>
92 static Evaluation dtau_dT(const Evaluation& temperature)
93 { return - 540.0 / (temperature*temperature); }
94
100 template <class Evaluation>
101 static Evaluation pi(const Evaluation& pressure)
102 { return pressure / 1e6; }
103
110 template <class Evaluation>
111 static Scalar dpi_dp(const Evaluation& /*pressure*/)
112 { return 1.0 / 1e6; }
113
120 template <class Evaluation>
121 static Evaluation dp_dpi(const Evaluation& /*pressure*/)
122 { return 1e6; }
123
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;
142
143 // ideal gas part
144 result = log(pi_);
145 for (int i = 0; i < 9; ++i)
146 result += n_g(i)*pow(tau_, J_g(i));
147
148 // residual part
149 for (int i = 0; i < 43; ++i)
150 result +=
151 n_r(i)*
152 pow(pi_, I_r(i))*
153 pow(tau_ - 0.5, J_r(i));
154 return result;
155 }
156
169 template <class Evaluation>
170 static Evaluation dgamma_dtau(const Evaluation& temperature, const Evaluation& pressure)
171 {
172 const Evaluation& tau_ = tau(temperature); /* reduced temperature */
173 const Evaluation& pi_ = pi(pressure); /* reduced pressure */
174
175 // ideal gas part
176 Evaluation result = 0.0;
177 for (int i = 0; i < 9; i++) {
178 result +=
179 n_g(i) *
180 J_g(i) *
181 pow(tau_, static_cast<Scalar>(J_g(i) - 1));
182 }
183
184 // residual part
185 for (int i = 0; i < 43; i++) {
186 result +=
187 n_r(i) *
188 pow(pi_, static_cast<Scalar>(I_r(i))) *
189 J_r(i) *
190 pow(tau_ - 0.5, static_cast<Scalar>(J_r(i) - 1));
191 }
192
193 return result;
194 }
195
208 template <class Evaluation>
209 static Evaluation dgamma_dpi(const Evaluation& temperature, const Evaluation& pressure)
210 {
211 const Evaluation& tau_ = tau(temperature); /* reduced temperature */
212 const Evaluation& pi_ = pi(pressure); /* reduced pressure */
213
214 // ideal gas part
215 Evaluation result = 1/pi_;
216
217 // residual part
218 for (int i = 0; i < 43; i++) {
219 result +=
220 n_r(i) *
221 I_r(i) *
222 pow(pi_, static_cast<Scalar>(I_r(i) - 1)) *
223 pow(tau_ - 0.5, static_cast<Scalar>(J_r(i)));
224 }
225
226 return result;
227 }
228
241 template <class Evaluation>
242 static Evaluation ddgamma_dtaudpi(const Evaluation& temperature, const Evaluation& pressure)
243 {
244 const Evaluation& tau_ = tau(temperature); /* reduced temperature */
245 const Evaluation& pi_ = pi(pressure); /* reduced pressure */
246
247 // ideal gas part
248 Evaluation result = 0.0;
249
250 // residual part
251 for (int i = 0; i < 43; i++) {
252 result +=
253 n_r(i) *
254 I_r(i) *
255 J_r(i) *
256 pow(pi_, static_cast<Scalar>(I_r(i) - 1)) *
257 pow(tau_ - 0.5, static_cast<Scalar>(J_r(i) - 1));
258 }
259
260 return result;
261 }
262
275 template <class Evaluation>
276 static Evaluation ddgamma_ddpi(const Evaluation& temperature, const Evaluation& pressure)
277 {
278 const Evaluation& tau_ = tau(temperature); /* reduced temperature */
279 const Evaluation& pi_ = pi(pressure); /* reduced pressure */
280
281 // ideal gas part
282 Evaluation result = -1/(pi_*pi_);
283
284 // residual part
285 for (int i = 0; i < 43; i++) {
286 result +=
287 n_r(i) *
288 I_r(i) *
289 (I_r(i) - 1) *
290 pow(pi_, static_cast<Scalar>(I_r(i) - 2)) *
291 pow(tau_ - 0.5, static_cast<Scalar>(J_r(i)));
292 }
293
294 return result;
295 }
296
309 template <class Evaluation>
310 static Evaluation ddgamma_ddtau(const Evaluation& temperature, const Evaluation& pressure)
311 {
312 const Evaluation& tau_ = tau(temperature); /* reduced temperature */
313 const Evaluation& pi_ = pi(pressure); /* reduced pressure */
314
315 // ideal gas part
316 Evaluation result = 0.0;
317 for (int i = 0; i < 9; i++) {
318 result +=
319 n_g(i) *
320 J_g(i) *
321 (J_g(i) - 1) *
322 pow(tau_, static_cast<Scalar>(J_g(i) - 2));
323 }
324
325 // residual part
326 for (int i = 0; i < 43; i++) {
327 result +=
328 n_r(i) *
329 pow(pi_, I_r(i)) *
330 J_r(i) *
331 (J_r(i) - 1.) *
332 pow(tau_ - 0.5, static_cast<Scalar>(J_r(i) - 2));
333 }
334
335 return result;
336 }
337
338
339private:
340 static Scalar n_g(int i)
341 {
342 static const Scalar n[9] = {
343 -0.96927686500217e1, 0.10086655968018e2, -0.56087911283020e-2,
344 0.71452738081455e-1, -0.40710498223928, 0.14240819171444e1,
345 -0.43839511319450e1, -0.28408632460772, 0.21268463753307e-1
346 };
347 return n[i];
348 }
349
350 static Scalar n_r(int i)
351 {
352 static const Scalar n[43] = {
353 -0.17731742473213e-2, -0.17834862292358e-1, -0.45996013696365e-1,
354 -0.57581259083432e-1, -0.50325278727930e-1, -0.33032641670203e-4,
355 -0.18948987516315e-3, -0.39392777243355e-2, -0.43797295650573e-1,
356 -0.26674547914087e-4, 0.20481737692309e-7, 0.43870667284435e-6,
357 -0.32277677238570e-4, -0.15033924542148e-2, -0.40668253562649e-1,
358 -0.78847309559367e-9, 0.12790717852285e-7, 0.48225372718507e-6,
359 0.22922076337661e-5, -0.16714766451061e-10, -0.21171472321355e-2,
360 -0.23895741934104e2, -0.59059564324270e-17, -0.12621808899101e-5,
361 -0.38946842435739e-1, 0.11256211360459e-10, -0.82311340897998e1,
362 0.19809712802088e-7, 0.10406965210174e-18, -0.10234747095929e-12,
363 -0.10018179379511e-8, -0.80882908646985e-10, 0.10693031879409,
364 -0.33662250574171, 0.89185845355421e-24, 0.30629316876232e-12,
365 -0.42002467698208e-5, -0.59056029685639e-25, 0.37826947613457e-5,
366 -0.12768608934681e-14, 0.73087610595061e-28, 0.55414715350778e-16,
367 -0.94369707241210e-6
368 };
369 return n[i];
370 }
371
372 static Scalar I_r(int i)
373 {
374 static const short int I[43] = {
375 1, 1, 1,
376 1, 1, 2,
377 2, 2, 2,
378 2, 3, 3,
379 3, 3, 3,
380 4, 4, 4,
381 5, 6, 6,
382 6, 7, 7,
383 7, 8, 8,
384 9, 10, 10,
385 10, 16, 16,
386 18, 20, 20,
387 20, 21, 22,
388 23, 24, 24,
389 24
390 };
391 return I[i];
392 }
393
394 static Scalar J_g(int i)
395 {
396 static const short int J[9] = {
397 0, 1, -5,
398 -4, -3, -2,
399 -1, 2, 3
400 };
401 return J[i];
402 }
403
404 static Scalar J_r(int i)
405 {
406 static const short int J[43] = {
407 0, 1, 2,
408 3, 6, 1,
409 2, 4, 7,
410 36, 0, 1,
411 3, 6, 35,
412 1, 2, 3,
413 7, 3, 16,
414 35, 0, 11,
415 25, 8, 36,
416 13, 4, 10,
417 14, 29, 50,
418 57, 20, 35,
419 48, 21, 53,
420 39, 26, 40,
421 58
422 };
423 return J[i];
424 }
425
426};
427
428} // namespace IAPWS
429} // namespace Opm
430
431#endif
A traits class which provides basic mathematical functions for arbitrary scalar floating point values...
Implements the equations for region 2 of the IAPWS '97 formulation.
Definition: Region2.hpp:52
static Scalar dpi_dp(const Evaluation &)
Returns the derivative of the reduced pressure to the pressure for IAPWS region 2 in .
Definition: Region2.hpp:111
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: Region2.hpp:310
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 2 ...
Definition: Region2.hpp:276
static Evaluation dtau_dT(const Evaluation &temperature)
Returns the derivative of the reduced temperature to the temperature for IAPWS region 2.
Definition: Region2.hpp:92
static Evaluation pi(const Evaluation &pressure)
Returns the reduced pressure (dimensionless) for IAPWS region 2.
Definition: Region2.hpp:101
static Evaluation dp_dpi(const Evaluation &)
Returns the derivative of the pressure to the reduced pressure for IAPWS region 2 (dimensionless).
Definition: Region2.hpp:121
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 2 (i....
Definition: Region2.hpp:170
static Evaluation tau(const Evaluation &temperature)
Returns the reduced temperature (dimensionless) for IAPWS region 2.
Definition: Region2.hpp:82
static Evaluation gamma(const Evaluation &temperature, const Evaluation &pressure)
The Gibbs free energy for IAPWS region 2 (i.e. sub-critical steam) (dimensionless).
Definition: Region2.hpp:136
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: Region2.hpp:242
static bool isValid(const Evaluation &temperature, const Evaluation &pressure)
Returns true if IAPWS region 2 applies for a (temperature, pressure) pair.
Definition: Region2.hpp:62
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 2 (i....
Definition: Region2.hpp:209
Definition: Air_Mesitylene.hpp:34
Evaluation log(const Evaluation &value)
Definition: MathToolbox.hpp:407
ReturnEval_< Evaluation1, Evaluation2 >::type pow(const Evaluation1 &base, const Evaluation2 &exp)
Definition: MathToolbox.hpp:416