Methods in the thermo class
ThermoPack - 2.1.0
The thermo
class, found in addon/pycThermopack/thermopack/thermo.py
, is the core of the ThermoPack Python interface. All equation of state classes inherit from thermo
. This is the class that contains the interface to all practical calculations that can be done from the python-side of ThermoPack. Derived classes only implement specific functions for parameter handling etc.
Table of contents
- TV-property interfaces
- Tp-property interfaces
- TVp-property interfaces
- Other property interfaces
- Flash interfaces
- Saturation interfaces
- Isolines
- Stability interfaces
- Virial interfaces
- Joule-Thompson interface
- Utility methods
- acentric_factor
- compmoleweight
- get_comp_name
- get_ideal_enthalpy_reference_value
- get_ideal_entropy_reference_value
- get_phase_flags
- get_phase_type
- get_pmax
- get_pmin
- get_tmax
- get_tmin
- getcompindex
- redefine_critical_parameters
- set_ideal_enthalpy_reference_value
- set_ideal_entropy_reference_value
- set_pmax
- set_pmin
- set_tmax
- set_tmin
- Internal methods
TV-property interfaces
Computing properties as a function of temperature and volume. Derivatives returned by methods in this section are computed as functions of (T, V, n).
Table of contents
chemical_potential_tv(self, temp, volume, n, dmudt=None, dmudv=None, dmudn=None, property_flag='IR')
Calculate chemical potential given temperature, volume and mol numbers.
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
dmudt (No type, optional):
Flag to activate calculation. Defaults to None.
dmudv (No type, optional):
Flag to activate calculation. Defaults to None.
dmudn (No type, optional):
Flag to activate calculation. Defaults to None.
property_flag (integer, optional):
Calculate residual (R) and/or ideal (I) entropy. Defaults to IR.
Returns:
float:
Chemical potential (J/mol)
Optionally chemical potential differentials
enthalpy_tv(self, temp, volume, n, dhdt=None, dhdv=None, dhdn=None, property_flag='IR')
Calculate enthalpy given temperature, volume and mol numbers.
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
dhdt (No type, optional):
Flag to activate calculation. Defaults to None.
dhdv (No type, optional):
Flag to activate calculation. Defaults to None.
dhdn (No type, optional):
Flag to activate calculation. Defaults to None.
property_flag (integer, optional):
Calculate residual (R) and/or ideal (I) entropy. Defaults to IR.
Returns:
float:
Enthalpy (J)
Optionally enthalpy differentials
entropy_tv(self, temp, volume, n, dsdt=None, dsdv=None, dsdn=None, property_flag='IR')
Calculate entropy given temperature, volume and mol numbers.
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
dsdt (No type, optional):
Flag to activate calculation. Defaults to None.
dsdv (No type, optional):
Flag to activate calculation. Defaults to None.
dsdn (No type, optional):
Flag to activate calculation. Defaults to None.
property_flag (integer, optional):
Calculate residual (R) and/or ideal (I) entropy. Defaults to IR.
Returns:
float:
Entropy (J/K)
Optionally entropy differentials
fugacity_tv(self, temp, volume, n, dlnphidt=None, dlnphidv=None, dlnphidn=None)
Calculate natural logarithm of fugacity given temperature, volume and mol numbers.
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
dlnphidt (No type, optional):
Flag to activate calculation. Defaults to None.
dlnphidv (No type, optional):
Flag to activate calculation. Defaults to None.
dlnphidn (No type, optional):
Flag to activate calculation. Defaults to None.
Returns:
ndarry:
Natural logarithm of fugacity
Optionally differentials
helmholtz_tv(self, temp, volume, n, dadt=None, dadv=None, dadn=None, property_flag='IR')
Calculate Helmholtz energy given temperature, volume and mol numbers.
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
dadt (No type, optional):
Flag to activate calculation. Defaults to None.
dadv (No type, optional):
Flag to activate calculation. Defaults to None.
dadn (No type, optional):
Flag to activate calculation. Defaults to None.
property_flag (integer, optional):
Calculate residual (R) and/or ideal (I) entropy. Defaults to IR.
Returns:
float:
Helmholtz energy (J)
Optionally energy differentials
internal_energy_tv(self, temp, volume, n, dedt=None, dedv=None, dedn=None, property_flag='IR')
Calculate internal energy given temperature, volume and mol numbers.
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
dedt (No type, optional):
Flag to activate calculation. Defaults to None.
dedv (No type, optional):
Flag to activate calculation. Defaults to None.
dedn (No type, optional):
Flag to activate calculation. Defaults to None.
property_flag (str, optional):
Calculate residual (‘R’), ideal (‘I’) or total (‘IR’) internal energy. Defaults to ‘IR’.
Returns:
float:
Energy (J)
Optionally energy differentials
pressure_tv(self, temp, volume, n, dpdt=None, dpdv=None, dpdn=None, property_flag='IR')
Calculate pressure given temperature, volume and mol numbers.
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
dpdt (No type, optional):
Flag to activate calculation. Defaults to None.
dpdv (No type, optional):
Flag to activate calculation. Defaults to None.
dpdn (No type, optional):
Flag to activate calculation. Defaults to None.
property_flag (str, optional):
Calculate residual (‘R’), ideal (‘I’) or total (‘IR’) pressure. Defaults to ‘IR’.
Returns:
float:
Pressure (Pa)
Optionally pressure differentials
speed_of_sound_tv(self, temp, volume, n)
Calculate speed of sound for single phase fluid
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
Returns:
float:
Speed of sound (m/s)
Tp-property interfaces
Computing properties as a function of temperature and pressure. Derivatives returned by methods in this section are computed as functions of (T, p, n).
Table of contents
enthalpy(self, temp, press, x, phase, dhdt=None, dhdp=None, dhdn=None, residual=False)
Calculate specific single-phase enthalpy Note that the order of the output match the default order of input for the differentials. Note further that dhdt, dhdp and dhdn only are flags to enable calculation.
Args:
temp (float):
Temperature (K)
press (float):
Pressure (Pa)
x (array_like):
Molar composition
phase (int):
Calculate root for specified phase
dhdt (logical, optional):
Calculate enthalpy differentials with respect to temperature while pressure and composition are held constant. Defaults to None.
dhdp (logical, optional):
Calculate enthalpy differentials with respect to pressure while temperature and composition are held constant. Defaults to None.
dhdn (logical, optional):
Calculate enthalpy differentials with respect to mol numbers while pressure and temperature are held constant. Defaults to None.
residual (logical, optional):
Calculate residual enthalpy. Defaults to False.
Returns:
float:
Specific enthalpy (J/mol), and optionally differentials
entropy(self, temp, press, x, phase, dsdt=None, dsdp=None, dsdn=None, residual=False)
Calculate specific single-phase entropy Note that the order of the output match the default order of input for the differentials. Note further that dsdt, dhsp and dsdn only are flags to enable calculation.
Args:
temp (float):
Temperature (K)
press (float):
Pressure (Pa)
x (array_like):
Molar composition
phase (int):
Calculate root for specified phase
dsdt (logical, optional):
Calculate entropy differentials with respect to temperature while pressure and composition are held constant. Defaults to None.
dsdp (logical, optional):
Calculate entropy differentials with respect to pressure while temperature and composition are held constant. Defaults to None.
dsdn (logical, optional):
Calculate entropy differentials with respect to mol numbers while pressure and temperature are held constant. Defaults to None.
residual (logical, optional):
Calculate residual entropy. Defaults to False.
Returns:
float:
Specific entropy (J/mol/K), and optionally differentials
idealenthalpysingle(self, temp, j, dhdt=None)
Calculate specific ideal enthalpy Note that the order of the output match the default order of input for the differentials. Note further that dhdt only are flags to enable calculation.
Args:
temp (float):
Temperature (K)
j (int):
Component index (FORTRAN)
dhdt (logical, optional):
Calculate ideal enthalpy differentials with respect to temperature while pressure and composition are held constant. Defaults to None.
Returns:
float:
Specific ideal enthalpy (J/mol), and optionally differentials
idealentropysingle(self, temp, press, j, dsdt=None, dsdp=None)
Calculate specific ideal entropy Note that the order of the output match the default order of input for the differentials. Note further that dhdt, and dhdp only are flags to enable calculation.
Args:
temp (float):
Temperature (K)
press (float):
Pressure (Pa)
j (int):
Component index (FORTRAN)
dsdt (logical, optional):
Calculate ideal entropy differentials with respect to temperature while pressure and composition are held constant. Defaults to None.
dsdp (logical, optional):
Calculate ideal entropy differentials with respect to pressure while temperature and composition are held constant. Defaults to None.
Returns:
float:
Specific ideal entropy (J/mol/K), and optionally differentials
specific_volume(self, temp, press, x, phase, dvdt=None, dvdp=None, dvdn=None)
Calculate single-phase specific volume Note that the order of the output match the default order of input for the differentials. Note further that dvdt, dvdp and dvdn only are flags to enable calculation.
Args:
temp (float):
Temperature (K)
press (float):
Pressure (Pa)
x (array_like):
Molar composition
phase (int):
Calculate root for specified phase
dvdt (logical, optional):
Calculate molar volume differentials with respect to temperature while pressure and composition are held constant. Defaults to None.
dvdp (logical, optional):
Calculate molar volume differentials with respect to pressure while temperature and composition are held constant. Defaults to None.
dvdn (logical, optional):
Calculate molar volume differentials with respect to mol numbers while pressure and temperature are held constant. Defaults to None.
Returns:
float:
Specific volume (m3/mol), and optionally differentials
speed_of_sound(self, temp, press, x, y, z, betaV, betaL, phase)
Calculate speed of sound for single phase or two phase mixture assuming mechanical, thermal and chemical equilibrium.
Args:
temp (float):
Temperature (K)
press (float):
Pressure (Pa)
x (array_like):
Liquid molar composition
y (array_like):
Gas molar composition
z (array_like):
Overall molar composition
betaV (float):
Molar gas phase fraction
betaL (float):
Molar liquid phase fraction
phase (int):
Calculate root for specified phase
Returns:
float:
Speed of sound (m/s)
thermo(self, temp, press, x, phase, dlnfugdt=None, dlnfugdp=None, dlnfugdn=None, ophase=None, v=None)
Calculate logarithm of fugacity coefficient given composition, temperature and pressure. Note that the order of the output match the default order of input for the differentials. Note further that dlnfugdt, dlnfugdp, dlnfugdn and ophase only are flags to enable calculation.
Args:
temp (float):
Temperature (K)
press (float):
Pressure (Pa)
x (array_like):
Molar composition (.)
phase (int):
Calculate root for specified phase
dlnfugdt (logical, optional):
Calculate fugacity coefficient differentials with respect to temperature while pressure and composition are held constant. Defaults to None.
dlnfugdp (logical, optional):
Calculate fugacity coefficient differentials with respect to pressure while temperature and composition are held constant. Defaults to None.
dlnfugdn (logical, optional):
Calculate fugacity coefficient differentials with respect to mol numbers while pressure and temperature are held constant. Defaults to None.
ophase (int, optional):
Phase flag. Only set when phase=MINGIBBSPH.
v (float, optional):
Specific volume (m3/mol)
Returns:
ndarray:
fugacity coefficient (-), and optionally differentials
zfac(self, temp, press, x, phase, dzdt=None, dzdp=None, dzdn=None)
Calculate single-phase compressibility Note that the order of the output match the default order of input for the differentials. Note further that dzdt, dzdp and dzdn only are flags to enable calculation.
Args:
temp (float):
Temperature (K)
press (float):
Pressure (Pa)
x (array_like):
Molar composition
phase (int):
Calculate root for specified phase
dzdt (logical, optional):
Calculate compressibility differentials with respect to temperature while pressure and composition are held constant. Defaults to None.
dzdp (logical, optional):
Calculate compressibility differentials with respect to pressure while temperature and composition are held constant. Defaults to None.
dzdn (logical, optional):
Calculate compressibility differentials with respect to mol numbers while pressure and temperature are held constant. Defaults to None.
Returns:
float:
Compressibility (-), and optionally differentials
TVp-property interfaces
Computing properties given Temperature, volume and mole numbers, but evaluate derivatives as functions of (T, p, n). See Advanced Usage => The different property interfaces for further explanation.
Table of contents
enthalpy_tvp(self, temp, volume, n, dhdt=None, dhdp=None, dhdn=None, property_flag='IR')
Calculate enthalpy given temperature, volume and mol numbers.
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
dhdt (No type, optional):
Flag to activate calculation. Defaults to None.
dhdp (No type, optional):
Flag to activate calculation. Defaults to None.
dhdn (No type, optional):
Flag to activate calculation. Defaults to None.
property_flag (integer, optional):
Calculate residual (R) and/or ideal (I) entropy. Defaults to IR.
Returns:
float:
Enthalpy (J)
Optionally enthalpy differentials
entropy_tvp(self, temp, volume, n, dsdt=None, dsdp=None, dsdn=None, property_flag='IR')
Calculate entropy given temperature, pressure and mol numbers.
Args:
temp (float):
Temperature (K)
volume (float):
Volume (m3)
n (array_like):
Mol numbers (mol)
dsdt (No type, optional):
Flag to activate calculation. Defaults to None.
dsdp (No type, optional):
Flag to activate calculation. Defaults to None.
dsdn (No type, optional):
Flag to activate calculation. Defaults to None.
property_flag (integer, optional):
Calculate residual (R) and/or ideal (I) entropy. Defaults to IR.
Returns:
float:
Entropy (J/K)
Optionally entropy differentials
thermo_tvp(self, temp, v, n, phase, dlnfugdt=None, dlnfugdp=None, dlnfugdn=None)
Calculate logarithm of fugacity coefficient given molar numbers, temperature and pressure. Note that the order of the output match the default order of input for the differentials. Note further that dlnfugdt, dlnfugdp, dlnfugdn and ophase only are flags to enable calculation.
Args:
temp (float):
Temperature (K)
v (float):
Volume (m3)
n (array_like):
Molar numbers (mol)
dlnfugdt (logical, optional):
Calculate fugacity coefficient differentials with respect to temperature while pressure and composition are held constant. Defaults to None.
dlnfugdp (logical, optional):
Calculate fugacity coefficient differentials with respect to pressure while temperature and composition are held constant. Defaults to None.
dlnfugdn (logical, optional):
Calculate fugacity coefficient differentials with respect to mol numbers while pressure and temperature are held constant. Defaults to None.
Returns:
ndarray:
fugacity coefficient (-), and optionally differentials
Other property interfaces
Property interfaces in other variables than $TV$ or $Tp$, for example computing density given $\mu - T$.
Table of contents
density_lnf_t(self, temp, lnf, rho_initial)
Solve densities (lnf=lnf(T,rho)) given temperature and fugcaity coefficients.
Args:
temp (float):
Temperature (K)
lnf (array_like):
Logaritm of fugacity coefficients.
rho_initial (array_like):
Initial guess for component densities (mol/m3).
Returns:
rho (array_like):
Array of component densities (mol/m3).
density_mu_t(self, temp, mu, rho_initial)
Solve for densities (mu=mu(T,rho)) given temperature and chemical potential.
Args:
temp (float):
Temperature (K)
mu (array_like):
Flag to activate calculation.
rho_initial (array_like):
Initial guess for component densities (mol/m3).
Returns:
rho (array_like):
Array of component densities (mol/m3).
Flash interfaces
Methods for flash calculations.
Table of contents
guess_phase(self, temp, press, z)
If only one root exsist for the equation of state the phase type can be determined from either the psedo-critical volume or a volume ratio to the co-volume
Args:
temp (float):
Temperature (K)
press (float):
Pressure (Pa)
Returns:
int:
Phase int (VAPPH or LIQPH)
set_ph_tolerance(self, tol)
Set tolerance of isobaric-isentalpic (PH) flash
Args:
tol (float):
Tolerance
two_phase_phflash(self, press, z, enthalpy, temp=None)
Do isenthalpic-isobaric (HP) flash
Args:
press (float):
Pressure (Pa)
z (array_like):
Overall molar composition
enthalpy (float):
Specific enthalpy (J/mol)
temp (float, optional):
Initial guess for temperature (K)
Returns:
tuple :
(Temperature, x, y, betaV, betaL, phase)
two_phase_psflash(self, press, z, entropy, temp=None)
Do isentropic-isobaric (SP) flash
Args:
press (float):
Pressure (Pa)
z (array_like):
Overall molar composition
entropy (float):
Specific entropy (J/mol/K)
temp (float, optional):
Initial guess for temperature (K)
Returns:
tuple :
(Temperature, x, y, betaV, betaL, phase)
two_phase_tpflash(self, temp, press, z)
Do isothermal-isobaric (TP) flash
Args:
temp (float):
Temperature (K)
press (float):
Pressure (Pa)
z (array_like):
Overall molar composition
Returns:
tuple :
(x, y, betaV, betaL, phase)
two_phase_uvflash(self, z, specific_energy, specific_volume, temp=None, press=None)
Do isoenergetic-isochoric (UV) flash
Args:
press (float):
Pressure (Pa)
z (array_like):
Overall molar composition
specific_energy (float):
Specific energy (J/mol)
specific_volume (float):
Specific volume (m3/mol)
temp (float, optional):
Initial guess for temperature (K)
press (float, optional):
Initial guess for pressure (Pa)
Returns:
tuple :
(T, p, x, y, betaV, betaL, phase)
Saturation interfaces
Bubble- and dew point calculations and phase envelopes.
Table of contents
binary_triple_point_pressure(self, temp, maximum_pressure=15000000.0, minimum_pressure=10000.0)
Calculate triple point for binary mixture at specified temperature
Args:
temp (float):
Temperature (K)
maximum_pressure (float, optional):
Exit on maximum pressure (Pa). Defaults to 1.5e7.
minimum_pressure (float, optional):
Exit on minimum pressure (Pa). Defaults to 1.0e4.
Returns:
has_triple_point (boolean):
Does the mixture have a triple point?
x (np.ndarray):
Liquid 1 composition
y (np.ndarray):
Gas composition
w (np.ndarray):
Liquid 2 composition
P (float):
Pressure (Pa)
bubble_pressure(self, temp, z)
Calculate bubble pressure given temperature and composition
Args:
temp (float):
Temperature (K)
z (array_like):
Composition (-)
Raises:
Exception:
Faild to calculate
Returns:
float:
Pressure (Pa)
ndarray:
Incipient phase composition
bubble_temperature(self, press, z)
Calculate bubble temperature given pressure and composition
Args:
press (float):
Pressure (Pa)
z (array_like):
Composition (-)
Raises:
Exception:
Faild to calculate
Returns:
float:
Temperature (K)
ndarray:
Incipient phase composition
dew_pressure(self, temp, z)
Calculate dew pressure given temperature and composition
Args:
temp (float):
Temperature (K)
z (float):
Composition (-)
Raises:
Exception:
Not able to solve for dew point
Returns:
float :
Pressure (Pa)
ndarray :
Incipient phase composition (-)
dew_temperature(self, press, z)
Calculate dew temperature given pressure and composition
Args:
press (float):
Pressure (Pa)
z (float):
Composition (-)
Raises:
Exception:
Not able to solve for dew point
Returns:
float :
Temperature (K)
ndarray :
Incipient phase composition (-)
get_binary_pxy(self, temp, maximum_pressure=15000000.0, minimum_pressure=100000.0, maximum_dz=0.003, maximum_dlns=0.01)
Calculate binary three phase envelope
Args:
temp (float):
Temperature (K)
maximum_pressure (float, optional):
Exit on maximum pressure (Pa). Defaults to 1.5e7.
minimum_pressure (float, optional):
Exit on minimum pressure (Pa). Defaults to 1.0e5.
maximum_dz (float, optional):
Maximum composition step. Defaults to 0.003.
maximum_dlns (float, optional):
Maximum step in most sensitive envelope variable (the specification variable), see doc/memo/binaryxy
for details on usage. Defaults to 0.01.
Returns:
tuple of arrays:
LLE, L1VE, L2VE
LLE :
Liquid 1 - Liquid 2 Equilibrium
LLE[0] -> Liquid 1 composition (mole fraction of component 1)
LLE[1] -> Liquid 2 composition (mole fraction of component 1)
LLE[2] -> Pressure [Pa]
L1VE :
Liquid 1 - Vapour Equilibrium
L1VE[0] -> Bubble line composition (mole fraction of component 1)
L1VE[1] -> Dew line composition (mole fraction of component 1)
L1VE[2] -> Pressure [Pa]
L2VE :
Liquid 2 - Vapour Equilibrium
L2VE[0] -> Bubble line composition (mole fraction of component 1)
L2VE[1] -> Dew line composition (mole fraction of component 1)
L2VE[2] -> Pressure [Pa]
If one or more of the equilibria are not found the corresponding tuple is (None, None, None)
get_bp_term(self, i_term)
Get error description for binary plot error
Args:
i_term (int):
binary plot error identifyer
Returns:
str:
Error message
get_envelope_twophase(self, initial_pressure, z, maximum_pressure=15000000.0, minimum_temperature=None, step_size_factor=1.0, step_size=None, calc_v=False, initial_temperature=None)
Get the phase-envelope at a given composition
Args:
initial_pressure (float):
Start mapping form dew point at initial pressure (Pa).
z (array_like):
Composition (-)
maximum_pressure (float , optional):
Exit on maximum pressure (Pa). Defaults to 1.5e7.
minimum_temperature (float , optional):
Exit on minimum pressure (Pa). Defaults to None.
step_size_factor (float , optional):
Scale default step size for envelope trace. Defaults to 1.0. Reducing step_size_factor will give a denser grid.
step_size (float , optional):
Set maximum step size for envelope trace. Overrides step_size_factor. Defaults to None.
calc_v (bool, optional):
Calculate specific volume of saturated phase? Defaults to False
initial_temperature (bool, optional):
Start mapping form dew point at initial temperature.
Overrides initial pressure. Defaults to None (K).
Returns:
ndarray:
Temperature values (K)
ndarray:
Pressure values (Pa)
ndarray (optional, if calc_v=True
):
Specific volume (m3/mol)
get_pure_fluid_saturation_curve(self, initial_pressure, initial_temperature=None, i=None, max_delta_press=20000.0, nmax=100, log_linear_grid=False)
Get the pure fluid saturation line
To start mapping from and initial temperature, use:
get_pure_fluid_saturation_curve(None, initial_temperature=
Args:
initial_pressure (float):
Start mapping form dew point at initial pressure (Pa).
initial_temperature (float, optional):
Start mapping form dew point at initial temperature (K). Default None.
i (int, optional):
FORTRAN component index. Default None. Must be given if self.nc > 1.
max_delta_press (float , optional):
Maximum delta pressure between points (Pa). Defaults to 0.2e5.
nmax (int, optional):
Maximum number of points on envelope. Defaults to 100.
log_linear_grid (logical, optional):
Use log-linear grid?. Defaults to False.
Returns:
ndarray:
Temperature values (K)
ndarray:
Pressure values (Pa)
ndarray:
Specific liquid volume (m3/mol)
ndarray:
Specific gas volume (m3/mol)
global_binary_plot(self, maximum_pressure=15000000.0, minimum_pressure=100000.0, minimum_temperature=150.0, maximum_temperature=500.0, include_azeotropes=False)
Calculate global binary phase envelope
Args:
maximum_pressure (float, optional):
Exit on maximum pressure (Pa). Defaults to 1.5e7.
minimum_pressure (float, optional):
Exit on minimum pressure (Pa). Defaults to 1.0e5.
minimum_temperature (float, optional):
Terminate phase line traceing at minimum temperature. Defaults to 150.0 K.
maximum_temperature (float, optional):
Terminate phase line traceing at maximum temperature. Defaults to 500.0 K.
include_azeotropes (bool, optional):
Include azeotropic lines. Defaults to False.
Returns:
tuple of arrays
melting_pressure_correlation(self, i, maximum_temperature=None, nmax=100, scale_to_eos=True)
Calculate melting line form correlation
Args:
i (int):
component FORTRAN index (first index is 1)
maximum_temperature (float, optional):
Get values up to maximum_temperature. Defaults to correlation limit.
nmax (int):
Number of points in equidistant grid. Defaults to 100.
scale_to_eos (bool, optional):
Scale pressures to match triple point pressure? Defaults to True
Returns:
T_melt (ndarray):
Melting temperature (K)
p_melt (ndarray):
Melting pressure (Pa)
solid_envelope_plot(self, initial_pressure, z, maximum_pressure=15000000.0, minimum_temperature=170.0, calc_esv=False)
Calculate phase envelope including solid lines
Args:
initial_pressure (float):
Start mapping from initial pressure (Pa).
z (array_like):
Composition (-)
maximum_pressure (float , optional):
Exit on maximum pressure (Pa). Defaults to 1.5e7.
calc_esv (bool, optional):
Calculate specific volume of saturated phase? Defaults to False
Returns:
tuple of arrays
sublimation_pressure_correlation(self, i, minimum_temperature=None, nmax=100, scale_to_eos=True)
Calculate melting line form correlation
Args:
i (int):
component FORTRAN index (first index is 1)
minimum_temperature (float, optional):
Get values from minimum_temperature. Defaults to correlation limit.
nmax (int):
Number of points in equidistant grid. Defaults to 100.
scale_to_eos (bool, optional):
Scale pressures to match triple point pressure? Defaults to True
Returns:
T_subl (ndarray):
Sublimation temperature (K)
p_subl (ndarray):
Sublimation pressure (Pa)
Isolines
Computing isolines.
Table of contents
get_isenthalp(self, enthalpy, z, minimum_pressure=100000.0, maximum_pressure=15000000.0, minimum_temperature=200.0, maximum_temperature=500.0, nmax=100)
Get isenthalpic line at specified enthalpy. Use as T, p, v, s = get_isenthalp(h, z)
, where (T, p, v, s)
is the temperature, pressure, specific volume and specific entropy along the isenthalp with specific enthalpy h
and molar composition z
.
Args:
enthalpy (float):
Enthalpy (J/mol)
z (array_like):
Composition (-)
minimum_pressure (float, optional):
Minimum pressure. Defaults to 1.0e5. (Pa)
maximum_pressure (float, optional):
Maximum pressure. Defaults to 1.5e7. (Pa)
minimum_temperature (float, optional):
Minimum temperature. Defaults to 200.0. (K)
maximum_temperature (float, optional):
Maximum temperature. Defaults to 500.0. (K)
nmax (int, optional):
Maximum number of points on isenthalp. Defaults to 100.
Returns:
(tuple of arrays) :
Corresponding to (temperature, pressure, specific volume, specific entropy) along the
isenthalp.
get_isentrope(self, entropy, z, minimum_pressure=100000.0, maximum_pressure=15000000.0, minimum_temperature=200.0, maximum_temperature=500.0, nmax=100)
Get isentrope at specified entropy. Use as T, p, v, h = get_isenthalp(s, z)
, where (T, p, v, h)
is the temperature, pressure, specific volume and specific enthalpy along the isentrope with specific entropy s
and molar composition z
.
Args:
entropy (float):
Entropy (J/mol/K)
z (array_like):
Composition (-)
minimum_pressure (float, optional):
Minimum pressure. Defaults to 1.0e5. (Pa)
maximum_pressure (float, optional):
Maximum pressure. Defaults to 1.5e7. (Pa)
minimum_temperature (float, optional):
Minimum temperature. Defaults to 200.0. (K)
maximum_temperature (float, optional):
Maximum temperature. Defaults to 500.0. (K)
nmax (int, optional):
Maximum number of points on isentrope. Defaults to 100.
Returns:
(tuple of arrays) :
Corresponding to (temperature, pressure, specific volume, specific enthalpy) along the
isentrope.
get_isobar(self, press, z, minimum_temperature=200.0, maximum_temperature=500.0, nmax=100)
Get isobar at specified pressure. Use as T, v, s, h = get_isobar(p, z)
, where (T, v, s, h)
is the temperature, specific volume, specific entropy and specific enthalpy along the isobar with pressure p
and molar composition z
.
Args:
press (float):
Pressure (Pa)
z (array_like):
Composition (-)
minimum_temperature (float, optional):
Minimum temperature. Defaults to 200.0. (K)
maximum_temperature (float, optional):
Maximum temperature. Defaults to 500.0. (K)
nmax (int, optional):
Maximum number of points on iso-bar. Defaults to 100.
Returns:
(tuple of arrays) :
Corresponding to (temperature, specific volume, specific entropy, specific enthalpy)
along the isobar.
get_isotherm(self, temp, z, minimum_pressure=100000.0, maximum_pressure=15000000.0, nmax=100)
Get iso-therm at specified temperature
Args:
temp (float):
Temperature (K)
z (array_like):
Composition (-)
minimum_pressure (float, optional):
Map to minimum pressure. Defaults to 1.0e5. (Pa)
maximum_pressure (float, optional):
Map to maximum pressure. Defaults to 1.5e7. (Pa)
nmax (int, optional):
Maximum number of points on iso-therm. Defaults to 100.
Returns:
Multiple numpy arrays.
Stability interfaces
Critical point solver.
Table of contents
critical(self, n, temp=0.0, v=0.0, tol=1e-07, v_min=None)
Calculate critical point in variables T and V
Args:
n (array_like):
Mol numbers (mol)
temp (float, optional):
Initial guess for temperature (K). Defaults to 0.0.
v (float, optional):
Initial guess for volume (m3/mol). Defaults to 0.0.
tol (float, optional):
Error tolerance (-). Defaults to 1.0e-8.
v_min (float, optional):
Minimum volume for search (m3/mol). Defaults to None.
Raises:
Exception:
Failure to solve for critical point
Returns:
float:
Temperature (K)
float:
Volume (m3/mol)
float:
Pressure (Pa)
Virial interfaces
Retrieve various virial coefficients.
Table of contents
binary_third_virial_matrix(self, temp)
Calculate composition-independent virial coefficients C, defined as P = RTrho + Brho2 + C*rho3 + O(rho**4) as rho->0. Including cross coefficients Currently the code only support binary mixtures
Args:
temp (float):
Temperature (K)
Returns:
ndarray:
C - Third virial coefficient matrix (m6/mol2)
second_virial_matrix(self, temp)
Calculate composition-independent virial coefficients B, defined as P = RTrho + Brho2 + C*rho3 + O(rho**4) as rho->0. Including cross coefficients.
Args:
temp (float):
Temperature (K)
Returns:
ndarray:
B - Second virial coefficient matrix (m3/mol)
virial_coeffcients(self, temp, n)
Calculate (composition-dependent) virial coefficients B and C, defined as P/RT = rho + B*rho2 + C*rho3 + O(rho**4) as rho->0.
Args:
temp (float):
Temperature
n (array_like):
Mol numbers (mol)
Returns:
float:
B (m3/mol)
float:
C (m6/mol2)
Joule-Thompson interface
Joule-Thompson inversion curves.
Table of contents
joule_thompson_inversion(self, z, nmax=1000)
Calculate Joule-Thompson inversion curve
Args:
z (array like):
Compozition
nmax (int):
Array size
Returns:
ndarray:
temp - Temperature (K)
ndarray:
press - Pressure (Pa)
ndarray:
vol - Volume (m3/mol)
Utility methods
Methods for setting … and getting …
Table of contents
- Utility methods
- acentric_factor
- compmoleweight
- get_comp_name
- get_ideal_enthalpy_reference_value
- get_ideal_entropy_reference_value
- get_phase_flags
- get_phase_type
- get_pmax
- get_pmin
- get_tmax
- get_tmin
- getcompindex
- redefine_critical_parameters
- set_ideal_enthalpy_reference_value
- set_ideal_entropy_reference_value
- set_pmax
- set_pmin
- set_tmax
- set_tmin
acentric_factor(self, i)
Get acentric factor of component i Args: i (int) component FORTRAN index returns: float: acentric factor
compmoleweight(self, comp)
Get component mole weight (g/mol)
Args:
comp (int):
Component FORTRAN index
Returns:
float:
Component mole weight (g/mol)
get_comp_name(self, index)
Get component name
Args:
index (int):
Component FORTRAN index
Returns:
comp (str):
Component name
get_ideal_enthalpy_reference_value(self, j)
Get specific ideal enthalpy reference value
Args:
j (integer):
Component index
Returns:
float:
Specific ideal enthalpy (J/mol)
get_ideal_entropy_reference_value(self, j)
Get specific ideal entropy reference value
Args:
j (integer):
Component index
Returns:
float:
Specific ideal entropy (J/mol/K)
get_phase_flags(self)
Get phase identifiers used by thermopack
Returns:
int:
Phase int identifiers
get_phase_type(self, i_phase)
Get phase type
Args:
i_phase (int):
Phase flag returned by thermopack
Returns:
str:
Phase type
get_pmax(self)
Get minimum pressure in Thermopack. Used to limit search domain for numerical solvers. Default value set on init is 100 MPa.
Args:
press (float):
Pressure (Pa)
get_pmin(self)
Get minimum pressure in Thermopack. Used to limit search domain for numerical solvers. Default value set on init is 10 Pa.
Args:
press (float):
Pressure (Pa)
get_tmax(self)
Get maximum temperature in Thermopack. Used to limit search domain for numerical solvers. Default value set on init is 999 K.
Returns:
float:
Temperature (K)
get_tmin(self)
Get minimum temperature in Thermopack. Used to limit search domain for numerical solvers. Default value set on init is 80 K.
Returns:
float:
Temperature (K)
getcompindex(self, comp)
Get component index
Args:
comp (str):
Component name
Returns:
int:
Component FORTRAN index
redefine_critical_parameters(self, silent=True, Tc_initials=None, vc_initials=None)
Recalculate critical properties of pure fluids
Args:
silent (bool):
Ignore warnings? Defaults to True
Tc_initials (array_like):
Initial value for pure fluid critical temperatures (K). Negative values will trigger use of SRK values from data base.
vc_initials (array_like):
Initial value for pure fluid critical volumes (m3/mol). Negative values will trigger use of SRK values from data base.
set_ideal_enthalpy_reference_value(self, j, h0)
Set specific ideal enthalpy reference value
Args:
j (integer):
Component index
h0 (float):
Ideal enthalpy reference (J/mol)
set_ideal_entropy_reference_value(self, j, s0)
Set specific ideal entropy reference value
Args:
j (integer):
Component index
s0 (float):
Ideal entropy reference (J/mol/K)
set_pmax(self, press)
Set minimum pressure in Thermopack. Used to limit search domain for numerical solvers. Default value set on init is 100 MPa.
Args:
press (float):
Pressure (Pa)
set_pmin(self, press)
Set minimum pressure in Thermopack. Used to limit search domain for numerical solvers. Default value set on init is 10 Pa.
Args:
press (float):
Pressure (Pa)
set_tmax(self, temp)
Set maximum temperature in Thermopack. Used to limit search domain for numerical solvers. Default value set on init is 999 K.
Args:
temp (float):
Temperature (K)
set_tmin(self, temp)
Set minimum temperature in Thermopack. Used to limit search domain for numerical solvers. Default value set on init is 80 K.
Args:
temp (float):
Temperature (K)
Internal methods
Methods for handling communication with the Fortran library.
Table of contents
__del__(self)
Delete FORTRAN memory allocated for this instance
__init__(self)
Load libthermopack.(so/dll) and initialize function pointers
activate(self)
Activate this instance of thermopack parameters for calculation
add_eos(self)
Allocate FORTRAN memory for this class instance
delete_eos(self)
de-allocate FORTRAN memory for this class instance
get_export_name(self, module, method)
Generate library export name based on module and method name
Args:
module (str):
Name of module
method (str):
Name of method
Returns:
str:
Library export name
get_model_id(self)
Get model identification
Returns:
str:
Eos name
init_peneloux_volume_translation(self, parameter_reference='Default')
Initialize Peneloux volume translations
Args:
parameter_reference (str):
String defining parameter set, Defaults to “Default”
init_solid(self, scomp)
Initialize pure solid
Args:
scomp (str):
Component name
init_thermo(self, eos, mixing, alpha, comps, nphases, liq_vap_discr_method=None, csp_eos=None, csp_ref_comp=None, kij_ref='Default', alpha_ref='Default', saft_ref='Default', b_exponent=None, TrendEosForCp=None, cptype=None, silent=None)
Initialize thermopack
Args:
eos (str):
Equation of state
mixing (str):
Mixture model for cubic eos
alpha (str):
Alpha formulations for cubic EOS
comps (string):
Comma separated list of components
nphases (int):
Maximum number of phases considered during multi-phase flash calculations
liq_vap_discr_method (int, optional):
Method to discriminate between liquid and vapor in case of an undefined single phase. Defaults to None.
csp_eos (str, optional):
Corresponding state equation. Defaults to None.
csp_ref_comp (str, optional):
CSP reference component. Defaults to None.
kij_ref (str, optional):
Data set identifiers. Defaults to “Default”.
alpha_ref (str, optional):
Data set identifiers. Defaults to “Default”.
saft_ref (str, optional):
Data set identifiers. Defaults to “Default”.
b_exponent (float, optional):
Exponent used in co-volume mixing. Defaults to None.
TrendEosForCp (str, optional):
Option to init trend for ideal gas properties. Defaults to None.
cptype (int array, optional):
Equation type number for Cp. Defaults to None.
silent (bool, optional):
Suppress messages during init?. Defaults to None.