General cubics¶
The reduced residual Helmholtz energy for the main cubic EOS (van der Waals, Peng-Robinson, and Soave-Redlich-Kwong) can be written in a common form (see https://www.ncbi.nlm.nih.gov/pmc/articles/PMC7365965/)
with the constants given by:
vdW: \(\Delta_1=0\), \(\Delta_2=0\)
SRK: \(\Delta_1=1\), \(\Delta_2=0\)
PR: \(\Delta_1=1+\sqrt{2}\), \(\Delta_2=1-\sqrt{2}\)
The quantities \(a_m\) and \(b_m\) are described (with exact solutions for the numerical coefficients) for each of these EOS in https://pubs.acs.org/doi/abs/10.1021/acs.iecr.1c00847.
The models in teqp are instantiated based on knowledge of the critical temperature, pressure, and acentric factor. Thereafter all quantities are obtained from derivatives of \(\alpha^r\).
The Python class is here: GeneralizedCubic
[1]:
import teqp
teqp.__version__
[1]:
'0.22.0'
[2]:
import json
import CoolProp.CoolProp as CP
# Values taken from http://dx.doi.org/10.6028/jres.121.011
Tc_K = [ 190.564, 154.581, 150.687 ]
pc_Pa = [ 4599200, 5042800, 4863000 ]
acentric = [ 0.011, 0.022, -0.002 ]
# Instantiate Peng-Robinson model
modelPR = teqp.canonical_PR(Tc_K, pc_Pa, acentric)
# Instantiate Soave-Redlich-Kwong model
modelSRK = teqp.canonical_SRK(Tc_K, pc_Pa, acentric)
[3]:
# And you can get information about the model in JSON format
# from the get_meta function
modelPR.get_meta()
[3]:
{'Delta1': 2.414213562373095,
'Delta2': -0.41421356237309515,
'OmegaA': 0.4572355289213822,
'OmegaB': 0.07779607390388846,
'R / J/mol/K': 8.31446261815324,
'kind': 'Peng-Robinson'}
Adjusting k_ij¶
Fine-tuned values of \(k_{ij}\) can be provided when instantiating the model, for Peng-Robinson and SRK. A complete matrix of all the \(k_{ij}\) values must be provided. This allows for asymmetric mixing models in which \(k_{ij}\neq k_{ji}\).
[4]:
k_12 = 0.01
kmat = [[0,k_12,0],[k_12,0,0],[0,0,0]]
teqp.canonical_PR(Tc_K, pc_Pa, acentric, kmat)
teqp.canonical_SRK(Tc_K, pc_Pa, acentric, kmat)
[4]:
<teqp.teqp.AbstractModel at 0x7f3cb6d7cbf0>
Superancillary¶
The superancillary equation gives the co-existing liquid and vapor (orthobaric) densities as a function of temperature. The set of Chebyshev expansions was developed in https://pubs.acs.org/doi/abs/10.1021/acs.iecr.1c00847 . These superancillary equations are more accurate than iterative calculations in double precision arithmetic and also at least 10 times faster to calculate, and cannot fail in iterative routines, even extremely close to the critical point.
The superancillary equation is only exposed for pure fluids to remove ambiguity when considering mixtures. The returned tuple is the liquid and vapor densities
[5]:
teqp.canonical_PR([Tc_K[0]], [pc_Pa[0]], [acentric[0]]).superanc_rhoLV(100)
[5]:
(30846.392909514052, 42.480231719002326)
a and b¶
For the cubic EOS, it can be useful to obtain the a and b parameters directly. The b parameter is particularly useful because 1/b is the maximum allowed density in the EOS
[6]:
import numpy as np
z = np.array([0.3, 0.4, 0.3])
modelPR.get_a(140, z), modelPR.get_b(140, z)
[6]:
(0.1874177858906821, 2.1984349667726406e-05)
alpha functions¶
It can be advantageous to modify the alpha function to allow for more accurate handling of the attractive interactions. Coefficients are tabulated for many species in https://pubs.acs.org/doi/10.1021/acs.jced.7b00967 for the Peng-Robinson EOS with Twu alpha function and the values from the SI of that paper are in the csv file next to this file.
[7]:
import pandas
dfTwu = pandas.read_csv('fitted_Twu_coeffs.csv')
def get_model(INCHIKey):
row = dfTwu.loc[dfTwu['inchikey']==INCHIKey]
if len(row) == 1:
row = row.iloc[0]
Tc_K = row['Tc_K']
pc_MPa = row['pc_MPa']
c = [row['c0'], row['c1'], row['c2']]
# The JSON definition of the EOS,
# here a generic cubic EOS to allow for
# specification of the alpha function(s)
j = {
'kind': 'cubic',
'model': {
'type': 'PR',
'Tcrit / K': [Tc_K],
'pcrit / Pa': [pc_MPa*1e6],
'acentric': [0.1],
'alpha': [{'type': 'Twu', 'c': c}]
}
}
model = teqp.make_model(j)
return model, j
# Hexane
model, j = get_model(INCHIKey='VLKZOEOYAKHREP-UHFFFAOYSA-N')
[8]:
# And how about we calculate the pressure and s^+ = -sr/R at NBP of water
model, j = get_model(INCHIKey='XLYOFNOQVPJJNP-UHFFFAOYSA-N') # WATER
T = 373.1242958476844 # K, NBP of water
rhoL, rhoV = model.superanc_rhoLV(T)
z = np.array([1.0])
pL = rhoL*model.get_R(z)*T*(1.0 + model.get_Ar01(T, rhoL, z))
splusL = model.get_splus(T, rhoL*z)
print(pL, splusL)
102739.27983424198 6.03697343297877
Also implemented in version 0.17 are the alpha functions of Mathias-Copeman.
with
Parameters are tabulated for many fluids in the paper of Horstmann (https://doi.org/10.1016/j.fluid.2004.11.002) for the SRK EOS (only)
[9]:
# Here is an example from Horstmann
j = {
"kind": "cubic",
"model": {
"type": "SRK",
"Tcrit / K": [647.30],
"pcrit / Pa": [22.048321e6],
"acentric": [0.3440],
"alpha": [
{"type": "Mathias-Copeman", "c": [1.07830, -0.58321, 0.54619]}
]
}
}
model = teqp.make_model(j)
T = 373.1242958476844 # K
rhoL, rhoV = model.superanc_rhoLV(T)
z = np.array([1.0])
pL = rhoL*model.get_R(z)*T*(1.0 + model.get_Ar01(T, rhoL, z))
print('And with SRK and Mathias-Copeman parameters:', pL, 'Pa')
And with SRK and Mathias-Copeman parameters: 101639.22259842217 Pa