Density/Volume (chemicals.volume)

This module contains various volume/density estimation routines, dataframes of fit coefficients, and mixing rules.

For reporting bugs, adding feature requests, or submitting pull requests, please use the GitHub issue tracker.

Pure Low Pressure Liquid Correlations

chemicals.volume.Rackett(T, Tc, Pc, Zc)[source]

Calculates saturation liquid volume, using Rackett CSP method and critical properties.

The molar volume of a liquid is given by:

Vs=RTcPcZc[1+(1T/Tc)2/7]V_s = \frac{RT_c}{P_c}{Z_c}^{[1+(1-{T/T_c})^{2/7} ]}

Units are all currently in m^3/mol - this can be changed to kg/m^3

Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

Pcfloat

Critical pressure of fluid [Pa]

Zcfloat

Critical compressibility of fluid, [-]

Returns
Vsfloat

Saturation liquid volume, [m^3/mol]

Notes

According to Reid et. al, underpredicts volume for compounds with Zc < 0.22

References

1

Rackett, Harold G. “Equation of State for Saturated Liquids.” Journal of Chemical & Engineering Data 15, no. 4 (1970): 514-517. doi:10.1021/je60047a012

Examples

Propane, example from the API Handbook

>>> from chemicals.utils import Vm_to_rho
>>> Vm_to_rho(Rackett(272.03889, 369.83, 4248000.0, 0.2763), 44.09562)
531.3221411755724
chemicals.volume.COSTALD(T, Tc, Vc, omega)[source]

Calculate saturation liquid density using the COSTALD CSP method.

A popular and accurate estimation method. If possible, fit parameters are used; alternatively critical properties work well.

The density of a liquid is given by:

Vs=VV(0)[1ωSRKV(δ)]V_s=V^*V^{(0)}[1-\omega_{SRK}V^{(\delta)}]
V(0)=11.52816(1Tr)1/3+1.43907(1Tr)2/30.81446(1Tr)+0.190454(1Tr)4/3V^{(0)}=1-1.52816(1-T_r)^{1/3}+1.43907(1-T_r)^{2/3} - 0.81446(1-T_r)+0.190454(1-T_r)^{4/3}
V(δ)=0.296123+0.386914Tr0.0427258Tr20.0480645Tr3Tr1.00001V^{(\delta)}=\frac{-0.296123+0.386914T_r-0.0427258T_r^2-0.0480645T_r^3} {T_r-1.00001}

Units are that of critical or fit constant volume.

Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

Vcfloat

Critical volume of fluid [m^3/mol]. This parameter is alternatively a fit parameter

omegafloat

(ideally SRK) Acentric factor for fluid, [-] This parameter is alternatively a fit parameter.

Returns
Vsfloat

Saturation liquid volume

Notes

196 constants are fit to this function in [1]. Range: 0.25 < Tr < 0.95, often said to be to 1.0

This function has been checked with the API handbook example problem.

References

1

Hankinson, Risdon W., and George H. Thomson. “A New Correlation for Saturated Densities of Liquids and Their Mixtures.” AIChE Journal 25, no. 4 (1979): 653-663. doi:10.1002/aic.690250412

Examples

Propane, from an example in the API Handbook:

>>> from chemicals.utils import Vm_to_rho
>>> Vm_to_rho(COSTALD(272.03889, 369.83333, 0.20008161E-3, 0.1532), 44.097)
530.3009967969844
chemicals.volume.Yen_Woods_saturation(T, Tc, Vc, Zc)[source]

Calculates saturation liquid volume, using the Yen and Woods [1] CSP method and a chemical’s critical properties.

The molar volume of a liquid is given by:

Vc/Vs=1+A(1Tr)1/3+B(1Tr)2/3+D(1Tr)4/3Vc/Vs = 1 + A(1-T_r)^{1/3} + B(1-T_r)^{2/3} + D(1-T_r)^{4/3}
D=0.93BD = 0.93-B
A=17.4425214.578Zc+989.625Zc21522.06Zc3A = 17.4425 - 214.578Z_c + 989.625Z_c^2 - 1522.06Z_c^3
B=3.28257+13.6377Zc+107.4844Zc2384.211Zc3 if Zc0.26B = -3.28257 + 13.6377Z_c + 107.4844Z_c^2-384.211Z_c^3 \text{ if } Zc \le 0.26
B=60.2091402.063Zc+501.0Zc2+641.0Zc3 if Zc0.26B = 60.2091 - 402.063Z_c + 501.0 Z_c^2 + 641.0 Z_c^3 \text{ if } Zc \ge 0.26
Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

Vcfloat

Critical volume of fluid [m^3/mol]

Zcfloat

Critical compressibility of fluid, [-]

Returns
Vsfloat

Saturation liquid volume, [m^3/mol]

Notes

Original equation was in terms of density, but it is converted here.

No example has been found, nor are there points in the article. However, it is believed correct. For compressed liquids with the Yen-Woods method, see the YenWoods_compressed function.

References

1

Yen, Lewis C., and S. S. Woods. “A Generalized Equation for Computer Calculation of Liquid Densities.” AIChE Journal 12, no. 1 (1966): 95-99. doi:10.1002/aic.690120119

Examples

>>> Yen_Woods_saturation(300, 647.14, 55.45E-6, 0.245)
1.769533076529574e-05
chemicals.volume.Yamada_Gunn(T, Tc, Pc, omega)[source]

Calculates saturation liquid volume, using Yamada and Gunn CSP method and a chemical’s critical properties and acentric factor.

The molar volume of a liquid is given by:

Vs=RTcPc(0.290560.08775ω)[1+(1T/Tc)2/7]V_s = \frac{RT_c}{P_c}{(0.29056-0.08775\omega)}^{[1+(1-{T/T_c})^{2/7}]}

Units are in m^3/mol.

Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

Pcfloat

Critical pressure of fluid [Pa]

omegafloat

Acentric factor for fluid, [-]

Returns
Vsfloat

saturation liquid volume, [m^3/mol]

Notes

This equation is an improvement on the Rackett equation. This is often presented as the Rackett equation. The acentric factor is used here, instead of the critical compressibility A variant using a reference fluid also exists

References

1

Gunn, R. D., and Tomoyoshi Yamada. “A Corresponding States Correlation of Saturated Liquid Volumes.” AIChE Journal 17, no. 6 (1971): 1341-45. doi:10.1002/aic.690170613

2

Yamada, Tomoyoshi, and Robert D. Gunn. “Saturated Liquid Molar Volumes. Rackett Equation.” Journal of Chemical & Engineering Data 18, no. 2 (1973): 234-36. doi:10.1021/je60057a006

Examples

>>> Yamada_Gunn(300, 647.14, 22048320.0, 0.245)
2.188284384699659e-05
chemicals.volume.Townsend_Hales(T, Tc, Vc, omega)[source]

Calculates saturation liquid density, using the Townsend and Hales CSP method as modified from the original Riedel equation. Uses chemical critical volume and temperature, as well as acentric factor

The density of a liquid is given by:

Vs=Vc/(1+0.85(1Tr)+(1.692+0.986ω)(1Tr)1/3)Vs = V_c/\left(1+0.85(1-T_r)+(1.692+0.986\omega)(1-T_r)^{1/3}\right)
Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

Vcfloat

Critical volume of fluid [m^3/mol]

omegafloat

Acentric factor for fluid, [-]

Returns
Vsfloat

Saturation liquid volume, [m^3/mol]

Notes

The requirement for critical volume and acentric factor requires all data.

References

1

Hales, J. L, and R Townsend. “Liquid Densities from 293 to 490 K of Nine Aromatic Hydrocarbons.” The Journal of Chemical Thermodynamics 4, no. 5 (1972): 763-72. doi:10.1016/0021-9614(72)90050-X

Examples

>>> Townsend_Hales(300, 647.14, 55.95E-6, 0.3449)
1.8007361992619923e-05
chemicals.volume.Bhirud_normal(T, Tc, Pc, omega)[source]

Calculates saturation liquid density using the Bhirud [1] CSP method. Uses Critical temperature and pressure and acentric factor.

The density of a liquid is given by:

lnPcρRT=lnU(0)+ωlnU(1)\ln \frac{P_c}{\rho RT} = \ln U^{(0)} + \omega\ln U^{(1)}
lnU(0)=1.3964424.076Tr+102.615Tr2255.719Tr3+355.805Tr4256.671Tr5+75.1088Tr6\ln U^{(0)} = 1.396 44 - 24.076T_r+ 102.615T_r^2 -255.719T_r^3+355.805T_r^4-256.671T_r^5 + 75.1088T_r^6
lnU(1)=13.4412135.7437Tr+533.380Tr21091.453Tr3+1231.43Tr4728.227Tr5+176.737Tr6\ln U^{(1)} = 13.4412 - 135.7437 T_r + 533.380T_r^2- 1091.453T_r^3+1231.43T_r^4 - 728.227T_r^5 + 176.737T_r^6
Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

Pcfloat

Critical pressure of fluid [Pa]

omegafloat

Acentric factor for fluid, [-]

Returns
Vmfloat

Saturated liquid molar volume, [mol/m^3]

Notes

Claimed inadequate by others.

An interpolation table for ln U values are used from Tr = 0.98 - 1.000. Has terrible behavior at low reduced temperatures.

References

1

Bhirud, Vasant L. “Saturated Liquid Densities of Normal Fluids.” AIChE Journal 24, no. 6 (November 1, 1978): 1127-31. doi:10.1002/aic.690240630

Examples

Pentane

>>> Bhirud_normal(280.0, 469.7, 33.7E5, 0.252)
0.00011249657842514176
chemicals.volume.Campbell_Thodos(T, Tb, Tc, Pc, MW, dipole=0.0, has_hydroxyl=False)[source]

Calculate saturation liquid density using the Campbell-Thodos [1] CSP method.

An old and uncommon estimation method.

Vs=RTcPcZRA[1+(1Tr)2/7]V_s = \frac{RT_c}{P_c}{Z_{RA}}^{[1+(1-T_r)^{2/7}]}
ZRA=α+β(1Tr)Z_{RA} = \alpha + \beta(1-T_r)
α=0.38830.0179s\alpha = 0.3883-0.0179s
s=TbrlnPc(1Tbr)s = T_{br} \frac{\ln P_c}{(1-T_{br})}
β=0.00318s0.0211+0.625Λ1.35\beta = 0.00318s-0.0211+0.625\Lambda^{1.35}
Λ=Pc1/3MW1/2Tc5/6\Lambda = \frac{P_c^{1/3}} { MW^{1/2} T_c^{5/6}}

For polar compounds:

θ=Pcμ2/Tc2\theta = P_c \mu^2/T_c^2
α=0.38830.0179s130540θ2.41\alpha = 0.3883 - 0.0179s - 130540\theta^{2.41}
β=0.00318s0.0211+0.625Λ1.35+9.74×106θ3.38\beta = 0.00318s - 0.0211 + 0.625\Lambda^{1.35} + 9.74\times 10^6 \theta^{3.38}

Polar Combounds with hydroxyl groups (water, alcohols)

α=[0.690Tbr0.3342+5.79×1010Tbr32.75]Pc0.145\alpha = \left[0.690T_{br} -0.3342 + \frac{5.79\times 10^{-10}} {T_{br}^{32.75}}\right] P_c^{0.145}
β=0.00318s0.0211+0.625Λ1.35+5.90Θ0.835\beta = 0.00318s - 0.0211 + 0.625 \Lambda^{1.35} + 5.90\Theta^{0.835}
Parameters
Tfloat

Temperature of fluid [K]

Tbfloat

Boiling temperature of the fluid [K]

Tcfloat

Critical temperature of fluid [K]

Pcfloat

Critical pressure of fluid [Pa]

MWfloat

Molecular weight of the fluid [g/mol]

dipolefloat, optional

Dipole moment of the fluid [debye]

has_hydroxylbool, optional

Swith to use the hydroxyl variant for polar fluids

Returns
Vsfloat

Saturation liquid volume, [m^3/mol]

Notes

If a dipole is provided, the polar chemical method is used. The paper is an excellent read. Pc is internally converted to atm.

References

1(1,2)

Campbell, Scott W., and George Thodos. “Prediction of Saturated Liquid Densities and Critical Volumes for Polar and Nonpolar Substances.” Journal of Chemical & Engineering Data 30, no. 1 (January 1, 1985): 102-11. doi:10.1021/je00039a032.

Examples

Ammonia, from [1].

>>> Campbell_Thodos(T=405.45, Tb=239.82, Tc=405.45, Pc=111.7*101325, MW=17.03, dipole=1.47)
7.347366126245e-05
chemicals.volume.SNM0(T, Tc, Vc, omega, delta_SRK=None)[source]

Calculates saturated liquid density using the Mchaweh, Moshfeghian model [1]. Designed for simple calculations.

Vs=Vc/(1+1.169τ1/3+1.818τ2/32.658τ+2.161τ4/3V_s = V_c/(1+1.169\tau^{1/3}+1.818\tau^{2/3}-2.658\tau+2.161\tau^{4/3}
τ=1(T/Tc)αSRK\tau = 1-\frac{(T/T_c)}{\alpha_{SRK}}
αSRK=[1+m(1T/TC]2\alpha_{SRK} = [1 + m(1-\sqrt{T/T_C}]^2
m=0.480+1.574ω0.176ω2m = 0.480+1.574\omega-0.176\omega^2

If the fit parameter delta_SRK is provided, the following is used:

Vs=VC/(1+1.169τ1/3+1.818τ2/32.658τ+2.161τ4/3)/[1+δSRK(αSRK1)1/3]V_s = V_C/(1+1.169\tau^{1/3}+1.818\tau^{2/3}-2.658\tau+2.161\tau^{4/3}) /\left[1+\delta_{SRK}(\alpha_{SRK}-1)^{1/3}\right]
Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

Vcfloat

Critical volume of fluid [m^3/mol]

omegafloat

Acentric factor for fluid, [-]

delta_SRKfloat, optional

Fitting parameter [-]

Returns
Vsfloat

Saturation liquid volume, [m^3/mol]

Notes

73 fit parameters have been gathered from the article.

References

1

Mchaweh, A., A. Alsaygh, Kh. Nasrifar, and M. Moshfeghian. “A Simplified Method for Calculating Saturated Liquid Densities.” Fluid Phase Equilibria 224, no. 2 (October 1, 2004): 157-67. doi:10.1016/j.fluid.2004.06.054

Examples

Argon, without the fit parameter and with it. Tabulated result in Perry’s is 3.4613e-05. The fit increases the error on this occasion.

>>> SNM0(121, 150.8, 7.49e-05, -0.004)
3.440225640273e-05
>>> SNM0(121, 150.8, 7.49e-05, -0.004, -0.03259620)
3.493288100008e-05

Pure High Pressure Liquid Correlations

chemicals.volume.COSTALD_compressed(T, P, Psat, Tc, Pc, omega, Vs)[source]

Calculates compressed-liquid volume, using the COSTALD [1] CSP method and a chemical’s critical properties.

The molar volume of a liquid is given by:

V=Vs(1ClnB+PB+Psat)V = V_s\left( 1 - C \ln \frac{B + P}{B + P^{sat}}\right)
BPc=1+aτ1/3+bτ2/3+dτ+eτ4/3\frac{B}{P_c} = -1 + a\tau^{1/3} + b\tau^{2/3} + d\tau + e\tau^{4/3}
e=exp(f+gωSRK+hωSRK2)e = \exp(f + g\omega_{SRK} + h \omega_{SRK}^2)
C=j+kωSRKC = j + k \omega_{SRK}
Parameters
Tfloat

Temperature of fluid [K]

Pfloat

Pressure of fluid [Pa]

Psatfloat

Saturation pressure of the fluid [Pa]

Tcfloat

Critical temperature of fluid [K]

Pcfloat

Critical pressure of fluid [Pa]

omegafloat

(ideally SRK) Acentric factor for fluid, [-] This parameter is alternatively a fit parameter.

Vsfloat

Saturation liquid volume, [m^3/mol]

Returns
V_densefloat

High-pressure liquid volume, [m^3/mol]

Notes

Original equation was in terms of density, but it is converted here.

The example is from DIPPR, and exactly correct. This is DIPPR Procedure 4C: Method for Estimating the Density of Pure Organic Liquids under Pressure.

References

1

Thomson, G. H., K. R. Brobst, and R. W. Hankinson. “An Improved Correlation for Densities of Compressed Liquids and Liquid Mixtures.” AIChE Journal 28, no. 4 (July 1, 1982): 671-76. doi:10.1002/aic.690280420

Examples

>>> COSTALD_compressed(303., 9.8E7, 85857.9, 466.7, 3640000.0, 0.281, 0.000105047)
9.287482879788505e-05

Liquid Mixing Rules

chemicals.volume.Amgat(xs, Vms)[source]

Calculate mixture liquid density using the Amgat mixing rule. Highly inacurate, but easy to use. Assumes idea liquids with no excess volume. Average molecular weight should be used with it to obtain density.

Vmix=ixiViV_{mix} = \sum_i x_i V_i

or in terms of density:

ρmix=xiρi\rho_{mix} = \sum\frac{x_i}{\rho_i}
Parameters
xsarray

Mole fractions of each component, []

Vmsarray

Molar volumes of each fluids at conditions [m^3/mol]

Returns
Vmfloat

Mixture liquid volume [m^3/mol]

Notes

Units are that of the given volumes. It has been suggested to use this equation with weight fractions, but the results have been less accurate.

Examples

>>> Amgat([0.5, 0.5], [4.057e-05, 5.861e-05])
4.9590000000000005e-05
chemicals.volume.Rackett_mixture(T, xs, MWs, Tcs, Pcs, Zrs)[source]

Calculate mixture liquid density using the Rackett-derived mixing rule as shown in [2].

Vm=ixiTciMWiPciZR,m(1+(1Tr)2/7)RixiMWiV_m = \sum_i\frac{x_i T_{ci}}{MW_i P_{ci}} Z_{R,m}^{(1 + (1 - T_r)^{2/7})} R \sum_i x_i MW_i
Parameters
Tfloat

Temperature of liquid [K]

xs: list

Mole fractions of each component, []

MWslist

Molecular weights of each component [g/mol]

Tcslist

Critical temperatures of each component [K]

Pcslist

Critical pressures of each component [Pa]

Zrslist

Rackett parameters of each component []

Returns
Vmfloat

Mixture liquid volume [m^3/mol]

Notes

Model for pure compounds in [1] forms the basis for this model, shown in [2]. Molecular weights are used as weighing by such has been found to provide higher accuracy in [2]. The model can also be used without molecular weights, but results are somewhat different.

As with the Rackett model, critical compressibilities may be used if Rackett parameters have not been regressed.

Critical mixture temperature, and compressibility are all obtained with simple mixing rules.

References

1

Rackett, Harold G. “Equation of State for Saturated Liquids.” Journal of Chemical & Engineering Data 15, no. 4 (1970): 514-517. doi:10.1021/je60047a012

2(1,2,3,4)

Danner, Ronald P, and Design Institute for Physical Property Data. Manual for Predicting Chemical Process Design Data. New York, N.Y, 1982.

Examples

Calculation in [2] for methanol and water mixture. Result matches example.

>>> Rackett_mixture(T=298., xs=[0.4576, 0.5424], MWs=[32.04, 18.01], Tcs=[512.58, 647.29], Pcs=[8.096E6, 2.209E7], Zrs=[0.2332, 0.2374])
2.6252894930056885e-05
chemicals.volume.COSTALD_mixture(xs, T, Tcs, Vcs, omegas)[source]

Calculate mixture liquid density using the COSTALD CSP method.

A popular and accurate estimation method. If possible, fit parameters are used; alternatively critical properties work well.

The mixing rules giving parameters for the pure component COSTALD equation are:

Tcm=ijxixj(VijTcij)VmT_{cm} = \frac{\sum_i\sum_j x_i x_j (V_{ij}T_{cij})}{V_m}
Vm=0.25[xiVi+3(xiVi2/3)(ixiVi1/3)]V_m = 0.25\left[ \sum x_i V_i + 3(\sum x_i V_i^{2/3})(\sum_i x_i V_i^{1/3})\right]
VijTcij=(ViTciVjTcj)0.5V_{ij}T_{cij} = (V_iT_{ci}V_{j}T_{cj})^{0.5}
ω=iziωi\omega = \sum_i z_i \omega_i
Parameters
xs: list

Mole fractions of each component

Tfloat

Temperature of fluid [K]

Tcslist

Critical temperature of fluids [K]

Vcslist

Critical volumes of fluids [m^3/mol]. This parameter is alternatively a fit parameter

omegaslist

(ideally SRK) Acentric factor of all fluids, [-] This parameter is alternatively a fit parameter.

Returns
Vsfloat

Saturation liquid mixture volume

Notes

Range: 0.25 < Tr < 0.95, often said to be to 1.0 No example has been found. Units are that of critical or fit constant volume.

References

1

Hankinson, Risdon W., and George H. Thomson. “A New Correlation for Saturated Densities of Liquids and Their Mixtures.” AIChE Journal 25, no. 4 (1979): 653-663. doi:10.1002/aic.690250412

Examples

>>> COSTALD_mixture([0.4576, 0.5424], 298.,  [512.58, 647.29], [0.000117, 5.6e-05], [0.559,0.344])
2.7065887732713534e-05

Gas Correlations

Gas volumes are predicted with one of:

  1. An equation of state

  2. A virial coefficient model

  3. The ideal gas law

Equations of state do much more than predict volume however. An implementation of many of them can be found in thermo.

Virial functions are implemented in chemicals.virial.

chemicals.volume.ideal_gas(T, P)[source]

Calculates ideal gas molar volume. The molar volume of an ideal gas is given by:

V=RTPV = \frac{RT}{P}
Parameters
Tfloat

Temperature of fluid [K]

Pfloat

Pressure of fluid [Pa]

Returns
Vfloat

Gas volume, [m^3/mol]

Examples

>>> ideal_gas(298.15, 101325.)
0.024465403697038125

Pure Solid Correlations

Solid density does not depend on pressure significantly, and unless operating in the geochemical or astronomical domain is normally neglected.

chemicals.volume.Goodman(T, Tt, Vml)[source]

Calculates solid density at T using the simple relationship by a member of the DIPPR.

The molar volume of a solid is given by:

1Vm=(1.280.16TTt)1VmL(Tt)\frac{1}{V_m} = \left( 1.28 - 0.16 \frac{T}{T_t}\right) \frac{1}{{Vm}_L(T_t)}
Parameters
Tfloat

Temperature of fluid [K]

Ttfloat

Triple temperature of fluid [K]

Vmlfloat

Liquid molar volume of the organic liquid at the triple point, [m^3/mol]

Returns
Vmsfloat

Solid molar volume, [m^3/mol]

Notes

Works to the next solid transition temperature or to approximately 0.3Tt.

References

1

Goodman, Benjamin T., W. Vincent Wilding, John L. Oscarson, and Richard L. Rowley. “A Note on the Relationship between Organic Solid Density and Liquid Density at the Triple Point.” Journal of Chemical & Engineering Data 49, no. 6 (2004): 1512-14. doi:10.1021/je034220e.

Examples

Decane at 200 K:

>>> Goodman(200, 243.225, 0.00023585)
0.0002053665090860923

Pure Component Liquid Fit Correlations

chemicals.volume.Rackett_fit(T, Tc, rhoc, b, n, MW=None)[source]

Calculates saturation liquid volume, using the Rackett equation form and a known or estimated critical temperature and density as well as fit parameters b and n.

The density of a liquid is given by:

ρsat=ρcb(1TTc)n\rho_{sat} = \rho_c b^{-\left(1 - \frac{T}{T_c}\right)^n}

The density is then converted to a specific volume by taking its inverse.

Note that the units of this equation in some sources are kg/m^3, g/mL in others, and m^3/mol in others. If the units for the coefficients are in molar units, do NOT provide MW or an incorrect value will be returned. If the units are mass units and MW is not provided, the output will have the same units as rhoc.

Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

rhocfloat

Critical density of fluid, often a fit parameter only [kg/m^3]

bfloat

Fit parameter, [-]

nfloat

Fit parameter, [-]

MWfloat, optional

Molecular weight, [g/mol]

Returns
Vsfloat

Saturation liquid volume, [m^3/mol if MW given; m^3/kg otherwise]

References

1

Frenkel, Michael, Robert D. Chirico, Vladimir Diky, Xinjian Yan, Qian Dong, and Chris Muzny. “ThermoData Engine (TDE): Software Implementation of the Dynamic Data Evaluation Concept.” Journal of Chemical Information and Modeling 45, no. 4 (July 1, 2005): 816-38. https://doi.org/10.1021/ci050067b.

2

Yaws, Carl L. “Liquid Density of the Elements: A Comprehensive Tabulation for All the Important Elements from Ag to Zr.” Chemical Engineering 114, no. 12 (2007): 44-47.

Examples

Input sample from NIST (naphthalene) (m^3/kg):

>>> Rackett_fit(T=400.0, Tc=748.402, rhoc=314.629, b=0.257033, n=0.280338)
0.00106174320755

Parameters in Yaws form (butane) (note the 1000 multiplier on rhoc, called A in Yaws) (m^3/kg):

>>> Rackett_fit(T=298.15, Tc=425.18, rhoc=0.2283*1000, b=0.2724, n=0.2863)
0.00174520519958

Same Yaws point, with MW provided:

>>> Rackett_fit(T=298.15, Tc=425.18, rhoc=0.2283*1000, b=0.2724, n=0.2863, MW=58.123)
0.00010143656181
chemicals.volume.volume_VDI_PPDS(T, Tc, rhoc, a, b, c, d, MW=None)[source]

Calculates saturation liquid volume, using the critical properties and fitted coefficients from [1]. This is also known as the PPDS equation 10 or PPDS10.

ρmass=ρc+aτ0.35+bτ2/3+cτ+dτ4/3\rho_{mass} = \rho_{c} + a\tau^{0.35} + b \tau^{2/3} + c\tau + d\tau^{4/3}
Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

rhocfloat

Critical density of fluid [kg/m^3]

a,b,c,dfloat

Fitted coefficients [-]

MWfloat, optional

Molecular weight of chemical [g/mol]

Returns
Vsfloat

Saturation liquid molar volume or density, [m^3/mol if MW given; kg/m^3 otherwise]

References

1

Gesellschaft, V. D. I., ed. VDI Heat Atlas. 2nd edition. Berlin; New York:: Springer, 2010.

Examples

Calculate density of nitrogen in kg/m3 at 300 K:

>>> volume_VDI_PPDS(300, 126.19, 313, 470.922, 493.251, -560.469, 389.611)
313.0

Calculate molar volume of nitrogen in m3/mol at 300 K:

>>> volume_VDI_PPDS(300, 126.19, 313, 470.922, 493.251, -560.469, 389.611, 28.01)
8.9488817891e-05
chemicals.volume.TDE_VDNS_rho(T, Tc, rhoc, a1, a2, a3, a4, MW=None)[source]

Calculates saturation liquid volume, using the critical properties and fitted coefficients in the TDE VDNW form from [1].

ρmass=ρc+aτ0.35+bτ+cτ2+dτ3\rho_{mass} = \rho_{c} + a\tau^{0.35} + b \tau + c\tau^2 + d\tau^3
Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

rhocfloat

Critical density of fluid [kg/m^3]

a1float

Regression parameter, [-]

a2float

Regression parameter, [-]

a3float

Regression parameter, [-]

a4float

Regression parameter, [-]

MWfloat, optional

Molecular weight of chemical [g/mol]

Returns
Vsfloat

Saturation liquid molar volume or density, [m^3/mol if MW given; kg/m^3 otherwise]

References

1

“ThermoData Engine (TDE103b V10.1) User`s Guide.” https://trc.nist.gov/TDE/Help/TDE103b/Eqns-Pure-DensityLG/VDNSExpansion.htm.

Examples

>>> TDE_VDNS_rho(T=400.0, Tc=772.999, rhoc=320.037, a1=795.092, a2=-169.132, a3=448.929, a4=-102.931)
947.4906064903
chemicals.volume.PPDS17(T, Tc, a0, a1, a2, MW=None)[source]

Calculates saturation liquid volume, using the critical temperature and fitted coefficients in the PPDS17 form in [1].

ρmass=1a0(a1+a2τ)(1+τ2/7)\rho_{mass} = \frac{1}{a_0(a_1 + a_2\tau)^{\left(1 + \tau^{2/7} \right)}}
Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

a0float

Regression parameter, [-]

a1float

Regression parameter, [-]

a2float

Regression parameter, [-]

MWfloat, optional

Molecular weight of chemical [g/mol]

Returns
Vsfloat

Saturation liquid molar volume or density, [m^3/mol if MW given; kg/m^3 otherwise]

References

1(1,2)

“ThermoData Engine (TDE103b V10.1) User`s Guide.” https://trc.nist.gov/TDE/TDE_Help/Eqns-Pure-DensityLG/PPDS17.htm.

Examples

Coefficients for the liquid density of benzene from [1] at 300 K:

>>> PPDS17(300, 562.05, a0=0.0115508, a1=0.281004, a2=-0.00635447)
871.520087707

Pure Component High-Pressure Liquid Fit Correlations

chemicals.volume.Tait(P, P_ref, rho_ref, B, C)[source]

Calculates compressed-liquid mass density using the Tait model [1] and fit coefficients B and C and the reference (usually saturation) liquid density. B and C are normally temperature dependent but it is assumed they are constant (or calculated earlier) in this function

The mass density of the compressed liquid is given by:

ρ=ρref1ClnB+PB+Pref\rho = \frac{\rho_{ref}}{1 - C \ln \frac{B+P}{B + P_{ref}}}
Parameters
Pfloat

Pressure of fluid [Pa]

P_reffloat

Pressure of the fluid at the reference density; normally saturation at higher pressures and either 1 atm or 1 MPa at low enough temperatures the saturation pressure stops being an important factor, [Pa]

rho_reffloat

The mass density of the fluid at the reference condition, [kg/m^3]

Bfloat

Fit coefficient, [Pa]

Cfloat

Fit coefficient, [-]

Returns
rhofloat

High-pressure liquid mass density, [kg/m^3]

Notes

If P is set to be lower than P_ref, it is adjusted to have the same value as P_ref (saturation condition).

If B becomes negative and higher than P_ref, the logarithm will become undefined.

References

1(1,2)

Haynes, W.M., Thomas J. Bruno, and David R. Lide. CRC Handbook of Chemistry and Physics. [Boca Raton, FL]: CRC press, 2014.

Examples

Coefficients for methanol from the CRC Handbook [1] at 300 K and 1E8 Pa.

>>> B_coeffs = [649718000.0, -4345830.0, 13072.2, -20.029200000000003, 0.0120566]
>>> C_coeffs = [0.115068, -5.322e-05]
>>> T = 300
>>> Tait(P=1e8, P_ref=101325, rho_ref=784.85, B=np.polyval(B_coeffs[::-1], T), C=np.polyval(C_coeffs[::-1], T))
853.744916
chemicals.volume.Tait_molar(P, P_ref, V_ref, B, C)[source]

Calculates compressed-liquid volume using the Tait model [1] and fit coefficients B and C and the reference (usually saturation) liquid density. B and C are normally temperature dependent but it is assumed they are constant (or calculated earlier) in this function

The molar volume of the compressed liquid is given by:

Vm=Vref(1ClnB+PB+Pref)V_m = V_{ref}\cdot \left({1 - C \ln \frac{B+P}{B + P_{ref}}}\right)
Parameters
Pfloat

Pressure of fluid [Pa]

P_reffloat

Pressure of the fluid at the reference density; normally saturation at higher pressures and either 1 atm or 1 MPa at low enough temperatures the saturation pressure stops being an important factor, [Pa]

V_reffloat

The molar volume of the fluid at the reference condition, [m^3/mol]

Bfloat

Fit coefficient, [Pa]

Cfloat

Fit coefficient, [-]

Returns
V_densefloat

High-pressure liquid volume, [m^3/mol]

Notes

If P is set to be lower than P_ref, it is adjusted to have the same value as P_ref (saturation condition).

References

1(1,2)

Haynes, W.M., Thomas J. Bruno, and David R. Lide. CRC Handbook of Chemistry and Physics. [Boca Raton, FL]: CRC press, 2014.

Examples

Coefficients for methanol from the CRC Handbook [1] at 300 K and 1E8 Pa.

>>> Tait_molar(P=1e8, P_ref=101325.0, V_ref=4.0825e-05, B=79337060.0, C=0.099102)
3.75305e-05

Pure Component Solid Fit Correlations

chemicals.volume.CRC_inorganic(T, rho0, k, Tm, MW=None)[source]

Calculates liquid density of a molten element or salt at temperature above the melting point. Some coefficients are given nearly up to the boiling point.

The mass density of the inorganic liquid is given by:

ρ=ρ0k(TTm)\rho = \rho_{0} - k(T-T_m)
Parameters
Tfloat

Temperature of the liquid, [K]

rho0float

Mass density of the liquid at Tm, [kg/m^3]

kfloat

Linear temperature dependence of the mass density, [kg/m^3/K]

Tmfloat

The normal melting point, used in the correlation [K]

MWfloat, optional

Molecular weight of chemical [g/mol]

Returns
rhofloat

Mass density of molten metal or salt, [m^3/mol if MW given; kg/m^3 otherwise]

Notes

[1] has units of g/mL. While the individual densities could have been converted to molar units, the temperature coefficient could only be converted by refitting to calculated data. To maintain compatibility with the form of the equations, this was not performed.

This linear form is useful only in small temperature ranges. Coefficients for one compound could be used to predict the temperature dependence of density of a similar compound.

References

1

Haynes, W.M., Thomas J. Bruno, and David R. Lide. CRC Handbook of Chemistry and Physics, 95E. [Boca Raton, FL]: CRC press, 2014.

Examples

>>> CRC_inorganic(300, 2370.0, 2.687, 239.08)
2206.30796

Fit Coefficients

All of these coefficients are lazy-loaded, so they must be accessed as an attribute of this module.

chemicals.volume.rho_data_COSTALD

Coefficients for the COSTALD method from [3]; 192 fluids have coefficients published.

chemicals.volume.rho_data_SNM0

Coefficients for the SNM0 method for 73 fluids from [2].

chemicals.volume.rho_data_Perry_8E_105_l

Coefficients for chemicals.dippr.EQ105 from [1] for 344 fluids. Note this is in terms of molar density; to obtain molar volume, invert the result!

chemicals.volume.rho_data_VDI_PPDS_2

Coefficients in [5] developed by the PPDS using chemicals.dippr.EQ116 but in terms of mass density [kg/m^3]; Valid up to the critical temperature, and extrapolates to very low temperatures well.

chemicals.volume.rho_data_CRC_inorg_l

Single-temperature coefficient linear model in terms of mass density for the density of inorganic liquids. Data is available for 177 fluids normally valid over a narrow range above the melting point, from [4]; described in CRC_inorganic.

chemicals.volume.rho_data_CRC_inorg_l_const

Constant inorganic liquid molar volumes published in [4].

chemicals.volume.rho_data_CRC_inorg_s_const

Constant solid densities molar volumes published in [4].

chemicals.volume.rho_data_CRC_virial

Coefficients for a tempereture polynomial (T in Kelvin) for the second B virial coefficient published in [4]. The form of the equation is B=(a1+t(a2+t(a3+t(a4+a5t))))×106B = (a_1 + t(a_2 + t(a_3 + t(a_4 + a_5 t)))) \times 10^{-6} with t=298.15T1t = \frac{298.15}{T} - 1 and then B will be in units of m^3/mol.

1

Green, Don, and Robert Perry. Perry’s Chemical Engineers’ Handbook, 8E. McGraw-Hill Professional, 2007.

2

Mchaweh, A., A. Alsaygh, Kh. Nasrifar, and M. Moshfeghian. “A Simplified Method for Calculating Saturated Liquid Densities.” Fluid Phase Equilibria 224, no. 2 (October 1, 2004): 157-67. doi:10.1016/j.fluid.2004.06.054

3

Hankinson, Risdon W., and George H. Thomson. “A New Correlation for Saturated Densities of Liquids and Their Mixtures.” AIChE Journal 25, no. 4 (1979): 653-663. doi:10.1002/aic.690250412

4(1,2,3,4)

Haynes, W.M., Thomas J. Bruno, and David R. Lide. CRC Handbook of Chemistry and Physics. [Boca Raton, FL]: CRC press, 2014.

5

Gesellschaft, V. D. I., ed. VDI Heat Atlas. 2nd edition. Berlin; New York:: Springer, 2010.

The structure of each dataframe is shown below:

In [1]: import chemicals

In [2]: chemicals.volume.rho_data_COSTALD
Out[2]: 
                              Chemical  omega_SRK     Vchar    Z_RA
CAS                                                                
60-29-7                    ethyl ether     0.2800  0.000281  0.2632
64-17-5                  ethyl alcohol     0.6378  0.000175  0.2502
67-56-1                 methyl alcohol     0.5536  0.000120  0.2334
67-63-0              isopropyl alcohol     0.6637  0.000231  0.2493
67-64-1                        acetone     0.3149  0.000208  0.2477
...                                ...        ...       ...     ...
14752-75-1           heptadecylbenzene     0.9404  0.001146     NaN
30453-31-7    ethyl n-propyl disulfide     0.3876  0.000440  0.2662
33672-51-4  propyl isopropyl disulfide     0.4059  0.000502  0.2680
53966-36-2   ethyl isopropyl disulfide     0.3556  0.000439  0.2711
61828-04-4             tricosylbenzene     1.1399  0.001995     NaN

[192 rows x 4 columns]

In [3]: chemicals.volume.rho_data_SNM0
Out[3]: 
                             Chemical  delta_SRK
CAS                                             
56-23-5    Tetrachlouromethane, R-10   -0.013152
60-29-7                   Ethylether    0.001062
64-19-7                  Acetic acid   -0.010347
65-85-0                 Benzoic acid    0.026866
67-56-1                     Methanol    0.007195
...                               ...        ...
7727-37-9                   Nitrogen   -0.007946
7782-39-0                  Deuterium   -0.053345
7782-41-4                   Flourine   -0.030398
7782-44-7                     Oxygen   -0.027049
7782-50-5                   Chlorine    0.013010

[73 rows x 2 columns]

In [4]: chemicals.volume.rho_data_Perry_8E_105_l
Out[4]: 
                          Chemical       C1       C2  ...       C4    Tmin    Tmax
CAS                                                   ...                         
50-00-0              Formaldehyde   1941.50  0.22309  ...  0.28571  181.15  408.00
55-21-0                 Benzamide    737.10  0.25487  ...  0.28571  403.00  824.00
56-23-5      Carbon tetrachloride    998.35  0.27400  ...  0.28700  250.33  556.35
57-55-6      1,2-Propylene glycol   1092.30  0.26106  ...  0.20459  213.15  626.00
60-29-7             Diethyl ether    955.40  0.26847  ...  0.28140  156.85  466.70
...                            ...      ...      ...  ...      ...     ...     ...
10028-15-6                  Ozone   3359.20  0.29884  ...  0.28523   80.15  261.00
10035-10-6       Hydrogen bromide   2832.00  0.28320  ...  0.28571  185.15  363.15
10102-43-9           Nitric oxide   5246.00  0.30440  ...  0.24200  109.50  180.15
13511-13-2    Propenylcyclohexene    612.55  0.26769  ...  0.28571  199.00  636.00
132259-10-0                   Air   2896.30  0.26733  ...  0.27341   59.15  132.45

[344 rows x 7 columns]

In [5]: chemicals.volume.rho_data_VDI_PPDS_2
Out[5]: 
                         Chemical      MW  ...          C          D
CAS                                        ...                      
50-00-0              Formaldehyde   30.03  ...   245.3425    43.9601
56-23-5      Carbon tetrachloride  153.82  ...   535.7568   -28.0071
56-81-5                  Glycerol   92.09  ...  1429.7634  -527.7710
60-29-7             Diethyl ether   74.12  ...  -489.2726   486.7458
62-53-3                   Aniline   93.13  ...   242.0930     0.7157
...                           ...     ...  ...        ...        ...
10097-32-2                Bromine  159.82  ...   676.7593    15.3973
10102-43-9           Nitric oxide   30.01  ...  2252.1437 -1031.3210
10102-44-0       Nitrogen dioxide   46.01  ...  2233.6217  -968.0655
10544-72-6    Dinitrogentetroxide   92.01  ...   604.1720  -135.9384
132259-10-0                   Air   28.96  ...  -841.3265   495.5129

[272 rows x 8 columns]

In [6]: chemicals.volume.rho_data_CRC_inorg_l
Out[6]: 
                          Chemical       MW      rho      k       Tm     Tmax
CAS                                                                          
497-19-8          Sodium carbonate  105.989   1972.0  0.448  1129.15  1277.15
584-09-8        Rubidium carbonate  230.945   2840.0  0.640  1110.15  1280.15
7429-90-5                 Aluminum   26.982   2377.0  0.311   933.47  1190.15
7429-91-6               Dysprosium  162.500   8370.0  1.430  1685.15  1813.15
7439-88-5                  Iridium  192.217  19000.0  0.000  2719.15  2739.15
...                            ...      ...      ...    ...      ...      ...
13572-98-0  Gadolinium(III) iodide  537.960   4120.0  0.908  1203.15  1305.15
13709-38-1      Lanthanum fluoride  195.900   4589.0  0.682  1766.15  2450.15
13709-59-6    Thorium(IV) fluoride  308.032   6058.0  0.759  1383.15  1651.15
13718-50-8           Barium iodide  391.136   4260.0  0.977   984.15  1248.15
13813-22-4        Lanthanum iodide  519.619   4290.0  1.110  1051.15  1180.15

[177 rows x 6 columns]

In [7]: chemicals.volume.rho_data_CRC_inorg_l_const
Out[7]: 
                              Chemical        Vm
CAS                                             
74-90-8               Hydrogen cyanide  0.000039
75-15-0               Carbon disulfide  0.000060
96-10-6          Chlorodiethylaluminum  0.000126
109-63-7    Boron trifluoride etherate  0.000126
289-22-5              Cyclopentasilane  0.000156
...                                ...       ...
19624-22-7              Pentaborane(9)  0.000105
20398-06-5      Thallium(I) ethanolate  0.000071
23777-80-2              Hexaborane(10)  0.000112
27218-16-2        Chlorine perchlorate  0.000075
52988-75-7          3-Silylpentasilane  0.000217

[116 rows x 2 columns]

In [8]: chemicals.volume.rho_data_CRC_inorg_s_const
Out[8]: 
                                        Chemical        Vm
CAS                                                       
62-54-4                          Calcium acetate  0.000105
62-76-0                           Sodium oxalate  0.000057
75-20-7                          Calcium carbide  0.000029
127-08-2                       Potassium acetate  0.000063
127-09-3                          Sodium acetate  0.000054
...                                          ...       ...
75926-28-2              Selenium sulfide [Se4S4]  0.000135
84359-31-9   Chromium(III) phosphate hexahydrate  0.000120
92141-86-1                     Cesium metaborate  0.000047
133578-89-9             Vanadyl selenite hydrate  0.000060
133863-98-6         Molybdenum(VI) metaphosphate  0.000174

[1872 rows x 2 columns]

In [9]: chemicals.volume.rho_data_CRC_virial
Out[9]: 
                         Chemical      a1       a2       a3       a4   a5
CAS                                                                      
56-23-5        Tetrachloromethane -1600.0  -4059.0  -4653.0      0.0  0.0
60-29-7             Diethyl ether -1226.0  -4458.0  -7746.0 -10005.0  0.0
64-17-5                   Ethanol -4475.0 -29719.0 -56716.0      0.0  0.0
67-56-1                  Methanol -1752.0  -4694.0      0.0      0.0  0.0
67-63-0                2-Propanol -3165.0 -16092.0 -24197.0      0.0  0.0
...                           ...     ...      ...      ...      ...  ...
7783-81-5    Uranium(VI) fluoride -1204.0  -2690.0  -2144.0      0.0  0.0
7783-82-6   Tungsten(VI) fluoride  -719.0  -1143.0      0.0      0.0  0.0
7803-51-2               Phosphine  -146.0   -733.0   1022.0  -1220.0  0.0
10024-97-2          Nitrous oxide  -130.0   -307.0   -248.0      0.0  0.0
10102-43-9           Nitric oxide   -12.0   -119.0     89.0    -73.0  0.0

[105 rows x 6 columns]