Bulk Phases (thermo.bulk)

This module contains a phase wrapper for obtaining properties of a pseudo-phase made of multiple other phases. This is useful in the context of multiple liquid phases; or multiple solid phases; or looking at all the phases together.

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

Bulk Class

class thermo.bulk.Bulk(T, P, zs, phases, phase_fractions, phase_bulk=None)[source]

Bases: thermo.phases.phase.Phase

Class to encapsulate multiple Phase objects and provide a unified interface for obtaining properties from a group of phases.

This class exists for three purposes:

  • Providing a common interface for obtaining properties like Cp - whether there is one phase or 100, calling Cp on the bulk will retrieve that value.

  • Retrieving “bulk” properties that do make sense to be calculated for a combination of phases together.

  • Allowing configurable estimations of non-bulk properties like isothermal compressibility or speed of sound for the group of phases together.

Parameters
Tfloat

Temperature of the bulk, [K]

Pfloat

Pressure of the bulk, [Pa]

zslist[float]

Mole fractions of the bulk, [-]

phaseslist[Phase]

Phase objects, [-]

phase_fractionslist[float]

Molar fractions of each phase, [-]

phase_bulkstr, optional

None to represent a bulk of all present phases; ‘l’ to represent a bulk of only liquid phases; s to represent a bulk of only solid phases, [-]

Notes

Please think carefully when retrieving a property of the bulk. If there are two liquid phases in a bulk, and a single viscosity value is retrieved, can that be used directly for a single phase pressure drop calculation? Not with any theoretical consistency, that’s for sure.

Attributes
beta

Phase fraction of the bulk phase.

beta_volume_liquid_ref

Method to return the standard liquid volume fraction of this phase.

betas_mass

Method to calculate and return the mass fraction of all of the phases in the bulk.

betas_volume

Method to calculate and return the volume fraction of all of the phases in the bulk.

Methods

Cp()

Method to calculate and return the constant-temperature and constant phase-fraction heat capacity of the bulk phase.

Cp_ideal_gas()

Method to calculate and return the ideal-gas heat capacity of the phase.

H()

Method to calculate and return the constant-temperature and constant phase-fraction enthalpy of the bulk phase.

H_ideal_gas()

Method to calculate and return the ideal-gas enthalpy of the phase.

H_reactive()

Method to calculate and return the constant-temperature and constant phase-fraction reactive enthalpy of the bulk phase.

Joule_Thomson()

Method to calculate and return the Joule-Thomson coefficient of the bulk according to the selected calculation methodology.

MW()

Method to calculate and return the molecular weight of the bulk phase.

Pmc()

Method to calculate and return the mechanical critical pressure of the phase.

S()

Method to calculate and return the constant-temperature and constant phase-fraction entropy of the bulk phase.

S_ideal_gas()

Method to calculate and return the ideal-gas entropy of the phase.

S_reactive()

Method to calculate and return the constant-temperature and constant phase-fraction reactive entropy of the bulk phase.

Tmc()

Method to calculate and return the mechanical critical temperature of the phase.

V()

Method to calculate and return the molar volume of the bulk phase.

V_iter([force])

Method to calculate and return the molar volume of the bulk phase, with precision suitable for a TV calculation to calculate a matching pressure.

Vmc()

Method to calculate and return the mechanical critical volume of the phase.

Zmc()

Method to calculate and return the mechanical critical compressibility of the phase.

d2P_dT2()

Method to calculate and return the second temperature derivative of pressure of the bulk according to the selected calculation methodology.

d2P_dT2_frozen()

Method to calculate and return the second constant-volume derivative of pressure with respect to temperature of the bulk phase, at constant phase fractions and phase compositions.

d2P_dTdV()

Method to calculate and return the second derivative of pressure with respect to temperature and volume of the bulk according to the selected calculation methodology.

d2P_dTdV_frozen()

Method to calculate and return the second derivative of pressure with respect to volume and temperature of the bulk phase, at constant phase fractions and phase compositions.

d2P_dV2()

Method to calculate and return the second volume derivative of pressure of the bulk according to the selected calculation methodology.

d2P_dV2_frozen()

Method to calculate and return the constant-temperature second derivative of pressure with respect to volume of the bulk phase, at constant phase fractions and phase compositions.

dA_dP()

Method to calculate and return the constant-temperature pressure derivative of Helmholtz energy.

dA_dT()

Method to calculate and return the constant-pressure temperature derivative of Helmholtz energy.

dG_dP()

Method to calculate and return the constant-temperature pressure derivative of Gibbs free energy.

dG_dT()

Method to calculate and return the constant-pressure temperature derivative of Gibbs free energy.

dP_dT()

Method to calculate and return the first temperature derivative of pressure of the bulk according to the selected calculation methodology.

dP_dT_frozen()

Method to calculate and return the constant-volume derivative of pressure with respect to temperature of the bulk phase, at constant phase fractions and phase compositions.

dP_dV()

Method to calculate and return the first volume derivative of pressure of the bulk according to the selected calculation methodology.

dP_dV_frozen()

Method to calculate and return the constant-temperature derivative of pressure with respect to volume of the bulk phase, at constant phase fractions and phase compositions.

dU_dP()

Method to calculate and return the constant-temperature pressure derivative of internal energy.

dU_dT()

Method to calculate and return the constant-pressure temperature derivative of internal energy.

isobaric_expansion()

Method to calculate and return the isobatic expansion coefficient of the bulk according to the selected calculation methodology.

k()

Calculate and return the thermal conductivity of the bulk according to the selected thermal conductivity settings in BulkSettings, the settings in ThermalConductivityGasMixture and ThermalConductivityLiquidMixture, and the configured pure-component settings in ThermalConductivityGas and ThermalConductivityLiquid.

kappa()

Method to calculate and return the isothermal compressibility of the bulk according to the selected calculation methodology.

mu()

Calculate and return the viscosity of the bulk according to the selected viscosity settings in BulkSettings, the settings in ViscosityGasMixture and ViscosityLiquidMixture, and the configured pure-component settings in ViscosityGas and ViscosityLiquid.

sigma()

Calculate and return the surface tension of the bulk according to the selected surface tension settings in BulkSettings, the settings in SurfaceTensionMixture and the configured pure-component settings in SurfaceTension.

speed_of_sound()

Method to calculate and return the molar speed of sound of the bulk according to the selected calculation methodology.

Cp()[source]

Method to calculate and return the constant-temperature and constant phase-fraction heat capacity of the bulk phase. This is a phase-fraction weighted calculation.

Cp=ipCp,iβiC_p = \sum_i^p C_{p,i} \beta_i
Returns
Cpfloat

Molar heat capacity, [J/(mol*K)]

Cp_ideal_gas()[source]

Method to calculate and return the ideal-gas heat capacity of the phase.

Cpig=iziCp,iigC_p^{ig} = \sum_i z_i {C_{p,i}^{ig}}
Returns
Cpfloat

Ideal gas heat capacity, [J/(mol*K)]

H()[source]

Method to calculate and return the constant-temperature and constant phase-fraction enthalpy of the bulk phase. This is a phase-fraction weighted calculation.

H=ipHiβiH = \sum_i^p H_{i} \beta_i
Returns
Hfloat

Molar enthalpy, [J/(mol)]

H_ideal_gas()[source]

Method to calculate and return the ideal-gas enthalpy of the phase.

Hig=iziHiigH^{ig} = \sum_i z_i {H_{i}^{ig}}
Returns
Hfloat

Ideal gas enthalpy, [J/(mol)]

H_reactive()[source]

Method to calculate and return the constant-temperature and constant phase-fraction reactive enthalpy of the bulk phase. This is a phase-fraction weighted calculation.

Hreactive=ipHreactive,iβiH_{\text{reactive}} = \sum_i^p H_{\text{reactive}, i} \beta_i
Returns
H_reactivefloat

Reactive molar enthalpy, [J/(mol)]

Joule_Thomson()[source]

Method to calculate and return the Joule-Thomson coefficient of the bulk according to the selected calculation methodology.

μJT=(TP)H\mu_{JT} = \left(\frac{\partial T}{\partial P}\right)_H
Returns
mu_JTfloat

Joule-Thomson coefficient [K/Pa]

MW()[source]

Method to calculate and return the molecular weight of the bulk phase. This is a phase-fraction weighted calculation.

MW=ipMWiβi\text{MW} = \sum_i^p \text{MW}_i \beta_i
Returns
MWfloat

Molecular weight, [g/mol]

Pmc()[source]

Method to calculate and return the mechanical critical pressure of the phase.

Returns
Pmcfloat

Mechanical critical pressure, [Pa]

S()[source]

Method to calculate and return the constant-temperature and constant phase-fraction entropy of the bulk phase. This is a phase-fraction weighted calculation.

S=ipSiβiS = \sum_i^p S_{i} \beta_i
Returns
Sfloat

Molar entropy, [J/(mol*K)]

S_ideal_gas()[source]

Method to calculate and return the ideal-gas entropy of the phase.

Sig=iziSiigRln(PPref)Riziln(zi)S^{ig} = \sum_i z_i S_{i}^{ig} - R\ln\left(\frac{P}{P_{ref}}\right) - R\sum_i z_i \ln(z_i)
Returns
Sfloat

Ideal gas molar entropy, [J/(mol*K)]

S_reactive()[source]

Method to calculate and return the constant-temperature and constant phase-fraction reactive entropy of the bulk phase. This is a phase-fraction weighted calculation.

Sreactive=ipSreactive,iβiS_{\text{reactive}} = \sum_i^p S_{\text{reactive}, i} \beta_i
Returns
S_reactivefloat

Reactive molar entropy, [J/(mol*K)]

Tmc()[source]

Method to calculate and return the mechanical critical temperature of the phase.

Returns
Tmcfloat

Mechanical critical temperature, [K]

V()[source]

Method to calculate and return the molar volume of the bulk phase. This is a phase-fraction weighted calculation.

V=ipViβiV = \sum_i^p V_i \beta_i
Returns
Vfloat

Molar volume, [m^3/mol]

V_iter(force=False)[source]

Method to calculate and return the molar volume of the bulk phase, with precision suitable for a TV calculation to calculate a matching pressure. This is a phase-fraction weighted calculation.

V=ipViβiV = \sum_i^p V_i \beta_i
Returns
Vfloat or mpf

Molar volume, [m^3/mol]

Vmc()[source]

Method to calculate and return the mechanical critical volume of the phase.

Returns
Vmcfloat

Mechanical critical volume, [m^3/mol]

Zmc()[source]

Method to calculate and return the mechanical critical compressibility of the phase.

Returns
Zmcfloat

Mechanical critical compressibility, [-]

property beta

Phase fraction of the bulk phase. Should always be 1 when representing all phases of a flash; but can be less than one if representing multiple solids or liquids as a single phase in a larger mixture.

Returns
betafloat

Phase fraction of bulk, [-]

property beta_volume_liquid_ref

Method to return the standard liquid volume fraction of this phase. This method is only available when the phase is linked to an EquilibriumState.

Returns
beta_volumefloat

Phase fraction on a volumetric basis, [-]

property betas_mass

Method to calculate and return the mass fraction of all of the phases in the bulk.

Returns
betas_masslist[float]

Mass phase fractions of all the phases in the bulk object, ordered vapor, liquid, then solid, [-]

property betas_volume

Method to calculate and return the volume fraction of all of the phases in the bulk.

Returns
betas_volumelist[float]

Volume phase fractions of all the phases in the bulk, ordered vapor, liquid, then solid , [-]

d2P_dT2()[source]

Method to calculate and return the second temperature derivative of pressure of the bulk according to the selected calculation methodology.

Returns
d2P_dT2float

Second temperature derivative of pressure, [Pa/K^2]

d2P_dT2_frozen()[source]

Method to calculate and return the second constant-volume derivative of pressure with respect to temperature of the bulk phase, at constant phase fractions and phase compositions. This is a molar phase-fraction weighted calculation.

(2PT2)V,β,zs=iphasesβi(2PT2)i,Vi,βi,zsi\left(\frac{\partial^2 P}{\partial T^2}\right)_{V, \beta, {zs}} = \sum_{i}^{\text{phases}} \beta_i \left(\frac{\partial^2 P} {\partial T^2}\right)_{i, V_i, \beta_i, {zs}_i}
Returns
d2P_dT2_frozenfloat

Frozen constant-volume second derivative of pressure with respect to temperature of the bulk phase, [Pa/K^2]

d2P_dTdV()[source]

Method to calculate and return the second derivative of pressure with respect to temperature and volume of the bulk according to the selected calculation methodology.

Returns
d2P_dTdVfloat

Second volume derivative of pressure, [mol*Pa^2/(J*K)]

d2P_dTdV_frozen()[source]

Method to calculate and return the second derivative of pressure with respect to volume and temperature of the bulk phase, at constant phase fractions and phase compositions. This is a molar phase-fraction weighted calculation.

(2PVT)β,zs=iphasesβi(2PVT)i,βi,zsi\left(\frac{\partial^2 P}{\partial V \partial T}\right)_{\beta, {zs}} = \sum_{i}^{\text{phases}} \beta_i \left(\frac{\partial^2 P} {\partial V \partial T}\right)_{i, \beta_i, {zs}_i}
Returns
d2P_dTdV_frozenfloat

Frozen second derivative of pressure with respect to volume and temperature of the bulk phase, [Pa*mol^2/m^6]

d2P_dV2()[source]

Method to calculate and return the second volume derivative of pressure of the bulk according to the selected calculation methodology.

Returns
d2P_dV2float

Second volume derivative of pressure, [Pa*mol^2/m^6]

d2P_dV2_frozen()[source]

Method to calculate and return the constant-temperature second derivative of pressure with respect to volume of the bulk phase, at constant phase fractions and phase compositions. This is a molar phase-fraction weighted calculation.

(2PV2)T,β,zs=iphasesβi(2PV2)i,T,βi,zsi\left(\frac{\partial^2 P}{\partial V^2}\right)_{T, \beta, {zs}} = \sum_{i}^{\text{phases}} \beta_i \left(\frac{\partial^2 P} {\partial V^2}\right)_{i, T, \beta_i, {zs}_i}
Returns
d2P_dV2_frozenfloat

Frozen constant-temperature second derivative of pressure with respect to volume of the bulk phase, [Pa*mol^2/m^6]

dA_dP()[source]

Method to calculate and return the constant-temperature pressure derivative of Helmholtz energy.

(AP)T=T(SP)T+(UP)T\left(\frac{\partial A}{\partial P}\right)_{T} = -T \left(\frac{\partial S}{\partial P}\right)_{T} + \left(\frac{\partial U}{\partial P}\right)_{T}
Returns
dA_dPfloat

Constant-temperature pressure derivative of Helmholtz energy, [J/(mol*Pa)]

dA_dT()[source]

Method to calculate and return the constant-pressure temperature derivative of Helmholtz energy.

(AT)P=T(ST)PS+(UT)P\left(\frac{\partial A}{\partial T}\right)_{P} = -T \left(\frac{\partial S}{\partial T}\right)_{P} - S + \left(\frac{\partial U}{\partial T}\right)_{P}
Returns
dA_dTfloat

Constant-pressure temperature derivative of Helmholtz energy, [J/(mol*K)]

dG_dP()[source]

Method to calculate and return the constant-temperature pressure derivative of Gibbs free energy.

(GP)T=T(SP)T+(HP)T\left(\frac{\partial G}{\partial P}\right)_{T} = -T\left(\frac{\partial S}{\partial P}\right)_{T} + \left(\frac{\partial H}{\partial P}\right)_{T}
Returns
dG_dPfloat

Constant-temperature pressure derivative of Gibbs free energy, [J/(mol*Pa)]

dG_dT()[source]

Method to calculate and return the constant-pressure temperature derivative of Gibbs free energy.

(GT)P=T(ST)PS+(HT)P\left(\frac{\partial G}{\partial T}\right)_{P} = -T\left(\frac{\partial S}{\partial T}\right)_{P} - S + \left(\frac{\partial H}{\partial T}\right)_{P}
Returns
dG_dTfloat

Constant-pressure temperature derivative of Gibbs free energy, [J/(mol*K)]

dP_dT()[source]

Method to calculate and return the first temperature derivative of pressure of the bulk according to the selected calculation methodology.

Returns
dP_dTfloat

First temperature derivative of pressure, [Pa/K]

dP_dT_frozen()[source]

Method to calculate and return the constant-volume derivative of pressure with respect to temperature of the bulk phase, at constant phase fractions and phase compositions. This is a molar phase-fraction weighted calculation.

(PT)V,β,zs=iphasesβi(PT)i,Vi,βi,zsi\left(\frac{\partial P}{\partial T}\right)_{V, \beta, {zs}} = \sum_{i}^{\text{phases}} \beta_i \left(\frac{\partial P} {\partial T}\right)_{i, V_i, \beta_i, {zs}_i}
Returns
dP_dT_frozenfloat

Frozen constant-volume derivative of pressure with respect to temperature of the bulk phase, [Pa/K]

dP_dV()[source]

Method to calculate and return the first volume derivative of pressure of the bulk according to the selected calculation methodology.

Returns
dP_dVfloat

First volume derivative of pressure, [Pa*mol/m^3]

dP_dV_frozen()[source]

Method to calculate and return the constant-temperature derivative of pressure with respect to volume of the bulk phase, at constant phase fractions and phase compositions. This is a molar phase-fraction weighted calculation.

(PV)T,β,zs=iphasesβi(PV)i,T,βi,zsi\left(\frac{\partial P}{\partial V}\right)_{T, \beta, {zs}} = \sum_{i}^{\text{phases}} \beta_i \left(\frac{\partial P} {\partial V}\right)_{i, T, \beta_i, {zs}_i}
Returns
dP_dV_frozenfloat

Frozen constant-temperature derivative of pressure with respect to volume of the bulk phase, [Pa*mol/m^3]

dU_dP()[source]

Method to calculate and return the constant-temperature pressure derivative of internal energy.

(UP)T=P(VP)TV+(HP)T\left(\frac{\partial U}{\partial P}\right)_{T} = -P \left(\frac{\partial V}{\partial P}\right)_{T} - V + \left(\frac{\partial H}{\partial P}\right)_{T}
Returns
dU_dPfloat

Constant-temperature pressure derivative of internal energy, [J/(mol*Pa)]

dU_dT()[source]

Method to calculate and return the constant-pressure temperature derivative of internal energy.

(UT)P=P(VT)P+(HT)P\left(\frac{\partial U}{\partial T}\right)_{P} = -P \left(\frac{\partial V}{\partial T}\right)_{P} + \left(\frac{\partial H}{\partial T}\right)_{P}
Returns
dU_dTfloat

Constant-pressure temperature derivative of internal energy, [J/(mol*K)]

isobaric_expansion()[source]

Method to calculate and return the isobatic expansion coefficient of the bulk according to the selected calculation methodology.

β=1V(VT)P\beta = \frac{1}{V}\left(\frac{\partial V}{\partial T} \right)_P
Returns
betafloat

Isobaric coefficient of a thermal expansion, [1/K]

k()[source]

Calculate and return the thermal conductivity of the bulk according to the selected thermal conductivity settings in BulkSettings, the settings in ThermalConductivityGasMixture and ThermalConductivityLiquidMixture, and the configured pure-component settings in ThermalConductivityGas and ThermalConductivityLiquid.

Returns
kfloat

Thermal Conductivity of bulk phase calculated with mixing rules, [Pa*s]

kappa()[source]

Method to calculate and return the isothermal compressibility of the bulk according to the selected calculation methodology.

κ=1V(VP)T\kappa = -\frac{1}{V}\left(\frac{\partial V}{\partial P} \right)_T
Returns
kappafloat

Isothermal coefficient of compressibility, [1/Pa]

mu()[source]

Calculate and return the viscosity of the bulk according to the selected viscosity settings in BulkSettings, the settings in ViscosityGasMixture and ViscosityLiquidMixture, and the configured pure-component settings in ViscosityGas and ViscosityLiquid.

Returns
mufloat

Viscosity of bulk phase calculated with mixing rules, [Pa*s]

sigma()[source]

Calculate and return the surface tension of the bulk according to the selected surface tension settings in BulkSettings, the settings in SurfaceTensionMixture and the configured pure-component settings in SurfaceTension.

Returns
sigmafloat

Surface tension of bulk phase calculated with mixing rules, [N/m]

Notes

A value is only returned if all phases in the bulk are liquids; this property is for a liquid-ideal gas calculation, not the interfacial tension between two liquid phases.

speed_of_sound()[source]

Method to calculate and return the molar speed of sound of the bulk according to the selected calculation methodology.

w=[V2(PV)TCpCv]1/2w = \left[-V^2 \left(\frac{\partial P}{\partial V}\right)_T \frac{C_p} {C_v}\right]^{1/2}

A similar expression based on molar density is:

w=[(Pρ)TCpCv]1/2w = \left[\left(\frac{\partial P}{\partial \rho}\right)_T \frac{C_p} {C_v}\right]^{1/2}
Returns
wfloat

Speed of sound for a real gas, [m*kg^0.5/(s*mol^0.5)]

Bulk Settings Class

class thermo.bulk.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', VL_ID_settings=None, S_ID='d2P_dVdT', S_ID_settings=None, 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: object

Class containing configuration methods for determining how properties of a Bulk phase made of different phases are handled. All parameters are also attributes.

Parameters
dP_dTstr, optional

The method used to calculate the constant-volume temperature derivative of pressure of the bulk. One of DP_DT_METHODS, [-]

dP_dVstr, optional

The method used to calculate the constant-temperature volume derivative of pressure of the bulk. One of DP_DV_METHODS, [-]

d2P_dV2str, optional

The method used to calculate the second constant-temperature volume derivative of pressure of the bulk. One of D2P_DV2_METHODS, [-]

d2P_dT2str, optional

The method used to calculate the second constant-volume temperature derivative of pressure of the bulk. One of D2P_DT2_METHODS, [-]

d2P_dTdVstr, optional

The method used to calculate the temperature and volume derivative of pressure of the bulk. One of D2P_DTDV_METHODS, [-]

T_liquid_volume_reffloat, optional

Liquid molar volume reference temperature; if this is 298.15 K exactly, the molar volumes in Vml_STPs will be used, and if it is 288.7055555555555 K exactly, Vml_60Fs will be used, and otherwise the molar liquid volumes will be obtained from the temperature-dependent correlations specified, [K]

T_gas_reffloat, optional

Reference temperature to use for the calculation of ideal-gas molar volume and flow rate, [K]

P_gas_reffloat, optional

Reference pressure to use for the calculation of ideal-gas molar volume and flow rate, [Pa]

T_normalfloat, optional

“Normal” gas reference temperature for the calculation of ideal-gas molar volume in the “normal” reference state; default 273.15 K (0 C) according to [1], [K]

P_normalfloat, optional

“Normal” gas reference pressure for the calculation of ideal-gas molar volume in the “normal” reference state; default 101325 Pa (1 atm) according to [1], [Pa]

T_standardfloat, optional

“Standard” gas reference temperature for the calculation of ideal-gas molar volume in the “standard” reference state; default 288.15 K (15° C) according to [2]; 288.7055555555555 is also often used (60° F), [K]

P_standardfloat, optional

“Standard” gas reference pressure for the calculation of ideal-gas molar volume in the “standard” reference state; default 101325 Pa (1 atm) according to [2], [Pa]

mu_LLstr, optional

Mixing rule for multiple liquid phase liquid viscosity calculations; see MU_LL_METHODS for available options, [-]

mu_LL_power_exponentfloat, optional

Liquid-liquid viscosity power-law mixing parameter, used only when a power law mixing rule is selected, [-]

mu_VLstr, optional

Mixing rule for vapor-liquid viscosity calculations; see MU_VL_METHODS for available options, [-]

mu_VL_power_exponentfloat, optional

Vapor-liquid viscosity power-law mixing parameter, used only when a power law mixing rule is selected, [-]

k_LLstr, optional

Mixing rule for multiple liquid phase liquid thermal conductivity calculations; see K_LL_METHODS for available options, [-]

k_LL_power_exponentfloat, optional

Liquid-liquid thermal conductivity power-law mixing parameter, used only when a power law mixing rule is selected, [-]

k_VLstr, optional

Mixing rule for vapor-liquid thermal conductivity calculations; see K_VL_METHODS for available options, [-]

k_VL_power_exponentfloat, optional

Vapor-liquid thermal conductivity power-law mixing parameter, used only when a power law mixing rule is selected, [-]

sigma_LLstr, optional

Mixing rule for multiple liquid phase, air-liquid surface tension calculations; see SIGMA_LL_METHODS for available options, [-]

sigma_LL_power_exponentfloat, optional

Air-liquid Liquid-liquid surface tension power-law mixing parameter, used only when a power law mixing rule is selected, [-]

equilibrium_perturbationfloat, optional

The relative perturbation to use when calculating equilibrium derivatives numerically; for example if this is 1e-3 and T is the perturbation variable and the statis is 500 K, the perturbation calculation temperature will be 500.5 K, [various]

isobaric_expansionstr, optional

Mixing rule for multiphase isobaric expansion calculations; see BETA_METHODS for available options, [-]

speed_of_soundstr, optional

Mixing rule for multiphase speed of sound calculations; see SPEED_OF_SOUND_METHODS for available options, [-]

kappastr, optional

Mixing rule for multiphase kappa calculations; see KAPPA_METHODS for available options, [-]

Joule_Thomsonstr, optional

Mixing rule for multiphase Joule-Thomson calculations; see JT_METHODS for available options, [-]

Notes

The linear mixing rules “MOLE_WEIGHTED”, “MASS_WEIGHTED”, and “VOLUME_WEIGHTED” have the following formula, with β\beta representing molar, mass, or volume phase fraction:

bulk property=(iphasesβiproperty)\text{bulk property} = \left(\sum_i^{phases} \beta_i \text{property} \right)

The power mixing rules “POWER_PROP_MOLE_WEIGHTED”, “POWER_PROP_MASS_WEIGHTED”, and “POWER_PROP_VOLUME_WEIGHTED” have the following formula, with β\beta representing molar, mass, or volume phase fraction:

bulk property=(iphasesβiproperty exponent)1/exponent\text{bulk property} = \left(\sum_i^{phases} \beta_i \text{property }^{\text{exponent}} \right)^{1/\text{exponent}}

The logarithmic mixing rules “LOG_PROP_MOLE_WEIGHTED”, “LOG_PROP_MASS_WEIGHTED”, and “LOG_PROP_VOLUME_WEIGHTED” have the following formula, with β\beta representing molar, mass, or volume phase fraction:

bulk property=exp(iphasesβiln(property ))\text{bulk property} = \exp\left(\sum_i^{phases} \beta_i \ln(\text{property })\right)

The mixing rule “MINIMUM_PHASE_PROP” selects the lowest phase value of the property, always. The mixing rule “MAXIMUM_PHASE_PROP” selects the highest phase value of the property, always.

The mixing rule “AS_ONE_LIQUID” calculates a property using the bulk composition but applied to the liquid model only. The mixing rule “AS_ONE_GAS” calculates a property using the bulk composition but applied to the gas model only.

The mixing rule “FROM_DERIVATIVE_SETTINGS” is used to indicate that the property depends on other configurable properties; and when this is the specified option, those configurations will be used in the calculation of this property.

The mixing rule “EQUILIBRIUM_DERIVATIVE” performs derivative calculations on flashes themselves. This is quite slow in comparison to other methods.

References

1(1,2)

14:00-17:00. “ISO 10780:1994.” ISO. Accessed March 29, 2021. https://www.iso.org/cms/render/live/en/sites/isoorg/contents/data/standard/01/88/18855.html.

2(1,2)

14:00-17:00. “ISO 13443:1996.” ISO. Accessed March 29, 2021. https://www.iso.org/cms/render/live/en/sites/isoorg/contents/data/standard/02/04/20461.html.

Methods

as_json

as_json()[source]
settings = ('dP_dT', 'dP_dV', 'd2P_dV2', 'd2P_dT2', 'd2P_dTdV', 'mu_LL', 'mu_LL_power_exponent', 'mu_VL', 'mu_VL_power_exponent', 'k_LL', 'k_LL_power_exponent', 'k_VL', 'k_VL_power_exponent', 'sigma_LL', 'sigma_LL_power_exponent', 'T_liquid_volume_ref', 'T_normal', 'P_normal', 'T_standard', 'P_standard', 'T_gas_ref', 'P_gas_ref', 'speed_of_sound', 'kappa', 'isobaric_expansion', 'Joule_Thomson', 'VL_ID', 'VL_ID_settings', 'S_ID', 'S_ID_settings', 'solid_sort_method', 'liquid_sort_method', 'liquid_sort_cmps', 'solid_sort_cmps', 'liquid_sort_cmps_neg', 'solid_sort_cmps_neg', 'liquid_sort_prop', 'solid_sort_prop', 'phase_sort_higher_first', 'water_sort', 'equilibrium_perturbation')
thermo.bulk.DP_DT_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'EQUILIBRIUM_DERIVATIVE', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP']

List of all valid and implemented calculation methods for the DP_DT bulk setting

thermo.bulk.DP_DV_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'EQUILIBRIUM_DERIVATIVE', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP']

List of all valid and implemented calculation methods for the DP_DV bulk setting

thermo.bulk.D2P_DV2_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP']

List of all valid and implemented calculation methods for the D2P_DV2 bulk setting

thermo.bulk.D2P_DT2_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP']

List of all valid and implemented calculation methods for the D2P_DT2 bulk setting

thermo.bulk.D2P_DTDV_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP']

List of all valid and implemented calculation methods for the D2P_DTDV bulk setting

thermo.bulk.MU_LL_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'AS_ONE_LIQUID', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'POWER_PROP_MOLE_WEIGHTED', 'POWER_PROP_MASS_WEIGHTED', 'POWER_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP']

List of all valid and implemented mixing rules for the MU_LL setting

thermo.bulk.MU_VL_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'AS_ONE_LIQUID', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'POWER_PROP_MOLE_WEIGHTED', 'POWER_PROP_MASS_WEIGHTED', 'POWER_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP', 'AS_ONE_GAS', 'Beattie Whalley', 'McAdams', 'Cicchitti', 'Lin Kwok', 'Fourar Bories', 'Duckler']

List of all valid and implemented mixing rules for the MU_VL setting

thermo.bulk.K_LL_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'AS_ONE_LIQUID', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'POWER_PROP_MOLE_WEIGHTED', 'POWER_PROP_MASS_WEIGHTED', 'POWER_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP']

List of all valid and implemented mixing rules for the K_LL setting

thermo.bulk.K_VL_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'AS_ONE_LIQUID', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'POWER_PROP_MOLE_WEIGHTED', 'POWER_PROP_MASS_WEIGHTED', 'POWER_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP', 'AS_ONE_GAS']

List of all valid and implemented mixing rules for the K_VL setting

thermo.bulk.SIGMA_LL_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'AS_ONE_LIQUID', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'POWER_PROP_MOLE_WEIGHTED', 'POWER_PROP_MASS_WEIGHTED', 'POWER_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP']

List of all valid and implemented mixing rules for the SIGMA_LL setting

thermo.bulk.BETA_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP', 'EQUILIBRIUM_DERIVATIVE', 'FROM_DERIVATIVE_SETTINGS']

List of all valid and implemented calculation methods for the isothermal_compressibility bulk setting

thermo.bulk.SPEED_OF_SOUND_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP', 'FROM_DERIVATIVE_SETTINGS', 'EQUILIBRIUM_DERIVATIVE', 'FIROOZABADI_PAN']

List of all valid and implemented calculation methods for the speed_of_sound bulk setting

thermo.bulk.KAPPA_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP', 'EQUILIBRIUM_DERIVATIVE', 'FROM_DERIVATIVE_SETTINGS']

List of all valid and implemented calculation methods for the kappa bulk setting

thermo.bulk.JT_METHODS = ['MOLE_WEIGHTED', 'MASS_WEIGHTED', 'VOLUME_WEIGHTED', 'LOG_PROP_MOLE_WEIGHTED', 'LOG_PROP_MASS_WEIGHTED', 'LOG_PROP_VOLUME_WEIGHTED', 'MINIMUM_PHASE_PROP', 'MAXIMUM_PHASE_PROP', 'EQUILIBRIUM_DERIVATIVE', 'FROM_DERIVATIVE_SETTINGS']

List of all valid and implemented calculation methods for the JT bulk setting