www.mooseframework.org
Water97FluidProperties.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "Water97FluidProperties.h"
11 #include "MathUtils.h"
12 #include "libmesh/utility.h"
13 
14 registerMooseObject("FluidPropertiesApp", Water97FluidProperties);
15 
18 {
20  params.addClassDescription("Fluid properties for water and steam (H2O) using IAPWS-IF97");
21  return params;
22 }
23 
25  : SinglePhaseFluidProperties(parameters),
26  _Mh2o(18.015e-3),
27  _Rw(461.526),
28  _p_critical(22.064e6),
29  _T_critical(647.096),
30  _rho_critical(322.0),
31  _p_triple(611.657),
32  _T_triple(273.16)
33 {
34 }
35 
37 
38 std::string
40 {
41  return "water";
42 }
43 
44 Real
46 {
47  return _Mh2o;
48 }
49 
50 Real
52 {
53  return _p_critical;
54 }
55 
56 Real
58 {
59  return _T_critical;
60 }
61 
62 Real
64 {
65  return _rho_critical;
66 }
67 
68 Real
70 {
71  return _p_triple;
72 }
73 
74 Real
76 {
77  return _T_triple;
78 }
79 
80 Real
82 {
84 }
85 
86 ADReal
88 {
90 }
91 
92 void
94  Real pressure, Real temperature, Real & rho, Real & drho_dp, Real & drho_dT) const
95 {
96  rho_from_p_T_template(pressure, temperature, rho, drho_dp, drho_dT);
97 }
98 
99 void
101  const ADReal & temperature,
102  ADReal & rho,
103  ADReal & drho_dp,
104  ADReal & drho_dT) const
105 {
106  rho_from_p_T_template(pressure, temperature, rho, drho_dp, drho_dT);
107 }
108 
109 Real
111 {
113 }
114 
115 ADReal
117 {
119 }
120 
121 void
123  Real pressure, Real temperature, Real & v, Real & dv_dp, Real & dv_dT) const
124 {
125  v_from_p_T_template(pressure, temperature, v, dv_dp, dv_dT);
126 }
127 
128 void
130  const ADReal & temperature,
131  ADReal & v,
132  ADReal & dv_dp,
133  ADReal & dv_dT) const
134 {
135  v_from_p_T_template(pressure, temperature, v, dv_dp, dv_dT);
136 }
137 
138 Real
139 Water97FluidProperties::p_from_v_e(const Real v, const Real e) const
140 {
141  return p_from_v_e_template(v, e);
142 }
143 
144 ADReal
146 {
147  return p_from_v_e_template(v, e);
148 }
149 
150 Real
152 {
153  return e_from_p_rho_template(p, rho);
154 }
155 
156 ADReal
158 {
159  return e_from_p_rho_template(p, rho);
160 }
161 
162 void
164  const Real p, const Real rho, Real & e, Real & de_dp, Real & de_drho) const
165 {
166  e_from_p_rho_template(p, rho, e, de_dp, de_drho);
167 }
168 
169 void
171  const ADReal & p, const ADReal & rho, ADReal & e, ADReal & de_dp, ADReal & de_drho) const
172 {
173  e_from_p_rho_template(p, rho, e, de_dp, de_drho);
174 }
175 
176 Real
178 {
180 }
181 
182 ADReal
184 {
186 }
187 
188 void
190  Real pressure, Real temperature, Real & e, Real & de_dp, Real & de_dT) const
191 {
192  e_from_p_T_template(pressure, temperature, e, de_dp, de_dT);
193 }
194 
195 void
197  const ADReal & temperature,
198  ADReal & e,
199  ADReal & de_dp,
200  ADReal & de_dT) const
201 {
202  e_from_p_T_template(pressure, temperature, e, de_dp, de_dT);
203 }
204 
205 ADReal
207 {
208  const auto [p, T] = p_T_from_v_h(v, h);
209  return e_from_p_T(p, T);
210 }
211 
212 Real
213 Water97FluidProperties::T_from_v_e(const Real v, const Real e) const
214 {
215  return p_T_from_v_e(v, e).second;
216 }
217 
218 ADReal
220 {
221  return p_T_from_v_e(v, e).second;
222 }
223 
224 ADReal
226 {
227  const auto [p, T] = p_T_from_v_e(v, e);
228  return c_from_p_T(p, T);
229 }
230 
231 Real
232 Water97FluidProperties::c_from_p_T(const Real p, const Real T) const
233 {
234  return c_from_p_T_template(p, T);
235 }
236 
237 ADReal
239 {
240  return c_from_p_T_template(p, T);
241 }
242 
243 Real
245 {
247 }
248 
249 ADReal
251 {
253 }
254 
255 ADReal
257 {
258  const auto [p, T] = p_T_from_v_e(v, e);
259  return cp_from_p_T(p, T);
260 }
261 
262 Real
264 {
266 }
267 
268 ADReal
270 {
272 }
273 
274 ADReal
276 {
277  const auto [p, T] = p_T_from_v_e(v, e);
278  return cv_from_p_T(p, T);
279 }
280 
281 Real
283 {
285 }
286 
287 ADReal
289 {
291 }
292 
293 void
295  Real pressure, Real temperature, Real & mu, Real & dmu_dp, Real & dmu_dT) const
296 {
297  Real rho, drho_dp, drho_dT;
298  this->rho_from_p_T(pressure, temperature, rho, drho_dp, drho_dT);
299  Real dmu_drho;
300  this->mu_from_rho_T(rho, temperature, drho_dT, mu, dmu_drho, dmu_dT);
301  dmu_dp = dmu_drho * drho_dp;
302 }
303 
304 Real
306 {
308 }
309 
310 void
312  Real temperature,
313  Real ddensity_dT,
314  Real & mu,
315  Real & dmu_drho,
316  Real & dmu_dT) const
317 {
318  const Real mu_star = 1.0e-6;
319  const Real rhobar = density / _rho_critical;
320  const Real Tbar = temperature / _T_critical;
321  const Real drhobar_drho = 1.0 / _rho_critical;
322  const Real dTbar_dT = 1.0 / _T_critical;
323 
324  // Limit of zero density. Derivative wrt rho is 0
325  Real sum0 = 0.0, dsum0_dTbar = 0.0;
326  for (std::size_t i = 0; i < _mu_H0.size(); ++i)
327  {
328  sum0 += _mu_H0[i] / MathUtils::pow(Tbar, i);
329  dsum0_dTbar -= i * _mu_H0[i] / MathUtils::pow(Tbar, i + 1);
330  }
331 
332  const Real mu0 = 100.0 * std::sqrt(Tbar) / sum0;
333  const Real dmu0_dTbar =
334  50.0 / std::sqrt(Tbar) / sum0 - 100.0 * std::sqrt(Tbar) * dsum0_dTbar / sum0 / sum0;
335 
336  // Residual component due to finite density
337  Real sum1 = 0.0, dsum1_drhobar = 0.0, dsum1_dTbar = 0.0;
338  for (unsigned int i = 0; i < 6; ++i)
339  {
340  const Real fact = MathUtils::pow(1.0 / Tbar - 1.0, i);
341  const Real dfact_dTbar = i * MathUtils::pow(1.0 / Tbar - 1.0, i - 1) / Tbar / Tbar;
342 
343  for (unsigned int j = 0; j < 7; ++j)
344  {
345  sum1 += fact * _mu_Hij[i][j] * MathUtils::pow(rhobar - 1.0, j);
346  dsum1_dTbar -= dfact_dTbar * _mu_Hij[i][j] * MathUtils::pow(rhobar - 1.0, j);
347  dsum1_drhobar += j * fact * _mu_Hij[i][j] * MathUtils::pow(rhobar - 1.0, j - 1);
348  }
349  }
350 
351  const Real mu1 = std::exp(rhobar * sum1);
352  const Real dmu1_drhobar = (sum1 + rhobar * dsum1_drhobar) * mu1;
353  const Real dmu1_dTbar = (rhobar * dsum1_dTbar) * mu1;
354 
355  // Viscosity and its derivatives are then
356  mu = mu_star * mu0 * mu1;
357  dmu_drho = mu_star * mu0 * dmu1_drhobar * drhobar_drho;
358  dmu_dT = mu_star * (dmu0_dTbar * mu1 + mu0 * dmu1_dTbar) * dTbar_dT + dmu_drho * ddensity_dT;
359 }
360 
361 ADReal
363 {
364  const auto [rho, T] = rho_T_from_v_e(v, e);
365  return mu_from_rho_T_template(rho, T);
366 }
367 
368 void
370  Real temperature,
371  Real & rho,
372  Real & mu) const
373 {
374  rho = this->rho_from_p_T(pressure, temperature);
375  mu = this->mu_from_rho_T(rho, temperature);
376 }
377 
378 void
380  Real temperature,
381  Real & rho,
382  Real & drho_dp,
383  Real & drho_dT,
384  Real & mu,
385  Real & dmu_dp,
386  Real & dmu_dT) const
387 {
388  this->rho_from_p_T(pressure, temperature, rho, drho_dp, drho_dT);
389  Real dmu_drho;
390  this->mu_from_rho_T(rho, temperature, drho_dT, mu, dmu_drho, dmu_dT);
391  dmu_dp = dmu_drho * drho_dp;
392 }
393 
394 Real
396 {
398 }
399 
400 ADReal
402 {
404 }
405 
406 void
407 Water97FluidProperties::k_from_p_T(Real, Real, Real &, Real &, Real &) const
408 {
409  mooseError("k_from_p_T() is not implemented.");
410 }
411 
412 Real
414 {
416 }
417 
418 Real
420 {
421  return k_from_v_e_template(v, e);
422 }
423 
424 ADReal
426 {
427  return k_from_v_e_template(v, e);
428 }
429 
430 Real
431 Water97FluidProperties::s_from_p_T(Real pressure, Real temperature) const
432 {
433  return s_from_p_T_template(pressure, temperature);
434 }
435 
436 ADReal
437 Water97FluidProperties::s_from_p_T(const ADReal & pressure, const ADReal & temperature) const
438 {
439  return s_from_p_T_template(pressure, temperature);
440 }
441 
442 void
443 Water97FluidProperties::s_from_p_T(
444  const Real pressure, const Real temperature, Real & s, Real & ds_dp, Real & ds_dT) const
445 {
446  s_from_p_T_template(pressure, temperature, s, ds_dp, ds_dT);
447 }
448 
449 void
450 Water97FluidProperties::s_from_p_T(const ADReal & pressure,
451  const ADReal & temperature,
452  ADReal & s,
453  ADReal & ds_dp,
454  ADReal & ds_dT) const
455 {
456  s_from_p_T_template(pressure, temperature, s, ds_dp, ds_dT);
457 }
458 
459 Real
461 {
463 }
464 
465 ADReal
467 {
469 }
470 
471 void
473  Real pressure, Real temperature, Real & h, Real & dh_dp, Real & dh_dT) const
474 {
475  h_from_p_T_template(pressure, temperature, h, dh_dp, dh_dT);
476 }
477 
478 void
480  const ADReal & temperature,
481  ADReal & h,
482  ADReal & dh_dp,
483  ADReal & dh_dT) const
484 {
485  h_from_p_T_template(pressure, temperature, h, dh_dp, dh_dT);
486 }
487 
488 Real
490 {
491  // Check whether the input temperature is within the region of validity of this equation.
492  // Valid for 273.15 K <= t <= 647.096 K
493  if (temperature < 273.15 || temperature > _T_critical)
494  mooseException(name(),
495  ": vaporPressure(): Temperature ",
496  temperature,
497  " is outside range 273.15 K <= T "
498  "<= 647.096 K");
499 
500  Real theta, theta2, a, b, c, p;
501  theta = temperature + _n4[8] / (temperature - _n4[9]);
502  theta2 = theta * theta;
503  a = theta2 + _n4[0] * theta + _n4[1];
504  b = _n4[2] * theta2 + _n4[3] * theta + _n4[4];
505  c = _n4[5] * theta2 + _n4[6] * theta + _n4[7];
506  p = Utility::pow<4>(2.0 * c / (-b + std::sqrt(b * b - 4.0 * a * c)));
507 
508  return p * 1.e6;
509 }
510 
513 {
514  // Check whether the input pressure is within the region of validity of this equation.
515  // Valid for 611.213 Pa <= p <= 22.064 MPa
516  if (pressure.value() < 611.23 || pressure.value() > _p_critical)
517  mooseException(name() + ": vaporTemperature(): Pressure ",
518  pressure.value(),
519  " is outside range 611.213 Pa <= p <= 22.064 MPa");
520 
521  const FPDualReal beta = std::pow(pressure / 1.e6, 0.25);
522  const FPDualReal beta2 = beta * beta;
523  const FPDualReal e = beta2 + _n4[2] * beta + _n4[5];
524  const FPDualReal f = _n4[0] * beta2 + _n4[3] * beta + _n4[6];
525  const FPDualReal g = _n4[1] * beta2 + _n4[4] * beta + _n4[7];
526  const FPDualReal d = 2.0 * g / (-f - std::sqrt(f * f - 4.0 * e * g));
527 
528  return (_n4[9] + d - std::sqrt((_n4[9] + d) * (_n4[9] + d) - 4.0 * (_n4[8] + _n4[9] * d))) / 2.0;
529 }
530 
531 Real
533 {
534  const FPDualReal p = pressure;
535 
536  return vaporTemperature_ad(p).value();
537 }
538 
539 void
540 Water97FluidProperties::vaporTemperature(Real pressure, Real & Tsat, Real & dTsat_dp) const
541 {
543  Moose::derivInsert(p.derivatives(), 0, 1.0);
544 
546 
547  Tsat = T.value();
548  dTsat_dp = T.derivatives()[0];
549 }
550 
551 Real
553 {
554  // Check whether the input temperature is within the region of validity of this equation.
555  // Valid for 623.15 K <= t <= 863.15 K
556  if (temperature < 623.15 || temperature > 863.15)
557  mooseException(name(),
558  ": b23p(): Temperature ",
559  temperature,
560  " is outside range of 623.15 K <= T <= 863.15 K");
561 
562  return (_n23[0] + _n23[1] * temperature + _n23[2] * temperature * temperature) * 1.e6;
563 }
564 
565 Real
567 {
568  // Check whether the input pressure is within the region of validity of this equation.
569  // Valid for 16.529 MPa <= p <= 100 MPa
570  if (pressure < 16.529e6 || pressure > 100.0e6)
571  mooseException(
572  name(), ": b23T(): Pressure ", pressure, "is outside range 16.529 MPa <= p <= 100 MPa");
573 
574  return _n23[3] + std::sqrt((pressure / 1.e6 - _n23[4]) / _n23[2]);
575 }
576 
577 unsigned int
579 {
580  unsigned int region;
581 
582  // Need to calculate enthalpies at the boundaries to delineate regions
583  Real p273 = vaporPressure(273.15);
584  Real p623 = vaporPressure(623.15);
585 
586  if (pressure >= p273 && pressure <= p623)
587  {
588  if (enthalpy >= h_from_p_T(pressure, 273.15) &&
590  region = 1;
591  else if (enthalpy > h_from_p_T(pressure, vaporTemperature(pressure)) &&
592  enthalpy <= h_from_p_T(pressure, 1073.15))
593  region = 2;
594  else if (enthalpy > h_from_p_T(pressure, 1073.15) && enthalpy <= h_from_p_T(pressure, 2273.15))
595  region = 5;
596  else
597  mooseException("Enthalpy ", enthalpy, " is out of range in ", name(), ": inRegionPH()");
598  }
599  else if (pressure > p623 && pressure <= 50.0e6)
600  {
601  if (enthalpy >= h_from_p_T(pressure, 273.15) && enthalpy <= h_from_p_T(pressure, 623.15))
602  region = 1;
603  else if (enthalpy > h_from_p_T(pressure, 623.15) &&
604  enthalpy <= h_from_p_T(pressure, b23T(pressure)))
605  region = 3;
606  else if (enthalpy > h_from_p_T(pressure, b23T(pressure)) &&
607  enthalpy <= h_from_p_T(pressure, 1073.15))
608  region = 2;
609  else if (enthalpy > h_from_p_T(pressure, 1073.15) && enthalpy <= h_from_p_T(pressure, 2273.15))
610  region = 5;
611  else
612  mooseException("Enthalpy ", enthalpy, " is out of range in ", name(), ": inRegionPH()");
613  }
614  else if (pressure > 50.0e6 && pressure <= 100.0e6)
615  {
616  if (enthalpy >= h_from_p_T(pressure, 273.15) && enthalpy <= h_from_p_T(pressure, 623.15))
617  region = 1;
618  else if (enthalpy > h_from_p_T(pressure, 623.15) &&
619  enthalpy <= h_from_p_T(pressure, b23T(pressure)))
620  region = 3;
621  else if (enthalpy > h_from_p_T(pressure, b23T(pressure)) &&
622  enthalpy <= h_from_p_T(pressure, 1073.15))
623  region = 2;
624  else
625  mooseException("Enthalpy ", enthalpy, " is out of range in ", name(), ": inRegionPH()");
626  }
627  else
628  mooseException("Pressure ", pressure, " is out of range in ", name(), ": inRegionPH()");
629 
630  return region;
631 }
632 
633 unsigned int
635 {
636  unsigned int subregion;
637 
638  if (pressure <= 4.0e6)
639  subregion = 1;
640  else if (pressure > 4.0e6 && pressure < 6.5467e6)
641  subregion = 2;
642  else
643  {
644  if (enthalpy >= b2bc(pressure))
645  subregion = 2;
646  else
647  subregion = 3;
648  }
649 
650  return subregion;
651 }
652 
653 unsigned int
655 {
656  unsigned int subregion;
657 
658  if (enthalpy <= b3ab(pressure))
659  subregion = 1;
660  else
661  subregion = 2;
662 
663  return subregion;
664 }
665 
666 Real
668 {
669  // Check whether the input pressure is within the region of validity of this equation.
670  if (pressure < 6.5467e6 || pressure > 100.0e6)
671  mooseException(
672  name(), ": b2bc(): Pressure ", pressure, " is outside range of 6.5467 MPa <= p <= 100 MPa");
673 
674  Real pi = pressure / 1.0e6;
675 
676  return (0.26526571908428e4 + std::sqrt((pi - 0.45257578905948e1) / 0.12809002730136e-3)) * 1.0e3;
677 }
678 
679 Real
681 {
682  const FPDualReal p = pressure;
683  const FPDualReal h = enthalpy;
684 
685  return T_from_p_h_ad(p, h).value();
686 }
687 
688 void
690  Real pressure, Real enthalpy, Real & temperature, Real & dT_dp, Real & dT_dh) const
691 {
693  Moose::derivInsert(p.derivatives(), 0, 1.0);
694  FPDualReal h = enthalpy;
695  Moose::derivInsert(h.derivatives(), 1, 1.0);
696 
697  const FPDualReal T = T_from_p_h_ad(p, h);
698 
699  temperature = T.value();
700  dT_dp = T.derivatives()[0];
701  dT_dh = T.derivatives()[1];
702 }
703 
706  const FPDualReal & enthalpy) const
707 {
708  FPDualReal temperature = 0.0;
709 
710  // Determine which region the point is in
711  const unsigned int region = inRegionPH(pressure.value(), enthalpy.value());
712 
713  switch (region)
714  {
715  case 1:
717  break;
718 
719  case 2:
720  {
721  // First, determine which subregion the point is in:
722  const unsigned int subregion = subregion2ph(pressure.value(), enthalpy.value());
723 
724  if (subregion == 1)
726  else if (subregion == 2)
728  else
730  break;
731  }
732 
733  case 3:
734  {
735  // First, determine which subregion the point is in:
736  const unsigned int subregion = subregion3ph(pressure.value(), enthalpy.value());
737 
738  if (subregion == 1)
740  else
742  break;
743  }
744 
745  case 5:
746  mooseError("temperature_from_ph() not implemented for region 5");
747  break;
748 
749  default:
750  mooseError("inRegionPH() has given an incorrect region");
751  }
752 
753  return temperature;
754 }
755 
758  const FPDualReal & enthalpy) const
759 {
760  const FPDualReal pi = pressure / 1.0e6;
761  const FPDualReal eta = enthalpy / 2500.0e3;
762  FPDualReal sum = 0.0;
763 
764  for (std::size_t i = 0; i < _nph1.size(); ++i)
765  sum += _nph1[i] * std::pow(pi, _Iph1[i]) * std::pow(eta + 1.0, _Jph1[i]);
766 
767  return sum;
768 }
769 
772  const FPDualReal & enthalpy) const
773 {
774  const FPDualReal pi = pressure / 1.0e6;
775  const FPDualReal eta = enthalpy / 2000.0e3;
776  FPDualReal sum = 0.0;
777 
778  // Factor out the negative in std::pow(eta - 2.1, _Jph2a[i]) to avoid fpe in dbg (see #13163)
779  const Real sgn = MathUtils::sign(eta.value() - 2.1);
780 
781  for (std::size_t i = 0; i < _nph2a.size(); ++i)
782  sum += _nph2a[i] * std::pow(pi, _Iph2a[i]) * std::pow(std::abs(eta - 2.1), _Jph2a[i]) *
783  std::pow(sgn, _Jph2a[i]);
784 
785  return sum;
786 }
787 
790  const FPDualReal & enthalpy) const
791 {
792  const FPDualReal pi = pressure / 1.0e6;
793  const FPDualReal eta = enthalpy / 2000.0e3;
794  FPDualReal sum = 0.0;
795 
796  // Factor out the negatives in std::pow(pi - 2.0, _Iph2b[i])* std::pow(eta - 2.6, _Jph2b[i])
797  // to avoid fpe in dbg (see #13163)
798  const Real sgn0 = MathUtils::sign(pi.value() - 2.0);
799  const Real sgn1 = MathUtils::sign(eta.value() - 2.6);
800 
801  for (std::size_t i = 0; i < _nph2b.size(); ++i)
802  sum += _nph2b[i] * std::pow(std::abs(pi - 2.0), _Iph2b[i]) * std::pow(sgn0, _Iph2b[i]) *
803  std::pow(std::abs(eta - 2.6), _Jph2b[i]) * std::pow(sgn1, _Jph2b[i]);
804 
805  return sum;
806 }
807 
810  const FPDualReal & enthalpy) const
811 {
812  const FPDualReal pi = pressure / 1.0e6;
813  const FPDualReal eta = enthalpy / 2000.0e3;
814  FPDualReal sum = 0.0;
815 
816  // Factor out the negative in std::pow(eta - 1.8, _Jph2c[i]) to avoid fpe in dbg (see #13163)
817  const Real sgn = MathUtils::sign(eta.value() - 1.8);
818 
819  for (std::size_t i = 0; i < _nph2c.size(); ++i)
820  sum += _nph2c[i] * std::pow(pi + 25.0, _Iph2c[i]) * std::pow(std::abs(eta - 1.8), _Jph2c[i]) *
821  std::pow(sgn, _Jph2c[i]);
822 
823  return sum;
824 }
825 
826 Real
828 {
829  // Check whether the input pressure is within the region of validity of this equation.
830  if (pressure < b23p(623.15) || pressure > 100.0e6)
831  mooseException(
832  name(), ": b3ab(): Pressure ", pressure, "is outside range of 16.529 MPa <= p <= 100 MPa");
833 
834  Real pi = pressure / 1.0e6;
835  Real eta = 0.201464004206875e4 + 0.374696550136983e1 * pi - 0.219921901054187e-1 * pi * pi +
836  0.87513168600995e-4 * pi * pi * pi;
837 
838  return eta * 1.0e3;
839 }
840 
843  const FPDualReal & enthalpy) const
844 {
845  const FPDualReal pi = pressure / 100.0e6;
846  const FPDualReal eta = enthalpy / 2300.0e3;
847  FPDualReal sum = 0.0;
848 
849  for (std::size_t i = 0; i < _nph3a.size(); ++i)
850  sum += _nph3a[i] * std::pow(pi + 0.24, _Iph3a[i]) * std::pow(eta - 0.615, _Jph3a[i]);
851 
852  return sum * 760.0;
853 }
854 
857  const FPDualReal & enthalpy) const
858 {
859  const FPDualReal pi = pressure / 100.0e6;
860  const FPDualReal eta = enthalpy / 2800.0e3;
861  FPDualReal sum = 0.0;
862 
863  for (std::size_t i = 0; i < _nph3b.size(); ++i)
864  sum += _nph3b[i] * std::pow(pi + 0.298, _Iph3b[i]) * std::pow(eta - 0.72, _Jph3b[i]);
865 
866  return sum * 860.0;
867 }
868 
869 Real
870 Water97FluidProperties::henryConstant(Real temperature, const std::vector<Real> & coeffs) const
871 {
872  const Real A = coeffs[0];
873  const Real B = coeffs[1];
874  const Real C = coeffs[2];
875 
876  const Real Tr = temperature / 647.096;
877  const Real tau = 1.0 - Tr;
878 
879  const Real lnkh =
880  A / Tr + B * std::pow(tau, 0.355) / Tr + C * std::pow(Tr, -0.41) * std::exp(tau);
881 
882  // The vapor pressure used in this formulation
883  const std::vector<Real> a{
884  -7.85951783, 1.84408259, -11.7866497, 22.6807411, -15.9618719, 1.80122502};
885  const std::vector<Real> b{1.0, 1.5, 3.0, 3.5, 4.0, 7.5};
886  Real sum = 0.0;
887 
888  for (std::size_t i = 0; i < a.size(); ++i)
889  sum += a[i] * std::pow(tau, b[i]);
890 
891  return 22.064e6 * std::exp(sum / Tr) * std::exp(lnkh);
892 }
893 
894 void
896  const std::vector<Real> & coeffs,
897  Real & Kh,
898  Real & dKh_dT) const
899 {
900  const Real A = coeffs[0];
901  const Real B = coeffs[1];
902  const Real C = coeffs[2];
903 
904  const Real pc = 22.064e6;
905  const Real Tc = 647.096;
906 
907  const Real Tr = temperature / Tc;
908  const Real tau = 1.0 - Tr;
909 
910  const Real lnkh =
911  A / Tr + B * std::pow(tau, 0.355) / Tr + C * std::pow(Tr, -0.41) * std::exp(tau);
912  const Real dlnkh_dT =
913  (-A / Tr / Tr - B * std::pow(tau, 0.355) / Tr / Tr - 0.355 * B * std::pow(tau, -0.645) / Tr -
914  0.41 * C * std::pow(Tr, -1.41) * std::exp(tau) - C * std::pow(Tr, -0.41) * std::exp(tau)) /
915  Tc;
916 
917  // The vapor pressure used in this formulation
918  const std::vector<Real> a{
919  -7.85951783, 1.84408259, -11.7866497, 22.6807411, -15.9618719, 1.80122502};
920  const std::vector<Real> b{1.0, 1.5, 3.0, 3.5, 4.0, 7.5};
921  Real sum = 0.0;
922  Real dsum = 0.0;
923 
924  for (std::size_t i = 0; i < a.size(); ++i)
925  {
926  sum += a[i] * std::pow(tau, b[i]);
927  dsum += a[i] * b[i] * std::pow(tau, b[i] - 1.0);
928  }
929 
930  const Real p = pc * std::exp(sum / Tr);
931  const Real dp_dT = -p / Tc / Tr * (sum / Tr + dsum);
932 
933  // Henry's constant and its derivative wrt temperature
934  Kh = p * std::exp(lnkh);
935  dKh_dT = (p * dlnkh_dT + dp_dT) * std::exp(lnkh);
936 }
937 
938 DualReal
940  const std::vector<Real> & coeffs) const
941 {
942  const Real T = temperature.value();
943  Real Kh_real = 0.0;
944  Real dKh_dT_real = 0.0;
945  henryConstant(T, coeffs, Kh_real, dKh_dT_real);
946 
947  DualReal Kh = Kh_real;
948  Kh.derivatives() = temperature.derivatives() * dKh_dT_real;
949 
950  return Kh;
951 }
952 
953 unsigned int
955 {
956  Real pMPa = pressure / 1.0e6;
957  const Real P3cd = 19.00881189173929;
958  unsigned int subregion = 0;
959 
960  if (pMPa > 40.0 && pMPa <= 100.0)
961  {
962  if (temperature <= tempXY(pressure, AB))
963  subregion = 0;
964  else // (temperature > tempXY(pressure, AB))
965  subregion = 1;
966  }
967  else if (pMPa > 25.0 && pMPa <= 40.0)
968  {
969  if (temperature <= tempXY(pressure, CD))
970  subregion = 2;
971  else if (temperature > tempXY(pressure, CD) && temperature <= tempXY(pressure, AB))
972  subregion = 3;
973  else if (temperature > tempXY(pressure, AB) && temperature <= tempXY(pressure, EF))
974  subregion = 4;
975  else // (temperature > tempXY(pressure, EF))
976  subregion = 5;
977  }
978  else if (pMPa > 23.5 && pMPa <= 25.0)
979  {
980  if (temperature <= tempXY(pressure, CD))
981  subregion = 2;
982  else if (temperature > tempXY(pressure, CD) && temperature <= tempXY(pressure, GH))
983  subregion = 6;
984  else if (temperature > tempXY(pressure, GH) && temperature <= tempXY(pressure, EF))
985  subregion = 7;
986  else if (temperature > tempXY(pressure, EF) && temperature <= tempXY(pressure, IJ))
987  subregion = 8;
988  else if (temperature > tempXY(pressure, IJ) && temperature <= tempXY(pressure, JK))
989  subregion = 9;
990  else // (temperature > tempXY(pressure, JK))
991  subregion = 10;
992  }
993  else if (pMPa > 23.0 && pMPa <= 23.5)
994  {
995  if (temperature <= tempXY(pressure, CD))
996  subregion = 2;
997  else if (temperature > tempXY(pressure, CD) && temperature <= tempXY(pressure, GH))
998  subregion = 11;
999  else if (temperature > tempXY(pressure, GH) && temperature <= tempXY(pressure, EF))
1000  subregion = 7;
1001  else if (temperature > tempXY(pressure, EF) && temperature <= tempXY(pressure, IJ))
1002  subregion = 8;
1003  else if (temperature > tempXY(pressure, IJ) && temperature <= tempXY(pressure, JK))
1004  subregion = 9;
1005  else // (temperature > tempXY(pressure, JK))
1006  subregion = 10;
1007  }
1008  else if (pMPa > 22.5 && pMPa <= 23.0)
1009  {
1010  if (temperature <= tempXY(pressure, CD))
1011  subregion = 2;
1012  else if (temperature > tempXY(pressure, CD) && temperature <= tempXY(pressure, GH))
1013  subregion = 11;
1014  else if (temperature > tempXY(pressure, GH) && temperature <= tempXY(pressure, MN))
1015  subregion = 12;
1016  else if (temperature > tempXY(pressure, MN) && temperature <= tempXY(pressure, EF))
1017  subregion = 13;
1018  else if (temperature > tempXY(pressure, EF) && temperature <= tempXY(pressure, OP))
1019  subregion = 14;
1020  else if (temperature > tempXY(pressure, OP) && temperature <= tempXY(pressure, IJ))
1021  subregion = 15;
1022  else if (temperature > tempXY(pressure, IJ) && temperature <= tempXY(pressure, JK))
1023  subregion = 9;
1024  else // (temperature > tempXY(pressure, JK))
1025  subregion = 10;
1026  }
1027  else if (pMPa > vaporPressure(643.15) * 1.0e-6 &&
1028  pMPa <= 22.5) // vaporPressure(643.15) = 21.04 MPa
1029  {
1030  if (temperature <= tempXY(pressure, CD))
1031  subregion = 2;
1032  else if (temperature > tempXY(pressure, CD) && temperature <= tempXY(pressure, QU))
1033  subregion = 16;
1034  else if (temperature > tempXY(pressure, QU) && temperature <= tempXY(pressure, RX))
1035  {
1036  if (pMPa > 22.11 && pMPa <= 22.5)
1037  {
1038  if (temperature <= tempXY(pressure, UV))
1039  subregion = 20;
1040  else if (temperature > tempXY(pressure, UV) && temperature <= tempXY(pressure, EF))
1041  subregion = 21;
1042  else if (temperature > tempXY(pressure, EF) && temperature <= tempXY(pressure, WX))
1043  subregion = 22;
1044  else // (temperature > tempXY(pressure, WX) && temperature <= tempXY(pressure, RX))
1045  subregion = 23;
1046  }
1047  else if (pMPa > 22.064 && pMPa <= 22.11)
1048  {
1049  if (temperature <= tempXY(pressure, UV))
1050  subregion = 20;
1051  else if (temperature > tempXY(pressure, UV) && temperature <= tempXY(pressure, EF))
1052  subregion = 24;
1053  else if (temperature > tempXY(pressure, EF) && temperature <= tempXY(pressure, WX))
1054  subregion = 25;
1055  else // (temperature > tempXY(pressure, WX) && temperature <= tempXY(pressure, RX))
1056  subregion = 23;
1057  }
1058  else if (temperature <= vaporTemperature(pressure))
1059  {
1060  if (pMPa > 21.93161551 && pMPa <= 22.064)
1062  subregion = 20;
1063  else
1064  subregion = 24;
1065  else // (pMPa > vaporPressure(643.15) * 1.0e-6 && pMPa <= 21.93161551)
1066  subregion = 20;
1067  }
1069  {
1070  if (pMPa > 21.90096265 && pMPa <= 22.064)
1071  {
1072  if (temperature <= tempXY(pressure, WX))
1073  subregion = 25;
1074  else
1075  subregion = 23;
1076  }
1077  else
1078  subregion = 23;
1079  }
1080  }
1081  else if (temperature > tempXY(pressure, RX) && temperature <= tempXY(pressure, JK))
1082  subregion = 17;
1083  else
1084  subregion = 10;
1085  }
1086  else if (pMPa > 20.5 &&
1087  pMPa <= vaporPressure(643.15) * 1.0e-6) // vaporPressure(643.15) = 21.04 MPa
1088  {
1089  if (temperature <= tempXY(pressure, CD))
1090  subregion = 2;
1092  subregion = 16;
1094  subregion = 17;
1095  else // (temperature > tempXY(pressure, JK))
1096  subregion = 10;
1097  }
1098  else if (pMPa > P3cd && pMPa <= 20.5) // P3cd = 19.00881189173929
1099  {
1100  if (temperature <= tempXY(pressure, CD))
1101  subregion = 2;
1103  subregion = 18;
1104  else
1105  subregion = 19;
1106  }
1107  else if (pMPa > vaporPressure(623.15) * 1.0e-6 && pMPa <= P3cd)
1108  {
1110  subregion = 2;
1111  else
1112  subregion = 19;
1113  }
1114  else if (pMPa > 22.11 && pMPa <= 22.5)
1115  {
1117  subregion = 20;
1118  else if (temperature > tempXY(pressure, UV) && temperature <= tempXY(pressure, EF))
1119  subregion = 21;
1120  else if (temperature > tempXY(pressure, EF) && temperature <= tempXY(pressure, WX))
1121  subregion = 22;
1122  else // (temperature > tempXY(pressure, WX) && temperature <= tempXY(pressure, RX))
1123  subregion = 23;
1124  }
1125  else if (pMPa > 22.064 && pMPa <= 22.11)
1126  {
1128  subregion = 20;
1129  else if (temperature > tempXY(pressure, UV) && temperature <= tempXY(pressure, EF))
1130  subregion = 24;
1131  else if (temperature > tempXY(pressure, EF) && temperature <= tempXY(pressure, WX))
1132  subregion = 25;
1133  else // (temperature > tempXY(pressure, WX) && temperature <= tempXY(pressure, RX))
1134  subregion = 23;
1135  }
1136  else
1137  mooseError("subregion3(): Shouldn't have got here!");
1138 
1139  return subregion;
1140 }
1141 
1142 unsigned int
1144 {
1145  // Valid for 273.15 K <= T <= 1073.15 K, p <= 100 MPa
1146  // 1073.15 K <= T <= 2273.15 K, p <= 50 Mpa
1147  if (temperature >= 273.15 && temperature <= 1073.15)
1148  {
1149  if (pressure < vaporPressure(273.15) || pressure > 100.0e6)
1150  mooseException("Pressure ", pressure, " is out of range in ", name(), ": inRegion()");
1151  }
1152  else if (temperature > 1073.15 && temperature <= 2273.15)
1153  {
1154  if (pressure < 0.0 || pressure > 50.0e6)
1155  mooseException("Pressure ", pressure, " is out of range in ", name(), ": inRegion()");
1156  }
1157  else
1158  mooseException("Temperature ", temperature, " is out of range in ", name(), ": inRegion()");
1159 
1160  // Determine the phase region that the (P, T) point lies in
1161  unsigned int region;
1162 
1163  if (temperature >= 273.15 && temperature <= 623.15)
1164  {
1165  if (pressure > vaporPressure(temperature) && pressure <= 100.0e6)
1166  region = 1;
1167  else
1168  region = 2;
1169  }
1170  else if (temperature > 623.15 && temperature <= 863.15)
1171  {
1172  if (pressure <= b23p(temperature))
1173  region = 2;
1174  else
1175  region = 3;
1176  }
1177  else if (temperature > 863.15 && temperature <= 1073.15)
1178  region = 2;
1179  else
1180  region = 5;
1181 
1182  return region;
1183 }
T v_from_p_T_template(const T &pressure, const T &temperature) const
const Real _T_critical
Critical temperature (K)
virtual Real mu_from_p_T(Real pressure, Real temperature) const override
static InputParameters validParams()
virtual void rho_mu_from_p_T(Real pressure, Real temperature, Real &rho, Real &mu) const override
Combined methods.
virtual Real triplePointTemperature() const override
Triple point temperature.
unsigned int inRegion(Real pressure, Real temperature) const
Determines the phase region that the given pressure and temperature values lie in.
const std::array< Real, 5 > _n23
Constants for the boundary between regions 2 and 3.
T tempXY(const T &pressure, subregionEnum xy) const
Boundaries between subregions in region 3.
std::pair< T, T > rho_T_from_v_e(const T &v, const T &e) const
Computes the density (first member of the pair) and temperature (second member of the pair) as functi...
Real vaporTemperature(Real pressure) const override
Saturation temperature as a function of pressure.
Real p_from_v_e(Real v, Real e) const override
const Real _p_critical
Critical pressure (Pa)
const std::array< int, 33 > _Jph3b
Real T_from_v_e(Real v, Real e) const override
const std::array< int, 20 > _Iph1
virtual Real e_from_p_T(Real pressure, Real temperature) const override
Real b2bc(Real pressure) const
Boundary between subregions b and c in region 2.
static InputParameters validParams()
virtual std::string fluidName() const override
const std::array< int, 23 > _Jph2c
const std::array< Real, 36 > _nph2a
T k_from_v_e_template(const T &v, const T &e) const
T e_from_p_T_template(const T &pressure, const T &temperature) const
registerMooseObject("FluidPropertiesApp", Water97FluidProperties)
virtual Real criticalDensity() const override
Critical density.
T k_from_p_T_template(const T &pressure, const T &temperature) const
Real b23T(Real pressure) const
Auxillary equation for the boundary between regions 2 and 3.
Water97FluidProperties(const InputParameters &parameters)
virtual ADReal cv_from_v_e(const ADReal &v, const ADReal &e) const override
T c_from_p_T_template(const T &pressure, const T &temperature) const
T rho_from_p_T_template(const T &pressure, const T &temperature) const
DualNumber< Real, DNDerivativeType, true > DualReal
unsigned int subregion3(Real pressure, Real temperature) const
Provides the correct subregion index for a (P,T) point in region 3.
const std::array< int, 36 > _Jph2a
static const std::string density
Definition: NS.h:33
FPDualReal temperature_from_ph1(const FPDualReal &pressure, const FPDualReal &enthalpy) const
Backwards equation T(p, h) in Region 1 Eq.
const std::array< std::array< Real, 7 >, 6 > _mu_Hij
int sgn(T val)
The sign function.
Definition: Numerics.h:39
T h_from_p_T_template(const T &pressure, const T &temperature) const
const std::array< Real, 33 > _nph3b
static const std::string temperature
Definition: NS.h:57
ADRealEigenVector< T, D, asd > sqrt(const ADRealEigenVector< T, D, asd > &)
const std::array< Real, 38 > _nph2b
std::pair< T, T > p_T_from_v_h(const T &v, const T &h) const
Computes the pressure (first member of the pair) and temperature (second member of the pair) as funct...
virtual Real h_from_p_T(Real pressure, Real temperature) const override
FPDualReal vaporTemperature_ad(const FPDualReal &pressure) const
AD version of saturation temperature as a function of pressure (used internally)
virtual Real criticalTemperature() const override
Critical temperature.
const std::array< int, 33 > _Iph3b
virtual Real k_from_p_T(Real pressure, Real temperature) const override
ADReal mu_from_v_e(const ADReal &v, const ADReal &e) const override
unsigned int subregion2ph(Real pressure, Real enthalpy) const
Provides the correct subregion index for a (P,h) point in region 2.
const Real _T_triple
Triple point temperature (K)
DualNumber< Real, DNDerivativeSize< 5 > > FPDualReal
virtual const std::string & name() const
std::pair< T, T > p_T_from_v_e(const T &v, const T &e) const
Computes the pressure (first member of the pair) and temperature (second member of the pair) as funct...
const std::array< int, 23 > _Iph2c
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
FPDualReal temperature_from_ph2a(const FPDualReal &pressure, const FPDualReal &enthalpy) const
Backwards equation T(p, h) in Region 2a Eq.
T cv_from_p_T_template(const T &pressure, const T &temperature) const
const std::array< int, 20 > _Jph1
const std::array< int, 36 > _Iph2a
Real k_from_v_e(Real v, Real e) const override
const std::array< int, 38 > _Jph2b
Real f(Real x)
Test function for Brents method.
const std::array< int, 31 > _Jph3a
FPDualReal temperature_from_ph2b(const FPDualReal &pressure, const FPDualReal &enthalpy) const
Backwards equation T(p, h) in Region 2b Eq.
static const std::string mu
Definition: NS.h:122
const std::array< int, 31 > _Iph3a
virtual Real v_from_p_T(Real pressure, Real temperature) const override
T sign(T x)
FPDualReal temperature_from_ph2c(const FPDualReal &pressure, const FPDualReal &enthalpy) const
Backwards equation T(p, h) in Region 2c Eq.
virtual Real triplePointPressure() const override
Triple point pressure.
Common class for single phase fluid properties.
const std::array< Real, 20 > _nph1
ADReal e_from_v_h(const ADReal &v, const ADReal &h) const override
DualReal ADReal
FPDualReal T_from_p_h_ad(const FPDualReal &pressure, const FPDualReal &enthalpy) const
AD version of backwards equation T(p, h) (used internally) From Revised Release on the IAPWS Industri...
e e e e s T T T T T rho v v T h
unsigned int subregion3ph(Real pressure, Real enthalpy) const
Provides the correct subregion index for a (P,h) point in region 3.
virtual Real T_from_p_h(Real pressure, Real enthalpy) const override
Backwards equation T(p, h) From Revised Release on the IAPWS Industrial Formulation 1997 for the Ther...
T mu_from_rho_T_template(const T &density, const T &temperature) const
virtual Real c_from_p_T(Real pressure, Real temperature) const override
T p_from_v_e_template(const T &v, const T &e) const
const std::array< Real, 31 > _nph3a
T k_from_rho_T_template(const T &density, const T &temperature) const
Real b23p(Real temperature) const
Auxillary equation for the boundary between regions 2 and 3.
const Real _rho_critical
Critical density (kg/m^3)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::string v
Definition: NS.h:82
const std::array< int, 38 > _Iph2b
virtual Real e_from_p_rho(Real p, Real rho) const override
virtual Real criticalPressure() const override
Critical pressure.
Real b3ab(Real pressure) const
Boundary between subregions a and b in region 3.
Water (H2O) fluid properties as a function of pressure (Pa) and temperature (K) from IAPWS-IF97: Revi...
static const std::string pressure
Definition: NS.h:56
virtual Real vaporPressure(Real temperature) const override
Vapor pressure.
Real henryConstant(Real temperature, const std::vector< Real > &coeffs) const
IAPWS formulation of Henry&#39;s law constant for dissolution in water From Guidelines on the Henry&#39;s con...
FPDualReal temperature_from_ph3b(const FPDualReal &pressure, const FPDualReal &enthalpy) const
Backwards equation T(p, h) in Region 3b Eq.
FPDualReal temperature_from_ph3a(const FPDualReal &pressure, const FPDualReal &enthalpy) const
Backwards equation T(p, h) in Region 3a Eq.
void mooseError(Args &&... args) const
unsigned int inRegionPH(Real pressure, Real enthalpy) const
Determines the phase region that the given pressure and enthaply values lie in.
void addClassDescription(const std::string &doc_string)
virtual Real k_from_rho_T(Real density, Real temperature) const override
T cp_from_p_T_template(const T &pressure, const T &temperature) const
void derivInsert(NumberArray< N, Real > &derivs, dof_id_type index, Real value)
virtual Real rho_from_p_T(Real pressure, Real temperature) const override
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const Real _p_triple
Triple point pressure (Pa)
T mu_from_p_T_template(const T &pressure, const T &temperature) const
T e_from_p_rho_template(const T &p, const T &rho) const
T pow(T x, int e)
virtual ADReal c_from_v_e(const ADReal &v, const ADReal &e) const override
virtual Real cv_from_p_T(Real pressure, Real temperature) const override
virtual Real mu_from_rho_T(Real density, Real temperature) const override
const std::array< Real, 23 > _nph2c
virtual Real cp_from_p_T(Real pressure, Real temperature) const override
const Real _Mh2o
Water molar mass (kg/mol)
MooseUnits pow(const MooseUnits &, int)
static const std::string C
Definition: NS.h:153
virtual ADReal cp_from_v_e(const ADReal &v, const ADReal &e) const override
virtual Real molarMass() const override
Fluid name.
const Real pi
const std::array< Real, 10 > _n4
Constants for region 4 (the saturation curve up to the critical point)
const std::array< Real, 4 > _mu_H0
Constants from Release on the IAPWS Formulation 2008 for the Viscosity of Ordinary Water Substance...