Viscosity (chemicals.viscosity)¶
This module contains various viscosity 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.viscosity.Letsou_Stiel(T, MW, Tc, Pc, omega)[source]¶
Calculates the viscosity of a liquid using an emperical model developed in [1]. However. the fitting parameters for tabulated values in the original article are found in ChemSep.
- Parameters
- Returns
- mu_l
float
Viscosity of liquid, [Pa*s]
- mu_l
Notes
The form of this equation is a polynomial fit to tabulated data. The fitting was performed by the DIPPR. This is DIPPR Procedure 8G: Method for the viscosity of pure, nonhydrocarbon liquids at high temperatures internal units are SI standard. [1]’s units were different. DIPPR test value for ethanol is used.
Average error 34%. Range of applicability is 0.76 < Tr < 0.98.
References
- 1(1,2)
Letsou, Athena, and Leonard I. Stiel. “Viscosity of Saturated Nonpolar Liquids at Elevated Pressures.” AIChE Journal 19, no. 2 (1973): 409-11. doi:10.1002/aic.690190241.
Examples
>>> Letsou_Stiel(400., 46.07, 516.25, 6.383E6, 0.6371) 0.0002036150875308
- chemicals.viscosity.Przedziecki_Sridhar(T, Tm, Tc, Pc, Vc, Vm, omega, MW)[source]¶
Calculates the viscosity of a liquid using an emperical formula developed in [1].
- Parameters
- T
float
Temperature of the fluid [K]
- Tm
float
Melting point of fluid [K]
- Tc
float
Critical temperature of the fluid [K]
- Pc
float
Critical pressure of the fluid [Pa]
- Vc
float
Critical volume of the fluid [m^3/mol]
- Vm
float
Molar volume of the fluid at temperature [m^3]
- omega
float
Acentric factor of compound
- MW
float
Molecular weight of fluid [g/mol]
- T
- Returns
- mu_l
float
Viscosity of liquid, [Pa*s]
- mu_l
Notes
A test by Reid (1983) is used, but only mostly correct. This function is not recommended. Internal units are bar and mL/mol.
References
- 1
Przedziecki, J. W., and T. Sridhar. “Prediction of Liquid Viscosities.” AIChE Journal 31, no. 2 (February 1, 1985): 333-35. doi:10.1002/aic.690310225.
Examples
>>> Przedziecki_Sridhar(383., 178., 591.8, 41E5, 316E-6, 95E-6, .263, 92.14) 0.00021981479956033846
Pure High Pressure Liquid Correlations¶
- chemicals.viscosity.Lucas(T, P, Tc, Pc, omega, Psat, mu_l)[source]¶
Adjustes for pressure the viscosity of a liquid using an emperical formula developed in [1], but as discussed in [2] as the original source is in German.
- Parameters
- Returns
- mu_l_dense
float
Viscosity of liquid, [Pa*s]
- mu_l_dense
Notes
This equation is entirely dimensionless; all dimensions cancel. The example is from Reid (1987); all results agree. Above several thousand bar, this equation does not represent true behavior. If Psat is larger than P, the fluid may not be liquid; dPr is set to 0.
References
- 1
Lucas, Klaus. “Ein Einfaches Verfahren Zur Berechnung Der Viskositat von Gasen Und Gasgemischen.” Chemie Ingenieur Technik 46, no. 4 (February 1, 1974): 157-157. doi:10.1002/cite.330460413.
- 2
Reid, Robert C.; Prausnitz, John M.; Poling, Bruce E. Properties of Gases and Liquids. McGraw-Hill Companies, 1987.
Examples
>>> Lucas(300., 500E5, 572.2, 34.7E5, 0.236, 0, 0.00068) # methylcyclohexane 0.0010683738499316494
Liquid Mixing Rules¶
No specific correlations are implemented but
chemicals.utils.mixing_logarithmic
with weight fractions is the
recommended form.
Pure Low Pressure Gas Correlations¶
- chemicals.viscosity.Yoon_Thodos(T, Tc, Pc, MW)[source]¶
Calculates the viscosity of a gas using an emperical formula developed in [1].
- Parameters
- Returns
- mu_g
float
Viscosity of gas, [Pa*s]
- mu_g
Notes
This equation has been tested. The equation uses SI units only internally. The constant 2173.424 is an adjustment factor for units. Average deviation within 3% for most compounds. Greatest accuracy with dipole moments close to 0. Hydrogen and helium have different coefficients, not implemented. This is DIPPR Procedure 8B: Method for the Viscosity of Pure, non hydrocarbon, nonpolar gases at low pressures
References
- 1
Yoon, Poong, and George Thodos. “Viscosity of Nonpolar Gaseous Mixtures at Normal Pressures.” AIChE Journal 16, no. 2 (1970): 300-304. doi:10.1002/aic.690160225.
Examples
>>> Yoon_Thodos(300., 556.35, 4.5596E6, 153.8) 1.019488572777e-05
- chemicals.viscosity.Stiel_Thodos(T, Tc, Pc, MW)[source]¶
Calculates the viscosity of a gas using an emperical formula developed in [1].
if :
else:
- Parameters
- Returns
- mu_g
float
Viscosity of gas, [Pa*s]
- mu_g
Notes
Claimed applicability from 0.2 to 5 atm. Developed with data from 52 nonpolar, and 53 polar gases. internal units are poise and atm. Seems to give reasonable results.
References
- 1
Stiel, Leonard I., and George Thodos. “The Viscosity of Nonpolar Gases at Normal Pressures.” AIChE Journal 7, no. 4 (1961): 611-15. doi:10.1002/aic.690070416.
Examples
>>> Stiel_Thodos(300., 556.35, 4.5596E6, 153.8) #CCl4 1.040892622360e-05
- chemicals.viscosity.Lucas_gas(T, Tc, Pc, Zc, MW, dipole=0.0, CASRN=None)[source]¶
Estimate the viscosity of a gas using an emperical formula developed in several sources, but as discussed in [1] as the original sources are in German or merely personal communications with the authors of [1].
- Parameters
- Returns
- mu_g
float
Viscosity of gas, [Pa*s]
- mu_g
Notes
The example is from [1]; all results agree. Viscosity is calculated in micropoise, and converted to SI internally (1E-7). Q for He = 1.38; Q for H2 = 0.76; Q for D2 = 0.52.
References
- 1(1,2,3)
Reid, Robert C.; Prausnitz, John M.; Poling, Bruce E. Properties of Gases and Liquids. McGraw-Hill Companies, 1987.
Examples
>>> Lucas_gas(T=550., Tc=512.6, Pc=80.9E5, Zc=0.224, MW=32.042, dipole=1.7) 1.7822676912698925e-05
- chemicals.viscosity.viscosity_gas_Gharagheizi(T, Tc, Pc, MW)[source]¶
Calculates the viscosity of a gas using an emperical formula developed in [1].
- Parameters
- Returns
- mu_g
float
Viscosity of gas, [Pa*s]
- mu_g
Notes
Example is first point in supporting information of article, for methane. This is the prefered function for gas viscosity. 7% average relative deviation. Deviation should never be above 30%. Developed with the DIPPR database. It is believed theoretically predicted values are included in the correlation.
Under 0.2Tc, this correlation has been modified to provide values at the limit.
References
- 1
Gharagheizi, Farhad, Ali Eslamimanesh, Mehdi Sattari, Amir H. Mohammadi, and Dominique Richon. “Corresponding States Method for Determination of the Viscosity of Gases at Atmospheric Pressure.” Industrial & Engineering Chemistry Research 51, no. 7 (February 22, 2012): 3179-85. doi:10.1021/ie202591f.
Examples
>>> viscosity_gas_Gharagheizi(120., 190.564, 45.99E5, 16.04246) 5.215761625399613e-06
Pure High Pressure Gas Correlations¶
No correlations are implemented yet.
Gas Mixing Rules¶
- chemicals.viscosity.Herning_Zipperer(zs, mus, MWs, MW_roots=None)[source]¶
Calculates viscosity of a gas mixture according to mixing rules in [1].
- Parameters
- Returns
- mug
float
Viscosity of gas mixture, [Pa*s]
- mug
Notes
This equation is entirely dimensionless; all dimensions cancel. The original source has not been reviewed.
Adding the square roots can speed up the calculation.
References
- 1
Herning, F. and Zipperer, L,: “Calculation of the Viscosity of Technical Gas Mixtures from the Viscosity of Individual Gases, german”, Gas u. Wasserfach (1936) 79, No. 49, 69.
Examples
>>> Herning_Zipperer([0.5, 0.25, 0.25], [1.78e-05, 1.12e-05, 9.35e-06], [28.0134, 16.043, 30.07]) 1.4174908599465168e-05
- chemicals.viscosity.Brokaw(T, ys, mus, MWs, molecular_diameters, Stockmayers)[source]¶
Calculates viscosity of a gas mixture according to mixing rules in [1].
- Parameters
- T
float
Temperature of fluid, [K]
- ys
float
Mole fractions of gas components, [-]
- mus
float
Gas viscosities of all components, [Pa*s]
- MWs
float
Molecular weights of all components, [g/mol]
- molecular_diameters
float
L-J molecular diameter of all components, [angstroms]
- Stockmayers
float
L-J Stockmayer energy parameters of all components, []
- T
- Returns
- mug
float
Viscosity of gas mixture, [Pa*s]
- mug
Notes
This equation is entirely dimensionless; all dimensions cancel. The original source has not been reviewed.
This is DIPPR Procedure 8D: Method for the Viscosity of Nonhydrocarbon Vapor Mixtures at Low Pressure (Polar and Nonpolar)
References
- 1
Brokaw, R. S. “Predicting Transport Properties of Dilute Gases.” Industrial & Engineering Chemistry Process Design and Development 8, no. 2 (April 1, 1969): 240-53. doi:10.1021/i260030a015.
- 2
Brokaw, R. S. Viscosity of Gas Mixtures, NASA-TN-D-4496, 1968.
- 3
Danner, Ronald P, and Design Institute for Physical Property Data. Manual for Predicting Chemical Process Design Data. New York, N.Y, 1982.
Examples
>>> Brokaw(308.2, [0.05, 0.95], [1.34E-5, 9.5029E-6], [64.06, 46.07], [0.42, 0.19], [347, 432]) 9.699085099801568e-06
- chemicals.viscosity.Wilke(ys, mus, MWs)[source]¶
Calculates viscosity of a gas mixture according to mixing rules in [1].
- Parameters
- Returns
- mug
float
Viscosity of gas mixture, [Pa*s]
- mug
See also
Notes
This equation is entirely dimensionless; all dimensions cancel. The original source has not been reviewed or found.
References
- 1
Wilke, C. R. “A Viscosity Equation for Gas Mixtures.” The Journal of Chemical Physics 18, no. 4 (April 1, 1950): 517-19. https://doi.org/10.1063/1.1747673.
Examples
>>> Wilke([0.05, 0.95], [1.34E-5, 9.5029E-6], [64.06, 46.07]) 9.701614885866193e-06
- chemicals.viscosity.Wilke_prefactors(MWs)[source]¶
The
Wilke
gas viscosity method can be sped up by precomputing several matrices. The memory used is proportional to N^2, so it can be significant, but is still a substantial performance increase even when they are so large they cannot fit into cached memory. These matrices are functions of molecular weights only. These are used by theWilke_prefactored
function.- Parameters
- Returns
Notes
These terms are derived as follows using SymPy. The viscosity terms are not known before hand so they are not included in the factors, but otherwise these parameters simplify the computation of the term to the following:
>>> from sympy import * >>> MWi, MWj, mui, muj = symbols('MW_i, MW_j, mu_i, mu_j') >>> f = (1 + sqrt(mui/muj)*(MWj/MWi)**Rational(1,4))**2 >>> denom = sqrt(8*(1+MWi/MWj)) >>> (expand(simplify(expand(f))/denom)) mu_i*sqrt(MW_j/MW_i)/(mu_j*sqrt(8*MW_i/MW_j + 8)) + 2*(MW_j/MW_i)**(1/4)*sqrt(mu_i/mu_j)/sqrt(8*MW_i/MW_j + 8) + 1/sqrt(8*MW_i/MW_j + 8)
Examples
>>> Wilke_prefactors([64.06, 46.07]) ([[0.25, 0.19392193320396522], [0.3179655106303118, 0.25]], [[0.5, 0.421161930934918], [0.5856226024677849, 0.5]], [[0.25, 0.22867110638055677], [0.2696470380083788, 0.25]]) >>> Wilke_prefactored([0.05, 0.95], [1.34E-5, 9.5029E-6], *Wilke_prefactors([64.06, 46.07])) 9.701614885866193e-06
- chemicals.viscosity.Wilke_prefactored(ys, mus, t0s, t1s, t2s)[source]¶
Calculates viscosity of a gas mixture according to mixing rules in [1], using precomputed parameters.
- Parameters
- Returns
- mug
float
Viscosity of gas mixture, [Pa*s]
- mug
See also
Notes
This equation is entirely dimensionless; all dimensions cancel.
References
- 1
Wilke, C. R. “A Viscosity Equation for Gas Mixtures.” The Journal of Chemical Physics 18, no. 4 (April 1, 1950): 517-19. https://doi.org/10.1063/1.1747673.
Examples
>>> Wilke_prefactored([0.05, 0.95], [1.34E-5, 9.5029E-6], *Wilke_prefactors([64.06, 46.07])) 9.701614885866193e-06
- chemicals.viscosity.Wilke_large(ys, mus, MWs)[source]¶
Calculates viscosity of a gas mixture according to mixing rules in [1].
This function is a slightly faster version of
Wilke
. It achieves its extra speed by avoiding some checks, some powers, and by allocating less memory during the computation. For very large component vectors, this function should be called instead.- Parameters
- Returns
- mug
float
Viscosity of gas mixture, [Pa*s]
- mug
See also
References
- 1
Wilke, C. R. “A Viscosity Equation for Gas Mixtures.” The Journal of Chemical Physics 18, no. 4 (April 1, 1950): 517-19. https://doi.org/10.1063/1.1747673.
Examples
>>> Wilke_large([0.05, 0.95], [1.34E-5, 9.5029E-6], [64.06, 46.07]) 9.701614885866193e-06
Correlations for Specific Substances¶
- chemicals.viscosity.mu_IAPWS(T, rho, drho_dP=None, drho_dP_Tr=None)[source]¶
Calculates and returns the viscosity of water according to the IAPWS (2008) release.
Viscosity is calculated as a function of three terms; the first is the dilute-gas limit; the second is the contribution due to finite density; and the third and most complex is a critical enhancement term.
- Parameters
- T
float
Temperature of water [K]
- rho
float
Density of water [kg/m^3]
- drho_dP
float
,optional
Partial derivative of density with respect to pressure at constant temperature (at the temperature and density of water), [kg/m^3/Pa]
- drho_dP_Tr
float
,optional
Partial derivative of density with respect to pressure at constant temperature (at the reference temperature (970.644 K) and the actual density of water), [kg/m^3/Pa]
- T
- Returns
- mu
float
Viscosity, [Pa*s]
- mu
Notes
There are three ways to use this formulation.
Compute the Industrial formulation value which does not include the critical enhacement, by leaving drho_dP and drho_dP_Tr None.
Compute the Scientific formulation value by accurately computing and providing drho_dP and drho_dP_Tr, both with IAPWS-95.
Get a non-standard but 8 decimal place matching result by providing drho_dP computed with either IAPWS-95 or IAPWS-97, but not providing drho_dP_Tr; which is calculated internally. There is a formulation for that term in the thermal conductivity IAPWS equation which is used.
xmu = 0.068
qc = (1.9E-9)**-1
qd = (1.1E-9)**-1
nu = 0.630
gamma = 1.239
xi0 = 0.13E-9
Gamma0 = 0.06
TRC = 1.5
This forulation is highly optimized, spending most of its time in the logarithm, power, and square root.
References
- 1
Huber, M. L., R. A. Perkins, A. Laesecke, D. G. Friend, J. V. Sengers, M. J. Assael, I. N. Metaxa, E. Vogel, R. Mares, and K. Miyagawa. “New International Formulation for the Viscosity of H2O.” Journal of Physical and Chemical Reference Data 38, no. 2 (June 1, 2009): 101-25. doi:10.1063/1.3088050.
Examples
>>> mu_IAPWS(298.15, 998.) 0.000889735100149808
>>> mu_IAPWS(1173.15, 400.) 6.415460784836147e-05
Point 4 of formulation, compared with MPEI and IAPWS, matches.
>>> mu_IAPWS(T=647.35, rho=322., drho_dP=1.213641949033E-2) 4.2961578738287e-05
Full scientific calculation:
>>> from chemicals.iapws import iapws95_properties, iapws95_P, iapws95_Tc >>> T, P = 298.15, 1e5 >>> rho, _, _, _, _, _, _, _, _, _, drho_dP = iapws95_properties(T, P) >>> P_ref = iapws95_P(1.5*iapws95_Tc, rho) >>> _, _, _, _, _, _, _, _, _, _, drho_dP_Tr = iapws95_properties(1.5*iapws95_Tc, P_ref) >>> mu_IAPWS(T, rho, drho_dP, drho_dP_Tr) 0.00089002267377
- chemicals.viscosity.mu_air_lemmon(T, rho)[source]¶
Calculates and returns the viscosity of air according to Lemmon and Jacobsen (2003) [1].
Viscosity is calculated as a function of two terms; the first is the dilute-gas limit; the second is the contribution due to finite density.
- Parameters
- Returns
- mu
float
Viscosity of air, [Pa*s]
- mu
Notes
The coefficients are:
Ni = [10.72, 1.122, 0.002019, -8.876, -0.02916]
ti = [0.2, 0.05, 2.4, 0.6, 3.6]
di = [1, 4, 9, 1, 8]
gammai = Ii = [0, 0, 0, 1, 1]
bi = [.431, -0.4623, 0.08406, 0.005341, -0.00331]
The reducing parameters are K and mol/m^3. Additional parameters used are nm, g/mol and K.
This is an implementation optimized for speed, spending its time in the calclulation of 1 log; 2 exp; 1 power; and 2 divisions.
References
- 1
Lemmon, E. W., and R. T. Jacobsen. “Viscosity and Thermal Conductivity Equations for Nitrogen, Oxygen, Argon, and Air.” International Journal of Thermophysics 25, no. 1 (January 1, 2004): 21-69. https://doi.org/10.1023/B:IJOT.0000022327.04529.f3.
Examples
Viscosity at 300 K and 1 bar:
>>> mu_air_lemmon(300.0, 40.10292351061862) 1.85371518556e-05
Calculate the density in-place:
>>> from chemicals.air import lemmon2000_rho >>> mu_air_lemmon(300.0, lemmon2000_rho(300.0, 1e5)) 1.85371518556e-05
Petroleum Correlations¶
- chemicals.viscosity.Twu_1985(T, Tb, rho)[source]¶
Calculate the viscosity of a petroleum liquid using the Twu (1985) correlation developed in [1]. Based on a fit to n-alkanes that used as a reference. Requires the boiling point and density of the system.
- Parameters
- Returns
- mu
float
Liquid viscosity, [Pa*s]
- mu
Notes
The formulas are as follows:
References
- 1
Twu, Chorng H. “Internally Consistent Correlation for Predicting Liquid Viscosities of Petroleum Fractions.” Industrial & Engineering Chemistry Process Design and Development 24, no. 4 (October 1, 1985): 1287-93. https://doi.org/10.1021/i200031a064.
Examples
Sample point from article:
>>> Twu_1985(T=338.7055, Tb=672.3166, rho=895.5189) 0.008235009644854494
- chemicals.viscosity.Lorentz_Bray_Clarke(T, P, Vm, zs, MWs, Tcs, Pcs, Vcs)[source]¶
Calculates the viscosity of a gas or a liquid using the method of Lorentz, Bray, and Clarke [1]. This method is not quite the same as the original, but rather the form commonly presented and used today. The original had a different formula for pressure correction for gases which was tabular and not presented entirely in [1]. However using that distinction introduces a discontinuity between the liquid and gas viscosity, so it is not normally used.
- Parameters
- T
float
Temperature of the fluid [K]
- P
float
Pressure of the fluid [Pa]
- Vm
float
Molar volume of the fluid at the actual conditions, [m^3/mol]
- zs
list
[float
] Mole fractions of chemicals in the fluid, [-]
- MWs
list
[float
] Molecular weights of chemicals in the fluid [g/mol]
- Tcs
float
Critical temperatures of chemicals in the fluid [K]
- Pcs
float
Critical pressures of chemicals in the fluid [Pa]
- Vcs
float
Critical molar volumes of chemicals in the fluid; these are often used as tuning parameters, fit to match a pure component experimental viscosity value [m^3/mol]
- T
- Returns
- mu
float
Viscosity of phase at actual conditions , [Pa*s]
- mu
Notes
An example from [2] was implemented and checked for validation. Somewhat different rounding is used in [2].
The mixing of the pure component Stiel-Thodos viscosities happens with the Herning-Zipperer mixing rule:
References
- 1(1,2)
Lohrenz, John, Bruce G. Bray, and Charles R. Clark. “Calculating Viscosities of Reservoir Fluids From Their Compositions.” Journal of Petroleum Technology 16, no. 10 (October 1, 1964): 1,171-1,176. https://doi.org/10.2118/915-PA.
- 2(1,2)
Whitson, Curtis H., and Michael R. Brulé. Phase Behavior. Henry L. Doherty Memorial Fund of AIME, Society of Petroleum Engineers, 2000.
Examples
>>> Lorentz_Bray_Clarke(T=300.0, P=1e6, Vm=0.0023025, zs=[.4, .3, .3], ... MWs=[16.04246, 30.06904, 44.09562], Tcs=[190.564, 305.32, 369.83], ... Pcs=[4599000.0, 4872000.0, 4248000.0], Vcs=[9.86e-05, 0.0001455, 0.0002]) 9.925488160761484e-06
Fit Correlations¶
- chemicals.viscosity.PPDS9(T, A, B, C, D, E)[source]¶
Calculate the viscosity of a liquid using the 5-term exponential power fit developed by the PPDS and named PPDS equation 9.
- Parameters
- Returns
- mu
float
Liquid viscosity, [Pa*s]
- mu
Notes
No other source for these coefficients has been found.
There can be a singularity in this equation when T approaches C or D; it may be helpful to take as a limit to this equation D + 5 K.
References
- 1
Gesellschaft, V. D. I., ed. VDI Heat Atlas. 2nd edition. Berlin; New York:: Springer, 2010.
Examples
>>> PPDS9(400.0, 1.74793, 1.33728, 482.347, 41.78, 9.963e-05) 0.00035091137378230684
- chemicals.viscosity.dPPDS9_dT(T, A, B, C, D, E)[source]¶
Calculate the temperature derivative of viscosity of a liquid using the 5-term exponential power fit developed by the PPDS and named PPDS equation 9.
Normally, the temperature derivative is:
For the low-temperature region:
- Parameters
- Returns
References
- 1
Gesellschaft, V. D. I., ed. VDI Heat Atlas. 2nd edition. Berlin; New York:: Springer, 2010.
Examples
>>> dPPDS9_dT(400.0, 1.74793, 1.33728, 482.347, 41.78, 9.963e-05) (-3.186540635882627e-06, 0.00035091137378230684)
- chemicals.viscosity.PPDS5(T, Tc, a0, a1, a2)[source]¶
Calculate the viscosity of a low-pressure gas using the 3-term exponential power fit developed by the PPDS and named PPDS equation 5.
- Parameters
- Returns
- mu
float
Low pressure gas viscosity, [Pa*s]
- mu
References
- 1
“ThermoData Engine (TDE103b V10.1) User`s Guide.” https://trc.nist.gov/TDE/Help/TDE103b/Eqns-Pure-ViscosityG/PPDS5-ViscosityGas.htm.
Examples
Sample coefficients for n-pentane in [1], at 350 K:
>>> PPDS5(T=350.0, Tc=470.008, a0=1.08003e-5, a1=0.19583, a2=0.811897) 8.096643275836e-06
- chemicals.viscosity.Viswanath_Natarajan_2(T, A, B)[source]¶
Calculate the viscosity of a liquid using the 2-term form representation developed in [1]. Requires input coefficients. The A coefficient is assumed to yield coefficients in Pa*s; if it yields values in 1E-3 Pa*s, remove log(100) for A.
- Parameters
- Returns
- mu
float
Liquid viscosity, [Pa*s]
- mu
Notes
No other source for these coefficients than [1] has been found.
References
- 1(1,2)
Viswanath, Dabir S., and G. Natarajan. Databook On The Viscosity Of Liquids. New York: Taylor & Francis, 1989
Examples
DDBST has 0.0004580 as a value at this temperature for 1-Butanol.
>>> Viswanath_Natarajan_2(348.15, -5.9719-log(100), 1007.0) 0.000459836869568295
- chemicals.viscosity.Viswanath_Natarajan_2_exponential(T, C, D)[source]¶
Calculate the viscosity of a liquid using the 2-term exponential form representation developed in [1]. Requires input coefficients. The A coefficient is assumed to yield coefficients in Pa*s, as all coefficients found so far have been.
- Parameters
- Returns
- mu
float
Liquid viscosity, [Pa*s]
- mu
Notes
No other source for these coefficients has been found.
References
- 1(1,2)
Viswanath, Dabir S., and G. Natarajan. Databook On The Viscosity Of Liquids. New York: Taylor & Francis, 1989
Examples
>>> Ts = [283.15, 288.15, 303.15, 349.65] >>> mus = [2.2173, 2.1530, 1.741, 1.0091] # in cP >>> Viswanath_Natarajan_2_exponential(288.15, 4900800, -3.8075) 0.002114798866203873
Calculation of the AARD of the fit (1% is the value stated in [1].:
>>> mu_calc = [Viswanath_Natarajan_2_exponential(T, 4900800, -3.8075) for T in Ts] >>> float(np.mean([abs((mu - mu_i*1000)/mu) for mu, mu_i in zip(mus, mu_calc)])) 0.010467928813061298
- chemicals.viscosity.Viswanath_Natarajan_3(T, A, B, C)[source]¶
Calculate the viscosity of a liquid using the 3-term Antoine form representation developed in [1]. Requires input coefficients. If the coefficients do not yield viscosity in Pa*s, but rather cP, remove log10(1000) from A.
- Parameters
- Returns
- mu
float
Liquid viscosity, [Pa*s]
- mu
Notes
No other source for these coefficients has been found.
References
- 1
Viswanath, Dabir S., and G. Natarajan. Databook On The Viscosity Of Liquids. New York: Taylor & Francis, 1989
Examples
>>> from math import log10 >>> Viswanath_Natarajan_3(298.15, -2.7173-log10(1000), -1071.18, -129.51) 0.0006129806445142113
- chemicals.viscosity.mu_Yaws(T, A, B, C=0.0, D=0.0)[source]¶
Calculate the viscosity of a liquid using the 4-term Yaws polynomial form. Requires input coefficients. If the coefficients do not yield viscosity in Pa*s, but rather cP, remove log10(1000) from A; this is required for the coefficients in [1].
- Parameters
- Returns
- mu
float
Liquid viscosity, [Pa*s]
- mu
References
- 1
Yaws, Carl L. Thermophysical Properties of Chemicals and Hydrocarbons, Second Edition. 2 edition. Amsterdam Boston: Gulf Professional Publishing, 2014.
Examples
>>> from math import log10 >>> mu_Yaws(300.0, -6.4406-log10(1000), 1117.6, 0.0137, -0.000015465) 0.0010066612081
- chemicals.viscosity.dmu_Yaws_dT(T, A, B, C=0.0, D=0.0)[source]¶
Calculate the temperature derivative of the viscosity of a liquid using the 4-term Yaws polynomial form. Requires input coefficients.
- Parameters
- Returns
- dmu_dT
float
First temperature derivative of liquid viscosity, [Pa*s/K]
- dmu_dT
Examples
>>> dmu_Yaws_dT(300.0, -9.4406, 1117.6, 0.0137, -0.000015465) -1.853591586963e-05
- chemicals.viscosity.mu_Yaws_fitting_jacobian(Ts, A, B, C, D)[source]¶
Compute and return the Jacobian of the property predicted by the Yaws viscosity equation with respect to all the coefficients. This is used in fitting parameters for chemicals.
- chemicals.viscosity.mu_TDE(T, A, B, C, D)[source]¶
Calculate the viscosity of a liquid using the 4-term exponential inverse-temperature fit equation used in NIST’s TDE.
- Parameters
- Returns
- mu
float
Liquid viscosity, [Pa*s]
- mu
References
- 1
“ThermoData Engine (TDE103b V10.1) User`s Guide.” https://trc.nist.gov/TDE/Help/TDE103b/Eqns-Pure-ViscositySatL/ViscosityL.htm.
Examples
Coefficients for isooctane at 400 K, as shown in [1].
>>> mu_TDE(400.0, -14.0878, 3500.26, -678132.0, 6.17706e7) 0.0001822175281438
Conversion functions¶
- chemicals.viscosity.viscosity_converter(val, old_scale, new_scale, extrapolate=False)[source]¶
Converts kinematic viscosity values from different scales which have historically been used. Though they may not be in use much, some standards still specify values in these scales.
- Parameters
- val
float
Viscosity value in the specified scale; [m^2/s] if ‘kinematic viscosity’; [degrees] if Engler or Barbey; [s] for the other scales.
- old_scale
str
String representing the scale that val is in originally.
- new_scale
str
String representing the scale that val should be converted to.
- extrapolatebool
If True, a conversion will be performed even if outside the limits of either scale; if False, and either value is outside a limit, an exception will be raised.
- val
- Returns
- result
float
Viscosity value in the specified scale; [m^2/s] if ‘kinematic viscosity’; [degrees] if Engler or Barbey; [s] for the other scales
- result
Notes
The valid scales for this function are any of the following:
[‘a&w b’, ‘a&w crucible’, ‘american can’, ‘astm 0.07’, ‘astm 0.10’, ‘astm 0.15’, ‘astm 0.20’, ‘astm 0.25’, ‘barbey’, ‘caspers tin plate’, ‘continental can’, ‘crown cork and seal’, ‘demmier #1’, ‘demmier #10’, ‘engler’, ‘ford cup #3’, ‘ford cup #4’, ‘kinematic viscosity’, ‘mac michael’, ‘murphy varnish’, ‘parlin cup #10’, ‘parlin cup #15’, ‘parlin cup #20’, ‘parlin cup #25’, ‘parlin cup #30’, ‘parlin cup #7’, ‘pratt lambert a’, ‘pratt lambert b’, ‘pratt lambert c’, ‘pratt lambert d’, ‘pratt lambert e’, ‘pratt lambert f’, ‘pratt lambert g’, ‘pratt lambert h’, ‘pratt lambert i’, ‘redwood admiralty’, ‘redwood standard’, ‘saybolt furol’, ‘saybolt universal’, ‘scott’, ‘stormer 100g load’, ‘westinghouse’, ‘zahn cup #1’, ‘zahn cup #2’, ‘zahn cup #3’, ‘zahn cup #4’, ‘zahn cup #5’]
Some of those scales are converted linearly; the rest use tabulated data and splines.
The conversion uses cubic spline interpolation and is reversible to nearly machine precision.
The method ‘Saybolt universal’ has a special formula implemented for its conversion, from [4]. It is designed for maximum backwards compatibility with prior experimental data. It is solved by newton’s method when kinematic viscosity is desired as an output.
References
- 1
Hydraulic Institute. Hydraulic Institute Engineering Data Book. Cleveland, Ohio: Hydraulic Institute, 1990.
- 2
Gardner/Sward. Paint Testing Manual. Physical and Chemical Examination of Paints, Varnishes, Lacquers, and Colors. 13th Edition. ASTM, 1972.
- 3
Euverard, M. R., The Efflux Type Viscosity Cup. National Paint, Varnish, and Lacquer Association, 1948.
- 4
API Technical Data Book: General Properties & Characterization. American Petroleum Institute, 7E, 2005.
- 5
ASTM. Standard Practice for Conversion of Kinematic Viscosity to Saybolt Universal Viscosity or to Saybolt Furol Viscosity. D 2161 - 93.
Examples
>>> viscosity_converter(8.79, 'engler', 'parlin cup #7') 52.5 >>> viscosity_converter(700, 'Saybolt Universal Seconds', 'kinematic viscosity') 0.00015108914751515542
- chemicals.viscosity.viscosity_index(nu_40, nu_100, rounding=False)[source]¶
Calculates the viscosity index of a liquid. Requires dynamic viscosity of a liquid at 40°C and 100°C. Value may either be returned with or without rounding. Rounding is performed per the standard.
if nu_100 < 70:
else:
if nu_40 > H:
else:
- Parameters
- Returns
- VI:
float
Viscosity index [-]
- VI:
Notes
VI is undefined for nu_100 under 2 mm^2/s. None is returned if this is the case. Internal units are mm^2/s. Higher values of viscosity index suggest a lesser decrease in kinematic viscosity as temperature increases.
Note that viscosity is a pressure-dependent property, and that the viscosity index is defined for a fluid at whatever pressure it is at. The viscosity index is thus also a function of pressure.
References
- 1
ASTM D2270-10(2016) Standard Practice for Calculating Viscosity Index from Kinematic Viscosity at 40 °C and 100 °C, ASTM International, West Conshohocken, PA, 2016, http://dx.doi.org/10.1520/D2270-10R16
Examples
>>> viscosity_index(73.3E-6, 8.86E-6, rounding=True) 92
Fit Coefficients¶
All of these coefficients are lazy-loaded, so they must be accessed as an attribute of this module.
- chemicals.viscosity.mu_data_Dutt_Prasad¶
Coefficient sfor
chemicals.viscosity.Viswanath_Natarajan_3
from [1] for 100 fluids.
- chemicals.viscosity.mu_data_VN3¶
Coefficients for
chemicals.viscosity.Viswanath_Natarajan_3
from [1] with data for 432 fluids.
- chemicals.viscosity.mu_data_VN2¶
Coefficients for
chemicals.viscosity.Viswanath_Natarajan_2
from [1] with data for 135 fluids.
- chemicals.viscosity.mu_data_VN2E¶
Coefficients for
chemicals.viscosity.Viswanath_Natarajan_2_exponential
from [1] with data for 14 fluids.
- chemicals.viscosity.mu_data_Perrys_8E_2_313¶
A collection of 337 coefficient sets for
chemicals.dippr.EQ101
from the DIPPR database published openly in [3].
- chemicals.viscosity.mu_data_Perrys_8E_2_312¶
A collection of 345 coefficient sets for
chemicals.dippr.EQ102
from the DIPPR database published openly in [3].
- chemicals.viscosity.mu_data_VDI_PPDS_7¶
Coefficients for the model equation
PPDS9
, published openly in [2]. Provides no temperature limits, but has been designed for extrapolation. Extrapolated to low temperatures it provides a smooth exponential increase. However, for some chemicals such as glycerol, extrapolated to higher temperatures viscosity is predicted to increase above a certain point.
- chemicals.viscosity.mu_data_VDI_PPDS_8¶
Coefficients for a tempereture polynomial (T in Kelvin) developed by the PPDS, published openly in [2]. .
- 1(1,2,3,4)
Viswanath, Dabir S., and G. Natarajan. Databook On The Viscosity Of Liquids. New York: Taylor & Francis, 1989
- 2(1,2)
Gesellschaft, V. D. I., ed. VDI Heat Atlas. 2nd edition. Berlin; New York:: Springer, 2010.
- 3(1,2)
Green, Don, and Robert Perry. Perry’s Chemical Engineers’ Handbook, Eighth Edition. McGraw-Hill Professional, 2007.
The structure of each dataframe is shown below:
In [1]: import chemicals
In [2]: chemicals.viscosity.mu_data_Dutt_Prasad
Out[2]:
Chemical A B C Tmin Tmax
CAS
56-23-5 Carbon tetrachloride -1.4708 -324.45 71.19 273.0 373.0
60-29-7 Ethyl ether -4.4735 -3623.26 -648.55 273.0 373.0
62-53-3 Aniline -1.1835 -224.31 170.82 268.0 393.0
64-17-5 Ethyl alcohol -2.8857 -1032.53 -55.95 248.0 348.0
64-18-6 Formic acid -1.4150 -297.43 114.74 281.0 373.0
... ... ... ... ... ... ...
629-59-4 Tetra decane -1.4424 -350.81 100.18 283.0 373.0
629-62-9 Penta decane -1.4073 -348.84 105.48 293.0 373.0
629-78-7 Hepta decane -1.7847 -577.32 51.72 303.0 553.0
693-02-7 1 - Hexyne -3.0941 -1404.92 -233.99 293.0 333.0
3744-21-6 2,2 - Dimethyl propane -0.9128 -30.15 202.98 258.0 283.0
[100 rows x 6 columns]
In [3]: chemicals.viscosity.mu_data_VN3
Out[3]:
Name ... Tmax
CAS ...
57-10-3 Palmitic acid ... 370.0
57-50-1 Sucrose ... 330.0
60-12-8 Phenethyl alcohol ... 380.0
60-35-5 Acetamide ... 500.0
62-53-3 Aniline ... 460.0
... ... ... ...
66538-96-3 1,2,3,4 - Tetrahydro -6 - butyl -hexyl naphtha... ... 380.0
87077-20-1 2-Methyl - 7 -heptanol ... 380.0
99332-99-7 Hexyl thiohexanoate ... 370.0
101433-18-5 Ethyl tetra decanol ... 380.0
109309-32-2 2,2-Di - p - toly butane ... 480.0
[432 rows x 7 columns]
In [4]: chemicals.viscosity.mu_data_VN2
Out[4]:
Name Formula ... Tmin Tmax
CAS ...
71-36-3 1-Butanol C4 H10O ... 220.0 390.0
74-87-3 Methyl chloride CH3Cl ... 250.0 310.0
74-88-4 Iodo methane CH3I ... 270.0 320.0
75-08-1 Ethane thiol C2H6S ... 270.0 300.0
75-18-3 Methyl sulfide C2 H6S ... 270.0 310.0
... ... ... ... ... ...
12200-64-5 Sodium hydroxide hydrate NaOH. H2O ... 330.0 360.0
13478-00-7 Nickle - nitrate hexa hydrate Ni(NO3)2. 6H2O ... 330.0 350.0
18358-66-2 3 - n - Propyl - 4 - methyl sydnone C6 H10N2 O2 ... 290.0 320.0
29136-19-4 Nona decyl benzene C25H44 ... 300.0 350.0
31304-44-6 Sodium acetate hydrate CH3COONa. 3H2O ... 330.0 360.0
[135 rows x 6 columns]
In [5]: chemicals.viscosity.mu_data_VN2E
Out[5]:
Substance Formula ... Tmin Tmax
CAS ...
60-29-7 Ether C4H10O ... 270.0 410.0
64-19-7 Acetic acid C2H4O2 ... 270.0 390.0
75-07-0 Acetaldehyde C2H4O2 ... 270.0 300.0
75-25-2 Bromoform CHBr3 ... 280.0 350.0
78-93-3 Methylketone ethyl C4H8O ... 240.0 360.0
109-73-9 Butyl amine C4H11N ... 270.0 360.0
110-58-7 Amyl amine C5H13N ... 270.0 360.0
111-26-2 n-Hexyl amine C6H15N ... 270.0 380.0
764-49-8 Allyl thiocynate C4H5NS ... 290.0 400.0
2307-17-7 Hexyl thio myrisate C20H40OS ... 300.0 370.0
10034-85-2 Hydrogen iodide HI ... 220.0 240.0
10035-10-6 Hydrogen bromide HBr ... 180.0 200.0
28488-34-8 Methylacetate C3H6O2 ... 270.0 420.0
37340-18-4 Perfluoro-1- isopropoxy hexane C9F20O ... 290.0 320.0
[14 rows x 6 columns]
In [6]: chemicals.viscosity.mu_data_Perrys_8E_2_313
Out[6]:
Chemical C1 C2 ... C5 Tmin Tmax
CAS ...
50-00-0 Formaldehyde -11.2400 751.69 ... 0.0 181.15 254.05
55-21-0 Benzamide -12.6320 2668.20 ... 0.0 403.00 563.15
56-23-5 Carbon tetrachloride -8.0738 1121.10 ... 0.0 250.00 455.00
57-55-6 1,2-Propylene glycol -804.5400 30487.00 ... 1.0 213.15 500.80
60-29-7 Diethyl ether 10.1970 -63.80 ... 0.0 200.00 373.15
... ... ... ... ... ... ... ...
10028-15-6 Ozone -10.9400 415.96 ... 0.0 77.55 208.80
10035-10-6 Hydrogen bromide -11.6330 316.38 ... 0.0 185.15 206.45
10102-43-9 Nitric oxide -246.6500 3150.30 ... 1.0 109.50 180.05
13511-13-2 Propenylcyclohexene -11.2080 1079.80 ... 0.0 199.00 508.80
132259-10-0 Air -20.0770 285.15 ... 10.0 59.15 130.00
[337 rows x 8 columns]
In [7]: chemicals.viscosity.mu_data_Perrys_8E_2_312
Out[7]:
Chemical C1 C2 ... C4 Tmin Tmax
CAS ...
50-00-0 Formaldehyde 4.758000e-07 0.64050 ... 0.0 181.15 1000.0
55-21-0 Benzamide 2.508200e-08 0.96663 ... 0.0 403.00 1000.0
56-23-5 Carbon tetrachloride 3.137000e-06 0.37420 ... 0.0 250.33 1000.0
57-55-6 1,2-Propylene glycol 4.543000e-08 0.91730 ... 0.0 213.15 1000.0
60-29-7 Diethyl ether 1.948000e-06 0.41000 ... 0.0 156.85 1000.0
... ... ... ... ... ... ... ...
10028-15-6 Ozone 1.196000e-07 0.84797 ... 0.0 80.15 1000.0
10035-10-6 Hydrogen bromide 9.170000e-08 0.92730 ... 0.0 206.45 800.0
10102-43-9 Nitric oxide 1.467000e-06 0.51230 ... 0.0 110.00 1500.0
13511-13-2 Propenylcyclohexene 5.474900e-07 0.53893 ... 0.0 199.00 1000.0
132259-10-0 Air 1.425000e-06 0.50390 ... 0.0 80.00 2000.0
[345 rows x 7 columns]
In [8]: chemicals.viscosity.mu_data_VDI_PPDS_7
Out[8]:
Chemical Formula A ... C D E
CAS ...
50-00-0 Formaldehyde CH2O 0.69796 ... 549.921 -44.110 0.000036
56-23-5 Carbon tetrachloride CC4l 0.83033 ... 562.119 -73.328 0.000099
56-81-5 Glycerol C3H8O3 -3.91153 ... 582.480 73.885 0.007996
60-29-7 Diethyl ether C4H10O 2.19245 ... 520.594 -370.873 0.000020
62-53-3 Aniline C6H7N 0.85750 ... 462.011 136.981 0.000282
... ... ... ... ... ... ... ...
10097-32-2 Bromine B2r 3.19074 ... 499.481 -209.817 0.000058
10102-43-9 Nitric oxide NO 7.22569 ... 202.500 -106.123 0.000002
10102-44-0 Nitrogen dioxide NO2 6.86768 ... 423.463 -446.706 0.000009
10544-72-6 Dinitrogentetroxide N2O4 -0.03739 ... 615.987 11.286 0.000139
132259-10-0 Air NaN 2.22755 ... 132.897 4.000 0.000016
[271 rows x 7 columns]
In [9]: chemicals.viscosity.mu_data_VDI_PPDS_8
Out[9]:
Chemical A ... D E
CAS ...
50-00-0 Formaldehyde -8.285000e-07 ... 0.000000e+00 0.000000e+00
56-23-5 Carbon tetrachloride -7.132000e-07 ... 0.000000e+00 0.000000e+00
56-81-5 Glycerol -1.460000e-08 ... 0.000000e+00 0.000000e+00
60-29-7 Diethyl ether -8.933000e-07 ... 0.000000e+00 0.000000e+00
62-53-3 Aniline -9.488000e-07 ... 0.000000e+00 0.000000e+00
... ... ... ... ... ...
10097-32-2 Bromine 1.948300e-06 ... 0.000000e+00 0.000000e+00
10102-43-9 Nitric oxide -9.105000e-07 ... 4.240000e-14 -1.020000e-17
10102-44-0 Nitrogen dioxide -2.285050e-05 ... 1.713400e-13 -4.920000e-17
10544-72-6 Dinitrogentetroxide -8.683000e-07 ... 0.000000e+00 0.000000e+00
132259-10-0 Air -1.702000e-07 ... 4.960000e-14 -1.388000e-17
[274 rows x 6 columns]