Flash Calculations (thermo.flash)

This module contains classes and functions for performing flash calculations.

For reporting bugs, adding feature requests, or submitting pull requests, please use the GitHub issue tracker.

Main Interfaces

Pure Components

class thermo.flash.FlashPureVLS(constants, correlations, gas, liquids, solids, settings=BulkSettings(dP_dT='MOLE_WEIGHTED', dP_dV='MOLE_WEIGHTED', d2P_dV2='MOLE_WEIGHTED', d2P_dT2='MOLE_WEIGHTED', d2P_dTdV='MOLE_WEIGHTED', mu_LL='LOG_PROP_MASS_WEIGHTED', mu_LL_power_exponent=0.4, mu_VL='McAdams', mu_VL_power_exponent=0.4, k_LL='MASS_WEIGHTED', k_LL_power_exponent=0.4, k_VL='MASS_WEIGHTED', k_VL_power_exponent=0.4, sigma_LL='MASS_WEIGHTED', sigma_LL_power_exponent=0.4, T_liquid_volume_ref=298.15, T_normal=273.15, P_normal=101325.0, T_standard=288.15, P_standard=101325.0, T_gas_ref=288.15, P_gas_ref=101325.0, speed_of_sound='MOLE_WEIGHTED', kappa='MOLE_WEIGHTED', isobaric_expansion='MOLE_WEIGHTED', Joule_Thomson='MOLE_WEIGHTED', VL_ID='PIP', S_ID='d2P_dVdT', solid_sort_method='prop', liquid_sort_method='prop', liquid_sort_cmps=[], solid_sort_cmps=[], liquid_sort_cmps_neg=[], solid_sort_cmps_neg=[], liquid_sort_prop='DENSITY_MASS', solid_sort_prop='DENSITY_MASS', phase_sort_higher_first=True, water_sort='water not special', equilibrium_perturbation=1e-07))[source]

Bases: thermo.flash.flash_base.Flash

Class for performing flash calculations on pure-component systems. This class is subtantially more robust than using multicomponent algorithms on pure species. It is also faster. All parameters are also attributes.

The minimum information that is needed in addition to the Phase objects is:

  • MW

  • Vapor pressure curve if including liquids

  • Sublimation pressure curve if including solids

  • Functioning enthalpy models for each phase

Parameters
constantsChemicalConstantsPackage object

Package of chemical constants; these are used as boundaries at times, initial guesses other times, and in all cases these properties are accessible as attributes of the resulting EquilibriumState object, [-]

correlationsPropertyCorrelationsPackage

Package of chemical T-dependent properties; these are used as boundaries at times, for initial guesses other times, and in all cases these properties are accessible as attributes of the resulting EquilibriumState object, [-]

gasPhase object

A single phase which can represent the gas phase, [-]

liquidslist[Phase]

A list of phases for representing the liquid phase; normally only one liquid phase is present for a pure-component system, but multiple liquids are allowed for the really weird cases like having both parahydrogen and orthohydrogen. The liquid phase which calculates a lower Gibbs free energy is always used. [-]

solidslist[Phase]

A list of phases for representing the solid phase; it is very common for multiple solid forms of a compound to exist. For water ice, the list is very long - normally ice is in phase Ih but other phases are Ic, II, III, IV, V, VI, VII, VIII, IX, X, XI, XII, XIII, XIV, XV, XVI, Square ice, and Amorphous ice. It is less common for there to be published, reliable, thermodynamic models for these different phases; for water there is the IAPWS-06 model for Ih, and another model here for phases Ih, Ic, II, III, IV, V, VI, IX, XI, XII. [-]

settingsBulkSettings object

Object containing settings for calculating bulk and transport properties, [-]

Notes

The algorithms in this object are mostly from [1] and [2]. They all boil down to newton methods with analytical derivatives. The phase with the lowest Gibbs energy is the most stable if there are multiple solutions.

Phase input combinations which have specific simplifying assumptions (and thus more speed) are:

Additional information that can be provided in the ChemicalConstantsPackage object and PropertyCorrelationsPackage object that may help convergence is:

  • Tc, Pc, omega, Tb, and atoms

  • Gas heat capacity correlations

  • Liquid molar volume correlations

  • Heat of vaporization correlations

References

1

Poling, Bruce E., John M. Prausnitz, and John P. O`Connell. The Properties of Gases and Liquids. 5th edition. New York: McGraw-Hill Professional, 2000.

2

Gmehling, Jürgen, Michael Kleiber, Bärbel Kolbe, and Jürgen Rarey. Chemical Thermodynamics for Process Simulation. John Wiley & Sons, 2019.

Examples

Create all the necessary objects using all of the default parameters for decane and do a flash at 300 K and 1 bar:

>>> from thermo import ChemicalConstantsPackage, PRMIX, CEOSLiquid, CEOSGas, FlashPureVLS
>>> constants, correlations = ChemicalConstantsPackage.from_IDs(['decane'])
>>> eos_kwargs = dict(Tcs=constants.Tcs, Pcs=constants.Pcs, omegas=constants.omegas)
>>> liquid = CEOSLiquid(PRMIX, HeatCapacityGases=correlations.HeatCapacityGases, eos_kwargs=eos_kwargs)
>>> gas = CEOSGas(PRMIX, HeatCapacityGases=correlations.HeatCapacityGases, eos_kwargs=eos_kwargs)
>>> flasher = FlashPureVLS(constants, correlations, gas=gas, liquids=[liquid], solids=[])
>>> print(flasher.flash(T=300, P=1e5))
<EquilibriumState, T=300.0000, P=100000.0000, zs=[1.0], betas=[1.0], phases=[<CEOSLiquid, T=300 K, P=100000 Pa>]>

Working with steam:

>>> from thermo import FlashPureVLS, IAPWS95Liquid, IAPWS95Gas, iapws_constants, iapws_correlations
>>> liquid = IAPWS95Liquid(T=300, P=1e5, zs=[1])
>>> gas = IAPWS95Gas(T=300, P=1e5, zs=[1])
>>> flasher = FlashPureVLS(iapws_constants, iapws_correlations, gas, [liquid], [])
>>> PT = flasher.flash(T=800.0, P=1e7)
>>> PT.rho_mass()
29.1071839176
>>> print(flasher.flash(T=600, VF=.5))
<EquilibriumState, T=600.0000, P=12344824.3572, zs=[1.0], betas=[0.5, 0.5], phases=[<IAPWS95Gas, T=600 K, P=1.23448e+07 Pa>, <IAPWS95Liquid, T=600 K, P=1.23448e+07 Pa>]>
>>> print(flasher.flash(T=600.0, H=50802))
<EquilibriumState, T=600.0000, P=10000469.1288, zs=[1.0], betas=[1.0], phases=[<IAPWS95Gas, T=600 K, P=1.00005e+07 Pa>]>
>>> print(flasher.flash(P=1e7, S=104.))
<EquilibriumState, T=599.6790, P=10000000.0000, zs=[1.0], betas=[1.0], phases=[<IAPWS95Gas, T=599.679 K, P=1e+07 Pa>]>
>>> print(flasher.flash(V=.00061, U=55850))
<EquilibriumState, T=800.5922, P=10144789.0899, zs=[1.0], betas=[1.0], phases=[<IAPWS95Gas, T=800.592 K, P=1.01448e+07 Pa>]>
Attributes
VL_IG_hackbool

Whether or not to trust the saturation curve of the liquid phase; applied automatically to the GibbsExcessLiquid phase if there is a single liquid only, [-]

VL_EOS_hacksbool

Whether or not to trust the saturation curve of the EOS liquid phase; applied automatically to the CEOSLiquid phase if there is a single liquid only, [-]

TPV_HSGUA_guess_maxiterint

Maximum number of iterations to try when converging a shortcut model for flashes with one (T, P, V) spec and one (H, S, G, U, A) spec, [-]

TPV_HSGUA_guess_xtolfloat

Convergence tolerance in the iteration variable when converging a shortcut model for flashes with one (T, P, V) spec and one (H, S, G, U, A) spec, [-]

TPV_HSGUA_maxiterint

Maximum number of iterations to try when converging a flashes with one (T, P, V) spec and one (H, S, G, U, A) spec; this is on a per-phase basis, so if there is a liquid and a gas phase, the maximum number of iterations that could end up being tried would be twice this, [-]

TPV_HSGUA_xtolfloat

Convergence tolerance in the iteration variable dimension when converging a flash with one (T, P, V) spec and one (H, S, G, U, A) spec, [-]

TVF_maxiterint

Maximum number of iterations to try when converging a flashes with a temperature and vapor fraction specification, [-]

TVF_xtolfloat

Convergence tolerance in the temperature dimension when converging a flashes with a temperature and vapor fraction specification, [-]

PVF_maxiterint

Maximum number of iterations to try when converging a flashes with a pressure and vapor fraction specification, [-]

PVF_xtolfloat

Convergence tolerance in the pressure dimension when converging a flashes with a pressure and vapor fraction specification, [-]

TSF_maxiterint

Maximum number of iterations to try when converging a flashes with a temperature and solid fraction specification, [-]

TSF_xtolfloat

Convergence tolerance in the temperature dimension when converging a flashes with a temperature and solid fraction specification, [-]

PSF_maxiterint

Maximum number of iterations to try when converging a flashes with a pressure and solid fraction specification, [-]

PSF_xtolfloat

Convergence tolerance in the pressure dimension when converging a flashes with a pressure and solid fraction specification, [-]

Vapor-Liquid Systems

class thermo.flash.FlashVL(constants, correlations, gas, liquid, settings=BulkSettings(dP_dT='MOLE_WEIGHTED', dP_dV='MOLE_WEIGHTED', d2P_dV2='MOLE_WEIGHTED', d2P_dT2='MOLE_WEIGHTED', d2P_dTdV='MOLE_WEIGHTED', mu_LL='LOG_PROP_MASS_WEIGHTED', mu_LL_power_exponent=0.4, mu_VL='McAdams', mu_VL_power_exponent=0.4, k_LL='MASS_WEIGHTED', k_LL_power_exponent=0.4, k_VL='MASS_WEIGHTED', k_VL_power_exponent=0.4, sigma_LL='MASS_WEIGHTED', sigma_LL_power_exponent=0.4, T_liquid_volume_ref=298.15, T_normal=273.15, P_normal=101325.0, T_standard=288.15, P_standard=101325.0, T_gas_ref=288.15, P_gas_ref=101325.0, speed_of_sound='MOLE_WEIGHTED', kappa='MOLE_WEIGHTED', isobaric_expansion='MOLE_WEIGHTED', Joule_Thomson='MOLE_WEIGHTED', VL_ID='PIP', S_ID='d2P_dVdT', solid_sort_method='prop', liquid_sort_method='prop', liquid_sort_cmps=[], solid_sort_cmps=[], liquid_sort_cmps_neg=[], solid_sort_cmps_neg=[], liquid_sort_prop='DENSITY_MASS', solid_sort_prop='DENSITY_MASS', phase_sort_higher_first=True, water_sort='water not special', equilibrium_perturbation=1e-07))[source]

Bases: thermo.flash.flash_base.Flash

Class for performing flash calculations on one and two phase vapor and liquid multicomponent systems. Use FlashVLN for systems which can have multiple liquid phases.

The minimum information that is needed in addition to the Phase objects is:

  • MWs

  • Vapor pressure curve

  • Functioning enthalpy models for each phase

Parameters
constantsChemicalConstantsPackage object

Package of chemical constants; these are used as boundaries at times, initial guesses other times, and in all cases these properties are accessible as attributes of the resulting EquilibriumState object, [-]

correlationsPropertyCorrelationsPackage

Package of chemical T-dependent properties; these are used as boundaries at times, for initial guesses other times, and in all cases these properties are accessible as attributes of the resulting EquilibriumState object, [-]

gasPhase object

A single phase which can represent the gas phase, [-]

liquidPhase

A single phase which can represent the liquid phase, [-]

settingsBulkSettings object

Object containing settings for calculating bulk and transport properties, [-]

Notes

The algorithms in this object are mostly from [1], [2] and [3]. Sequential substitution without acceleration is used by default to converge two-phase systems.

Quasi-newton methods are used by default to converge bubble and dew point calculations.

Flashes with one (T, P, V) spec and one (H, S, G, U, A) spec are solved by a 1D search over PT flashes.

Additional information that can be provided in the ChemicalConstantsPackage object and PropertyCorrelationsPackage object that may help convergence is:

  • Tc, Pc, omega, Tb, and atoms

  • Gas heat capacity correlations

  • Liquid molar volume correlations

  • Heat of vaporization correlations

Warning

If this flasher is used on systems that can form two or more liquid phases, and the flash specs are in that region, there is no guarantee which solution is returned. Sometimes it is almost random, jumping back and forth and providing nasty discontinuities.

References

1

Michelsen, Michael L., and Jørgen M. Mollerup. Thermodynamic Models: Fundamentals & Computational Aspects. Tie-Line Publications, 2007.

2

Poling, Bruce E., John M. Prausnitz, and John P. O`Connell. The Properties of Gases and Liquids. 5th edition. New York: McGraw-Hill Professional, 2000.

3

Gmehling, Jürgen, Michael Kleiber, Bärbel Kolbe, and Jürgen Rarey. Chemical Thermodynamics for Process Simulation. John Wiley & Sons, 2019.

Examples

For the system methane-ethane-nitrogen with a composition [0.965, 0.018, 0.017], calculate the vapor fraction of the system and equilibrium phase compositions at 110 K and 1 bar. Use the Peng-Robinson equation of state and the chemsep sample interaction parameter database.

>>> from thermo import ChemicalConstantsPackage, CEOSGas, CEOSLiquid, PRMIX, FlashVL
>>> from thermo.interaction_parameters import IPDB
>>> constants, properties = ChemicalConstantsPackage.from_IDs(['methane', 'ethane', 'nitrogen'])
>>> kijs = IPDB.get_ip_asymmetric_matrix('ChemSep PR', constants.CASs, 'kij')
>>> kijs
[[0.0, -0.0059, 0.0289], [-0.0059, 0.0, 0.0533], [0.0289, 0.0533, 0.0]]
>>> eos_kwargs = {'Pcs': constants.Pcs, 'Tcs': constants.Tcs, 'omegas': constants.omegas, 'kijs': kijs}
>>> gas = CEOSGas(PRMIX, eos_kwargs=eos_kwargs, HeatCapacityGases=properties.HeatCapacityGases)
>>> liquid = CEOSLiquid(PRMIX, eos_kwargs=eos_kwargs, HeatCapacityGases=properties.HeatCapacityGases)
>>> flasher = FlashVL(constants, properties, liquid=liquid, gas=gas)
>>> zs = [0.965, 0.018, 0.017]
>>> PT = flasher.flash(T=110.0, P=1e5, zs=zs)
>>> PT.VF, PT.gas.zs, PT.liquid0.zs
(0.0890, [0.8688, 2.5765e-05, 0.13115], [0.9744, 0.01975, 0.00584])

A few more flashes with the same system to showcase the functionality of the flash interface:

>>> flasher.flash(P=1e5, VF=1, zs=zs).T
133.8
>>> flasher.flash(T=133, VF=0, zs=zs).P
515029.6
>>> flasher.flash(P=PT.P, H=PT.H(), zs=zs).T
110.0
>>> flasher.flash(P=PT.P, S=PT.S(), zs=zs).T
110.0
>>> flasher.flash(T=PT.T, H=PT.H(), zs=zs).T
110.0
>>> flasher.flash(T=PT.T, S=PT.S(), zs=zs).T
110.0
Attributes
PT_SS_MAXITERint

Maximum number of sequential substitution iterations to try when converging a two-phase solution, [-]

PT_SS_TOLfloat

Convergence tolerance in sequential substitution [-]

PT_SS_POLISHbool

When set to True, flashes which are very near a vapor fraction of 0 or 1 are converged to a higher tolerance to ensure the solution is correct; without this, a flash might converge to a vapor fraction of -1e-7 and be called single phase, but with this the correct solution may be found to be 1e-8 and will be correctly returned as two phase.[-]

PT_SS_POLISH_VFfloat

What tolerance to a vapor fraction of 0 or 1; this is an absolute vapor fraction value, [-]

PT_SS_POLISH_MAXITERint

Maximum number of sequential substitution iterations to try when converging a two-phase solution that has been detected to be very sensitive, with a vapor fraction near 0 or 1 [-]

PT_SS_POLISH_TOLfloat

Convergence tolerance in sequential substitution when converging a two-phase solution that has been detected to be very sensitive, with a vapor fraction near 0 or 1 [-]

PT_STABILITY_MAXITERint

Maximum number of iterations to try when converging a stability test, [-]

PT_STABILITY_XTOLfloat

Convergence tolerance in the stability test [-]

DEW_BUBBLE_VF_K_COMPOSITION_INDEPENDENT_XTOLfloat

Convergence tolerance in Newton solver for bubble, dew, and vapor fraction spec flashes when both the liquid and gas model’s K values do not dependent on composition, [-]

DEW_BUBBLE_QUASI_NEWTON_XTOLfloat

Convergence tolerance in quasi-Newton bubble and dew point flashes, [-]

DEW_BUBBLE_QUASI_NEWTON_MAXITERint

Maximum number of iterations to use in quasi-Newton bubble and dew point flashes, [-]

DEW_BUBBLE_NEWTON_XTOLfloat

Convergence tolerance in Newton bubble and dew point flashes, [-]

DEW_BUBBLE_NEWTON_MAXITERint

Maximum number of iterations to use in Newton bubble and dew point flashes, [-]

TPV_HSGUA_BISECT_XTOLfloat

Tolerance in the iteration variable when converging a flash with one (T, P, V) spec and one (H, S, G, U, A) spec using a bisection-type solver, [-]

TPV_HSGUA_BISECT_YTOLfloat

Absolute tolerance in the (H, S, G, U, A) spec when converging a flash with one (T, P, V) spec and one (H, S, G, U, A) spec using a bisection-type solver, [-]

TPV_HSGUA_BISECT_YTOL_ONLYbool

When True, the TPV_HSGUA_BISECT_XTOL setting is ignored and the flash is considered converged once TPV_HSGUA_BISECT_YTOL is satisfied, [-]

TPV_HSGUA_NEWTON_XTOLfloat

Tolerance in the iteration variable when converging a flash with one (T, P, V) spec and one (H, S, G, U, A) spec using a full newton solver, [-]

TPV_HSGUA_NEWTON_MAXITERfloat

Maximum number of iterations when converging a flash with one (T, P, V) spec and one (H, S, G, U, A) spec using full newton solver, [-]

TPV_HSGUA_SECANT_MAXITERfloat

Maximum number of iterations when converging a flash with one (T, P, V) spec and one (H, S, G, U, A) spec using a secant solver, [-]

HSGUA_NEWTON_ANALYTICAL_JACbool

Whether or not to calculate the full newton jacobian analytically or numerically; this would need to be set to False if the phase objects used in the flash do not have complete analytical derivatives implemented, [-]

Vapor and Multiple Liquid Systems

class thermo.flash.FlashVLN(constants, correlations, liquids, gas, solids=None, settings=BulkSettings(dP_dT='MOLE_WEIGHTED', dP_dV='MOLE_WEIGHTED', d2P_dV2='MOLE_WEIGHTED', d2P_dT2='MOLE_WEIGHTED', d2P_dTdV='MOLE_WEIGHTED', mu_LL='LOG_PROP_MASS_WEIGHTED', mu_LL_power_exponent=0.4, mu_VL='McAdams', mu_VL_power_exponent=0.4, k_LL='MASS_WEIGHTED', k_LL_power_exponent=0.4, k_VL='MASS_WEIGHTED', k_VL_power_exponent=0.4, sigma_LL='MASS_WEIGHTED', sigma_LL_power_exponent=0.4, T_liquid_volume_ref=298.15, T_normal=273.15, P_normal=101325.0, T_standard=288.15, P_standard=101325.0, T_gas_ref=288.15, P_gas_ref=101325.0, speed_of_sound='MOLE_WEIGHTED', kappa='MOLE_WEIGHTED', isobaric_expansion='MOLE_WEIGHTED', Joule_Thomson='MOLE_WEIGHTED', VL_ID='PIP', S_ID='d2P_dVdT', solid_sort_method='prop', liquid_sort_method='prop', liquid_sort_cmps=[], solid_sort_cmps=[], liquid_sort_cmps_neg=[], solid_sort_cmps_neg=[], liquid_sort_prop='DENSITY_MASS', solid_sort_prop='DENSITY_MASS', phase_sort_higher_first=True, water_sort='water not special', equilibrium_perturbation=1e-07))[source]

Bases: thermo.flash.flash_vl.FlashVL

Class for performing flash calculations on multiphase vapor-liquid systems. This rigorous class does not make any assumptions and will search for up to the maximum amount of liquid phases specified by the user. Vapor and each liquid phase do not need to use a consistent thermodynamic model.

The minimum information that is needed in addition to the Phase objects is:

  • MWs

  • Vapor pressure curve

  • Functioning enthalpy models for each phase

Parameters
constantsChemicalConstantsPackage object

Package of chemical constants; these are used as boundaries at times, initial guesses other times, and in all cases these properties are accessible as attributes of the resulting EquilibriumState object, [-]

correlationsPropertyCorrelationsPackage

Package of chemical T-dependent properties; these are used as boundaries at times, for initial guesses other times, and in all cases these properties are accessible as attributes of the resulting EquilibriumState object, [-]

gasPhase object

A single phase which can represent the gas phase, [-]

liquidslist[Phase]

A list of phase objects that can represent the liquid phases; if working with a VLL system with a consistent model, specify the same liquid phase twice; the length of this list is the maximum number of liquid phases that will be searched for, [-]

solidslist[Phase]

Not used, [-]

settingsBulkSettings object

Object containing settings for calculating bulk and transport properties, [-]

Notes

The algorithms in this object are mostly from [1], [2] and [3]. Sequential substitution without acceleration is used by default to converge multiphase systems.

Additional information that can be provided in the ChemicalConstantsPackage object and PropertyCorrelationsPackage object that may help convergence is:

  • Tc, Pc, omega, Tb, and atoms

  • Gas heat capacity correlations

  • Liquid molar volume correlations

  • Heat of vaporization correlations

References

1

Michelsen, Michael L., and Jørgen M. Mollerup. Thermodynamic Models: Fundamentals & Computational Aspects. Tie-Line Publications, 2007.

2

Poling, Bruce E., John M. Prausnitz, and John P. O`Connell. The Properties of Gases and Liquids. 5th edition. New York: McGraw-Hill Professional, 2000.

3

Gmehling, Jürgen, Michael Kleiber, Bärbel Kolbe, and Jürgen Rarey. Chemical Thermodynamics for Process Simulation. John Wiley & Sons, 2019.

Examples

A three-phase flash of butanol, water, and ethanol with the SRK EOS without BIPs:

>>> from thermo import ChemicalConstantsPackage, CEOSGas, CEOSLiquid, SRKMIX, FlashVLN, PropertyCorrelationsPackage, HeatCapacityGas
>>> constants = ChemicalConstantsPackage(Tcs=[563.0, 647.14, 514.0], Pcs=[4414000.0, 22048320.0, 6137000.0], omegas=[0.59, 0.344, 0.635], MWs=[74.1216, 18.01528, 46.06844], CASs=['71-36-3', '7732-18-5', '64-17-5'])
>>> properties = PropertyCorrelationsPackage(constants=constants,
...                                     HeatCapacityGases=[HeatCapacityGas(poly_fit=(50.0, 1000.0, [-3.787200194613107e-20, 1.7692887427654656e-16, -3.445247207129205e-13, 3.612771874320634e-10, -2.1953250181084466e-07, 7.707135849197655e-05, -0.014658388538054169, 1.5642629364740657, -7.614560475001724])),
...                                     HeatCapacityGas(poly_fit=(50.0, 1000.0, [5.543665000518528e-22, -2.403756749600872e-18, 4.2166477594350336e-15, -3.7965208514613565e-12, 1.823547122838406e-09, -4.3747690853614695e-07, 5.437938301211039e-05, -0.003220061088723078, 33.32731489750759])),
...                                     HeatCapacityGas(poly_fit=(50.0, 1000.0, [-1.162767978165682e-20, 5.4975285700787494e-17, -1.0861242757337942e-13, 1.1582703354362728e-10, -7.160627710867427e-08, 2.5392014654765875e-05, -0.004732593693568646, 0.5072291035198603, 20.037826650765965])),], )
>>> eos_kwargs = dict(Tcs=constants.Tcs, Pcs=constants.Pcs, omegas=constants.omegas)
>>> gas = CEOSGas(SRKMIX, eos_kwargs, HeatCapacityGases=properties.HeatCapacityGases)
>>> liq = CEOSLiquid(SRKMIX, eos_kwargs, HeatCapacityGases=properties.HeatCapacityGases)
>>> flashN = FlashVLN(constants, properties, liquids=[liq, liq], gas=gas)
>>> res = flashN.flash(T=361, P=1e5, zs=[.25, 0.7, .05])
>>> res.phase_count
3
Attributes
SS_NP_MAXITERint

Maximum number of sequential substitution iterations to try when converging a three or more phase solution, [-]

SS_NP_TOLfloat

Convergence tolerance in sequential substitution for a three or more phase solution [-]

SS_NP_TRIVIAL_TOLfloat

Tolerance at which to quick a three-phase flash because it is converging to the trivial solution, [-]

SS_STAB_AQUEOUS_CHECKbool

If True, the first three-phase stability check will be on water (if it is present) as it forms a three-phase solution more than any other component, [-]

DOUBLE_CHECK_2Pbool

This parameter should be set to True if any issues in the solution are noticed. It can slow down two-phase solution. It ensures that all potential vapor-liquid and liquid-liquid phase pairs are searched for stability, instead of testing first for a vapor-liquid solution and then moving on to a three phase flash if an instability is detected, [-]

Base Flash Class

class thermo.flash.Flash[source]

Bases: object

Base class for performing flash calculations. All Flash objects need to inherit from this, and common methods can be added to it.

Methods

flash([zs, T, P, VF, SF, V, H, S, G, U, A, ...])

Method to perform a flash calculation and return the result as an EquilibriumState object.

plot_PT(zs[, Pmin, Pmax, pts, branches, ...])

Method to create a plot of the phase envelope as can be calculated from a series of pressure & vapor fraction spec flashes.

plot_Pxy(T[, pts, ignore_errors, values, show])

Method to create a Pxy plot for a binary system (holding temperature constant); the mole fraction of the first species is varied.

plot_TP(zs[, Tmin, Tmax, pts, branches, ...])

Method to create a plot of the phase envelope as can be calculated from a series of temperature & vapor fraction spec flashes.

plot_Txy(P[, pts, ignore_errors, values, show])

Method to create a Txy plot for a binary system (holding pressure constant); the mole fraction of the first species is varied.

plot_ternary([T, P, scale])

Method to create a ternary plot of the system at either a specified temperature or pressure.

plot_xy([P, T, pts, ignore_errors, values, ...])

Method to create a xy diagram for a binary system.

flash(zs=None, T=None, P=None, VF=None, SF=None, V=None, H=None, S=None, G=None, U=None, A=None, solution=None, hot_start=None, retry=False, dest=None, rho=None, rho_mass=None, H_mass=None, S_mass=None, G_mass=None, U_mass=None, A_mass=None, spec_fun=None, H_reactive=None)[source]

Method to perform a flash calculation and return the result as an EquilibriumState object. This generic interface allows flashes with any combination of valid specifications; if a flash is unimplemented and error will be raised.

Parameters
zslist[float], optional

Mole fractions of each component, required unless there is only one component, [-]

Tfloat, optional

Temperature, [K]

Pfloat, optional

Pressure, [Pa]

VFfloat, optional

Vapor fraction, [-]

SFfloat, optional

Solid fraction, [-]

Vfloat, optional

Molar volume of the overall bulk, [m^3/mol]

Hfloat, optional

Molar enthalpy of the overall bulk, [J/mol]

Sfloat, optional

Molar entropy of the overall bulk, [J/(mol*K)]

Gfloat, optional

Molar Gibbs free energy of the overall bulk, [J/mol]

Ufloat, optional

Molar internal energy of the overall bulk, [J/mol]

Afloat, optional

Molar Helmholtz energy of the overall bulk, [J/mol]

solutionstr or int, optional

When multiple solutions exist, if more than one is found they will be sorted by T (and then P) increasingly; this number will index into the multiple solution array. Negative indexing is supported. ‘high’ is an alias for 0, and ‘low’ an alias for -1. Setting this parameter may make a flash slower because in some cases more checks are performed. [-]

hot_startEquilibriumState

A previously converged flash or initial guessed state from which the flash can begin; this parameter can save time in some cases, [-]

retrybool

Usually for flashes like UV or PH, there are multiple sets of possible iteration variables. For the UV case, the prefered iteration variable is P, so each iteration a PV solve is done on the phase; but equally the flash can be done iterating on T, where a TV solve is done on the phase each iteration. Depending on the tolerances, the flash type, the thermodynamic consistency of the phase, and other factors, it is possible the flash can fail. If retry is set to True, the alternate variable set will be iterated as a backup if the first flash fails. [-]

destNone or EquilibriumState or EquilibriumStream

What type of object the flash result is set into; leave as None to obtain the normal EquilibriumState results, [-]

rhofloat, optional

Molar density of the overall bulk; this is trivially converted to a V spec, [mol/m^3]

rho_massfloat, optional

Mass density of the overall bulk; this is trivially converted to a rho spec, [kg/m^3]

H_massfloat, optional

Mass enthalpy of the overall bulk; this is trivially converted to a H spec, [J/kg]

S_massfloat, optional

Mass entropy of the overall bulk; this is trivially converted to a S spec, [J/(kg*K)]

G_massfloat, optional

Mass Gibbs free energy of the overall bulk; this is trivially converted to a G spec, [J/kg]

U_massfloat, optional

Mass internal energy of the overall bulk; this is trivially converted to a U spec, [J/kg]

A_massfloat, optional

Mass Helmholtz energy of the overall bulk; this is trivially converted to a A spec, [J/kg]

H_reactivefloat, optional

Molar reactive enthalpy, [J/mol]

Returns
resultsEquilibriumState

Equilibrium object containing the state of the phases after the flash calculation [-]

Notes

Warning

Not all flash specifications have a unique solution. Not all flash specifications will converge, whether from a bad model, bad inputs, or simply a lack of convergence by the implemented algorithms. You are welcome to submit these cases to the author but the library is provided AS IS, with NO SUPPORT.

Warning

Convergence of a flash may be impaired by providing hot_start. If reliability is desired, do not use this parameter.

Warning

The most likely thermodynamic methods to converge are thermodynamically consistent ones. This means e.g. an ideal liquid and an ideal gas; or an equation of state for both phases. Mixing thermodynamic models increases the possibility of multiple solutions, discontinuities, and other not-fun issues for the algorithms.

plot_PT(zs, Pmin=None, Pmax=None, pts=50, branches=[], ignore_errors=True, values=False, show=True, hot=False)[source]

Method to create a plot of the phase envelope as can be calculated from a series of pressure & vapor fraction spec flashes. By default vapor fractions of 0 and 1 are plotted; additional vapor fraction specifications can be specified in the branches argument as a list.

Parameters
zslist[float]

Mole fractions of the feed, [-]

Pminfloat, optional

Minimum pressure to begin the plot, [Pa]

Pmaxfloat, optional

Maximum pressure to end the plot, [Pa]

ptsint, optional

The number of points to calculated for each vapor fraction value, [-]

brancheslist[float], optional

Extra vapor fraction values to plot, [-]

ignore_errorsbool, optional

Whether to fail on a calculation failure, or to ignore the bad point, [-]

valuesbool, optional

If True, the calculated values will be returned instead of plotted, [-]

showbool, optional

If True, the plot will be created and displayed; if False and values is False, the plot Figure object will be returned but not displayed; and if False and values is False, no plot will be created or shown [-]

hotbool, optional

Whether to restart the next flash from the previous flash or not (intended to speed the call when True), [-]

Returns
Pslist[float]

Pressures, [Pa]

T_dewslist[float]

Bubble point temperatures at the evaluated points, [K]

T_bubbleslist[float]

Dew point temperatures at the evaluated points, [K]

branch_TsNone or list[list[float]]

Temperatures which yield the equilibrium vapor fractions specified; formatted as [[T1_VFx, T2_VFx, … Tn_VFx], …, [T1_VFy, T2_VFy, … Tn_VFy]], [k]

plot_Pxy(T, pts=30, ignore_errors=True, values=False, show=True)[source]

Method to create a Pxy plot for a binary system (holding temperature constant); the mole fraction of the first species is varied.

Parameters
Tfloat

Temperature for the plot, [K]

ptsint, optional

The number of points to calculated for each vapor fraction value, [-]

ignore_errorsbool, optional

Whether to fail on a calculation failure, or to ignore the bad point, [-]

valuesbool, optional

If True, the calculated values will be returned instead of plotted, [-]

showbool, optional

If True, the plot will be created and displayed; if False and values is False, the plot Figure object will be returned but not displayed; and if False and values is False, no plot will be created or shown [-]

Returns
figFigure

Plot object, [-]

z1list[float]

Mole fractions of the first component at each point, [-]

z2list[float]

Mole fractions of the second component at each point, [-]

P_dewslist[float]

Bubble point pressures at the evaluated points, [Pa]

P_bubbleslist[float]

Dew point pressures at the evaluated points, [Pa]

plot_TP(zs, Tmin=None, Tmax=None, pts=50, branches=None, ignore_errors=True, values=False, show=True, hot=False)[source]

Method to create a plot of the phase envelope as can be calculated from a series of temperature & vapor fraction spec flashes. By default vapor fractions of 0 and 1 are plotted; additional vapor fraction specifications can be specified in the branches argument as a list.

Parameters
zslist[float]

Mole fractions of the feed, [-]

Tminfloat, optional

Minimum temperature to begin the plot, [K]

Tmaxfloat, optional

Maximum temperature to end the plot, [K]

ptsint, optional

The number of points to calculated for each vapor fraction value, [-]

brancheslist[float], optional

Extra vapor fraction values to plot, [-]

ignore_errorsbool, optional

Whether to fail on a calculation failure, or to ignore the bad point, [-]

valuesbool, optional

If True, the calculated values will be returned instead of plotted, [-]

showbool, optional

If True, the plot will be created and displayed; if False and values is False, the plot Figure object will be returned but not displayed; and if False and values is False, no plot will be created or shown [-]

hotbool, optional

Whether to restart the next flash from the previous flash or not (intended to speed the call when True), [-]

Returns
Tslist[float]

Temperatures, [K]

P_dewslist[float]

Bubble point pressures at the evaluated points, [Pa]

P_bubbleslist[float]

Dew point pressures at the evaluated points, [Pa]

branch_PsNone or list[list[float]]

Pressures which yield the equilibrium vapor fractions specified; formatted as [[P1_VFx, P2_VFx, … Pn_VFx], …, [P1_VFy, P2_VFy, … Pn_VFy]], [Pa]

plot_Txy(P, pts=30, ignore_errors=True, values=False, show=True)[source]

Method to create a Txy plot for a binary system (holding pressure constant); the mole fraction of the first species is varied.

Parameters
Pfloat

Pressure for the plot, [Pa]

ptsint, optional

The number of points to calculated for each vapor fraction value, [-]

ignore_errorsbool, optional

Whether to fail on a calculation failure, or to ignore the bad point, [-]

valuesbool, optional

If True, the calculated values will be returned instead of plotted, [-]

showbool, optional

If True, the plot will be created and displayed; if False and values is False, the plot Figure object will be returned but not displayed; and if False and values is False, no plot will be created or shown [-]

Returns
figFigure

Plot object, [-]

z1list[float]

Mole fractions of the first component at each point, [-]

z2list[float]

Mole fractions of the second component at each point, [-]

T_dewslist[float]

Bubble point temperatures at the evaluated points, [K]

T_bubbleslist[float]

Dew point temperatures at the evaluated points, [K]

plot_ternary(T=None, P=None, scale=10)[source]

Method to create a ternary plot of the system at either a specified temperature or pressure.

Parameters
Tfloat, optional

Temperature, [K]

Pfloat, optional

Pressure, [Pa]

plot_xy(P=None, T=None, pts=30, ignore_errors=True, values=False, show=True, VF=0.0)[source]

Method to create a xy diagram for a binary system. Either a temperature or pressure can be specified. By default, bubble point flashes are performed; this can be varied by changing VF.

Parameters
Pfloat, optional

The specified pressure, [Pa]

Tfloat, optional

The specified temperature, [K]

ptsint, optional

The number of points in the plot [-]

ignore_errorsbool, optional

Whether to fail on a calculation failure, or to ignore the bad point, [-]

valuesbool, optional

If True, the calculated values will be returned instead of plotted, [-]

showbool, optional

If True, the plot will be created and displayed; if False and values is False, the plot Figure object will be returned but not displayed; and if False and values is False, no plot will be created or shown [-]

Returns
figFigure

Plot object, [-]

z1list[float]

Overall mole fractions of the first component at each point, [-]

z2list[float]

Overall mole fractions of the second component at each point, [-]

x1list[float]

Liquid mole fractions of component 1 at each point, [-]

y1list[float]

Vapor mole fractions of component 1 at each point, [-]

Specific Flash Algorithms

It is recommended to use the Flash classes, which are designed to have generic interfaces. The implemented specific flash algorithms may be changed in the future, but reading their source code may be helpful for instructive purposes.