Cantera  3.1.0a1
WaterPropsIAPWS Class Reference

Class for calculating the equation of state of water. More...

#include <WaterPropsIAPWS.h>

Detailed Description

Class for calculating the equation of state of water.

This is a helper class for WaterSSTP and PDSS_Water and does not constitute a complete implementation of a thermo phase by itself (see Thermodynamic Properties and classes WaterSSTP and PDSS_Water).

The reference is W. Wagner, A. Pruss, "The IAPWS Formulation 1995 for the Thermodynamic Properties of Ordinary Water Substance for General and Scientific Use," J. Phys. Chem. Ref. Dat, 31, 387, 2002.

This class provides a very complicated polynomial for the specific Helmholtz free energy of water, as a function of temperature and density.

\[ \frac{M\hat{f}(\rho,T)}{R T} = \phi(\delta, \tau) = \phi^o(\delta, \tau) + \phi^r(\delta, \tau) \]

where

\[ \delta = \rho / \rho_c \quad \mathrm{and} \quad \tau = T_c / T \]

The following constants are assumed

\[ T_c = 647.096\mathrm{\;K} \]

\[ \rho_c = 322 \mathrm{\;kg\,m^{-3}} \]

\[ R/M = 0.46151805 \mathrm{\;kJ\,kg^{-1}\,K^{-1}} \]

The free energy is a unique single-valued function of the temperature and density over its entire range.

Note, the base thermodynamic state for this class is the one used in the steam tables, such that the liquid at the triple point for water has the following properties:

  • u(273.16, rho) = 0.0
  • s(273.16, rho) = 0.0
  • psat(273.16) = 611.655 Pascal
  • rho(273.16, psat) = 999.793 kg m-3

Therefore, to use this class within Cantera, offsets to u() and s() must be used to put the water class onto the same basis as other thermodynamic quantities. For example, in the WaterSSTP class, these offsets are calculated in the following way. The thermodynamic base state for water is set to the NIST basis here by specifying constants EW_Offset and SW_Offset. These offsets are calculated on the fly so that the following properties hold:

  • Delta_Hfo_idealGas(298.15, 1bar) = -241.826 kJ/gmol
  • So_idealGas(298.15, 1bar) = 188.835 J/gmolK

The offsets are calculated by actually computing the above quantities and then calculating the correction factor.

This class provides an interface to the WaterPropsIAPWSphi class, which actually calculates the \( \phi^o(\delta, \tau) \) and the \( \phi^r(\delta, \tau) \) polynomials in dimensionless form.

All thermodynamic results from this class are returned in dimensional form. This is because the gas constant (and molecular weight) used within this class is allowed to be potentially different than that used elsewhere in Cantera. Therefore, everything has to be in dimensional units. Note, however, the thermodynamic basis is set to that used in the steam tables. (u = s = 0 for liquid water at the triple point).

This class is not a ThermoPhase. However, it does maintain an internal state of the object that is dependent on temperature and density. The internal state is characterized by an internally stored \( \tau \) and a \( \delta \) value, and an iState value, which indicates whether the point is a liquid, a gas, or a supercritical fluid. Along with that the \( \tau \) and a \( \delta \) values are polynomials of \( \tau \) and a \( \delta \) that are kept by the WaterPropsIAPWSphi class. Therefore, whenever \( \tau \) or \( \delta \) is changed, the function setState() must be called in order for the internal state to be kept up to date.

The class is pretty straightforward. However, one function deserves mention. The density() function calculates the density that is consistent with a particular value of the temperature and pressure. It may therefore be multivalued or potentially there may be no answer from this function. It therefore takes a phase guess and a density guess as optional parameters. If no guesses are supplied to density(), a gas phase guess is assumed. This may or may not be what is wanted. Therefore, density() should usually at least be supplied with a phase guess so that it may manufacture an appropriate density guess. density() manufactures the initial density guess, nondimensionalizes everything, and then calls WaterPropsIAPWSphi::dfind(), which does the iterative calculation to find the density condition that matches the desired input pressure.

The phase guess defines are located in the .h file. they are

  • WATER_GAS
  • WATER_LIQUID
  • WATER_SUPERCRIT

There are only three functions which actually change the value of the internal state of this object after it's been instantiated

  • setState_TD(temperature, rho)
  • density(temperature, pressure, phase, rhoguess)
  • psat(temperature, waterState);

The setState_TD() is the main function that sets the temperature and rho value. The density() function serves as a setState_TP() function, in that it sets internal state to a temperature and pressure. However, note that this is potentially multivalued. Therefore, we need to supply in addition a phase guess and a rho guess to the input temperature and pressure. The psat() function sets the internal state to the saturated liquid or saturated gas state, depending on the waterState parameter.

Because the underlying object WaterPropsIAPWSphi is privately held, you can be sure that the underlying state of this object doesn't change except due to the three function calls listed above.

Definition at line 161 of file WaterPropsIAPWS.h.

Public Member Functions

 WaterPropsIAPWS ()=default
 Base constructor. More...
 
 WaterPropsIAPWS (const WaterPropsIAPWS &right)=delete
 
WaterPropsIAPWSoperator= (const WaterPropsIAPWS &right)=delete
 
void setState_TD (double temperature, double rho)
 Set the internal state of the object wrt temperature and density. More...
 
double gibbs_mass () const
 Get the Gibbs free energy (J/kg) at the current temperature and density. More...
 
double enthalpy_mass () const
 Get the enthalpy (J/kg) at the current temperature and density. More...
 
double intEnergy_mass () const
 Get the internal energy (J/kg) at the current temperature and density. More...
 
double entropy_mass () const
 Get the entropy (J/kg/K) at the current temperature and density. More...
 
double cv_mass () const
 Get the constant volume heat capacity (J/kg/K) at the current temperature and density. More...
 
double cp_mass () const
 Get the constant pressure heat capacity (J/kg/K) at the current temperature and density. More...
 
double pressure () const
 Calculates the pressure (Pascals), given the current value of the temperature and density. More...
 
double density (double temperature, double pressure, int phase=-1, double rhoguess=-1.0)
 Calculates the density given the temperature and the pressure, and a guess at the density. More...
 
double density_const (double pressure, int phase=-1, double rhoguess=-1.0) const
 Calculates the density given the temperature and the pressure, and a guess at the density, while not changing the internal state. More...
 
double density () const
 Returns the density (kg m-3) More...
 
double temperature () const
 Returns the temperature (Kelvin) More...
 
double coeffThermExp () const
 Returns the coefficient of thermal expansion. More...
 
double coeffPresExp () const
 Returns the isochoric pressure derivative wrt temperature. More...
 
double isothermalCompressibility () const
 Returns the coefficient of isothermal compressibility for the state of the object. More...
 
double dpdrho () const
 Returns the value of dp / drho at constant T for the state of the object. More...
 
double psat_est (double temperature) const
 This function returns an estimated value for the saturation pressure. More...
 
double psat (double temperature, int waterState=WATER_LIQUID)
 This function returns the saturation pressure given the temperature as an input parameter, and sets the internal state to the saturated conditions. More...
 
double densSpinodalWater () const
 Return the value of the density at the water spinodal point (on the liquid side) for the current temperature. More...
 
double densSpinodalSteam () const
 Return the value of the density at the water spinodal point (on the gas side) for the current temperature. More...
 
int phaseState (bool checkState=false) const
 Returns the Phase State flag for the current state of the object. More...
 
double Tcrit () const
 Returns the critical temperature of water (Kelvin) More...
 
double Pcrit () const
 Returns the critical pressure of water (22.064E6 Pa) More...
 
double Rhocrit () const
 Return the critical density of water (kg m-3) More...
 

Private Member Functions

void calcDim (double temperature, double rho)
 Calculate the dimensionless temp and rho and store internally. More...
 
void corr (double temperature, double pressure, double &densLiq, double &densGas, double &delGRT)
 Utility routine in the calculation of the saturation pressure. More...
 
void corr1 (double temperature, double pressure, double &densLiq, double &densGas, double &pcorr)
 Utility routine in the calculation of the saturation pressure. More...
 

Private Attributes

WaterPropsIAPWSphi m_phi
 pointer to the underlying object that does the calculations. More...
 
double tau = -1.0
 Dimensionless temperature, tau = T_C / T. More...
 
double delta = -1.0
 Dimensionless density, delta = rho / rho_c. More...
 
int iState = -30000
 Current state of the system. More...
 

Constructor & Destructor Documentation

◆ WaterPropsIAPWS()

WaterPropsIAPWS ( )
default

Base constructor.

Member Function Documentation

◆ setState_TD()

void setState_TD ( double  temperature,
double  rho 
)

Set the internal state of the object wrt temperature and density.

Parameters
temperaturetemperature (kelvin)
rhodensity (kg m-3)
Since
New in Cantera 3.0.

Definition at line 548 of file WaterPropsIAPWS.cpp.

◆ gibbs_mass()

double gibbs_mass ( ) const

Get the Gibbs free energy (J/kg) at the current temperature and density.

Definition at line 554 of file WaterPropsIAPWS.cpp.

◆ enthalpy_mass()

double enthalpy_mass ( ) const

Get the enthalpy (J/kg) at the current temperature and density.

Definition at line 559 of file WaterPropsIAPWS.cpp.

◆ intEnergy_mass()

double intEnergy_mass ( ) const

Get the internal energy (J/kg) at the current temperature and density.

Definition at line 564 of file WaterPropsIAPWS.cpp.

◆ entropy_mass()

double entropy_mass ( ) const

Get the entropy (J/kg/K) at the current temperature and density.

Definition at line 569 of file WaterPropsIAPWS.cpp.

◆ cv_mass()

double cv_mass ( ) const

Get the constant volume heat capacity (J/kg/K) at the current temperature and density.

Definition at line 574 of file WaterPropsIAPWS.cpp.

◆ cp_mass()

double cp_mass ( ) const

Get the constant pressure heat capacity (J/kg/K) at the current temperature and density.

Definition at line 579 of file WaterPropsIAPWS.cpp.

◆ pressure()

double pressure ( ) const

Calculates the pressure (Pascals), given the current value of the temperature and density.

The density is an independent variable in the underlying equation of state

Returns
the pressure (Pascal)

Definition at line 46 of file WaterPropsIAPWS.cpp.

◆ density() [1/2]

double density ( double  temperature,
double  pressure,
int  phase = -1,
double  rhoguess = -1.0 
)

Calculates the density given the temperature and the pressure, and a guess at the density.

Sets the internal state.

Note, below T_c, this is a multivalued function.

The density() function calculates the density that is consistent with a particular value of the temperature and pressure. It may therefore be multivalued or potentially there may be no answer from this function. It therefore takes a phase guess and a density guess as optional parameters. If no guesses are supplied to density(), a gas phase guess is assumed. This may or may not be what is wanted. Therefore, density() should usually at least be supplied with a phase guess so that it may manufacture an appropriate density guess. density() manufactures the initial density guess, nondimensionalizes everything, and then calls WaterPropsIAPWSphi::dfind(), which does the iterative calculation to find the density condition that matches the desired input pressure.

Parameters
temperatureKelvin
pressurePressure in Pascals (Newton/m**2)
phaseguessed phase of water; -1: no guessed phase
rhoguessguessed density of the water; -1.0 no guessed density
Returns
the density. If an error is encountered in the calculation the value of -1.0 is returned.

Definition at line 54 of file WaterPropsIAPWS.cpp.

◆ density_const()

double density_const ( double  pressure,
int  phase = -1,
double  rhoguess = -1.0 
) const

Calculates the density given the temperature and the pressure, and a guess at the density, while not changing the internal state.

Note, below T_c, this is a multivalued function.

The density() function calculates the density that is consistent with a particular value of the temperature and pressure. It may therefore be multivalued or potentially there may be no answer from this function. It therefore takes a phase guess and a density guess as optional parameters. If no guesses are supplied to density(), a gas phase guess is assumed. This may or may not be what is wanted. Therefore, density() should usually at least be supplied with a phase guess so that it may manufacture an appropriate density guess. density() manufactures the initial density guess, nondimensionalizes everything, and then calls WaterPropsIAPWSphi::dfind(), which does the iterative calculation to find the density condition that matches the desired input pressure.

Parameters
pressurePressure in Pascals (Newton/m**2)
phaseguessed phase of water; -1: no guessed phase
rhoguessguessed density of the water; -1.0: no guessed density
Returns
the density. If an error is encountered in the calculation the value of -1.0 is returned.

Definition at line 114 of file WaterPropsIAPWS.cpp.

◆ density() [2/2]

double density ( ) const

Returns the density (kg m-3)

The density is an independent variable in the underlying equation of state

Returns
the density (kg m-3)

Definition at line 167 of file WaterPropsIAPWS.cpp.

◆ temperature()

double temperature ( ) const

Returns the temperature (Kelvin)

Returns
s the internally stored temperature

Definition at line 172 of file WaterPropsIAPWS.cpp.

◆ coeffThermExp()

double coeffThermExp ( ) const

Returns the coefficient of thermal expansion.

alpha = d (ln V) / dT at constant P.

Returns
the coefficient of thermal expansion

Definition at line 234 of file WaterPropsIAPWS.cpp.

◆ coeffPresExp()

double coeffPresExp ( ) const

Returns the isochoric pressure derivative wrt temperature.

beta = M / (rho * Rgas) (d (pressure) / dT) at constant rho

Note for ideal gases this is equal to one.

beta = delta (phi0_d() + phiR_d()) - tau delta (phi0_dt() + phiR_dt())

Definition at line 229 of file WaterPropsIAPWS.cpp.

◆ isothermalCompressibility()

double isothermalCompressibility ( ) const

Returns the coefficient of isothermal compressibility for the state of the object.

kappa = - d (ln V) / dP at constant T.

units - 1/Pascal

Returns
the isothermal compressibility

Definition at line 215 of file WaterPropsIAPWS.cpp.

◆ dpdrho()

double dpdrho ( ) const

Returns the value of dp / drho at constant T for the state of the object.

units - Joules / kg

Returns
dpdrho

Definition at line 222 of file WaterPropsIAPWS.cpp.

◆ psat_est()

double psat_est ( double  temperature) const

This function returns an estimated value for the saturation pressure.

It does this via a polynomial fit of the vapor pressure curve. units = (Pascals)

Parameters
temperatureInput temperature (Kelvin)
Returns
the estimated saturation pressure

Definition at line 177 of file WaterPropsIAPWS.cpp.

◆ psat()

double psat ( double  temperature,
int  waterState = WATER_LIQUID 
)

This function returns the saturation pressure given the temperature as an input parameter, and sets the internal state to the saturated conditions.

Note this function will return the saturation pressure, given the temperature. It will then set the state of the system to the saturation condition. The input parameter waterState is used to either specify the liquid state or the gas state at the desired temperature and saturated pressure.

If the input temperature, T, is above T_c, this routine will set the internal state to T and the pressure to P_c. Then, return P_c.

Parameters
temperatureinput temperature (kelvin)
waterStateinteger specifying the water state
Returns
the saturation pressure. units = Pascal

Definition at line 291 of file WaterPropsIAPWS.cpp.

◆ densSpinodalWater()

double densSpinodalWater ( ) const

Return the value of the density at the water spinodal point (on the liquid side) for the current temperature.

Returns
the density with units of kg m-3

Definition at line 373 of file WaterPropsIAPWS.cpp.

◆ densSpinodalSteam()

double densSpinodalSteam ( ) const

Return the value of the density at the water spinodal point (on the gas side) for the current temperature.

Returns
the density with units of kg m-3

Definition at line 461 of file WaterPropsIAPWS.cpp.

◆ phaseState()

int phaseState ( bool  checkState = false) const

Returns the Phase State flag for the current state of the object.

Parameters
checkStateIf true, this function does a complete check to see where in parameters space we are

There are three values:

  • WATER_GAS below the critical temperature but below the critical density
  • WATER_LIQUID below the critical temperature but above the critical density
  • WATER_SUPERCRIT above the critical temperature

Definition at line 333 of file WaterPropsIAPWS.cpp.

◆ Tcrit()

double Tcrit ( ) const
inline

Returns the critical temperature of water (Kelvin)

This is hard coded to the value 647.096 Kelvin

Definition at line 370 of file WaterPropsIAPWS.h.

◆ Pcrit()

double Pcrit ( ) const
inline

Returns the critical pressure of water (22.064E6 Pa)

This is hard coded to the value of 22.064E6 pascals

Definition at line 378 of file WaterPropsIAPWS.h.

◆ Rhocrit()

double Rhocrit ( ) const
inline

Return the critical density of water (kg m-3)

This is equal to 322 kg m-3.

Definition at line 386 of file WaterPropsIAPWS.h.

◆ calcDim()

void calcDim ( double  temperature,
double  rho 
)
private

Calculate the dimensionless temp and rho and store internally.

Parameters
temperatureinput temperature (kelvin)
rhodensity in kg m-3

Definition at line 29 of file WaterPropsIAPWS.cpp.

◆ corr()

void corr ( double  temperature,
double  pressure,
double &  densLiq,
double &  densGas,
double &  delGRT 
)
private

Utility routine in the calculation of the saturation pressure.

Calculate the Gibbs free energy in mks units of J kmol-1 K-1.

Parameters
temperaturetemperature (kelvin)
pressurepressure (Pascal)
densLiqOutput density of liquid
densGasoutput Density of gas
delGRToutput delGRT

Definition at line 242 of file WaterPropsIAPWS.cpp.

◆ corr1()

void corr1 ( double  temperature,
double  pressure,
double &  densLiq,
double &  densGas,
double &  pcorr 
)
private

Utility routine in the calculation of the saturation pressure.

Parameters
temperaturetemperature (kelvin)
pressurepressure (Pascal)
densLiqOutput density of liquid
densGasoutput Density of gas
pcorroutput corrected pressure

Definition at line 266 of file WaterPropsIAPWS.cpp.

Member Data Documentation

◆ m_phi

WaterPropsIAPWSphi m_phi
mutableprivate

pointer to the underlying object that does the calculations.

Definition at line 423 of file WaterPropsIAPWS.h.

◆ tau

double tau = -1.0
private

Dimensionless temperature, tau = T_C / T.

Definition at line 426 of file WaterPropsIAPWS.h.

◆ delta

double delta = -1.0
mutableprivate

Dimensionless density, delta = rho / rho_c.

Definition at line 429 of file WaterPropsIAPWS.h.

◆ iState

int iState = -30000
mutableprivate

Current state of the system.

Definition at line 432 of file WaterPropsIAPWS.h.


The documentation for this class was generated from the following files: