Cubic Equations of State for Mixtures (thermo.eos_mix)¶
This module contains implementations of most cubic equations of state for mixtures. This includes Peng-Robinson, SRK, Van der Waals, PRSV, TWU and many other variants.
For reporting bugs, adding feature requests, or submitting pull requests, please use the GitHub issue tracker.
Base Class¶
- class thermo.eos_mix.GCEOSMIX[source]¶
Bases:
GCEOSClass for solving a generic pressure-explicit three-parameter cubic equation of state for a mixture. Does not implement any parameters itself; must be subclassed by a mixture equation of state class which subclasses it.
- Attributes:
A_dep_gDeparture molar Helmholtz energy from ideal gas behavior for the gas phase, [J/mol].
A_dep_lDeparture molar Helmholtz energy from ideal gas behavior for the liquid phase, [J/mol].
Cp_minus_Cv_gCp - Cv for the gas phase, [J/mol/K].
Cp_minus_Cv_lCp - Cv for the liquid phase, [J/mol/K].
Jaubert_r1Parameter used in kij interconversion / calculation of Eij G_res^(E, inf)
Jaubert_r2Parameter used in kij interconversion / calculation of Eij G_res^(E, inf)
U_dep_gDeparture molar internal energy from ideal gas behavior for the gas phase, [J/mol].
U_dep_lDeparture molar internal energy from ideal gas behavior for the liquid phase, [J/mol].
V_dep_gDeparture molar volume from ideal gas behavior for the gas phase, [m^3/mol].
V_dep_lDeparture molar volume from ideal gas behavior for the liquid phase, [m^3/mol].
V_g_mpmathThe molar volume of the gas phase calculated with mpmath to a higher precision, [m^3/mol].
V_l_mpmathThe molar volume of the liquid phase calculated with mpmath to a higher precision, [m^3/mol].
VcCritical volume, [m^3/mol].
a_alpha_ijsCalculate and return the matrix .
beta_gIsobaric (constant-pressure) expansion coefficient for the gas phase, [1/K].
beta_lIsobaric (constant-pressure) expansion coefficient for the liquid phase, [1/K].
- c1
- c2
d2H_dep_dT2_gSecond temperature derivative of departure enthalpy with respect to temperature for the gas phase, [(J/mol)/K^2].
d2H_dep_dT2_g_PSecond temperature derivative of departure enthalpy with respect to temperature for the gas phase, [(J/mol)/K^2].
d2H_dep_dT2_g_VSecond temperature derivative of departure enthalpy with respect to temperature at constant volume for the gas phase, [(J/mol)/K^2].
d2H_dep_dT2_lSecond temperature derivative of departure enthalpy with respect to temperature for the liquid phase, [(J/mol)/K^2].
d2H_dep_dT2_l_PSecond temperature derivative of departure enthalpy with respect to temperature for the liquid phase, [(J/mol)/K^2].
d2H_dep_dT2_l_VSecond temperature derivative of departure enthalpy with respect to temperature at constant volume for the liquid phase, [(J/mol)/K^2].
d2H_dep_dTdP_gTemperature and pressure derivative of departure enthalpy at constant pressure then temperature for the gas phase, [(J/mol)/K/Pa].
d2H_dep_dTdP_lTemperature and pressure derivative of departure enthalpy at constant pressure then temperature for the liquid phase, [(J/mol)/K/Pa].
d2P_dT2_PV_gSecond derivative of pressure with respect to temperature twice, but with pressure held constant the first time and volume held constant the second time for the gas phase, [Pa/K^2].
d2P_dT2_PV_lSecond derivative of pressure with respect to temperature twice, but with pressure held constant the first time and volume held constant the second time for the liquid phase, [Pa/K^2].
d2P_dTdP_gSecond derivative of pressure with respect to temperature and, then pressure; and with volume held constant at first, then temperature, for the gas phase, [1/K].
d2P_dTdP_lSecond derivative of pressure with respect to temperature and, then pressure; and with volume held constant at first, then temperature, for the liquid phase, [1/K].
d2P_dTdrho_gDerivative of pressure with respect to molar density, and temperature for the gas phase, [Pa/(K*mol/m^3)].
d2P_dTdrho_lDerivative of pressure with respect to molar density, and temperature for the liquid phase, [Pa/(K*mol/m^3)].
d2P_dVdP_gSecond derivative of pressure with respect to molar volume and then pressure for the gas phase, [mol/m^3].
d2P_dVdP_lSecond derivative of pressure with respect to molar volume and then pressure for the liquid phase, [mol/m^3].
d2P_dVdT_TP_gSecond derivative of pressure with respect to molar volume and then temperature at constant temperature then pressure for the gas phase, [Pa*mol/m^3/K].
d2P_dVdT_TP_lSecond derivative of pressure with respect to molar volume and then temperature at constant temperature then pressure for the liquid phase, [Pa*mol/m^3/K].
d2P_dVdT_gAlias of
GCEOS.d2P_dTdV_gd2P_dVdT_lAlias of
GCEOS.d2P_dTdV_ld2P_drho2_gSecond derivative of pressure with respect to molar density for the gas phase, [Pa/(mol/m^3)^2].
d2P_drho2_lSecond derivative of pressure with respect to molar density for the liquid phase, [Pa/(mol/m^3)^2].
d2S_dep_dT2_gSecond temperature derivative of departure entropy with respect to temperature for the gas phase, [(J/mol)/K^3].
d2S_dep_dT2_g_VSecond temperature derivative of departure entropy with respect to temperature at constant volume for the gas phase, [(J/mol)/K^3].
d2S_dep_dT2_lSecond temperature derivative of departure entropy with respect to temperature for the liquid phase, [(J/mol)/K^3].
d2S_dep_dT2_l_VSecond temperature derivative of departure entropy with respect to temperature at constant volume for the liquid phase, [(J/mol)/K^3].
d2S_dep_dTdP_gTemperature and pressure derivative of departure entropy at constant pressure then temperature for the gas phase, [(J/mol)/K^2/Pa].
d2S_dep_dTdP_lTemperature and pressure derivative of departure entropy at constant pressure then temperature for the liquid phase, [(J/mol)/K^2/Pa].
d2T_dP2_gSecond partial derivative of temperature with respect to pressure (constant volume) for the gas phase, [K/Pa^2].
d2T_dP2_lSecond partial derivative of temperature with respect to pressure (constant temperature) for the liquid phase, [K/Pa^2].
d2T_dPdV_gSecond partial derivative of temperature with respect to pressure (constant volume) and then volume (constant pressure) for the gas phase, [K*mol/(Pa*m^3)].
d2T_dPdV_lSecond partial derivative of temperature with respect to pressure (constant volume) and then volume (constant pressure) for the liquid phase, [K*mol/(Pa*m^3)].
d2T_dPdrho_gDerivative of temperature with respect to molar density, and pressure for the gas phase, [K/(Pa*mol/m^3)].
d2T_dPdrho_lDerivative of temperature with respect to molar density, and pressure for the liquid phase, [K/(Pa*mol/m^3)].
d2T_dV2_gSecond partial derivative of temperature with respect to volume (constant pressure) for the gas phase, [K*mol^2/m^6].
d2T_dV2_lSecond partial derivative of temperature with respect to volume (constant pressure) for the liquid phase, [K*mol^2/m^6].
d2T_dVdP_gSecond partial derivative of temperature with respect to pressure (constant volume) and then volume (constant pressure) for the gas phase, [K*mol/(Pa*m^3)].
d2T_dVdP_lSecond partial derivative of temperature with respect to pressure (constant volume) and then volume (constant pressure) for the liquid phase, [K*mol/(Pa*m^3)].
d2T_drho2_gSecond derivative of temperature with respect to molar density for the gas phase, [K/(mol/m^3)^2].
d2T_drho2_lSecond derivative of temperature with respect to molar density for the liquid phase, [K/(mol/m^3)^2].
d2V_dP2_gSecond partial derivative of volume with respect to pressure (constant temperature) for the gas phase, [m^3/(Pa^2*mol)].
d2V_dP2_lSecond partial derivative of volume with respect to pressure (constant temperature) for the liquid phase, [m^3/(Pa^2*mol)].
d2V_dPdT_gSecond partial derivative of volume with respect to pressure (constant temperature) and then presssure (constant temperature) for the gas phase, [m^3/(K*Pa*mol)].
d2V_dPdT_lSecond partial derivative of volume with respect to pressure (constant temperature) and then presssure (constant temperature) for the liquid phase, [m^3/(K*Pa*mol)].
d2V_dT2_gSecond partial derivative of volume with respect to temperature (constant pressure) for the gas phase, [m^3/(mol*K^2)].
d2V_dT2_lSecond partial derivative of volume with respect to temperature (constant pressure) for the liquid phase, [m^3/(mol*K^2)].
d2V_dTdP_gSecond partial derivative of volume with respect to pressure (constant temperature) and then presssure (constant temperature) for the gas phase, [m^3/(K*Pa*mol)].
d2V_dTdP_lSecond partial derivative of volume with respect to pressure (constant temperature) and then presssure (constant temperature) for the liquid phase, [m^3/(K*Pa*mol)].
d2a_alpha_dT2_dnsHelper method for calculating the mole number derivatives of d2a_alpha_dT2.
d2a_alpha_dT2_dzsHelper method for calculating the mole number derivatives of d2a_alpha_dT2.
d2a_alpha_dT2_ijsCalculate and return the matrix of the second temperature derivatives of the alpha terms.
d2a_alpha_dTdP_g_VDerivative of the temperature derivative of a_alpha with respect to pressure at constant volume (varying T) for the gas phase, [J^2/mol^2/Pa^2/K].
d2a_alpha_dTdP_l_VDerivative of the temperature derivative of a_alpha with respect to pressure at constant volume (varying T) for the liquid phase, [J^2/mol^2/Pa^2/K].
d2a_alpha_dninjsHelper method for calculating the second partial molar derivatives of a_alpha (hessian).
d2a_alpha_dzizjsHelper method for calculating the second composition derivatives of a_alpha (hessian).
d2b_dninjsHelper method for calculating the second partial mole number derivatives of b.
d2b_dzizjsHelper method for calculating the second partial mole fraction derivatives of b.
d2rho_dP2_gSecond derivative of molar density with respect to pressure for the gas phase, [(mol/m^3)/Pa^2].
d2rho_dP2_lSecond derivative of molar density with respect to pressure for the liquid phase, [(mol/m^3)/Pa^2].
d2rho_dPdT_gSecond derivative of molar density with respect to pressure and temperature for the gas phase, [(mol/m^3)/(K*Pa)].
d2rho_dPdT_lSecond derivative of molar density with respect to pressure and temperature for the liquid phase, [(mol/m^3)/(K*Pa)].
d2rho_dT2_gSecond derivative of molar density with respect to temperature for the gas phase, [(mol/m^3)/K^2].
d2rho_dT2_lSecond derivative of molar density with respect to temperature for the liquid phase, [(mol/m^3)/K^2].
d3a_alpha_dT3Method to calculate the third temperature derivative of , [J^2/mol^2/Pa/K^3].
d3a_alpha_dninjnksHelper method for calculating the third mole number derivatives of a_alpha.
d3a_alpha_dzizjzksHelper method for calculating the third composition derivatives of a_alpha.
d3b_dninjnksHelper method for calculating the third partial mole number derivatives of b.
d3b_dzizjzksHelper method for calculating the third partial mole fraction derivatives of b.
d3delta_dzizjzksHelper method for calculating the third composition derivatives of delta.
d3epsilon_dzizjzksHelper method for calculating the third composition derivatives of epsilon.
dH_dep_dP_gDerivative of departure enthalpy with respect to pressure for the gas phase, [(J/mol)/Pa].
dH_dep_dP_g_VDerivative of departure enthalpy with respect to pressure at constant volume for the liquid phase, [(J/mol)/Pa].
dH_dep_dP_lDerivative of departure enthalpy with respect to pressure for the liquid phase, [(J/mol)/Pa].
dH_dep_dP_l_VDerivative of departure enthalpy with respect to pressure at constant volume for the gas phase, [(J/mol)/Pa].
dH_dep_dT_gDerivative of departure enthalpy with respect to temperature for the gas phase, [(J/mol)/K].
dH_dep_dT_g_VDerivative of departure enthalpy with respect to temperature at constant volume for the gas phase, [(J/mol)/K].
dH_dep_dT_lDerivative of departure enthalpy with respect to temperature for the liquid phase, [(J/mol)/K].
dH_dep_dT_l_VDerivative of departure enthalpy with respect to temperature at constant volume for the liquid phase, [(J/mol)/K].
dH_dep_dV_g_PDerivative of departure enthalpy with respect to volume at constant pressure for the gas phase, [J/m^3].
dH_dep_dV_g_TDerivative of departure enthalpy with respect to volume at constant temperature for the gas phase, [J/m^3].
dH_dep_dV_l_PDerivative of departure enthalpy with respect to volume at constant pressure for the liquid phase, [J/m^3].
dH_dep_dV_l_TDerivative of departure enthalpy with respect to volume at constant temperature for the gas phase, [J/m^3].
dP_drho_gDerivative of pressure with respect to molar density for the gas phase, [Pa/(mol/m^3)].
dP_drho_lDerivative of pressure with respect to molar density for the liquid phase, [Pa/(mol/m^3)].
dS_dep_dP_gDerivative of departure entropy with respect to pressure for the gas phase, [(J/mol)/K/Pa].
dS_dep_dP_g_VDerivative of departure entropy with respect to pressure at constant volume for the gas phase, [(J/mol)/K/Pa].
dS_dep_dP_lDerivative of departure entropy with respect to pressure for the liquid phase, [(J/mol)/K/Pa].
dS_dep_dP_l_VDerivative of departure entropy with respect to pressure at constant volume for the liquid phase, [(J/mol)/K/Pa].
dS_dep_dT_gDerivative of departure entropy with respect to temperature for the gas phase, [(J/mol)/K^2].
dS_dep_dT_g_VDerivative of departure entropy with respect to temperature at constant volume for the gas phase, [(J/mol)/K^2].
dS_dep_dT_lDerivative of departure entropy with respect to temperature for the liquid phase, [(J/mol)/K^2].
dS_dep_dT_l_VDerivative of departure entropy with respect to temperature at constant volume for the liquid phase, [(J/mol)/K^2].
dS_dep_dV_g_PDerivative of departure entropy with respect to volume at constant pressure for the gas phase, [J/K/m^3].
dS_dep_dV_g_TDerivative of departure entropy with respect to volume at constant temperature for the gas phase, [J/K/m^3].
dS_dep_dV_l_PDerivative of departure entropy with respect to volume at constant pressure for the liquid phase, [J/K/m^3].
dS_dep_dV_l_TDerivative of departure entropy with respect to volume at constant temperature for the gas phase, [J/K/m^3].
dT_drho_gDerivative of temperature with respect to molar density for the gas phase, [K/(mol/m^3)].
dT_drho_lDerivative of temperature with respect to molar density for the liquid phase, [K/(mol/m^3)].
dZ_dP_gDerivative of compressibility factor with respect to pressure for the gas phase, [1/Pa].
dZ_dP_lDerivative of compressibility factor with respect to pressure for the liquid phase, [1/Pa].
dZ_dT_gDerivative of compressibility factor with respect to temperature for the gas phase, [1/K].
dZ_dT_lDerivative of compressibility factor with respect to temperature for the liquid phase, [1/K].
da_alpha_dP_g_VDerivative of the a_alpha with respect to pressure at constant volume (varying T) for the gas phase, [J^2/mol^2/Pa^2].
da_alpha_dP_l_VDerivative of the a_alpha with respect to pressure at constant volume (varying T) for the liquid phase, [J^2/mol^2/Pa^2].
da_alpha_dT_dnsHelper method for calculating the mole number derivatives of da_alpha_dT.
da_alpha_dT_dzsHelper method for calculating the composition derivatives of da_alpha_dT.
da_alpha_dT_ijsCalculate and return the matrix for the temperature derivatives of the alpha terms.
da_alpha_dnsHelper method for calculating the mole number derivatives of a_alpha.
da_alpha_dzsHelper method for calculating the composition derivatives of a_alpha.
db_dnsHelper method for calculating the mole number derivatives of b.
db_dzsHelper method for calculating the composition derivatives of b.
dbeta_dP_gDerivative of isobaric expansion coefficient with respect to pressure for the gas phase, [1/(Pa*K)].
dbeta_dP_lDerivative of isobaric expansion coefficient with respect to pressure for the liquid phase, [1/(Pa*K)].
dbeta_dT_gDerivative of isobaric expansion coefficient with respect to temperature for the gas phase, [1/K^2].
dbeta_dT_lDerivative of isobaric expansion coefficient with respect to temperature for the liquid phase, [1/K^2].
dfugacity_dP_gDerivative of fugacity with respect to pressure for the gas phase, [-].
dfugacity_dP_lDerivative of fugacity with respect to pressure for the liquid phase, [-].
dfugacity_dT_gDerivative of fugacity with respect to temperature for the gas phase, [Pa/K].
dfugacity_dT_lDerivative of fugacity with respect to temperature for the liquid phase, [Pa/K].
dna_alpha_dT_dnsHelper method for calculating the mole number derivatives of da_alpha_dT.
dna_alpha_dnsHelper method for calculating the partial molar derivatives of a_alpha.
dnb_dnsHelper method for calculating the partial molar derivative of b.
dphi_dP_gDerivative of fugacity coefficient with respect to pressure for the gas phase, [1/Pa].
dphi_dP_lDerivative of fugacity coefficient with respect to pressure for the liquid phase, [1/Pa].
dphi_dT_gDerivative of fugacity coefficient with respect to temperature for the gas phase, [1/K].
dphi_dT_lDerivative of fugacity coefficient with respect to temperature for the liquid phase, [1/K].
drho_dP_gDerivative of molar density with respect to pressure for the gas phase, [(mol/m^3)/Pa].
drho_dP_lDerivative of molar density with respect to pressure for the liquid phase, [(mol/m^3)/Pa].
drho_dT_gDerivative of molar density with respect to temperature for the gas phase, [(mol/m^3)/K].
drho_dT_lDerivative of molar density with respect to temperature for the liquid phase, [(mol/m^3)/K].
fugacity_gFugacity for the gas phase, [Pa].
fugacity_lFugacity for the liquid phase, [Pa].
isothermal_compressibility_gIsothermal (constant-temperature) expansion coefficient for the gas phase, [1/Pa].
isothermal_compressibility_lIsothermal (constant-temperature) expansion coefficient for the liquid phase, [1/Pa].
kappa_gIsothermal (constant-temperature) expansion coefficient for the gas phase, [1/Pa].
kappa_lIsothermal (constant-temperature) expansion coefficient for the liquid phase, [1/Pa].
lnphi_gThe natural logarithm of the fugacity coefficient for the gas phase, [-].
lnphi_lThe natural logarithm of the fugacity coefficient for the liquid phase, [-].
more_stable_phaseChecks the Gibbs energy of each possible phase, and returns ‘l’ if the liquid-like phase is more stable, and ‘g’ if the vapor-like phase is more stable.
mpmath_volume_ratiosMethod to compare, as ratios, the volumes of the implemented cubic solver versus those calculated using mpmath.
mpmath_volumesMethod to calculate to a high precision the exact roots to the cubic equation, using mpmath.
mpmath_volumes_floatMethod to calculate real roots of a cubic equation, using mpmath, but returned as floats.
phi_gFugacity coefficient for the gas phase, [Pa].
phi_lFugacity coefficient for the liquid phase, [Pa].
pseudo_PcApply a linear mole-fraction mixing rule to compute the average critical pressure, [Pa].
pseudo_TcApply a linear mole-fraction mixing rule to compute the average critical temperature, [K].
pseudo_aApply a linear mole-fraction mixing rule to compute the average a coefficient, [-].
pseudo_omegaApply a linear mole-fraction mixing rule to compute the average omega, [-].
rho_gGas molar density, [mol/m^3].
rho_lLiquid molar density, [mol/m^3].
sorted_volumesList of lexicographically-sorted molar volumes available from the root finding algorithm used to solve the PT point.
state_specsConvenience method to return the two specified state specs (T, P, or V) as a dictionary.
Methods
Hvap(T)Method to calculate enthalpy of vaporization for a pure fluid from an equation of state, without iteration.
PT_surface_special([Tmin, Tmax, Pmin, Pmax, ...])Method to create a plot of the special curves of a pure fluid - vapor pressure, determinant zeros, pseudo critical point, and mechanical critical point.
P_PIP_transition(T[, low_P_limit])Method to calculate the pressure which makes the phase identification parameter exactly 1.
P_discriminant_zero_g()Method to calculate the pressure which zero the discriminant function of the general cubic eos, and is likely to sit on a boundary between not having a vapor-like volume; and having a vapor-like volume.
P_discriminant_zero_l()Method to calculate the pressure which zero the discriminant function of the general cubic eos, and is likely to sit on a boundary between not having a liquid-like volume; and having a liquid-like volume.
P_discriminant_zeros()Method to calculate the pressures which zero the discriminant function of the general cubic eos, at the current temperature.
P_discriminant_zeros_analytical(T, b, delta, ...)Method to calculate the pressures which zero the discriminant function of the general cubic eos.
P_max_at_V(V)Dummy method.
Psat(T[, polish])Generic method to calculate vapor pressure of a pure-component equation of state for a specified T.
Psat_errors([Tmin, Tmax, pts, plot, show, ...])Method to create a plot of vapor pressure and the relative error of its calculation vs. the iterative polish approach.
T_discriminant_zero_g([T_guess])Method to calculate the temperature which zeros the discriminant function of the general cubic eos, and is likely to sit on a boundary between not having a vapor-like volume; and having a vapor-like volume.
T_discriminant_zero_l([T_guess])Method to calculate the temperature which zeros the discriminant function of the general cubic eos, and is likely to sit on a boundary between not having a liquid-like volume; and having a liquid-like volume.
T_max_at_V(V[, Pmax])Method to calculate the maximum temperature the EOS can create at a constant volume, if one exists; returns None otherwise.
T_min_at_V(V[, Pmin])Returns the minimum temperature for the EOS to have the volume as specified.
Tsat(P[, polish])Generic method to calculate the temperature for a specified vapor pressure of the pure fluid.
V_g_sat(T)Method to calculate molar volume of the vapor phase along the saturation line.
V_l_sat(T)Method to calculate molar volume of the liquid phase along the saturation line.
Vs_mpmath()Method to calculate real roots of a cubic equation, using mpmath.
a_alpha_and_derivatives(T[, full, quick, ...])Method to calculate a_alpha and its first and second derivatives for an EOS with the Van der Waals mixing rules.
a_alpha_and_derivatives_pure(T)Dummy method to calculate and its first and second derivatives.
a_alpha_for_Psat(T, Psat[, a_alpha_guess])Method to calculate which value of is required for a given T, Psat pair.
a_alpha_for_V(T, P, V)Method to calculate which value of is required for a given T, P pair to match a specified V.
a_alpha_plot([Tmin, Tmax, pts, plot, show])Method to create a plot of the parameter and its first two derivatives.
as_json([cache, option])Method to create a JSON-friendly serialization of the eos which can be stored, and reloaded later.
check_sufficient_inputs()Method to an exception if none of the pairs (T, P), (T, V), or (P, V) are given.
Calculates the molar departure Gibbs energy mole number derivatives (where the mole fractions sum to 1).
Calculates the molar departure Gibbs energy second composition derivative (where the mole fractions do not sum to 1).
d2Scomp_dninjs() diff/P_ref), n1, n2)d2V_dninjs(Z)Calculates the molar volume second mole number derivatives (where the mole fractions sum to 1).
d2V_dzizjs(Z)Calculates the molar volume second composition derivative (where the mole fractions do not sum to 1).
Calculates the mixture log fugacity coefficient second mole number derivatives (where the mole fraction sum to 1).
Calculates the mixture log fugacity coefficient second mole fraction derivatives (where the mole fractions do not sum to 1).
d2phi_sat_dT2(T[, polish])Method to calculate the second temperature derivative of saturation fugacity coefficient of the compound.
dA_dep_dns_Vt(phase)from sympy import * Vt, P, T, R, n1, n2, n3 = symbols('Vt, P, T, R, n1, n2, n3') # doctest:+SKIP P, V, a_alpha, delta, epsilon, b = symbols('P, V, a\alpha, delta, epsilon, b', cls=Function) # doctest:+SKIP da_alpha_dT, d2a_alpha_dT2 = symbols('da_alpha_dT, d2a_alpha_dT2', cls=Function) # doctest:+SKIP ns = [n1, n2, n3]
dG_dep_dns(Z)Calculates the molar departure Gibbs energy mole number derivatives (where the mole fractions sum to 1).
dG_dep_dzs(Z)Calculates the molar departure Gibbs energy composition derivative (where the mole fractions do not sum to 1).
dH_dep_dT_sat_g(T[, polish])Method to calculate and return the temperature derivative of saturation vapor excess enthalpy.
dH_dep_dT_sat_l(T[, polish])Method to calculate and return the temperature derivative of saturation liquid excess enthalpy.
dH_dep_dns(Z)Calculates the molar departure enthalpy mole number derivatives (where the mole fractions sum to 1).
dH_dep_dzs(Z)Calculates the molar departure enthalpy composition derivative (where the mole fractions do not sum to 1).
dP_dns_Vt(phase)from sympy import * Vt, P, T, R, n1, n2, n3, no = symbols('Vt, P, T, R, n1, n2, n3, no') # doctest:+SKIP n, P, V, a_alpha, delta, epsilon, b = symbols('n, P, V, a\alpha, delta, epsilon, b', cls=Function) # doctest:+SKIP da_alpha_dT, d2a_alpha_dT2 = symbols('da_alpha_dT, d2a_alpha_dT2', cls=Function) # doctest:+SKIP n = no + n1 + n2 + n3
dPsat_dT(T[, polish, also_Psat])Generic method to calculate the temperature derivative of vapor pressure for a specified T.
dS_dep_dT_sat_g(T[, polish])Method to calculate and return the temperature derivative of saturation vapor excess entropy.
dS_dep_dT_sat_l(T[, polish])Method to calculate and return the temperature derivative of saturation liquid excess entropy.
dS_dep_dns(Z)Calculates the molar departure entropy mole number derivatives (where the mole fractions sum to 1).
dS_dep_dzs(Z)Calculates the molar departure entropy composition derivative (where the mole fractions do not sum to 1).
dV_dns(Z)Calculates the molar volume mole number derivatives (where the mole fractions sum to 1).
dV_dzs(Z)Calculates the molar volume composition derivative (where the mole fractions do not sum to 1).
dZ_dns(Z)Calculates the compressibility mole number derivatives (where the mole fractions sum to 1).
dZ_dzs(Z)Calculates the compressibility composition derivatives (where the mole fractions do not sum to 1).
dfugacities_dns(phase)Generic formula for calculating the mole number derivaitves of fugacities for each species in a mixture.
discriminant([T, P])Method to compute the discriminant of the cubic volume solution with the current EOS parameters, optionally at the same (assumed) T, and P or at different ones, if values are specified.
dlnfugacities_dns(phase)Generic formula for calculating the mole number derivaitves of log fugacities for each species in a mixture.
dlnphi_dns(Z)Calculates the mixture log fugacity coefficient mole number derivatives (where the mole fractions sum to 1).
dlnphi_dzs(Z)Calculates the mixture log fugacity coefficient mole fraction derivatives (where the mole fractions do not sum to 1).
dlnphis_dP(phase)Generic formula for calculating the pressure derivaitve of log fugacity coefficients for each species in a mixture.
dlnphis_dT(phase)Generic formula for calculating the temperature derivaitve of log fugacity coefficients for each species in a mixture.
dlnphis_dns(Z)Generic formula for calculating the mole number derivaitves of log fugacity coefficients for each species in a mixture.
dlnphis_dzs(Z)Generic formula for calculating the mole fraction derivaitves of log fugacity coefficients for each species in a mixture.
dnG_dep_dns(Z)Calculates the partial molar departure Gibbs energy.
dnH_dep_dns(Z)Calculates the partial molar departure enthalpy.
dnV_dns(Z)Calculates the partial molar volume of the specified phase No specific formula is implemented for this property - it is calculated from the molar volume mole fraction derivative.
dnZ_dns(Z)Calculates the partial compressibility of the specified phase No specific formula is implemented for this property - it is calculated from the compressibility mole fraction derivative.
dphi_sat_dT(T[, polish])Method to calculate the temperature derivative of saturation fugacity coefficient of the compound.
from_json(json_repr[, cache])Method to create a mixture cubic equation of state from a JSON friendly serialization of another mixture cubic equation of state.
fugacities([only_l, only_g])Helper method for calculating fugacity coefficients for any phases present, using either the overall mole fractions for both phases or using specified mole fractions for each phase.
Generic formula for calculating log fugacity coefficients for each species in a mixture.
Method to calculate the mechanical critical point of a mixture of defined composition.
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.
phi_sat(T[, polish])Method to calculate the saturation fugacity coefficient of the compound.
pures()Helper method which returns a list of pure EOSs at the same T and P and base EOS as the mixture.
resolve_full_alphas()Generic method to resolve the eos with fully calculated alpha derviatives.
saturation_prop_plot(prop[, Tmin, Tmax, ...])Method to create a plot of a specified property of the EOS along the (pure component) saturation line.
set_dnzs_derivatives_and_departures([n, x, ...])Sets a number of mole number and/or composition partial derivatives of thermodynamic partial derivatives.
set_from_PT(Vs[, only_l, only_g])Counts the number of real volumes in Vs, and determines what to do.
set_properties_from_solution(T, P, V, b, ...)Sets all interesting properties which can be calculated from an EOS alone.
solve([pure_a_alphas, only_l, only_g, ...])First EOS-generic method; should be called by all specific EOSs.
solve_T(P, V[, quick, solution])Generic method to calculate T from a specified P and V.
solve_missing_volumes()Generic method to ensure both volumes, if solutions are physical, have calculated properties.
state_hash()Basic method to calculate a hash of the state of the model and its model parameters.
subset(idxs, **state_specs)Method to construct a new
GCEOSMIXthat removes all components not specified in the idxs argument.to([zs, T, P, V, fugacities])Method to construct a new
GCEOSMIXobject at two of T, P or V with the specified composition.to_PV(P, V)Method to construct a new
GCEOSMIXobject at the spcified P and V with the current composition.to_PV_zs(P, V, zs[, fugacities, only_l, only_g])Method to construct a new
GCEOSMIXinstance at P, V, and zs with the same parameters as the existing object.to_TP(T, P)Method to construct a new
GCEOSMIXobject at the spcified T and P with the current composition.to_TPV_pure(i[, T, P, V])Helper method which returns a pure EOSs at the specs (two of T, P and V) and base EOS as the mixture for a particular index.
to_TP_zs(T, P, zs[, fugacities, only_l, only_g])Method to construct a new
GCEOSMIXinstance at T, P, and zs with the same parameters as the existing object.to_TP_zs_fast(T, P, zs[, only_l, only_g, ...])Method to construct a new
GCEOSMIXinstance with the same parameters as the existing object.to_TV(T, V)Method to construct a new
GCEOSMIXobject at the spcified T and V with the current composition.Method to construct a new
GCEOSMIXobject at the current object's properties and composition, but which is at the mechanical critical point.volume_error()Method to calculate the relative absolute error in the calculated molar volumes.
volume_errors([Tmin, Tmax, Pmin, Pmax, pts, ...])Method to create a plot of the relative absolute error in the cubic volume solution as compared to a higher-precision calculation.
volume_solutions(T, P, b, delta, epsilon, ...)Halley's method based solver for cubic EOS volumes based on the idea of initializing from a single liquid-like guess which is solved precisely, deflating the cubic analytically, solving the quadratic equation for the next two volumes, and then performing two halley steps on each of them to obtain the final solutions.
volume_solutions_full(T, P, b, delta, ...[, ...])Newton-Raphson based solver for cubic EOS volumes based on the idea of initializing from an analytical solver.
volume_solutions_mp(T, P, b, delta, epsilon, ...)Solution of this form of the cubic EOS in terms of volumes, using the mpmath arbitrary precision library.
a_alpha_and_derivatives_numpy
a_alpha_and_derivatives_py
broyden2_VL
d2A_dep_dninjs
d2A_dep_dninjs_Vt
d2A_dninjs_Vt
d2A_dninjs_Vt_another
d2P_dninjs_Vt
d2nA_dninjs_Vt
d3P_dninjnks_Vt
dScomp_dns
derivatives_and_departures
main_derivatives_and_departures
newton_VL
sequential_substitution_VL
stability_Michelsen
stabiliy_iteration_Michelsen
- Psat(T, polish=False)[source]¶
Generic method to calculate vapor pressure of a pure-component equation of state for a specified T. An explicit solution is used unless polish is True.
The result of this function has no physical meaning for multicomponent mixtures, and does not represent either a dew point or a bubble point!
- Parameters:
- Returns:
- Psat
float Vapor pressure using the pure-component approach, [Pa]
- Psat
Notes
For multicomponent mixtures this may serve as a useful guess for the dew and the bubble pressure.
- a_alpha_and_derivatives(T, full=True, quick=True, pure_a_alphas=True)[source]¶
Method to calculate a_alpha and its first and second derivatives for an EOS with the Van der Waals mixing rules. Uses the parent class’s interface to compute pure component values. Returns a_alpha, da_alpha_dT, and d2a_alpha_dT2.
For use in
solve_Tthis returns only a_alpha if full is False.- Parameters:
- T
float Temperature, [K]
- fullbool,
optional If False, calculates and returns only a_alpha
- quickbool,
optional Only the quick variant is implemented; it is little faster anyhow
- pure_a_alphasbool,
optional Whether or not to recalculate the a_alpha terms of pure components (for the case of mixtures only) which stay the same as the composition changes (i.e in a PT flash), [-]
- T
- Returns:
Notes
The exact expressions can be obtained with the following SymPy expression below, commented out for brevity.
>>> from sympy import * >>> kij, T = symbols('kij, T ') >>> a_alpha_i, a_alpha_j = symbols('a_alpha_i, a_alpha_j', cls=Function) >>> a_alpha_ij = (1-kij)*sqrt(a_alpha_i(T)*a_alpha_j(T)) >>> diff(a_alpha_ij, T) >>> diff(a_alpha_ij, T, T)
- property a_alpha_ijs¶
Calculate and return the matrix .
- Returns:
Notes
In an earlier implementation this matrix was stored each EOS solve; however, allocating that much memory becomes quite expensive for large number of component cases and this is now calculated on-demand only.
- d2G_dep_dninjs(Z)[source]¶
Calculates the molar departure Gibbs energy mole number derivatives (where the mole fractions sum to 1). No specific formula is implemented for this property - it is calculated from the mole fraction derivative.
- d2G_dep_dzizjs(Z)[source]¶
Calculates the molar departure Gibbs energy second composition derivative (where the mole fractions do not sum to 1). Verified numerically. Useful in solving for gibbs minimization calculations or for solving for the true critical point. Also forms the basis for the molar departure Gibbs energy mole second number derivative.
- Parameters:
- Z
float Compressibility of the mixture for a desired phase, [-]
- Z
- Returns:
- d2G_dep_dzizjs
float Departure Gibbs free energy second composition derivatives, [J/mol]
- d2G_dep_dzizjs
Notes
The derivation for the derivative is performed as follows using SymPy. The function source code is an optimized variant created with the cse SymPy function, and hand optimized further.
>>> from sympy import * >>> P, T, R, x1, x2 = symbols('P, T, R, x1, x2') >>> a_alpha, delta, epsilon, V, b = symbols('a\ \\alpha, delta, epsilon, V, b', cls=Function) >>> da_alpha_dT, d2a_alpha_dT2 = symbols('da_alpha_dT, d2a_alpha_dT2', cls=Function) >>> S_dep = R*log(P*V(x1, x2)/(R*T)) + R*log(V(x1, x2)-b(x1, x2))+2*da_alpha_dT(x1, x2)*atanh((2*V(x1, x2)+delta(x1, x2))/sqrt(delta(x1, x2)**2-4*epsilon(x1, x2)))/sqrt(delta(x1, x2)**2-4*epsilon(x1, x2))-R*log(V(x1, x2)) >>> H_dep = P*V(x1, x2) - R*T + 2*atanh((2*V(x1, x2)+delta(x1, x2))/sqrt(delta(x1, x2)**2-4*epsilon(x1, x2)))*(da_alpha_dT(x1, x2)*T-a_alpha(x1, x2))/sqrt(delta(x1, x2)**2-4*epsilon(x1, x2)) >>> G_dep = simplify(H_dep - T*S_dep) >>> diff(G_dep, x1, x2)
- d2V_dninjs(Z)[source]¶
Calculates the molar volume second mole number derivatives (where the mole fractions sum to 1). No specific formula is implemented for this property - it is calculated from the second mole fraction derivatives.
- d2V_dzizjs(Z)[source]¶
Calculates the molar volume second composition derivative (where the mole fractions do not sum to 1). Verified numerically. Used in many other derivatives, and for the molar volume second mole number derivative.
- Parameters:
- Z
float Compressibility of the mixture for a desired phase, [-]
- Z
- Returns:
- d2V_dzizjs
float Molar volume second composition derivatives, [m^3/mol]
- d2V_dzizjs
Notes
The derivation for the derivative is performed as follows using SymPy. The function source code is an optimized variant created with the cse SymPy function, and hand optimized further.
>>> from sympy import * >>> P, T, R, x1, x2 = symbols('P, T, R, x1, x2') >>> V, delta, epsilon, a_alpha, b = symbols('V, delta, epsilon, a\ \\alpha, b', cls=Function) >>> CUBIC = R*T/(V(x1, x2) - b(x1, x2)) - a_alpha(x1, x2)/(V(x1, x2)*V(x1, x2) + delta(x1, x2)*V(x1, x2) + epsilon(x1, x2)) - P >>> solve(diff(CUBIC, x1, x2), Derivative(V(x1, x2), x1, x2))
- property d2a_alpha_dT2_dns¶
Helper method for calculating the mole number derivatives of d2a_alpha_dT2. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2a_alpha_dT2_dzs¶
Helper method for calculating the mole number derivatives of d2a_alpha_dT2. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2a_alpha_dT2_ijs¶
Calculate and return the matrix of the second temperature derivatives of the alpha terms.
- Returns:
Notes
In an earlier implementation this matrix was stored each EOS solve; however, allocating that much memory becomes quite expensive for large number of component cases and this is now calculated on-demand only.
- property d2a_alpha_dninjs¶
Helper method for calculating the second partial molar derivatives of a_alpha (hessian). Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2a_alpha_dzizjs¶
Helper method for calculating the second composition derivatives of a_alpha (hessian). Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2b_dninjs¶
Helper method for calculating the second partial mole number derivatives of b. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2b_dzizjs¶
Helper method for calculating the second partial mole fraction derivatives of b. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- d2lnphi_dninjs(Z)[source]¶
Calculates the mixture log fugacity coefficient second mole number derivatives (where the mole fraction sum to 1). No specific formula is implemented for this property - it is calculated from the second mole fraction derivative of Gibbs free energy.
- d2lnphi_dzizjs(Z)[source]¶
Calculates the mixture log fugacity coefficient second mole fraction derivatives (where the mole fractions do not sum to 1). No specific formula is implemented for this property - it is calculated from the second mole fraction derivative of Gibbs free energy.
- property d3a_alpha_dninjnks¶
Helper method for calculating the third mole number derivatives of a_alpha. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3a_alpha_dzizjzks¶
Helper method for calculating the third composition derivatives of a_alpha. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3b_dninjnks¶
Helper method for calculating the third partial mole number derivatives of b. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3b_dzizjzks¶
Helper method for calculating the third partial mole fraction derivatives of b. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3delta_dzizjzks¶
Helper method for calculating the third composition derivatives of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3epsilon_dzizjzks¶
Helper method for calculating the third composition derivatives of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- dG_dep_dns(Z)[source]¶
Calculates the molar departure Gibbs energy mole number derivatives (where the mole fractions sum to 1). No specific formula is implemented for this property - it is calculated from the mole fraction derivative.
Apart from the ideal term, this is the formulation for chemical potential.
- dG_dep_dzs(Z)[source]¶
Calculates the molar departure Gibbs energy composition derivative (where the mole fractions do not sum to 1). Verified numerically. Useful in solving for gibbs minimization calculations or for solving for the true critical point. Also forms the basis for the molar departure Gibbs energy mole number derivative and molar partial departure Gibbs energy.
- Parameters:
- Z
float Compressibility of the mixture for a desired phase, [-]
- Z
- Returns:
- dG_dep_dzs
float Departure Gibbs free energy composition derivatives, [J/mol]
- dG_dep_dzs
Notes
The derivation for the derivative is performed as follows using SymPy. The function source code is an optimized variant created with the cse SymPy function, and hand optimized further.
>>> from sympy import * >>> P, T, R, x = symbols('P, T, R, x') >>> a_alpha, a, delta, epsilon, V, b, da_alpha_dT = symbols('a\ \\alpha, a, delta, epsilon, V, b, da_alpha_dT', cls=Function) >>> S_dep = R*log(P*V(x)/(R*T)) + R*log(V(x)-b(x))+2*da_alpha_dT(x)*atanh((2*V(x)+delta(x))/sqrt(delta(x)**2-4*epsilon(x)))/sqrt(delta(x)**2-4*epsilon(x))-R*log(V(x)) >>> H_dep = P*V(x) - R*T + 2*atanh((2*V(x)+delta(x))/sqrt(delta(x)**2-4*epsilon(x)))*(da_alpha_dT(x)*T-a_alpha(x))/sqrt(delta(x)**2-4*epsilon(x)) >>> G_dep = simplify(H_dep - T*S_dep) >>> diff(G_dep, x) P*Derivative(V(x), x) - R*T*(Derivative(V(x), x) - Derivative(b(x), x))/(V(x) - b(x)) - 2*(-delta(x)*Derivative(delta(x), x) + 2*Derivative(epsilon(x), x))*a \alpha(x)*atanh(2*V(x)/sqrt(delta(x)**2 - 4*epsilon(x)) + delta(x)/sqrt(delta(x)**2 - 4*epsilon(x)))/(delta(x)**2 - 4*epsilon(x))**(3/2) - 2*atanh(2*V(x)/sqrt(delta(x)**2 - 4*epsilon(x)) + delta(x)/sqrt(delta(x)**2 - 4*epsilon(x)))*Derivative(a \alpha(x), x)/sqrt(delta(x)**2 - 4*epsilon(x)) - 2*(2*(-delta(x)*Derivative(delta(x), x) + 2*Derivative(epsilon(x), x))*V(x)/(delta(x)**2 - 4*epsilon(x))**(3/2) + (-delta(x)*Derivative(delta(x), x) + 2*Derivative(epsilon(x), x))*delta(x)/(delta(x)**2 - 4*epsilon(x))**(3/2) + 2*Derivative(V(x), x)/sqrt(delta(x)**2 - 4*epsilon(x)) + Derivative(delta(x), x)/sqrt(delta(x)**2 - 4*epsilon(x)))*a \alpha(x)/((1 - (2*V(x)/sqrt(delta(x)**2 - 4*epsilon(x)) + delta(x)/sqrt(delta(x)**2 - 4*epsilon(x)))**2)*sqrt(delta(x)**2 - 4*epsilon(x)))
- dH_dep_dns(Z)[source]¶
Calculates the molar departure enthalpy mole number derivatives (where the mole fractions sum to 1). No specific formula is implemented for this property - it is calculated from the mole fraction derivative.
- dH_dep_dzs(Z)[source]¶
Calculates the molar departure enthalpy composition derivative (where the mole fractions do not sum to 1). Verified numerically. Useful in solving for enthalpy specifications in newton-type methods, and forms the basis for the molar departure enthalpy mole number derivative and molar partial departure enthalpy.
- Parameters:
- Z
float Compressibility of the mixture for a desired phase, [-]
- Z
- Returns:
- dH_dep_dzs
float Departure enthalpy composition derivatives, [J/mol]
- dH_dep_dzs
Notes
The derivation for the derivative is performed as follows using SymPy. The function source code is an optimized variant created with the cse SymPy function, and hand optimized further.
>>> from sympy import * >>> P, T, V, R, b, a, delta, epsilon, x = symbols('P, T, V, R, b, a, delta, epsilon, x') >>> V, delta, epsilon, a_alpha, b = symbols('V, delta, epsilon, a_alpha, b', cls=Function) >>> H_dep = (P*V(x) - R*T + 2/sqrt(delta(x)**2 - 4*epsilon(x))*(T*Derivative(a_alpha(T, x), T) ... - a_alpha(x))*atanh((2*V(x)+delta(x))/sqrt(delta(x)**2-4*epsilon(x)))) >>> diff(H_dep, x) P*Derivative(V(x), x) + 2*(T*Derivative(a \alpha(T, x), T) - a \alpha(x))*(-delta(x)*Derivative(delta(x), x) + 2*Derivative(epsilon(x), x))*atanh((2*V(x) + delta(x))/sqrt(delta(x)**2 - 4*epsilon(x)))/(delta(x)**2 - 4*epsilon(x))**(3/2) + 2*(T*Derivative(a \alpha(T, x), T) - a \alpha(x))*((-delta(x)*Derivative(delta(x), x) + 2*Derivative(epsilon(x), x))*(2*V(x) + delta(x))/(delta(x)**2 - 4*epsilon(x))**(3/2) + (2*Derivative(V(x), x) + Derivative(delta(x), x))/sqrt(delta(x)**2 - 4*epsilon(x)))/((-(2*V(x) + delta(x))**2/(delta(x)**2 - 4*epsilon(x)) + 1)*sqrt(delta(x)**2 - 4*epsilon(x))) + 2*(T*Derivative(a \alpha(T, x), T, x) - Derivative(a \alpha(x), x))*atanh((2*V(x) + delta(x))/sqrt(delta(x)**2 - 4*epsilon(x)))/sqrt(delta(x)**2 - 4*epsilon(x))
- dS_dep_dns(Z)[source]¶
Calculates the molar departure entropy mole number derivatives (where the mole fractions sum to 1). No specific formula is implemented for this property - it is calculated from the mole fraction derivative.
- dS_dep_dzs(Z)[source]¶
Calculates the molar departure entropy composition derivative (where the mole fractions do not sum to 1). Verified numerically. Useful in solving for entropy specifications in newton-type methods, and forms the basis for the molar departure entropy mole number derivative and molar partial departure entropy.
- dV_dns(Z)[source]¶
Calculates the molar volume mole number derivatives (where the mole fractions sum to 1). No specific formula is implemented for this property - it is calculated from the mole fraction derivative.
- dV_dzs(Z)[source]¶
Calculates the molar volume composition derivative (where the mole fractions do not sum to 1). Verified numerically. Used in many other derivatives, and for the molar volume mole number derivative and partial molar volume calculation.
- Parameters:
- Z
float Compressibility of the mixture for a desired phase, [-]
- Z
- Returns:
- dV_dzs
float Molar volume composition derivatives, [m^3/mol]
- dV_dzs
Notes
The derivation for the derivative is performed as follows using SymPy. The function source code is an optimized variant created with the cse SymPy function, and hand optimized further.
>>> from sympy import * >>> P, T, R, x = symbols('P, T, R, x') >>> V, delta, epsilon, a_alpha, b = symbols('V, delta, epsilon, a\ \\alpha, b', cls=Function) >>> CUBIC = R*T/(V(x) - b(x)) - a_alpha(x)/(V(x)*V(x) + delta(x)*V(x) + epsilon(x)) - P >>> solve(diff(CUBIC, x), Derivative(V(x), x)) [(-R*T*(V(x)**2 + V(x)*delta(x) + epsilon(x))**3*Derivative(b(x), x) + (V(x) - b(x))**2*(V(x)**2 + V(x)*delta(x) + epsilon(x))**2*Derivative(a \alpha(x), x) - (V(x) - b(x))**2*V(x)**3*a \alpha(x)*Derivative(delta(x), x) - (V(x) - b(x))**2*V(x)**2*a \alpha(x)*delta(x)*Derivative(delta(x), x) - (V(x) - b(x))**2*V(x)**2*a \alpha(x)*Derivative(epsilon(x), x) - (V(x) - b(x))**2*V(x)*a \alpha(x)*delta(x)*Derivative(epsilon(x), x) - (V(x) - b(x))**2*V(x)*a \alpha(x)*epsilon(x)*Derivative(delta(x), x) - (V(x) - b(x))**2*a \alpha(x)*epsilon(x)*Derivative(epsilon(x), x))/(-R*T*(V(x)**2 + V(x)*delta(x) + epsilon(x))**3 + 2*(V(x) - b(x))**2*V(x)**3*a \alpha(x) + 3*(V(x) - b(x))**2*V(x)**2*a \alpha(x)*delta(x) + (V(x) - b(x))**2*V(x)*a \alpha(x)*delta(x)**2 + 2*(V(x) - b(x))**2*V(x)*a \alpha(x)*epsilon(x) + (V(x) - b(x))**2*a \alpha(x)*delta(x)*epsilon(x))]
- dZ_dns(Z)[source]¶
Calculates the compressibility mole number derivatives (where the mole fractions sum to 1). No specific formula is implemented for this property - it is calculated from the mole fraction derivative.
- dZ_dzs(Z)[source]¶
Calculates the compressibility composition derivatives (where the mole fractions do not sum to 1). No specific formula is implemented for this property - it is calculated from the composition derivative of molar volume, which does have its formula implemented.
- property da_alpha_dT_dns¶
Helper method for calculating the mole number derivatives of da_alpha_dT. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property da_alpha_dT_dzs¶
Helper method for calculating the composition derivatives of da_alpha_dT. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property da_alpha_dT_ijs¶
Calculate and return the matrix for the temperature derivatives of the alpha terms.
- Returns:
Notes
In an earlier implementation this matrix was stored each EOS solve; however, allocating that much memory becomes quite expensive for large number of component cases and this is now calculated on-demand only.
- property da_alpha_dns¶
Helper method for calculating the mole number derivatives of a_alpha. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property da_alpha_dzs¶
Helper method for calculating the composition derivatives of a_alpha. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property db_dns¶
Helper method for calculating the mole number derivatives of b. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- property db_dzs¶
Helper method for calculating the composition derivatives of b. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- dfugacities_dns(phase)[source]¶
Generic formula for calculating the mole number derivaitves of fugacities for each species in a mixture. Verified numerically. Applicable to all cubic equations of state which can be cast in the form used here.
- dlnfugacities_dns(phase)[source]¶
Generic formula for calculating the mole number derivaitves of log fugacities for each species in a mixture. Verified numerically. Applicable to all cubic equations of state which can be cast in the form used here.
- dlnphi_dns(Z)[source]¶
Calculates the mixture log fugacity coefficient mole number derivatives (where the mole fractions sum to 1). No specific formula is implemented for this property - it is calculated from the mole fraction derivative of Gibbs free energy.
This property can be converted into a partial molar property to obtain the individual fugacity coefficients.
- dlnphi_dzs(Z)[source]¶
Calculates the mixture log fugacity coefficient mole fraction derivatives (where the mole fractions do not sum to 1). No specific formula is implemented for this property - it is calculated from the mole fraction derivative of Gibbs free energy.
- dlnphis_dP(phase)[source]¶
Generic formula for calculating the pressure derivaitve of log fugacity coefficients for each species in a mixture. Verified numerically. Applicable to all cubic equations of state which can be cast in the form used here.
Normally this routine is slower than EOS-specific ones, as it does not make assumptions that certain parameters are zero or equal to other parameters.
- Parameters:
- phase
str One of ‘l’ or ‘g’, [-]
- phase
- Returns:
- dlnphis_dP
float Pressure derivatives of log fugacity coefficient for each species, [1/Pa]
- dlnphis_dP
Notes
This expression for the partial derivative of the mixture lnphi with respect to pressure and mole number can be derived as follows; to convert to the partial molar lnphi pressure and temperature derivative, add ::math::frac{G_{dep}/(RT)}{partial P}_{T, n}.
>>> from sympy import * >>> P, T, R, n = symbols('P, T, R, n') >>> a_alpha, a, delta, epsilon, V, b, da_alpha_dT, d2a_alpha_dT2 = symbols('a_alpha, a, delta, epsilon, V, b, da_alpha_dT, d2a_alpha_dT2', cls=Function) >>> S_dep = R*log(P*V(n, P)/(R*T)) + R*log(V(n, P)-b(n))+2*da_alpha_dT(n, T)*atanh((2*V(n, P)+delta(n))/sqrt(delta(n)**2-4*epsilon(n)))/sqrt(delta(n)**2-4*epsilon(n))-R*log(V(n, P)) >>> H_dep = P*V(n, P) - R*T + 2*atanh((2*V(n, P)+delta(n))/sqrt(delta(n)**2-4*epsilon(n)))*(da_alpha_dT(n, T)*T-a_alpha(n, T))/sqrt(delta(n)**2-4*epsilon(n)) >>> G_dep = H_dep - T*S_dep >>> lnphi = simplify(G_dep/(R*T)) >>> diff(diff(lnphi, P), n) P*Derivative(V(n, P), P, n)/(R*T) + Derivative(V(n, P), P, n)/V(n, P) - Derivative(V(n, P), P)*Derivative(V(n, P), n)/V(n, P)**2 - Derivative(V(n, P), P, n)/(V(n, P) - b(n)) - (-Derivative(V(n, P), n) + Derivative(b(n), n))*Derivative(V(n, P), P)/(V(n, P) - b(n))**2 + Derivative(V(n, P), n)/(R*T) - 4*(-2*delta(n)*Derivative(delta(n), n) + 4*Derivative(epsilon(n), n))*a_alpha(n, T)*Derivative(V(n, P), P)/(R*T*(1 - (2*V(n, P)/sqrt(delta(n)**2 - 4*epsilon(n)) + delta(n)/sqrt(delta(n)**2 - 4*epsilon(n)))**2)*(delta(n)**2 - 4*epsilon(n))**2) - 4*a_alpha(n, T)*Derivative(V(n, P), P, n)/(R*T*(1 - (2*V(n, P)/sqrt(delta(n)**2 - 4*epsilon(n)) + delta(n)/sqrt(delta(n)**2 - 4*epsilon(n)))**2)*(delta(n)**2 - 4*epsilon(n))) - 4*Derivative(V(n, P), P)*Derivative(a_alpha(n, T), n)/(R*T*(1 - (2*V(n, P)/sqrt(delta(n)**2 - 4*epsilon(n)) + delta(n)/sqrt(delta(n)**2 - 4*epsilon(n)))**2)*(delta(n)**2 - 4*epsilon(n))) - 4*(2*V(n, P)/sqrt(delta(n)**2 - 4*epsilon(n)) + delta(n)/sqrt(delta(n)**2 - 4*epsilon(n)))*(4*(-delta(n)*Derivative(delta(n), n) + 2*Derivative(epsilon(n), n))*V(n, P)/(delta(n)**2 - 4*epsilon(n))**(3/2) + 2*(-delta(n)*Derivative(delta(n), n) + 2*Derivative(epsilon(n), n))*delta(n)/(delta(n)**2 - 4*epsilon(n))**(3/2) + 4*Derivative(V(n, P), n)/sqrt(delta(n)**2 - 4*epsilon(n)) + 2*Derivative(delta(n), n)/sqrt(delta(n)**2 - 4*epsilon(n)))*a_alpha(n, T)*Derivative(V(n, P), P)/(R*T*(1 - (2*V(n, P)/sqrt(delta(n)**2 - 4*epsilon(n)) + delta(n)/sqrt(delta(n)**2 - 4*epsilon(n)))**2)**2*(delta(n)**2 - 4*epsilon(n))) + R*T*(P*Derivative(V(n, P), P)/(R*T) + V(n, P)/(R*T))*Derivative(V(n, P), n)/(P*V(n, P)**2) - R*T*(P*Derivative(V(n, P), P, n)/(R*T) + Derivative(V(n, P), n)/(R*T))/(P*V(n, P))
- dlnphis_dT(phase)[source]¶
Generic formula for calculating the temperature derivaitve of log fugacity coefficients for each species in a mixture. Verified numerically. Applicable to all cubic equations of state which can be cast in the form used here.
Normally this routine is slower than EOS-specific ones, as it does not make assumptions that certain parameters are zero or equal to other parameters.
- Parameters:
- phase
str One of ‘l’ or ‘g’, [-]
- phase
- Returns:
- dlnphis_dT
float Temperature derivatives of log fugacity coefficient for each species, [1/K]
- dlnphis_dT
Notes
This expression for the partial derivative of the mixture lnphi with respect to pressure and mole number can be derived as follows; to convert to the partial molar lnphi pressure and temperature derivative, add ::math::frac{G_{dep}/(RT)}{partial T}_{P, n}.
>>> from sympy import * >>> P, T, R, n = symbols('P, T, R, n') >>> a_alpha, a, delta, epsilon, V, b, da_alpha_dT, d2a_alpha_dT2 = symbols('a_alpha, a, delta, epsilon, V, b, da_alpha_dT, d2a_alpha_dT2', cls=Function) >>> S_dep = R*log(P*V(n, T)/(R*T)) + R*log(V(n, T)-b(n))+2*da_alpha_dT(n, T)*atanh((2*V(n, T)+delta(n))/sqrt(delta(n)**2-4*epsilon(n)))/sqrt(delta(n)**2-4*epsilon(n))-R*log(V(n, T)) >>> H_dep = P*V(n, T) - R*T + 2*atanh((2*V(n, T)+delta(n))/sqrt(delta(n)**2-4*epsilon(n)))*(da_alpha_dT(n, T)*T-a_alpha(n, T))/sqrt(delta(n)**2-4*epsilon(n)) >>> G_dep = H_dep - T*S_dep >>> lnphi = simplify(G_dep/(R*T)) >>> diff(diff(lnphi, T), n)
- dlnphis_dns(Z)[source]¶
Generic formula for calculating the mole number derivaitves of log fugacity coefficients for each species in a mixture. Verified numerically. Applicable to all cubic equations of state which can be cast in the form used here.
- dlnphis_dzs(Z)[source]¶
Generic formula for calculating the mole fraction derivaitves of log fugacity coefficients for each species in a mixture. Verified numerically. Applicable to all cubic equations of state which can be cast in the form used here.
- dnG_dep_dns(Z)[source]¶
Calculates the partial molar departure Gibbs energy. No specific formula is implemented for this property - it is calculated from the mole fraction derivative.
- dnH_dep_dns(Z)[source]¶
Calculates the partial molar departure enthalpy. No specific formula is implemented for this property - it is calculated from the mole fraction derivative.
- dnV_dns(Z)[source]¶
Calculates the partial molar volume of the specified phase No specific formula is implemented for this property - it is calculated from the molar volume mole fraction derivative.
- dnZ_dns(Z)[source]¶
Calculates the partial compressibility of the specified phase No specific formula is implemented for this property - it is calculated from the compressibility mole fraction derivative.
- property dna_alpha_dT_dns¶
Helper method for calculating the mole number derivatives of da_alpha_dT. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property dna_alpha_dns¶
Helper method for calculating the partial molar derivatives of a_alpha. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property dnb_dns¶
Helper method for calculating the partial molar derivative of b. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- classmethod from_json(json_repr, cache=None)[source]¶
Method to create a mixture cubic equation of state from a JSON friendly serialization of another mixture cubic equation of state.
- Parameters:
- json_repr
dict Json representation, [-]
- json_repr
- Returns:
- eos_mix
GCEOSMIX Newly created object from the json serialization, [-]
- eos_mix
Notes
It is important that the input string be in the same format as that created by
GCEOS.as_json.Examples
>>> import pickle >>> eos = PRSV2MIX(Tcs=[507.6], Pcs=[3025000], omegas=[0.2975], zs=[1], T=299., P=1E6, kappa1s=[0.05104], kappa2s=[0.8634], kappa3s=[0.460]) >>> json_stuff = pickle.dumps(eos.as_json()) >>> new_eos = PRSV2MIX.from_json(pickle.loads(json_stuff)) >>> assert new_eos == eos
- fugacities(only_l=False, only_g=False)[source]¶
Helper method for calculating fugacity coefficients for any phases present, using either the overall mole fractions for both phases or using specified mole fractions for each phase.
Requires
fugacity_coefficientsto be implemented by each subclassing EOS.In addition to setting fugacities_l and/or fugacities_g, this also sets the fugacity coefficients phis_l and/or phis_g.
Note that in a flash calculation, each phase requires their own EOS object.
- Parameters:
Notes
It is helpful to check that
fugacity_coefficientshas been implemented correctly using the following expression, from [1].For reference, several expressions for fugacity of a component are as follows, shown in [1] and [2].
References
[1] (1,2)Hu, Jiawen, Rong Wang, and Shide Mao. “Some Useful Expressions for Deriving Component Fugacity Coefficients from Mixture Fugacity Coefficient.” Fluid Phase Equilibria 268, no. 1-2 (June 25, 2008): 7-13. doi:10.1016/j.fluid.2008.03.007.
[2]Walas, Stanley M. Phase Equilibria in Chemical Engineering. Butterworth-Heinemann, 1985.
- fugacity_coefficients(Z)[source]¶
Generic formula for calculating log fugacity coefficients for each species in a mixture. Verified numerically. Applicable to all cubic equations of state which can be cast in the form used here. Normally this routine is slower than EOS-specific ones, as it does not make assumptions that certain parameters are zero or equal to other parameters.
- kwargs_linear = ()¶
Tuple of 1D arguments used by the specific EOS in addition to the conventional ones.
- kwargs_square = ('kijs',)¶
Tuple of 2D arguments used by the specific EOS.
- mechanical_critical_point()[source]¶
Method to calculate the mechanical critical point of a mixture of defined composition.
The mechanical critical point is where:
Notes
One useful application of the mechanical critical temperature is that the phase identification approach of Venkatarathnam is valid only up to it.
Note that the equation of state, when solved at these conditions, will have fairly large (1e-3 - 1e-6) results for the derivatives; but they are the minimum. This is just from floating point precision.
It can also be checked looking at the calculated molar volumes - all three (available with
sorted_volumes) will be very close (1e-5 difference in practice), again differing because of floating point error.The algorithm here is a custom implementation, using Newton-Raphson’s method with the initial guesses described in [1] (mole-weighted critical pressure average, critical temperature average using a quadratic mixing rule). Normally ~4 iterations are needed to solve the system. It is relatively fast, as only one evaluation of a_alpha and da_alpha_dT are needed per call to function and its jacobian.
References
[1]Watson, Harry A. J., and Paul I. Barton. “Reliable Flash Calculations: Part 3. A Nonsmooth Approach to Density Extrapolation and Pseudoproperty Evaluation.” Industrial & Engineering Chemistry Research, November 11, 2017. https://doi.org/10.1021/acs.iecr.7b03233.
[2]Mathias P. M., Boston J. F., and Watanasiri S. “Effective Utilization of Equations of State for Thermodynamic Properties in Process Simulation.” AIChE Journal 30, no. 2 (June 17, 2004): 182-86. https://doi.org/10.1002/aic.690300203.
- mix_kwargs_to_pure = {}¶
- multicomponent = True¶
All inherited classes of GCEOSMIX are multicomponent.
- nonstate_constants = ('N', 'cmps', 'Tcs', 'Pcs', 'omegas', 'kijs', 'kwargs', 'ais', 'bs')¶
- property pseudo_Pc¶
Apply a linear mole-fraction mixing rule to compute the average critical pressure, [Pa].
Examples
>>> base = RKMIX(T=150.0, P=4e6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> base.pseudo_Pc 3878000.0
- property pseudo_Tc¶
Apply a linear mole-fraction mixing rule to compute the average critical temperature, [K].
Examples
>>> base = RKMIX(T=150.0, P=4e6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> base.pseudo_Tc 151.9
- property pseudo_a¶
Apply a linear mole-fraction mixing rule to compute the average a coefficient, [-].
Examples
>>> base = RKMIX(T=150.0, P=4e6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> base.pseudo_a 0.17634464184
- property pseudo_omega¶
Apply a linear mole-fraction mixing rule to compute the average omega, [-].
Examples
>>> base = RKMIX(T=150.0, P=4e6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> base.pseudo_omega 0.0284
- pures()[source]¶
Helper method which returns a list of pure EOSs at the same T and P and base EOS as the mixture.
- Returns:
- eos_pures
list[eos] A list of pure-species EOSs at the same T and P as the system, [-]
- eos_pures
Notes
This is useful for i.e. comparing mixture fugacities with the Lewis-Randall rule or when using an activity coefficient model which require pure component fugacities.
- set_dnzs_derivatives_and_departures(n=True, x=True, only_l=False, only_g=False)[source]¶
Sets a number of mole number and/or composition partial derivatives of thermodynamic partial derivatives.
The list of properties set is as follows, with all properties suffixed with ‘_l’ or ‘_g’
if n is True: d2P_dTdns, d2P_dVdns, d2V_dTdns, d2V_dPdns, d2T_dVdns, d2T_dPdns, d3P_dT2dns, d3P_dV2dns, d3V_dT2dns, d3V_dP2dns, d3T_dV2dns, d3T_dP2dns, d3V_dPdTdns, d3P_dTdVdns, d3T_dPdVdns, dV_dep_dns, dG_dep_dns, dH_dep_dns, dU_dep_dns, dS_dep_dns, dA_dep_dns
if x is True: d2P_dTdzs, d2P_dVdzs, d2V_dTdzs, d2V_dPdzs, d2T_dVdzs, d2T_dPdzs, d3P_dT2dzs, d3P_dV2dzs, d3V_dT2dzs, d3V_dP2dzs, d3T_dV2dzs, d3T_dP2dzs, d3V_dPdTdzs, d3P_dTdVdzs, d3T_dPdVdzs, dV_dep_dzs, dG_dep_dzs, dH_dep_dzs, dU_dep_dzs, dS_dep_dzs, dA_dep_dzs
- Parameters:
- nbool,
optional Whether or not to set the mole number derivatives (sums up to one), [-]
- xbool,
optional Whether or not to set the composition derivatives (does not sum up to one), [-]
- only_lbool,
optional Whether or not to set only the liquid-like phase properties (if there are two phases), [-]
- only_gbool,
optional Whether or not to set only the gas-like phase properties (if there are two phases), [-]
- nbool,
- solve_T(P, V, quick=True, solution=None)[source]¶
Generic method to calculate T from a specified P and V. Provides SciPy’s newton solver, and iterates to solve the general equation for P, recalculating a_alpha as a function of temperature using
a_alpha_and_derivativeseach iteration.- Parameters:
- P
float Pressure, [Pa]
- V
float Molar volume, [m^3/mol]
- quickbool,
optional Unimplemented, although it may be possible to derive explicit expressions as done for many pure-component EOS
- solution
strorNone,optional ‘l’ or ‘g’ to specify a liquid of vapor solution (if one exists); if None, will select a solution more likely to be real (closer to STP, attempting to avoid temperatures like 60000 K or 0.0001 K).
- P
- Returns:
- T
float Temperature, [K]
- T
- subset(idxs, **state_specs)[source]¶
Method to construct a new
GCEOSMIXthat removes all components not specified in the idxs argument.- Parameters:
- Returns:
Notes
Subclassing equations of state require their
kwargs_linearandkwargs_squareattributes to be correct for this to work. Tcs, Pcs, and omegas are always assumed to be used.Examples
>>> kijs = [[0.0, 0.00076, 0.00171], [0.00076, 0.0, 0.00061], [0.00171, 0.00061, 0.0]] >>> PR3 = PRMIX(Tcs=[469.7, 507.4, 540.3], zs=[0.8168, 0.1501, 0.0331], omegas=[0.249, 0.305, 0.349], Pcs=[3.369E6, 3.012E6, 2.736E6], T=322.29, P=101325.0, kijs=kijs) >>> PR3.subset([1,2]) PRMIX(Tcs=[507.4, 540.3], Pcs=[3012000.0, 2736000.0], omegas=[0.305, 0.349], kijs=[[0.0, 0.00061], [0.00061, 0.0]], zs=[0.8193231441048036, 0.1806768558951965], T=322.29, P=101325.0) >>> PR3.subset([1,2], T=500.0, P=1e5, zs=[.2, .8]) PRMIX(Tcs=[507.4, 540.3], Pcs=[3012000.0, 2736000.0], omegas=[0.305, 0.349], kijs=[[0.0, 0.00061], [0.00061, 0.0]], zs=[0.2, 0.8], T=500.0, P=100000.0) >>> PR3.subset([1,2], zs=[.2, .8]) PRMIX(Tcs=[507.4, 540.3], Pcs=[3012000.0, 2736000.0], omegas=[0.305, 0.349], kijs=[[0.0, 0.00061], [0.00061, 0.0]], zs=[0.2, 0.8], T=322.29, P=101325.0)
- to(zs=None, T=None, P=None, V=None, fugacities=True)[source]¶
Method to construct a new
GCEOSMIXobject at two of T, P or V with the specified composition. In the event the specs match those of the current object, it will be returned unchanged.- Parameters:
- Returns:
Notes
Constructs the object with parameters Tcs, Pcs, omegas, and kwargs.
Examples
>>> base = PRMIX(T=500.0, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> base.to(T=300.0, P=1e9).state_specs {'T': 300.0, 'P': 1000000000.0} >>> base.to(T=300.0, V=1.0).state_specs {'T': 300.0, 'V': 1.0} >>> base.to(P=1e5, V=1.0).state_specs {'P': 100000.0, 'V': 1.0}
- to_PV(P, V)[source]¶
Method to construct a new
GCEOSMIXobject at the spcified P and V with the current composition. In the event the P and V match the current object’s P and V, it will be returned unchanged.- Parameters:
- Returns:
Notes
Constructs the object with parameters Tcs, Pcs, omegas, and kwargs.
Examples
>>> base = RKMIX(T=500.0, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> new = base.to_PV(P=1000000.0, V=1.0) >>> base.state_specs, new.state_specs ({'T': 500.0, 'P': 1000000.0}, {'P': 1000000.0, 'V': 1.0})
- to_PV_zs(P, V, zs, fugacities=True, only_l=False, only_g=False)[source]¶
Method to construct a new
GCEOSMIXinstance at P, V, and zs with the same parameters as the existing object. Optionally, only one set of phase properties can be solved for, increasing speed. The fugacities calculation can be be skipped by by setting fugacities to False.- Parameters:
- P
float Pressure, [Pa]
- V
float Molar volume, [m^3/mol]
- zs
list[float] Mole fractions of each component, [-]
- fugacitiesbool
Whether or not to calculate and set the fugacities of each component, [-]
- only_lbool
When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set.
- only_gbool
When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set.
- P
- Returns:
Notes
A check for whether or not P, V, and zs are the same as the existing instance is performed; if it is, the existing object is returned.
Examples
>>> base = RKMIX(T=500.0, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> base.to_PV_zs(V=0.004162, P=1e5, zs=[.1, 0.9]) RKMIX(Tcs=[126.1, 190.6], Pcs=[3394000.0, 4604000.0], omegas=[0.04, 0.011], kijs=[[0.0, 0.0], [0.0, 0.0]], zs=[0.1, 0.9], P=100000.0, V=0.004162)
- to_TP(T, P)[source]¶
Method to construct a new
GCEOSMIXobject at the spcified T and P with the current composition. In the event the T and P match the current object’s T and P, it will be returned unchanged.- Parameters:
- Returns:
Notes
Constructs the object with parameters Tcs, Pcs, omegas, and kwargs.
Examples
>>> base = RKMIX(T=500.0, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> new = base.to_TP(T=10.0, P=2000.0) >>> base.state_specs, new.state_specs ({'T': 500.0, 'P': 1000000.0}, {'T': 10.0, 'P': 2000.0})
- to_TPV_pure(i, T=None, P=None, V=None)[source]¶
Helper method which returns a pure EOSs at the specs (two of T, P and V) and base EOS as the mixture for a particular index.
- to_TP_zs(T, P, zs, fugacities=True, only_l=False, only_g=False)[source]¶
Method to construct a new
GCEOSMIXinstance at T, P, and zs with the same parameters as the existing object. Optionally, only one set of phase properties can be solved for, increasing speed. The fugacities calculation can be be skipped by by setting fugacities to False.- Parameters:
- T
float Temperature, [K]
- P
float Pressure, [Pa]
- zs
list[float] Mole fractions of each component, [-]
- fugacitiesbool
Whether or not to calculate and set the fugacities of each component, [-]
- only_lbool
When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set.
- only_gbool
When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set.
- T
- Returns:
Notes
A check for whether or not T, P, and zs are the same as the existing instance is performed; if it is, the existing object is returned.
Examples
>>> base = RKMIX(T=500.0, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> base.to_TP_zs(T=300.0, P=1e5, zs=[.1, 0.9]) RKMIX(Tcs=[126.1, 190.6], Pcs=[3394000.0, 4604000.0], omegas=[0.04, 0.011], kijs=[[0.0, 0.0], [0.0, 0.0]], zs=[0.1, 0.9], T=300.0, P=100000.0)
- to_TP_zs_fast(T, P, zs, only_l=False, only_g=False, full_alphas=True)[source]¶
Method to construct a new
GCEOSMIXinstance with the same parameters as the existing object. If both instances are at the same temperature, a_alphas and da_alpha_dTs and d2a_alpha_dT2s are shared between the instances. It is always assumed the new object has a differet composition. Optionally, only one set of phase properties can be solved for, increasing speed. Additionally, if full_alphas is set to False no temperature derivatives of a_alpha will be computed. Those derivatives are not needed in the context of a PT or PVF flash.- Parameters:
- T
float Temperature, [K]
- P
float Pressure, [Pa]
- zs
list[float] Mole fractions of each component, [-]
- only_lbool
When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set.
- only_gbool
When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set.
- T
- Returns:
Examples
>>> base = RKMIX(T=500.0, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> base.to_TP_zs_fast(T=300.0, P=1e5, zs=base.zs) RKMIX(Tcs=[126.1, 190.6], Pcs=[3394000.0, 4604000.0], omegas=[0.04, 0.011], kijs=[[0.0, 0.0], [0.0, 0.0]], zs=[0.6, 0.4], T=300.0, P=100000.0)
- to_TV(T, V)[source]¶
Method to construct a new
GCEOSMIXobject at the spcified T and V with the current composition. In the event the T and V match the current object’s T and V, it will be returned unchanged.- Parameters:
- Returns:
Notes
Constructs the object with parameters Tcs, Pcs, omegas, and kwargs.
Examples
>>> base = RKMIX(T=500.0, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> new = base.to_TV(T=1000000.0, V=1.0) >>> base.state_specs, new.state_specs ({'T': 500.0, 'P': 1000000.0}, {'T': 1000000.0, 'V': 1.0})
- to_mechanical_critical_point()[source]¶
Method to construct a new
GCEOSMIXobject at the current object’s properties and composition, but which is at the mechanical critical point.Examples
>>> base = RKMIX(T=500.0, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.6, 0.4]) >>> base.to_mechanical_critical_point() RKMIX(Tcs=[126.1, 190.6], Pcs=[3394000.0, 4604000.0], omegas=[0.04, 0.011], kijs=[[0.0, 0.0], [0.0, 0.0]], zs=[0.6, 0.4], T=151.861, P=3908737.9)
- translated = False¶
Whether or not the model implements volume translation.
- vectorized = False¶
Whether the model is implemented using numpy arrays of float64, or pure-Python lists of floats.
Peng-Robinson Family EOSs¶
Standard Peng Robinson¶
- class thermo.eos_mix.PRMIX(Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
-
Class for solving the Peng-Robinson [1] [2] cubic equation of state for a mixture of any number of compounds. Subclasses PR. Solves the EOS on initialization and calculates fugacities for all components in all phases.
Two of T, P, and V are needed to solve the EOS.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
- Attributes:
d2delta_dninjsHelper method for calculating the second mole number derivatives (hessian) of delta.
d2delta_dzizjsHelper method for calculating the second composition derivatives (hessian) of delta.
d2epsilon_dninjsHelper method for calculating the second mole number derivatives (hessian) of epsilon.
d2epsilon_dzizjsHelper method for calculating the second composition derivatives (hessian) of epsilon.
d3a_alpha_dT3Method to calculate approximately the third temperature derivative of a_alpha for the PR EOS.
d3delta_dninjnksHelper method for calculating the third partial mole number derivatives of delta.
d3epsilon_dninjnksHelper method for calculating the third partial mole number derivatives of epsilon.
ddelta_dnsHelper method for calculating the mole number derivatives of delta.
ddelta_dzsHelper method for calculating the composition derivatives of delta.
depsilon_dnsHelper method for calculating the mole number derivatives of epsilon.
depsilon_dzsHelper method for calculating the composition derivatives of epsilon.
Methods
Method to calculate the pure-component a_alphas and their first and second derivatives for the PR EOS.
Method to calculate the pure-component a_alphas for the PR EOS.
dlnphis_dP(phase)Generic formula for calculating the pressure derivaitve of log fugacity coefficients for each species in a mixture for the Peng-Robinson EOS.
dlnphis_dT(phase)Formula for calculating the temperature derivaitve of log fugacity coefficients for each species in a mixture for the Peng-Robinson equation of state.
dlnphis_dzs(Z)Calculate and return the mole fraction derivaitves of log fugacity coefficients for each species in a mixture.
Literature formula for calculating fugacity coefficients for each species in a mixture.
Notes
For P-V initializations, a numerical solver is used to find T.
References
[1]Peng, Ding-Yu, and Donald B. Robinson. “A New Two-Constant Equation of State.” Industrial & Engineering Chemistry Fundamentals 15, no. 1 (February 1, 1976): 59-64. doi:10.1021/i160057a011.
[2]Robinson, Donald B., Ding-Yu Peng, and Samuel Y-K Chung. “The Development of the Peng - Robinson Equation and Its Application to Phase Equilibrium in a System Containing Methanol.” Fluid Phase Equilibria 24, no. 1 (January 1, 1985): 25-41. doi:10.1016/0378-3812(85)87035-7.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = PRMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.V_l, eos.V_g (3.6257362939e-05, 0.00070066592313) >>> eos.fugacities_l, eos.fugacities_g ([793860.83821, 73468.552253], [436530.92470, 358114.63827])
- a_alpha_and_derivatives_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas and their first and second derivatives for the PR EOS. This vectorized implementation is added for extra speed.
- Parameters:
- T
float Temperature, [K]
- T
- Returns:
- a_alphas
list[float] Coefficient calculated by EOS-specific method, [J^2/mol^2/Pa]
- da_alpha_dTs
list[float] Temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K]
- d2a_alpha_dT2s
list[float] Second temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K**2]
- a_alphas
- a_alphas_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas for the PR EOS. This vectorized implementation is added for extra speed.
- property d2delta_dninjs¶
Helper method for calculating the second mole number derivatives (hessian) of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2delta_dzizjs¶
Helper method for calculating the second composition derivatives (hessian) of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2epsilon_dninjs¶
Helper method for calculating the second mole number derivatives (hessian) of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2epsilon_dzizjs¶
Helper method for calculating the second composition derivatives (hessian) of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3a_alpha_dT3¶
Method to calculate approximately the third temperature derivative of a_alpha for the PR EOS. A rigorous calculation has not been implemented.
- property d3delta_dninjnks¶
Helper method for calculating the third partial mole number derivatives of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3epsilon_dninjnks¶
Helper method for calculating the third partial mole number derivatives of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property ddelta_dns¶
Helper method for calculating the mole number derivatives of delta. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- property ddelta_dzs¶
Helper method for calculating the composition derivatives of delta. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- property depsilon_dns¶
Helper method for calculating the mole number derivatives of epsilon. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- property depsilon_dzs¶
Helper method for calculating the composition derivatives of epsilon. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- dlnphis_dP(phase)[source]¶
Generic formula for calculating the pressure derivaitve of log fugacity coefficients for each species in a mixture for the Peng-Robinson EOS. Verified numerically.
- Parameters:
- phase
str One of ‘l’ or ‘g’, [-]
- phase
- Returns:
- dlnphis_dP
float Pressure derivatives of log fugacity coefficient for each species, [1/Pa]
- dlnphis_dP
Notes
This expression was derived using SymPy and optimized with the cse technique.
- dlnphis_dT(phase)[source]¶
Formula for calculating the temperature derivaitve of log fugacity coefficients for each species in a mixture for the Peng-Robinson equation of state. Verified numerically.
- Parameters:
- phase
str One of ‘l’ or ‘g’, [-]
- phase
- Returns:
- dlnphis_dT
float Temperature derivatives of log fugacity coefficient for each species, [1/K]
- dlnphis_dT
Notes
This expression was derived using SymPy and optimized with the cse technique.
- dlnphis_dzs(Z)[source]¶
Calculate and return the mole fraction derivaitves of log fugacity coefficients for each species in a mixture. This formula is specific to the Peng-Robinson equation of state.
- Parameters:
- Z
float Compressibility of the mixture for a desired phase, [-]
- Z
- Returns:
Notes
This formula is from [1] but is validated to match the generic implementation.
References
[1]Chang, Yih-Bor. “Development and Application of an Equation of State Compositional Simulator” 1990. https://repositories.lib.utexas.edu/handle/2152/80585.
Examples
>>> kijs = [[0, 0.00076, 0.00171], [0.00076, 0, 0.00061], [0.00171, 0.00061, 0]] >>> eos = PRMIX(Tcs=[469.7, 507.4, 540.3], zs=[0.8168, 0.1501, 0.0331], omegas=[0.249, 0.305, 0.349], Pcs=[3.369E6, 3.012E6, 2.736E6], T=322.29, P=101325, kijs=kijs) >>> eos.dlnphis_dzs(eos.Z_l) [[0.009938069276, 0.0151503498382, 0.018297235797], [-0.038517738793, -0.05958926042, -0.068438990795], [-0.07057106923, -0.10363920720, -0.14116283024]]
- fugacity_coefficients(Z)[source]¶
Literature formula for calculating fugacity coefficients for each species in a mixture. Verified numerically. Applicable to most derivatives of the Peng-Robinson equation of state as well. Called by
fugacitieson initialization, or by a solver routine which is performing a flash calculation.
Peng Robinson (1978)¶
- class thermo.eos_mix.PR78MIX(Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
PRMIXClass for solving the Peng-Robinson cubic equation of state for a mixture of any number of compounds according to the 1978 variant. Subclasses PR. Solves the EOS on initialization and calculates fugacities for all components in all phases.
Two of T, P, and V are needed to solve the EOS.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Notes
This variant is recommended over the original.
References
[1]Peng, Ding-Yu, and Donald B. Robinson. “A New Two-Constant Equation of State.” Industrial & Engineering Chemistry Fundamentals 15, no. 1 (February 1, 1976): 59-64. doi:10.1021/i160057a011.
[2]Robinson, Donald B., Ding-Yu Peng, and Samuel Y-K Chung. “The Development of the Peng - Robinson Equation and Its Application to Phase Equilibrium in a System Containing Methanol.” Fluid Phase Equilibria 24, no. 1 (January 1, 1985): 25-41. doi:10.1016/0378-3812(85)87035-7.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa, with modified acentric factors to show the difference between
PRMIX>>> eos = PR78MIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.6, 0.7], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.V_l, eos.V_g (3.2396438915e-05, 0.00050433802024) >>> eos.fugacities_l, eos.fugacities_g ([833048.45119, 6160.9088153], [460717.27767, 279598.90103])
Peng Robinson Stryjek-Vera¶
- class thermo.eos_mix.PRSVMIX(Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None, kappa1s=None, fugacities=True, only_l=False, only_g=False)[source]¶
-
Class for solving the Peng-Robinson-Stryjek-Vera equations of state for a mixture as given in [1]. Subclasses
PRMIXandPRSV. Solves the EOS on initialization and calculates fugacities for all components in all phases.Inherits the method of calculating fugacity coefficients from
PRMIX. Two of T, P, and V are needed to solve the EOS.- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- kappa1s
list[float],optional Fit parameter; available in [1] for over 90 compounds, SRKMIXTranslated[-]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Method to calculate the pure-component a_alphas and their first and second derivatives for the PRSV EOS.
Method to calculate the pure-component a_alphas for the PRSV EOS.
Notes
[1] recommends that kappa1 be set to 0 for Tr > 0.7. This is not done by default; the class boolean kappa1_Tr_limit may be set to True and the problem re-solved with that specified if desired. kappa1_Tr_limit is not supported for P-V inputs.
For P-V initializations, a numerical solver is used to find T.
[2] and [3] are two more resources documenting the PRSV EOS. [4] lists kappa values for 69 additional compounds. See also
PRSV2MIX. Note that tabulated kappa values should be used with the critical parameters used in their fits. Both [1] and [4] only considered vapor pressure in fitting the parameter.References
[1] (1,2,3,4)Stryjek, R., and J. H. Vera. “PRSV: An Improved Peng-Robinson Equation of State for Pure Compounds and Mixtures.” The Canadian Journal of Chemical Engineering 64, no. 2 (April 1, 1986): 323-33. doi:10.1002/cjce.5450640224.
[2]Stryjek, R., and J. H. Vera. “PRSV - An Improved Peng-Robinson Equation of State with New Mixing Rules for Strongly Nonideal Mixtures.” The Canadian Journal of Chemical Engineering 64, no. 2 (April 1, 1986): 334-40. doi:10.1002/cjce.5450640225.
[3]Stryjek, R., and J. H. Vera. “Vapor-liquid Equilibrium of Hydrochloric Acid Solutions with the PRSV Equation of State.” Fluid Phase Equilibria 25, no. 3 (January 1, 1986): 279-90. doi:10.1016/0378-3812(86)80004-8.
Examples
P-T initialization, two-phase, nitrogen and methane
>>> eos = PRSVMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.phase, eos.V_l, eos.H_dep_l, eos.S_dep_l ('l/g', 3.6235536165e-05, -6349.0055583, -49.1240502472)
- a_alpha_and_derivatives_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas and their first and second derivatives for the PRSV EOS. This vectorized implementation is added for extra speed.
- Parameters:
- T
float Temperature, [K]
- T
- Returns:
- a_alphas
list[float] Coefficient calculated by EOS-specific method, [J^2/mol^2/Pa]
- da_alpha_dTs
list[float] Temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K]
- d2a_alpha_dT2s
list[float] Second temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K**2]
- a_alphas
Peng Robinson Stryjek-Vera 2¶
- class thermo.eos_mix.PRSV2MIX(Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None, kappa1s=None, kappa2s=None, kappa3s=None, fugacities=True, only_l=False, only_g=False)[source]¶
-
Class for solving the Peng-Robinson-Stryjek-Vera 2 equations of state for a Mixture as given in [1]. Subclasses
PRMIXand PRSV2 <thermo.eos.PRSV2>. Solves the EOS on initialization and calculates fugacities for all components in all phases.Inherits the method of calculating fugacity coefficients from
PRMIX. Two of T, P, and V are needed to solve the EOS.- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- kappa1s
list[float],optional Fit parameter; available in [1] for over 90 compounds, [-]
- kappa2s
list[float],optional Fit parameter; available in [1] for over 90 compounds, [-]
- kappa3s
list[float],optional Fit parameter; available in [1] for over 90 compounds, [-]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Method to calculate the pure-component a_alphas and their first and second derivatives for the PRSV2 EOS.
Method to calculate the pure-component a_alphas for the PRSV2 EOS.
Notes
For P-V initializations, a numerical solver is used to find T.
Note that tabulated kappa values should be used with the critical parameters used in their fits. [1] considered only vapor pressure in fitting the parameter.
References
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = PRSV2MIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.V_l, eos.V_g (3.6235536165e-05, 0.00070024238654) >>> eos.fugacities_l, eos.fugacities_g ([794057.58318, 72851.22327], [436553.65618, 357878.11066])
- a_alpha_and_derivatives_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas and their first and second derivatives for the PRSV2 EOS. This vectorized implementation is added for extra speed.
- Parameters:
- T
float Temperature, [K]
- T
- Returns:
- a_alphas
list[float] Coefficient calculated by EOS-specific method, [J^2/mol^2/Pa]
- da_alpha_dTs
list[float] Temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K]
- d2a_alpha_dT2s
list[float] Second temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K**2]
- a_alphas
- a_alphas_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas for the PRSV2 EOS. This vectorized implementation is added for extra speed.
- Parameters:
- T
float Temperature, [K]
- T
- Returns:
Examples
>>> eos = PRSV2MIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.a_alphas_vectorized(300) [0.0860568595, 0.20174345803]
Peng Robinson Twu (1995)¶
- class thermo.eos_mix.TWUPRMIX(Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
TwuPR95_a_alpha,PRMIXClass for solving the Twu [1] variant of the Peng-Robinson cubic equation of state for a mixture. Solves the EOS on initialization and calculates fugacities for all components in all phases.
Two of T, P, and V are needed to solve the EOS.
For sub-critical conditions:
L0, M0, N0 = 0.125283, 0.911807, 1.948150;
L1, M1, N1 = 0.511614, 0.784054, 2.812520
For supercritical conditions:
L0, M0, N0 = 0.401219, 4.963070, -0.2;
L1, M1, N1 = 0.024955, 1.248089, -8.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Notes
For P-V initializations, a numerical solver is used to find T. Claimed to be more accurate than the PR, PR78 and PRSV equations.
References
[1]Twu, Chorng H., John E. Coon, and John R. Cunningham. “A New Generalized Alpha Function for a Cubic Equation of State Part 1. Peng-Robinson Equation.” Fluid Phase Equilibria 105, no. 1 (March 15, 1995): 49-59. doi:10.1016/0378-3812(94)02601-V.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = TWUPRMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.V_l, eos.V_g (3.624571041e-05, 0.0007004401318) >>> eos.fugacities_l, eos.fugacities_g ([792155.022163, 73305.88829], [436468.967764, 358049.2495573])
Peng Robinson Translated¶
- class thermo.eos_mix.PRMIXTranslated(Tcs, Pcs, omegas, zs, kijs=None, cs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
PRMIXClass for solving the Peng-Robinson [1] [2] translated cubic equation of state for a mixture of any number of compounds. Solves the EOS on initialization and calculates fugacities for all components in all phases.
Two of T, P, and V are needed to solve the EOS.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- cs
list[float],optional Volume translation parameters; always zero in the original implementation, [m^3/mol]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
- Attributes:
d2delta_dninjsHelper method for calculating the second mole number derivatives (hessian) of delta.
d2delta_dzizjsHelper method for calculating the second composition derivatives (hessian) of delta.
d2epsilon_dninjsHelper method for calculating the second mole number derivatives (hessian) of epsilon.
d2epsilon_dzizjsHelper method for calculating the second composition derivatives (hessian) of epsilon.
d3delta_dninjnksHelper method for calculating the third partial mole number derivatives of delta.
d3delta_dzizjzksHelper method for calculating the third composition derivatives of delta.
d3epsilon_dninjnksHelper method for calculating the third partial mole number derivatives of epsilon.
d3epsilon_dzizjzksHelper method for calculating the third composition derivatives of epsilon.
ddelta_dnsHelper method for calculating the mole number derivatives of delta.
ddelta_dzsHelper method for calculating the composition derivatives of delta.
depsilon_dnsHelper method for calculating the mole number derivatives of epsilon.
depsilon_dzsHelper method for calculating the composition derivatives of epsilon.
Methods
Notes
For P-V initializations, a numerical solver is used to find T.
References
[1]Peng, Ding-Yu, and Donald B. Robinson. “A New Two-Constant Equation of State.” Industrial & Engineering Chemistry Fundamentals 15, no. 1 (February 1, 1976): 59-64. doi:10.1021/i160057a011.
[2]Robinson, Donald B., Ding-Yu Peng, and Samuel Y-K Chung. “The Development of the Peng - Robinson Equation and Its Application to Phase Equilibrium in a System Containing Methanol.” Fluid Phase Equilibria 24, no. 1 (January 1, 1985): 25-41. doi:10.1016/0378-3812(85)87035-7.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = PRMIXTranslated(T=115, P=1E6, cs=[-4.4e-6, -4.35e-6], Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.2, 0.8], kijs=[[0,0.03],[0.03,0]]) >>> eos.V_l, eos.V_g (3.9079056337e-05, 0.00060231393016) >>> eos.fugacities_l, eos.fugacities_g ([442838.8615, 108854.48589], [184396.972, 565531.7709])
- property d2delta_dninjs¶
Helper method for calculating the second mole number derivatives (hessian) of delta. Note this is independent of the phase. refers to the original b parameter not involving any translation.
- Returns:
Notes
This derivative is checked numerically.
- property d2delta_dzizjs¶
Helper method for calculating the second composition derivatives (hessian) of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2epsilon_dninjs¶
Helper method for calculating the second mole number derivatives (hessian) of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2epsilon_dzizjs¶
Helper method for calculating the second composition derivatives (hessian) of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3delta_dninjnks¶
Helper method for calculating the third partial mole number derivatives of delta. Note this is independent of the phase. refers to the original b parameter not involving any translation.
- Returns:
Notes
This derivative is checked numerically.
- property d3delta_dzizjzks¶
Helper method for calculating the third composition derivatives of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3epsilon_dninjnks¶
Helper method for calculating the third partial mole number derivatives of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3epsilon_dzizjzks¶
Helper method for calculating the third composition derivatives of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property ddelta_dns¶
Helper method for calculating the mole number derivatives of delta. Note this is independent of the phase. refers to the original b parameter not involving any translation.
Notes
This derivative is checked numerically.
- property ddelta_dzs¶
Helper method for calculating the composition derivatives of delta. Note this is independent of the phase. refers to the original b parameter not involving any translation.
Notes
This derivative is checked numerically.
- property depsilon_dns¶
Helper method for calculating the mole number derivatives of epsilon. Note this is independent of the phase. refers to the original b parameter not involving any translation.
Notes
This derivative is checked numerically.
- property depsilon_dzs¶
Helper method for calculating the composition derivatives of epsilon. Note this is independent of the phase. refers to the original b parameter not involving any translation.
Notes
This derivative is checked numerically.
- eos_pure¶
alias of
PRTranslated
Peng Robinson Translated-Consistent¶
- class thermo.eos_mix.PRMIXTranslatedConsistent(Tcs, Pcs, omegas, zs, kijs=None, cs=None, alpha_coeffs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
Twu91_a_alpha,PRMIXTranslatedClass for solving the volume translated Le Guennec, Privat, and Jaubert revision of the Peng-Robinson equation of state according to [1].
Two of T, P, and V are needed to solve the EOS.
If c is not provided, they are estimated as:
If alpha_coeffs is not provided, the parameters L and M are estimated from the acentric factor as follows:
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- cs
list[float],optional Volume translation parameters, [m^3/mol]
- alpha_coeffs
list[tuple(float[3])],optional Coefficients L, M, N (also called C1, C2, C3) of TWU 1991 form, [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Notes
For P-V initializations, a numerical solver is used to find T.
References
[1]Le Guennec, Yohann, Romain Privat, and Jean-Noël Jaubert. “Development of the Translated-Consistent Tc-PR and Tc-RK Cubic Equations of State for a Safe and Accurate Prediction of Volumetric, Energetic and Saturation Properties of Pure Compounds in the Sub- and Super-Critical Domains.” Fluid Phase Equilibria 429 (December 15, 2016): 301-12. https://doi.org/10.1016/j.fluid.2016.09.003.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = PRMIXTranslatedConsistent(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.2, 0.8], kijs=[[0,0.03],[0.03,0]]) >>> eos.V_l, eos.V_g (3.675235812e-05, 0.00059709319879) >>> eos.fugacities_l, eos.fugacities_g ([443454.9336, 106184.004057], [184122.74082, 563037.785])
- eos_pure¶
alias of
PRTranslatedConsistent
Peng Robinson Translated (Pina-Martinez, Privat, and Jaubert Variant)¶
- class thermo.eos_mix.PRMIXTranslatedPPJP(Tcs, Pcs, omegas, zs, kijs=None, cs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
PRMIXTranslatedClass for solving the Pina-Martinez, Privat, Jaubert, and Peng revision of the Peng-Robinson equation of state.
Two of T, P, and V are needed to solve the EOS.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- cs
list[float],optional Volume translation parameters, [m^3/mol]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Notes
For P-V initializations, a numerical solver is used to find T.
References
[1]Pina-Martinez, Andrés, Romain Privat, Jean-Noël Jaubert, and Ding-Yu Peng. “Updated Versions of the Generalized Soave a-Function Suitable for the Redlich-Kwong and Peng-Robinson Equations of State.” Fluid Phase Equilibria, December 7, 2018. https://doi.org/10.1016/j.fluid.2018.12.007.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = PRMIXTranslatedPPJP(T=115, P=1E6, cs=[-4.4e-6, -4.35e-6], Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.2, 0.8], kijs=[[0,0.03],[0.03,0]]) >>> eos.V_l, eos.V_g (3.8989032701e-05, 0.00059686183724) >>> eos.fugacities_l, eos.fugacities_g ([444791.13707, 104520.280997], [184782.600238, 563352.147])
- eos_pure¶
alias of
PRTranslatedPPJP
SRK Family EOSs¶
Standard SRK¶
- class thermo.eos_mix.SRKMIX(Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
EpsilonZeroMixingRules,GCEOSMIX,SRKClass for solving the Soave-Redlich-Kwong cubic equation of state for a mixture of any number of compounds. Solves the EOS on initialization and calculates fugacities for all components in all phases.
The implemented method here is
fugacity_coefficients, which implements the formula for fugacity coefficients in a mixture as given in [1]. Two of T, P, and V are needed to solve the EOS.- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Method to calculate the pure-component a_alphas and their first and second derivatives for the SRK EOS.
Method to calculate the pure-component a_alphas for the SRK EOS.
dlnphis_dP(phase)Generic formula for calculating the pressure derivaitve of log fugacity coefficients for each species in a mixture for the SRK EOS.
dlnphis_dT(phase)Formula for calculating the temperature derivaitve of log fugacity coefficients for each species in a mixture for the SRK equation of state.
Literature formula for calculating fugacity coefficients for each species in a mixture.
Notes
For P-V initializations, a numerical solver is used to find T.
References
[1]Soave, Giorgio. “Equilibrium Constants from a Modified Redlich-Kwong Equation of State.” Chemical Engineering Science 27, no. 6 (June 1972): 1197-1203. doi:10.1016/0009-2509(72)80096-4.
[2]Poling, Bruce E. The Properties of Gases and Liquids. 5th edition. New York: McGraw-Hill Professional, 2000.
[3]Walas, Stanley M. Phase Equilibria in Chemical Engineering. Butterworth-Heinemann, 1985.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> SRK_mix = SRKMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> SRK_mix.V_l, SRK_mix.V_g (4.1047569614e-05, 0.0007110158049)
- a_alpha_and_derivatives_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas and their first and second derivatives for the SRK EOS. This vectorized implementation is added for extra speed.
- Parameters:
- T
float Temperature, [K]
- T
- Returns:
- a_alphas
list[float] Coefficient calculated by EOS-specific method, [J^2/mol^2/Pa]
- da_alpha_dTs
list[float] Temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K]
- d2a_alpha_dT2s
list[float] Second temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K**2]
- a_alphas
- a_alphas_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas for the SRK EOS. This vectorized implementation is added for extra speed.
- dlnphis_dP(phase)[source]¶
Generic formula for calculating the pressure derivaitve of log fugacity coefficients for each species in a mixture for the SRK EOS. Verified numerically.
- Parameters:
- phase
str One of ‘l’ or ‘g’, [-]
- phase
- Returns:
- dlnphis_dP
float Pressure derivatives of log fugacity coefficient for each species, [1/Pa]
- dlnphis_dP
Notes
This expression was derived using SymPy and optimized with the cse technique.
- dlnphis_dT(phase)[source]¶
Formula for calculating the temperature derivaitve of log fugacity coefficients for each species in a mixture for the SRK equation of state. Verified numerically.
- Parameters:
- phase
str One of ‘l’ or ‘g’, [-]
- phase
- Returns:
- dlnphis_dT
float Temperature derivatives of log fugacity coefficient for each species, [1/K]
- dlnphis_dT
Notes
This expression was derived using SymPy and optimized with the cse technique.
- fugacity_coefficients(Z)[source]¶
Literature formula for calculating fugacity coefficients for each species in a mixture. Verified numerically. Applicable to most derivatives of the SRK equation of state as well. Called by
fugacitieson initialization, or by a solver routine which is performing a flash calculation.
Twu SRK (1995)¶
- class thermo.eos_mix.TWUSRKMIX(Tcs, Pcs, omegas, zs, kijs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
TwuSRK95_a_alpha,SRKMIXClass for solving the Twu variant of the Soave-Redlich-Kwong cubic equation of state for a mixture. Solves the EOS on initialization and calculates fugacities for all components in all phases.
Two of T, P, and V are needed to solve the EOS.
For sub-critical conditions:
L0, M0, N0 = 0.141599, 0.919422, 2.496441
L1, M1, N1 = 0.500315, 0.799457, 3.291790
For supercritical conditions:
L0, M0, N0 = 0.441411, 6.500018, -0.20
L1, M1, N1 = 0.032580, 1.289098, -8.0
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Notes
For P-V initializations, a numerical solver is used to find T. Claimed to be more accurate than the SRK equation.
References
[1]Twu, Chorng H., John E. Coon, and John R. Cunningham. “A New Generalized Alpha Function for a Cubic Equation of State Part 2. Redlich-Kwong Equation.” Fluid Phase Equilibria 105, no. 1 (March 15, 1995): 61-69. doi:10.1016/0378-3812(94)02602-W.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = TWUSRKMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.V_l, eos.V_g (4.1087927542e-05, 0.00071170732525) >>> eos.fugacities_l, eos.fugacities_g ([809692.830826, 74093.6388157], [441783.431489, 362470.3174107])
API SRK¶
- class thermo.eos_mix.APISRKMIX(Tcs, Pcs, zs, omegas=None, kijs=None, T=None, P=None, V=None, S1s=None, S2s=None, fugacities=True, only_l=False, only_g=False)[source]¶
-
Class for solving the Refinery Soave-Redlich-Kwong cubic equation of state for a mixture of any number of compounds, as shown in the API Databook [1]. Subclasses
APISRK. Solves the EOS on initialization and calculates fugacities for all components in all phases.Two of T, P, and V are needed to solve the EOS.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- S1s
float,optional Fit constant or estimated from acentric factor if not provided [-]
- S2s
float,optional Fit constant or 0 if not provided [-]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Notes
For P-V initializations, a numerical solver is used to find T.
References
[1]API Technical Data Book: General Properties & Characterization. American Petroleum Institute, 7E, 2005.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = APISRKMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.V_l, eos.V_g (4.101592310e-05, 0.00071046883030) >>> eos.fugacities_l, eos.fugacities_g ([817882.3033, 71620.4823812], [442158.29113, 361519.79877])
SRK Translated¶
- class thermo.eos_mix.SRKMIXTranslated(Tcs, Pcs, omegas, zs, kijs=None, cs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
SRKMIXClass for solving the volume translated Soave-Redlich-Kwong cubic equation of state for a mixture of any number of compounds. Subclasses
SRKMIX. Solves the EOS on initialization and calculates fugacities for all components in all phases.Two of T, P, and V are needed to solve the EOS.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- cs
list[float],optional Volume translation parameters; always zero in the original implementation, [m^3/mol]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
- Attributes:
d2delta_dninjsHelper method for calculating the second mole number derivatives (hessian) of delta.
d2delta_dzizjsHelper method for calculating the second composition derivatives (hessian) of delta.
d2epsilon_dninjsHelper method for calculating the second mole number derivatives (hessian) of epsilon.
d2epsilon_dzizjsHelper method for calculating the second composition derivatives (hessian) of epsilon.
d3delta_dninjnksHelper method for calculating the third partial mole number derivatives of delta.
d3delta_dzizjzksHelper method for calculating the third composition derivatives of delta.
d3epsilon_dninjnksHelper method for calculating the third partial mole number derivatives of epsilon.
d3epsilon_dzizjzksHelper method for calculating the third composition derivatives of epsilon.
ddelta_dnsHelper method for calculating the mole number derivatives of delta.
ddelta_dzsHelper method for calculating the composition derivatives of delta.
depsilon_dnsHelper method for calculating the mole number derivatives of epsilon.
depsilon_dzsHelper method for calculating the composition derivatives of epsilon.
Methods
Notes
For P-V initializations, a numerical solver is used to find T.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = SRKMIXTranslated(T=115, P=1E6, cs=[-4.4e-6, -4.35e-6], Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.2, 0.8], kijs=[[0,0.03],[0.03,0]]) >>> eos.V_l, eos.V_g (4.35928920e-05, 0.00060927202)
- property d2delta_dninjs¶
Helper method for calculating the second mole number derivatives (hessian) of delta. Note this is independent of the phase. refers to the original b parameter not involving any translation.
- Returns:
Notes
This derivative is checked numerically.
- property d2delta_dzizjs¶
Helper method for calculating the second composition derivatives (hessian) of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2epsilon_dninjs¶
Helper method for calculating the second mole number derivatives (hessian) of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2epsilon_dzizjs¶
Helper method for calculating the second composition derivatives (hessian) of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3delta_dninjnks¶
Helper method for calculating the third partial mole number derivatives of delta. Note this is independent of the phase. refers to the original b parameter not involving any translation.
- Returns:
Notes
This derivative is checked numerically.
- property d3delta_dzizjzks¶
Helper method for calculating the third composition derivatives of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3epsilon_dninjnks¶
Helper method for calculating the third partial mole number derivatives of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3epsilon_dzizjzks¶
Helper method for calculating the third composition derivatives of epsilon. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property ddelta_dns¶
Helper method for calculating the mole number derivatives of delta. Note this is independent of the phase. refers to the original b parameter not involving any translation.
Notes
This derivative is checked numerically.
- property ddelta_dzs¶
Helper method for calculating the composition derivatives of delta. Note this is independent of the phase. refers to the original b parameter not involving any translation.
Notes
This derivative is checked numerically.
- property depsilon_dns¶
Helper method for calculating the mole number derivatives of epsilon. Note this is independent of the phase. refers to the original b parameter not involving any translation.
Notes
This derivative is checked numerically.
- property depsilon_dzs¶
Helper method for calculating the composition derivatives of epsilon. Note this is independent of the phase. refers to the original b parameter not involving any translation.
Notes
This derivative is checked numerically.
- eos_pure¶
alias of
SRKTranslated
SRK Translated-Consistent¶
- class thermo.eos_mix.SRKMIXTranslatedConsistent(Tcs, Pcs, omegas, zs, kijs=None, cs=None, alpha_coeffs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
Twu91_a_alpha,SRKMIXTranslatedClass for solving the volume translated Le Guennec, Privat, and Jaubert revision of the SRK equation of state according to [1].
Two of T, P, and V are needed to solve the EOS.
If cs is not provided, they are estimated as:
If alpha_coeffs is not provided, the parameters L and M are estimated from each of the acentric factors as follows:
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- cs
list[float],optional Volume translation parameters, [m^3/mol]
- alpha_coeffs
list[list[float]] Coefficients for
thermo.eos_alpha_functions.Twu91_a_alpha, [-]- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Notes
For P-V initializations, a numerical solver is used to find T.
References
[1]Le Guennec, Yohann, Romain Privat, and Jean-Noël Jaubert. “Development of the Translated-Consistent Tc-PR and Tc-RK Cubic Equations of State for a Safe and Accurate Prediction of Volumetric, Energetic and Saturation Properties of Pure Compounds in the Sub- and Super-Critical Domains.” Fluid Phase Equilibria 429 (December 15, 2016): 301-12. https://doi.org/10.1016/j.fluid.2016.09.003.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = SRKMIXTranslatedConsistent(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.2, 0.8], kijs=[[0,0.03],[0.03,0]]) >>> eos.V_l, eos.V_g (3.591044498e-05, 0.0006020501621)
- eos_pure¶
alias of
SRKTranslatedConsistent
MSRK Translated¶
- class thermo.eos_mix.MSRKMIXTranslated(Tcs, Pcs, omegas, zs, kijs=None, cs=None, alpha_coeffs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
Soave_1979_a_alpha,SRKMIXTranslatedConsistentClass for solving the volume translated Soave (1980) alpha function, revision of the Soave-Redlich-Kwong equation of state for a pure compound according to [1]. Uses two fitting parameters N and M to more accurately fit the vapor pressure of pure species.
Two of T, P, and V are needed to solve the EOS.
This is an older correlation that offers lower accuracy on many properties which were sacrificed to obtain the vapor pressure accuracy. The alpha function of this EOS does not meet any of the consistency requriements for alpha functions.
Coefficients can be found in [2], or estimated with the method in [3]. The estimation method in [3] works as follows, using the acentric factor and true critical compressibility:
An alternate estimation scheme is provided in [1], which provides analytical solutions to calculate the parameters M and N from two points on the vapor pressure curve, suggested as 10 mmHg and 1 atm. This is used as an estimation method here if the parameters are not provided, and the two vapor pressure points are obtained from the original SRK equation of state.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- cs
list[float],optional Volume translation parameters, [m^3/mol]
- alpha_coeffs
list[list[float]] Coefficients for
thermo.eos_alpha_functions.Soave_1979_a_alpha, [-]- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
Notes
For P-V initializations, a numerical solver is used to find T.
References
[1] (1,2)Soave, G. “Rigorous and Simplified Procedures for Determining the Pure-Component Parameters in the Redlich—Kwong—Soave Equation of State.” Chemical Engineering Science 35, no. 8 (January 1, 1980): 1725-30. https://doi.org/10.1016/0009-2509(80)85007-X.
[2]Sandarusi, Jamal A., Arthur J. Kidnay, and Victor F. Yesavage. “Compilation of Parameters for a Polar Fluid Soave-Redlich-Kwong Equation of State.” Industrial & Engineering Chemistry Process Design and Development 25, no. 4 (October 1, 1986): 957-63. https://doi.org/10.1021/i200035a020.
[3] (1,2)Valderrama, Jose O., Héctor De la Puente, and Ahmed A. Ibrahim. “Generalization of a Polar-Fluid Soave-Redlich-Kwong Equation of State.” Fluid Phase Equilibria 93 (February 11, 1994): 377-83. https://doi.org/10.1016/0378-3812(94)87021-7.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = MSRKMIXTranslated(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.2, 0.8], kijs=[[0,0.03],[0.03,0]]) >>> eos.V_l, eos.V_g (3.9222990198e-05, 0.00060438075638)
- eos_pure¶
alias of
MSRKTranslated
Cubic Equation of State with Activity Coefficients¶
- class thermo.eos_mix.PSRK(Tcs, Pcs, omegas, zs, alpha_coeffs, ge_model, kijs=None, cs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
Mathias_Copeman_poly_a_alpha,PSRKMixingRules,SRKMIXTranslatedClass for solving the Predictive Soave-Redlich-Kwong [1] equation of state for a mixture of any number of compounds. Solves the EOS on initialization.
Two of T, P, and V are needed to solve the EOS.
Warning
This class is not complete! Fugacities and their derivatives among others are not yet implemented.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- omegas
list[float] Acentric factors of all compounds, [-]
- zs
list[float] Overall mole fractions of all species, [-]
- alpha_coeffs
list[list[float]] Coefficients for
thermo.eos_alpha_functions.Mathias_Copeman_poly_a_alpha, [-]- ge_model
thermo.activity.GibbsExcessobject Excess Gibbs free energy model; to match the PSRK model, this is a
thermo.unifac.UNIFACobject, [-]- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- cs
list[float],optional Volume translation parameters; always zero in the original implementation, [m^3/mol]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
Methods
References
[1]Holderbaum, T., and J. Gmehling. “PSRK: A Group Contribution Equation of State Based on UNIFAC.” Fluid Phase Equilibria 70, no. 2-3 (December 30, 1991): 251-65. https://doi.org/10.1016/0378-3812(91)85038-V.
Examples
T-P initialization, equimolar CO2, n-hexane:
>>> from thermo.unifac import UNIFAC, PSRKIP, PSRKSG >>> Tcs = [304.2, 507.4] >>> Pcs = [7.37646e6, 3.014419e6] >>> omegas = [0.2252, 0.2975] >>> zs = [0.5, 0.5] >>> Mathias_Copeman_coeffs = [[-1.7039, 0.2515, 0.8252, 1.0], [2.9173, -1.4411, 1.1061, 1.0]] >>> T = 313. >>> P = 1E6 >>> ge_model = UNIFAC.from_subgroups(T=T, xs=zs, chemgroups=[{117: 1}, {1:2, 2:4}], subgroups=PSRKSG, interaction_data=PSRKIP, version=0) >>> eos = PSRK(Tcs=Tcs, Pcs=Pcs, omegas=omegas, zs=zs, ge_model=ge_model, alpha_coeffs=Mathias_Copeman_coeffs, T=T, P=P) >>> eos PSRK(Tcs=[304.2, 507.4], Pcs=[7376460.0, 3014419.0], omegas=[0.2252, 0.2975], kijs=[[0.0, 0.0], [0.0, 0.0]], alpha_coeffs=[[-1.7039, 0.2515, 0.8252, 1.0], [2.9173, -1.4411, 1.1061, 1.0]], cs=[0.0, 0.0], ge_model=UNIFAC(T=313.0, xs=[0.5, 0.5], rs=[1.3, 4.4998000000000005], qs=[0.982, 3.856], Qs=[0.848, 0.54, 0.982], vs=[[0.0, 2.0], [0.0, 4.0], [1.0, 0.0]], psi_a=[[0.0, 0.0, 919.8], [0.0, 0.0, 919.8], [-38.672, -38.672, 0.0]], psi_b=[[0.0, 0.0, -3.9132], [0.0, 0.0, -3.9132], [0.8615, 0.8615, 0.0]], psi_c=[[0.0, 0.0, 0.0046309], [0.0, 0.0, 0.0046309], [-0.0017906, -0.0017906, 0.0]], version=0), zs=[0.5, 0.5], T=313.0, P=1000000.0) >>> eos.phase, eos.V_l, eos.V_g ('l/g', 0.000110889753959, 0.00197520225546)
- eos_pure¶
alias of
SRKTranslated
Van der Waals Equation of State¶
- class thermo.eos_mix.VDWMIX(Tcs, Pcs, zs, kijs=None, T=None, P=None, V=None, omegas=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
EpsilonZeroMixingRules,GCEOSMIX,VDWClass for solving the Van der Waals [1] [2] cubic equation of state for a mixture of any number of compounds. Solves the EOS on initialization and calculates fugacities for all components in all phases.
Two of T, P, and V are needed to solve the EOS.
- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- omegas
list[float],optional Acentric factors of all compounds - Not used in equation of state!, [-]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
- Attributes:
d2delta_dninjsHelper method for calculating the second mole number derivatives (hessian) of delta.
d2delta_dzizjsHelper method for calculating the second composition derivatives (hessian) of delta.
d3delta_dninjnksHelper method for calculating the third partial mole number derivatives of delta.
ddelta_dnsHelper method for calculating the mole number derivatives of delta.
ddelta_dzsHelper method for calculating the composition derivatives of delta.
Methods
Method to calculate the pure-component a_alphas and their first and second derivatives for the VDW EOS.
Method to calculate the pure-component a_alphas for the VDW EOS.
dlnphis_dP(phase)Generic formula for calculating the pressure derivaitve of log fugacity coefficients for each species in a mixture for the VDW EOS.
dlnphis_dT(phase)Formula for calculating the temperature derivaitve of log fugacity coefficients for each species in a mixture for the VDW equation of state.
Literature formula for calculating fugacity coefficients for each species in a mixture.
Notes
For P-V initializations, a numerical solver is used to find T.
References
[1]Walas, Stanley M. Phase Equilibria in Chemical Engineering. Butterworth-Heinemann, 1985.
[2]Poling, Bruce E. The Properties of Gases and Liquids. 5th edition. New York: McGraw-Hill Professional, 2000.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = VDWMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.V_l, eos.V_g (5.881369844883e-05, 0.00077708723758) >>> eos.fugacities_l, eos.fugacities_g ([854533.266920, 207126.8497276], [448470.736338, 397826.543999])
- a_alpha_and_derivatives_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas and their first and second derivatives for the VDW EOS. This vectorized implementation is added for extra speed.
- Parameters:
- T
float Temperature, [K]
- T
- Returns:
- a_alphas
list[float] Coefficient calculated by EOS-specific method, [J^2/mol^2/Pa]
- da_alpha_dTs
list[float] Temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K]
- d2a_alpha_dT2s
list[float] Second temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K**2]
- a_alphas
- a_alphas_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas for the VDW EOS. This vectorized implementation is added for extra speed.
- property d2delta_dninjs¶
Helper method for calculating the second mole number derivatives (hessian) of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2delta_dzizjs¶
Helper method for calculating the second composition derivatives (hessian) of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3delta_dninjnks¶
Helper method for calculating the third partial mole number derivatives of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property ddelta_dns¶
Helper method for calculating the mole number derivatives of delta. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- property ddelta_dzs¶
Helper method for calculating the composition derivatives of delta. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- dlnphis_dP(phase)[source]¶
Generic formula for calculating the pressure derivaitve of log fugacity coefficients for each species in a mixture for the VDW EOS. Verified numerically.
- Parameters:
- phase
str One of ‘l’ or ‘g’, [-]
- phase
- Returns:
- dlnphis_dP
float Pressure derivatives of log fugacity coefficient for each species, [1/Pa]
- dlnphis_dP
Notes
This expression was derived using SymPy and optimized with the cse technique.
- dlnphis_dT(phase)[source]¶
Formula for calculating the temperature derivaitve of log fugacity coefficients for each species in a mixture for the VDW equation of state. Verified numerically.
- Parameters:
- phase
str One of ‘l’ or ‘g’, [-]
- phase
- Returns:
- dlnphis_dT
float Temperature derivatives of log fugacity coefficient for each species, [1/K]
- dlnphis_dT
Notes
This expression was derived using SymPy and optimized with the cse technique.
- fugacity_coefficients(Z)[source]¶
Literature formula for calculating fugacity coefficients for each species in a mixture. Verified numerically. Called by fugacities on initialization, or by a solver routine which is performing a flash calculation.
- Parameters:
- Z
float Compressibility of the mixture for a desired phase, [-]
- Z
- Returns:
- log_phis
float Log fugacity coefficient for each species, [-]
- log_phis
References
[1]Walas, Stanley M. Phase Equilibria in Chemical Engineering. Butterworth-Heinemann, 1985.
Redlich-Kwong Equation of State¶
- class thermo.eos_mix.RKMIX(Tcs, Pcs, zs, omegas=None, kijs=None, T=None, P=None, V=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
EpsilonZeroMixingRules,GCEOSMIX,RKClass for solving the Redlich Kwong [1] [2] cubic equation of state for a mixture of any number of compounds. Subclasses
thermo.eos.RK. Solves the EOS on initialization and calculates fugacities for all components in all phases. Two of T, P, and V are needed to solve the EOS.- Parameters:
- Tcs
list[float] Critical temperatures of all compounds, [K]
- Pcs
list[float] Critical pressures of all compounds, [Pa]
- zs
list[float] Overall mole fractions of all species, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- omegas
list[float],optional Acentric factors of all compounds - Not used in this equation of state!, [-]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- Tcs
- Attributes:
d2delta_dninjsHelper method for calculating the second mole number derivatives (hessian) of delta.
d2delta_dzizjsHelper method for calculating the second composition derivatives (hessian) of delta.
d3delta_dninjnksHelper method for calculating the third partial mole number derivatives of delta.
ddelta_dnsHelper method for calculating the mole number derivatives of delta.
ddelta_dzsHelper method for calculating the composition derivatives of delta.
Methods
Method to calculate the pure-component a_alphas and their first and second derivatives for the RK EOS.
Method to calculate the pure-component a_alphas for the RK EOS.
Notes
The PV solution for T is iterative.
References
[1]Walas, Stanley M. Phase Equilibria in Chemical Engineering. Butterworth-Heinemann, 1985.
[2]Poling, Bruce E. The Properties of Gases and Liquids. 5th edition. New York: McGraw-Hill Professional, 2000.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = RKMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.V_l, eos.V_g (4.048414781e-05, 0.00070060605863)
- a_alpha_and_derivatives_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas and their first and second derivatives for the RK EOS. This vectorized implementation is added for extra speed.
- Parameters:
- T
float Temperature, [K]
- T
- Returns:
- a_alphas
list[float] Coefficient calculated by EOS-specific method, [J^2/mol^2/Pa]
- da_alpha_dTs
list[float] Temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K]
- d2a_alpha_dT2s
list[float] Second temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K**2]
- a_alphas
Examples
>>> eos = RKMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.a_alpha_and_derivatives_vectorized(115) ([0.1449810919468, 0.30019773677], [-0.000630352573681, -0.00130520755121], [8.2219900915e-06, 1.7024446320e-05])
- a_alphas_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas for the RK EOS. This vectorized implementation is added for extra speed.
- Parameters:
- T
float Temperature, [K]
- T
- Returns:
Examples
>>> eos = RKMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, 0.011], zs=[0.5, 0.5], kijs=[[0,0],[0,0]]) >>> eos.a_alphas_vectorized(115) [0.1449810919468, 0.30019773677]
- property d2delta_dninjs¶
Helper method for calculating the second mole number derivatives (hessian) of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d2delta_dzizjs¶
Helper method for calculating the second composition derivatives (hessian) of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property d3delta_dninjnks¶
Helper method for calculating the third partial mole number derivatives of delta. Note this is independent of the phase.
- Returns:
Notes
This derivative is checked numerically.
- property ddelta_dns¶
Helper method for calculating the mole number derivatives of delta. Note this is independent of the phase.
Notes
This derivative is checked numerically.
- property ddelta_dzs¶
Helper method for calculating the composition derivatives of delta. Note this is independent of the phase.
Notes
This derivative is checked numerically.
Ideal Gas Equation of State¶
- class thermo.eos_mix.IGMIX(zs, T=None, P=None, V=None, Tcs=None, Pcs=None, omegas=None, kijs=None, fugacities=True, only_l=False, only_g=False)[source]¶
Bases:
EpsilonZeroMixingRules,GCEOSMIX,IGClass for solving the ideal gas [1] [2] equation of state for a mixture of any number of compounds. Subclasses
thermo.eos.IG. Solves the EOS on initialization. Two of T, P, and V are needed to solve the EOS.- Parameters:
- zs
list[float] Overall mole fractions of all species, [-]
- T
float,optional Temperature, [K]
- P
float,optional Pressure, [Pa]
- V
float,optional Molar volume, [m^3/mol]
- Tcs
list[float],optional Critical temperatures of all compounds, [K]
- Pcs
list[float],optional Critical pressures of all compounds, [Pa]
- omegas
list[float],optional Acentric factors of all compounds - Not used in this equation of state!, [-]
- kijs
list[list[float]],optional n*n size list of lists with binary interaction parameters for the Van der Waals mixing rules, default all 0 and not used[-]
- fugacitiesbool,
optional Whether or not to calculate fugacity related values (phis, log phis, and fugacities); default True, [-]
- only_lbool,
optional When true, if there is a liquid and a vapor root, only the liquid root (and properties) will be set; default False, [-]
- only_gbool,
optional When true, if there is a liquid and a vapor root, only the vapor root (and properties) will be set; default False, [-]
- zs
Methods
Method to calculate the pure-component a_alphas and their first and second derivatives for the Ideal Gas EOS.
Method to calculate the pure-component a_alphas for the Ideal Gas EOS.
Notes
Many properties of this object are zero. Many of the arguments are not used and are provided for consistency only.
References
[1]Walas, Stanley M. Phase Equilibria in Chemical Engineering. Butterworth-Heinemann, 1985.
[2]Poling, Bruce E. The Properties of Gases and Liquids. 5th edition. New York: McGraw-Hill Professional, 2000.
Examples
T-P initialization, nitrogen-methane at 115 K and 1 MPa:
>>> eos = IGMIX(T=115, P=1E6, Tcs=[126.1, 190.6], Pcs=[33.94E5, 46.04E5], omegas=[0.04, .008], zs=[0.5, 0.5]) >>> eos.phase, eos.V_g ('g', 0.0009561632010876225)
- a_alpha_and_derivatives_vectorized(T)[source]¶
Method to calculate the pure-component a_alphas and their first and second derivatives for the Ideal Gas EOS. This vectorized implementation is added for extra speed.
- Parameters:
- T
float Temperature, [K]
- T
- Returns:
- a_alphas
list[float] Coefficient calculated by EOS-specific method, [J^2/mol^2/Pa]
- da_alpha_dTs
list[float] Temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K]
- d2a_alpha_dT2s
list[float] Second temperature derivative of coefficient calculated by EOS-specific method, [J^2/mol^2/Pa/K**2]
- a_alphas
Different Mixing Rules¶
- class thermo.eos_mix.EpsilonZeroMixingRules[source]¶
- Attributes:
d2epsilon_dninjsHelper method for calculating the second mole number derivatives (hessian) of epsilon.
d2epsilon_dzizjsHelper method for calculating the second composition derivatives (hessian) of epsilon.
d3epsilon_dninjnksHelper method for calculating the third partial mole number derivatives of epsilon.
depsilon_dnsHelper method for calculating the mole number derivatives of epsilon.
depsilon_dzsHelper method for calculating the composition derivatives of epsilon.
- class thermo.eos_mix.PSRKMixingRules[source]¶
Bases:
objectMethods
a_alpha_and_derivatives(T[, full, quick, ...])Method to calculate a_alpha and its first and second derivatives for an EOS with the PSRK mixing rules.
- A = -0.6466271649250525¶
- a_alpha_and_derivatives(T, full=True, quick=True, pure_a_alphas=True)[source]¶
Method to calculate a_alpha and its first and second derivatives for an EOS with the PSRK mixing rules. Returns a_alpha, da_alpha_dT, and d2a_alpha_dT2.
For use in some methods, this returns only a_alpha if full is False.
- Parameters:
- T
float Temperature, [K]
- fullbool,
optional If False, calculates and returns only a_alpha
- quickbool,
optional Only the quick variant is implemented; it is little faster anyhow
- pure_a_alphasbool,
optional Whether or not to recalculate the a_alpha terms of pure components (for the case of mixtures only) which stay the same as the composition changes (i.e in a PT flash), [-]
- T
- Returns:
- a_alpha
float Coefficient calculated by PSRK-specific method, [J^2/mol^2/Pa]
- da_alpha_dT
float Temperature derivative of coefficient calculated by PSRK-specific method, [J^2/mol^2/Pa/K]
- d2a_alpha_dT2
float Second temperature derivative of coefficient calculated by PSRK-specific method, [J^2/mol^2/Pa/K**2]
- a_alpha
- u = 1.1¶
Lists of Equations of State¶
- thermo.eos_mix.eos_mix_list = [<class 'thermo.eos_mix.PRMIX'>, <class 'thermo.eos_mix.SRKMIX'>, <class 'thermo.eos_mix.PR78MIX'>, <class 'thermo.eos_mix.VDWMIX'>, <class 'thermo.eos_mix.PRSVMIX'>, <class 'thermo.eos_mix.PRSV2MIX'>, <class 'thermo.eos_mix.TWUPRMIX'>, <class 'thermo.eos_mix.TWUSRKMIX'>, <class 'thermo.eos_mix.APISRKMIX'>, <class 'thermo.eos_mix.IGMIX'>, <class 'thermo.eos_mix.RKMIX'>, <class 'thermo.eos_mix.PRMIXTranslatedConsistent'>, <class 'thermo.eos_mix.PRMIXTranslatedPPJP'>, <class 'thermo.eos_mix.SRKMIXTranslatedConsistent'>, <class 'thermo.eos_mix.PRMIXTranslated'>, <class 'thermo.eos_mix.SRKMIXTranslated'>]¶
List of all exported EOS classes.
- thermo.eos_mix.eos_mix_no_coeffs_list = [<class 'thermo.eos_mix.PRMIX'>, <class 'thermo.eos_mix.SRKMIX'>, <class 'thermo.eos_mix.PR78MIX'>, <class 'thermo.eos_mix.VDWMIX'>, <class 'thermo.eos_mix.TWUPRMIX'>, <class 'thermo.eos_mix.TWUSRKMIX'>, <class 'thermo.eos_mix.IGMIX'>, <class 'thermo.eos_mix.RKMIX'>, <class 'thermo.eos_mix.PRMIXTranslatedConsistent'>, <class 'thermo.eos_mix.PRMIXTranslated'>, <class 'thermo.eos_mix.SRKMIXTranslated'>, <class 'thermo.eos_mix.PRMIXTranslatedPPJP'>, <class 'thermo.eos_mix.SRKMIXTranslatedConsistent'>]¶
List of all exported EOS classes that do not require special parameters or can fill in their special parameters from other specified parameters.