3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
h2o.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 *****************************************************************************/
24#ifndef DUMUX_H2O_HH
25#define DUMUX_H2O_HH
26
27#include <cmath>
28#include <cassert>
29
33
34#include "iapws/common.hh"
35#include "iapws/region1.hh"
36#include "iapws/region2.hh"
37#include "iapws/region4.hh"
38
42
43namespace Dumux {
44namespace Components {
45
56template <class Scalar>
57class H2O
58: public Components::Base<Scalar, H2O<Scalar> >
59, public Components::Liquid<Scalar, H2O<Scalar> >
60, public Components::Gas<Scalar, H2O<Scalar> >
61{
62
67
68 // specific gas constant of water
69 static constexpr Scalar Rs = Common::Rs;
70
71public:
75 static std::string name()
76 { return "H2O"; }
77
81 static constexpr Scalar molarMass()
82 { return Common::molarMass; }
83
87 static constexpr Scalar acentricFactor()
88 { return Common::acentricFactor; }
89
93 static constexpr Scalar criticalTemperature()
95
99 static constexpr Scalar criticalPressure()
100 { return Common::criticalPressure; }
101
105 static constexpr Scalar criticalMolarVolume()
107
111 static constexpr Scalar tripleTemperature()
112 { return Common::tripleTemperature; }
113
117 static constexpr Scalar triplePressure()
118 { return Common::triplePressure; }
119
132 {
133 using std::min;
134 using std::max;
135 T = min(T, criticalTemperature());
136 T = max(T,tripleTemperature());
137
139 }
140
153 {
154 using std::min;
155 using std::max;
158
160 }
161
175 {
177
178 // regularization
179 if (pressure < triplePressure() - 100) {
180 // We assume an ideal gas for low pressures to avoid the
181 // 0/0 for the gas enthalpy at very low pressures. The
182 // enthalpy of an ideal gas does not exhibit any
183 // dependence on pressure, so we can just return the
184 // specific enthalpy at the point of regularization, i.e.
185 // the triple pressure - 100Pa
186 return enthalpyRegion2_(temperature, triplePressure() - 100);
187 }
189 if (pressure > pv) {
190 // the pressure is too high, in this case we use the slope
191 // of the enthalpy at the vapor pressure to regularize
192 Scalar dh_dp =
193 Rs*temperature*
195 Region2::dPi_dp(pv)*
197
198 return
199 enthalpyRegion2_(temperature, pv) +
200 (pressure - pv)*dh_dp;
201 }
202
203 return enthalpyRegion2_(temperature, pressure);
204 }
205
219 {
221
222 // regularization
224 if (pressure < pv) {
225 // the pressure is too low, in this case we use the slope
226 // of the enthalpy at the vapor pressure to regularize
227 Scalar dh_dp =
228 Rs * temperature*
230 Region1::dPi_dp(pv)*
232
233 return
234 enthalpyRegion1_(temperature, pv) +
235 (pressure - pv)*dh_dp;
236 }
237
238 return enthalpyRegion1_(temperature, pressure);
239 }
240
254 {
256
257 // regularization
258 if (pressure < triplePressure() - 100) {
259 return heatCap_p_Region2_(temperature, triplePressure() - 100);
260 }
262 if (pressure > pv) {
263 // the pressure is too high, in this case we use the heat
264 // cap at the vapor pressure to regularize
265 return heatCap_p_Region2_(temperature, pv);
266 }
267 return heatCap_p_Region2_(temperature, pressure);
268 }
269
283 {
285
286 // regularization
288 if (pressure < pv) {
289 // the pressure is too low, in this case we use the heat cap at the vapor pressure to regularize
290 return heatCap_p_Region1_(temperature, pv);
291 }
292
293 return heatCap_p_Region1_(temperature, pressure);
294 }
295
309 {
311
312 // regularization
314 if (pressure < pv) {
315 // the pressure is too low, in this case we use the slope
316 // of the internal energy at the vapor pressure to
317 // regularize
318
319 /*
320 // calculate the partial derivative of the internal energy
321 // to the pressure at the vapor pressure.
322 Scalar tau = Region1::tau(temperature);
323 Scalar dGamma_dPi = Region1::dGamma_dPi(temperature, pv);
324 Scalar ddGamma_dTaudPi = Region1::ddGamma_dTaudPi(temperature, pv);
325 Scalar ddGamma_ddPi = Region1::ddGamma_ddPi(temperature, pv);
326 Scalar pi = Region1::pi(pv);
327 Scalar dPi_dp = Region1::dPi_dp(pv);
328 Scalar du_dp =
329 Rs*temperature*
330 (tau*dPi_dp*ddGamma_dTaudPi + dPi_dp*dPi_dp*dGamma_dPi + pi*dPi_dp*ddGamma_ddPi);
331 */
332
333 // use a straight line for extrapolation. use forward
334 // differences to calculate the partial derivative to the
335 // pressure at the vapor pressure
336 static const Scalar eps = 1e-7;
337 Scalar uv = internalEnergyRegion1_(temperature, pv);
338 Scalar uvPEps = internalEnergyRegion1_(temperature, pv + eps);
339 Scalar du_dp = (uvPEps - uv)/eps;
340 return uv + du_dp*(pressure - pv);
341 }
342
343 return internalEnergyRegion1_(temperature, pressure);
344 }
345
358 {
360
361 // regularization
362 if (pressure < triplePressure() - 100) {
363 // We assume an ideal gas for low pressures to avoid the
364 // 0/0 for the internal energy of gas at very low
365 // pressures. The enthalpy of an ideal gas does not
366 // exhibit any dependence on pressure, so we can just
367 // return the specific enthalpy at the point of
368 // regularization, i.e. the triple pressure - 100Pa, and
369 // subtract the work required to change the volume for an
370 // ideal gas.
371 return
372 enthalpyRegion2_(temperature, triplePressure() - 100)
373 -
374 Rs*temperature; // = p*v for an ideal gas!
375 }
377 if (pressure > pv) {
378 // the pressure is too high, in this case we use the slope
379 // of the internal energy at the vapor pressure to
380 // regularize
381
382 /*
383 // calculate the partial derivative of the internal energy
384 // to the pressure at the vapor pressure.
385 Scalar tau = Region2::tau(temperature);
386 Scalar dGamma_dPi = Region2::dGamma_dPi(temperature, pv);
387 Scalar ddGamma_dTaudPi = Region2::ddGamma_dTaudPi(temperature, pv);
388 Scalar ddGamma_ddPi = Region2::ddGamma_ddPi(temperature, pv);
389 Scalar pi = Region2::pi(pv);
390 Scalar dPi_dp = Region2::dPi_dp(pv);
391 Scalar du_dp =
392 Rs*temperature*
393 (tau*dPi_dp*ddGamma_dTaudPi + dPi_dp*dPi_dp*dGamma_dPi + pi*dPi_dp*ddGamma_ddPi);
394
395 // use a straight line for extrapolation
396 Scalar uv = internalEnergyRegion2_(temperature, pv);
397 return uv + du_dp*(pressure - pv);
398 */
399
400 // use a straight line for extrapolation. use backward
401 // differences to calculate the partial derivative to the
402 // pressure at the vapor pressure
403 static const Scalar eps = 1e-7;
404 Scalar uv = internalEnergyRegion2_(temperature, pv);
405 Scalar uvMEps = internalEnergyRegion2_(temperature, pv - eps);
406 Scalar du_dp = (uv - uvMEps)/eps;
407 return uv + du_dp*(pressure - pv);
408 }
409
410 return internalEnergyRegion2_(temperature, pressure);
411 }
412
426 {
427 Region1::checkValidityRange(temperature, pressure, "Heat capacity for a constant volume");
428
429 // regularization
431 if (pressure < pv) {
432 // the pressure is too low, in this case we use the heat cap at the vapor pressure to regularize
433
434 return heatCap_v_Region1_(temperature, pv);
435 }
436
437 return heatCap_v_Region1_(temperature, pressure);
438 }
439
453 {
454 Region2::checkValidityRange(temperature, pressure, "Heat capacity for a constant volume");
455
456 // regularization
457 if (pressure < triplePressure() - 100) {
458 return
459 heatCap_v_Region2_(temperature, triplePressure() - 100);
460 }
462 if (pressure > pv) {
463 return heatCap_v_Region2_(temperature, pv);
464 }
465
466 return heatCap_v_Region2_(temperature, pressure);
467 }
468
472 static constexpr bool gasIsCompressible()
473 { return true; }
474
478 static constexpr bool liquidIsCompressible()
479 { return true; }
480
494 {
496
497 // regularization
498 if (pressure < triplePressure() - 100) {
499 // We assume an ideal gas for low pressures to avoid the
500 // 0/0 for the internal energy and enthalpy.
501 Scalar rho0IAPWS = 1.0/volumeRegion2_(temperature,
502 triplePressure() - 100);
505 triplePressure() - 100);
506 return
507 rho0IAPWS/rho0Id *
510 pressure);
511 }
513 if (pressure > pv) {
514 // the pressure is too high, in this case we use the slope
515 // of the density energy at the vapor pressure to
516 // regularize
517
518 // calculate the partial derivative of the specific volume
519 // to the pressure at the vapor pressure.
520 const Scalar eps = pv*1e-8;
521 Scalar v0 = volumeRegion2_(temperature, pv);
522 Scalar v1 = volumeRegion2_(temperature, pv + eps);
523 Scalar dv_dp = (v1 - v0)/eps;
524 /*
525 Scalar pi = Region2::pi(pv);
526 Scalar dp_dPi = Region2::dp_dPi(pv);
527 Scalar dGamma_dPi = Region2::dGamma_dPi(temperature, pv);
528 Scalar ddGamma_ddPi = Region2::ddGamma_ddPi(temperature, pv);
529
530 Scalar RT = Rs*temperature;
531 Scalar dv_dp =
532 RT/(dp_dPi*pv)
533 *
534 (dGamma_dPi + pi*ddGamma_ddPi - v0*dp_dPi/RT);
535 */
536
537 // calculate the partial derivative of the density to the
538 // pressure at vapor pressure
539 Scalar drho_dp = - 1/(v0*v0)*dv_dp;
540
541 // use a straight line for extrapolation
542 return 1.0/v0 + (pressure - pv)*drho_dp;
543 }
544
545 return 1.0/volumeRegion2_(temperature, pressure);
546 }
547
557
561 static constexpr bool gasIsIdeal()
562 { return false; }
563
577 {
578 // We use the newton method for this. For the initial value we
579 // assume steam to be an ideal gas
581 Scalar eps = pressure*1e-7;
582
583 Scalar deltaP = pressure*2;
584 using std::abs;
585 for (int i = 0; i < 5 && abs(pressure*1e-9) < abs(deltaP); ++i)
586 {
588
589 Scalar df_dp;
590 df_dp = gasDensity(temperature, pressure + eps);
591 df_dp -= gasDensity(temperature, pressure - eps);
592 df_dp /= 2*eps;
593
594 deltaP = - f/df_dp;
595
596 pressure += deltaP;
597 }
598
599 return pressure;
600 }
601
615 {
617
618 // regularization
620 if (pressure < pv) {
621 // the pressure is too low, in this case we use the slope
622 // of the density at the vapor pressure to regularize
623
624 // calculate the partial derivative of the specific volume
625 // to the pressure at the vapor pressure.
626 const Scalar eps = pv*1e-8;
627 Scalar v0 = volumeRegion1_(temperature, pv);
628 Scalar v1 = volumeRegion1_(temperature, pv + eps);
629 Scalar dv_dp = (v1 - v0)/eps;
630
631 /*
632 Scalar v0 = volumeRegion1_(temperature, pv);
633 Scalar pi = Region1::pi(pv);
634 Scalar dp_dPi = Region1::dp_dPi(pv);
635 Scalar dGamma_dPi = Region1::dGamma_dPi(temperature, pv);
636 Scalar ddGamma_ddPi = Region1::ddGamma_ddPi(temperature, pv);
637
638 Scalar RT = Rs*temperature;
639 Scalar dv_dp =
640 RT/(dp_dPi*pv)
641 *
642 (dGamma_dPi + pi*ddGamma_ddPi - v0*dp_dPi/RT);
643 */
644
645 // calculate the partial derivative of the density to the
646 // pressure at vapor pressure
647 Scalar drho_dp = - 1/(v0*v0)*dv_dp;
648
649 // use a straight line for extrapolation
650 return 1.0/v0 + (pressure - pv)*drho_dp;
651 }
652
653 return 1/volumeRegion1_(temperature, pressure);
654 }
655
665
680 {
681 // We use Brent's method for this, with a pressure range including the surroundings of the
682 // vapor pressure line
683 Scalar minPressure = vaporPressure(temperature)/1.11;
684 Scalar maxPressure = 100e6;
685 const auto residualFunction = [&] (const Scalar pressure) {
687 };
688 try
689 {
690 return findScalarRootBrent(minPressure, maxPressure, residualFunction);
691 }
692 catch (const NumericalProblem& e)
693 {
694 DUNE_THROW(NumericalProblem,
695 "searched for pressure(T=" << temperature << ",rho=" << density
696 <<") in [" << minPressure << ", " << maxPressure << "]: "
697 << e.what());
698 }
699 }
700
717 {
719
721 return Common::viscosity(temperature, rho);
722 }
723
724
736 {
738
740 return Common::viscosity(temperature, rho);
741 }
742
759 {
760 // Thermal conductivity of water is empirically fit.
761 // Evaluating that fitting-function outside the area of validity does not make sense.
762 if ( !( (pressure <= 400e6 && (273.15 <= temperature) && (temperature <= 398.15))
763 || (pressure <= 200e6 && (398.15 < temperature) && (temperature <= 523.15))
764 || (pressure <= 150e6 && (523.15 < temperature) && (temperature <= 673.15))
765 || (pressure <= 100e6 && (673.15 < temperature) && (temperature <= 1073.15)) ))
766 {
767 DUNE_THROW(NumericalProblem,
768 "Evaluating the IAPWS fit function for thermal conductivity outside range of applicability."
769 "(T=" << temperature << ", p=" << pressure << ")");
770 }
771
774 }
775
792 {
793 // Thermal conductivity of water is empirically fit.
794 // Evaluating that fitting-function outside the area of validity does not make sense.
795 if ( !( (pressure <= 400e6 && (273.15 <= temperature) && (temperature <= 398.15))
796 || (pressure <= 200e6 && (398.15 < temperature) && (temperature <= 523.15))
797 || (pressure <= 150e6 && (523.15 < temperature) && (temperature <= 673.15))
798 || (pressure <= 100e6 && (673.15 < temperature) && (temperature <= 1073.15)) ))
799 {
800 DUNE_THROW(NumericalProblem,
801 "Evaluating the IAPWS fit function for thermal conductivity outside range of applicability."
802 "(T=" << temperature << ", p=" << pressure << ")");
803 }
804
807 }
808
809private:
810 // the unregularized specific enthalpy for liquid water
811 static constexpr Scalar enthalpyRegion1_(Scalar temperature, Scalar pressure)
812 {
813 return
816 Rs*temperature;
817 }
818
819 // the unregularized specific isobaric heat capacity
820 static constexpr Scalar heatCap_p_Region1_(Scalar temperature, Scalar pressure)
821 {
822 return
825 Rs;
826 }
827
828 // the unregularized specific isochoric heat capacity
829 static Scalar heatCap_v_Region1_(Scalar temperature, Scalar pressure)
830 {
834
835 return
836 - tau * tau *
838 diff;
839 }
840
841 // the unregularized specific internal energy for liquid water
842 static constexpr Scalar internalEnergyRegion1_(Scalar temperature, Scalar pressure)
843 {
844 return
845 Rs * temperature *
848 }
849
850 // the unregularized specific volume for liquid water
851 static constexpr Scalar volumeRegion1_(Scalar temperature, Scalar pressure)
852 {
853 return
856 Rs * temperature / pressure;
857 }
858
859 // the unregularized specific enthalpy for steam
860 static constexpr Scalar enthalpyRegion2_(Scalar temperature, Scalar pressure)
861 {
862 return
865 Rs*temperature;
866 }
867
868 // the unregularized specific internal energy for steam
869 static constexpr Scalar internalEnergyRegion2_(Scalar temperature, Scalar pressure)
870 {
871 return
872 Rs * temperature *
875 }
876
877 // the unregularized specific isobaric heat capacity
878 static constexpr Scalar heatCap_p_Region2_(Scalar temperature, Scalar pressure)
879 {
880 return
883 Rs;
884 }
885
886 // the unregularized specific isochoric heat capacity
887 static Scalar heatCap_v_Region2_(Scalar temperature, Scalar pressure)
888 {
892 Scalar diff = num * num / (1 - pi * pi * Region2::ddGamma_ddPi(temperature, pressure));
893 return
894 - tau * tau *
896 - diff;
897 }
898
899 // the unregularized specific volume for steam
900 static constexpr Scalar volumeRegion2_(Scalar temperature, Scalar pressure)
901 {
902 return
905 Rs * temperature / pressure;
906 }
907};
908
909template <class Scalar>
910struct IsAqueous<H2O<Scalar>> : public std::true_type {};
911
912} // end namespace Components
913
914namespace FluidSystems::Detail {
915// viscosity according to Reid, R.C.
916template <class Scalar>
918{
919 constexpr Scalar tc = 647.3;
920 using std::max;
921 const Scalar tr = max(temperature/tc, 1e-8);
922
923 const Scalar fp0 = 1.0 + 0.221*(0.96 + 0.1*(tr - 0.7));
924 constexpr Scalar xi = 3.334e-3;
925 using std::pow;
926 using std::exp;
927 const Scalar eta_xi = (0.807*pow(tr, 0.618) - 0.357*exp((-0.449)*tr)
928 + 0.34*exp((-4.058)*tr) + 0.018)*fp0;
929
930 return 1.0e-7*eta_xi/xi;
931}
932
933// viscosity according to Nagel, T. et al.
934template <class Scalar>
936{
937 constexpr Scalar a1 = -4.4189440e-6;
938 constexpr Scalar a2 = 4.6876380e-8;
939 constexpr Scalar a3 = -5.3894310e-12;
940 constexpr Scalar a4 = 3.2028560e-16;
941 constexpr Scalar a5 = 4.9191790e-22;
942
943 return a1 + a2*temperature + a3*temperature*temperature
946}
947} // end FluidSystems::Detail
948namespace FluidSystems {
972template <class Scalar>
974{
975 if (temperature < 480.0)
977
978 else if (temperature > 500.0)
980
981 else // interpolate
982 {
983 const Scalar op = (500.0 - temperature)/20.0;
984
987 }
988}
989} // end namespace FluidSystems
990} // end namespace Dumux
991
992#endif
Some exceptions thrown in DuMux
Implements the equations for region 1 of the IAPWS '97 formulation. See:
Implements the equations for region 2 of the IAPWS '97 formulation. See:
Implements relations common for all regions of the IAPWS '97 formulation. See:
Implements the equations for region 4 of the IAPWS '97 formulation.
Interface for components that have a gas state.
Interface for components that have a liquid state.
Relations valid for an ideal gas.
Root finding algorithms for scalar functions.
Scalar findScalarRootBrent(Scalar a, Scalar b, const ResFunc &residual, const Scalar tol=1e-13, const int maxIter=200)
Brent's root finding algorithm for scalar functions.
Definition: findscalarroot.hh:111
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Scalar h2oGasViscosityInMixture(Scalar temperature, Scalar pressure)
The dynamic viscosity of steam in a gas mixture.
Definition: h2o.hh:973
Scalar viscosityNagel_(Scalar temperature)
Definition: h2o.hh:935
Scalar viscosityReid_(Scalar temperature)
Definition: h2o.hh:917
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
IsAqueous struct.
Definition: components/base.hh:47
Base class for all components Components provide the thermodynamic relations for the liquid,...
Definition: components/base.hh:59
Scalar Scalar
export the scalar type used by the component
Definition: components/base.hh:63
Interface for components that have a gas state.
Definition: gas.hh:41
Material properties of pure water .
Definition: h2o.hh:61
static Scalar liquidViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of pure water.
Definition: h2o.hh:735
static Scalar vaporTemperature(Scalar pressure)
The vapor temperature in of pure water at a given pressure.
Definition: h2o.hh:152
static Scalar liquidPressure(Scalar temperature, Scalar density)
The pressure of liquid water in at a given density and temperature.
Definition: h2o.hh:679
static Scalar liquidThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of water (IAPWS) .
Definition: h2o.hh:758
static const Scalar liquidInternalEnergy(Scalar temperature, Scalar pressure)
Specific internal energy of liquid water .
Definition: h2o.hh:307
static Scalar gasHeatCapacityConstVolume(Scalar temperature, Scalar pressure)
Specific isochoric heat capacity of steam and water vapor .
Definition: h2o.hh:452
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: h2o.hh:561
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of steam.
Definition: h2o.hh:716
static const Scalar liquidHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of liquid water .
Definition: h2o.hh:281
static const Scalar liquidEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of liquid water .
Definition: h2o.hh:217
static Scalar liquidDensity(Scalar temperature, Scalar pressure)
The density of pure water in at a given pressure and temperature.
Definition: h2o.hh:614
static const Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of water steam .
Definition: h2o.hh:173
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of steam in at a given pressure and temperature.
Definition: h2o.hh:493
static Scalar liquidMolarDensity(Scalar temperature, Scalar pressure)
The molar density of water in at a given pressure and temperature.
Definition: h2o.hh:663
static constexpr Scalar acentricFactor()
The acentric factor of water.
Definition: h2o.hh:87
static Scalar vaporPressure(Scalar T)
The vapor pressure in of pure water at a given temperature.
Definition: h2o.hh:131
static constexpr Scalar criticalTemperature()
Returns the critical temperature of water.
Definition: h2o.hh:93
static Scalar gasThermalConductivity(const Scalar temperature, const Scalar pressure)
Thermal conductivity of steam (IAPWS) .
Definition: h2o.hh:791
static constexpr Scalar triplePressure()
Returns the pressure at water's triple point.
Definition: h2o.hh:117
static Scalar gasInternalEnergy(Scalar temperature, Scalar pressure)
Specific internal energy of steam and water vapor .
Definition: h2o.hh:357
static const Scalar gasHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of water steam .
Definition: h2o.hh:252
static constexpr bool gasIsCompressible()
Returns true if the gas phase is assumed to be compressible.
Definition: h2o.hh:472
static const Scalar liquidHeatCapacityConstVolume(Scalar temperature, Scalar pressure)
Specific isochoric heat capacity of liquid water .
Definition: h2o.hh:424
static std::string name()
A human readable name for the water.
Definition: h2o.hh:75
static constexpr Scalar molarMass()
The molar mass in of water.
Definition: h2o.hh:81
static constexpr Scalar criticalPressure()
Returns the critical pressure of water.
Definition: h2o.hh:99
static constexpr Scalar criticalMolarVolume()
Returns the molar volume of water at the critical point.
Definition: h2o.hh:105
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: h2o.hh:478
static Scalar gasPressure(Scalar temperature, Scalar density)
The pressure of steam in at a given density and temperature.
Definition: h2o.hh:576
static constexpr Scalar tripleTemperature()
Returns the temperature at water's triple point.
Definition: h2o.hh:111
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of steam in at a given pressure and temperature.
Definition: h2o.hh:555
Implements relations which are common for all regions of the IAPWS '97 formulation.
Definition: common.hh:57
static constexpr Scalar criticalMolarVolume
Critical molar volume of water .
Definition: common.hh:72
static constexpr Scalar Rs
Specific gas constant of water .
Definition: common.hh:63
static constexpr Scalar triplePressure
Triple pressure of water .
Definition: common.hh:84
static constexpr Scalar molarMass
The molar mass of water .
Definition: common.hh:60
static Scalar viscosity(Scalar temperature, Scalar rho)
The dynamic viscosity of pure water.
Definition: common.hh:100
static constexpr Scalar acentricFactor
The acentric factor of water .
Definition: common.hh:75
static constexpr Scalar criticalTemperature
Critical temperature of water .
Definition: common.hh:66
static constexpr Scalar tripleTemperature
Triple temperature of water .
Definition: common.hh:81
static constexpr Scalar criticalPressure
Critical pressure of water .
Definition: common.hh:69
static Scalar thermalConductivityIAPWS(const Scalar T, const Scalar rho)
Thermal conductivity water (IAPWS) .
Definition: common.hh:161
Implements the equations for region 1 of the IAPWS '97 formulation.
Definition: region1.hh:50
static Scalar dGamma_dPi(Scalar temperature, Scalar pressure)
The partial derivative of the Gibbs free energy to the normalized pressure for IAPWS region 1 (i....
Definition: region1.hh:182
static Scalar ddGamma_ddPi(Scalar temperature, Scalar pressure)
The second partial derivative of the Gibbs free energy to the normalized pressure for IAPWS region 1 ...
Definition: region1.hh:241
static constexpr Scalar dPi_dp(Scalar pressure)
Returns the derivative of the reduced pressure to the pressure for IAPWS region 1 in .
Definition: region1.hh:107
static constexpr Scalar tau(Scalar temperature)
Returns the reduced temperature for IAPWS region 1.
Definition: region1.hh:81
static void checkValidityRange(Scalar temperature, Scalar pressure, const std::string &propertyName="This property")
Returns true if IAPWS region 1 applies for a (temperature in , pressure in ) pair.
Definition: region1.hh:60
static Scalar dGamma_dTau(Scalar temperature, Scalar pressure)
The partial derivative of the Gibbs free energy to the normalized temperature for IAPWS region 1 (i....
Definition: region1.hh:154
static constexpr Scalar pi(Scalar pressure)
Returns the reduced pressure for IAPWS region 1.
Definition: region1.hh:98
static Scalar ddGamma_ddTau(Scalar temperature, Scalar pressure)
The second partial derivative of the Gibbs free energy to the normalized temperature for IAPWS region...
Definition: region1.hh:270
static Scalar ddGamma_dTaudPi(Scalar temperature, Scalar pressure)
The partial derivative of the Gibbs free energy to the normalized pressure and to the normalized temp...
Definition: region1.hh:211
Implements the equations for region 2 of the IAPWS '97 formulation.
Definition: region2.hh:52
static Scalar ddGamma_ddTau(Scalar temperature, Scalar pressure)
The second partial derivative of the Gibbs free energy to the normalized temperature for IAPWS region...
Definition: region2.hh:300
static Scalar dGamma_dTau(Scalar temperature, Scalar pressure)
The partial derivative of the Gibbs free energy to the normalized temperature for IAPWS region 2 (i....
Definition: region2.hh:165
static Scalar ddGamma_dTaudPi(Scalar temperature, Scalar pressure)
The partial derivative of the Gibbs free energy to the normalized pressure and to the normalized temp...
Definition: region2.hh:234
static Scalar ddGamma_ddPi(Scalar temperature, Scalar pressure)
The second partial derivative of the Gibbs free energy to the normalized pressure for IAPWS region 2 ...
Definition: region2.hh:267
static void checkValidityRange(Scalar temperature, Scalar pressure, const std::string &propertyName="This property")
Returns true if IAPWS region 2 applies for a (temperature, pressure) pair.
Definition: region2.hh:62
static constexpr Scalar pi(Scalar pressure)
Returns the reduced pressure (dimensionless) for IAPWS region 2.
Definition: region2.hh:100
static constexpr Scalar dPi_dp(Scalar pressure)
Returns the derivative of the reduced pressure to the pressure for IAPWS region 2 in .
Definition: region2.hh:109
static constexpr Scalar tau(Scalar temperature)
Returns the reduced temperature (dimensionless) for IAPWS region 2.
Definition: region2.hh:83
static Scalar dGamma_dPi(Scalar temperature, Scalar pressure)
The partial derivative of the Gibbs free energy to the normalized pressure for IAPWS region 2 (i....
Definition: region2.hh:202
Implements the equations for region 4 of the IAPWS '97 formulation.
Definition: region4.hh:53
static Scalar saturationPressure(Scalar temperature)
Returns the saturation pressure in of pure water at a given temperature.
Definition: region4.hh:64
static Scalar vaporTemperature(Scalar pressure)
Returns the saturation temperature in of pure water at a given pressure.
Definition: region4.hh:95
Interface for components that have a liquid state.
Definition: liquid.hh:41
static constexpr Scalar pressure(Scalar temperature, Scalar rhoMolar)
The pressure of the gas in , depending on the molar density and temperature.
Definition: idealgas.hh:60
static constexpr Scalar density(Scalar avgMolarMass, Scalar temperature, Scalar pressure)
The density of the gas in , depending on pressure, temperature and average molar mass of the gas.
Definition: idealgas.hh:49
Base class for all components Components provide the thermodynamic relations for the liquid,...