FPROPS
Fluids supported by FPROPS
- models/johnpye/fprops/ammonia.c
- models/johnpye/fprops/hydrogen.c
- models/johnpye/fprops/nitrogen.c
- models/johnpye/fprops/water.c
- models/johnpye/fprops/carbondioxide.c
- Some more, yet to be merged: hongke:models/johnpye/fprops/
FPROPS is a free open-source C-based library for evaluating thermodynamic properties of a number of pure substances via published data for the Helmholtz fundamental equation for those substances. It has been developed by John Pye and will happily function as standalone code, but is also provided with external library code for ASCEND so that it can be used to access these accurate property correlations from within a MODEL.
Currently FPROPS supports calculation of the properties of the pure substances shown at right. The properties that can be calculated are internal energy, entropy, pressure, enthalpy and Helmholtz energy, as well as various partial derivatives of these with respect to temperature and density.
FPROPS reproduces a limited subset of the functionality of commercial programs such as REFPROP, PROPATH, EES, FLUIDCAL, freesteam, SteamTab, and others, but is open source software and licensed under the GPL. Contributions are most welcome!
IMPORTANT NOTE: FPROPS does not yet implement calculation of the saturation curve, so it gives incorrect property evaluations anywhere within the saturation region. This is because the code for evaluating the 'Maxwell criterion' is not yet completely implemented.
FPROPS is included with ASCEND as of version 0.9.5.116, or you can browse the code at models/johnpye/fprops/ or access it via subversion (there have been changes since the 0.9.5.116 release).
Implementation
The calculation routines implement calculation of reduced Helmholtz energy φ by dividing it into ideal (<math>\phi^0</math>) and residual (<math>\phi^r</math>) parts, as described in IAPWS-95[1].
Ideal part
Calculation of the ideal part of the reduced Helmholtz energy is derived from the expression for the zero-pressure limit of reduced specific heat capacity at constant pressure, <math>c_p^0 / R</math>, where <math>R</math> is the specific gas constant. This quantity <math>{c_p^0}(T)/R</math> can be expressed as a series composed of two different types of terms:
- power terms, <math>c_i T^{t_i}</math>
- exponential terms, <math>\frac{b_i x_i^2 exp(-x_i)}{[1-exp(-x_i)]^2}</math> where <math>x_i = \frac{\beta_i}{T}</math>
An expression for calculating the ideal component of Helmholtz energy <math>a^0 \left(T , \rho \right)</math> from <math>c_p^0</math> is given by Span et al[2]:
- <math>\frac{a^0}{RT} = \phi^0 = \frac{h_0^0 \tau}{R T_c} - \frac{s_0^0}{R} - 1 + \ln {\frac{\delta / \delta_0}{\tau / \tau_0}} - \frac{\tau}{R} \int_{\tau_0}^{\tau}{\frac{c_p^0}{\tau^2}d\tau} + \frac{1}{R} \int_{\tau_0}^{\tau}{\frac{c_p^0}{\tau} d \tau}</math>
where the reduced temperature is <math>\tau = \frac{T^{*}}{T}</math> and the reduced density is <math>\delta = \frac{\rho}{\rho^{*}}</math>.
On the other hand, using the relations for the Helmholtz fundamental equation, one can see that
- <math>\frac{c_p^0}{R} = 1 -\tau^2 \phi^0_{\tau \tau}</math>
so coefficients of <math>c_p^0</math> can be calculated from coefficients of <math>\phi^0</math> (providing, we assume, that the fundamental equation has smooth first and second (and possibly third) derivatives, and that when terms of <math>{c_p^0}(T)</math> are integrated they don't 'interact'). One can go the other way, and determine <math>\phi^0</math> from <math>c_p^0</math>, with the exception of the integration constants, which would in that case be determined using specified values of enthalpy <math>h_ref</math> and entropy <math>s_ref</math> at a reference state <math>\left(T_{ref},\rho_{ref} \right)</math>
Residual part
For the residual part of the reduced Helmholtz energy, FPROPS currently supports three different kinds of terms:
| type | equation | notes | correlation constants |
|---|---|---|---|
| power terms | <math>
\begin{cases} a_i \tau^{t_i} \delta^{d_i} exp(-\delta^{l_i}), & \mbox{if }l \ne 0 \\
a_i \tau^{t_i} \delta^{d_i}, & \mbox{otherwise}
\end{cases}
| ||
| gaussian terms | <math>n_i \tau^{t_i} \delta^{d_i} exp[-\alpha_i(\delta - \epsilon_i)^2 - \beta_i(\tau - \gamma_i)^2]</math> | <math> n, t, d, \alpha, \beta, \gamma, \epsilon</math> | |
| critical terms | <math>n_i \Delta^{b_i} \delta \psi</math> | where
|
<math>n, a, b, \beta, A, B, C, D</math> |
These different terms are used in various combinations across a number of high-accuracy publications of thermodynamic property data. Overall, the residual series permits calculation of 'real' reduced Helmholtz energy <math>\phi \left(\tau,\delta \right) = \phi^0 \left(\tau,\delta \right) + \phi^r \left( \tau,\delta \right)</math> and hence the specific Helmholtz energy <math>\frac{a}{R T} = \phi \left(\delta,\tau\right)</math>.
From the expression for <math>\phi</math>, we can then derive expressions for <math>p</math>, <math>h</math>, <math>u</math>, <math>s</math> and their partial derivatives with respect to density and temperature. These expressions are summarised in IAPWS-95 as well as by Tillner Roth et al[3].
Maxwell criterion
The Helmholtz function <math>a \left(\rho,T \right)</math>, given above, applies for all regions outside the saturation region. During saturation conditions, stable equilibrium of phases occurs when the liquid phase and vapour phase have equal Gibbs free energy <math>g = a + p / \rho</math>. This requirement is equivalent to the Maxwell equal area rule[4][1][5]. The condition can be expressed as a set of three simultaneous equations in the variables <math>p_{sat}, \rho_f, \rho_g</math>:
<math> \begin{align} p_{sat} & = p \left(\rho_f, T \right) \\ p_{sat} & = p \left(\rho_g, T \right) \\ g\left(\rho_f,T \right) &= g\left(\rho_g,T \right) \end{align} </math>
This can be further reduced to :
<math> \begin{align} p_{sat} & = R T \rho_f \left(1 + \delta_f \phi^r_\delta \left(\delta_f,\tau \right) \right) \\ p_{sat} & = R T \rho_g \left(1 + \delta_g \phi^r_\delta \left(\delta_g,\tau \right) \right) \\ \frac{p_{sat}}{R T} \left( \frac{1}{\rho_g} - \frac{1}{\rho_f} \right) & = \phi \left(\delta_f, \tau \right) - \phi \left(\delta_g, \tau \right)
= \ln \left( \frac{\rho_f}{\rho_g} \right) + \phi^r \left(\delta_f, \tau \right) - \phi^r \left(\delta_g, \tau \right)
\end{align} </math>
Calculation of the phase equilibrium conditions requires iterative solution of three equations; one suitable approach[6][7] for example, is:
- make an initial guess for <math>p_{sat} \left(T\right)</math>. This can be derived using the definition of the acentric factor <math>\omega</math>, together with an assumed form <math>\log_{10} \left(p \right) = A + \frac{B}{T}</math>.
- use first-order (or second-order) Newton's method to solve <math>\rho_f \left( p_{sat}, T \right)</math>, then again for <math>\rho_g \left( p_{sat}, T \right)</math>, taking care to avoid guesses in the region where <math>\left. \frac{\partial p}{\partial \rho} \right|_T < 0</math>. For the iteration, <math>\ln v</math> is used in preference to the native variable <math>\rho</math>.
- calculate <math>\Delta a = a \left(\rho_f, T \right) - a \left(\rho_g, T \right)</math>, the difference in the values of Helmholtz energy.
- replace <math>p_{sat} = \frac{\Delta a}{\frac{1}{\rho_g} - \frac{1}{\rho_f}}</math>
- once <math>p_{sat}</math> stops changing, solution has converged.
This approach (see models/johnpye/fprops/sat2.c) converges well except for in the vicinity of the critical point. New code will be required for values of temperature near the critical point. It is proposed to use the method of Akasaka[6], but this is not yet implemented.
Testing the code
FPROPS includes a suite of self tests for each fluid that has been implemented. Some of the test data comes from the publications from which the correlations have been taken; other test data comes from using the commercial package REFPROP to calculate some sample data points, and comparing the FPROPS output with that data.
To run the FPROPS tests for a given fluid, cd to the FPROPS directory (eg cd ~/ascend/models/johnpye/fprops</tdd>) then run the test script <tt>./test.py fluidname, for example:
./test.py hydrogen ./test.py nitrogen
Requires GNU scientific library(libgsl0-dev) to be installed in the system.
The output will tell you if all calculated values have conformed to expected error margins.
Test code for each fluid is embedded at the end of the corresponding fluid file, for example models/johnpye/fprops/hydrogen.c.
Adding new fluids
Currently, new fluids must be added by declaring constant data structures in C code.
An example of a fluid declaration is shown in models/johnpye/fprops/hydrogen.c.
After adding the data structures in "component".c file, a corresponding "component".h is required.
An example of which is shown in models/johnpye/fprops/hydrogen.h
The new component has to be declared in SConstruct and asc_helmholtz.c.
The data structures which are to be filled are declared in models/johnpye/fprops/ideal.h and models/johnpye/fprops/helmholtz.h.
Note: when adding new fluids it is important to observe that different papers give different 'zero points' for enthalpy and entropy for the substances in question, and the zero points vary from paper to paper. These Helmholtz equations of state can be adjusted for different 'zero points' by adjusting the constant and linear parameters ('c' and 'm') in the IdealData object that is provided in your data structures. The equation to fix these values is
- <math>\frac{\Delta h}{RT} = \tau \Delta \phi_{\tau}^0 = T^{*} \Delta a_2^0</math>
- <math>\frac{\Delta s}{R} = \tau \Delta \phi_{\tau}^0 - \Delta \phi^0 = - \Delta a_1^0</math>
where we assume
- <math>\phi^0(\tau,\delta) = a_1^0 + a_2^0 \tau + \mathrm{other\ terms\ in\ \tau, \delta}</math>
In other words, to correct an offset in h and s, these formulae tell you the offset to apply to the constant and linear terms in φ0(τ,δ).
Thermodynamic property derivatives
The IAPWS has made a release[8] with information about how partial derivatives of thermodynamic properties with respect to other properties can be calculated in general. It is intended that we will add support for this in FPROPS using code already developed for the IAPWS-IF97 code freesteam.
See also derivs.c from freesteam.
Using within ASCEND
To use FPROPS from ASCEND, see the example code in models/johnpye/fprops/fprops_test.a4c. Essentially you must declare a DATA instance into which you must place the name of the substance for which you wish to calculate properties. Then you can use the functions helmholtz_p, helmholtz_u, helmholtz_h, helmholtz_p, helmholtz_s, helmholtz_a to calculate the property you need. A function helmholtz_phsx_vT is also available which returns all of the above properties, including correct saturation properties.
Note. When solving properties in 'reverse' mode (eg solving for ρ,T when p,h are know), you will get significantly better results by switching the QRSlv solver parameter for 'convergence test' (convopt) to 'RELNOM_SCALE'. This significantly improves ASCEND's ability to solve these problems. Also, use OPTION iterationlimit 200 for larger models.
A first system model using FPROPS is available in models/johnpye/rankine_fprops.a4c.
Some details on the implementation of the 'wrapper' code for connecting FPROPS with ASCEND is given in writing ASCEND external relations in C.
Using from C
Using the FPROPS code from C should also be very straightforward; example code is shown in models/johnpye/fprops/test.c.
Further work
FPROPS is still under development! We aim to implement a few more features including
- add support for calculation of the saturation curve using Maxwell phase-equilibrium condition (Assigned to Shrikanth)
support for Water using IAPWS-95 correlation (partially complete, but needs support for another type of residual function term).refactor the 'HelmholtzExpTerm' which is actually a double of 'HelmholtzGausTerm' with epsilon = 1.add support for calculation of partial derivatives dp/dT, dp/drho, dh/dT, dh/drho, du/dT, du/drho,ds/dT, ds/drho, da/dT, da/drhoto improve convergence of ASCEND models that calculate p, h, u,s, a etc.- investigate apparent disagreement between cp0 for ammonia in FPROPS and values from REFPROP (is it consistent with Tillner-Roth?)
- investigate apparent very small discrepancy between REFPROP and IAPWS95 in value of speed of sound near critical point.
- move specification of the zero point out of IdealData and into HelmholtzData.
- implement a pre-calculation routine where for example coefficients of <math>\phi^0</math> can be prepared?
- implement support for thermophysical properties conductivity, viscosity, surface tension, and possibly others.
implement support for calculation of the speed of sound.- resolve the confusion between T* and T_crit in the data structures. Probably all correlations will be using T_crit -- but is that sure?
- remove p_c from HelmholtzData structure, because that value can be calculated using p(T_c,rho_c).
Further down the track we would like to support:
- reading property correlations from a database (SQLite?) or text file.
- some simpler equations of state such as MBWR or Peng-Robinson.(Assign to Ankit)
- more fluids of general interest.
- mixing models for multi-phase, multi-component systems (or else incorporate FPROPS into exist code that supports this).
We'd be very pleased to be able to collaborate on this work with others; the code has been made deliberately able to stand alone; it has no dependency on any other code from ASCEND, and is suitable for incorporation into other open source programs.
References
- ↑ 1.0 1.1 http://www.iapws.org/relguide/IAPWS95-Rev.pdf
- ↑ Span & Wagner (1996), "A new equation of state for carbon dioxide covering from the triple point temperature to 1100 K at pressures up to 800 MPa" J. Phys. Chem. Ref. Data, 25 p 1509 ff
- ↑ Tillner-Roth, Harms-Watzenberg and Baehr (1993), "Eine neue Fundamentalgleichung für Ammoniak", DKV-Tagungsbericht, 20, 167-181
- ↑ http://en.wikipedia.org/wiki/Van_der_Waals_equation
- ↑ Eubank & Hall (2004), "Equal area rule and algorithm for determining phase compositions", AIChE Journal, 41, http://dx.doi.org/10.1002/aic.690410419
- ↑ 6.0 6.1 Akasaka (2008) "A Reliable and Useful Method to Determine the Saturation State from Helmholtz Energy Equations of State", Journal of Thermal Science and Technology, 3, 442-451, http://dx.doi.org/10.1299/jtst.3.442
- ↑ Soave (1986), "Direct calculation of pure-compound vapour pressures through cubic equations of state ", Fluid Phase Equilibria, 31, 203-207, http://dx.doi.org/10.1016/0378-3812(86)90013-0
- ↑ IAPWS, 2008, Revised Advisory Note No. 3 Thermodynamic Derivatives from IAPWS Formulation, http://www.iapws.org
Contributors
The people who have contributed code and/or implementable advice for FPROPS: