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.

ξ=2173.424Tc1/6MWPc2/3\xi = \frac{2173.424 T_c^{1/6}}{\sqrt{MW} P_c^{2/3}}
ξ(0)=(1.51742.135Tr+0.75Tr2)105\xi^{(0)} = (1.5174 - 2.135T_r + 0.75T_r^2)\cdot 10^{-5}
ξ(1)=(4.25527.674Tr+3.4Tr2)105\xi^{(1)} = (4.2552 - 7.674 T_r + 3.4 T_r^2)\cdot 10^{-5}
μ=(ξ(0)+ωξ(1))/ξ\mu = (\xi^{(0)} + \omega \xi^{(1)})/\xi
Parameters
Tfloat

Temperature of fluid [K]

MWfloat

Molecular weight of fluid [g/mol]

Tcfloat

Critical temperature of the fluid [K]

Pcfloat

Critical pressure of the fluid [Pa]

omegafloat

Acentric factor of compound

Returns
mu_lfloat

Viscosity of liquid, [Pa*s]

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].

μ=VoE(VVo)\mu=\frac{V_o}{E(V-V_o)}
E=1.12+Vc12.94+0.10MW0.23Pc+0.0424Tm11.58(Tm/Tc)E=-1.12+\frac{V_c}{12.94+0.10MW-0.23P_c+0.0424T_{m}-11.58(T_{m}/T_c)}
Vo=0.0085ωTc2.02+Vm0.342(Tm/Tc)+0.894V_o = 0.0085\omega T_c-2.02+\frac{V_{m}}{0.342(T_m/T_c)+0.894}
Parameters
Tfloat

Temperature of the fluid [K]

Tmfloat

Melting point of fluid [K]

Tcfloat

Critical temperature of the fluid [K]

Pcfloat

Critical pressure of the fluid [Pa]

Vcfloat

Critical volume of the fluid [m^3/mol]

Vmfloat

Molar volume of the fluid at temperature [m^3]

omegafloat

Acentric factor of compound

MWfloat

Molecular weight of fluid [g/mol]

Returns
mu_lfloat

Viscosity of liquid, [Pa*s]

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.

μμsat=1+D(ΔPr/2.118)A1+CωΔPr\frac{\mu}{\mu_{sat}}=\frac{1+D(\Delta P_r/2.118)^A}{1+C\omega \Delta P_r}
ΔPr=PPsatPc\Delta P_r = \frac{P-P^{sat}}{P_c}
A=0.99914.674×1041.0523Tr0.038771.0513A=0.9991-\frac{4.674\times 10^{-4}}{1.0523T_r^{-0.03877}-1.0513}
D=0.3257(1.0039Tr2.573)0.29060.2086D = \frac{0.3257}{(1.0039-T_r^{2.573})^{0.2906}}-0.2086
C=0.07921+2.1616Tr13.4040Tr2+44.1706Tr384.8291Tr4+96.1209Tr559.8127Tr6+15.6719Tr7C = -0.07921+2.1616T_r-13.4040T_r^2+44.1706T_r^3-84.8291T_r^4+ 96.1209T_r^5-59.8127T_r^6+15.6719T_r^7
Parameters
Tfloat

Temperature of fluid [K]

Pfloat

Pressure of fluid [Pa]

Tc: float

Critical point of fluid [K]

Pcfloat

Critical pressure of the fluid [Pa]

omegafloat

Acentric factor of compound

Psatfloat

Saturation pressure of the fluid [Pa]

mu_lfloat

Viscosity of liquid at 1 atm or saturation, [Pa*s]

Returns
mu_l_densefloat

Viscosity of liquid, [Pa*s]

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].

ηξ×108=46.10Tr0.61820.40exp(0.449Tr)+19.40exp(4.058Tr)+1\eta \xi \times 10^8 = 46.10 T_r^{0.618} - 20.40 \exp(-0.449T_r) + 1 9.40\exp(-4.058T_r)+1
ξ=2173.424Tc1/6MW1/2Pc2/3\xi = 2173.424 T_c^{1/6} MW^{-1/2} P_c^{-2/3}
Parameters
Tfloat

Temperature of the fluid [K]

Tcfloat

Critical temperature of the fluid [K]

Pcfloat

Critical pressure of the fluid [Pa]

MWfloat

Molecular weight of fluid [g/mol]

Returns
mu_gfloat

Viscosity of gas, [Pa*s]

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 Tr>1.5T_r > 1.5:

μg=17.78×105(4.58Tr1.67)0.625/ξ\mu_g = 17.78\times 10^{-5} (4.58T_r - 1.67)^{0.625}/\xi

else:

μg=34×105Tr0.94/ξ\mu_g = 34\times 10^{-5} T_r^{0.94}/\xi
ξ=Tc(1/6)MWPc2/3\xi = \frac{T_c^{(1/6)}}{\sqrt{MW} P_c^{2/3}}
Parameters
Tfloat

Temperature of the fluid [K]

Tcfloat

Critical temperature of the fluid [K]

Pcfloat

Critical pressure of the fluid [Pa]

MWfloat

Molecular weight of fluid [g/mol]

Returns
mu_gfloat

Viscosity of gas, [Pa*s]

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].

η=[0.807Tr0.6180.357exp(0.449Tr)+0.340exp(4.058Tr)+0.018]FpFQ/ξ\eta = \left[0.807T_r^{0.618}-0.357\exp(-0.449T_r) + 0.340\exp(-4.058 T_r) + 0.018\right]F_p^\circ F_Q^\circ /\xi
Fp=1,0μr<0.022F_p^\circ=1, 0 \le \mu_{r} < 0.022
Fp=1+30.55(0.292Zc)1.72,0.022μr<0.075F_p^\circ = 1+30.55(0.292-Z_c)^{1.72}, 0.022 \le \mu_{r} < 0.075
Fp=1+30.55(0.292Zc)1.720.96+0.1(Tr0.7)0.075<μrF_p^\circ = 1+30.55(0.292-Z_c)^{1.72}|0.96+0.1(T_r-0.7)| 0.075 < \mu_{r}
FQ=1.22Q0.15{1+0.00385[(Tr12)2]1/Msign(Tr12)}F_Q^\circ = 1.22Q^{0.15}\left\{ 1+0.00385[(T_r-12)^2]^{1/M}\text{sign} (T_r-12)\right\}
μr=52.46μ2PcTc2\mu_r = 52.46 \frac{\mu^2 P_c}{T_c^2}
ξ=0.176(TcMW3Pc4)1/6\xi=0.176\left(\frac{T_c}{MW^3 P_c^4}\right)^{1/6}
Parameters
Tfloat

Temperature of fluid [K]

Tc: float

Critical point of fluid [K]

Pcfloat

Critical pressure of the fluid [Pa]

Zcfloat

Critical compressibility of the fluid [Pa]

dipolefloat

Dipole moment of fluid [debye]

CASRNstr, optional

CAS of the fluid

Returns
mu_gfloat

Viscosity of gas, [Pa*s]

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].

μ=107105PcTr+(0.0910.477M)T+M(105Pc8M2T2)(10.7639Tc4.1929T)\mu = 10^{-7} | 10^{-5} P_cT_r + \left(0.091-\frac{0.477}{M}\right)T + M \left(10^{-5}P_c-\frac{8M^2}{T^2}\right) \left(\frac{10.7639}{T_c}-\frac{4.1929}{T}\right)|
Parameters
Tfloat

Temperature of the fluid [K]

Tcfloat

Critical temperature of the fluid [K]

Pcfloat

Critical pressure of the fluid [Pa]

MWfloat

Molecular weight of fluid [g/mol]

Returns
mu_gfloat

Viscosity of gas, [Pa*s]

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].

μ=xiμiMWixiMWi\mu = \frac{\sum x_i \mu_i \sqrt{MW_i}} {\sum x_i \sqrt{MW_i}}
Parameters
zsfloat

Mole fractions of components, [-]

musfloat

Gas viscosities of all components, [Pa*s]

MWsfloat

Molecular weights of all components, [g/mol]

MW_rootsfloat, optional

Square roots of molecular weights of all components, [g^0.5/mol^0.5]

Returns
mugfloat

Viscosity of gas mixture, [Pa*s]

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].

ηmix=i=1nyiηij=1nyjϕij\eta_{mix} = \sum_{i=1}^n \frac{y_i \eta_i}{\sum_{j=1}^n y_j \phi_{ij}}
ϕij=(ηiηj)0.5SijAij\phi_{ij} = \left( \frac{\eta_i}{\eta_j} \right)^{0.5} S_{ij} A_{ij}
Aij=mijMij0.5[1+MijMij0.452(1+Mij)+(1+Mij0.45)mij0.51+mij]A_{ij} = m_{ij} M_{ij}^{-0.5} \left[1 + \frac{M_{ij} - M_{ij}^{0.45}} {2(1+M_{ij}) + \frac{(1 + M_{ij}^{0.45}) m_{ij}^{-0.5}}{1 + m_{ij}}} \right]
mij=[4(1+Mij1)(1+Mij)]0.25m_{ij} = \left[ \frac{4}{(1+M_{ij}^{-1})(1+M_{ij})}\right]^{0.25}
Mij=MiMjM_{ij} = \frac{M_i}{M_j}
Sij=1+(TiTj)0.5+(δiδj/4)[1+Ti+(δi2/4)]0.5[1+Tj+(δj2/4)]0.5S_{ij} = \frac{1 + (T_i^* T_j^*)^{0.5} + (\delta_i \delta_j/4)} {[1+T_i^* + (\delta_i^2/4)]^{0.5}[1+T_j^*+(\delta_j^2/4)]^{0.5}}
T=kT/ϵT^* = kT/\epsilon
Parameters
Tfloat

Temperature of fluid, [K]

ysfloat

Mole fractions of gas components, [-]

musfloat

Gas viscosities of all components, [Pa*s]

MWsfloat

Molecular weights of all components, [g/mol]

molecular_diametersfloat

L-J molecular diameter of all components, [angstroms]

Stockmayersfloat

L-J Stockmayer energy parameters of all components, []

Returns
mugfloat

Viscosity of gas mixture, [Pa*s]

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].

ηmix=i=1nyiηij=1nyjϕij\eta_{mix} = \sum_{i=1}^n \frac{y_i \eta_i}{\sum_{j=1}^n y_j \phi_{ij}}
ϕij=(1+ηi/ηj(MWj/MWi)0.25)28(1+MWi/MWj)\phi_{ij} = \frac{(1 + \sqrt{\eta_i/\eta_j}(MW_j/MW_i)^{0.25})^2} {\sqrt{8(1+MW_i/MW_j)}}
Parameters
ysfloat

Mole fractions of gas components, [-]

musfloat

Gas viscosities of all components, [Pa*s]

MWsfloat

Molecular weights of all components, [g/mol]

Returns
mugfloat

Viscosity of gas mixture, [Pa*s]

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 the Wilke_prefactored function.

t0i,j=MWjMWi8MWiMWj+8t0_{i,j} = \frac{ \sqrt{\frac{MW_{j}}{MW_{i}}}}{ \sqrt{\frac{8 MW_{i}}{MW_{j}} + 8}}
t1i,j=2MWjMWi48MWiMWj+8t1_{i,j} = \frac{2 \sqrt[4]{\frac{MW_{j}}{MW_{i}}} }{\sqrt{\frac{8 MW_{i}}{MW_{j}} + 8}}
t2i,j=18MWiMWj+8t2_{i,j} = \frac{1}{\sqrt{\frac{8 MW_{i}}{MW_{j}} + 8}}
Parameters
MWslist[float]

Molecular weights of all components, [g/mol]

Returns
t0slist[list[float]]

First terms, [-]

t1slist[list[float]]

Second terms, [-]

t2slist[list[float]]

Third terms, [-]

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 ϕij\phi_{ij} term to the following:

ϕij=μiμjt0i,j+μiμjt1i,j+t2i,j\phi_{ij} = \frac{\mu_i}{\mu_j}t0_{i,j} + \sqrt{\frac{\mu_i}{\mu_j}}t1_{i,j} + t2_{i,j}
>>> 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.

ηmix=i=1nyiηij=1nyjϕij\eta_{mix} = \sum_{i=1}^n \frac{y_i \eta_i}{\sum_{j=1}^n y_j \phi_{ij}}
ϕij=μiμjt0i,j+μiμjt1i,j+t2i,j\phi_{ij} = \frac{\mu_i}{\mu_j}t0_{i,j} + \sqrt{\frac{\mu_i}{\mu_j}} t1_{i,j} + t2_{i,j}
Parameters
ysfloat

Mole fractions of gas components, [-]

musfloat

Gas viscosities of all components, [Pa*s]

t0slist[list[float]]

First terms, [-]

t1slist[list[float]]

Second terms, [-]

t2slist[list[float]]

Third terms, [-]

Returns
mugfloat

Viscosity of gas mixture, [Pa*s]

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
ysfloat

Mole fractions of gas components, [-]

musfloat

Gas viscosities of all components, [Pa*s]

MWsfloat

Molecular weights of all components, [g/mol]

Returns
mugfloat

Viscosity of gas mixture, [Pa*s]

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.

μ=μ0μ1(T,ρ)μ2(T,ρ)\mu = \mu_0 \cdot \mu_1(T, \rho) \cdot \mu_2(T, \rho)
μ0(T)=100Ti=03HiTi\mu_0(T) = \frac{100\sqrt{T}}{\sum_{i=0}^3 \frac{H_i}{T^i}}
μ1(T,ρ)=exp[ρi=05((1T1)ij=06Hij(ρ1)j)]\mu_1(T, \rho) = \exp\left[\rho \sum_{i=0}^5 \left(\left(\frac{1}{T} - 1 \right)^i \sum_{j=0}^6 H_{ij}(\rho - 1)^j\right)\right]
if ξ<0.3817016416 nm:\text{if }\xi < 0.3817016416 \text{ nm:}
Y=0.2qcξ(qDξ)5(1qcξ+(qcξ)2765504(qDξ)2)Y = 0.2 q_c \xi(q_D \xi)^5 \left(1 - q_c\xi + (q_c\xi)^2 - \frac{765}{504}(q_D\xi)^2\right)
else:\text{else:}
Y=112sin(3ψD)14qcξsin(2ψD)+1(qcξ)2[11.25(qcξ)2]sin(ψD)1(qcξ)3{[11.5(qcξ)2]ψD(qcξ)211.5L(w)}Y = \frac{1}{12}\sin(3\psi_D) - \frac{1}{4q_c \xi}\sin(2\psi_D) + \frac{1}{(q_c\xi)^2}\left[1 - 1.25(q_c\xi)^2\right]\sin(\psi_D) - \frac{1}{(q_c\xi)^3}\left\{\left[1 - 1.5(q_c\xi)^2\right]\psi_D - \left|(q_c\xi)^2 - 1\right|^{1.5}L(w)\right\}
w=qcξ1qcξ+10.5tan(ψD2)w = \left| \frac{q_c \xi -1}{q_c \xi +1}\right|^{0.5} \tan\left( \frac{\psi_D}{2}\right)
L(w)=ln1+w1w if qcξ>1L(w) = \ln\frac{1 + w}{1 - w} \text{ if }q_c \xi > 1
L(w)=2arctanw if qcξ1L(w) = 2\arctan|w| \text{ if }q_c \xi \le 1
ψD=arccos[(1+qD2ξ2)0.5]\psi_D = \arccos\left[\left(1 + q_D^2 \xi^2\right)^{-0.5}\right]
Δχˉ(Tˉ,ρˉ)=ρˉ[ζ(Tˉ,ρˉ)ζ(TˉR,ρˉ)TˉRTˉ]\Delta \bar\chi(\bar T, \bar \rho) = \bar\rho\left[\zeta(\bar T, \bar \rho) - \zeta(\bar T_R, \bar \rho)\frac{\bar T_R}{\bar T}\right]
ξ=ξ0(ΔχˉΓ0)ν/γ\xi = \xi_0 \left(\frac{\Delta \bar\chi}{\Gamma_0}\right)^{\nu/\gamma}
ζ=(ρˉpˉ)Tˉ\zeta = \left(\frac{\partial\bar\rho}{\partial \bar p}\right)_{\bar T}
Parameters
Tfloat

Temperature of water [K]

rhofloat

Density of water [kg/m^3]

drho_dPfloat, 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_Trfloat, 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]

Returns
mufloat

Viscosity, [Pa*s]

Notes

There are three ways to use this formulation.

  1. Compute the Industrial formulation value which does not include the critical enhacement, by leaving drho_dP and drho_dP_Tr None.

  2. Compute the Scientific formulation value by accurately computing and providing drho_dP and drho_dP_Tr, both with IAPWS-95.

  3. 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.

μ=μ0(T)+μr(T,ρ)\mu = \mu^0(T) + \mu^r(T, \rho)
μ0(T)=0.9266958MTσ2Ω(T)\mu^0(T) = \frac{0.9266958\sqrt{MT}}{\sigma^2 \Omega(T^*)}
Ω(T)=exp(i=04bi[ln(T)]i)\Omega(T^*) = \exp\left( \sum_{i=0}^4 b_i [\ln(T^*)]^i \right)
μr=i=1nNiτtiδdiexp(γiδli)\mu^r = \sum_{i=1}^n N_i \tau^{t_i} \delta^{d_i} \exp\left( -\gamma_i \delta^{l_i}\right)
Parameters
Tfloat

Temperature of air [K]

rhofloat

Molar density of air [mol/m^3]

Returns
mufloat

Viscosity of air, [Pa*s]

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 Tc=132.6312T_c = 132.6312 K and ρc=10447.7\rho_c = 10447.7 mol/m^3. Additional parameters used are σ=0.36\sigma = 0.36 nm, M=28.9586M = 28.9586 g/mol and ek=103.3\frac{e}{k} = 103.3 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
Tfloat

Temperature of fluid [K]

Tbfloat

Normal boiling point, [K]

rhofloat

Liquid density liquid as measured at 60 deg F, [kg/m^3]

Returns
mufloat

Liquid viscosity, [Pa*s]

Notes

The formulas are as follows:

Tc=Tb(0.533272+0.191017×103Tb+0.779681×107Tb20.284376×1010Tb3+0.959468×1028/Tb13)1T_{c}^{\circ}=T_{b}\left(0.533272+0.191017 \times 10^{-3} T_{b} +0.779681 \times 10^{-7} T_{b}^{2} -0.284376 \times 10^{-10} T_{b}^{3}+0.959468 \times 10^{28}/T_{b}^{13}\right)^{-1}
α=1Tb/Tc\alpha=1-T_{b} / T_{c}^{\circ}
ln(ν2+1.5)=4.7322727.0975α+49.4491α250.4706α4\ln \left(\nu_2^{\circ}+1.5\right)=4.73227-27.0975 \alpha +49.4491 \alpha^{2}-50.4706 \alpha^{4}
ln(ν1)=0.801621+1.37179ln(ν2)\ln \left(\nu_1^{\circ}\right)=0.801621+1.37179 \ln \left(\nu_2^{\circ}\right)
SG=0.8435930.128624α3.36159α313749.5α12{SG}^{\circ}=0.843593-0.128624 \alpha-3.36159 \alpha^{3}-13749.5 \alpha^{12}
ΔSG=SGSG\Delta {SG} = {SG} - {SG}^\circ
x=1.9987356.7394/Tb|x|=\left|1.99873-56.7394 / \sqrt{T_{b}}\right|
f1=1.33932xΔSG21.1141ΔSG2/Tbf_{1}=1.33932|x| \Delta {SG} - 21.1141 \Delta {SG}^{2} / \sqrt{T_{b}}
f2=xΔSG21.1141ΔSG2/Tbf_{2}=|x| \Delta {SG}-21.1141 \Delta {SG}^{2} / \sqrt{T_{b}}
ln(ν1+450Tb)=ln(ν1+450Tb)(1+2f112f1)2\ln \left(\nu_{1}+\frac{450}{T_{b}}\right)=\ln \left(\nu_{1}^{\circ} +\frac{450}{T_{b}}\right)\left(\frac{1+2 f_{1}}{1-2 f_{1}}\right)^{2}
ln(ν2+450Tb)=ln(ν2+450Tb)(1+2f212f2)2\ln \left(\nu_{2}+\frac{450}{T_{b}}\right)=\ln \left(\nu_{2}^{\circ} +\frac{450}{T_{b}}\right)\left(\frac{1+2 f_{2}}{1-2 f_{2}}\right)^{2}
Z=ν+0.7+exp(1.471.84ν0.51ν2)Z = \nu+0.7+\exp \left(-1.47-1.84 \nu-0.51 \nu^{2}\right)
B=lnlnZ1lnlnZ2lnT1lnT2B=\frac{\ln \ln Z_{1}-\ln \ln Z_{2}}{\ln T_1-\ln T_2}
lnlnZ=lnlnZ1+B(lnTlnT1)\ln \ln Z=\ln \ln Z_{1}+B(\ln T-\ln T_1)
ν=(Z0.7)exp(0.74873.295Z0.7)+0.6119Z0.7)20.3193Z0.7)3)\nu=(Z-0.7)-\exp \left(-0.7487-3.295Z-0.7)+0.6119Z-0.7)^{2}-0.3193Z-0.7)^{3}\right)

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.

μ[centipoise]=μP low, Stiel-hThodos[centipoise]+poly40.0001ξ\mu [\text{centipoise}] = \mu_{\text{P low, Stiel-hThodos}} [\text{centipoise}] + \frac{\text{poly}^4 - 0.0001}{\xi}
poly=(0.1023+0.023364ρr+0.058533ρr20.040758ρr3+0.0093724ρr4)\text{poly} = (0.1023 + 0.023364 \rho_r + 0.058533\rho_r^2 - 0.040758\rho_r^3 + 0.0093724\rho_r^4)
ξ=Tc1/6MW1/2(Pc[atm])2/3\xi = T_c^{1/6} MW^{-1/2} (P_c\text{[atm]})^{-2/3}
Parameters
Tfloat

Temperature of the fluid [K]

Pfloat

Pressure of the fluid [Pa]

Vmfloat

Molar volume of the fluid at the actual conditions, [m^3/mol]

zslist[float]

Mole fractions of chemicals in the fluid, [-]

MWslist[float]

Molecular weights of chemicals in the fluid [g/mol]

Tcsfloat

Critical temperatures of chemicals in the fluid [K]

Pcsfloat

Critical pressures of chemicals in the fluid [Pa]

Vcsfloat

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]

Returns
mufloat

Viscosity of phase at actual conditions , [Pa*s]

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:

μ=xiμiMWixiMWi\mu = \frac{\sum x_i \mu_i \sqrt{MW_i}}{\sum x_i \sqrt{MW_i}}

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.

μ=Eexp[A(CTTD)1/3+B(CTTD)4/3]\mu = E \exp\left[A \left(\frac{C-T}{T-D}\right)^{1/3} + B \left(\frac{C-T}{T-D}\right)^{4/3} \right]
Parameters
Tfloat

Temperature of fluid [K]

Afloat

Coefficient, [-]

Bfloat

Coefficient, [-]

Cfloat

Coefficient, [K]

Dfloat

Coefficient, [K]

Efloat

Coefficient, [Pa*s]

Returns
mufloat

Liquid viscosity, [Pa*s]

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:

μT=E(ACTD+T3(D+T)(CT3(D+T)213(D+T))CTBCTD+T3(CT)(D+T)2+BCTD+T3(CT3(D+T)213(D+T))BCTD+T3D+T)eACTD+T3+BCTD+T3(CT)D+T\frac{\partial \mu}{\partial T} = E \left(\frac{A \sqrt[3]{\frac{C - T} {- D + T}} \left(- D + T\right) \left(- \frac{C - T}{3 \left(- D + T\right)^{2}} - \frac{1}{3 \left(- D + T\right)}\right)}{C - T} - \frac{B \sqrt[3]{\frac{C - T}{- D + T}} \left(C - T\right)}{\left( - D + T\right)^{2}} + B \sqrt[3]{\frac{C - T}{- D + T}} \left(- \frac{ C - T}{3 \left(- D + T\right)^{2}} - \frac{1}{3 \left(- D + T\right)} \right) - \frac{B \sqrt[3]{\frac{C - T}{- D + T}}}{- D + T}\right) e^{A \sqrt[3]{\frac{C - T}{- D + T}} + \frac{B \sqrt[3]{\frac{C - T} {- D + T}} \left(C - T\right)}{- D + T}}

For the low-temperature region:

μT=E(AC+TD+T3(D+T)(C+T3(D+T)2+13(D+T))C+T+BC+TD+T3(CT)(D+T)2+BC+TD+T3D+TBC+TD+T3(CT)(C+T3(D+T)2+13(D+T))C+T)eAC+TD+T3BC+TD+T3(CT)D+T\frac{\partial \mu}{\partial T} = E \left(- \frac{A \sqrt[3]{\frac{ - C + T}{- D + T}} \left(- D + T\right) \left(- \frac{- C + T}{3 \left(- D + T\right)^{2}} + \frac{1}{3 \left(- D + T\right)}\right) }{- C + T} + \frac{B \sqrt[3]{\frac{- C + T}{- D + T}} \left(C - T\right)}{\left(- D + T\right)^{2}} + \frac{B \sqrt[3]{\frac{ - C + T}{- D + T}}}{- D + T} - \frac{B \sqrt[3]{\frac{- C + T}{ - D + T}} \left(C - T\right) \left(- \frac{- C + T}{3 \left(- D + T\right)^{2}} + \frac{1}{3 \left(- D + T\right)}\right)}{- C + T}\right) e^{- A \sqrt[3]{\frac{- C + T}{- D + T}} - \frac{B \sqrt[3]{\frac{- C + T}{- D + T}} \left(C - T\right)}{- D + T}}
Parameters
Tfloat

Temperature of fluid [K]

Afloat

Coefficient, [-]

Bfloat

Coefficient, [-]

Cfloat

Coefficient, [K]

Dfloat

Coefficient, [K]

Efloat

Coefficient, [Pa*s]

Returns
dmu_dTfloat

First temperature derivative of liquid viscosity, [Pa*s]

mufloat

Liquid viscosity, [Pa*s]

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.

μ=a0Tr(1+a1Tra2(Tr1))1/6\mu = \frac{a_0 T_r}{\left( 1 + a_1 T_r^{a_2}(T_r - 1) \right)^{1/6}}
Parameters
Tfloat

Temperature of fluid [K]

Tcfloat

Critical temperature of fluid [K]

a0float

Coefficient, [-]

a1float

Coefficient, [-]

a2float

Coefficient, [-]

Returns
mufloat

Low pressure gas viscosity, [Pa*s]

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.

μ=exp(A+BT)\mu = \exp\left(A + \frac{B}{T}\right)
Parameters
Tfloat

Temperature of fluid [K]

Afloat

Coefficient, [-]

Bfloat

Coefficient, [K]

Returns
mufloat

Liquid viscosity, [Pa*s]

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.

μ=CTD\mu = C T^D
Parameters
Tfloat

Temperature of fluid [K]

Cfloat

Linear coefficient, [Pa*s]

Dfloat

Exponential coefficient, [-]

Returns
mufloat

Liquid viscosity, [Pa*s]

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]
>>> 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.

log10μ=A+B/(T+C)\log_{10} \mu = A + B/(T + C)
Parameters
Tfloat

Temperature of fluid [K]

Afloat

Coefficient, [-]

Bfloat

Coefficient, [K]

Cfloat

Coefficient, [K]

Returns
mufloat

Liquid viscosity, [Pa*s]

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].

log10μ=A+B/T+CT+DT2\log_{10} \mu = A + B/T + CT + DT^2
Parameters
Tfloat

Temperature of fluid [K]

Afloat

Coefficient, [-]

Bfloat

Coefficient, [K]

Cfloat

Coefficient, [1/K]

Dfloat

Coefficient, [1/K^2]

Returns
mufloat

Liquid viscosity, [Pa*s]

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.

μT=10A+BT+T(C+DT)(BT2+C+2DT)log(10)\frac{\partial \mu}{\partial T} = 10^{A + \frac{B}{T} + T \left(C + D T\right)} \left(- \frac{B}{T^{2}} + C + 2 D T\right) \log{\left(10 \right)}
Parameters
Tfloat

Temperature of fluid [K]

Afloat

Coefficient, [-]

Bfloat

Coefficient, [K]

Cfloat

Coefficient, [1/K]

Dfloat

Coefficient, [1/K^2]

Returns
dmu_dTfloat

First temperature derivative of liquid viscosity, [Pa*s/K]

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.

Parameters
Tslist[float]

Temperatures of the experimental data points, [K]

Afloat

Coefficient, [-]

Bfloat

Coefficient, [K]

Cfloat

Coefficient, [1/K]

Dfloat

Coefficient, [1/K^2]

Returns
jaclist[list[float, 4], len(Ts)]

Matrix of derivatives of the equation with respect to the fitting parameters, [various]

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.

μ=exp[A+BT+CT2+DT3]\mu = \exp\left[A + \frac{B}{T} + \frac{C}{T^2} + \frac{D}{T^3}\right]
Parameters
Tfloat

Temperature of fluid [K]

Afloat

Coefficient, [-]

Bfloat

Coefficient, [K]

Cfloat

Coefficient, [K^2]

Dfloat

Coefficient, [K^3]

Returns
mufloat

Liquid viscosity, [Pa*s]

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
valfloat

Viscosity value in the specified scale; [m^2/s] if ‘kinematic viscosity’; [degrees] if Engler or Barbey; [s] for the other scales.

old_scalestr

String representing the scale that val is in originally.

new_scalestr

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.

Returns
resultfloat

Viscosity value in the specified scale; [m^2/s] if ‘kinematic viscosity’; [degrees] if Engler or Barbey; [s] for the other scales

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.

Because the conversion is performed by spline functions, a re-conversion of a value will not yield exactly the original value. However, it is quite close.

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.

SUSeq=4.6324νt+[1.0+0.03264νt][(3930.2+262.7νt+23.97νt2+1.646νt3)×105)]SUS_{eq} = 4.6324\nu_t + \frac{[1.0 + 0.03264\nu_t]} {[(3930.2 + 262.7\nu_t + 23.97\nu_t^2 + 1.646\nu_t^3)\times10^{-5})]}

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:

L,H=interp(nu100)L, H = \text{interp}(nu_100)

else:

L=0.8353ν1002+14.67ν100216L = 0.8353\nu_{100}^2 + 14.67\nu_{100} - 216
H=0.1684ν1002+11.85ν10097H = 0.1684\nu_{100}^2 + 11.85\nu_{100} - 97

if nu_40 > H:

VI=Lnu40LH100VI = \frac{L-nu_{40}}{L-H}\cdot 100

else:

N=ln(H)ln(ν40)ln(ν100)N = \frac{\ln(H) - \ln(\nu_{40})}{\ln (\nu_{100})}
VI=10N10.00715+100VI = \frac{10^N-1}{0.00715} + 100
Parameters
nu_40float

Dynamic viscosity of fluid at 40°C, [m^2/s]

nu_100float

Dynamic viscosity of fluid at 100°C, [m^2/s]

roundingbool, optional

Whether to round the value or not.

Returns
VI: float

Viscosity index [-]

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]. μ=A+BT+CT2+DT3+ET4\mu = A + BT + CT^2 + DT^3 + ET^4.

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]