Legacy Mixtures (thermo.mixture)

class thermo.mixture.Mixture(IDs=None, zs=None, ws=None, Vfls=None, Vfgs=None, T=None, P=None, VF=None, H=None, Hm=None, S=None, Sm=None, pkg=None, Vf_TP=(None, None))[source]

Bases: object

Creates a Mixture object which contains basic information such as molecular weight and the structure of the species, as well as thermodynamic and transport properties as a function of two of the variables temperature, pressure, vapor fraction, enthalpy, or entropy.

The components of the mixture must be specified by specifying the names of the chemicals; the composition can be specified by providing any one of the following parameters:

  • Mass fractions ws

  • Mole fractions zs

  • Liquid volume fractions (based on pure component densities) Vfls

  • Gas volume fractions (based on pure component densities) Vfgs

If volume fractions are provided, by default the pure component volumes are calculated at the specified T and P. To use another reference temperature and pressure specify it as a tuple for the argument Vf_TP.

If no thermodynamic conditions are specified, or if only one of T and P are specifed without another thermodynamic variable as well, the T and P 298.15 K and/or 101325 Pa will be set instead of the missing variables.

Parameters
IDslist, optional

List of chemical identifiers - names, CAS numbers, SMILES or InChi strings can all be recognized and may be mixed [-]

zslist or dict, optional

Mole fractions of all components in the mixture [-]

wslist or dict, optional

Mass fractions of all components in the mixture [-]

Vflslist or dict, optional

Volume fractions of all components as a hypothetical liquid phase based on pure component densities [-]

Vfgslist, or dict optional

Volume fractions of all components as a hypothetical gas phase based on pure component densities [-]

Tfloat, optional

Temperature of the mixture (default 298.15 K), [K]

Pfloat, optional

Pressure of the mixture (default 101325 Pa) [Pa]

VFfloat, optional

Vapor fraction (mole basis) of the mixture, [-]

Hmfloat, optional

Molar enthalpy of the mixture, [J/mol]

Hfloat, optional

Mass enthalpy of the mixture, [J/kg]

Smfloat, optional

Molar entropy of the mixture, [J/mol/K]

Sfloat, optional

Mass entropy of the mixture, [J/kg/K]

pkgobject

The thermodynamic property package to use for flash calculations; one of the caloric packages in thermo.property_package; defaults to the ideal model [-]

Vf_TPtuple(2, float), optional

The (T, P) at which the volume fractions are specified to be at, [K] and [Pa]

Notes

Warning

The Mixture class is not designed for high-performance or the ability to use different thermodynamic models. It is especially limited in its multiphase support and the ability to solve with specifications other than temperature and pressure. It is impossible to change constant properties such as a compound’s critical temperature in this interface.

It is recommended to switch over to the thermo.flash interface which solves those problems and is better positioned to grow. That interface also requires users to be responsible for their chemical constants and pure component correlations; while default values can easily be loaded for most compounds, the user is ultimately responsible for them.

Examples

Creating Mixture objects:

>>> Mixture(['water', 'ethanol'], Vfls=[.6, .4], T=300, P=1E5)
<Mixture, components=['water', 'ethanol'], mole fractions=[0.8299, 0.1701], T=300.00 K, P=100000 Pa>

For mixtures with large numbers of components, it may be confusing to enter the composition separate from the names of the chemicals. For that case, the syntax using dictionaries as follows is supported with any composition specification:

>>> comp = OrderedDict([('methane', 0.96522),
...                     ('nitrogen', 0.00259),
...                     ('carbon dioxide', 0.00596),
...                     ('ethane', 0.01819),
...                     ('propane', 0.0046),
...                     ('isobutane', 0.00098),
...                     ('butane', 0.00101),
...                     ('2-methylbutane', 0.00047),
...                     ('pentane', 0.00032),
...                     ('hexane', 0.00066)])
>>> m = Mixture(zs=comp)
Attributes
MWfloat

Mole-weighted average molecular weight all chemicals in the mixture, [g/mol]

IDslist of str

Names of all the species in the mixture as given in the input, [-]

nameslist of str

Names of all the species in the mixture, [-]

CASslist of str

CAS numbers of all species in the mixture, [-]

MWslist of float

Molecular weights of all chemicals in the mixture, [g/mol]

Tmslist of float

Melting temperatures of all chemicals in the mixture, [K]

Tbslist of float

Boiling temperatures of all chemicals in the mixture, [K]

Tcslist of float

Critical temperatures of all chemicals in the mixture, [K]

Pcslist of float

Critical pressures of all chemicals in the mixture, [Pa]

Vcslist of float

Critical volumes of all chemicals in the mixture, [m^3/mol]

Zcslist of float

Critical compressibilities of all chemicals in the mixture, [-]

rhocslist of float

Critical densities of all chemicals in the mixture, [kg/m^3]

rhocmslist of float

Critical molar densities of all chemicals in the mixture, [mol/m^3]

omegaslist of float

Acentric factors of all chemicals in the mixture, [-]

StielPolarslist of float

Stiel Polar factors of all chemicals in the mixture, see chemicals.acentric.Stiel_polar_factor for the definition, [-]

Ttslist of float

Triple temperatures of all chemicals in the mixture, [K]

Ptslist of float

Triple pressures of all chemicals in the mixture, [Pa]

Hfusslist of float

Enthalpy of fusions of all chemicals in the mixture, [J/kg]

Hfusmslist of float

Molar enthalpy of fusions of all chemicals in the mixture, [J/mol]

Hsubslist of float

Enthalpy of sublimations of all chemicals in the mixture, [J/kg]

Hsubmslist of float

Molar enthalpy of sublimations of all chemicals in the mixture, [J/mol]

Hfmslist of float

Molar enthalpy of formations of all chemicals in the mixture, [J/mol]

Hfslist of float

Enthalpy of formations of all chemicals in the mixture, [J/kg]

Gfmslist of float

Molar Gibbs free energies of formation of all chemicals in the mixture, [J/mol]

Gfslist of float

Gibbs free energies of formation of all chemicals in the mixture, [J/kg]

Sfmslist of float

Molar entropy of formation of all chemicals in the mixture, [J/mol/K]

Sfslist of float

Entropy of formation of all chemicals in the mixture, [J/kg/K]

S0mslist of float

Standard absolute entropies of all chemicals in the mixture, [J/mol/K]

S0slist of float

Standard absolute entropies of all chemicals in the mixture, [J/kg/K]

Hcmslist of float

Molar higher heats of combustions of all chemicals in the mixture, [J/mol]

Hcslist of float

Higher heats of combustions of all chemicals in the mixture, [J/kg]

Hcms_lowerlist of float

Molar lower heats of combustions of all chemicals in the mixture, [J/mol]

Hcs_lowerlist of float

Higher lower of combustions of all chemicals in the mixture, [J/kg]

Tflashslist of float

Flash points of all chemicals in the mixture, [K]

Tautoignitionslist of float

Autoignition points of all chemicals in the mixture, [K]

LFLslist of float

Lower flammability limits of the gases in an atmosphere at STP, mole fractions, [-]

UFLslist of float

Upper flammability limit of the gases in an atmosphere at STP, mole fractions, [-]

TWAslist of list of tuple(quantity, unit)

Time-Weighted Average limits on worker exposure to dangerous chemicals.

STELslist of tuple(quantity, unit)

Short-term Exposure limits on worker exposure to dangerous chemicals.

Ceilingslist of tuple(quantity, unit)

Ceiling limits on worker exposure to dangerous chemicals.

Skinslist of bool

Whether or not each of the chemicals can be absorbed through the skin.

Carcinogenslist of str or dict

Carcinogen status information for each chemical in the mixture.

Chemicalslist of Chemical instances

Chemical instances used in calculating mixture properties, [-]

dipoleslist of float

Dipole moments of all chemicals in the mixture in debye, [3.33564095198e-30 ampere*second^2]

Stockmayerslist of float

Lennard-Jones depth of potential-energy minimum over k for all chemicals in the mixture, [K]

molecular_diameterslist of float

Lennard-Jones molecular diameters of all chemicals in the mixture, [angstrom]

GWPslist of float

Global warming potentials (default 100-year outlook) (impact/mass chemical)/(impact/mass CO2) of all chemicals in the mixture, [-]

ODPslist of float

Ozone Depletion potentials (impact/mass chemical)/(impact/mass CFC-11), of all chemicals in the mixture, [-]

logPslist of float

Octanol-water partition coefficients of all chemicals in the mixture, [-]

Psat_298slist of float

Vapor pressure of the chemicals in the mixture at 298.15 K, [Pa]

phase_STPslist of str

Phase of the chemicals in the mixture at 298.15 K and 101325 Pa; one of ‘s’, ‘l’, ‘g’, or ‘l/g’.

Vml_Tbslist of float

Molar volumes of the chemicals in the mixture as liquids at their normal boiling points, [m^3/mol]

Vml_Tmslist of float

Molar volumes of the chemicals in the mixture as liquids at their melting points, [m^3/mol]

Vml_STPslist of float

Molar volume of the chemicals in the mixture as liquids at 298.15 K and 101325 Pa, [m^3/mol]

rhoml_STPslist of float

Molar densities of the chemicals in the mixture as liquids at 298.15 K and 101325 Pa, [mol/m^3]

Vmg_STPslist of float

Molar volume of the chemicals in the mixture as gases at 298.15 K and 101325 Pa, [m^3/mol]

Vms_Tmslist of float

Molar volumes of solid phase at the melting point [m^3/mol]

rhos_Tmslist of float

Mass densities of solid phase at the melting point [kg/m^3]

Hvap_Tbmslist of float

Molar enthalpies of vaporization of the chemicals in the mixture at their normal boiling points, [J/mol]

Hvap_Tbslist of float

Mass enthalpies of vaporization of the chemicals in the mixture at their normal boiling points, [J/kg]

alpha

Thermal diffusivity of the mixture at its current temperature, pressure, and phase in units of [m^2/s].

alphag

Thermal diffusivity of the gas phase of the mixture if one exists at its current temperature and pressure, in units of [m^2/s].

alphags

Pure component thermal diffusivities of the chemicals in the mixture in the gas phase at the current temperature and pressure, in units of [m^2/s].

alphal

Thermal diffusivity of the liquid phase of the mixture if one exists at its current temperature and pressure, in units of [m^2/s].

alphals

Pure component thermal diffusivities of the chemicals in the mixture in the liquid phase at the current temperature and pressure, in units of [m^2/s].

A

Helmholtz energy of the mixture at its current state, in units of [J/kg].

Am

Helmholtz energy of the mixture at its current state, in units of [J/mol].

atom_fractions

Dictionary of atomic fractions for each atom in the mixture.

atom_fractionss

List of dictionaries of atomic fractions for all chemicals in the mixture.

atomss

List of dictionaries of atom counts for all chemicals in the mixture.

Bvirial

Second virial coefficient of the gas phase of the mixture at its current temperature, pressure, and composition in units of [mol/m^3].

charges

Charges for all chemicals in the mixture, [faraday].

Cp

Mass heat capacity of the mixture at its current phase and temperature, in units of [J/kg/K].

Cpg

Gas-phase heat capacity of the mixture at its current temperature , and composition in units of [J/kg/K].

Cpgm

Gas-phase heat capacity of the mixture at its current temperature and composition, in units of [J/mol/K].

Cpgms

Gas-phase ideal gas heat capacity of the chemicals at its current temperature, in units of [J/mol/K].

Cpgs

Gas-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/kg/K].

Cpl

Liquid-phase heat capacity of the mixture at its current temperature and composition, in units of [J/kg/K].

Cplm

Liquid-phase heat capacity of the mixture at its current temperature and composition, in units of [J/mol/K].

Cplms

Liquid-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/mol/K].

Cpls

Liquid-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/kg/K].

Cpm

Molar heat capacity of the mixture at its current phase and temperature, in units of [J/mol/K].

Cps

Solid-phase heat capacity of the mixture at its current temperature and composition, in units of [J/kg/K].

Cpsm

Solid-phase heat capacity of the mixture at its current temperature and composition, in units of [J/mol/K].

Cpsms

Solid-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/mol/K].

Cpss

Solid-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/kg/K].

Cvg

Gas-phase ideal-gas contant-volume heat capacity of the mixture at its current temperature, in units of [J/kg/K].

Cvgm

Gas-phase ideal-gas contant-volume heat capacity of the mixture at its current temperature and composition, in units of [J/mol/K].

Cvgms

Gas-phase pure component ideal-gas contant-volume heat capacities of the chemicals in the mixture at its current temperature, in units of [J/mol/K].

Cvgs

Gas-phase pure component ideal-gas contant-volume heat capacities of the chemicals in the mixture at its current temperature, in units of [J/kg/K].

economic_statuses

List of dictionaries of the economic status for all chemicals in the mixture.

eos

Equation of state object held by the mixture.

formulas

Chemical formulas for all chemicals in the mixture.

Hvapms

Pure component enthalpies of vaporization of the chemicals in the mixture at its current temperature, in units of [J/mol].

Hvaps

Enthalpy of vaporization of the chemicals in the mixture at its current temperature, in units of [J/kg].

InChI_Keys

InChI keys for all chemicals in the mixture.

InChIs

InChI strings for all chemicals in the mixture.

isentropic_exponent

Gas-phase ideal-gas isentropic exponent of the mixture at its current temperature, [dimensionless].

isentropic_exponents

Gas-phase pure component ideal-gas isentropic exponent of the chemicals in the mixture at its current temperature, [dimensionless].

isobaric_expansion

Isobaric (constant-pressure) expansion of the mixture at its current phase, temperature, and pressure in units of [1/K].

isobaric_expansion_g

Isobaric (constant-pressure) expansion of the gas phase of the mixture at its current temperature and pressure, in units of [1/K].

isobaric_expansion_gs

Pure component isobaric (constant-pressure) expansions of the chemicals in the mixture in the gas phase at its current temperature and pressure, in units of [1/K].

isobaric_expansion_l

Isobaric (constant-pressure) expansion of the liquid phase of the mixture at its current temperature and pressure, in units of [1/K].

isobaric_expansion_ls

Pure component isobaric (constant-pressure) expansions of the chemicals in the mixture in the liquid phase at its current temperature and pressure, in units of [1/K].

IUPAC_names

IUPAC names for all chemicals in the mixture.

JT

Joule Thomson coefficient of the mixture at its current phase, temperature, and pressure in units of [K/Pa].

JTg

Joule Thomson coefficient of the gas phase of the mixture if one exists at its current temperature and pressure, in units of [K/Pa].

JTgs

Pure component Joule Thomson coefficients of the chemicals in the mixture in the gas phase at its current temperature and pressure, in units of [K/Pa].

JTl

Joule Thomson coefficient of the liquid phase of the mixture if one exists at its current temperature and pressure, in units of [K/Pa].

JTls

Pure component Joule Thomson coefficients of the chemicals in the mixture in the liquid phase at its current temperature and pressure, in units of [K/Pa].

k

Thermal conductivity of the mixture at its current phase, temperature, and pressure in units of [W/m/K].

kg

Thermal conductivity of the mixture in the gas phase at its current temperature, pressure, and composition in units of [Pa*s].

kgs

Pure component thermal conductivies of the chemicals in the mixture in the gas phase at its current temperature and pressure, in units of [W/m/K].

kl

Thermal conductivity of the mixture in the liquid phase at its current temperature, pressure, and composition in units of [Pa*s].

kls

Pure component thermal conductivities of the chemicals in the mixture in the liquid phase at its current temperature and pressure, in units of [W/m/K].

legal_statuses

List of dictionaries of the legal status for all chemicals in the mixture.

mass_fractions

Dictionary of mass fractions for each atom in the mixture.

mass_fractionss

List of dictionaries of mass fractions for all chemicals in the mixture.

mu

Viscosity of the mixture at its current phase, temperature, and pressure in units of [Pa*s].

mug

Viscosity of the mixture in the gas phase at its current temperature, pressure, and composition in units of [Pa*s].

mugs

Pure component viscosities of the chemicals in the mixture in the gas phase at its current temperature and pressure, in units of [Pa*s].

mul

Viscosity of the mixture in the liquid phase at its current temperature, pressure, and composition in units of [Pa*s].

muls

Pure component viscosities of the chemicals in the mixture in the liquid phase at its current temperature and pressure, in units of [Pa*s].

nu

Kinematic viscosity of the the mixture at its current temperature, pressure, and phase in units of [m^2/s].

nug

Kinematic viscosity of the gas phase of the mixture if one exists at its current temperature and pressure, in units of [m^2/s].

nugs

Pure component kinematic viscosities of the gas phase of the chemicals in the mixture at its current temperature and pressure, in units of [m^2/s].

nul

Kinematic viscosity of the liquid phase of the mixture if one exists at its current temperature and pressure, in units of [m^2/s].

nuls

Pure component kinematic viscosities of the liquid phase of the chemicals in the mixture at its current temperature and pressure, in units of [m^2/s].

permittivites

Pure component relative permittivities of the chemicals in the mixture at its current temperature, [dimensionless].

Pr

Prandtl number of the mixture at its current temperature, pressure, and phase; [dimensionless].

Prg

Prandtl number of the gas phase of the mixture if one exists at its current temperature and pressure, [dimensionless].

Prgs

Pure component Prandtl numbers of the gas phase of the chemicals in the mixture at its current temperature and pressure, [dimensionless].

Prl

Prandtl number of the liquid phase of the mixture if one exists at its current temperature and pressure, [dimensionless].

Prls

Pure component Prandtl numbers of the liquid phase of the chemicals in the mixture at its current temperature and pressure, [dimensionless].

Psats

Pure component vapor pressures of the chemicals in the mixture at its current temperature, in units of [Pa].

PSRK_groups

List of dictionaries of PSRK subgroup: count groups for each chemical in the mixture.

PubChems

PubChem Component ID numbers for all chemicals in the mixture.

rho

Mass density of the mixture at its current phase and temperature and pressure, in units of [kg/m^3].

rhog

Gas-phase mass density of the mixture at its current temperature, pressure, and composition in units of [kg/m^3].

rhogm

Molar density of the mixture in the gas phase at the current temperature, pressure, and composition in units of [mol/m^3].

rhogms

Pure component molar densities of the chemicals in the gas phase at the current temperature and pressure, in units of [mol/m^3].

rhogm_STP

Molar density of the mixture in the gas phase at 298.15 K and 101.325 kPa, and the current composition, in units of [mol/m^3].

rhogs

Pure-component gas-phase mass densities of the chemicals in the mixture at its current temperature and pressure, in units of [kg/m^3].

rhog_STP

Gas-phase mass density of the mixture at 298.15 K and 101.325 kPa, and the current composition in units of [kg/m^3].

rhol

Liquid-phase mass density of the mixture at its current temperature, pressure, and composition in units of [kg/m^3].

rholm

Molar density of the mixture in the liquid phase at the current temperature, pressure, and composition in units of [mol/m^3].

rholms

Pure component molar densities of the chemicals in the mixture in the liquid phase at the current temperature and pressure, in units of [mol/m^3].

rholm_STP

Molar density of the mixture in the liquid phase at 298.15 K and 101.325 kPa, and the current composition, in units of [mol/m^3].

rhols

Pure-component liquid-phase mass density of the chemicals in the mixture at its current temperature and pressure, in units of [kg/m^3].

rhol_STP

Liquid-phase mass density of the mixture at 298.15 K and 101.325 kPa, and the current composition in units of [kg/m^3].

rhom

Molar density of the mixture at its current phase and temperature and pressure, in units of [mol/m^3].

rhosms

Pure component molar densities of the chemicals in the solid phase at the current temperature and pressure, in units of [mol/m^3].

rhoss

Pure component solid-phase mass density of the chemicals in the mixture at its current temperature, in units of [kg/m^3].

ringss

List of ring counts for all chemicals in the mixture.

sigma

Surface tension of the mixture at its current temperature and composition, in units of [N/m].

sigmas

Pure component surface tensions of the chemicals in the mixture at its current temperature, in units of [N/m].

smiless

SMILES strings for all chemicals in the mixture.

solubility_parameters

Pure component solubility parameters of the chemicals in the mixture at its current temperature and pressure, in units of [Pa^0.5].

synonymss

Lists of synonyms for all chemicals in the mixture.

U

Internal energy of the mixture at its current state, in units of [J/kg].

Um

Internal energy of the mixture at its current state, in units of [J/mol].

UNIFAC_Dortmund_groups

List of dictionaries of Dortmund UNIFAC subgroup: count groups for each chemcial in the mixture.

UNIFAC_groups

List of dictionaries of UNIFAC subgroup: count groups for each chemical in the mixture.

Vm

Molar volume of the mixture at its current phase and temperature and pressure, in units of [m^3/mol].

Vmg

Gas-phase molar volume of the mixture at its current temperature, pressure, and composition in units of [m^3/mol].

Vmgs

Pure component gas-phase molar volumes of the chemicals in the mixture at its current temperature and pressure, in units of [m^3/mol].

Vmg_STP

Gas-phase molar volume of the mixture at 298.15 K and 101.325 kPa, and the current composition in units of [m^3/mol].

Vml

Liquid-phase molar volume of the mixture at its current temperature, pressure, and composition in units of [m^3/mol].

Vmls

Pure component liquid-phase molar volumes of the chemicals in the mixture at its current temperature and pressure, in units of [m^3/mol].

Vml_STP

Liquid-phase molar volume of the mixture at 298.15 K and 101.325 kPa, and the current composition in units of [m^3/mol].

Vmss

Pure component solid-phase molar volumes of the chemicals in the mixture at its current temperature, in units of [m^3/mol].

Z

Compressibility factor of the mixture at its current phase and temperature and pressure, [dimensionless].

Zg

Compressibility factor of the mixture in the gas phase at the current temperature, pressure, and composition, [dimensionless].

Zgs

Pure component compressibility factors of the chemicals in the mixture in the gas phase at the current temperature and pressure, [dimensionless].

Zg_STP

Gas-phase compressibility factor of the mixture at 298.15 K and 101.325 kPa, and the current composition, [dimensionless].

Zl

Compressibility factor of the mixture in the liquid phase at the current temperature, pressure, and composition, [dimensionless].

Zls

Pure component compressibility factors of the chemicals in the liquid phase at the current temperature and pressure, [dimensionless].

Zl_STP

Liquid-phase compressibility factor of the mixture at 298.15 K and 101.325 kPa, and the current composition, [dimensionless].

Zss

Pure component compressibility factors of the chemicals in the mixture in the solid phase at the current temperature and pressure, [dimensionless].

Methods

Hc_volumetric_g([T, P])

Standard higher molar heat of combustion of the mixture, in units of [J/m^3] at the specified T and P in the gas phase.

Hc_volumetric_g_lower([T, P])

Standard lower molar heat of combustion of the mixture, in units of [J/m^3] at the specified T and P in the gas phase.

Vfgs([T, P])

Volume fractions of all species in a hypothetical pure-gas phase at the current or specified temperature and pressure.

Vfls([T, P])

Volume fractions of all species in a hypothetical pure-liquid phase at the current or specified temperature and pressure.

draw_2d([Hs])

Interface for drawing a 2D image of all the molecules in the mixture.

set_chemical_TP([T, P])

Basic method to change all chemical instances to be at the T and P specified.

set_chemical_constants()

Basic method which retrieves and sets constants of chemicals to be accessible as lists from a Mixture object.

Bond

Capillary

Grashof

Jakob

Peclet_heat

Reynolds

Weber

compound_index

eos_pures

flash_caloric

properties

set_Chemical_property_objects

set_TP_sources

set_constant_sources

set_constants

set_eos

set_property_package

property A

Helmholtz energy of the mixture at its current state, in units of [J/kg].

This property requires that the property package of the mixture found a solution to the given state variables. It also depends on the molar volume of the mixture at its current conditions.

property API

API gravity of the hypothetical liquid phase of the mixture, [degrees]. The reference condition is water at 15.6 °C (60 °F) and 1 atm (rho=999.016 kg/m^3, standardized).

Examples

>>> Mixture(['hexane', 'decane'], ws=[0.5, 0.5]).API
71.34707841728181
property Am

Helmholtz energy of the mixture at its current state, in units of [J/mol].

This property requires that the property package of the mixture found a solution to the given state variables. It also depends on the molar volume of the mixture at its current conditions.

Bond(L=None)[source]
property Bvirial

Second virial coefficient of the gas phase of the mixture at its current temperature, pressure, and composition in units of [mol/m^3].

This property uses the object-oriented interface thermo.volume.VolumeGasMixture, converting its result with thermo.utils.B_from_Z.

Examples

>>> Mixture(['hexane'], ws=[1], T=300, P=1E5).Bvirial
-0.001486976173801296
Capillary(V=None)[source]
property Cp

Mass heat capacity of the mixture at its current phase and temperature, in units of [J/kg/K].

Examples

>>> w = Mixture(['water'], ws=[1])
>>> w.Cp, w.phase
(4180.597021827336, 'l')
>>> Pd = Mixture(['palladium'], ws=[1])
>>> Pd.Cp, Pd.phase
(234.26767209171211, 's')
property Cpg

Gas-phase heat capacity of the mixture at its current temperature , and composition in units of [J/kg/K]. For calculation of this property at other temperatures or compositions, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.heat_capacity.HeatCapacityGasMixture; each Mixture instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> Mixture(['oxygen', 'nitrogen'], ws=[.4, .6], T=350, P=1E6).Cpg
995.8911053614883
property Cpgm

Gas-phase heat capacity of the mixture at its current temperature and composition, in units of [J/mol/K]. For calculation of this property at other temperatures or compositions, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.heat_capacity.HeatCapacityGasMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['oxygen', 'nitrogen'], ws=[.4, .6], T=350, P=1E6).Cpgm
29.361044582498046
property Cpgms

Gas-phase ideal gas heat capacity of the chemicals at its current temperature, in units of [J/mol/K].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Cpgms
[89.55804092586159, 111.70390334788907]
property Cpgs

Gas-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/kg/K].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Cpgs
[1146.5360555565146, 1212.3488046342566]
property Cpl

Liquid-phase heat capacity of the mixture at its current temperature and composition, in units of [J/kg/K]. For calculation of this property at other temperatures or compositions, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.heat_capacity.HeatCapacityLiquidMixture; each Mixture instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> Mixture(['water', 'sodium chloride'], ws=[.9, .1], T=301.5).Cpl
3735.4604049449786
property Cplm

Liquid-phase heat capacity of the mixture at its current temperature and composition, in units of [J/mol/K]. For calculation of this property at other temperatures or compositions, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.heat_capacity.HeatCapacityLiquidMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['toluene', 'decane'], ws=[.9, .1], T=300).Cplm
168.29127923518843
property Cplms

Liquid-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/mol/K].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Cplms
[140.9113971170526, 163.62584810669068]
property Cpls

Liquid-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/kg/K].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Cpls
[1803.9697581961016, 1775.869915141704]
property Cpm

Molar heat capacity of the mixture at its current phase and temperature, in units of [J/mol/K]. Available only if single phase.

Examples

>>> Mixture(['ethylbenzene'], ws=[1], T=550, P=3E6).Cpm
294.18449553310046
property Cps

Solid-phase heat capacity of the mixture at its current temperature and composition, in units of [J/kg/K]. For calculation of this property at other temperatures or compositions, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.heat_capacity.HeatCapacitySolidMixture; each Mixture instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> Mixture(['silver', 'platinum'], ws=[0.95, 0.05]).Cps
229.55166388430328
property Cpsm

Solid-phase heat capacity of the mixture at its current temperature and composition, in units of [J/mol/K]. For calculation of this property at other temperatures or compositions, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.heat_capacity.HeatCapacitySolidMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['silver', 'platinum'], ws=[0.95, 0.05]).Cpsm
25.32745796347474
property Cpsms

Solid-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/mol/K].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Cpsms
[109.77384365511931, 135.22614707678474]
property Cpss

Solid-phase pure component heat capacity of the chemicals in the mixture at its current temperature, in units of [J/kg/K].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Cpss
[1405.341925822248, 1467.6412627521154]
property Cvg

Gas-phase ideal-gas contant-volume heat capacity of the mixture at its current temperature, in units of [J/kg/K]. Subtracts R from the ideal-gas heat capacity; does not include pressure-compensation from an equation of state.

Examples

>>> Mixture(['water'], ws=[1], T=520).Cvg
1506.1471795798861
property Cvgm

Gas-phase ideal-gas contant-volume heat capacity of the mixture at its current temperature and composition, in units of [J/mol/K]. Subtracts R from the ideal-gas heat capacity; does not include pressure-compensation from an equation of state.

Examples

>>> Mixture(['water'], ws=[1], T=520).Cvgm
27.13366316134193
property Cvgms

Gas-phase pure component ideal-gas contant-volume heat capacities of the chemicals in the mixture at its current temperature, in units of [J/mol/K]. Subtracts R from the ideal-gas heat capacities; does not include pressure-compensation from an equation of state.

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Cvgms
[81.2435811258616, 103.38944354788907]
property Cvgs

Gas-phase pure component ideal-gas contant-volume heat capacities of the chemicals in the mixture at its current temperature, in units of [J/kg/K]. Subtracts R from the ideal-gas heat capacity; does not include pressure-compensation from an equation of state.

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Cvgs
[1040.093040003431, 1122.1100117398266]
Grashof(Tw=None, L=None)[source]
H = None
property Hc

Standard higher heat of combustion of the mixture, in units of [J/kg].

This property depends on the bulk composition only.

property Hc_lower

Standard lower heat of combustion of the mixture, in units of [J/kg].

This property depends on the bulk composition only.

Hc_volumetric_g(T=288.7055555555555, P=101325.0)[source]

Standard higher molar heat of combustion of the mixture, in units of [J/m^3] at the specified T and P in the gas phase.

This property depends on the bulk composition only.

Parameters
Tfloat, optional

Reference temperature, [K]

Pfloat, optional

Reference pressure, [Pa]

Returns
Hc_volumetric_gfloat, optional

Higher heat of combustion on a volumetric basis, [J/m^3]

Hc_volumetric_g_lower(T=288.7055555555555, P=101325.0)[source]

Standard lower molar heat of combustion of the mixture, in units of [J/m^3] at the specified T and P in the gas phase.

This property depends on the bulk composition only.

Parameters
Tfloat, optional

Reference temperature, [K]

Pfloat, optional

Reference pressure, [Pa]

Returns
Hc_volumetric_gfloat, optional

Lower heat of combustion on a volumetric basis, [J/m^3]

property Hcm

Standard higher molar heat of combustion of the mixture, in units of [J/mol].

This property depends on the bulk composition only.

property Hcm_lower

Standard lower molar heat of combustion of the mixture, in units of [J/mol].

This property depends on the bulk composition only.

Hm = None
property Hvapms

Pure component enthalpies of vaporization of the chemicals in the mixture at its current temperature, in units of [J/mol].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Hvapms
[32639.806783391632, 36851.7902195611]
property Hvaps

Enthalpy of vaporization of the chemicals in the mixture at its current temperature, in units of [J/kg].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Hvaps
[417859.9144942896, 399961.16950519773]
property IUPAC_names

IUPAC names for all chemicals in the mixture.

Examples

>>> Mixture(['1-hexene', '1-nonene'], zs=[.7, .3]).IUPAC_names
['hex-1-ene', 'non-1-ene']
property InChI_Keys

InChI keys for all chemicals in the mixture.

Examples

>>> Mixture(['1-nonene'], zs=[1]).InChI_Keys
['JRZJOMJEPLMPRA-UHFFFAOYSA-N']
property InChIs

InChI strings for all chemicals in the mixture.

Examples

>>> Mixture(['methane', 'ethane', 'propane', 'butane'],
... zs=[0.25, 0.25, 0.25, 0.25]).InChIs
['CH4/h1H4', 'C2H6/c1-2/h1-2H3', 'C3H8/c1-3-2/h3H2,1-2H3', 'C4H10/c1-3-4-2/h3-4H2,1-2H3']
property JT

Joule Thomson coefficient of the mixture at its current phase, temperature, and pressure in units of [K/Pa]. Available only if single phase.

μJT=(TP)H=1Cp[T(VT)PV]=VCp(βT1)\mu_{JT} = \left(\frac{\partial T}{\partial P}\right)_H = \frac{1}{C_p} \left[T \left(\frac{\partial V}{\partial T}\right)_P - V\right] = \frac{V}{C_p}\left(\beta T-1\right)

Examples

>>> Mixture(['water'], ws=[1]).JT
-2.2150394958666412e-07
property JTg

Joule Thomson coefficient of the gas phase of the mixture if one exists at its current temperature and pressure, in units of [K/Pa].

μJT=(TP)H=1Cp[T(VT)PV]=VCp(βT1)\mu_{JT} = \left(\frac{\partial T}{\partial P}\right)_H = \frac{1}{C_p} \left[T \left(\frac{\partial V}{\partial T}\right)_P - V\right] = \frac{V}{C_p}\left(\beta T-1\right)

Examples

>>> Mixture(['dodecane'], ws=[1], T=400, P=1000).JTg
5.4089897835384913e-05
property JTgs

Pure component Joule Thomson coefficients of the chemicals in the mixture in the gas phase at its current temperature and pressure, in units of [K/Pa].

μJT=(TP)H=1Cp[T(VT)PV]=VCp(βT1)\mu_{JT} = \left(\frac{\partial T}{\partial P}\right)_H = \frac{1}{C_p} \left[T \left(\frac{\partial V}{\partial T}\right)_P - V\right] = \frac{V}{C_p}\left(\beta T-1\right)

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).JTgs
[6.0940046688790938e-05, 4.1290005523287549e-05]
property JTl

Joule Thomson coefficient of the liquid phase of the mixture if one exists at its current temperature and pressure, in units of [K/Pa].

μJT=(TP)H=1Cp[T(VT)PV]=VCp(βT1)\mu_{JT} = \left(\frac{\partial T}{\partial P}\right)_H = \frac{1}{C_p} \left[T \left(\frac{\partial V}{\partial T}\right)_P - V\right] = \frac{V}{C_p}\left(\beta T-1\right)

Examples

>>> Mixture(['dodecane'], ws=[1], T=400).JTl
-3.193910574559279e-07
property JTls

Pure component Joule Thomson coefficients of the chemicals in the mixture in the liquid phase at its current temperature and pressure, in units of [K/Pa].

μJT=(TP)H=1Cp[T(VT)PV]=VCp(βT1)\mu_{JT} = \left(\frac{\partial T}{\partial P}\right)_H = \frac{1}{C_p} \left[T \left(\frac{\partial V}{\partial T}\right)_P - V\right] = \frac{V}{C_p}\left(\beta T-1\right)

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).JTls
[-3.8633730709853161e-07, -3.464395792560331e-07]
Jakob(Tw=None)[source]
property PSRK_groups

List of dictionaries of PSRK subgroup: count groups for each chemical in the mixture. Uses the PSRK subgroups, as determined by DDBST’s online service.

Examples

>>> Mixture(['1-pentanol', 'decane'], ws=[0.5, 0.5]).PSRK_groups
[{1: 1, 2: 4, 14: 1}, {1: 2, 2: 8}]
P_default = 101325.0
property Parachor

Parachor of the mixture at its current temperature and pressure, in units of [N^0.25*m^2.75/mol].

P=σ0.25MWρLρVP = \frac{\sigma^{0.25} MW}{\rho_L - \rho_V}

Calculated based on surface tension, density of the liquid and gas phase, and molecular weight. For uses of this property, see thermo.utils.Parachor.

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).Parachor
4.233407085050756e-05
property Parachors

Pure component Parachor parameters of the chemicals in the mixture at its current temperature and pressure, in units of [N^0.25*m^2.75/mol].

P=σ0.25MWρLρVP = \frac{\sigma^{0.25} MW}{\rho_L - \rho_V}

Calculated based on surface tension, density of the liquid and gas phase, and molecular weight. For uses of this property, see thermo.utils.Parachor.

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).Parachors
[3.6795616000855504e-05, 4.82947303150274e-05]
property Pbubble

Bubble point pressure of the mixture at its current temperature and composition, in units of [Pa].

This property requires that the property package of the mixture found a solution to the given state variables.

property Pdew

Dew point pressure of the mixture at its current temperature and composition, in units of [Pa].

This property requires that the property package of the mixture found a solution to the given state variables.

Peclet_heat(V=None, D=None)[source]
property Pr

Prandtl number of the mixture at its current temperature, pressure, and phase; [dimensionless]. Available only if single phase.

Pr=CpμkPr = \frac{C_p \mu}{k}

Examples

>>> Mixture(['acetone'], ws=[1]).Pr
4.183039103542711
property Prg

Prandtl number of the gas phase of the mixture if one exists at its current temperature and pressure, [dimensionless].

Pr=CpμkPr = \frac{C_p \mu}{k}

Examples

>>> Mixture(['NH3'], ws=[1]).Prg
0.8472637319330079
property Prgs

Pure component Prandtl numbers of the gas phase of the chemicals in the mixture at its current temperature and pressure, [dimensionless].

Pr=CpμkPr = \frac{C_p \mu}{k}

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).Prgs
[0.7810364900059606, 0.784358381123896]
property Prl

Prandtl number of the liquid phase of the mixture if one exists at its current temperature and pressure, [dimensionless].

Pr=CpμkPr = \frac{C_p \mu}{k}

Examples

>>> Mixture(['nitrogen'], ws=[1], T=70).Prl
2.782821450148889
property Prls

Pure component Prandtl numbers of the liquid phase of the chemicals in the mixture at its current temperature and pressure, [dimensionless].

Pr=CpμkPr = \frac{C_p \mu}{k}

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).Prls
[6.13542244155373, 5.034355147908088]
property Psats

Pure component vapor pressures of the chemicals in the mixture at its current temperature, in units of [Pa].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Psats
[32029.25774454549, 10724.419010511821]
property PubChems

PubChem Component ID numbers for all chemicals in the mixture.

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5]).PubChems
[241, 1140]
property R_specific

Specific gas constant of the mixture, in units of [J/kg/K].

Examples

>>> Mixture(['N2', 'O2'], zs=[0.79, .21]).R_specific
288.1928437986195
Reynolds(V=None, D=None)[source]
property SG

Specific gravity of the mixture, [dimensionless].

For gas-phase conditions, this is calculated at 15.6 °C (60 °F) and 1 atm for the mixture and the reference fluid, air. For liquid and solid phase conditions, this is calculated based on a reference fluid of water at 4°C at 1 atm, but the with the liquid or solid mixture’s density at the currently specified conditions.

Examples

>>> Mixture('MTBE').SG
0.7428160596603596
property SGg

Specific gravity of a hypothetical gas phase of the mixture, . [dimensionless]. The reference condition is air at 15.6 °C (60 °F) and 1 atm (rho=1.223 kg/m^3). The definition for gases uses the compressibility factor of the reference gas and the mixture both at the reference conditions, not the conditions of the mixture.

Examples

>>> Mixture('argon').SGg
1.3800407778218216
property SGl

Specific gravity of a hypothetical liquid phase of the mixture at the specified temperature and pressure, [dimensionless]. The reference condition is water at 4 °C and 1 atm (rho=999.017 kg/m^3). For liquids, SG is defined that the reference chemical’s T and P are fixed, but the chemical itself varies with the specified T and P.

Examples

>>> Mixture('water', ws=[1], T=365).SGl
0.9650065522428539
property SGs

Specific gravity of a hypothetical solid phase of the mixture at the specified temperature and pressure, [dimensionless]. The reference condition is water at 4 °C and 1 atm (rho=999.017 kg/m^3). The SG varries with temperature and pressure but only very slightly.

T_default = 298.15
property Tbubble

Bubble point temperature of the mixture at its current pressure and composition, in units of [K].

This property requires that the property package of the mixture found a solution to the given state variables.

property Tdew

Dew point temperature of the mixture at its current pressure and composition, in units of [K].

This property requires that the property package of the mixture found a solution to the given state variables.

property U

Internal energy of the mixture at its current state, in units of [J/kg].

This property requires that the property package of the mixture found a solution to the given state variables. It also depends on the molar volume of the mixture at its current conditions.

property UNIFAC_Dortmund_groups

List of dictionaries of Dortmund UNIFAC subgroup: count groups for each chemcial in the mixture. Uses the Dortmund UNIFAC subgroups, as determined by DDBST’s online service.

Examples

>>> Mixture(['1-pentanol', 'decane'], ws=[0.5, 0.5]).UNIFAC_Dortmund_groups
[{1: 1, 2: 4, 14: 1}, {1: 2, 2: 8}]
property UNIFAC_Qs

UNIFAC Q (normalized Van der Waals area) values, dimensionless. Used in the UNIFAC model.

Examples

>>> Mixture(['o-xylene', 'decane'], zs=[.5, .5]).UNIFAC_Qs
[3.536, 6.016]
property UNIFAC_Rs

UNIFAC R (normalized Van der Waals volume) values, dimensionless. Used in the UNIFAC model.

Examples

>>> Mixture(['o-xylene', 'm-xylene'], zs=[.5, .5]).UNIFAC_Rs
[4.6578, 4.6578]
property UNIFAC_groups

List of dictionaries of UNIFAC subgroup: count groups for each chemical in the mixture. Uses the original UNIFAC subgroups, as determined by DDBST’s online service.

Examples

>>> Mixture(['1-pentanol', 'decane'], ws=[0.5, 0.5]).UNIFAC_groups
[{1: 1, 2: 4, 14: 1}, {1: 2, 2: 8}]
property Um

Internal energy of the mixture at its current state, in units of [J/mol].

This property requires that the property package of the mixture found a solution to the given state variables. It also depends on the molar volume of the mixture at its current conditions.

V_over_F = None
property Van_der_Waals_areas

List of unnormalized Van der Waals areas of all the chemicals in the mixture, in units of [m^2/mol].

Examples

>>> Mixture(['1-pentanol', 'decane'], ws=[0.5, 0.5]).Van_der_Waals_areas
[1052000.0, 1504000.0]
property Van_der_Waals_volumes

List of unnormalized Van der Waals volumes of all the chemicals in the mixture, in units of [m^3/mol].

Examples

>>> Mixture(['1-pentanol', 'decane'], ws=[0.5, 0.5]).Van_der_Waals_volumes
[6.9762279e-05, 0.00010918455800000001]
Vfgs(T=None, P=None)[source]

Volume fractions of all species in a hypothetical pure-gas phase at the current or specified temperature and pressure. If temperature or pressure are specified, the non-specified property is assumed to be that of the mixture. Note this is a method, not a property. Volume fractions are calculated based on pure species volumes only.

Examples

>>> Mixture(['sulfur hexafluoride', 'methane'], zs=[.2, .9], T=315).Vfgs()
[0.18062059238682632, 0.8193794076131737]
>>> S = Mixture(['sulfur hexafluoride', 'methane'], zs=[.1, .9])
>>> S.Vfgs(P=1E2)
[0.0999987466608421, 0.9000012533391578]
Vfls(T=None, P=None)[source]

Volume fractions of all species in a hypothetical pure-liquid phase at the current or specified temperature and pressure. If temperature or pressure are specified, the non-specified property is assumed to be that of the mixture. Note this is a method, not a property. Volume fractions are calculated based on pure species volumes only.

Examples

>>> Mixture(['hexane', 'pentane'], zs=[.5, .5], T=315).Vfls()
[0.5299671144566751, 0.47003288554332484]
>>> S = Mixture(['hexane', 'decane'], zs=[0.25, 0.75])
>>> S.Vfls(298.16, 101326)
[0.18301434895886864, 0.8169856510411313]
property Vm

Molar volume of the mixture at its current phase and temperature and pressure, in units of [m^3/mol]. Available only if single phase.

Examples

>>> Mixture(['ethylbenzene'], ws=[1], T=550, P=3E6).Vm
0.00017758024401627633
property Vmg

Gas-phase molar volume of the mixture at its current temperature, pressure, and composition in units of [m^3/mol]. For calculation of this property at other temperatures or pressures or compositions, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.volume.VolumeGasMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['hexane'], ws=[1], T=300, P=2E5).Vmg
0.010888694235142216
property Vmg_STP

Gas-phase molar volume of the mixture at 298.15 K and 101.325 kPa, and the current composition in units of [m^3/mol].

Examples

>>> Mixture(['nitrogen'], ws=[1]).Vmg_STP
0.02445443688838904
property Vmgs

Pure component gas-phase molar volumes of the chemicals in the mixture at its current temperature and pressure, in units of [m^3/mol].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Vmgs
[0.024929001982294974, 0.024150186467130488]
property Vml

Liquid-phase molar volume of the mixture at its current temperature, pressure, and composition in units of [m^3/mol]. For calculation of this property at other temperatures or pressures or compositions, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.volume.VolumeLiquidMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['cyclobutane'], ws=[1], T=225).Vml
7.42395423425395e-05
property Vml_STP

Liquid-phase molar volume of the mixture at 298.15 K and 101.325 kPa, and the current composition in units of [m^3/mol].

Examples

>>> Mixture(['cyclobutane'], ws=[1]).Vml_STP
8.143327329133706e-05
property Vmls

Pure component liquid-phase molar volumes of the chemicals in the mixture at its current temperature and pressure, in units of [m^3/mol].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Vmls
[9.188896727673715e-05, 0.00010946199496993461]
Vms = None
property Vmss

Pure component solid-phase molar volumes of the chemicals in the mixture at its current temperature, in units of [m^3/mol].

Examples

>>> Mixture(['iron'], ws=[1], T=320).Vmss
[7.09593392630242e-06]
Weber(V=None, D=None)[source]
property Z

Compressibility factor of the mixture at its current phase and temperature and pressure, [dimensionless]. Available only if single phase.

Examples

>>> Mixture(['MTBE'], ws=[1], T=900, P=1E-2).Z
0.9999999999056374
property Zg

Compressibility factor of the mixture in the gas phase at the current temperature, pressure, and composition, [dimensionless].

Utilizes the object oriented interface and thermo.volume.VolumeGasMixture to perform the actual calculation of molar volume.

Examples

>>> Mixture(['hexane'], ws=[1], T=300, P=1E5).Zg
0.9403859376888885
property Zg_STP

Gas-phase compressibility factor of the mixture at 298.15 K and 101.325 kPa, and the current composition, [dimensionless].

Examples

>>> Mixture(['nitrogen'], ws=[1]).Zg_STP
0.9995520809691023
property Zgs

Pure component compressibility factors of the chemicals in the mixture in the gas phase at the current temperature and pressure, [dimensionless].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Zgs
[0.9493743379816593, 0.9197146081359057]
property Zl

Compressibility factor of the mixture in the liquid phase at the current temperature, pressure, and composition, [dimensionless].

Utilizes the object oriented interface and thermo.volume.VolumeLiquidMixture to perform the actual calculation of molar volume.

Examples

>>> Mixture(['water'], ws=[1]).Zl
0.0007385375470263454
property Zl_STP

Liquid-phase compressibility factor of the mixture at 298.15 K and 101.325 kPa, and the current composition, [dimensionless].

Examples

>>> Mixture(['cyclobutane'], ws=[1]).Zl_STP
0.0033285083663950068
property Zls

Pure component compressibility factors of the chemicals in the liquid phase at the current temperature and pressure, [dimensionless].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).Zls
[0.0034994191720201235, 0.004168655010037687]
property Zss

Pure component compressibility factors of the chemicals in the mixture in the solid phase at the current temperature and pressure, [dimensionless].

Examples

>>> Mixture(['palladium'], ws=[1]).Zss
[0.00036248477437931853]
property alpha

Thermal diffusivity of the mixture at its current temperature, pressure, and phase in units of [m^2/s]. Available only if single phase.

α=kρCp\alpha = \frac{k}{\rho Cp}

Examples

>>> Mixture(['furfural'], ws=[1]).alpha
8.696537158635412e-08
property alphag

Thermal diffusivity of the gas phase of the mixture if one exists at its current temperature and pressure, in units of [m^2/s].

α=kρCp\alpha = \frac{k}{\rho Cp}

Examples

>>> Mixture(['ammonia'], ws=[1]).alphag
1.6968517002221566e-05
property alphags

Pure component thermal diffusivities of the chemicals in the mixture in the gas phase at the current temperature and pressure, in units of [m^2/s].

α=kρCp\alpha = \frac{k}{\rho Cp}

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).alphags
[3.3028044028118324e-06, 2.4412958544059014e-06]
property alphal

Thermal diffusivity of the liquid phase of the mixture if one exists at its current temperature and pressure, in units of [m^2/s].

α=kρCp\alpha = \frac{k}{\rho Cp}

Examples

>>> Mixture(['nitrogen'], ws=[1], T=70).alphal
9.444949636299626e-08
property alphals

Pure component thermal diffusivities of the chemicals in the mixture in the liquid phase at the current temperature and pressure, in units of [m^2/s].

α=kρCp\alpha = \frac{k}{\rho Cp}

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).alphals
[8.732683564481583e-08, 7.57355434073289e-08]
property atom_fractions

Dictionary of atomic fractions for each atom in the mixture.

Examples

>>> Mixture(['CO2', 'O2'], zs=[0.5, 0.5]).atom_fractions
{'C': 0.2, 'O': 0.8}
property atom_fractionss

List of dictionaries of atomic fractions for all chemicals in the mixture.

Examples

>>> Mixture(['oxygen', 'nitrogen'], zs=[.5, .5]).atom_fractionss
[{'O': 1.0}, {'N': 1.0}]
property atoms

Mole-averaged dictionary of atom counts for all atoms of the chemicals in the mixture.

Examples

>>> Mixture(['nitrogen', 'oxygen'], zs=[.01, .99]).atoms
{'O': 1.98, 'N': 0.02}
property atomss

List of dictionaries of atom counts for all chemicals in the mixture.

Examples

>>> Mixture(['nitrogen', 'oxygen'], zs=[.01, .99]).atomss
[{'N': 2}, {'O': 2}]
autoflash = True
property charge_balance

Charge imbalance of the mixture, in units of [faraday]. Mixtures meeting the electroneutrality condition will have an imbalance of 0.

Examples

>>> Mixture(['Na+', 'Cl-', 'water'], zs=[.01, .01, .98]).charge_balance
0.0
property charges

Charges for all chemicals in the mixture, [faraday].

Examples

>>> Mixture(['water', 'sodium ion', 'chloride ion'], zs=[.9, .05, .05]).charges
[0, 1, -1]
compound_index(CAS)[source]
conductivity = None
property constants

Returns a :obj:`thermo.chemical_package.ChemicalConstantsPackage instance with constants from the mixture, [-].

draw_2d(Hs=False)[source]

Interface for drawing a 2D image of all the molecules in the mixture. Requires an HTML5 browser, and the libraries RDKit and IPython. An exception is raised if either of these libraries is absent.

Parameters
Hsbool

Whether or not to show hydrogen

Examples

Mixture([‘natural gas’]).draw_2d()

property economic_statuses

List of dictionaries of the economic status for all chemicals in the mixture.

Examples

>>> Mixture(['o-xylene', 'm-xylene'], zs=[.5, .5]).economic_statuses
[["US public: {'Manufactured': 0.0, 'Imported': 0.0, 'Exported': 0.0}",
  u'100,000 - 1,000,000 tonnes per annum',
  'OECD HPV Chemicals'],
 ["US public: {'Manufactured': 39.805, 'Imported': 0.0, 'Exported': 0.0}",
  u'100,000 - 1,000,000 tonnes per annum',
  'OECD HPV Chemicals']]
property eos

Equation of state object held by the mixture. See : obj:thermo.eos_mix for a full listing.

eos_in_a_box = []
eos_pures(eos=<class 'thermo.eos.PR'>, T=None, P=None)[source]
flash_caloric(T=None, P=None, VF=None, Hm=None, Sm=None, H=None, S=None)[source]
flashed = True
property formulas

Chemical formulas for all chemicals in the mixture.

Examples

>>> Mixture(['ethanol', 'trichloroethylene', 'furfuryl alcohol'],
... ws=[0.5, 0.2, 0.3]).formulas
['C2H6O', 'C2HCl3', 'C5H6O2']
property isentropic_exponent

Gas-phase ideal-gas isentropic exponent of the mixture at its current temperature, [dimensionless]. Does not include pressure-compensation from an equation of state.

Examples

>>> Mixture(['hydrogen'], ws=[1]).isentropic_exponent
1.405237786321222
property isentropic_exponents

Gas-phase pure component ideal-gas isentropic exponent of the chemicals in the mixture at its current temperature, [dimensionless].

Does not include pressure-compensation from an equation of state.

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).isentropic_exponents
[1.1023398979313739, 1.080418846592871]
property isobaric_expansion

Isobaric (constant-pressure) expansion of the mixture at its current phase, temperature, and pressure in units of [1/K]. Available only if single phase.

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

Examples

>>> Mixture(['water'], ws=[1], T=647.1, P=22048320.0).isobaric_expansion
0.34074205839222449
property isobaric_expansion_g

Isobaric (constant-pressure) expansion of the gas phase of the mixture at its current temperature and pressure, in units of [1/K]. Available only if single phase.

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

Examples

>>> Mixture(['argon'], ws=[1], T=647.1, P=22048320.0).isobaric_expansion_g
0.0015661100323025273
property isobaric_expansion_gs

Pure component isobaric (constant-pressure) expansions of the chemicals in the mixture in the gas phase at its current temperature and pressure, in units of [1/K].

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

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).isobaric_expansion_gs
[0.0038091518363900499, 0.0043556759306508453]
property isobaric_expansion_l

Isobaric (constant-pressure) expansion of the liquid phase of the mixture at its current temperature and pressure, in units of [1/K]. Available only if single phase.

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

Examples

>>> Mixture(['argon'], ws=[1], T=647.1, P=22048320.0).isobaric_expansion_l
0.001859152875154442
property isobaric_expansion_ls

Pure component isobaric (constant-pressure) expansions of the chemicals in the mixture in the liquid phase at its current temperature and pressure, in units of [1/K].

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

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).isobaric_expansion_ls
[0.0012736035771253886, 0.0011234157437069571]
property k

Thermal conductivity of the mixture at its current phase, temperature, and pressure in units of [W/m/K]. Available only if single phase.

Examples

>>> Mixture(['ethanol'], ws=[1], T=300).kl
0.16313594741877802
property kg

Thermal conductivity of the mixture in the gas phase at its current temperature, pressure, and composition in units of [Pa*s].

For calculation of this property at other temperatures and pressures, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.thermal_conductivity.ThermalConductivityGasMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['water'], ws=[1], T=500).kg
0.036035173297862676
property kgs

Pure component thermal conductivies of the chemicals in the mixture in the gas phase at its current temperature and pressure, in units of [W/m/K].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).kgs
[0.011865404482987936, 0.010981336502491088]
property kl

Thermal conductivity of the mixture in the liquid phase at its current temperature, pressure, and composition in units of [Pa*s].

For calculation of this property at other temperatures and pressures, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.thermal_conductivity.ThermalConductivityLiquidMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['water'], ws=[1], T=320).kl
0.6369957248212118
property kls

Pure component thermal conductivities of the chemicals in the mixture in the liquid phase at its current temperature and pressure, in units of [W/m/K].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).kls
[0.13391538485205587, 0.12429339088930591]
ks = None
property legal_statuses

List of dictionaries of the legal status for all chemicals in the mixture.

Examples

>>> Mixture(['oxygen', 'nitrogen'], zs=[.5, .5]).legal_statuses
[{'DSL': 'LISTED',
  'EINECS': 'LISTED',
  'NLP': 'UNLISTED',
  'SPIN': 'LISTED',
  'TSCA': 'LISTED'},
 {'DSL': 'LISTED',
  'EINECS': 'LISTED',
  'NLP': 'UNLISTED',
  'SPIN': 'LISTED',
  'TSCA': 'LISTED'}]
property mass_fractions

Dictionary of mass fractions for each atom in the mixture.

Examples

>>> Mixture(['CO2', 'O2'], zs=[0.5, 0.5]).mass_fractions
{'C': 0.15801826905745822, 'O': 0.8419817309425419}
property mass_fractionss

List of dictionaries of mass fractions for all chemicals in the mixture.

Examples

>>> Mixture(['oxygen', 'nitrogen'], zs=[.5, .5]).mass_fractionss
[{'O': 1.0}, {'N': 1.0}]
property mu

Viscosity of the mixture at its current phase, temperature, and pressure in units of [Pa*s]. Available only if single phase.

Examples

>>> Mixture(['ethanol'], ws=[1], T=400).mu
1.1853097849748213e-05
property mug

Viscosity of the mixture in the gas phase at its current temperature, pressure, and composition in units of [Pa*s].

For calculation of this property at other temperatures and pressures, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.viscosity.ViscosityGasMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['water'], ws=[1], T=500).mug
1.7298722343367148e-05
property mugs

Pure component viscosities of the chemicals in the mixture in the gas phase at its current temperature and pressure, in units of [Pa*s].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).mugs
[8.082880451060605e-06, 7.442602145854158e-06]
property mul

Viscosity of the mixture in the liquid phase at its current temperature, pressure, and composition in units of [Pa*s].

For calculation of this property at other temperatures and pressures, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.viscosity.ViscosityLiquidMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['water'], ws=[1], T=320).mul
0.0005767262693751547
property muls

Pure component viscosities of the chemicals in the mixture in the liquid phase at its current temperature and pressure, in units of [Pa*s].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).muls
[0.00045545522798131764, 0.00043274394349114754]
property nu

Kinematic viscosity of the the mixture at its current temperature, pressure, and phase in units of [m^2/s]. Available only if single phase.

ν=μρ\nu = \frac{\mu}{\rho}

Examples

>>> Mixture(['argon'], ws=[1]).nu
1.3842643382482236e-05
property nug

Kinematic viscosity of the gas phase of the mixture if one exists at its current temperature and pressure, in units of [m^2/s].

ν=μρ\nu = \frac{\mu}{\rho}

Examples

>>> Mixture(['methane'], ws=[1], T=115).nug
2.5118460023343146e-06
property nugs

Pure component kinematic viscosities of the gas phase of the chemicals in the mixture at its current temperature and pressure, in units of [m^2/s].

ν=μρ\nu = \frac{\mu}{\rho}

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).nugs
[5.357870271650772e-07, 3.8127962283230277e-07]
property nul

Kinematic viscosity of the liquid phase of the mixture if one exists at its current temperature and pressure, in units of [m^2/s].

ν=μρ\nu = \frac{\mu}{\rho}

Examples

>>> Mixture(['methane'], ws=[1], T=110).nul
2.858088468937333e-07
property nuls

Pure component kinematic viscosities of the liquid phase of the chemicals in the mixture at its current temperature and pressure, in units of [m^2/s].

ν=μρ\nu = \frac{\mu}{\rho}

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).nuls
[5.357870271650772e-07, 3.8127962283230277e-07]
property permittivites

Pure component relative permittivities of the chemicals in the mixture at its current temperature, [dimensionless].

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).permittivites
[2.23133472, 1.8508128]
phase = None
properties(copy_pures=True, copy_mixtures=True)[source]
property rho

Mass density of the mixture at its current phase and temperature and pressure, in units of [kg/m^3]. Available only if single phase.

Examples

>>> Mixture(['decane'], ws=[1], T=550, P=2E6).rho
498.67008448640604
property rhog

Gas-phase mass density of the mixture at its current temperature, pressure, and composition in units of [kg/m^3]. For calculation of this property at other temperatures, pressures, or compositions or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.volume.VolumeGasMixture; each Mixture instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> Mixture(['hexane'], ws=[1], T=300, P=2E5).rhog
7.914447603999089
property rhog_STP

Gas-phase mass density of the mixture at 298.15 K and 101.325 kPa, and the current composition in units of [kg/m^3].

Examples

>>> Mixture(['nitrogen'], ws=[1]).rhog_STP
1.145534453639403
property rhogm

Molar density of the mixture in the gas phase at the current temperature, pressure, and composition in units of [mol/m^3].

Utilizes the object oriented interface and thermo.volume.VolumeGasMixture to perform the actual calculation of molar volume.

Examples

>>> Mixture(['water'], ws=[1], T=500).rhogm
24.467426039789093
property rhogm_STP

Molar density of the mixture in the gas phase at 298.15 K and 101.325 kPa, and the current composition, in units of [mol/m^3].

Examples

>>> Mixture(['nitrogen'], ws=[1]).rhogm_STP
40.892374850585895
property rhogms

Pure component molar densities of the chemicals in the gas phase at the current temperature and pressure, in units of [mol/m^3].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).rhogms
[40.11392035309789, 41.407547778608084]
property rhogs

Pure-component gas-phase mass densities of the chemicals in the mixture at its current temperature and pressure, in units of [kg/m^3].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).rhogs
[3.1333721283939258, 3.8152260283954584]
property rhol

Liquid-phase mass density of the mixture at its current temperature, pressure, and composition in units of [kg/m^3]. For calculation of this property at other temperatures, pressures, compositions or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.volume.VolumeLiquidMixture; each Mixture instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> Mixture(['o-xylene'], ws=[1], T=297).rhol
876.9946785618097
property rhol_STP

Liquid-phase mass density of the mixture at 298.15 K and 101.325 kPa, and the current composition in units of [kg/m^3].

Examples

>>> Mixture(['cyclobutane'], ws=[1]).rhol_STP
688.9851989526821
property rholm

Molar density of the mixture in the liquid phase at the current temperature, pressure, and composition in units of [mol/m^3].

Utilizes the object oriented interface and thermo.volume.VolumeLiquidMixture to perform the actual calculation of molar volume.

Examples

>>> Mixture(['water'], ws=[1], T=300).rholm
55317.352773503124
property rholm_STP

Molar density of the mixture in the liquid phase at 298.15 K and 101.325 kPa, and the current composition, in units of [mol/m^3].

Examples

>>> Mixture(['water'], ws=[1]).rholm_STP
55344.59086372442
property rholms

Pure component molar densities of the chemicals in the mixture in the liquid phase at the current temperature and pressure, in units of [mol/m^3].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).rholms
[10882.699301520635, 9135.590853014008]
property rhols

Pure-component liquid-phase mass density of the chemicals in the mixture at its current temperature and pressure, in units of [kg/m^3].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).rhols
[850.0676666084917, 841.7389069631628]
property rhom

Molar density of the mixture at its current phase and temperature and pressure, in units of [mol/m^3]. Available only if single phase.

Examples

>>> Mixture(['1-hexanol'], ws=[1]).rhom
7983.414573003429
rhos = None
property rhosms

Pure component molar densities of the chemicals in the solid phase at the current temperature and pressure, in units of [mol/m^3].

Examples

>>> Mixture(['iron'], ws=[1], T=320).rhosms
[140925.7767033753]
property rhoss

Pure component solid-phase mass density of the chemicals in the mixture at its current temperature, in units of [kg/m^3].

Examples

>>> Mixture(['iron'], ws=[1], T=320).rhoss
[7869.999999999994]
property ringss

List of ring counts for all chemicals in the mixture.

Examples

>>> Mixture(['Docetaxel', 'Paclitaxel'], zs=[.5, .5]).ringss
[6, 7]
set_Chemical_property_objects()[source]
set_TP_sources()[source]
set_chemical_TP(T=None, P=None)[source]

Basic method to change all chemical instances to be at the T and P specified. If they are not specified, the the values of the mixture will be used. This is not necessary for using the Mixture instance unless values specified to chemicals are required.

set_chemical_constants()[source]

Basic method which retrieves and sets constants of chemicals to be accessible as lists from a Mixture object. This gets called automatically on the instantiation of a new Mixture instance.

set_constant_sources()[source]
set_constants()[source]
set_eos(T, P, eos=<class 'thermo.eos_mix.PRMIX'>)[source]
set_property_package(pkg=None)[source]
property sigma

Surface tension of the mixture at its current temperature and composition, in units of [N/m].

For calculation of this property at other temperatures, or specifying manually the method used to calculate it, and more - see the object oriented interface thermo.interface.SurfaceTensionMixture; each Mixture instance creates one to actually perform the calculations.

Examples

>>> Mixture(['water'], ws=[1], T=300, P=1E5).sigma
0.07176932405246211
property sigmas

Pure component surface tensions of the chemicals in the mixture at its current temperature, in units of [N/m].

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5], T=320).sigmas
[0.02533469712937521, 0.025254723406585546]
property similarity_variables

Similarity variables for all chemicals in the mixture, see chemicals.elements.similarity_variable for the definition, [mol/g]

Examples

>>> Mixture(['benzene', 'toluene'], ws=[0.5, 0.5]).similarity_variables
[0.15362587797189262, 0.16279853724428964]
property smiless

SMILES strings for all chemicals in the mixture.

Examples

>>> Mixture(['methane', 'ethane', 'propane', 'butane'],
... zs=[0.25, 0.25, 0.25, 0.25]).smiless
['C', 'CC', 'CCC', 'CCCC']
property solubility_parameters

Pure component solubility parameters of the chemicals in the mixture at its current temperature and pressure, in units of [Pa^0.5].

δ=ΔHvapRTVm\delta = \sqrt{\frac{\Delta H_{vap} - RT}{V_m}}

Examples

>>> Mixture(['benzene', 'hexane'], ws=[0.5, 0.5], T=320).solubility_parameters
[18062.51359608708, 14244.12852702228]
property speed_of_sound

Bulk speed of sound of the mixture at its current temperature, [m/s].

Examples

>>> Mixture(['toluene'], P=1E5, VF=0.5, ws=[1]).speed_of_sound
478.99527258140211
property speed_of_sound_g

Gas-phase speed of sound of the mixture at its current temperature, [m/s].

Examples

>>> Mixture(['nitrogen'], ws=[1]).speed_of_sound_g
351.77445481641661
property speed_of_sound_l

Liquid-phase speed of sound of the mixture at its current temperature, [m/s].

Examples

>>> Mixture(['toluene'], P=1E5, T=300, ws=[1]).speed_of_sound_l
1116.0852487852942
property synonymss

Lists of synonyms for all chemicals in the mixture.

Examples

>>> Mixture(['Tetradecene', 'Pentadecene'], zs=[.1, .9]).synonymss
[['tetradec-2-ene', 'tetradecene', '2-tetradecene', 'tetradec-2-ene', '26952-13-6', '35953-53-8', '1652-97-7'], ['pentadec-1-ene', '1-pentadecene', 'pentadecene,1-', 'pentadec-1-ene', '13360-61-7', 'pentadecene']]
xs = None
ys = None