opm-common
Loading...
Searching...
No Matches
H2.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*/
31#ifndef OPM_H2_HPP
32#define OPM_H2_HPP
33
38
39#include <cmath>
40
41namespace Opm {
42
43// this class collects all the H2tabulated quantities in one convenient place
44struct H2Tables {
45 static const Opm::UniformTabulated2DFunction<double> tabulatedEnthalpy;
46 static const Opm::UniformTabulated2DFunction<double> tabulatedDensity;
47 static constexpr double brineSalinity = 1.000000000000000e-01;
48};
49
51 typedef double Scalar;
52 static const char *name;
53 static const int numX = 200;
54 static const Scalar xMin;
55 static const Scalar xMax;
56 static const int numY = 500;
57 static const Scalar yMin;
58 static const Scalar yMax;
59
60 static const Scalar vals[200][500];
61};
62
64 typedef double Scalar;
65 static const char *name;
66 static const int numX = 200;
67 static const Scalar xMin;
68 static const Scalar xMax;
69 static const int numY = 500;
70 static const Scalar yMin;
71 static const Scalar yMax;
72 static const Scalar vals[200][500];
73};
74
88template <class Scalar>
89class H2 : public Component<Scalar, H2<Scalar> >
90{
91 using IdealGas = Opm::IdealGas<Scalar>;
92
93 static const UniformTabulated2DFunction<double>& tabulatedEnthalpy;
94 static const UniformTabulated2DFunction<double>& tabulatedDensity;
95
96public:
97 // For H2Tables class
98 static const Scalar brineSalinity;
99
103 static std::string name()
104 { return "H2"; }
105
109 static constexpr Scalar molarMass()
110 { return 2.01588e-3; }
111
115 static Scalar criticalTemperature()
116 { return 33.145; /* [K] */ }
117
121 static Scalar criticalPressure()
122 { return 1.2964e6; /* [N/m^2] */ }
123
127 static Scalar criticalDensity()
128 { return 15.508e-3; /* [mol/cm^3] */ }
129
133 static Scalar tripleTemperature()
134 { return 13.957; /* [K] */ }
135
139 static Scalar triplePressure()
140 { return 0.00736e6; /* [N/m^2] */ }
141
145 static Scalar tripleDensity()
146 { return 38.2e-3; /* [mol/cm^3] */ }
147
151 static Scalar criticalVolume() {return 6.45e-2; }
152
156 static Scalar acentricFactor() { return -0.22; }
157
165 template <class Evaluation>
166 static Evaluation vaporPressure(Evaluation temperature)
167 {
168 if (temperature > criticalTemperature())
169 return criticalPressure();
170 if (temperature < tripleTemperature())
171 return 0; // H2 is solid: We don't take sublimation into
172 // account
173
174 // Intermediate calculations involving temperature
175 Evaluation sigma = 1 - temperature/criticalTemperature();
176 Evaluation T_recp = criticalTemperature() / temperature;
177
178 // Parameters for normal hydrogen in Table 8
179 static const Scalar N[4] = {-4.89789, 0.988558, 0.349689, 0.499356};
180 static const Scalar k[4] = {1.0, 1.5, 2.0, 2.85};
181
182 // Eq. (33)
183 Evaluation s = 0.0; // sum calculation
184 for (int i = 0; i < 4; ++i) {
185 s += N[i] * pow(sigma, k[i]);
186 }
187 Evaluation lnPsigmaPc = T_recp * s;
188
189 return exp(lnPsigmaPc) * criticalPressure();
190 }
191
199 template <class Evaluation>
200 static Evaluation gasDensity(Evaluation temperature, Evaluation pressure, bool extrapolate = false)
201 {
202 return tabulatedDensity.eval(temperature, pressure, extrapolate);
203 }
204
212 template <class Evaluation>
213 static Evaluation gasMolarDensity(Evaluation temperature, Evaluation pressure, bool extrapolate = false)
214 { return gasDensity(temperature, pressure, extrapolate) / molarMass(); }
215
219 static constexpr bool gasIsCompressible()
220 { return true; }
221
225 static constexpr bool gasIsIdeal()
226 { return false; }
227
234 template <class Evaluation>
235 static Evaluation gasPressure(Evaluation temperature, Evaluation density)
236 {
237 // Eq. (56) in Span et al. (2000)
238 Scalar R = IdealGas::R;
239 Evaluation rho_red = density / (molarMass() * criticalDensity() * 1e6);
240 Evaluation T_red = H2::criticalTemperature() / temperature;
241 return rho_red * criticalDensity() * R * temperature
242 * (1 + rho_red * derivResHelmholtzWrtRedRho(T_red, rho_red));
243 }
244
248 template <class Evaluation>
249 static Evaluation gasInternalEnergy(const Evaluation& temperature,
250 const Evaluation& pressure,
251 bool extrapolate = false)
252 {
253 const Evaluation h = gasEnthalpy(temperature, pressure, extrapolate);
254 const Evaluation rho = gasDensity(temperature, pressure, extrapolate);
255
256 return h - (pressure / rho);
257 }
258
266 template <class Evaluation>
267 static const Evaluation gasEnthalpy(Evaluation temperature,
268 Evaluation pressure,
269 bool extrapolate = false)
270 {
271 return tabulatedEnthalpy.eval(temperature, pressure, extrapolate);
272 }
273
284 template <class Evaluation>
285 static Evaluation gasViscosity(const Evaluation& temperature,
286 const Evaluation& pressure,
287 bool extrapolate = false)
288 {
289 // Some needed parameters and variables
290 const Scalar M = molarMass() * 1e3; // g/mol
291 const Scalar epsilon_div_kb = 30.41;
292 const Scalar sigma = 0.297; // nm
293 const Scalar Na = 6.022137e23; // Avogadro's number
294
295 Evaluation T_star = temperature / epsilon_div_kb;
296 Evaluation ln_T_star = log(T_star);
297 Evaluation T_r = temperature / criticalTemperature();
298 Evaluation rho = gasDensity(temperature, pressure, extrapolate);
299 Evaluation rho_r = rho / 90.909090909; // see corrections
300
301 //
302 // Eta_0 terms (zero-density viscocity)
303 //
304 // Parameters in Table 2 for Eq. (4)
305 static constexpr Scalar a[5] =
306 {2.0963e-1, -4.55274e-1, 1.43602e-1, -3.35325e-2, 2.76981e-3};
307
308 // Eq. (4)
309 Evaluation ln_S_star = 0.0;
310 for (int i = 0; i < 5; ++i) {
311 ln_S_star += a[i] * pow(ln_T_star, i);
312 }
313
314 // Eq. (3)
315 Evaluation eta_0 = 0.021357 * sqrt(M * temperature) / (sigma * sigma * exp(ln_S_star));
316
317 //
318 // Eta_1 terms (excess contribution)
319 //
320 //
321 // Parameters in Table 3 for Eq. (7)
322 static constexpr Scalar b[7] =
323 {-0.187, 2.4871, 3.7151, -11.0972, 9.0965, -3.8292, 0.5166};
324
325 // Eq. (7) with corrections
326 Evaluation B_star = 0.0;
327 for (int i = 0; i < 7; ++i) {
328 B_star += b[i] * pow(T_star, -i);
329 }
330
331 // Eq. (9) (eta_1 part) using Eq. (5-7) with corrections and sigma in m (instead of nm). Note that Na * sigma_m
332 // ^ 3 have unit [m3/mol], so we need to multiply with molar density (= rho / molarMass) and not mass density in
333 // Eq. (9)
334 const Scalar sigma_m = sigma * 1e-9;
335 Evaluation eta_1 = B_star * Na * sigma_m * sigma_m * sigma_m * eta_0 * rho / M;
336
337 //
338 // delta eta_h terms (higher-order effects)
339 //
340 // Parameters in Table 4 for Eq. (9)
341 static constexpr Scalar c[6] =
342 {6.43449673, 4.56334068e-2, 2.32797868e-1, 9.58326120e-1, 1.27941189e-1, 3.63576595e-1};
343
344 // delta eta_h terms in Eq. (9)
345 Evaluation delta_eta_h = c[0] * rho_r * rho_r * exp(c[1] * T_r + c[2] / T_r +
346 (c[3] * rho_r * rho_r) / (c[4] + T_r) + c[5] * pow(rho_r, 6));
347
348 // return all terms in Eq. (9) converted from muPa*s to Pa*s
349 return (eta_0 + eta_1 + delta_eta_h) * 1e-6;
350 }
351
360 template <class Evaluation>
361 static const Evaluation gasHeatCapacity(Evaluation temperature,
362 Evaluation pressure)
363 {
364 // Reduced variables
365 Evaluation rho_red = reducedMolarDensity(temperature, pressure);
366 Evaluation T_red = criticalTemperature() / temperature;
367
368 // Need Eq. (62) in Span et al. (2000)
369 Evaluation cv = gasIsochoricHeatCapacity(temperature, pressure); // [J/(kg*K)]
370
371 // Some intermediate calculations
372 Evaluation numerator = pow(1 + rho_red * derivResHelmholtzWrtRedRho(T_red, rho_red)
373 - rho_red * T_red * secDerivResHelmholtzWrtRecipRedTempAndRedRho(T_red, rho_red), 2);
374
375 Evaluation denominator = 1 + 2 * rho_red * derivResHelmholtzWrtRedRho(T_red, rho_red)
376 + pow(rho_red, 2) * secDerivResHelmholtzWrtRedRho(T_red, rho_red);
377
378 // Eq. (63) in Span et al. (2000).
379 Scalar R = IdealGas::R;
380 Evaluation cp = cv + R * (numerator / denominator) / molarMass(); // divide by M to get [J/(kg*K)]
381
382 // Return
383 return cp;
384 }
385
393 template <class Evaluation>
394 static const Evaluation gasIsochoricHeatCapacity(Evaluation temperature,
395 Evaluation pressure)
396 {
397 // Reduced variables
398 Evaluation rho_red = reducedMolarDensity(temperature, pressure);
399 Evaluation T_red = criticalTemperature() / temperature;
400
401 // Eq. (62) in Span et al. (2000)
402 Scalar R = IdealGas::R;
403 Evaluation cv = R * (-pow(T_red, 2) * (secDerivIdealHelmholtzWrtRecipRedTemp(T_red)
404 + secDerivResHelmholtzWrtRecipRedTemp(T_red, rho_red))); // [J/(mol*K)]
405
406 return cv / molarMass();
407 }
408
416 template <class Evaluation>
417 static Evaluation reducedMolarDensity(const Evaluation& temperature,
418 const Evaluation& pg,
419 bool extrapolate = false)
420 {
421 return gasDensity(temperature, pg, extrapolate) / (molarMass() * criticalDensity() * 1e6);
422 }
423
430 template <class Evaluation>
431 static Evaluation idealGasPartHelmholtz(const Evaluation& T_red, const Evaluation& rho_red)
432 {
433 // Eq. (31), which can be compared with Eq. (53) in Span et al. (2000)
434 // Terms not in sum
435 Evaluation s1 = log(rho_red) + 1.5*log(T_red) + a_[0] + a_[1] * T_red;
436
437 // Sum term
438 Evaluation s2 = 0.0;
439 for (int i = 2; i < 7; ++i) {
440 s1 += a_[i] * log(1 - exp(b_[i-2] * T_red));
441 }
442
443 // Return total
444 Evaluation s = s1 + s2;
445 return s;
446 }
447
453 template <class Evaluation>
454 static Evaluation derivIdealHelmholtzWrtRecipRedTemp(const Evaluation& T_red)
455 {
456 // Derivative of Eq. (31) wrt. reciprocal reduced temperature, which can be compared with Eq. (79) in Span et
457 // al. (2000)
458 // Terms not in sum
459 Evaluation s1 = (1.5 / T_red) + a_[1];
460
461 // Sum term
462 Evaluation s2 = 0.0;
463 for (int i = 2; i < 7; ++i) {
464 s2 += (-a_[i] * b_[i-2] * exp(b_[i-2] * T_red)) / (1 - exp(b_[i-2] * T_red));
465 }
466
467 // Return total
468 Evaluation s = s1 + s2;
469 return s;
470 }
471
478 template <class Evaluation>
479 static Evaluation secDerivIdealHelmholtzWrtRecipRedTemp(const Evaluation& T_red)
480 {
481 // Second derivative of Eq. (31) wrt. reciprocal reduced temperature, which can be compared with Eq. (80) in
482 // Span et al. (2000)
483 // Sum term
484 Evaluation s1 = 0.0;
485 for (int i = 2; i < 7; ++i) {
486 s1 += (-a_[i] * pow(b_[i-2], 2) * exp(b_[i-2] * T_red)) / pow(1 - exp(b_[i-2] * T_red), 2);
487 }
488
489 // Return total
490 Evaluation s = (-1.5 / pow(T_red, 2)) + s1;
491 return s;
492 }
493
500 template <class Evaluation>
501 static Evaluation residualPartHelmholtz(const Evaluation& T_red, const Evaluation& rho_red)
502 {
503 // Eq. (32), which can be compared with Eq. (55) in Span et al. (2000)
504 // First sum term
505 Evaluation s1 = 0.0;
506 for (int i = 0; i < 7; ++i) {
507 s1 += N_[i] * pow(rho_red, d_[i]) * pow(T_red, t_[i]);
508 }
509
510 // Second sum term
511 Evaluation s2 = 0.0;
512 for (int i = 7; i < 9; ++i) {
513 s2 += N_[i] * pow(T_red, t_[i]) * pow(rho_red, d_[i]) * exp(-pow(rho_red, p_[i-7]));
514 }
515
516 // Third, and last, sum term
517 Evaluation s3 = 0.0;
518 for (int i = 9; i < 14; ++i) {
519 s3 += N_[i] * pow(T_red, t_[i]) * pow(rho_red, d_[i]) *
520 exp(phi_[i-9] * pow(rho_red - D_[i-9], 2) + beta_[i-9] * pow(T_red - gamma_[i-9], 2));
521 }
522
523 // Return total sum
524 Evaluation s = s1 + s2 + s3;
525 return s;
526 }
527
534 template <class Evaluation>
535 static Evaluation derivResHelmholtzWrtRedRho(const Evaluation& T_red, const Evaluation& rho_red)
536 {
537 // Derivative of Eq. (32) wrt to reduced density, which can be compared with Eq. (81) in Span et al. (2000)
538 // First sum term
539 Evaluation s1 = 0.0;
540 for (int i = 0; i < 7; ++i) {
541 s1 += d_[i] * N_[i] * pow(rho_red, d_[i]-1) * pow(T_red, t_[i]);
542 }
543
544 // Second sum term
545 Evaluation s2 = 0.0;
546 for (int i = 7; i < 9; ++i) {
547 s2 += N_[i] * pow(T_red, t_[i]) * pow(rho_red, d_[i]-1) * exp(-pow(rho_red, p_[i-7])) *
548 (d_[i] - p_[i-7]*pow(rho_red, p_[i-7]));
549 }
550
551 // Third, and last, sum term
552 Evaluation s3 = 0.0;
553 for (int i = 9; i < 14; ++i) {
554 s3 += N_[i] * pow(T_red, t_[i]) * pow(rho_red, d_[i]-1) *
555 exp(phi_[i-9] * pow(rho_red - D_[i-9], 2) + beta_[i-9] * pow(T_red - gamma_[i-9], 2)) *
556 (d_[i] + 2 * phi_[i-9] * rho_red * (rho_red - D_[i-9]));
557 }
558
559 // Return total sum
560 Evaluation s = s1 + s2 + s3;
561 return s;
562 }
563
570 template <class Evaluation>
571 static Evaluation secDerivResHelmholtzWrtRedRho(const Evaluation& T_red, const Evaluation& rho_red)
572 {
573 // Second derivative of Eq. (32) wrt to reduced density, which can be compared with Eq. (82) in Span et al.
574 // (2000)
575 // First sum term
576 Evaluation s1 = 0.0;
577 for (int i = 0; i < 7; ++i) {
578 s1 += d_[i] * (d_[i] - 1) * N_[i] * pow(rho_red, d_[i]-2) * pow(T_red, t_[i]);
579 }
580
581 // Second sum term
582 Evaluation s2 = 0.0;
583 for (int i = 7; i < 9; ++i) {
584 s2 += N_[i] * pow(T_red, t_[i]) * pow(rho_red, d_[i]-2) * exp(-pow(rho_red, p_[i-7])) *
585 ((d_[i] - p_[i-7] * pow(rho_red, p_[i-7])) * (d_[i] - p_[i-7] * pow(rho_red, p_[i-7]) - 1.0)
586 - pow(p_[i-7], 2) * pow(rho_red, p_[i-7]));
587 }
588
589 // Third, and last, sum term
590 Evaluation s3 = 0.0;
591 for (int i = 9; i < 14; ++i) {
592 s3 += N_[i] * pow(T_red, t_[i]) * pow(rho_red, d_[i]-2) *
593 exp(phi_[i-9] * pow(rho_red - D_[i-9], 2) + beta_[i-9] * pow(T_red - gamma_[i-9], 2)) *
594 (pow(d_[i] + 2 * phi_[i-9] * rho_red * (rho_red - D_[i-9]), 2)
595 - d_[i] + 2 * phi_[i-9] * pow(rho_red, 2));
596 }
597
598 // Return total sum
599 Evaluation s = s1 + s2 + s3;
600 return s;
601 }
602
609 template <class Evaluation>
610 static Evaluation derivResHelmholtzWrtRecipRedTemp(const Evaluation& T_red, const Evaluation& rho_red)
611 {
612 // Derivative of Eq. (32) wrt to reciprocal reduced temperature, which can be compared with Eq. (84) in Span et
613 // al. (2000).
614 // First sum term
615 Evaluation s1 = 0.0;
616 for (int i = 0; i < 7; ++i) {
617 s1 += t_[i] * N_[i] * pow(rho_red, d_[i]) * pow(T_red, t_[i]-1);
618 }
619
620 // Second sum term
621 Evaluation s2 = 0.0;
622 for (int i = 7; i < 9; ++i) {
623 s2 += t_[i] * N_[i] * pow(T_red, t_[i]-1) * pow(rho_red, d_[i]) * exp(-pow(rho_red, p_[i-7]));
624 }
625
626 // Third, and last, sum term
627 Evaluation s3 = 0.0;
628 for (int i = 9; i < 14; ++i) {
629 s3 += N_[i] * pow(T_red, t_[i]-1) * pow(rho_red, d_[i]) *
630 exp(phi_[i-9] * pow(rho_red - D_[i-9], 2) + beta_[i-9] * pow(T_red - gamma_[i-9], 2)) *
631 (t_[i] + 2 * beta_[i-9] * T_red * (T_red - gamma_[i-9]));
632 }
633
634 // Return total sum
635 Evaluation s = s1 + s2 + s3;
636 return s;
637 }
638
645 template <class Evaluation>
646 static Evaluation secDerivResHelmholtzWrtRecipRedTemp(const Evaluation& T_red, const Evaluation& rho_red)
647 {
648 // Second derivative of Eq. (32) wrt to reciprocal reduced temperature, which can be compared with Eq. (85) in
649 // Span et al. (2000).
650 // First sum term
651 Evaluation s1 = 0.0;
652 for (int i = 0; i < 7; ++i) {
653 s1 += t_[i] * (t_[i] - 1) * N_[i] * pow(rho_red, d_[i]) * pow(T_red, t_[i]-2);
654 }
655
656 // Second sum term
657 Evaluation s2 = 0.0;
658 for (int i = 7; i < 9; ++i) {
659 s2 += t_[i] * (t_[i] - 1) * N_[i] * pow(T_red, t_[i]-2) * pow(rho_red, d_[i]) * exp(-pow(rho_red, p_[i-7]));
660 }
661
662 // Third, and last, sum term
663 Evaluation s3 = 0.0;
664 for (int i = 9; i < 14; ++i) {
665 s3 += N_[i] * pow(T_red, t_[i]-2) * pow(rho_red, d_[i]) *
666 exp(phi_[i-9] * pow(rho_red - D_[i-9], 2) + beta_[i-9] * pow(T_red - gamma_[i-9], 2)) *
667 (pow(t_[i] + 2 * beta_[i-9] * T_red * (T_red - gamma_[i-9]), 2)
668 - t_[i] + 2 * beta_[i-9] * pow(T_red, 2));
669 }
670
671 // Return total sum
672 Evaluation s = s1 + s2 + s3;
673 return s;
674 }
675
683 template <class Evaluation>
684 static Evaluation secDerivResHelmholtzWrtRecipRedTempAndRedRho(const Evaluation& T_red, const Evaluation& rho_red)
685 {
686 // Second derivative of Eq. (32) wrt to reciprocal reduced temperature and reduced density, which can be
687 // compared with Eq. (86) in Span et al. (2000).
688 // First sum term
689 Evaluation s1 = 0.0;
690 for (int i = 0; i < 7; ++i) {
691 s1 += t_[i] * d_[i] * N_[i] * pow(rho_red, d_[i]-1) * pow(T_red, t_[i]-1);
692 }
693
694 // Second sum term
695 Evaluation s2 = 0.0;
696 for (int i = 7; i < 9; ++i) {
697 s2 += t_[i] * N_[i] * pow(T_red, t_[i]-1) * pow(rho_red, d_[i]-1) * exp(-pow(rho_red, p_[i-7]))
698 * (d_[i] - p_[i-7] * pow(rho_red, p_[i-7]));
699 }
700
701 // Third, and last, sum term
702 Evaluation s3 = 0.0;
703 for (int i = 9; i < 14; ++i) {
704 s3 += N_[i] * pow(T_red, t_[i]-1) * pow(rho_red, d_[i]-1) *
705 exp(phi_[i-9] * pow(rho_red - D_[i-9], 2) + beta_[i-9] * pow(T_red - gamma_[i-9], 2)) *
706 (t_[i] + 2 * beta_[i-9] * T_red * (T_red - gamma_[i-9]))
707 * (d_[i] + 2 * phi_[i-9] * rho_red * (rho_red - D_[i-9]));
708 }
709
710 // Return total sum
711 Evaluation s = s1 + s2 + s3;
712 return s;
713 }
714
715private:
716
717 // Parameter values need in the ideal-gas contribution to the reduced Helmholtz free energy given in Table 4
718 static constexpr Scalar a_[7] = {-1.4579856475, 1.888076782, 1.616, -0.4117, -0.792, 0.758, 1.217};
719 static constexpr Scalar b_[5] = {-16.0205159149, -22.6580178006, -60.0090511389, -74.9434303817, -206.9392065168};
720
721 // Parameter values needed in the residual contribution to the reduced Helmholtz free energy given in Table 5.
722 static constexpr Scalar N_[14] = {-6.93643, 0.01, 2.1101, 4.52059, 0.732564, -1.34086, 0.130985, -0.777414,
723 0.351944, -0.0211716, 0.0226312, 0.032187, -0.0231752, 0.0557346};
724 static constexpr Scalar t_[14] = {0.6844, 1.0, 0.989, 0.489, 0.803, 1.1444, 1.409, 1.754, 1.311, 4.187, 5.646,
725 0.791, 7.249, 2.986};
726 static constexpr Scalar d_[14] = {1, 4, 1, 1, 2, 2, 3, 1, 3, 2, 1, 3, 1, 1};
727 static constexpr Scalar p_[2] = {1, 1};
728 static constexpr Scalar phi_[5] = {-1.685, -0.489, -0.103, -2.506, -1.607};
729 static constexpr Scalar beta_[5] = {-0.1710, -0.2245, -0.1304, -0.2785, -0.3967};
730 static constexpr Scalar gamma_[5] = {0.7164, 1.3444, 1.4517, 0.7204, 1.5445};
731 static constexpr Scalar D_[5] = {1.506, 0.156, 1.736, 0.670, 1.662};
732
733#if 0
741 template <class Evaluation>
742 static Evaluation rootFindingObj_(const Evaluation& rho_red, const Evaluation& temperature, const Evaluation& pg)
743 {
744 // Temporary calculations
745 Evaluation T_red = criticalTemperature() / temperature; // reciprocal reduced temp.
746 Evaluation p_MPa = pg / 1.0e6; // Pa --> MPa
747 Scalar R = IdealGas::R;
748 Evaluation rho_cRT = criticalDensity() * R * temperature;
749
750 // Eq. (56) in Span et al. (2000)
751 Evaluation dResHelm_dRedRho = derivResHelmholtzWrtRedRho(T_red, rho_red);
752 Evaluation obj = rho_red * rho_cRT * (1 + rho_red * dResHelm_dRedRho) - p_MPa;
753 return obj;
754 }
755#endif
756};
757
758} // end namespace Opm
759
760#endif
Abstract base class of a pure chemical species.
Relations valid for an ideal gas.
A number of commonly used algebraic functions for the localized OPM automatic differentiation (AD) fr...
Implements a scalar function that depends on two variables and which is sampled on an uniform X-Y gri...
Abstract base class of a pure chemical species.
Definition Component.hpp:44
Properties of pure molecular hydrogen .
Definition H2.hpp:90
static Evaluation secDerivResHelmholtzWrtRedRho(const Evaluation &T_red, const Evaluation &rho_red)
Second derivative of the residual part of Helmholtz energy wrt.
Definition H2.hpp:571
static std::string name()
A human readable name for the .
Definition H2.hpp:103
static Evaluation secDerivResHelmholtzWrtRecipRedTempAndRedRho(const Evaluation &T_red, const Evaluation &rho_red)
Second derivative of the residual part of Helmholtz energy first wrt.
Definition H2.hpp:684
static Evaluation vaporPressure(Evaluation temperature)
The vapor pressure in of pure molecular hydrogen at a given temperature.
Definition H2.hpp:166
static Evaluation residualPartHelmholtz(const Evaluation &T_red, const Evaluation &rho_red)
The residual part of Helmholtz energy.
Definition H2.hpp:501
static Scalar criticalTemperature()
Returns the critical temperature of molecular hydrogen.
Definition H2.hpp:115
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition H2.hpp:225
static Scalar tripleDensity()
Returns the density of molecular hydrogen's triple point.
Definition H2.hpp:145
static Scalar acentricFactor()
Acentric factor of .
Definition H2.hpp:156
static Scalar criticalDensity()
Returns the critical density of molecular hydrogen.
Definition H2.hpp:127
static const Evaluation gasIsochoricHeatCapacity(Evaluation temperature, Evaluation pressure)
Specific isochoric heat capacity of pure hydrogen gas.
Definition H2.hpp:394
static Evaluation secDerivResHelmholtzWrtRecipRedTemp(const Evaluation &T_red, const Evaluation &rho_red)
Second derivative of the residual part of Helmholtz energy wrt.
Definition H2.hpp:646
static Scalar criticalPressure()
Returns the critical pressure of molecular hydrogen.
Definition H2.hpp:121
static Evaluation gasPressure(Evaluation temperature, Evaluation density)
The pressure of gaseous in at a given density and temperature.
Definition H2.hpp:235
static Scalar tripleTemperature()
Returns the temperature at molecular hydrogen's triple point.
Definition H2.hpp:133
static Evaluation idealGasPartHelmholtz(const Evaluation &T_red, const Evaluation &rho_red)
The ideal-gas part of Helmholtz energy.
Definition H2.hpp:431
static Evaluation reducedMolarDensity(const Evaluation &temperature, const Evaluation &pg, bool extrapolate=false)
Calculate reduced density (rho/rho_crit) from pressure and temperature.
Definition H2.hpp:417
static Evaluation derivResHelmholtzWrtRedRho(const Evaluation &T_red, const Evaluation &rho_red)
Derivative of the residual part of Helmholtz energy wrt.
Definition H2.hpp:535
static Evaluation gasInternalEnergy(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate=false)
Specific internal energy of H2 [J/kg].
Definition H2.hpp:249
static const Evaluation gasHeatCapacity(Evaluation temperature, Evaluation pressure)
Specific isobaric heat capacity of pure hydrogen gas.
Definition H2.hpp:361
static Scalar criticalVolume()
Critical volume of [m2/kmol].
Definition H2.hpp:151
static Evaluation gasViscosity(const Evaluation &temperature, const Evaluation &pressure, bool extrapolate=false)
The dynamic viscosity of at a given pressure and temperature.
Definition H2.hpp:285
static Evaluation derivResHelmholtzWrtRecipRedTemp(const Evaluation &T_red, const Evaluation &rho_red)
Derivative of the residual part of Helmholtz energy wrt.
Definition H2.hpp:610
static constexpr bool gasIsCompressible()
Returns true if the gas phase is assumed to be compressible.
Definition H2.hpp:219
static Evaluation secDerivIdealHelmholtzWrtRecipRedTemp(const Evaluation &T_red)
Second derivative of the ideal-gas part of Helmholtz energy wrt to reciprocal reduced temperature.
Definition H2.hpp:479
static constexpr Scalar molarMass()
The molar mass in of molecular hydrogen.
Definition H2.hpp:109
static Scalar triplePressure()
Returns the pressure of molecular hydrogen's triple point.
Definition H2.hpp:139
static Evaluation gasDensity(Evaluation temperature, Evaluation pressure, bool extrapolate=false)
The density of at a given pressure and temperature.
Definition H2.hpp:200
static Evaluation gasMolarDensity(Evaluation temperature, Evaluation pressure, bool extrapolate=false)
The molar density of in , depending on pressure and temperature.
Definition H2.hpp:213
static Evaluation derivIdealHelmholtzWrtRecipRedTemp(const Evaluation &T_red)
Derivative of the ideal-gas part of Helmholtz energy wrt to reciprocal reduced temperature.
Definition H2.hpp:454
static const Evaluation gasEnthalpy(Evaluation temperature, Evaluation pressure, bool extrapolate=false)
Specific enthalpy of pure hydrogen gas.
Definition H2.hpp:267
Relations valid for an ideal gas.
Definition IdealGas.hpp:39
static constexpr Scalar R
The ideal gas constant .
Definition IdealGas.hpp:42
Implements a scalar function that depends on two variables and which is sampled on an uniform X-Y gri...
Definition UniformTabulated2DFunction.hpp:68
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition H2.hpp:44
Definition H2.hpp:50
Definition H2.hpp:63