UNIFAC Gibbs Excess Model (thermo.unifac)¶
This module contains functions and classes related to the UNIFAC and its many variants. The bulk of the code relates to calculating derivativies, or is tables of data.
For reporting bugs, adding feature requests, or submitting pull requests, please use the GitHub issue tracker or contact the author at Caleb.Andrew.Bell@gmail.com.
Main Model (Object-Oriented)¶
- class thermo.unifac.UNIFAC(T, xs, rs, qs, Qs, vs, psi_coeffs=None, psi_abc=None, version=0)[source]¶
Class for representing an a liquid with excess gibbs energy represented by the UNIFAC equation. This model is capable of representing VL and LL behavior, provided the correct interaction parameters are used. [1] and [2] are good references on this model.
- Parameters
- T
float
Temperature, [K]
- xs
list
[float
] Mole fractions, [-]
- rs
list
[float
] r parameters , [-]
- qs
list
[float
] q parameters , [-]
- Qs
list
[float
] Q parameter for each subgroup; subgroups are not required to but are suggested to be sorted from lowest number to highest number, [-]
- vs
list
[list
[float
]] Indexed by [subgroup][count], this variable is the count of each subgroups in each compound, [-]
- psi_abc
tuple
(list
[list
[float
]], 3),optional
psi interaction parameters between each subgroup; indexed [subgroup][subgroup], not symmetrical; first arg is the matrix for a, then b, and then c. Only one of psi_abc or psi_coeffs is required, [-]
- psi_coeffs
list
[list
[tuple
(float
, 3)]],optional
psi interaction parameters between each subgroup; indexed [subgroup][subgroup][letter], not symmetrical. Only one of psi_abc or psi_coeffs is required, [-]
- version
int
,optional
Which version of the model to use [-]
0 - original UNIFAC, OR UNIFAC LLE
1 - Dortmund UNIFAC (adds T dept, 3/4 power)
2 - PSRK (original with T dept function)
3 - VTPR (drops combinatorial term, Dortmund UNIFAC otherwise)
4 - Lyngby/Larsen has different combinatorial, 2/3 power
5 - UNIFAC KT (2 params for psi, Lyngby/Larsen formulation; otherwise same as original)
- T
Notes
In addition to the methods presented here, the methods of its base class
thermo.activity.GibbsExcess
are available as well.References
- 1
Poling, Bruce E., John M. Prausnitz, and John P. O`Connell. The Properties of Gases and Liquids. 5th edition. New York: McGraw-Hill Professional, 2000.
- 2(1,2)
Gmehling, Jürgen, Michael Kleiber, Bärbel Kolbe, and Jürgen Rarey. Chemical Thermodynamics for Process Simulation. John Wiley & Sons, 2019.
Examples
The DDBST has published numerous sample problems using UNIFAC; a simple binary system from example P05.22a in [2] with n-hexane and butanone-2 is shown below:
>>> from thermo.unifac import UFIP, UFSG >>> GE = UNIFAC.from_subgroups(chemgroups=[{1:2, 2:4}, {1:1, 2:1, 18:1}], T=60+273.15, xs=[0.5, 0.5], version=0, interaction_data=UFIP, subgroups=UFSG) >>> GE.gammas() [1.4276025835, 1.3646545010] >>> GE.GE(), GE.dGE_dT(), GE.d2GE_dT2() (923.641197, 0.206721488, -0.00380070204) >>> GE.HE(), GE.SE(), GE.dHE_dT(), GE.dSE_dT() (854.77193363, -0.2067214889, 1.266203886, 0.0038007020460)
The solution given by the DDBST has the same values [1.428, 1.365], and can be found here: http://chemthermo.ddbst.com/Problems_Solutions/Mathcad_Files/05.22a%20VLE%20of%20Hexane-Butanone-2%20Via%20UNIFAC%20-%20Step%20by%20Step.xps
Methods
CpE
()Calculate and return the first temperature derivative of excess enthalpy of a liquid phase using an activity coefficient model.
Fis
()Calculate the terms used in calculating the combinatorial part.
GE
()Calculate the excess Gibbs energy with the UNIFAC model.
HE
()Calculate and return the excess entropy of a liquid phase using an activity coefficient model.
SE
()Calculates the excess entropy of a liquid phase using an activity coefficient model.
Thetas
()Calculate the parameters used in calculating the residual part.
Calculate the parameters for each chemical in the mixture as a pure species, used in calculating the residual part.
Vis
()Calculate the terms used in calculating the combinatorial part.
Calculate the terms used in calculating the combinatorial part.
Xs
()Calculate the parameters used in calculating the residual part.
Xs_pure
()Calculate the parameters for each chemical in the mixture as a pure species, used in calculating the residual part.
as_json
()Method to create a JSON-friendly representation of the Gibbs Excess model which can be stored, and reloaded later.
Calculate the second mole fraction derivative of the terms used in calculating the combinatorial part.
d2GE_dT2
()Calculate the second temperature derivative of excess Gibbs energy with the UNIFAC model.
d2GE_dTdns
()Calculate and return the mole number derivative of the first temperature derivative of excess Gibbs energy of a liquid phase using an activity coefficient model.
Calculate the first composition derivative and temperature derivative of excess Gibbs energy with the UNIFAC model.
Calculate the second composition derivative of excess Gibbs energy with the UNIFAC model.
Calculate the mole fraction derivatives of the parameters.
Calculate the second mole fraction derivative of the terms used in calculating the combinatorial part.
Calculate the second mole fraction derivative of the terms used in calculating the combinatorial part.
Calculate the second temperature derivative of the parameters for the phase; depends on the phases's composition and temperature.
Calculate the temperature and mole fraction derivatives of the parameters for the phase; depends on the phases's composition and temperature.
Calculate the second mole fraction derivatives of the parameters for the phase; depends on the phases's composition and temperature.
Calculate the second temperature derivative of pure component parameters for the phase; depends on the phases's temperature only.
Second temperature derivatives of the combinatorial part of the UNIFAC model.
Second temperature derivative and first mole fraction derivative of the combinatorial part of the UNIFAC model.
Second composition derivative of the combinatorial part of the UNIFAC model.
Calculates the second temperature derivative of the residual part of the UNIFAC model.
Calculates the second temperature derivative of the residual part of the UNIFAC model.
Calculates the first mole fraction derivative of the temperature derivative of the residual part of the UNIFAC model.
Calculates the second mole fraction derivative of the residual part of the UNIFAC model.
d2nGE_dTdns
()Calculate and return the partial mole number derivative of the first temperature derivative of excess Gibbs energy of a liquid phase using an activity coefficient model.
d2nGE_dninjs
()Calculate and return the second partial mole number derivative of excess Gibbs energy of a liquid phase using an activity coefficient model.
Calculate the term second temperature derivative matrix for all groups interacting with all other groups.
Calculate the third mole fraction derivative of the terms used in calculating the combinatorial part.
d3GE_dT3
()Calculate the third temperature derivative of excess Gibbs energy with the UNIFAC model.
Calculate the third mole fraction derivative of the terms used in calculating the combinatorial part.
Calculate the third mole fraction derivative of the terms used in calculating the combinatorial part.
Calculate the third temperature derivative of the parameters for the phase; depends on the phases's composition and temperature.
Calculate the third temperature derivative of pure component parameters for the phase; depends on the phases's temperature only.
Third temperature derivatives of the combinatorial part of the UNIFAC model.
Third composition derivative of the combinatorial part of the UNIFAC model.
Calculates the third temperature derivative of the residual part of the UNIFAC model.
Calculates the third temperature derivative of the residual part of the UNIFAC model.
Calculate the term third temperature derivative matrix for all groups interacting with all other groups.
dFis_dxs
()Calculate the mole fraction derivative of the terms used in calculating the combinatorial part.
dGE_dT
()Calculate the first temperature derivative of excess Gibbs energy with the UNIFAC model.
dGE_dns
()Calculate and return the mole number derivative of excess Gibbs energy of a liquid phase using an activity coefficient model.
dGE_dxs
()Calculate the first composition derivative of excess Gibbs energy with the UNIFAC model.
dHE_dT
()Calculate and return the first temperature derivative of excess enthalpy of a liquid phase using an activity coefficient model.
dHE_dns
()Calculate and return the mole number derivative of excess enthalpy of a liquid phase using an activity coefficient model.
dHE_dxs
()Calculate and return the mole fraction derivative of excess enthalpy of a liquid phase using an activity coefficient model.
dSE_dT
()Calculate and return the first temperature derivative of excess entropy of a liquid phase using an activity coefficient model.
dSE_dns
()Calculate and return the mole number derivative of excess entropy of a liquid phase using an activity coefficient model.
dSE_dxs
()Calculate and return the mole fraction derivative of excess entropy of a liquid phase using an activity coefficient model.
Calculate the mole fraction derivatives of the parameters.
dVis_dxs
()Calculate the mole fraction derivative of the terms used in calculating the combinatorial part.
Calculate the mole fraction derivative of the terms used in calculating the combinatorial part.
Calculates the first temperature derivative of activity coefficients with the UNIFAC model.
Calculate and return the mole number derivative of activity coefficients of a liquid phase using an activity coefficient model.
Calculates the first mole fraction derivative of activity coefficients with the UNIFAC model.
Calculate the first temperature derivative of the parameters for the phase; depends on the phases's composition and temperature.
Calculate the mole fraction derivatives of the parameters for the phase; depends on the phases's composition and temperature.
Calculate the first temperature derivative of pure component parameters for the phase; depends on the phases's temperature only.
Temperature derivatives of the combinatorial part of the UNIFAC model.
First composition derivative of the combinatorial part of the UNIFAC model.
Calculates the first temperature derivative of the residual part of the UNIFAC model.
Calculates the first temperature derivative of the residual part of the UNIFAC model.
Calculates the first mole fraction derivative of the residual part of the UNIFAC model.
dnGE_dns
()Calculate and return the partial mole number derivative of excess Gibbs energy of a liquid phase using an activity coefficient model.
dnHE_dns
()Calculate and return the partial mole number derivative of excess enthalpy of a liquid phase using an activity coefficient model.
dnSE_dns
()Calculate and return the partial mole number derivative of excess entropy of a liquid phase using an activity coefficient model.
dpsis_dT
()Calculate the term first temperature derivative matrix for all groups interacting with all other groups.
from_json
(json_repr)Method to create a Gibbs Excess model from a JSON-friendly serialization of another Gibbs Excess model.
from_subgroups
(T, xs, chemgroups[, ...])Method to construct a UNIFAC object from a dictionary of interaction parameters parameters and a list of dictionaries of UNIFAC keys.
gammas
()Calculates the activity coefficients with the UNIFAC model.
gammas_args
([T])Return a tuple of arguments at the specified tempearture that can be used to efficiently compute gammas at the specified temperature but with varying compositions.
gammas_infinite_dilution
()Calculate and return the infinite dilution activity coefficients of each component.
Calculate the parameters for the phase; depends on the phases's composition and temperature.
Calculate the pure component parameters for the phase; depends on the phases's temperature only.
Calculates the combinatorial part of the UNIFAC model.
Calculates the residual part of the UNIFAC model.
model_hash
()Basic method to calculate a hash of the non-state parts of the model This is useful for comparing to models to determine if they are the same, i.e. in a VLL flash it is important to know if both liquids have the same model.
psis
()Calculate the term matrix for all groups interacting with all other groups.
state_hash
()Basic method to calculate a hash of the state of the model and its model parameters.
to_T_xs
(T, xs)Method to construct a new
UNIFAC
instance at temperature T, and mole fractions xs with the same parameters as the existing object.gammas_from_args
- Fis()[source]¶
Calculate the terms used in calculating the combinatorial part. A function of mole fractions and the parameters q only.
This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.
- GE()[source]¶
Calculate the excess Gibbs energy with the UNIFAC model.
For the VTPR model, the combinatorial component is set to zero.
- Returns
- GE
float
Excess Gibbs energy, [J/mol]
- GE
- Thetas()[source]¶
Calculate the parameters used in calculating the residual part. A function of mole fractions and group counts only.
- Thetas_pure()[source]¶
Calculate the parameters for each chemical in the mixture as a pure species, used in calculating the residual part. A function of group counts only.
- Vis()[source]¶
Calculate the terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.
This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.
- Vis_modified()[source]¶
Calculate the terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.
This is used in the UNIFAC Dortmund and UNIFAC-NIST model with n=0.75, and the Lyngby model with n=2/3.
- Xs()[source]¶
Calculate the parameters used in calculating the residual part. A function of mole fractions and group counts only.
- Xs_pure()[source]¶
Calculate the parameters for each chemical in the mixture as a pure species, used in calculating the residual part. A function of group counts only, not even mole fractions or temperature.
- d2Fis_dxixjs()[source]¶
Calculate the second mole fraction derivative of the terms used in calculating the combinatorial part. A function of mole fractions and the parameters q only.
This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.
- d2GE_dT2()[source]¶
Calculate the second temperature derivative of excess Gibbs energy with the UNIFAC model.
- Returns
- d2GE_dT2
float
Second temperature derivative of excess Gibbs energy, [J/mol/K^2]
- d2GE_dT2
- d2GE_dTdxs()[source]¶
Calculate the first composition derivative and temperature derivative of excess Gibbs energy with the UNIFAC model.
- d2GE_dxixjs()[source]¶
Calculate the second composition derivative of excess Gibbs energy with the UNIFAC model.
- d2Thetas_dxixjs()[source]¶
Calculate the mole fraction derivatives of the parameters. A function of mole fractions and group counts only.
- d2Vis_dxixjs()[source]¶
Calculate the second mole fraction derivative of the terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.
This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.
- d2Vis_modified_dxixjs()[source]¶
Calculate the second mole fraction derivative of the terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.
This is used in the UNIFAC Dortmund and UNIFAC-NIST model with n=0.75, and the Lyngby model with n=2/3.
- d2lnGammas_subgroups_dT2()[source]¶
Calculate the second temperature derivative of the parameters for the phase; depends on the phases’s composition and temperature.
- d2lnGammas_subgroups_dTdxs()[source]¶
Calculate the temperature and mole fraction derivatives of the parameters for the phase; depends on the phases’s composition and temperature.
The following groups are used as follows to simplfy the number of evaluations:
In the below expression, k` refers to a group, and i refers to a component.
- d2lnGammas_subgroups_dxixjs()[source]¶
Calculate the second mole fraction derivatives of the parameters for the phase; depends on the phases’s composition and temperature.
The following groups are used as follows to simplfy the number of evaluations:
- d2lnGammas_subgroups_pure_dT2()[source]¶
Calculate the second temperature derivative of pure component parameters for the phase; depends on the phases’s temperature only.
In this model, the values come from the
UNIFAC.Thetas_pure
method, where each compound is assumed to be pure.
- d2lngammas_c_dT2()[source]¶
Second temperature derivatives of the combinatorial part of the UNIFAC model. Zero in all variations.
- d2lngammas_c_dTdx()[source]¶
Second temperature derivative and first mole fraction derivative of the combinatorial part of the UNIFAC model. Zero in all variations.
- d2lngammas_c_dxixjs()[source]¶
Second composition derivative of the combinatorial part of the UNIFAC model. For the modified UNIFAC model, the equation is as follows; for the original UNIFAC and UNIFAC LLE, replace with .
For the Lyngby model, the following equations are used:
- d2lngammas_dT2()¶
Calculates the second temperature derivative of the residual part of the UNIFAC model.
where the second Gamma is the pure-component Gamma of group k in component i.
- d2lngammas_r_dT2()[source]¶
Calculates the second temperature derivative of the residual part of the UNIFAC model.
where the second Gamma is the pure-component Gamma of group k in component i.
- d2lngammas_r_dTdxs()[source]¶
Calculates the first mole fraction derivative of the temperature derivative of the residual part of the UNIFAC model.
- d2lngammas_r_dxixjs()[source]¶
Calculates the second mole fraction derivative of the residual part of the UNIFAC model.
- d2psis_dT2()[source]¶
Calculate the term second temperature derivative matrix for all groups interacting with all other groups.
The main model calculates the derivative as a function of three coefficients;
Only the first, a coefficient, is used in the original UNIFAC model as well as the UNIFAC-LLE model, so the expression simplifies to:
For the Lyngby model, the second temperature derivative is:
with K and the a coefficients are specific to each pair of main groups, and they are asymmetric, so .
- d3Fis_dxixjxks()[source]¶
Calculate the third mole fraction derivative of the terms used in calculating the combinatorial part. A function of mole fractions and the parameters q only.
This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.
- d3GE_dT3()[source]¶
Calculate the third temperature derivative of excess Gibbs energy with the UNIFAC model.
- Returns
- d3GE_dT3
float
Third temperature derivative of excess Gibbs energy, [J/mol/K^3]
- d3GE_dT3
- d3Vis_dxixjxks()[source]¶
Calculate the third mole fraction derivative of the terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.
This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.
- d3Vis_modified_dxixjxks()[source]¶
Calculate the third mole fraction derivative of the terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.
This is used in the UNIFAC Dortmund and UNIFAC-NIST model with n=0.75, and the Lyngby model with n=2/3.
- d3lnGammas_subgroups_dT3()[source]¶
Calculate the third temperature derivative of the parameters for the phase; depends on the phases’s composition and temperature.
- d3lnGammas_subgroups_pure_dT3()[source]¶
Calculate the third temperature derivative of pure component parameters for the phase; depends on the phases’s temperature only.
In this model, the values come from the
UNIFAC.Thetas_pure
method, where each compound is assumed to be pure.
- d3lngammas_c_dT3()[source]¶
Third temperature derivatives of the combinatorial part of the UNIFAC model. Zero in all variations.
- d3lngammas_c_dxixjxks()[source]¶
Third composition derivative of the combinatorial part of the UNIFAC model. For the modified UNIFAC model, the equation is as follows; for the original UNIFAC and UNIFAC LLE, replace with .
For the Lyngby model, the following equations are used:
- d3lngammas_dT3()¶
Calculates the third temperature derivative of the residual part of the UNIFAC model.
where the second Gamma is the pure-component Gamma of group k in component i.
- d3lngammas_r_dT3()[source]¶
Calculates the third temperature derivative of the residual part of the UNIFAC model.
where the second Gamma is the pure-component Gamma of group k in component i.
- d3psis_dT3()[source]¶
Calculate the term third temperature derivative matrix for all groups interacting with all other groups.
The main model calculates the derivative as a function of three coefficients;
Only the first, a coefficient, is used in the original UNIFAC model as well as the UNIFAC-LLE model, so the expression simplifies to:
For the Lyngby model, the third temperature derivative is:
with K and the a coefficients are specific to each pair of main groups, and they are asymmetric, so .
- dFis_dxs()[source]¶
Calculate the mole fraction derivative of the terms used in calculating the combinatorial part. A function of mole fractions and the parameters q only.
This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.
- dGE_dT()[source]¶
Calculate the first temperature derivative of excess Gibbs energy with the UNIFAC model.
- Returns
- dGE_dT
float
First temperature derivative of excess Gibbs energy, [J/mol/K]
- dGE_dT
- dGE_dxs()[source]¶
Calculate the first composition derivative of excess Gibbs energy with the UNIFAC model.
- dThetas_dxs()[source]¶
Calculate the mole fraction derivatives of the parameters. A function of mole fractions and group counts only.
- dVis_dxs()[source]¶
Calculate the mole fraction derivative of the terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.
This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.
- dVis_modified_dxs()[source]¶
Calculate the mole fraction derivative of the terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.
This is used in the UNIFAC Dortmund and UNIFAC-NIST model with n=0.75, and the Lyngby model with n=2/3.
- dgammas_dT()[source]¶
Calculates the first temperature derivative of activity coefficients with the UNIFAC model.
- dgammas_dns()[source]¶
Calculate and return the mole number derivative of activity coefficients of a liquid phase using an activity coefficient model.
- dgammas_dxs()[source]¶
Calculates the first mole fraction derivative of activity coefficients with the UNIFAC model.
For the VTPR variant, the combinatorial part is skipped:
- dlnGammas_subgroups_dT()[source]¶
Calculate the first temperature derivative of the parameters for the phase; depends on the phases’s composition and temperature.
- dlnGammas_subgroups_dxs()[source]¶
Calculate the mole fraction derivatives of the parameters for the phase; depends on the phases’s composition and temperature.
The group W is used internally as follows to simplfy the number of evaluations.
- dlnGammas_subgroups_pure_dT()[source]¶
Calculate the first temperature derivative of pure component parameters for the phase; depends on the phases’s temperature only.
In this model, the values come from the
UNIFAC.Thetas_pure
method, where each compound is assumed to be pure.
- dlngammas_c_dT()[source]¶
Temperature derivatives of the combinatorial part of the UNIFAC model. Zero in all variations.
- dlngammas_c_dxs()[source]¶
First composition derivative of the combinatorial part of the UNIFAC model. For the modified UNIFAC model, the equation is as follows; for the original UNIFAC and UNIFAC LLE, replace with .
For the Lyngby model, the following equations are used:
- dlngammas_dT()¶
Calculates the first temperature derivative of the residual part of the UNIFAC model.
where the second Gamma is the pure-component Gamma of group k in component i.
- dlngammas_r_dT()[source]¶
Calculates the first temperature derivative of the residual part of the UNIFAC model.
where the second Gamma is the pure-component Gamma of group k in component i.
- dlngammas_r_dxs()[source]¶
Calculates the first mole fraction derivative of the residual part of the UNIFAC model.
- dpsis_dT()[source]¶
Calculate the term first temperature derivative matrix for all groups interacting with all other groups.
The main model calculates the derivative as a function of three coefficients;
Only the first, a coefficient, is used in the original UNIFAC model as well as the UNIFAC-LLE model, so the expression simplifies to:
For the Lyngby model, the first temperature derivative is:
with K and the a coefficients are specific to each pair of main groups, and they are asymmetric, so .
- static from_subgroups(T, xs, chemgroups, subgroups=None, interaction_data=None, version=0)[source]¶
Method to construct a UNIFAC object from a dictionary of interaction parameters parameters and a list of dictionaries of UNIFAC keys. As the actual implementation is matrix based not dictionary based, this method can be quite convenient.
- Parameters
- T
float
Temperature, [K]
- xs
list
[float
] Mole fractions, [-]
- chemgroups
list
[dict
] List of dictionaries of subgroup IDs and their counts for all species in the mixture, [-]
- subgroups
dict
[int:UNIFAC_subgroup
],optional
UNIFAC subgroup data; available dictionaries in this module include UFSG (original), DOUFSG (Dortmund), or NISTUFSG. The default depends on the given version, [-]
- interaction_data
dict
[int:dict
[int:tuple
(a_mn
,b_mn
,c_mn
)]],optional
UNIFAC interaction parameter data; available dictionaries in this module include UFIP (original), DOUFIP2006 (Dortmund parameters published in 2006), DOUFIP2016 (Dortmund parameters published in 2016), and NISTUFIP. The default depends on the given version, [-]
- version
int
,optional
Which version of the model to use. Defaults to 0, [-]
0 - original UNIFAC, OR UNIFAC LLE
1 - Dortmund UNIFAC (adds T dept, 3/4 power)
2 - PSRK (original with T dept function)
3 - VTPR (drops combinatorial term, Dortmund UNIFAC otherwise)
4 - Lyngby/Larsen has different combinatorial, 2/3 power
5 - UNIFAC KT (2 params for psi, Lyngby/Larsen formulation; otherwise same as original)
- T
- Returns
- UNIFAC
UNIFAC
Object for performing calculations with the UNIFAC activity coefficient model, [-]
- UNIFAC
Notes
Warning
For version 0, the interaction data and subgroups default to the original UNIFAC model (not LLE).
For version 1, the interaction data defaults to the Dortmund parameters publshed in 2016 (not 2006).
Examples
Mixture of [‘benzene’, ‘cyclohexane’, ‘acetone’, ‘ethanol’] according to the Dortmund UNIFAC model:
>>> from thermo.unifac import DOUFIP2006, DOUFSG >>> T = 373.15 >>> xs = [0.2, 0.3, 0.1, 0.4] >>> chemgroups = [{9: 6}, {78: 6}, {1: 1, 18: 1}, {1: 1, 2: 1, 14: 1}] >>> GE = UNIFAC.from_subgroups(T=T, xs=xs, chemgroups=chemgroups, version=1, interaction_data=DOUFIP2006, subgroups=DOUFSG) >>> GE UNIFAC(T=373.15, xs=[0.2, 0.3, 0.1, 0.4], rs=[2.2578, 4.2816, 2.3373, 2.4951999999999996], qs=[2.5926, 5.181, 2.7308, 2.6616], Qs=[1.0608, 0.7081, 0.4321, 0.8927, 1.67, 0.8635], vs=[[0, 0, 1, 1], [0, 0, 0, 1], [6, 0, 0, 0], [0, 0, 0, 1], [0, 0, 1, 0], [0, 6, 0, 0]], psi_abc=([[0.0, 0.0, 114.2, 2777.0, 433.6, -117.1], [0.0, 0.0, 114.2, 2777.0, 433.6, -117.1], [16.07, 16.07, 0.0, 3972.0, 146.2, 134.6], [1606.0, 1606.0, 3049.0, 0.0, -250.0, 3121.0], [199.0, 199.0, -57.53, 653.3, 0.0, 168.2], [170.9, 170.9, -2.619, 2601.0, 464.5, 0.0]], [[0.0, 0.0, 0.0933, -4.674, 0.1473, 0.5481], [0.0, 0.0, 0.0933, -4.674, 0.1473, 0.5481], [-0.2998, -0.2998, 0.0, -13.16, -1.237, -1.231], [-4.746, -4.746, -12.77, 0.0, 2.857, -13.69], [-0.8709, -0.8709, 1.212, -1.412, 0.0, -0.8197], [-0.8062, -0.8062, 1.094, -1.25, 0.1542, 0.0]], [[0.0, 0.0, 0.0, 0.001551, 0.0, -0.00098], [0.0, 0.0, 0.0, 0.001551, 0.0, -0.00098], [0.0, 0.0, 0.0, 0.01208, 0.004237, 0.001488], [0.0009181, 0.0009181, 0.01435, 0.0, -0.006022, 0.01446], [0.0, 0.0, -0.003715, 0.000954, 0.0, 0.0], [0.001291, 0.001291, -0.001557, -0.006309, 0.0, 0.0]]), version=1)
- gammas()[source]¶
Calculates the activity coefficients with the UNIFAC model.
For the VTPR variant, the combinatorial part is skipped:
- gammas_args(T=None)[source]¶
Return a tuple of arguments at the specified tempearture that can be used to efficiently compute gammas at the specified temperature but with varying compositions. This is useful in the context of a TP flash.
- lnGammas_subgroups()[source]¶
Calculate the parameters for the phase; depends on the phases’s composition and temperature.
- lnGammas_subgroups_pure()[source]¶
Calculate the pure component parameters for the phase; depends on the phases’s temperature only.
In this model, the values come from the
UNIFAC.Thetas_pure
method, where each compound is assumed to be pure.
- lngammas_c()[source]¶
Calculates the combinatorial part of the UNIFAC model. For the modified UNIFAC model, the equation is as follows; for the original UNIFAC and UNIFAC LLE, replace with .
For the Lyngby model:
- lngammas_r()[source]¶
Calculates the residual part of the UNIFAC model.
where the second Gamma is the pure-component Gamma of group k in component i.
- property model_id¶
A unique numerical identifier refering to the thermodynamic model being implemented. For internal use.
- psis()[source]¶
Calculate the term matrix for all groups interacting with all other groups.
The main model calculates it as a function of three coefficients;
Only the first, a coefficient, is used in the original UNIFAC model as well as the UNIFAC-LLE model, so the expression simplifies to:
For the Lyngby model, the temperature dependence is modified slightly, as follows:
with K and the a coefficients are specific to each pair of main groups, and they are asymmetric, so .
- to_T_xs(T, xs)[source]¶
Method to construct a new
UNIFAC
instance at temperature T, and mole fractions xs with the same parameters as the existing object.- Parameters
- Returns
Notes
If the new temperature is the same temperature as the existing temperature, if the psi terms or their derivatives have been calculated, they will be set to the new object as well. If the mole fractions are the same, various subgroup terms are also kept.
Main Model (Functional)¶
- thermo.unifac.UNIFAC_gammas(T, xs, chemgroups, cached=None, subgroup_data=None, interaction_data=None, modified=False)[source]¶
Calculates activity coefficients using the UNIFAC model (optionally modified), given a mixture’s temperature, liquid mole fractions, and optionally the subgroup data and interaction parameter data of your choice. The default is to use the original UNIFAC model, with the latest parameters published by DDBST. The model supports modified forms (Dortmund, NIST) when the modified parameter is True.
- Parameters
- T
float
Temperature of the system, [K]
- xs
list
[float
] Mole fractions of all species in the system in the liquid phase, [-]
- chemgroups
list
[dict
] List of dictionaries of subgroup IDs and their counts for all species in the mixture, [-]
- subgroup_data
dict
[UNIFAC_subgroup
] UNIFAC subgroup data; available dictionaries in this module are UFSG (original), DOUFSG (Dortmund), or NISTUFSG ([4]).
- interaction_data
dict
[dict
[tuple
(a_mn
,b_mn
,c_mn
)]] UNIFAC interaction parameter data; available dictionaries in this module are UFIP (original), DOUFIP2006 (Dortmund parameters as published by 2006), DOUFIP2016 (Dortmund parameters as published by 2016), and NISTUFIP ([4]).
- modifiedbool
True if using the modified form and temperature dependence, otherwise False.
- T
- Returns
Notes
The actual implementation of UNIFAC is formulated slightly different than the formulas above for computational efficiency. DDBST switched to using the more efficient forms in their publication, but the numerical results are identical.
The model is as follows:
Combinatorial component
Residual component
R and Q
The newer forms of UNIFAC (Dortmund, NIST) calculate the combinatorial part slightly differently:
Although this form looks substantially different than the original, it infact reverts to the original form if only is replaced by . This is more clear when looking at the full rearranged form as in [3].
In some publications such as [5], the nomenclature is such that and do not contain the top , making and . [5] is also notable for having supporting information containing very nice sets of analytical derivatives.
UNIFAC LLE uses the original formulation of UNIFAC, and otherwise only different interaction parameters.
References
- 1
Gmehling, Jurgen. Chemical Thermodynamics: For Process Simulation. Weinheim, Germany: Wiley-VCH, 2012.
- 2
Fredenslund, Aage, Russell L. Jones, and John M. Prausnitz. “Group Contribution Estimation of Activity Coefficients in Nonideal Liquid Mixtures.” AIChE Journal 21, no. 6 (November 1, 1975): 1086-99. doi:10.1002/aic.690210607.
- 3
Jakob, Antje, Hans Grensemann, Jürgen Lohmann, and Jürgen Gmehling. “Further Development of Modified UNIFAC (Dortmund): Revision and Extension 5.” Industrial & Engineering Chemistry Research 45, no. 23 (November 1, 2006): 7924-33. doi:10.1021/ie060355c.
- 4(1,2)
Kang, Jeong Won, Vladimir Diky, and Michael Frenkel. “New Modified UNIFAC Parameters Using Critically Evaluated Phase Equilibrium Data.” Fluid Phase Equilibria 388 (February 25, 2015): 128-41. doi:10.1016/j.fluid.2014.12.042.
- 5(1,2)
Jäger, Andreas, Ian H. Bell, and Cornelia Breitkopf. “A Theoretically Based Departure Function for Multi-Fluid Mixture Models.” Fluid Phase Equilibria 469 (August 15, 2018): 56-69. https://doi.org/10.1016/j.fluid.2018.04.015.
Examples
>>> UNIFAC_gammas(T=333.15, xs=[0.5, 0.5], chemgroups=[{1:2, 2:4}, {1:1, 2:1, 18:1}]) [1.427602583562, 1.364654501010]
>>> from thermo.unifac import DOUFIP2006 >>> UNIFAC_gammas(373.15, [0.2, 0.3, 0.2, 0.2], ... [{9:6}, {78:6}, {1:1, 18:1}, {1:1, 2:1, 14:1}], ... subgroup_data=DOUFSG, interaction_data=DOUFIP2006, modified=True) [1.1864311137, 1.44028013391, 1.20447983349, 1.972070609029]
- thermo.unifac.UNIFAC_psi(T, subgroup1, subgroup2, subgroup_data, interaction_data, modified=False)[source]¶
Calculates the interaction parameter psi(m, n) for two UNIFAC subgroups, given the system temperature, the UNIFAC subgroups considered for the variant of UNIFAC used, the interaction parameters for the variant of UNIFAC used, and whether or not the temperature dependence is modified from the original form, as shown below.
Original temperature dependence:
Modified temperature dependence:
- Parameters
- T
float
Temperature of the system, [K]
- subgroup1
int
First UNIFAC subgroup for identifier, [-]
- subgroup2
int
Second UNIFAC subgroup for identifier, [-]
- subgroup_data
dict
[UNIFAC_subgroup
] Normally provided as inputs to UNIFAC.
- interaction_data
dict
[dict
[tuple
(a_mn
,b_mn
,c_mn
)]] Normally provided as inputs to UNIFAC.
- modifiedbool
True if the modified temperature dependence is used by the interaction parameters, otherwise False
- T
- Returns
- psi
float
UNIFAC interaction parameter term, [-]
- psi
Notes
UNIFAC interaction parameters are asymmetric. No warning is raised if an interaction parameter is missing.
References
- 1
Gmehling, Jurgen. Chemical Thermodynamics: For Process Simulation. Weinheim, Germany: Wiley-VCH, 2012.
- 2
Fredenslund, Aage, Russell L. Jones, and John M. Prausnitz. “Group Contribution Estimation of Activity Coefficients in Nonideal Liquid Mixtures.” AIChE Journal 21, no. 6 (November 1, 1975): 1086-99. doi:10.1002/aic.690210607.
Examples
>>> from thermo.unifac import UFSG, UFIP, DOUFSG, DOUFIP2006
>>> UNIFAC_psi(307, 18, 1, UFSG, UFIP) 0.9165248264184787
>>> UNIFAC_psi(373.15, 9, 78, DOUFSG, DOUFIP2006, modified=True) 1.3703140538273264
Misc Functions¶
- thermo.unifac.UNIFAC_RQ(groups, subgroup_data=None)[source]¶
Calculates UNIFAC parameters R and Q for a chemical, given a dictionary of its groups, as shown in [1]. Most UNIFAC methods use the same subgroup values; however, a dictionary of UNIFAC_subgroup instances may be specified as an optional second parameter.
- Parameters
- Returns
Notes
These parameters have some predictive value for other chemical properties.
References
- 1
Gmehling, Jurgen. Chemical Thermodynamics: For Process Simulation. Weinheim, Germany: Wiley-VCH, 2012.
Examples
Hexane
>>> UNIFAC_RQ({1:2, 2:4}) (4.4998000000000005, 3.856)
- thermo.unifac.Van_der_Waals_volume(R)[source]¶
Calculates a species Van der Waals molar volume with the UNIFAC method, given a species’s R parameter.
- Parameters
- R
float
R UNIFAC parameter (normalized Van der Waals Volume) [-]
- R
- Returns
- V_vdw
float
Unnormalized Van der Waals volume, [m^3/mol]
- V_vdw
Notes
The volume was originally given in cm^3/mol, but is converted to SI here.
References
- 1
Wei, James, Morton M. Denn, John H. Seinfeld, Arup Chakraborty, Jackie Ying, Nicholas Peppas, and George Stephanopoulos. Molecular Modeling and Theory in Chemical Engineering. Academic Press, 2001.
Examples
>>> Van_der_Waals_volume(4.4998) 6.826196599999999e-05
- thermo.unifac.Van_der_Waals_area(Q)[source]¶
Calculates a species Van der Waals molar surface area with the UNIFAC method, given a species’s Q parameter.
- Parameters
- Q
float
Q UNIFAC parameter (normalized Van der Waals Area) [-]
- Q
- Returns
- A_vdw
float
Unnormalized Van der Waals surface area, [m^2/mol]
- A_vdw
Notes
The volume was originally given in cm^2/mol, but is converted to SI here.
References
- 1
Wei, James, Morton M. Denn, John H. Seinfeld, Arup Chakraborty, Jackie Ying, Nicholas Peppas, and George Stephanopoulos. Molecular Modeling and Theory in Chemical Engineering. Academic Press, 2001.
Examples
>>> Van_der_Waals_area(3.856) 964000.0
- thermo.unifac.chemgroups_to_matrix(chemgroups)[source]¶
Index by [group index][compound index]
>>> chemgroups_to_matrix([{9: 6}, {2: 6}, {1: 1, 18: 1}, {1: 1, 2: 1, 14: 1}]) [[0, 0, 1, 1], [0, 6, 0, 1], [6, 0, 0, 0], [0, 0, 0, 1], [0, 0, 1, 0]]
- thermo.unifac.UNIFAC_group_assignment_DDBST(CAS, model)[source]¶
Lookup the group assignment of a compound in either the ‘UNIFAC’ the ‘MODIFIED_UNIFAC’, or the ‘PSRK’ model. These values are read from a sqlite database on demand.
- Parameters
- Returns
- asssignments
dict
The group assignments and their counts; note that an empty dictionary indicates the fragmentation is not available, [-]
- asssignments
Examples
>>> UNIFAC_group_assignment_DDBST('50-14-6', 'UNIFAC') {1: 5, 2: 8, 3: 6, 4: 1, 6: 1, 7: 1, 8: 2, 14: 1}
Data for Original UNIFAC¶
- thermo.unifac.UFSG = {1: <CH3>, 2: <CH2>, 3: <CH>, 4: <C>, 5: <CH2=CH>, 6: <CH=CH>, 7: <CH2=C>, 8: <CH=C>, 9: <ACH>, 10: <AC>, 11: <ACCH3>, 12: <ACCH2>, 13: <ACCH>, 14: <OH>, 15: <CH3OH>, 16: <H2O>, 17: <ACOH>, 18: <CH3CO>, 19: <CH2CO>, 20: <CHO>, 21: <CH3COO>, 22: <CH2COO>, 23: <HCOO>, 24: <CH3O>, 25: <CH2O>, 26: <CHO>, 27: <THF>, 28: <CH3NH2>, 29: <CH2NH2>, 30: <CHNH2>, 31: <CH3NH>, 32: <CH2NH>, 33: <CHNH>, 34: <CH3N>, 35: <CH2N>, 36: <ACNH2>, 37: <C5H5N>, 38: <C5H4N>, 39: <C5H3N>, 40: <CH3CN>, 41: <CH2CN>, 42: <COOH>, 43: <HCOOH>, 44: <CH2CL>, 45: <CHCL>, 46: <CCL>, 47: <CH2CL2>, 48: <CHCL2>, 49: <CCL2>, 50: <CHCL3>, 51: <CCL3>, 52: <CCL4>, 53: <ACCL>, 54: <CH3NO2>, 55: <CH2NO2>, 56: <CHNO2>, 57: <ACNO2>, 58: <CS2>, 59: <CH3SH>, 60: <CH2SH>, 61: <FURFURAL>, 62: <DOH>, 63: <I>, 64: <BR>, 65: <CH=-C>, 66: <C=-C>, 67: <DMSO>, 68: <ACRY>, 69: <CL-(C=C)>, 70: <C=C>, 71: <ACF>, 72: <DMF>, 73: <HCON(CH2)2>, 74: <CF3>, 75: <CF2>, 76: <CF>, 77: <COO>, 78: <SIH3>, 79: <SIH2>, 80: <SIH>, 81: <SI>, 82: <SIH2O>, 83: <SIHO>, 84: <SIO>, 85: <NMP>, 86: <CCL3F>, 87: <CCL2F>, 88: <HCCL2F>, 89: <HCCLF>, 90: <CCLF2>, 91: <HCCLF2>, 92: <CCLF3>, 93: <CCL2F2>, 94: <AMH2>, 95: <AMHCH3>, 96: <AMHCH2>, 97: <AM(CH3)2>, 98: <AMCH3CH2>, 99: <AM(CH2)2>, 100: <C2H5O2>, 101: <C2H4O2>, 102: <CH3S>, 103: <CH2S>, 104: <CHS>, 105: <MORPH>, 106: <C4H4S>, 107: <C4H3S>, 108: <C4H2S>, 109: <NCO>, 118: <(CH2)2SU>, 119: <CH2CHSU>, 178: <IMIDAZOL>, 179: <BTI>}¶
- thermo.unifac.UFMG = {1: ('CH2', [1, 2, 3, 4]), 2: ('C=C', [5, 6, 7, 8, 70]), 3: ('ACH', [9, 10]), 4: ('ACCH2', [11, 12, 13]), 5: ('OH', [14]), 6: ('CH3OH', [15]), 7: ('H2O', [16]), 8: ('ACOH', [17]), 9: ('CH2CO', [18, 19]), 10: ('CHO', [20]), 11: ('CCOO', [21, 22]), 12: ('HCOO', [23]), 13: ('CH2O', [24, 25, 26, 27]), 14: ('CNH2', [28, 29, 30]), 15: ('CNH', [31, 32, 33]), 16: ('(C)3N', [34, 35]), 17: ('ACNH2', [36]), 18: ('PYRIDINE', [37, 38, 39]), 19: ('CCN', [40, 41]), 20: ('COOH', [42, 43]), 21: ('CCL', [44, 45, 46]), 22: ('CCL2', [47, 48, 49]), 23: ('CCL3', [50, 51]), 24: ('CCL4', [52]), 25: ('ACCL', [53]), 26: ('CNO2', [54, 55, 56]), 27: ('ACNO2', [57]), 28: ('CS2', [58]), 29: ('CH3SH', [59, 60]), 30: ('FURFURAL', [61]), 31: ('DOH', [62]), 32: ('I', [63]), 33: ('BR', [64]), 34: ('C=-C', [65, 66]), 35: ('DMSO', [67]), 36: ('ACRY', [68]), 37: ('CLCC', [69]), 38: ('ACF', [71]), 39: ('DMF', [72, 73]), 40: ('CF2', [74, 75, 76]), 41: ('COO', [77]), 42: ('SIH2', [78, 79, 80, 81]), 43: ('SIO', [82, 83, 84]), 44: ('NMP', [85]), 45: ('CCLF', [86, 87, 88, 89, 90, 91, 92, 93]), 46: ('CON(AM)', [94, 95, 96, 97, 98, 99]), 47: ('OCCOH', [100, 101]), 48: ('CH2S', [102, 103, 104]), 49: ('MORPH', [105]), 50: ('THIOPHEN', [106, 107, 108]), 51: ('NCO', [109]), 55: ('SULFONES', [118, 119]), 84: ('IMIDAZOL', [178]), 85: ('BTI', [179])}¶
Data for Dortmund UNIFAC¶
- thermo.unifac.DOUFSG = {1: <CH3>, 2: <CH2>, 3: <CH>, 4: <C>, 5: <CH2=CH>, 6: <CH=CH>, 7: <CH2=C>, 8: <CH=C>, 9: <ACH>, 10: <AC>, 11: <ACCH3>, 12: <ACCH2>, 13: <ACCH>, 14: <OH(P)>, 15: <CH3OH>, 16: <H2O>, 17: <ACOH>, 18: <CH3CO>, 19: <CH2CO>, 20: <CHO>, 21: <CH3COO>, 22: <CH2COO>, 23: <HCOO>, 24: <CH3O>, 25: <CH2O>, 26: <CHO>, 27: <THF>, 28: <CH3NH2>, 29: <CH2NH2>, 30: <CHNH2>, 31: <CH3NH>, 32: <CH2NH>, 33: <CHNH>, 34: <CH3N>, 35: <CH2N>, 36: <ACNH2>, 37: <AC2H2N>, 38: <AC2HN>, 39: <AC2N>, 40: <CH3CN>, 41: <CH2CN>, 42: <COOH>, 43: <HCOOH>, 44: <CH2CL>, 45: <CHCL>, 46: <CCL>, 47: <CH2CL2>, 48: <CHCL2>, 49: <CCL2>, 50: <CHCL3>, 51: <CCL3>, 52: <CCL4>, 53: <ACCL>, 54: <CH3NO2>, 55: <CH2NO2>, 56: <CHNO2>, 57: <ACNO2>, 58: <CS2>, 59: <CH3SH>, 60: <CH2SH>, 61: <FURFURAL>, 62: <DOH>, 63: <I>, 64: <BR>, 65: <CH=-C>, 66: <C=-C>, 67: <DMSO>, 68: <ACRY>, 69: <CL-(C=C)>, 70: <C=C>, 71: <ACF>, 72: <DMF>, 73: <HCON(CH2)2>, 74: <CF3>, 75: <CF2>, 76: <CF>, 77: <COO>, 78: <CY-CH2>, 79: <CY-CH>, 80: <CY-C>, 81: <OH(S)>, 82: <OH(T)>, 83: <CY-CH2O>, 84: <TRIOXAN>, 85: <CNH2>, 86: <NMP>, 87: <NEP>, 88: <NIPP>, 89: <NTBP>, 91: <CONH2>, 92: <CONHCH3>, 100: <CONHCH2>, 101: <AM(CH3)2>, 102: <AMCH3CH2>, 103: <AM(CH2)2>, 104: <AC2H2S>, 105: <AC2HS>, 106: <AC2S>, 107: <H2COCH>, 108: <COCH>, 109: <HCOCH>, 110: <(CH2)2SU>, 111: <CH2SUCH>, 112: <(CH3)2CB>, 113: <(CH2)2CB>, 114: <CH2CH3CB>, 119: <H2COCH2>, 122: <CH3S>, 123: <CH2S>, 124: <CHS>, 153: <H2COC>, 178: <C3H2N2+>, 179: <BTI->, 184: <C3H3N2+>, 189: <C4H8N+>, 195: <BF4->, 196: <C5H5N+>, 197: <OTF->, 201: <-S-S->}¶
- thermo.unifac.DOUFMG = {1: ('CH2', [1, 2, 3, 4]), 2: ('C=C', [5, 6, 7, 8, 70]), 3: ('ACH', [9, 10]), 4: ('ACCH2', [11, 12, 13]), 5: ('OH', [14, 81, 82]), 6: ('CH3OH', [15]), 7: ('H2O', [16]), 8: ('ACOH', [17]), 9: ('CH2CO', [18, 19]), 10: ('CHO', [20]), 11: ('CCOO', [21, 22]), 12: ('HCOO', [23]), 13: ('CH2O', [24, 25, 26]), 14: ('CH2NH2', [28, 29, 30, 85]), 15: ('CH2NH', [31, 32, 33]), 16: ('(C)3N', [34, 35]), 17: ('ACNH2', [36]), 18: ('PYRIDINE', [37, 38, 39]), 19: ('CH2CN', [40, 41]), 20: ('COOH', [42]), 21: ('CCL', [44, 45, 46]), 22: ('CCL2', [47, 48, 49]), 23: ('CCL3', [51]), 24: ('CCL4', [52]), 25: ('ACCL', [53]), 26: ('CNO2', [54, 55, 56]), 27: ('ACNO2', [57]), 28: ('CS2', [58]), 29: ('CH3SH', [59, 60]), 30: ('FURFURAL', [61]), 31: ('DOH', [62]), 32: ('I', [63]), 33: ('BR', [64]), 34: ('C=-C', [65, 66]), 35: ('DMSO', [67]), 36: ('ACRY', [68]), 37: ('CLCC', [69]), 38: ('ACF', [71]), 39: ('DMF', [72, 73]), 40: ('CF2', [74, 75, 76]), 41: ('COO', [77]), 42: ('CY-CH2', [78, 79, 80]), 43: ('CY-CH2O', [27, 83, 84]), 44: ('HCOOH', [43]), 45: ('CHCL3', [50]), 46: ('CY-CONC', [86, 87, 88, 89]), 47: ('CONR', [91, 92, 100]), 48: ('CONR2', [101, 102, 103]), 49: ('HCONR', [93, 94]), 52: ('ACS', [104, 105, 106]), 53: ('EPOXIDES', [107, 108, 109, 119, 153]), 55: ('CARBONAT', [112, 113, 114]), 56: ('SULFONE', [110, 111]), 61: ('SULFIDES', [122, 123, 124]), 84: ('IMIDAZOL', [178, 184]), 85: ('BTI', [179]), 87: ('PYRROL', [189]), 89: ('BF4', [195]), 90: ('PYRIDIN', [196]), 91: ('OTF', [197]), 93: ('DISULFIDES', [201])}¶
Data for NIST UNIFAC (2015)¶
- thermo.unifac.NISTUFSG = {1: <CH3>, 2: <CH2>, 3: <CH>, 4: <C>, 5: <CH2=CH>, 6: <CH=CH>, 7: <CH2=C>, 8: <CH=C>, 9: <ACH>, 10: <AC>, 11: <ACCH3>, 12: <ACCH2>, 13: <ACCH>, 14: <OH prim>, 15: <CH3OH>, 16: <H2O>, 17: <ACOH>, 18: <CH3CO>, 19: <CH2CO>, 20: <CHO>, 21: <CH3COO>, 22: <CH2COO>, 23: <HCOO>, 24: <CH3O>, 25: <CH2O>, 26: <CHO>, 27: <CH2-O-CH2>, 28: <CH3NH2>, 29: <CH2NH2>, 30: <CHNH2>, 31: <CH3NH>, 32: <CH2NH>, 33: <CHNH>, 34: <CH3N>, 35: <CH2N>, 36: <ACNH2>, 37: <AC2H2N>, 38: <AC2HN>, 39: <AC2N>, 40: <CH3CN>, 41: <CH2CN>, 42: <COOH>, 43: <HCOOH>, 44: <CH2Cl>, 45: <CHCl>, 46: <CCl>, 47: <CH2Cl2>, 48: <CHCl2>, 49: <CCl2>, 50: <CHCl3>, 51: <CCl3>, 52: <CCl4>, 53: <ACCl>, 54: <CH3NO2>, 55: <CH2NO2>, 56: <CHNO2>, 57: <ACNO2>, 58: <CS2>, 59: <CH3SH>, 60: <CH2SH>, 61: <Furfural>, 62: <CH2(OH)-CH2(OH)>, 63: <I>, 64: <Br>, 65: <CH#C>, 66: <C#C>, 67: <DMSO>, 68: <Acrylonitrile>, 69: <Cl-(C=C)>, 70: <C=C>, 71: <ACF>, 72: <DMF>, 73: <HCON(CH2)2>, 74: <CF3>, 75: <CF2>, 76: <CF>, 77: <COO>, 78: <c-CH2>, 79: <c-CH>, 80: <c-C>, 81: <OH sec>, 82: <OH tert>, 83: <CH2-O-[CH2-O]1/2>, 84: <[O-CH2]1/2-O-[CH2-O]1/2>, 85: <CNH2>, 86: <c-CON-CH3>, 87: <c-CON-CH2>, 88: <c-CON-CH>, 89: <c-CON-C>, 92: <CONHCH3>, 93: <HCONHCH3>, 94: <HCONHCH2>, 100: <CONHCH2>, 101: <CON(CH3)2>, 102: <CON(CH3)CH2>, 103: <CON(CH2)2>, 104: <AC2H2S>, 105: <AC2HS>, 106: <AC2S>, 107: <H2COCH>, 109: <HCOCH>, 110: <CH2SuCH2>, 111: <CH2SuCH >, 112: <(CH3O)2CO>, 113: <(CH2O)2CO>, 114: <(CH3O)COOCH2>, 116: <ACCN>, 117: <CH3NCO>, 118: <CH2NCO>, 119: <CHNCO>, 120: <ACNCO>, 121: <COOCO>, 122: <ACSO2>, 123: <ACCHO>, 124: <ACCOOH>, 125: <c-CO-NH>, 126: <c-CO-O>, 127: <AC-O-CO-CH3 >, 128: <AC-O-CO-CH2>, 129: <AC-O-CO-CH>, 130: <AC-O-CO-C>, 131: <-O-CH2-CH2-OH>, 132: <-O-CH-CH2-OH>, 133: <-O-CH2-CH-OH>, 134: <CH3-S->, 135: <-CH2-S->, 136: <>CH-S->, 137: <->C-S->, 138: <CH3O-(O)>, 139: <CH2O-(O)>, 140: <CHO-(O)>, 141: <CO-(O)>, 142: <ACO-(O)>, 143: <CFH>, 144: <CFCl>, 145: <CFCl2>, 146: <CF2H>, 147: <CF2ClH>, 148: <CF2Cl2>, 149: <CF3H>, 150: <CF3Cl>, 151: <CF4>, 152: <C(O)2>, 153: <ACN(CH3)2>, 154: <ACN(CH3)CH2>, 155: <ACN(CH2)2>, 156: <ACNHCH3>, 157: <ACNHCH2>, 158: <ACNHCH>, 159: <AC2H2O>, 160: <AC2HO>, 161: <AC2O>, 162: <c-CH-NH>, 163: <c-C-NH>, 164: <c-CH-NCH3>, 165: <c-CH-NCH2>, 166: <c-CH-NCH>, 170: <SiH3->, 171: <-SiH2->, 172: <>SiH->, 173: <>Si<>, 174: <-SiH2-O->, 175: <>SiH-O->, 176: <->Si-O->, 177: <C=NOH>, 178: <ACCO>, 179: <C2Cl4>, 180: <c-CHH2>, 186: <CH(O)2>, 187: <ACS>, 188: <c-CH2-NH>, 189: <c-CH2-NCH3>, 190: <c-CH2-NCH2>, 191: <c-CH2-NCH>, 192: <CHSH>, 193: <CSH>, 194: <ACSH>, 195: <ACC>, 196: <AC2H2NH>, 197: <AC2HNH>, 198: <AC2NH>, 199: <(ACO)COOCH2>, 200: <(ACO)CO(OAC)>, 201: <c-CH=CH>, 202: <c-CH=C>, 203: <c-C=C>, 204: <Glycerol>, 205: <-CH(OH)-CH2(OH)>, 206: <-CH(OH)-CH(OH)->, 207: <>C(OH)-CH2(OH)>, 208: <>C(OH)-CH(OH)->, 209: <>C(OH)-C(OH)<>, 301: <CHCO>, 302: <CCO>, 303: <CHCN>, 304: <CCN>, 305: <CNO2>, 306: <ACNH>, 307: <ACN>, 308: <HCHO>, 309: <CH=NOH>}¶
- thermo.unifac.NISTUFMG = {1: ('CH2', [1, 2, 3, 4], 'Alkyl chains'), 2: ('C=C', [5, 6, 7, 8, 9], 'Double bonded alkyl chains'), 3: ('ACH', [15, 16, 17], 'Aromatic carbon'), 4: ('ACCH2', [18, 19, 20, 21], 'Aromatic carbon plus alkyl chain'), 5: ('OH', [34, 204, 205], 'Alcohols'), 6: ('CH3OH', [35], 'Methanol'), 7: ('H2O', [36], 'Water'), 8: ('ACOH', [37], 'Phenolic –OH groups '), 9: ('CH2CO', [42, 43, 44, 45], 'Ketones'), 10: ('CHO', [48], 'Aldehydes'), 11: ('CCOO', [51, 52, 53, 54], 'Esters'), 12: ('HCOO', [55], 'Formates'), 13: ('CH2O', [59, 60, 61, 62, 63], 'Ethers'), 14: ('CNH2', [66, 67, 68, 69], 'Amines with 1-alkyl group'), 15: ('(C)2NH', [71, 72, 73], 'Amines with 2-alkyl groups'), 16: ('(C)3N', [74, 75], 'Amines with 3-alkyl groups'), 17: ('ACNH2', [79, 80, 81], 'Anilines'), 18: ('PYRIDINE', [76, 77, 78], 'Pyridines'), 19: ('CCN', [85, 86, 87, 88], 'Nitriles'), 20: ('COOH', [94, 95], 'Acids'), 21: ('CCl', [99, 100, 101], 'Chlorocarbons'), 22: ('CCl2', [102, 103, 104], 'Dichlorocarbons'), 23: ('CCl3', [105, 106], 'Trichlorocarbons'), 24: ('CCl4', [107], 'Tetrachlorocarbons'), 25: ('ACCl', [109], 'Chloroaromatics'), 26: ('CNO2', [132, 133, 134, 135], 'Nitro alkanes'), 27: ('ACNO2', [136], 'Nitroaromatics'), 28: ('CS2', [146], 'Carbon disulfide'), 29: ('CH3SH', [138, 139, 140, 141], 'Mercaptans'), 30: ('FURFURAL', [50], 'Furfural'), 31: ('DOH', [38], 'Ethylene Glycol'), 32: ('I', [128], 'Iodides'), 33: ('BR', [130], 'Bromides'), 34: ('C≡C', [13, 14], 'Triplebonded alkyl chains'), 35: ('DMSO', [153], 'Dimethylsulfoxide'), 36: ('ACRY', [90], 'Acrylic'), 37: ('ClC=C', [108], 'Chlorine attached to double bonded alkyl chain'), 38: ('ACF', [118], 'Fluoroaromatics'), 39: ('DMF', [161, 162, 163, 164, 165], 'Amides'), 40: ('CF2', [111, 112, 113, 114, 115, 116, 117], 'Fluorines'), 41: ('COO', [58], 'Esters'), 42: ('SiH2', [197, 198, 199, 200], 'Silanes'), 43: ('SiO', [201, 202, 203], 'Siloxanes'), 44: ('NMP', [195], 'N-Methyl-2-pyrrolidone'), 45: ('CClF', [120, 121, 122, 123, 124, 125, 126, 127], 'Chloro-Fluorides'), 46: ('CONCH2', [166, 167, 168, 169], 'Amides'), 47: ('OCCOH', [39, 40, 41], 'Oxygenated Alcohols'), 48: ('CH2S', [142, 143, 144, 145], 'Sulfides'), 49: ('MORPHOLIN', [196], 'Morpholine'), 50: ('THIOPHENE', [147, 148, 149], 'Thiophene'), 51: ('CH2(cy)', [27, 28, 29], 'Cyclic hydrocarbon chains'), 52: ('C=C(cy)', [30, 31, 32], 'Cyclic unsaturated hydrocarbon chains')}¶
Data for NIST KT UNIFAC (2011)¶
- thermo.unifac.NISTKTUFSG = {1: <CH3->, 2: <-CH2->, 3: <-CH<>, 4: <>C<>, 5: <CH2=CH->, 6: <-CH=CH->, 7: <CH2=C<>, 8: <-CH=C<>, 9: <>C=C<>, 13: <CH≡C->, 14: <-C≡C->, 15: <-ACH->, 16: <>AC- (link)>, 17: <>AC- (cond)>, 18: <>AC-CH3>, 19: <>AC-CH2->, 20: <>AC-CH<>, 21: <>AC-C<->, 27: <-CH2- (cy)>, 28: <>CH- (cy)>, 29: <>C< (cy)>, 30: <-CH=CH- (cy)>, 31: <CH2=C< (cy)>, 32: <-CH=C< (cy)>, 34: <-OH(primary)>, 35: <CH3OH>, 36: <H2O>, 37: <>AC-OH>, 38: <(CH2OH)2>, 39: <-O-CH2-CH2-OH>, 40: <-O-CH-CH2-OH>, 41: <-O-CH2-CH-OH>, 42: <CH3-CO->, 43: <-CH2-CO->, 44: <>CH-CO->, 45: <->C-CO->, 48: <-CHO>, 50: <C5H4O2>, 51: <CH3-COO->, 52: <-CH2-COO->, 53: <>CH-COO->, 54: <->C-COO->, 55: <HCOO->, 58: <-COO->, 59: <CH3-O->, 60: <-CH2-O->, 61: <>CH-O->, 62: <->CO->, 63: <-CH2-O- (cy)>, 66: <CH3-NH2>, 67: <-CH2-NH2>, 68: <>CH-NH2>, 69: <->C-NH2>, 71: <CH3-NH->, 72: <-CH2-NH->, 73: <>CH-NH->, 74: <CH3-N<>, 75: <-CH2-N<>, 76: <C5H5N>, 77: <C5H4N->, 78: <C5H3N<>, 79: <>AC-NH2>, 80: <>AC-NH->, 81: <>AC-N<>, 85: <CH3-CN>, 86: <-CH2-CN>, 87: <>CH-CN>, 88: <->C-CN>, 90: <CH2=CH-CN>, 94: <-COOH>, 95: <HCOOH>, 99: <-CH2-Cl>, 100: <>CH-Cl>, 101: <->CCl>, 102: <CH2Cl2>, 103: <-CHCl2>, 104: <>CCl2>, 105: <CHCl3>, 106: <-CCl3>, 107: <CCl4>, 108: <Cl(C=C)>, 109: <>AC-Cl>, 111: <CHF3>, 112: <-CF3>, 113: <-CHF2>, 114: <>CF2>, 115: <-CH2F>, 116: <>CH-F>, 117: <->CF>, 118: <>AC-F>, 120: <CCl3F>, 121: <-CCl2F>, 122: <HCCl2F>, 123: <-HCClF>, 124: <-CClF2>, 125: <HCClF2>, 126: <CClF3>, 127: <CCl2F2>, 128: <-I>, 130: <-Br>, 132: <CH3-NO2>, 133: <-CH2-NO2>, 134: <>CH-NO2>, 135: <->C-NO2>, 136: <>AC-NO2>, 138: <CH3-SH>, 139: <-CH2-SH>, 140: <>CH-SH>, 141: <->C-SH>, 142: <CH3-S->, 143: <-CH2-S->, 144: <>CH-S->, 145: <->C-S->, 146: <CS2>, 147: <THIOPHENE>, 148: <C4H3S->, 149: <C4H2S<>, 153: <DMSO>, 161: <DMF>, 162: <-CON(CH3)2>, 163: <-CON(CH2)(CH3)->, 164: <HCON(CH2)2<>, 165: <-CON(CH2)2<>, 166: <-CONH(CH3)>, 167: <HCONH(CH2)->, 168: <-CONH(CH2)->, 169: <-CONH2>, 195: <NMP>, 196: <MORPHOLIN>, 197: <SiH3->, 198: <-SiH2->, 199: <>SiH->, 200: <>Si<>, 201: <-SiH2-O->, 202: <>SiH-O->, 203: <->Si-O->, 204: <-OH(secondary)>, 205: <-OH(tertiary)>}¶
- thermo.unifac.NISTKTUFMG = {1: ('C', [1, 2, 3, 4]), 2: ('C=C', [5, 6, 7, 8, 9]), 3: ('ACH', [15, 16, 17]), 4: ('ACCH2', [18, 19, 20, 21]), 5: ('OH', [34, 204, 205]), 6: ('CH2OH', [35]), 7: ('H2O', [36]), 8: ('ACOH', [37]), 9: ('CH2CO', [42, 43, 44, 45]), 10: ('CHO', [48]), 11: ('CCOO', [51, 52, 53, 54]), 12: ('HCOO', [55]), 13: ('CH2O', [59, 60, 61, 62]), 14: ('CNH2', [66, 67, 68, 69]), 15: ('(C)2NH', [71, 72, 73]), 16: ('(C)3N', [74, 75]), 17: ('ACNH2', [79, 80, 81]), 18: ('Pyridine', [76, 77, 78]), 19: ('CCN', [85, 86, 87, 88]), 20: ('COOH', [94, 95]), 21: ('CCl', [99, 100, 101]), 22: ('CCl2', [102, 103, 104]), 23: ('CCl3', [105, 106]), 24: ('CCl4', [107]), 25: ('ACCl', [109]), 26: ('CNO2', [132, 133, 134, 135]), 27: ('ACNO2', [136]), 28: ('CS2', [146]), 29: ('CH3SH', [138, 139, 140, 141]), 30: ('Furfural', [50]), 31: ('DOH', [38]), 32: ('I', [128]), 33: ('Br', [130]), 34: ('C=-C', [13, 14]), 35: ('DMSO', [153]), 36: ('ACRY', [90]), 37: ('Cl(C=C)', [108]), 38: ('ACF', [118]), 39: ('DMF', [161, 162, 163, 164, 165]), 40: ('CF2', [111, 112, 113, 114, 115, 116, 117]), 41: ('COO', [58]), 42: ('SiH2', [197, 198, 199, 200]), 43: ('SiO', [201, 202, 203]), 44: ('NMP', [195]), 45: ('CClF', [120, 121, 122, 123, 124, 125, 126, 127]), 46: ('CONCH2', [166, 167, 168, 169]), 47: ('OCCOH', [39, 40, 41]), 48: ('CH2S', [142, 143, 144, 145]), 49: ('Morpholin', [196]), 50: ('THIOPHENE', [147, 148, 149]), 51: ('CH2(cyc)', [27, 28, 29]), 52: ('C=C(cyc)', [30, 31, 32])}¶
Compared to storing the values in dict[(int1, int2)] = (values), the dict-in-dict structure is found emperically to take 111608 bytes vs. 79096 bytes, or 30% less memory.
Data for UNIFAC LLE¶
- thermo.unifac.LLEUFSG = {1: <CH3>, 2: <CH2>, 3: <CH>, 4: <C>, 5: <CH2=CH>, 6: <CH=CH>, 7: <CH=C>, 8: <CH2=C>, 9: <ACH>, 10: <AC>, 11: <ACCH3>, 12: <ACCH2>, 13: <ACCH>, 14: <OH>, 15: <P1>, 16: <P2>, 17: <H2O>, 18: <ACOH>, 19: <CH3CO>, 20: <CH2CO>, 21: <CHO>, 22: <Furfural>, 23: <COOH>, 24: <HCOOH>, 25: <CH3COO>, 26: <CH2COO>, 27: <CH3O>, 28: <CH2O>, 29: <CHO>, 30: <FCH2O>, 31: <CH2CL>, 32: <CHCL>, 33: <CCL>, 34: <CH2CL2>, 35: <CHCL2>, 36: <CCL2>, 37: <CHCL3>, 38: <CCL3>, 39: <CCL4>, 40: <ACCL>, 41: <CH3CN>, 42: <CH2CN>, 43: <ACNH2>, 44: <CH3NO2>, 45: <CH2NO2>, 46: <CHNO2>, 47: <ACNO2>, 48: <DOH>, 49: <(HOCH2CH2)2O>, 50: <C5H5N>, 51: <C5H4N>, 52: <C5H3N>, 53: <CCl2=CHCl>, 54: <HCONHCH3>, 55: <DMF>, 56: <(CH2)4SO2>, 57: <DMSO>}¶
- thermo.unifac.LLEMG = {1: ('CH2', [1, 2, 3, 4]), 2: ('C=C', [5, 6, 7, 8]), 3: ('ACH', [9, 10]), 4: ('ACCH2', [11, 12, 13]), 5: ('OH', [14]), 6: ('P1', [15]), 7: ('P2', [16]), 8: ('H2O', [17]), 9: ('ACOH', [18]), 10: ('CH2CO', [19, 20]), 11: ('CHO', [21]), 12: ('Furfural', [22]), 13: ('COOH', [23, 24]), 14: ('CCOO', [25, 26]), 15: ('CH2O', [27, 28, 29, 30]), 16: ('CCL', [31, 32, 33]), 17: ('CCL2', [34, 35, 36]), 18: ('CCL3', [37, 38]), 19: ('CCL4', [39]), 20: ('ACCL', [40]), 21: ('CCN', [41, 42]), 22: ('ACNH2', [43]), 23: ('CNO2', [44, 45, 46]), 24: ('ACNO2', [47]), 25: ('DOH', [48]), 26: ('DEOH', [49]), 27: ('PYRIDINE', [50, 51, 52]), 28: ('TCE', [53]), 29: ('MFA', [54]), 30: ('DMFA', [55]), 31: ('TMS', [56]), 32: ('DMSO', [57])}¶
Larsen, Bent L., Peter Rasmussen, and Aage Fredenslund. “A Modified UNIFAC Group-Contribution Model for Prediction of Phase Equilibria and Heats of Mixing.” Industrial & Engineering Chemistry Research 26, no. 11 (November 1, 1987): 2274-86. https://doi.org/10.1021/ie00071a018.
Data for Lyngby UNIFAC¶
- thermo.unifac.LUFSG = {1: <CH3>, 2: <CH2>, 3: <CH>, 4: <C>, 5: <CH2=CH>, 6: <CH=CH>, 7: <CH2=C>, 8: <CH=C>, 9: <C=C>, 10: <ACH>, 11: <AC>, 12: <OH>, 13: <CH3OH>, 14: <H2O>, 15: <CH3CO>, 16: <CH2CO>, 17: <CHO>, 18: <CH3COO>, 19: <CH2COO>, 20: <CH3O>, 21: <CH2O>, 22: <CHO>, 23: <THF>, 24: <NH2>, 25: <CH3NH>, 26: <CH2NH>, 27: <CHNH>, 28: <CH3N>, 29: <CH2N>, 30: <ANH2>, 31: <C5H5N>, 32: <C5H4N>, 33: <C5H3N>, 34: <CH3CN>, 35: <CH2CN>, 36: <COOH>, 37: <CH2CL>, 38: <CHCL>, 39: <CCL>, 40: <CH2CL2>, 41: <CHCL2>, 42: <CCL2>, 43: <CHCL3>, 44: <CCL3>, 45: <CCL4>}¶
- thermo.unifac.LUFMG = {1: ('CH2', [1, 2, 3, 4]), 2: ('C=C', [5, 6, 7, 8, 9]), 3: ('ACH', [10, 11]), 4: ('OH', [12]), 5: ('CH3OH', [13]), 6: ('H2O', [14]), 7: ('CH2CO', [15, 16]), 8: ('CHO', [17]), 9: ('CCOO', [18, 19]), 10: ('CH2O', [20, 21, 22, 23]), 11: ('NH2', [24]), 12: ('CNH2NG', [25, 26, 27]), 13: ('CH2N', [28, 29]), 14: ('ANH2', [30]), 15: ('PYRIDINE', [31, 32, 33]), 16: ('CCN', [34, 35]), 17: ('COOH', [36]), 18: ('CCL', [37, 38, 39]), 19: ('CCL2', [40, 41, 42]), 20: ('CCL3', [43, 44]), 21: ('CCL4', [45])}¶
Data for PSRK UNIFAC¶
- thermo.unifac.PSRKSG = {1: <CH3>, 2: <CH2>, 3: <CH>, 4: <C>, 5: <CH2=CH>, 6: <CH=CH>, 7: <CH2=C>, 8: <CH=C>, 9: <ACH>, 10: <AC>, 11: <ACCH3>, 12: <ACCH2>, 13: <ACCH>, 14: <OH>, 15: <CH3OH>, 16: <H2O>, 17: <ACOH>, 18: <CH3CO>, 19: <CH2CO>, 20: <CHO>, 21: <CH3COO>, 22: <CH2COO>, 23: <HCOO>, 24: <CH3O>, 25: <CH2O>, 26: <CHO>, 27: <THF>, 28: <CH3NH2>, 29: <CH2NH2>, 30: <CHNH2>, 31: <CH3NH>, 32: <CH2NH>, 33: <CHNH>, 34: <CH3N>, 35: <CH2N>, 36: <ACNH2>, 37: <C5H5N>, 38: <C5H4N>, 39: <C5H3N>, 40: <CH3CN>, 41: <CH2CN>, 42: <COOH>, 43: <HCOOH>, 44: <CH2CL>, 45: <CHCL>, 46: <CCL>, 47: <CH2CL2>, 48: <CHCL2>, 49: <CCL2>, 50: <CHCL3>, 51: <CCL3>, 52: <CCL4>, 53: <ACCL>, 54: <CH3NO2>, 55: <CH2NO2>, 56: <CHNO2>, 57: <ACNO2>, 58: <CS2>, 59: <CH3SH>, 60: <CH2SH>, 61: <FURFURAL>, 62: <DOH>, 63: <I>, 64: <BR>, 65: <CH=-C>, 66: <C=-C>, 67: <DMSO>, 68: <ACRY>, 69: <CL-(C=C)>, 70: <C=C>, 71: <ACF>, 72: <DMF>, 73: <HCON(CH2)2>, 74: <CF3>, 75: <CF2>, 76: <CF>, 77: <COO>, 78: <SIH3>, 79: <SIH2>, 80: <SIH>, 81: <SI>, 82: <SIH2O>, 83: <SIHO>, 84: <SIO>, 85: <NMP>, 86: <CCL3F>, 87: <CCL2F>, 88: <HCCL2F>, 89: <HCCLF>, 90: <CCLF2>, 91: <HCCLF2>, 92: <CCLF3>, 93: <CCL2F2>, 94: <AMH2>, 95: <AMHCH3>, 96: <AMHCH2>, 97: <AM(CH3)2>, 98: <AMCH3CH2>, 99: <AM(CH2)2>, 100: <C2H5O2>, 101: <C2H4O2>, 102: <CH3S>, 103: <CH2S>, 104: <CHS>, 105: <MORPH>, 106: <C4H4S>, 107: <C4H3S>, 108: <C4H2S>, 109: <H2C=CH2>, 110: <CH=-CH>, 111: <NH3>, 112: <CO>, 113: <H2>, 114: <H2S>, 115: <N2>, 116: <AR>, 117: <CO2>, 118: <CH4>, 119: <O2>, 120: <D2>, 121: <SO2>, 122: <NO>, 123: <N2O>, 124: <SF6>, 125: <HE>, 126: <NE>, 127: <KR>, 128: <XE>, 129: <HF>, 130: <HCL>, 131: <HBR>, 132: <HI>, 133: <COS>, 134: <CHSH>, 135: <CSH>, 136: <H2COCH>, 137: <HCOCH>, 138: <HCOC>, 139: <H2COCH2>, 140: <H2COC>, 141: <COC>, 142: <F2>, 143: <CL2>, 144: <BR2>, 145: <HCN>, 146: <NO2>, 147: <CF4>, 148: <O3>, 149: <CLNO>, 152: <CNH2>}¶
- thermo.unifac.PSRKMG = {1: ('CH2', [1, 2, 3, 4]), 2: ('C=C', [5, 6, 7, 8, 70, 109]), 3: ('ACH', [9, 10]), 4: ('ACCH2', [11, 12, 13]), 5: ('OH', [14]), 6: ('CH3OH', [15]), 7: ('H2O', [16]), 8: ('ACOH', [17]), 9: ('CH2CO', [18, 19]), 10: ('CHO', [20]), 11: ('CCOO', [21, 22]), 12: ('HCOO', [23]), 13: ('CH2O', [24, 25, 26, 27]), 14: ('CNH2', [28, 29, 30, 152]), 15: ('CNH', [31, 32, 33]), 16: ('(C)3N', [34, 35]), 17: ('ACNH2', [36]), 18: ('PYRIDINE', [37, 38, 39]), 19: ('CCN', [40, 41]), 20: ('COOH', [42, 43]), 21: ('CCL', [44, 45, 46]), 22: ('CCL2', [47, 48, 49]), 23: ('CCL3', [50, 51]), 24: ('CCL4', [52]), 25: ('ACCL', [53]), 26: ('CNO2', [54, 55, 56]), 27: ('ACNO2', [57]), 28: ('CS2', [58]), 29: ('CH3SH', [59, 60, 134, 135]), 30: ('FURFURAL', [61]), 31: ('DOH', [62]), 32: ('I', [63]), 33: ('BR', [64]), 34: ('C=-C', [65, 66, 110]), 35: ('DMSO', [67]), 36: ('ACRY', [68]), 37: ('CLCC', [69]), 38: ('ACF', [71]), 39: ('DMF', [72, 73]), 40: ('CF2', [74, 75, 76]), 41: ('COO', [77]), 42: ('SIH2', [78, 79, 80, 81]), 43: ('SIO', [82, 83, 84]), 44: ('NMP', [85]), 45: ('CCLF', [86, 87, 88, 89, 90, 91, 92, 93]), 46: ('CON (AM)', [94, 95, 96, 97, 98, 99]), 47: ('OCCOH', [100, 101]), 48: ('CH2S', [102, 103, 104]), 49: ('MORPH', [105]), 50: ('THIOPHEN', [106, 107, 108]), 51: ('EPOXY', [136, 137, 138, 139, 140, 141]), 55: ('NH3', [111]), 56: ('CO2', [117]), 57: ('CH4', [118]), 58: ('O2', [119]), 59: ('AR', [116]), 60: ('N2', [115]), 61: ('H2S', [114]), 62: ('H2', [113, 120]), 63: ('CO', [112]), 65: ('SO2', [121]), 66: ('NO', [122]), 67: ('N2O', [123]), 68: ('SF6', [124]), 69: ('HE', [125]), 70: ('NE', [126]), 71: ('KR', [127]), 72: ('XE', [128]), 73: ('HF', [129]), 74: ('HCL', [130]), 75: ('HBR', [131]), 76: ('HI', [132]), 77: ('COS', [133]), 78: ('F2', [142]), 79: ('CL2', [143]), 80: ('BR2', [144]), 81: ('HCN', [145]), 82: ('NO2', [146]), 83: ('CF4', [147]), 84: ('O3', [148]), 85: ('CLNO', [149])}¶
- Magnussen, Thomas, Peter Rasmussen, and Aage Fredenslund. “UNIFAC Parameter Table for Prediction of Liquid-Liquid Equilibriums.”
Industrial & Engineering Chemistry Process Design and Development 20, no. 2 (April 1, 1981): 331-39. https://doi.org/10.1021/i200013a024.
Data for VTPR UNIFAC¶
- thermo.unifac.VTPRSG = {1: <CH3>, 2: <CH2>, 3: <CH>, 4: <C>, 5: <CH2=CH>, 6: <CH=CH>, 7: <CH2=C>, 8: <CH=C>, 9: <ACH>, 10: <AC>, 11: <ACCH3>, 12: <ACCH2>, 13: <ACCH>, 14: <OH(P)>, 15: <CH3OH>, 16: <H2O>, 17: <ACOH>, 18: <CH3CO>, 19: <CH2CO>, 20: <CHO>, 21: <CH3COO>, 22: <CH2COO>, 23: <HCOO>, 24: <CH3O>, 25: <CH2O>, 26: <CHO>, 27: <THF>, 28: <CH3NH2>, 29: <CH2NH2>, 30: <CHNH2>, 31: <CH3NH>, 32: <CH2NH>, 33: <CHNH>, 34: <CH3N>, 35: <CH2N>, 36: <ACNH2>, 40: <CH3CN>, 41: <CH2CN>, 44: <CH2CL>, 45: <CHCL>, 46: <CCL>, 47: <CH2CL2>, 48: <CHCL2>, 49: <CCL2>, 50: <CHCL3>, 51: <CCL3>, 52: <CCL4>, 53: <ACCL>, 54: <CH3NO2>, 55: <CH2NO2>, 56: <CHNO2>, 58: <CS2>, 59: <CH3SH>, 60: <CH2SH>, 61: <FURFURAL>, 62: <DOH>, 63: <I>, 64: <BR>, 67: <DMSO>, 70: <C=C>, 72: <DMF>, 73: <HCON(..>, 78: <CY-CH2>, 79: <CY-CH>, 80: <CY-C>, 81: <OH(S)>, 82: <OH(T)>, 83: <CY-CH2O>, 84: <TRIOXAN>, 85: <CNH2>, 86: <NMP>, 87: <NEP>, 88: <NIPP>, 89: <NTBP>, 97: <Allene>, 98: <=CHCH=>, 99: <=CCH=>, 107: <H2COCH>, 108: <COCH>, 109: <HCOCH>, 116: <AC-CHO>, 119: <H2COCH2>, 129: <CHCOO>, 139: <CF2H>, 140: <CF2H2>, 142: <CF2Cl>, 143: <CF2Cl2>, 146: <CF4>, 148: <CF3Br>, 153: <H2COC>, 180: <CHCOO>, 250: <H2C=CH2>, 300: <NH3>, 301: <CO>, 302: <H2>, 303: <H2S>, 304: <N2>, 305: <Ar>, 306: <CO2>, 307: <CH4>, 308: <O2>, 309: <D2>, 310: <SO2>, 312: <N2O>, 314: <He>, 315: <Ne>, 319: <HCl>, 345: <Hg>}¶
- thermo.unifac.VTPRMG = {1: ('CH2', [1, 2, 3, 4]), 2: ('H2C=CH2', [5, 6, 7, 8, 70, 97, 98, 99, 250]), 3: ('ACH', [9, 10]), 4: ('ACCH2', [11, 12, 13]), 5: ('OH', [14, 81, 82]), 6: ('CH3OH', [15]), 7: ('H2O', [16]), 8: ('ACOH', [17]), 9: ('CH2CO', [18, 19]), 10: ('CHO', [20]), 11: ('CCOO', [21, 22, 129, 180]), 12: ('HCOO', [23]), 13: ('CH2O', [24, 25, 26]), 14: ('CH2NH2', [28, 29, 30, 85]), 15: ('CH2NH', [31, 32, 33]), 16: ('(C)3N', [34, 35]), 17: ('ACNH2', [36]), 19: ('CH2CN', [40, 41]), 21: ('CCL', [44, 45, 46]), 22: ('CCL2', [47, 48, 49]), 23: ('CCL3', [51]), 24: ('CCL4', [52]), 25: ('ACCL', [53]), 26: ('CNO2', [54, 55, 56]), 28: ('CS2', [58]), 29: ('CH3SH', [59, 60]), 30: ('FURFURAL', [61]), 31: ('DOH', [62]), 32: ('I', [63]), 33: ('BR', [64]), 35: ('DMSO', [67]), 39: ('DMF', [72, 73]), 42: ('CY-CH2', [78, 79, 80]), 43: ('CY-CH2O', [27, 83, 84]), 45: ('CHCL3', [50]), 46: ('CY-CONC', [86, 87, 88, 89]), 53: ('EPOXIDES', [107, 108, 109, 119, 153]), 57: ('AC-CHO', [116]), 68: ('CF2H', [139, 140]), 70: ('CF2Cl2', [142, 143, 148]), 73: ('CF4', [146]), 150: ('NH3', [300]), 151: ('CO2', [306]), 152: ('CH4', [307]), 153: ('O2', [308]), 154: ('Ar', [305]), 155: ('N2', [304]), 156: ('H2S', [303]), 157: ('D2', [302, 309]), 158: ('CO', [301]), 160: ('SO2', [310]), 162: ('N2O', [312]), 164: ('He', [314]), 165: ('Ne', [315]), 169: ('HCl', [319]), 185: ('Hg', [345])}¶