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(*, xs, rs, qs, Qs, vs, T=298.15, psi_coeffs=None, psi_abc=None, psi_a=None, psi_b=None, psi_c=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
Tfloat

Temperature, [K]

xslist[float]

Mole fractions, [-]

rslist[float]

r parameters ri=k=1nνkRkr_i = \sum_{k=1}^{n} \nu_k R_k, [-]

qslist[float]

q parameters qi=k=1nνkQkq_i = \sum_{k=1}^{n}\nu_k Q_k, [-]

Qslist[float]

Q parameter for each subgroup; subgroups are not required to but are suggested to be sorted from lowest number to highest number, [-]

vslist[list[float]]

Indexed by [subgroup][count], this variable is the count of each subgroups in each compound, [-]

psi_abctuple(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 or psi_a and psi_b and psi_c is required, [-]

psi_coeffslist[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 or psi_a and psi_b and psi_c is required, [-]

psi_alist[list[float]], optional

psi a term interaction parameters between each subgroup; indexed [subgroup][subgroup]. Only one of psi_abc or psi_coeffs or psi_a and psi_b and psi_c is required, [-]

psi_blist[list[float]], optional

psi b term interaction parameters between each subgroup; indexed [subgroup][subgroup]. Only one of psi_abc or psi_coeffs or psi_a and psi_b and psi_c is required, [-]

psi_clist[list[float]], optional

psi c term interaction parameters between each subgroup; indexed [subgroup][subgroup]. Only one of psi_abc or psi_coeffs or psi_a and psi_b and psi_c is required, [-]

versionint, 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)

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

Attributes
Tfloat

Temperature, [K]

xslist[float]

Mole fractions, [-]

Methods

CpE()

Calculate and return the first temperature derivative of excess enthalpy of a liquid phase using an activity coefficient model.

Fis()

Calculate the FiF_i terms used in calculating the combinatorial part.

Fs()

Computes the following:

Fs_pure()

Computes the following:

GE()

Calculate the excess Gibbs energy with the UNIFAC model.

Gs()

Computes the following:

Gs_pure()

Computes the following:

HE()

Calculate and return the excess entropy of a liquid phase using an activity coefficient model.

Hs()

Computes the following:

Hs_pure()

Computes the following:

SE()

Calculates the excess entropy of a liquid phase using an activity coefficient model.

Theta_Psi_sum_invs()

Computes the following term for each group k, size number of groups.

Theta_Psi_sums()

Computes the following term for each group k, size number of groups.

Theta_pure_Psi_sum_invs()

Computes the following term for each group k, size number of groups.

Thetas()

Calculate the Θm\Theta_m parameters used in calculating the residual part.

Thetas_pure()

Calculate the Θm\Theta_m parameters for each chemical in the mixture as a pure species, used in calculating the residual part.

Vis()

Calculate the ViV_i terms used in calculating the combinatorial part.

Vis_modified()

Calculate the ViV_i' terms used in calculating the combinatorial part.

Ws()

Computes the following for each k and each i, indexed by [k][i] k is in groups, and i is in components.

Xs()

Calculate the XmX_m parameters used in calculating the residual part.

Xs_pure()

Calculate the XmX_m parameters for each chemical in the mixture as a pure species, used in calculating the residual part.

as_json([cache, option])

Method to create a JSON-friendly representation of the Gibbs Excess model which can be stored, and reloaded later.

d2Fis_dxixjs()

Calculate the second mole fraction derivative of the FiF_i 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.

d2GE_dTdxs()

Calculate the first composition derivative and temperature derivative of excess Gibbs energy with the UNIFAC model.

d2GE_dxixjs()

Calculate the second composition derivative of excess Gibbs energy with the UNIFAC model.

d2Thetas_dxixjs()

Calculate the mole fraction derivatives of the Θm\Theta_m parameters.

d2Vis_dxixjs()

Calculate the second mole fraction derivative of the ViV_i terms used in calculating the combinatorial part.

d2Vis_modified_dxixjs()

Calculate the second mole fraction derivative of the ViV_i' terms used in calculating the combinatorial part.

d2lnGammas_subgroups_dT2()

Calculate the second temperature derivative of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases's composition and temperature.

d2lnGammas_subgroups_dTdxs()

Calculate the temperature and mole fraction derivatives of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases's composition and temperature.

d2lnGammas_subgroups_dxixjs()

Calculate the second mole fraction derivatives of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases's composition and temperature.

d2lnGammas_subgroups_pure_dT2()

Calculate the second temperature derivative of lnΓk\ln \Gamma_k pure component parameters for the phase; depends on the phases's temperature only.

d2lngammas_c_dT2()

Second temperature derivatives of the combinatorial part of the UNIFAC model.

d2lngammas_c_dTdx()

Second temperature derivative and first mole fraction derivative of the combinatorial part of the UNIFAC model.

d2lngammas_c_dxixjs()

Second composition derivative of the combinatorial part of the UNIFAC model.

d2lngammas_dT2()

Calculates the second temperature derivative of the residual part of the UNIFAC model.

d2lngammas_r_dT2()

Calculates the second temperature derivative of the residual part of the UNIFAC model.

d2lngammas_r_dTdxs()

Calculates the first mole fraction derivative of the temperature derivative of the residual part of the UNIFAC model.

d2lngammas_r_dxixjs()

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.

d2psis_dT2()

Calculate the Ψ\Psi term second temperature derivative matrix for all groups interacting with all other groups.

d3Fis_dxixjxks()

Calculate the third mole fraction derivative of the FiF_i terms used in calculating the combinatorial part.

d3GE_dT3()

Calculate the third temperature derivative of excess Gibbs energy with the UNIFAC model.

d3Vis_dxixjxks()

Calculate the third mole fraction derivative of the ViV_i terms used in calculating the combinatorial part.

d3Vis_modified_dxixjxks()

Calculate the third mole fraction derivative of the ViV_i' terms used in calculating the combinatorial part.

d3lnGammas_subgroups_dT3()

Calculate the third temperature derivative of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases's composition and temperature.

d3lnGammas_subgroups_pure_dT3()

Calculate the third temperature derivative of lnΓk\ln \Gamma_k pure component parameters for the phase; depends on the phases's temperature only.

d3lngammas_c_dT3()

Third temperature derivatives of the combinatorial part of the UNIFAC model.

d3lngammas_c_dxixjxks()

Third composition derivative of the combinatorial part of the UNIFAC model.

d3lngammas_dT3()

Calculates the third temperature derivative of the residual part of the UNIFAC model.

d3lngammas_r_dT3()

Calculates the third temperature derivative of the residual part of the UNIFAC model.

d3psis_dT3()

Calculate the Ψ\Psi term third temperature derivative matrix for all groups interacting with all other groups.

dFis_dxs()

Calculate the mole fraction derivative of the FiF_i 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.

dThetas_dxs()

Calculate the mole fraction derivatives of the Θm\Theta_m parameters.

dVis_dxs()

Calculate the mole fraction derivative of the ViV_i terms used in calculating the combinatorial part.

dVis_modified_dxs()

Calculate the mole fraction derivative of the ViV_i' terms used in calculating the combinatorial part.

dgammas_dT()

Calculates the first temperature derivative of activity coefficients with the UNIFAC model.

dgammas_dns()

Calculate and return the mole number derivative of activity coefficients of a liquid phase using an activity coefficient model.

dgammas_dxs()

Calculates the first mole fraction derivative of activity coefficients with the UNIFAC model.

dlnGammas_subgroups_dT()

Calculate the first temperature derivative of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases's composition and temperature.

dlnGammas_subgroups_dxs()

Calculate the mole fraction derivatives of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases's composition and temperature.

dlnGammas_subgroups_pure_dT()

Calculate the first temperature derivative of lnΓk\ln \Gamma_k pure component parameters for the phase; depends on the phases's temperature only.

dlngammas_c_dT()

Temperature derivatives of the combinatorial part of the UNIFAC model.

dlngammas_c_dxs()

First composition derivative of the combinatorial part of the UNIFAC model.

dlngammas_dT()

Calculates the first temperature derivative of the residual part of the UNIFAC model.

dlngammas_r_dT()

Calculates the first temperature derivative of the residual part of the UNIFAC model.

dlngammas_r_dxs()

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 Ψ\Psi term first temperature derivative matrix for all groups interacting with all other groups.

exact_hash()

Method to calculate and return a hash representing the exact state of the object.

from_json(json_repr[, cache])

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.

lnGammas_subgroups()

Calculate the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases's composition and temperature.

lnGammas_subgroups_pure()

Calculate the lnΓk\ln \Gamma_k pure component parameters for the phase; depends on the phases's temperature only.

lngammas()

Calculate and return the natural logarithm of the activity coefficients of a liquid phase using an activity coefficient model.

lngammas_c()

Calculates the combinatorial part of the UNIFAC model.

lngammas_r()

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 Ψ\Psi 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.

Theta_pure_Psi_sums

Thetas_sum_inv

VSXS

Xs_sum_inv

d2GE_dT2_numerical

d2GE_dTdxs_numerical

d2GE_dxixjs_numerical

d3GE_dT2dxs_numerical

d3GE_dT3_numerical

d3GE_dTdxixjs_numerical

d4GE_dT2dxixjs_numerical

d4GE_dT3dxs_numerical

d4GE_dT4_numerical

d5GE_dT3dxixjs_numerical

dGE_dT_numerical

dGE_dxs_numerical

gammas_dGE_dxs

gammas_from_args

gammas_numerical

Fis()[source]

Calculate the FiF_i terms used in calculating the combinatorial part. A function of mole fractions and the parameters q only.

Fi=qijqjxjF_i = \frac{q_i}{\sum_j q_j x_j}

This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.

Returns
Fislist[float]

F terms size number of components, [-]

Fs()[source]

Computes the following:

F(k)=mgrθmψm,kTF(k) = \sum_m^{gr} \theta_m \frac{\partial \psi_{m,k}}{\partial T}
Fs_pure()[source]

Computes the following:

F(k)=mgrθmψm,kTF(k) = \sum_m^{gr} \theta_m \frac{\partial \psi_{m,k}}{\partial T}
GE()[source]

Calculate the excess Gibbs energy with the UNIFAC model.

GE=RTixi(lnγic+lnγir)G^E = RT\sum_i x_i \left(\ln \gamma_i^c + \ln \gamma_i^r \right)

For the VTPR model, the combinatorial component is set to zero.

Returns
GEfloat

Excess Gibbs energy, [J/mol]

Gs()[source]

Computes the following:

G(k)=mgrθm2ψm,kT2G(k) = \sum_m^{gr} \theta_m \frac{\partial^2 \psi_{m,k}}{\partial T^2}
Gs_pure()[source]

Computes the following:

G(k)=mgrθm2ψm,kT2G(k) = \sum_m^{gr} \theta_m \frac{\partial^2 \psi_{m,k}}{\partial T^2}
Hs()[source]

Computes the following:

H(k)=mgrθm3ψm,kT3H(k) = \sum_m^{gr} \theta_m \frac{\partial^3 \psi_{m,k}}{\partial T^3}
Hs_pure()[source]

Computes the following:

H(k)=mgrθm3ψm,kT3H(k) = \sum_m^{gr} \theta_m \frac{\partial^3 \psi_{m,k}}{\partial T^3}
Theta_Psi_sum_invs()[source]

Computes the following term for each group k, size number of groups.

U(k)=1mΘmΨm,kU(k) = \frac{1}{\sum_m \Theta_m \Psi_{m,k}}
Theta_Psi_sums()[source]

Computes the following term for each group k, size number of groups.

mΘmΨmk\sum_m \Theta_m \Psi_{mk}
Theta_pure_Psi_sum_invs()[source]

Computes the following term for each group k, size number of groups.

U(k)=1mΘmΨm,kU(k) = \frac{1}{\sum_m \Theta_m \Psi_{m,k}}
Thetas()[source]

Calculate the Θm\Theta_m parameters used in calculating the residual part. A function of mole fractions and group counts only.

Θm=QmXmnQnXn\Theta_m = \frac{Q_m X_m}{\sum_{n} Q_n X_n}
Returns
Thetaslist[float]

Θm\Theta_m terms, size number of subgroups, [-]

Thetas_pure()[source]

Calculate the Θm\Theta_m parameters for each chemical in the mixture as a pure species, used in calculating the residual part. A function of group counts only.

Θm=QmXmnQnXn\Theta_m = \frac{Q_m X_m}{\sum_{n} Q_n X_n}
Returns
Thetas_purelist[list[float]]

Θm\Theta_m terms, size number of components by number of subgroups and indexed in that order, [-]

Vis()[source]

Calculate the ViV_i terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.

Vi=rijrjxjV_i = \frac{r_i}{\sum_j r_j x_j}

This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.

Returns
Vislist[float]

V terms size number of components, [-]

Vis_modified()[source]

Calculate the ViV_i' terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.

Vi=rinjrjnxjV_i' = \frac{r_i^n}{\sum_j r_j^n x_j}

This is used in the UNIFAC Dortmund and UNIFAC-NIST model with n=0.75, and the Lyngby model with n=2/3.

Returns
Vis_modifiedlist[float]

Modified V terms size number of components, [-]

Ws()[source]

Computes the following for each k and each i, indexed by [k][i] k is in groups, and i is in components.

W(k,i)=mgrψm,kθmxiW(k,i) = \sum_m^{gr} \psi_{m,k} \frac{\partial \theta_m}{\partial x_i}
Xs()[source]

Calculate the XmX_m parameters used in calculating the residual part. A function of mole fractions and group counts only.

Xm=jνmjxjjnνnjxjX_m = \frac{ \sum_j \nu^j_m x_j}{\sum_j \sum_n \nu_n^j x_j}
Returns
Xslist[float]

XmX_m terms, size number of subgroups, [-]

Xs_pure()[source]

Calculate the XmX_m 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.

Xm=νmngrνnX_m = \frac{\nu_m}{\sum^{gr}_n \nu_n}
Returns
Xs_purelist[list[float]]

XmX_m terms, size number of subgroups by number of components and indexed in that order, [-]

d2Fis_dxixjs()[source]

Calculate the second mole fraction derivative of the FiF_i terms used in calculating the combinatorial part. A function of mole fractions and the parameters q only.

Fixjxk=2qiqjqkGsum3\frac{\partial F_i}{\partial x_j \partial x_k} = 2 q_i q_j q_k G_{sum}^3
Gsum=1jqjxjG_{sum} = \frac{1}{\sum_j q_j x_j}

This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.

Returns
d2Fis_dxixjslist[list[list[float]]]

F terms size number of components by number of components by number of components, [-]

d2GE_dT2()[source]

Calculate the second temperature derivative of excess Gibbs energy with the UNIFAC model.

2GET2=RTixi2lnγirT2+2RixilnγirT\frac{\partial^2 G^E}{\partial T^2} = RT\sum_i x_i \frac{\partial^2 \ln \gamma_i^r}{\partial T^2} + 2R\sum_i x_i \frac{\partial \ln \gamma_i^r}{\partial T}
Returns
d2GE_dT2float

Second temperature derivative of excess Gibbs energy, [J/mol/K^2]

d2GE_dTdxs()[source]

Calculate the first composition derivative and temperature derivative of excess Gibbs energy with the UNIFAC model.

2GETxi=RT(lnγirT+jxjlnγjrxi)+R[lnγicxi+lnγirxi+jxj(lnγjcxi+lnγjrxi)]\frac{\partial^2 G^E}{\partial T\partial x_i} = RT\left(\frac{\partial \ln \gamma_i^r}{\partial T} + \sum_j x_j \frac{\partial \ln \gamma_j^r}{\partial x_i} \right) + R\left[ \frac{\partial \ln \gamma_i^c}{\partial x_i} + \frac{\partial \ln \gamma_i^r}{\partial x_i} + \sum_j x_j \left( \frac{\partial \ln \gamma_j^c}{\partial x_i} + \frac{\partial \ln \gamma_j^r}{\partial x_i}\right)\right]
Returns
dGE_dxslist[float]

First composition derivative and first temperature derivative of excess Gibbs energy, [J/mol/K]

d2GE_dxixjs()[source]

Calculate the second composition derivative of excess Gibbs energy with the UNIFAC model.

2GExjxk=RT[i(lnγicxjxk+lnγirxjxk)+lnγjcxk+lnγjrxk+lnγkcxj+lnγkrxj]\frac{\partial^2 G^E}{\partial x_j \partial x_k} = RT \left[\sum_i \left( \frac{\partial \ln \gamma_i^c}{\partial x_j \partial x_k} + \frac{\partial \ln \gamma_i^r}{\partial x_j \partial x_k} \right) + \frac{\partial \ln \gamma_j^c}{\partial x_k} + \frac{\partial \ln \gamma_j^r}{\partial x_k} + \frac{\partial \ln \gamma_k^c}{\partial x_j} + \frac{\partial \ln \gamma_k^r}{\partial x_j}\right]
Returns
d2GE_dxixjslist[list[float]]

Second composition derivative of excess Gibbs energy, [J/mol]

d2Thetas_dxixjs()[source]

Calculate the mole fraction derivatives of the Θm\Theta_m parameters. A function of mole fractions and group counts only.

2Θixjxk=QinQn(νx)sum,n[F(ν)sum,jνi,kF(ν)sum,kνi,j+2F2(ν)sum,j(ν)sum,k(νx)sum,i+F(νx)sum,i[n(2FQn(ν)sum,j(ν)sum,k(νx)sum,n+Qn(ν)sum,jνn,k+Qn(ν)sum,kνn,j)]ngrQn(νx)sum,n+2(νx)sum,i(ngr[FQn(ν)sum,j(νx)sum,n+Qnνn,j])(ngr[FQn(ν)sum,k(νx)sum,n+Qnνn,k])(ngrQn(νx)sum,n)2νi,j(ngrFQn(ν)sum,k(νx)sum,n+Qnνn,k)(ngrQn(νx)sum,n)νi,k(ngrFQn(ν)sum,j(νx)sum,n+Qnνn,j)(ngrQn(νx)sum,n)+F(ν)sum,j(νx)sum,i(ngrFQn(ν)sum,k(νx)sum,n+Qnνn,k)(ngrQn(νx)sum,n)+F(ν)sum,k(νx)sum,i(ngrFQn(ν)sum,j(νx)sum,n+Qnνn,j)(ngrQn(νx)sum,n)]\frac{\partial^2 \Theta_i}{\partial x_j \partial x_k} = \frac{Q_i}{\sum_n Q_n (\nu x)_{sum,n}}\left[ -F(\nu)_{sum,j} \nu_{i,k} - F (\nu)_{sum,k}\nu_{i,j} + 2F^2(\nu)_{sum,j} (\nu)_{sum,k} (\nu x)_{sum,i} + \frac{F (\nu x)_{sum,i}\left[ \sum_n(-2 F Q_n (\nu)_{sum,j} (\nu)_{sum,k} (\nu x)_{sum,n} + Q_n (\nu)_{sum,j} \nu_{n,k} + Q_n (\nu)_{sum,k}\nu_{n,j} )\right] } {\sum_n^{gr} Q_n (\nu x)_{sum,n} } + \frac{2(\nu x)_{sum,i}(\sum_n^{gr}[-FQ_n (\nu)_{sum,j} (\nu x)_{sum,n} + Q_n \nu_{n,j}]) (\sum_n^{gr}[-FQ_n (\nu)_{sum,k} (\nu x)_{sum,n} + Q_n \nu_{n,k}]) } {\left( \sum_n^{gr} Q_n (\nu x)_{sum,n} \right)^2} - \frac{\nu_{i,j}(\sum_n^{gr} -FQ_n (\nu)_{sum,k} (\nu x)_{sum,n} + Q_n \nu_{n,k} )} {\left( \sum_n^{gr} Q_n (\nu x)_{sum,n} \right)} - \frac{\nu_{i,k}(\sum_n^{gr} -FQ_n (\nu)_{sum,j} (\nu x)_{sum,n} + Q_n \nu_{n,j} )} {\left( \sum_n^{gr} Q_n (\nu x)_{sum,n} \right)} + \frac{F(\nu)_{sum,j} (\nu x)_{sum,i} (\sum_n^{gr} -FQ_n (\nu)_{sum,k} (\nu x)_{sum,n} + Q_n \nu_{n,k})} {\left(\sum_n^{gr} Q_n (\nu x)_{sum,n} \right)} + \frac{F(\nu)_{sum,k} (\nu x)_{sum,i} (\sum_n^{gr} -FQ_n (\nu)_{sum,j} (\nu x)_{sum,n} + Q_n \nu_{n,j})} {\left(\sum_n^{gr} Q_n (\nu x)_{sum,n} \right)} \right]
G=1jQjXjG = \frac{1}{\sum_j Q_j X_j}
F=1jnνnjxjF = \frac{1}{\sum_j \sum_n \nu_n^j x_j}
(ν)sum,i=jνj,i(\nu)_{sum,i} = \sum_j \nu_{j,i}
(νx)sum,i=jνi,jxj(\nu x)_{sum,i} = \sum_j \nu_{i,j}x_j
Returns
d2Thetas_dxixjslist[list[list[float]]]

Θm\Theta_m terms, size number of subgroups by mole fractions and indexed in that order, [-]

d2Vis_dxixjs()[source]

Calculate the second mole fraction derivative of the ViV_i terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.

Vixjxk=2rirjrkVsum3\frac{\partial V_i}{\partial x_j \partial x_k} = 2 r_i r_j r_k V_{sum}^3
Vsum=1jrjxjV_{sum} = \frac{1}{\sum_j r_j x_j}

This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.

Returns
d2Vis_dxixjslist[list[list[float]]]

V terms size number of components by number of components by number of components, [-]

d2Vis_modified_dxixjs()[source]

Calculate the second mole fraction derivative of the ViV_i' terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.

Vixjxk=2rinrjnrknVsum3\frac{\partial V_i'}{\partial x_j \partial x_k} = 2 r_i^n r_j^n r_k^n V_{sum}^3
Vsum=1jrjnxjV_{sum} = \frac{1}{\sum_j r_j^n x_j}

This is used in the UNIFAC Dortmund and UNIFAC-NIST model with n=0.75, and the Lyngby model with n=2/3.

Returns
d2Vis_modified_dxixjslist[list[list[float]]]

V’ terms size number of components by number of components by number of components, [-]

d2lnGammas_subgroups_dT2()[source]

Calculate the second temperature derivative of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases’s composition and temperature.

2lnΓiT2=Qi[Z(i)G(i)F(i)2Z(i)2+j(θjZ(j)2ψi,jTZ(j)2(G(j)θjψi,j+2Fjθjψi,jT)+2Z(j)3F(j)2θjψi,j)]\frac{\partial^2 \ln \Gamma_i}{\partial T^2} = -Q_i\left[ Z(i)G(i) - F(i)^2 Z(i)^2 + \sum_j\left( \theta_j Z(j)\frac{\partial^2 \psi_{i,j}}{\partial T} - Z(j)^2 \left(G(j)\theta_j \psi_{i,j} + 2 F_j \theta_j \frac{\partial \psi_{i,j}}{\partial T}\right) + 2Z(j)^3F(j)^2 \theta_j \psi_{i,j} \right)\right]
F(k)=mgrθmψm,kTF(k) = \sum_m^{gr} \theta_m \frac{\partial \psi_{m,k}}{\partial T}
G(k)=mgrθm2ψm,kT2G(k) = \sum_m^{gr} \theta_m \frac{\partial^2 \psi_{m,k}}{\partial T^2}
Z(k)=1mΘmΨm,kZ(k) = \frac{1}{\sum_m \Theta_m \Psi_{m,k}}
Returns
d2lnGammas_subgroups_dT2list[float]

Second temperature derivative of ln Gamma parameters for each subgroup, size number of subgroups, [1/K^2]

d2lnGammas_subgroups_dTdxs()[source]

Calculate the temperature and mole fraction derivatives of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases’s composition and temperature.

2lnΓkxiT=Qk(D(k,i)Z(k)B(k)W(k,i)Z(k)2+mgr(Z(m)θmxiψk,mT)mgr(B(m)Z(m)2ψk,mθmxi)mgr(D(m,i)Z(m)2θmψk,m)mgr(W(m,i)Z(m)2θmψk,mT)+mgr2B(m)W(m,i)Z(m)3θmψk,m)\frac{\partial^2 \ln \Gamma_k}{\partial x_i \partial T} = -Q_k\left( D(k,i) Z(k) - B(k)W(k,i) Z(k)^2 + \sum_m^{gr} (Z(m) \frac{\partial \theta_m}{\partial x_i}\frac{\partial \psi_{k,m}}{\partial T}) -\sum_m^{gr} (B(m) Z(m)^2 \psi_{k,m} \frac{\partial \theta_m}{\partial x_i}) -\sum_m^{gr}(D(m,i) Z(m)^2 \theta_m \psi_{k,m}) - \sum_m^{gr} (W(m,i) Z(m)^2 \theta_m \frac{\partial \psi_{k,m}}{\partial T}) + \sum_m^{gr} 2 B(m) W(m,i) Z(m)^3 \theta_m \psi_{k,m} \right)

The following groups are used as follows to simplfy the number of evaluations:

W(k,i)=mgrψm,kθmxiW(k,i) = \sum_m^{gr} \psi_{m,k} \frac{\partial \theta_m}{\partial x_i}
Z(k)=1mΘmΨmkZ(k) = \frac{1}{\sum_m \Theta_m \Psi_{mk}}
F(k)=mgrθmψm,kTF(k) = \sum_m^{gr} \theta_m \frac{\partial \psi_{m,k}}{\partial T}

In the below expression, k` refers to a group, and i refers to a component.

D(k,i)=mgrθmxiψm,kTD(k,i) = \sum_m^{gr} \frac{\partial \theta_m}{\partial x_i} \frac{\partial \psi_{m,k}}{\partial T}
Returns
d2lnGammas_subgroups_dTdxslist[list[float]]

Temperature and mole fraction derivatives of Gamma parameters for each subgroup, size number of subgroups by number of components and indexed in that order, [1/K]

d2lnGammas_subgroups_dxixjs()[source]

Calculate the second mole fraction derivatives of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases’s composition and temperature.

2lnΓkxixj=Qk(Z(k)K(k,i,j)mgrZ(m)2K(m,i,j)θmψk,mW(k,i)W(k,j)Z(k)2+mgrZmψk,m2θmxixjm(W(m,j)Z(m)2ψk,mθmxi+W(m,i)Z(m)2ψ(k,m)θmxj)+mgr2W(m,i)W(m,j)Z(m)3θmψk,m)\frac{\partial^2 \ln \Gamma_k}{\partial x_i \partial x_j} = -Q_k\left( -Z(k) K(k,i,j) - \sum_m^{gr} Z(m)^2 K(m,i,j)\theta_m \psi_{k,m} -W(k,i) W(k,j) Z(k)^2 + \sum_m^{gr} Z_m \psi_{k,m} \frac{\partial^2 \theta_m}{\partial x_i \partial x_j} - \sum_m \left(W(m,j) Z(m)^2 \psi_{k,m} \frac{\partial \theta_m}{\partial x_i} + W(m,i) Z(m)^2 \psi(k,m) \frac{\partial \theta_m}{\partial x_j}\right) + \sum_m^{gr} 2 W(m,i) W(m,j) Z(m)^3 \theta_m \psi_{k,m}\right)

The following groups are used as follows to simplfy the number of evaluations:

W(k,i)=mgrψm,kθmxiW(k,i) = \sum_m^{gr} \psi_{m,k} \frac{\partial \theta_m}{\partial x_i}
Z(k)=1mΘmΨmkZ(k) = \frac{1}{\sum_m \Theta_m \Psi_{mk}}
K(k,i,j)=mgrψm,k2θmxixjK(k, i, j) = \sum_m^{gr} \psi_{m,k} \frac{\partial^2 \theta_m}{\partial x_i \partial x_j}
Returns
d2lnGammas_subgroups_dxixjslist[list[list[float]]]

Second mole fraction derivatives of Gamma parameters for each subgroup, size number of components by number of components by number of subgroups and indexed in that order, [-]

d2lnGammas_subgroups_pure_dT2()[source]

Calculate the second temperature derivative of lnΓk\ln \Gamma_k pure component parameters for the phase; depends on the phases’s temperature only.

2lnΓiT2=Qi[Z(i)G(i)F(i)2Z(i)2+j(θjZ(j)2ψi,jTZ(j)2(G(j)θjψi,j+2Fjθjψi,jT)+2Z(j)3F(j)2θjψi,j)]\frac{\partial^2 \ln \Gamma_i}{\partial T^2} = -Q_i\left[ Z(i)G(i) - F(i)^2 Z(i)^2 + \sum_j\left( \theta_j Z(j)\frac{\partial^2 \psi_{i,j}}{\partial T} - Z(j)^2 \left(G(j)\theta_j \psi_{i,j} + 2 F_j \theta_j \frac{\partial \psi_{i,j}}{\partial T}\right) + 2Z(j)^3F(j)^2 \theta_j \psi_{i,j} \right)\right]
F(k)=mgrθmψm,kTF(k) = \sum_m^{gr} \theta_m \frac{\partial \psi_{m,k}}{\partial T}
G(k)=mgrθm2ψm,kT2G(k) = \sum_m^{gr} \theta_m \frac{\partial^2 \psi_{m,k}}{\partial T^2}
Z(k)=1mΘmΨm,kZ(k) = \frac{1}{\sum_m \Theta_m \Psi_{m,k}}

In this model, the Θ\Theta values come from the UNIFAC.Thetas_pure method, where each compound is assumed to be pure.

Returns
d2lnGammas_subgroups_pure_dT2list[list[float]]

Second temperature derivative of ln Gamma parameters for each subgroup, size number of subgroups by number of components and indexed in that order, [1/K^2]

d2lngammas_c_dT2()[source]

Second temperature derivatives of the combinatorial part of the UNIFAC model. Zero in all variations.

2lnγicT2=0\frac{\partial^2 \ln \gamma_i^c}{\partial T^2} = 0
Returns
d2lngammas_c_dT2list[float]

Combinatorial lngammas term second temperature derivatives, size number of components, [-]

d2lngammas_c_dTdx()[source]

Second temperature derivative and first mole fraction derivative of the combinatorial part of the UNIFAC model. Zero in all variations.

3lnγicT2xj=0\frac{\partial^3 \ln \gamma_i^c}{\partial T^2 \partial x_j} = 0
Returns
d2lngammas_c_dTdxlist[list[float]]

Combinatorial lngammas term second temperature derivatives, size number of components by number of components, [-]

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 ViV_i' with ViV_i.

lnγicxjxk=5qi(d2dxkdxjVi+Vid2dxkdxjFiFi+ddxjFiddxkViFi+ddxkFiddxjViFi2ViddxjFiddxkFiFi2Vi+(ddxjViViddxjFiFi)ddxkViVi2+d2dxkdxjViFi(ddxjViViddxjFiFi)ddxkFiFiViVid2dxkdxjFiFi2ddxjFiddxkViFi2ddxkFiddxjViFi2+2ViddxjFiddxkFiFi3)d2dxkdxjVi+d2dxkdxjViViddxjViddxkViVi2\frac{\partial \ln \gamma^c_i}{\partial x_j \partial x_k} = 5 q_{i} \left(\frac{- \frac{d^{2}}{d x_{k}d x_{j}} V_{i} + \frac{V_{i} \frac{d^{2}}{d x_{k}d x_{j}} F_{i}}{F_{i}} + \frac{\frac{d}{d x_{j}} F_{i} \frac{d}{d x_{k}} V_{i}}{F_{i}} + \frac{\frac{d}{d x_{k}} F_{i} \frac{d}{d x_{j}} V_{i}}{F_{i}} - \frac{2 V_{i} \frac{d}{d x_{j}} F_{i} \frac{d}{d x_{k}} F_{i}}{F_{i}^{2}}}{V_{i}} + \frac{\left( \frac{d}{d x_{j}} V_{i} - \frac{V_{i} \frac{d}{d x_{j}} F_{i}} {F_{i}}\right) \frac{d}{d x_{k}} V_{i}}{V_{i}^{2}} + \frac{\frac{d^{2}}{d x_{k}d x_{j}} V_{i}}{F_{i}} - \frac{\left( \frac{d}{d x_{j}} V_{i} - \frac{V_{i} \frac{d}{d x_{j}} F_{i}}{ F_{i}}\right) \frac{d}{d x_{k}} F_{i}}{F_{i} V_{i}} - \frac{V_{i} \frac{d^{2}}{d x_{k}d x_{j}} F_{i}}{F_{i}^{2}} - \frac{\frac{d} {d x_{j}} F_{i} \frac{d}{d x_{k}} V_{i}}{F_{i}^{2}} - \frac{\frac{d}{d x_{k}} F_{i} \frac{d}{d x_{j}} V_{i}}{F_{i}^{2}} + \frac{2 V_{i} \frac{d}{d x_{j}} F_{i} \frac{d}{d x_{k}} F_{i}} {F_{i}^{3}}\right) - \frac{d^{2}}{d x_{k}d x_{j}} Vi' + \frac{\frac{d^{2}}{d x_{k}d x_{j}} Vi'}{Vi'} - \frac{\frac{d} {d x_{j}} Vi' \frac{d}{d x_{k}} Vi'}{Vi'^{2}}

For the Lyngby model, the following equations are used:

2lnγicxjxk=2Vixjxk+1Vi2Vixjxk1(Vi)2VixjVixk\frac{\partial^2 \ln \gamma^c_i}{\partial x_j \partial x_k} = -\frac{\partial^2 V_i'}{\partial x_j \partial x_k} + \frac{1}{V_i'} \frac{\partial^2 V_i'}{\partial x_j \partial x_k} - \frac{1}{\left(V_i'\right)^2} \frac{\partial V_i'}{\partial x_j} \frac{\partial V_i'}{\partial x_k}
Returns
d2lngammas_c_dxixjslist[list[list[float]]]

Combinatorial lngammas term second composition derivative, size number of components by number of components by number of components, [-]

d2lngammas_dT2()

Calculates the second temperature derivative of the residual part of the UNIFAC model.

2lnγirT2=kgrνk(i)[2lnΓkT22lnΓk(i)T2]\frac{\partial^2 \ln \gamma_i^r}{\partial T^2} = \sum_{k}^{gr} \nu_k^{(i)} \left[ \frac{\partial^2 \ln \Gamma_k}{\partial T^2} - \frac{\partial^2 \ln \Gamma_k^{(i)}}{\partial T^2} \right]

where the second Gamma is the pure-component Gamma of group k in component i.

Returns
d2lngammas_r_dT2list[float]

Residual lngammas terms second temperature derivative, size number of components [1/K^2]

d2lngammas_r_dT2()[source]

Calculates the second temperature derivative of the residual part of the UNIFAC model.

2lnγirT2=kgrνk(i)[2lnΓkT22lnΓk(i)T2]\frac{\partial^2 \ln \gamma_i^r}{\partial T^2} = \sum_{k}^{gr} \nu_k^{(i)} \left[ \frac{\partial^2 \ln \Gamma_k}{\partial T^2} - \frac{\partial^2 \ln \Gamma_k^{(i)}}{\partial T^2} \right]

where the second Gamma is the pure-component Gamma of group k in component i.

Returns
d2lngammas_r_dT2list[float]

Residual lngammas terms second temperature derivative, size number of components [1/K^2]

d2lngammas_r_dTdxs()[source]

Calculates the first mole fraction derivative of the temperature derivative of the residual part of the UNIFAC model.

2lnγirxjT=mgrνm(i)2lnΓmxjT\frac{\partial^2 \ln \gamma_i^r}{\partial x_j \partial T} = \sum_{m}^{gr} \nu_m^{(i)} \frac{\partial^2 \ln \Gamma_m} {\partial x_j \partial T}
Returns
d2lngammas_r_dTdxslist[list[float]]

First mole fraction derivative and temperature derivative of residual lngammas terms, size number of components by number of components [-]

d2lngammas_r_dxixjs()[source]

Calculates the second mole fraction derivative of the residual part of the UNIFAC model.

2lnγirxj2=mgrνm(i)2lnΓmxj2\frac{\partial^2 \ln \gamma_i^r}{\partial x_j^2} = \sum_{m}^{gr} \nu_m^{(i)} \frac{\partial^2 \ln \Gamma_m}{\partial x_j^2}
Returns
d2lngammas_r_dxixjslist[list[list[float]]]

Second mole fraction derivative of the residual lngammas terms, size number of components by number of components by number of components [-]

d2psis_dT2()[source]

Calculate the Ψ\Psi 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;

2ΨmnT2=(2cmn+2(2Tcmn+bmn)T+(2Tcmn+bmnT2cmn+Tbmn+amnT)2T2(T2cmn+Tbmn+amn)T2)eT2cmn+Tbmn+amnTT\frac{\partial^2 \Psi_{mn}}{\partial T^2} = \frac{\left(- 2 c_{mn} + \frac{2 \left(2 T c_{mn} + b_{mn}\right)}{T} + \frac{\left(2 T c_{mn} + b_{mn} - \frac{T^{2} c_{mn} + T b_{mn} + a_{mn}}{T} \right)^{2}}{T} - \frac{2 \left(T^{2} c_{mn} + T b_{mn} + a_{mn} \right)}{T^{2}}\right) e^{- \frac{T^{2} c_{mn} + T b_{mn} + a_{mn}} {T}}}{T}

Only the first, a coefficient, is used in the original UNIFAC model as well as the UNIFAC-LLE model, so the expression simplifies to:

2ΨmnT2=amn(2+amnT)eamnTT3\frac{\partial^2 \Psi_{mn}}{\partial T^2} = \frac{a_{mn} \left(-2 + \frac{a_{mn}}{T}\right) e^{- \frac{a_{mn}}{T}}}{T^{3}}

For the Lyngby model, the second temperature derivative is:

2ΨmkT2=(2a2+2a3ln(T0T)+a3+(a2+a3ln(T0T)a1+a2(TT0)+a3(Tln(T0T)+TT0)T)22(a1+a2(TT0)+a3(Tln(T0T)+TT0))T)ea1+a2(TT0)+a3(Tln(T0T)+TT0)TT2\frac{\partial^2 \Psi_{mk}}{\partial T^2} = \frac{\left(2 a_{2} + 2 a_{3} \ln{\left(\frac{T_{0}}{T} \right)} + a_{3} + \left(a_{2} + a_{3} \ln{\left(\frac{T_{0}}{T} \right)} - \frac{a_{1} + a_{2} \left(T - T_{0}\right) + a_{3} \left(T \ln{\left(\frac{T_{0}}{T} \right)} + T - T_{0}\right)}{T}\right)^{2} - \frac{2 \left(a_{1} + a_{2} \left(T - T_{0}\right) + a_{3} \left(T \ln{\left( \frac{T_{0}}{T} \right)} + T - T_{0}\right)\right)}{T}\right) e^{- \frac{a_{1} + a_{2} \left(T - T_{0}\right) + a_{3} \left( T \ln{\left(\frac{T_{0}}{T} \right)} + T - T_{0}\right)}{T}}} {T^{2}}

with T0=298.15T_0 = 298.15 K and the a coefficients are specific to each pair of main groups, and they are asymmetric, so a0,mka0,kma_{0,mk} \ne a_{0,km}.

Returns
d2psis_dT2list[list[float]]

Second temperature derivative of`psi` terms, size subgroups x subgroups [-]

d3Fis_dxixjxks()[source]

Calculate the third mole fraction derivative of the FiF_i terms used in calculating the combinatorial part. A function of mole fractions and the parameters q only.

Fixjxkxm=6qiqjqkqmGsum4\frac{\partial F_i}{\partial x_j \partial x_k \partial x_m} = -6 q_i q_j q_k q_m G_{sum}^4
Gsum=1jqjxjG_{sum} = \frac{1}{\sum_j q_j x_j}

This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.

Returns
d3Fis_dxixjxkslist[list[list[list[float]]]]

F terms size number of components by number of components by number of components by number of components, [-]

d3GE_dT3()[source]

Calculate the third temperature derivative of excess Gibbs energy with the UNIFAC model.

3GET3=RTixi3lnγirT3+3Rixi2lnγirT2\frac{\partial^3 G^E}{\partial T^3} = RT\sum_i x_i \frac{\partial^3 \ln \gamma_i^r}{\partial T^3} + 3R\sum_i x_i \frac{\partial^2 \ln \gamma_i^r}{\partial T^2}
Returns
d3GE_dT3float

Third temperature derivative of excess Gibbs energy, [J/mol/K^3]

d3Vis_dxixjxks()[source]

Calculate the third mole fraction derivative of the ViV_i terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.

Vixjxkxm=6rirjrkrmVsum4\frac{\partial V_i}{\partial x_j \partial x_k \partial x_m} = -6 r_i r_j r_k r_m V_{sum}^4
Vsum=1jrjxjV_{sum} = \frac{1}{\sum_j r_j x_j}

This is used in the UNIFAC, UNIFAC-LLE, UNIFAC Dortmund, UNIFAC-NIST, and PSRK models.

Returns
d3Vis_dxixjxkslist[list[list[list[float]]]]

V terms size number of components by number of components by number of components by number of components, [-]

d3Vis_modified_dxixjxks()[source]

Calculate the third mole fraction derivative of the ViV_i' terms used in calculating the combinatorial part. A function of mole fractions and the parameters r only.

Vixjxkxm=6rinrjnrknrmnVsum4\frac{\partial V_i'}{\partial x_j \partial x_k \partial x_m} = -6 r_i^n r_j^n r_k^n r_m^n V_{sum}^4
Vsum=1jrjxjV_{sum} = \frac{1}{\sum_j r_j x_j}

This is used in the UNIFAC Dortmund and UNIFAC-NIST model with n=0.75, and the Lyngby model with n=2/3.

Returns
d3Vis_modified_dxixjxkslist[list[list[list[float]]]]

V’ terms size number of components by number of components by number of components by number of components, [-]

d3lnGammas_subgroups_dT3()[source]

Calculate the third temperature derivative of the lnΓk\ln \Gamma_k parameters for the phase; depends on the phases’s composition and temperature.

3lnΓiT3=Qi[H(i)Z(i)2F(i)3Z(i)3+3F(i)G(i)Z(i)2+(θjZ(j)3ψT3+H(j)Z(j)2θ(j)ψi,j6F(j)2Z(j)3θjψi,jT+3F(j)Z(j)2θ(j)2ψi,jT2++3G(j)θ(j)Z(j)2ψi,jT+6F(j)3θ(j)Z(j)4ψi,j6F(j)G(j)θ(j)Z(j)3ψi,j)]\frac{\partial^3 \ln \Gamma_i}{\partial T^3} =Q_i\left[-H(i) Z(i) - 2F(i)^3 Z(i)^3 + 3F(i) G(i) Z(i)^2+ \left( -\theta_j Z(j) \frac{\partial^3 \psi}{\partial T^3} + H(j) Z(j)^2 \theta(j)\psi_{i,j} - 6F(j)^2 Z(j)^3 \theta_j \frac{\partial \psi_{i,j}}{\partial T} + 3 F(j) Z(j)^2 \theta(j) \frac{\partial^2 \psi_{i,j}}{\partial T^2} ++ 3G(j) \theta(j) Z(j)^2 \frac{\partial \psi_{i,j}}{\partial T} + 6F(j)^3 \theta(j) Z(j)^4 \psi_{i,j} - 6F(j) G(j) \theta(j) Z(j)^3 \psi_{i,j} \right) \right]
F(k)=mgrθmψm,kTF(k) = \sum_m^{gr} \theta_m \frac{\partial \psi_{m,k}}{\partial T}
G(k)=mgrθm2ψm,kT2G(k) = \sum_m^{gr} \theta_m \frac{\partial^2 \psi_{m,k}}{\partial T^2}
H(k)=mgrθm3ψm,kT3H(k) = \sum_m^{gr} \theta_m \frac{\partial^3 \psi_{m,k}}{\partial T^3}
Z(k)=1mΘmΨm,kZ(k) = \frac{1}{\sum_m \Theta_m \Psi_{m,k}}
Returns
d3lnGammas_subgroups_dT3list[float]

Third temperature derivative of ln Gamma parameters for each subgroup, size number of subgroups, [1/K^3]

d3lnGammas_subgroups_pure_dT3()[source]

Calculate the third temperature derivative of lnΓk\ln \Gamma_k pure component parameters for the phase; depends on the phases’s temperature only.

3lnΓiT3=Qi[H(i)Z(i)2F(i)3Z(i)3+3F(i)G(i)Z(i)2+(θjZ(j)3ψT3+H(j)Z(j)2θ(j)ψi,j6F(j)2Z(j)3θjψi,jT+3F(j)Z(j)2θ(j)2ψi,jT2++3G(j)θ(j)Z(j)2ψi,jT+6F(j)3θ(j)Z(j)4ψi,j6F(j)G(j)θ(j)Z(j)3ψi,j)]\frac{\partial^3 \ln \Gamma_i}{\partial T^3} =Q_i\left[-H(i) Z(i) - 2F(i)^3 Z(i)^3 + 3F(i) G(i) Z(i)^2+ \left( -\theta_j Z(j) \frac{\partial^3 \psi}{\partial T^3} + H(j) Z(j)^2 \theta(j)\psi_{i,j} - 6F(j)^2 Z(j)^3 \theta_j \frac{\partial \psi_{i,j}}{\partial T} + 3 F(j) Z(j)^2 \theta(j) \frac{\partial^2 \psi_{i,j}}{\partial T^2} ++ 3G(j) \theta(j) Z(j)^2 \frac{\partial \psi_{i,j}}{\partial T} + 6F(j)^3 \theta(j) Z(j)^4 \psi_{i,j} - 6F(j) G(j) \theta(j) Z(j)^3 \psi_{i,j} \right) \right]
F(k)=mgrθmψm,kTF(k) = \sum_m^{gr} \theta_m \frac{\partial \psi_{m,k}}{\partial T}
G(k)=mgrθm2ψm,kT2G(k) = \sum_m^{gr} \theta_m \frac{\partial^2 \psi_{m,k}}{\partial T^2}
H(k)=mgrθm3ψm,kT3H(k) = \sum_m^{gr} \theta_m \frac{\partial^3 \psi_{m,k}}{\partial T^3}
Z(k)=1mΘmΨm,kZ(k) = \frac{1}{\sum_m \Theta_m \Psi_{m,k}}