FPROPS

From ASCEND
Jump to: navigation, search

Fluids supported by FPROPS

38 high-accuracy Helmholtz fluids: acetone, ammonia, butane, butene, carbondioxide, carbonmonoxide, carbonylsulfide, cisbutene, decane, dimethylether, ethane, ethanol, hydrogen, hydrogensulfide, isobutane, isobutene, isohexane, isopentane, krypton, methane, neopentane, nitrogen, nitrousoxide, nonane, parahydrogen, propane, r116, r134a, r141b, r142b, r218, r245fa, r41, sulfurdioxide, toluene, transbutene, water, xenon

408 cubic EOS fluids: _rpp

FPROPS is a free open-source C-based library for high-accuracy evaluation of thermodynamic properties for a number of pure substances. It makes use of published data for the Helmholtz fundamental equation for those substances. It has been developed by John Pye and others, 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 on the right. The properties that can be calculated are internal energy u, entropy s, pressure p, enthalpy h and Helmholtz energy a, 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, SteamTab and others, but is free open-source software, licensed under the GPL. Similar open-source libraries include freesteam and CoolProp. Contributions are most welcome!

Enhancement to FPROPS feature set and list of supported fluids has been going on as part of ASCEND's involvement in GSOC.

FPROPS is included in current versions of ASCEND, but you can also browse the source code at models/johnpye/fprops/ (current trunk) . Currently (Mar 2014), viscosity and thermal conductivity correlations are being implemented in FPROPS.

Implementation

The calculation routines implement calculation of reduced Helmholtz energy \phi = a / R T as a function of reduced density \delta = \rho / \rho_c and inverse reduced temperature \tau = T_c / T, by dividing it into ideal \phi^0 \left(\tau,\delta\right) and 'real' or 'residual' \phi^r \left(\tau,\delta \right) 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, c_p^0 / R, where R is the specific gas constant.

An expression for c_p^0 can then be used to calculate ideal enthalpy, ideal entropy and ideal helmholtz energy as follows [2]:

h^0 = h_0^0 + \int_{T_0}^{T}{c_p^0 dT}
s^0 = s_0^0 + \int_{T_0}^{T}{\frac{c_p^0}{T}dT} - R \ln \left(\frac{\rho T}{\rho_0 T_0}\right)
a^0 = h^0  - RT - T s^0

Combining these, an expression for calculating the nondimensionalised ideal component of Helmholtz energy a^0 \left(T , \rho \right) from c_p^0 can be found, which incorporates additional correlation constants relating to the reference state \left(T_0, \rho_0, h_0^0, s_0^0 \right): [3]

\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}

where the reduced temperature is \tau = \frac{T^{*}}{T} and the reduced density is \delta = \frac{\rho}{\rho^{*}}.

In the reverse direction, using the relations for the Helmholtz fundamental equation, an important result is

\frac{c_p^0}{R} = 1 -\tau^2 \phi^0_{\tau \tau}

so coefficients of c_p^0 can be calculated from coefficients of \phi^0 (providing, we assume, that the fundamental equation has smooth first and second (and possibly third) derivatives). One can go the other way, and determine \phi^0 from c_p^0, with the exception of the integration constants, which are, in that case, determined using specified values of enthalpy h^0_{ref} and entropy s^0_{ref} at a reference state \left(T_{ref},\rho_{ref} \right)

In FPROPS, then, we permit \phi^0 expressions of this form:

\phi^0 = \ln \left(\delta \right) + c + m \tau + \sum_{i=0}^{\mathrm{np-1}} {\begin{cases}a_i \ln \tau &\mbox{if }p_i =0 \\ a_i \tau^p &\mbox{otherwise}\end{cases}} + \sum_0^{\mathrm{ne}-1} {n_i \ln \left[1 - \exp\left(-\gamma_i \tau\right)\right]}

This corresponds, after the above transformation, to:

\frac{c_p^0}{R} = 1 + \sum_{i=0}^{\mathrm{np}-1} {c_i T^{t_i}} + \sum_{i=0}^{\mathrm{ne}-1} \frac{b_i x_i^2 \exp\left(-x_i\right)}{[1-\exp\left(-x_i\right)]^2} \mbox{  where }x_i = \frac{\beta_i}{T}

The conversion from \phi^0\left(\tau,\delta\right) terms to c_p^0 terms is as follows (\tau = T^{*} / T):

type term in \phi^0 term in c_p^0 conversion
power terms \begin{cases} a \ln \tau, & \mbox{if } p = 0 \\
a \tau^p, & \mbox {otherwise} \end{cases} c T^{t} c = \begin{cases}a, & \mbox{if } p = 0 \\ - a p (p-1) \left(T^{*}\right)^p & \mbox{otherwise}\end{cases}
t = -p
Planck-Einstein terms n_i \ln \left[1 - e^{-\gamma_i \tau}\right] \frac{b_i x_i^2 \exp\left(-x_i\right)}{[1-\exp\left(-x_i\right)]^2} where x_i = \frac{\beta_i}{T} b = n
\beta = T^{*}\gamma
'lead' terms \ln\left(\delta\right) + c + m \tau 1 none required

Note the factor '1' in the lead terms is not usually specified separately in publications; our code compensates for that when reading in c_p^0 data; you don't need to subtract 1 from the provided constant term.

If a publication provides \phi^0, then the above allows calculation of c_p^0.

If on the other hand a publication provides c_p^0 then we need to use the reference state properties \left(T_0,\rho_0,h_0^0, s_0^0\right), to define c and m, as follows:

m = \frac{h_0^0}{R T_c}
c = -\frac{s_0^0}{R} - 1 - \ln\delta_0 + ln\tau_0

where \tau_0 =  T^{*}/T_0 and \delta_0 = \rho_0 / \rho^{*}.

However, if we are using this ideal correlation together with a 'real' part, then these constants c and m will instead take values determined from reference state properties \left(T_0,\rho_0,h_0, s_0\right), meaning that these expressions will be offset as a result of h_0^0 - h_0 and s_0^0 - s_0 variations.

Residual part

For the residual part of the reduced Helmholtz energy, FPROPS currently supports three different kinds of terms:

type expression correlation constants
power terms 

\begin{cases} 
  a_i \tau^{t_i} \delta^{d_i} \exp\left(-\delta^{l_i}\right), & \mbox{if }l \ne 0 \\
  a_i \tau^{t_i} \delta^{d_i},  & \mbox{otherwise} 
\end{cases}


a, t, d, l
gaussian terms n_i \tau^{t_i} \delta^{d_i} \exp\left[-\alpha_i(\delta - \epsilon_i)^2 - \beta_i(\tau - \gamma_i)^2\right]  n, t, d, \alpha, \beta, \gamma, \epsilon
critical terms n_i \Delta^{b_i} \delta \psi \mbox{  where } \begin{cases}
  \Delta = \theta^2 + B_i \left[ \left( \delta - 1 \right)^2 \right]^{a_i} \\
  \theta = \left( 1- \tau \right) + A_i \left[ \left( \delta - 1 \right)^2 \right]^{\frac{1}{2 \beta_i}} \\
  \psi = \exp \left[-C_i \left(\delta-1 \right)^2 - D_i \left( \tau - 1 \right)^2 \right]
\end{cases} n, a, b, \beta, A, B, C, D

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 \phi \left(\tau,\delta \right) = \phi^0 \left(\tau,\delta \right) + \phi^r \left( \tau,\delta \right) and hence the specific Helmholtz energy a = R T \phi \left(\tau,\delta\right).

From the expression for \phi, we can then derive expressions for p, h, u, s and their partial derivatives with respect to density and temperature. These expressions are summarised in IAPWS-95[1] as well as by Tillner Roth et al[4].

Saturated mixture properties

The Helmholtz function a \left(\rho,T \right), given above, applies for all regions outside the saturation region. The saturation region, on the other hand, is not calculated directly from the raw Helmholtz correlation. Instead within the saturation region, the stable equilibrium of phases must be calculated such that the liquid phase vapour phases have equal Gibbs free energy g = a + p / \rho. This requirement is equivalent to the Maxwell equal area rule[5][1][6]. The condition can be expressed as a set of three simultaneous equations in the variables p_{sat}, \rho_f, \rho_g:


\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}

This can be further reduced to :


\begin{align}
p_{sat} & = R T \rho_f \left(1 + \delta_f \phi^r_\delta \left(\tau,\delta_f \right) \right) \\
p_{sat} & = R T \rho_g \left(1 + \delta_g \phi^r_\delta \left(\tau,\delta_g, \right) \right) \\
\frac{p_{sat}}{R T} \left( \frac{1}{\rho_g} - \frac{1}{\rho_f} \right) & = \phi \left(\tau, \delta_f \right) - \phi \left(\tau, \delta_g \right)
 = \ln \left( \frac{\rho_f}{\rho_g} \right) + \phi^r \left(\tau, \delta_f \right) - \phi^r \left(\tau, \delta_g  \right)
\end{align}

Calculation of the phase equilibrium conditions requires iterative solution of these three equations; one suitable approach is that of Akasaka[7], which first eliminates p_{sat} from the above equations, then iterates only on \rho_f and \rho_g. The approach seems to be robust, with tests currently performed with a sample of five fluids including water, carbon-dioxide and nitrogen. Our code is in models/johnpye/fprops/sat.c.

The Akasaka approach is markedly simpler to implement that the complicated algorithm used in REFPROP (which is explained by Akasaka, and also by Soave[8]). A limitation of the Soave approach is that it does not converge for temperatures T \gtrapprox 0.99 T_c, and hence requires yet another, even more complicated, algorithm to be implemented for that region (including need for the very tedious-to-calculate third partial derivatives \phi^r_{\delta \delta \delta}!).

Saturation Tsat(p)

To calculation the saturation temperature as a function of pressure, a different iteration scheme from the above for p_{sat}(T) is needed. [9][10]

At present, we have implemented fprops_sat_p in models/johnpye/fprops/sat.c which provides a solution via simple Newton iteration that makes repeated calls to fprops_sat_T. The iteration makes use of the Clapeyron equation to calculate the gradients,

\left(\frac{\partial p}{\partial T} \right)_{sat} = \frac{h_{fg}}{T v_{fg}}

There is probably a more efficient approach, but this is adequate for the moment.

Thermodynamic property derivatives

This page documents an experimental feature. Please tell us if you experience any problems.

The IAPWS has made a release[11] with information about how partial derivatives of thermodynamic properties with respect to other properties can be calculated in general.

FPROPS includes an initial implementation of these derivative calculations in models/johnpye/fprops/derivs.c, but they have not yet been thoroughly tested.

Properties in terms of (p,h)

This page documents an experimental feature. Please tell us if you experience any problems.

It seems to be very helpful to have property correlations calculable in terms of the independent variables (p,h). This seems to aid stability of calculation in larger systems, compared to the (v,T) property pair. This would appear to be in part due to the very small, steep subcooled region when using (v,T) coordinates.

An initial implementation of an iterative solver for properties in terms of (p,h) has been implemented for FPROPS in models/johnpye/fprops/solve_ph.c. Currently, it solves for all regions except for pressures below the triple-point pressure. We hope to overcome this problem soon.

Properties in terms of (T,x)

For convenience some simple functions for calculating properties in terms of (T,x) have been written. The implementation is in models/johnpye/fprops/solve_Tx.c. Any input with x outside [0,1] or T outside \left[ T_t, T_c \right] will return an error code. Otherwise, the corresponding rho (\rho) value will be returned such that (T,\rho) is equivalent to the desired state.

A python demonstration of this code is in models/johnpye/fprops/python/solve_Tx.py, which gives the following output:

john@thunder:~/ascend/models/johnpye/fprops/python$ python solve_Tx.py 
solved state
T = 300.000000, x = 0.500000
--> rho = 0.309802, h = 395476.135716, s = 381991.898474, u = 1491.770280

checking reverse solve:
--> x = 0.500000

Transport properties

Initial implementations of dynamic viscosity and thermal conductivity have been created for carbon dioxide. Providing the fluid data supports this, you should be able to access these properties from Python, C and ASCEND. In Python,

from fprops import *
D = fluid("carbondioxide")
S = D.set_Trho(220,2.5)
print "mu =",S.mu
print "lam =",S.lam

Both properties are returned in base SI units (Pa*s for viscosity, and W/m·K for thermal conductvity).

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 (should be) embedded at the end of the corresponding fluid file, for example models/johnpye/fprops/fluids/hydrogen.c.

During development, some additional testing routines are in the python subdirectory. These include

  • satcvgc.py to test converge of saturation region routines fprops_sat_T and fprops_sat_p are working.
  • solve_ph_array to test convergence of the fprops_solve_ph routine, including plotting of regions where the routine is failing.

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/fluids/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/fluids/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

\frac{\Delta h}{RT} = \tau \Delta \phi_{\tau}^0 = T^{*} \Delta a_2^0
\frac{\Delta s}{R} = \tau \Delta \phi_{\tau}^0 -  \Delta \phi^0 = - \Delta a_1^0

where we assume

\phi^0(\tau,\delta) = a_1^0 + a_2^0 \tau + \mathrm{other\ terms\ in\ \tau, \delta}

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 \phi^0(\tau,\delta).

Using within ASCEND

Note that the interface for accessing FPROPS is being changed significantly in 2012, in the 'fprops' branch code and work by Sean Muratet.

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, to calculate p(T,rho)
  • helmholtz_u, to calculate u(T,rho)
  • helmholtz_h, to calculate h(T,rho)
  • helmholtz_s, to calculate s(T,rho)
  • helmholtz_a, to calculate a(T,rho)
  • helmholtz_cp, to calculate cp(T,rho)
  • helmholtz_cv, to calculate cv(T,rho)
  • helmholtz_w, to calculate speed of sound for given (T,rho)
  • helmholtz_phsx_vT is also available which returns p, h, s and x in for specified (v,T).
  • helmholtz_Tvsx_ph returns T, v, s and x for specified (p,h). Recommended for energy system models.

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/fprops/rankine_fprops.a4c. The key part is

cd IS_A fluid;
cd.component :== 'water';
(* ...definitions of p, h, T, v, s, x omitted here... *)

calc_ph: helmholtz_Tvsx_ph(
    p, h : INPUT;
    T, v, s, x : OUTPUT;
    cd : DATA
);

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 once the relevant header files have been examined. Warning! We are still working on the API and it is subject to change, especially in the light of current GSOC2011 work.

A very simple test of selecting a fluid, then solving for (p,h) values and then returning results, is given in models/johnpye/fprops/test_ph.c (Note: this file is currently missing and needs to be restored.). An extract of that file is shown below:

double p = 1e5;
double h = 3500e3;const char *fluid = "water";
double T = 0,rho = 0;
const HelmholtzData *D = fprops_fluid(fluid);
int res = fprops_solve_ph(p,h,&T, &rho, 0, D);
if(res)printf("ERROR %d!)\n",res);
printf("Result: T = %f K, rho = %f kg/m3\n",T,rho);
double p1, h1;
p1 = helmholtz_p(T,rho, D);
h1 = helmholtz_h(T,rho, D);
printf("Double check: h = %f kJ/kg\n",h1/1e3);

Other test code written in C is shown in models/johnpye/fprops/test.c, which may be helpful for seeing some other features in use.

At this stage no external documentation for the use of the C code has been written, other that what is on this page.

Using from Python

FPROPS can also be used from the Python language. There are a number of sample programs in models/johnpye/fprops/python. The API is still evolving as we seek the cleanest ways to report errors back to the user, and catch floating-point errors, etc.

Further work

This article is about planned development or proposed functionality. Comments welcome.

At the conclusion of GSOC2012 with Sean Muratet, we completed implementation of a cubic EOS (Peng Robinson) as well as heavily-rewritten internal fluids database and data structures. Work continues in Sean's branch as well as the fprops2:models/johnpye/fprops branch.

At the end of GSOC2011 with Richard Towers, we concluded the following steps are required using his richard: branch as the basis. Further work is now (Aug 2011) goin on in the fprops2:models/johnpye/fprops branch.

  1. checking/evaluatingthe proposed data structure (is the design there right?) and modifying if required
  2. fixing up broken ASCEND bindings
  3. adding more fluid data, eg from critical point/omega data from RPP
  4. adding test code for PR EOS-based evaluations
  5. fixing up the PR EOS evaluation code until the test code works
  6. merging to trunk
  7. extending python bindings to raise exceptions instead of returning error flags (nicer API design)
  8. performance testing, profiling, etc.

FPROPS is still under development! We aim to implement a few more features including

  • Adding incompressible fluids EOS
  • Adding additional cubic EOS
  • Adding option for interpolation of saturation curve (as seen in CoolProp)
  • Extend test suite with more comparison of RPP data and cubic EOS code
  • Revisit derivatives calculation for new EOSes.
  • Expose derivatives code to ASCEND and Python; extend ASCEND to support optimisation with external relations.
  • add support for calculation of the saturation curve using Maxwell phase-equilibrium condition
  • add iterative solver for properties in terms of (p,h)
  • 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/drho to 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. Third law of thermodynamics? How to avoid correlations with negative h or s values?
  • implement a pre-calculation routine where for example coefficients of \phi^0 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).
  • add support for solving in terms of (p,h) where p < p_t, the triple point pressure.
  • possibly some additional coefficients needs in order to support correlations such as for Methanol[12]
  • add properties of R-134a (see doi:10.1063/1.555958)

Further down the track we would like to support:

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, subject to the restrictions of the GPL license.

References

  1. 1.0 1.1 1.2 http://www.iapws.org/relguide/IAPWS95-Rev.pdf
  2. E W Lemmon and R Span, 2006. Short Fundamental Equations of State for 20 Industrial Fluids, J. Chem. Eng. Data 51, doi:10.1021/je050186n
  3. 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, http://link.aip.org/link/?jpr/25/1509%26agg=NIST, doi:10.1063/1.555991
  4. Tillner-Roth, Harms-Watzenberg and Baehr (1993), "Eine neue Fundamentalgleichung für Ammoniak", DKV-Tagungsbericht, 20, 167-181
  5. http://en.wikipedia.org/wiki/Van_der_Waals_equation
  6. Eubank & Hall (2004), "Equal area rule and algorithm for determining phase compositions", AIChE Journal, 41, doi:10.1002/aic.690410419
  7. 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, doi:10.1299/jtst.3.442
  8. Soave (1986), "Direct calculation of pure-compound vapour pressures through cubic equations of state ", Fluid Phase Equilibria, 31, 203-207, doi:10.1016/0378-3812(86)90013-0
  9. Long X. Nghiem and Yau-Kun Li, 1985, Application of the tangent plane criterion to saturation pressure and temperature computations, Fluid Phase Equilibria, doi:10.1016/0378-3812(85)90059-7
  10. M Michelson, 1984, Saturation point calculations, Fluid Phase Equilibria, doi:10.1016/0378-3812(85)90005-6
  11. IAPWS, 2008, Revised Advisory Note No. 3 Thermodynamic Derivatives from IAPWS Formulation, http://www.iapws.org
  12. K de Reuck, 1993, Methanol, International Thermodynamic Tables of the Fluid State Vol. 12 (1993), Blackwell/IUPAC, http://www.iupac.org/objID/Source/sou65427568843800267934799

Contributors

The people who have contributed code and/or implementable advice for FPROPS:

See Also