Joback Group Contribution Method (thermo.group_contribution.joback)

This module contains an implementation of the Joback group-contribution method. This functionality requires the RDKit library to work.

For submitting pull requests, please use the GitHub issue tracker.

Warning

The Joback class method does not contain all the groups for every chemical. There are often multiple ways of fragmenting a chemical. Other times, the fragmentation algorithm will fail. These limitations are present in both the implementation and the method itself. You are welcome to seek to improve this code but no to little help can be offered.

class thermo.group_contribution.joback.Joback(mol, atom_count=None, MW=None, Tb=None)[source]

Bases: object

Class for performing chemical property estimations with the Joback group contribution method as described in [1] and [2]. This is a very common method with low accuracy but wide applicability. This routine can be used with either its own automatic fragmentation routine, or user specified groups. It is applicable to organic compounds only, and has only 41 groups with no interactions between them. Each method’s documentation describes its accuracy. The automatic fragmentation routine is possible only because of the development of SMARTS expressions to match the Joback groups by Dr. Jason Biggs. The list of SMARTS expressions was posted publically on the RDKit mailing list.

Parameters
molrdkitmol or smiles str

Input molecule for the analysis, [-]

atom_countint, optional

The total number of atoms including hydrogen in the molecule; this will be counted by rdkit if it not provided, [-]

MWfloat, optional

Molecular weight of the molecule; this will be calculated by rdkit if not provided, [g/mol]

Tbfloat, optional

An experimentally known boiling temperature for the chemical; this increases the accuracy of the calculated critical point if provided. [K]

Notes

Be sure to check the status of the automatic fragmentation; not all chemicals with the Joback method are applicable.

Approximately 68% of chemcials in the thermo database seem to be able to be estimated with the Joback method.

If a group which was identified is missign a regressed contribution, the estimated property will be None. However, if not all atoms of a molecule are identified as particular groups, property estimation will go ahead with heavily reduced accuracy. Check the status attribute to be sure a molecule was properly fragmented.

References

1

Joback, Kevin G. “A Unified Approach to Physical Property Estimation Using Multivariate Statistical Techniques.” Thesis, Massachusetts Institute of Technology, 1984.

2

Joback, K.G., and R.C. Reid. “Estimation of Pure-Component Properties from Group-Contributions.” Chemical Engineering Communications 57, no. 1-6 (July 1, 1987): 233-43. doi:10.1080/00986448708960487.

Examples

Analysis of Acetone:

>>> J = Joback('CC(=O)C') 
>>> J.Hfus(J.counts) 
5125.0
>>> J.Cpig(350) 
84.69109750000001
>>> J.status 
'OK'

All properties can be obtained in one go with the estimate method:

>>> J.estimate(callables=False) 
{'Tb': 322.11, 'Tm': 173.5, 'Tc': 500.5590049525365, 'Pc': 4802499.604994407, 'Vc': 0.0002095, 'Hf': -217829.99999999997, 'Gf': -154540.00000000003, 'Hfus': 5125.0, 'Hvap': 29018.0, 'mul_coeffs': [839.1099999999998, -14.99], 'Cpig_coeffs': [7.520000000000003, 0.26084, -0.0001207, 1.545999999999998e-08]}

The results for propionic anhydride (if the status is not OK) should not be used.

>>> J = Joback('CCC(=O)OC(=O)CC') 
>>> J.status 
'Matched some atoms repeatedly: [4]'
>>> J.Cpig(300) 
175.85999999999999

None of the routines need to use the automatic routine; they can be used manually too:

>>> Joback.Tb({1: 2, 24: 1})
322.11
Attributes
calculated_Cpig_coeffs
calculated_mul_coeffs

Methods

Cpig(T)

Computes ideal-gas heat capacity at a specified temperature of an organic compound using the Joback method as a function of chemical structure only.

Cpig_coeffs(counts)

Computes the ideal-gas polynomial heat capacity coefficients of an organic compound using the Joback method as a function of chemical structure only.

Gf(counts)

Estimates the ideal-gas Gibbs energy of formation at 298.15 K of an organic compound using the Joback method as a function of chemical structure only.

Hf(counts)

Estimates the ideal-gas enthalpy of formation at 298.15 K of an organic compound using the Joback method as a function of chemical structure only.

Hfus(counts)

Estimates the enthalpy of fusion of an organic compound at its melting point using the Joback method as a function of chemical structure only.

Hvap(counts)

Estimates the enthalpy of vaporization of an organic compound at its normal boiling point using the Joback method as a function of chemical structure only.

Pc(counts, atom_count)

Estimates the critcal pressure of an organic compound using the Joback method as a function of chemical structure only.

Tb(counts)

Estimates the normal boiling temperature of an organic compound using the Joback method as a function of chemical structure only.

Tc(counts[, Tb])

Estimates the critcal temperature of an organic compound using the Joback method as a function of chemical structure only, or optionally improved by using an experimental boiling point.

Tm(counts)

Estimates the melting temperature of an organic compound using the Joback method as a function of chemical structure only.

Vc(counts)

Estimates the critcal volume of an organic compound using the Joback method as a function of chemical structure only.

estimate([callables])

Method to compute all available properties with the Joback method; returns their results as a dict.

mul(T)

Computes liquid viscosity at a specified temperature of an organic compound using the Joback method as a function of chemical structure only.

mul_coeffs(counts)

Computes the liquid phase viscosity Joback coefficients of an organic compound using the Joback method as a function of chemical structure only.

Cpig(T)[source]

Computes ideal-gas heat capacity at a specified temperature of an organic compound using the Joback method as a function of chemical structure only.

Cpig=iai37.93+[ibi+0.210]T+[ici3.91104]T2+[idi+2.06107]T3C_p^{ig} = \sum_i a_i - 37.93 + \left[ \sum_i b_i + 0.210 \right] T + \left[ \sum_i c_i - 3.91 \cdot 10^{-4} \right] T^2 + \left[\sum_i d_i + 2.06 \cdot 10^{-7}\right] T^3
Parameters
Tfloat

Temperature, [K]

Returns
Cpigfloat

Ideal-gas heat capacity, [J/mol/K]

Examples

>>> J = Joback('CC(=O)C') 
>>> J.Cpig(300) 
75.32642000000001
static Cpig_coeffs(counts)[source]

Computes the ideal-gas polynomial heat capacity coefficients of an organic compound using the Joback method as a function of chemical structure only.

Cpig=iai37.93+[ibi+0.210]T+[ici3.91104]T2+[idi+2.06107]T3C_p^{ig} = \sum_i a_i - 37.93 + \left[ \sum_i b_i + 0.210 \right] T + \left[ \sum_i c_i - 3.91 \cdot 10^{-4} \right] T^2 + \left[\sum_i d_i + 2.06 \cdot 10^{-7}\right] T^3

288 compounds were used by Joback in this determination. No overall error was reported.

The ideal gas heat capacity values used in developing the heat capacity polynomials used 9 data points between 298 K and 1000 K.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Returns
coefficientslist[float]

Coefficients which will result in a calculated heat capacity in in units of J/mol/K, [-]

Examples

>>> c = Joback.Cpig_coeffs({1: 2, 24: 1})
>>> c
[7.520000000000003, 0.26084, -0.0001207, 1.545999999999998e-08]
>>> Cp = lambda T : c[0] + c[1]*T + c[2]*T**2 + c[3]*T**3
>>> Cp(300)
75.32642000000001
static Gf(counts)[source]

Estimates the ideal-gas Gibbs energy of formation at 298.15 K of an organic compound using the Joback method as a function of chemical structure only.

Gformation=53.88+Gf,iG_{formation} = 53.88 + \sum {G_{f,i}}

In the above equation, Gibbs energy of formation is calculated in kJ/mol; it is converted to J/mol here.

328 compounds were used by Joback in this determination, with an absolute average error of 2.0 kcal/mol, standard devaition 4.37 kcal/mol, and AARE of 15.7%.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Returns
Gffloat

Estimated ideal-gas Gibbs energy of formation at 298.15 K, [J/mol]

Examples

>>> Joback.Gf({1: 2, 24: 1})
-154540.00000000003
static Hf(counts)[source]

Estimates the ideal-gas enthalpy of formation at 298.15 K of an organic compound using the Joback method as a function of chemical structure only.

Hformation=68.29+iHf,iH_{formation} = 68.29 + \sum_i {H_{f,i}}

In the above equation, enthalpy of formation is calculated in kJ/mol; it is converted to J/mol here.

370 compounds were used by Joback in this determination, with an absolute average error of 2.2 kcal/mol, standard devaition 2.0 kcal/mol, and AARE of 15.2%.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Returns
Hffloat

Estimated ideal-gas enthalpy of formation at 298.15 K, [J/mol]

Examples

>>> Joback.Hf({1: 2, 24: 1})
-217829.99999999997
static Hfus(counts)[source]

Estimates the enthalpy of fusion of an organic compound at its melting point using the Joback method as a function of chemical structure only.

ΔHfus=0.88+iHfus,i\Delta H_{fus} = -0.88 + \sum_i H_{fus,i}

In the above equation, enthalpy of fusion is calculated in kJ/mol; it is converted to J/mol here.

For 155 compounds tested by Joback, the absolute average error was 485.2 cal/mol and standard deviation was 661.4 cal/mol; the average relative error was 38.7%.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Returns
Hfusfloat

Estimated enthalpy of fusion of the compound at its melting point, [J/mol]

Examples

>>> Joback.Hfus({1: 2, 24: 1})
5125.0
static Hvap(counts)[source]

Estimates the enthalpy of vaporization of an organic compound at its normal boiling point using the Joback method as a function of chemical structure only.

ΔHvap=15.30+iHvap,i\Delta H_{vap} = 15.30 + \sum_i H_{vap,i}

In the above equation, enthalpy of fusion is calculated in kJ/mol; it is converted to J/mol here.

For 368 compounds tested by Joback, the absolute average error was 303.5 cal/mol and standard deviation was 429 cal/mol; the average relative error was 3.88%.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Returns
Hvapfloat

Estimated enthalpy of vaporization of the compound at its normal boiling point, [J/mol]

Examples

>>> Joback.Hvap({1: 2, 24: 1})
29018.0
static Pc(counts, atom_count)[source]

Estimates the critcal pressure of an organic compound using the Joback method as a function of chemical structure only. This correlation was developed using the actual number of atoms forming the molecule as well.

Pc=[0.113+0.0032NAiPc,i]2P_c = \left [0.113 + 0.0032N_A - \sum_i {P_{c,i}}\right ]^{-2}

In the above equation, critical pressure is calculated in bar; it is converted to Pa here.

392 compounds were used by Joback in this determination, with an absolute average error of 2.06 bar, standard devaition 3.2 bar, and AARE of 5.2%.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

atom_countint

Total number of atoms (including hydrogens) in the molecule, [-]

Returns
Pcfloat

Estimated critical pressure, [Pa]

Examples

>>> Joback.Pc({1: 2, 24: 1}, 10)
4802499.604994407
static Tb(counts)[source]

Estimates the normal boiling temperature of an organic compound using the Joback method as a function of chemical structure only.

Tb=198.2+iTb,iT_b = 198.2 + \sum_i {T_{b,i}}

For 438 compounds tested by Joback, the absolute average error was 12.91 K and standard deviation was 17.85 K; the average relative error was 3.6%.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Returns
Tbfloat

Estimated normal boiling temperature, [K]

Examples

>>> Joback.Tb({1: 2, 24: 1})
322.11
static Tc(counts, Tb=None)[source]

Estimates the critcal temperature of an organic compound using the Joback method as a function of chemical structure only, or optionally improved by using an experimental boiling point. If the experimental boiling point is not provided it will be estimated with the Joback method as well.

Tc=Tb[0.584+0.965iTc,i(iTc,i)2]1T_c = T_b \left[0.584 + 0.965 \sum_i {T_{c,i}} - \left(\sum_i {T_{c,i}}\right)^2 \right]^{-1}

For 409 compounds tested by Joback, the absolute average error was 4.76 K and standard deviation was 6.94 K; the average relative error was 0.81%.

Appendix BI of Joback’s work lists 409 estimated critical temperatures.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Tbfloat, optional

Experimental normal boiling temperature, [K]

Returns
Tcfloat

Estimated critical temperature, [K]

Examples

>>> Joback.Tc({1: 2, 24: 1}, Tb=322.11)
500.5590049525365
static Tm(counts)[source]

Estimates the melting temperature of an organic compound using the Joback method as a function of chemical structure only.

Tm=122.5+iTm,iT_m = 122.5 + \sum_i {T_{m,i}}

For 388 compounds tested by Joback, the absolute average error was 22.6 K and standard deviation was 24.68 K; the average relative error was 11.2%.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Returns
Tmfloat

Estimated melting temperature, [K]

Examples

>>> Joback.Tm({1: 2, 24: 1})
173.5
static Vc(counts)[source]

Estimates the critcal volume of an organic compound using the Joback method as a function of chemical structure only.

Vc=17.5+iVc,iV_c = 17.5 + \sum_i {V_{c,i}}

In the above equation, critical volume is calculated in cm^3/mol; it is converted to m^3/mol here.

310 compounds were used by Joback in this determination, with an absolute average error of 7.54 cm^3/mol, standard devaition 13.16 cm^3/mol, and AARE of 2.27%.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Returns
Vcfloat

Estimated critical volume, [m^3/mol]

Examples

>>> Joback.Vc({1: 2, 24: 1})
0.0002095
calculated_Cpig_coeffs = None
calculated_mul_coeffs = None
estimate(callables=True)[source]

Method to compute all available properties with the Joback method; returns their results as a dict. For the tempearture dependent values Cpig and mul, both the coefficients and objects to perform calculations are returned.

mul(T)[source]

Computes liquid viscosity at a specified temperature of an organic compound using the Joback method as a function of chemical structure only.

μliq=MWexp(iμa597.82T+iμb11.202)\mu_{liq} = \text{MW} \exp\left( \frac{ \sum_i \mu_a - 597.82}{T} + \sum_i \mu_b - 11.202 \right)
Parameters
Tfloat

Temperature, [K]

Returns
mulfloat

Liquid viscosity, [Pa*s]

Examples

>>> J = Joback('CC(=O)C') 
>>> J.mul(300) 
0.0002940378347162687
static mul_coeffs(counts)[source]

Computes the liquid phase viscosity Joback coefficients of an organic compound using the Joback method as a function of chemical structure only.

μliq=MWexp(iμa597.82T+iμb11.202)\mu_{liq} = \text{MW} \exp\left( \frac{ \sum_i \mu_a - 597.82}{T} + \sum_i \mu_b - 11.202 \right)

288 compounds were used by Joback in this determination. No overall error was reported.

The liquid viscosity data used was specified to be at “several temperatures for each compound” only. A small and unspecified number of compounds were used in this estimation.

Parameters
countsdict

Dictionary of Joback groups present (numerically indexed) and their counts, [-]

Returns
coefficientslist[float]

Coefficients which will result in a liquid viscosity in in units of Pa*s, [-]

Examples

>>> mu_ab = Joback.mul_coeffs({1: 2, 24: 1})
>>> mu_ab
[839.1099999999998, -14.99]
>>> MW = 58.041864812
>>> mul = lambda T : MW*exp(mu_ab[0]/T + mu_ab[1])
>>> mul(300)
0.0002940378347162687
thermo.group_contribution.joback.J_BIGGS_JOBACK_SMARTS = [['Methyl', '-CH3', '[CX4H3]'], ['Secondary acyclic', '-CH2-', '[!R;CX4H2]'], ['Tertiary acyclic', '>CH-', '[!R;CX4H]'], ['Quaternary acyclic', '>C<', '[!R;CX4H0]'], ['Primary alkene', '=CH2', '[CX3H2]'], ['Secondary alkene acyclic', '=CH-', '[!R;CX3H1;!$([CX3H1](=O))]'], ['Tertiary alkene acyclic', '=C<', '[$([!R;CX3H0]);!$([!R;CX3H0]=[#8])]'], ['Cumulative alkene', '=C=', '[$([CX2H0](=*)=*)]'], ['Terminal alkyne', '≡CH', '[$([CX2H1]#[!#7])]'], ['Internal alkyne', '≡C-', '[$([CX2H0]#[!#7])]'], ['Secondary cyclic', '-CH2- (ring)', '[R;CX4H2]'], ['Tertiary cyclic', '>CH- (ring)', '[R;CX4H]'], ['Quaternary cyclic', '>C< (ring)', '[R;CX4H0]'], ['Secondary alkene cyclic', '=CH- (ring)', '[R;CX3H1,cX3H1]'], ['Tertiary alkene cyclic', '=C< (ring)', '[$([R;CX3H0]);!$([R;CX3H0]=[#8])]'], ['Fluoro', '-F', '[F]'], ['Chloro', '-Cl', '[Cl]'], ['Bromo', '-Br', '[Br]'], ['Iodo', '-I', '[I]'], ['Alcohol', '-OH (alcohol)', '[OX2H;!$([OX2H]-[#6]=[O]);!$([OX2H]-a)]'], ['Phenol', '-OH (phenol)', '[$([OX2H]-a)]'], ['Ether acyclic', '-O- (nonring)', '[OX2H0;!R;!$([OX2H0]-[#6]=[#8])]'], ['Ether cyclic', '-O- (ring)', '[#8X2H0;R;!$([#8X2H0]~[#6]=[#8])]'], ['Carbonyl acyclic', '>C=O (nonring)', '[$([CX3H0](=[OX1]));!$([CX3](=[OX1])-[OX2]);!R]=O'], ['Carbonyl cyclic', '>C=O (ring)', '[$([#6X3H0](=[OX1]));!$([#6X3](=[#8X1])~[#8X2]);R]=O'], ['Aldehyde', 'O=CH- (aldehyde)', '[CX3H1](=O)'], ['Carboxylic acid', '-COOH (acid)', '[OX2H]-[C]=O'], ['Ester', '-COO- (ester)', '[#6X3H0;!$([#6X3H0](~O)(~O)(~O))](=[#8X1])[#8X2H0]'], ['Oxygen double bond other', '=O (other than above)', '[OX1H0;!$([OX1H0]~[#6X3]);!$([OX1H0]~[#7X3]~[#8])]'], ['Primary amino', '-NH2', '[NX3H2]'], ['Secondary amino acyclic', '>NH (nonring)', '[NX3H1;!R]'], ['Secondary amino cyclic', '>NH (ring)', '[#7X3H1;R]'], ['Tertiary amino', '>N- (nonring)', '[#7X3H0;!$([#7](~O)~O)]'], ['Imine acyclic', '-N= (nonring)', '[#7X2H0;!R]'], ['Imine cyclic', '-N= (ring)', '[#7X2H0;R]'], ['Aldimine', '=NH', '[#7X2H1]'], ['Cyano', '-CN', '[#6X2]#[#7X1H0]'], ['Nitro', '-NO2', '[$([#7X3,#7X3+][!#8])](=[O])~[O-]'], ['Thiol', '-SH', '[SX2H]'], ['Thioether acyclic', '-S- (nonring)', '[#16X2H0;!R]'], ['Thioether cyclic', '-S- (ring)', '[#16X2H0;R]']]

Metadata for the Joback groups. The first element is the group name; the second is the group symbol; and the third is the SMARTS matching string.

thermo.group_contribution.joback.J_BIGGS_JOBACK_SMARTS_id_dict = {1: '[CX4H3]', 2: '[!R;CX4H2]', 3: '[!R;CX4H]', 4: '[!R;CX4H0]', 5: '[CX3H2]', 6: '[!R;CX3H1;!$([CX3H1](=O))]', 7: '[$([!R;CX3H0]);!$([!R;CX3H0]=[#8])]', 8: '[$([CX2H0](=*)=*)]', 9: '[$([CX2H1]#[!#7])]', 10: '[$([CX2H0]#[!#7])]', 11: '[R;CX4H2]', 12: '[R;CX4H]', 13: '[R;CX4H0]', 14: '[R;CX3H1,cX3H1]', 15: '[$([R;CX3H0]);!$([R;CX3H0]=[#8])]', 16: '[F]', 17: '[Cl]', 18: '[Br]', 19: '[I]', 20: '[OX2H;!$([OX2H]-[#6]=[O]);!$([OX2H]-a)]', 21: '[$([OX2H]-a)]', 22: '[OX2H0;!R;!$([OX2H0]-[#6]=[#8])]', 23: '[#8X2H0;R;!$([#8X2H0]~[#6]=[#8])]', 24: '[$([CX3H0](=[OX1]));!$([CX3](=[OX1])-[OX2]);!R]=O', 25: '[$([#6X3H0](=[OX1]));!$([#6X3](=[#8X1])~[#8X2]);R]=O', 26: '[CX3H1](=O)', 27: '[OX2H]-[C]=O', 28: '[#6X3H0;!$([#6X3H0](~O)(~O)(~O))](=[#8X1])[#8X2H0]', 29: '[OX1H0;!$([OX1H0]~[#6X3]);!$([OX1H0]~[#7X3]~[#8])]', 30: '[NX3H2]', 31: '[NX3H1;!R]', 32: '[#7X3H1;R]', 33: '[#7X3H0;!$([#7](~O)~O)]', 34: '[#7X2H0;!R]', 35: '[#7X2H0;R]', 36: '[#7X2H1]', 37: '[#6X2]#[#7X1H0]', 38: '[$([#7X3,#7X3+][!#8])](=[O])~[O-]', 39: '[SX2H]', 40: '[#16X2H0;!R]', 41: '[#16X2H0;R]'}