Legacy Chemicals (thermo.chemical)

class thermo.chemical.Chemical(ID, T=298.15, P=101325, autocalc=True)[source]

Bases: object

Creates a Chemical 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 temperature and pressure.

Parameters
IDstr
One of the following [-]:
  • Name, in IUPAC form or common form or a synonym registered in PubChem

  • InChI name, prefixed by ‘InChI=1S/’ or ‘InChI=1/’

  • InChI key, prefixed by ‘InChIKey=’

  • PubChem CID, prefixed by ‘PubChem=’

  • SMILES (prefix with ‘SMILES=’ to ensure smiles parsing)

  • CAS number

Tfloat, optional

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

Pfloat, optional

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

Notes

Warning

The Chemical 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 chemical objects:

>>> Chemical('hexane')
<Chemical [hexane], T=298.15 K, P=101325 Pa>
>>> Chemical('CCCCCCCC', T=500, P=1E7)
<Chemical [octane], T=500.00 K, P=10000000 Pa>
>>> Chemical('7440-36-0', P=1000)
<Chemical [antimony], T=298.15 K, P=1000 Pa>

Getting basic properties:

>>> N2 = Chemical('Nitrogen')
>>> N2.Tm, N2.Tb, N2.Tc # melting, boiling, and critical points [K]
(63.15, 77.355, 126.2)
>>> N2.Pt, N2.Pc # sublimation and critical pressure [Pa]
(12526.9697368421, 3394387.5)
>>> N2.CAS, N2.formula, N2.InChI, N2.smiles, N2.atoms # CAS number, formula, InChI string, smiles string, dictionary of atomic elements and their count
('7727-37-9', 'N2', 'N2/c1-2', 'N#N', {'N': 2})

Changing the T/P of the chemical, and gettign temperature-dependent properties:

>>> N2.Cp, N2.rho, N2.mu # Heat capacity [J/kg/K], density [kg/m^3], viscosity [Pa*s]
(1039., 1.14, 1.78e-05)
>>> N2.calculate(T=65, P=1E6) # set it to a liquid at 65 K and 1 MPa
>>> N2.phase
'l'
>>> N2.Cp, N2.rho, N2.mu # properties are now of the liquid phase
(2002., 861., 0.000285)

Molar units are also available for properties:

>>> N2.Cpm, N2.Vm, N2.Hvapm # heat capacity [J/mol/K], molar volume [m^3/mol], enthalpy of vaporization [J/mol]
(56., 3.25e-05, 5982.)

A great deal of properties are available; for a complete list look at the attributes list.

>>> N2.alpha, N2.JT # thermal diffusivity [m^2/s], Joule-Thompson coefficient [K/Pa]
(9.87e-08, -4.0e-07)
>>> N2.isentropic_exponent, N2.isobaric_expansion
(1.4, 0.0047)

For pure species, the phase is easily identified, allowing for properties to be obtained without needing to specify the phase. However, the properties are also available in the hypothetical gas phase (when under the boiling point) and in the hypothetical liquid phase (when above the boiling point) as these properties are needed to evaluate mixture properties. Specify the phase of a property to be retrieved by appending ‘l’ or ‘g’ or ‘s’ to the property.

>>> C50 = Chemical('pentacontane')
>>> C50.rhog, C50.Cpg, C50.kg, C50.mug
(4.241646701894199, 1126.5533755283168, 0.00941385692301755, 6.973325939594919e-06)

Temperature dependent properties are calculated by objects which provide many useful features related to the properties. To determine the temperature at which nitrogen has a saturation pressure of 1 MPa:

>>> N2.VaporPressure.solve_property(1E6)
103.73528598652341

To compute an integral of the ideal-gas heat capacity of nitrogen to determine the enthalpy required for a given change in temperature. Note the thermodynamic objects calculate values in molar units always.

>>> N2.HeatCapacityGas.T_dependent_property_integral(100, 120) # J/mol/K
582.0121860897898

Derivatives of properties can be calculated as well, as may be needed by for example heat transfer calculations:

>>> N2.SurfaceTension.T_dependent_property_derivative(77)
-0.00022695346296730534

If a property is needed at multiple temperatures or pressures, it is faster to use the object directly to perform the calculation rather than setting the conditions for the chemical.

>>> [N2.VaporPressure(T) for T in range(80, 120, 10)]
[136979.4840843189, 360712.5746603142, 778846.276691705, 1466996.7208525643]

These objects are also how the methods by which the properties are calculated can be changed. To see the available methods for a property:

>>> N2.VaporPressure.all_methods
set(['VDI_PPDS', 'BOILING_CRITICAL', 'WAGNER_MCGARRY', 'AMBROSE_WALTON', 'COOLPROP', 'LEE_KESLER_PSAT', 'EOS', 'ANTOINE_POLING', 'SANJARI', 'DIPPR_PERRY_8E', 'Edalat', 'WAGNER_POLING'])

To specify the method which should be used for calculations of a property. In the example below, the Lee-kesler correlation for vapor pressure is specified.

>>> N2.calculate(80)
>>> N2.Psat
136979.4840843189
>>> N2.VaporPressure.method = 'LEE_KESLER_PSAT'
>>> N2.Psat
134987.76815364443

Properties may also be plotted via these objects:

>>> N2.VaporPressure.plot_T_dependent_property() 
>>> N2.VolumeLiquid.plot_isotherm(T=77, Pmin=1E5, Pmax=1E7) 
>>> N2.VolumeLiquid.plot_isobar(P=1E6,  Tmin=66, Tmax=120) 
>>> N2.VolumeLiquid.plot_TP_dependent_property(Tmin=60, Tmax=100,  Pmin=1E5, Pmax=1E7) 
Attributes
Tfloat

Temperature of the chemical, [K]

Pfloat

Pressure of the chemical, [Pa]

phasestr

Phase of the chemical; one of ‘s’, ‘l’, ‘g’, or ‘l/g’.

IDstr

User specified string by which the chemical’s CAS was looked up.

CASstr

The CAS number of the chemical.

PubChemint

PubChem Compound identifier (CID) of the chemical; all chemicals are sourced from their database. Chemicals can be looked at online at https://pubchem.ncbi.nlm.nih.gov.

MWfloat

Molecular weight of the compound, [g/mol]

formulastr

Molecular formula of the compound.

atomsdict

dictionary of counts of individual atoms, indexed by symbol with proper capitalization, [-]

similarity_variablefloat

Similarity variable, see chemicals.elements.similarity_variable for the definition, [mol/g]

smilesstr

Simplified molecular-input line-entry system representation of the compound.

InChIstr

IUPAC International Chemical Identifier of the compound.

InChI_Keystr

25-character hash of the compound’s InChI.

IUPAC_namestr

Preferred IUPAC name for a compound.

synonymslist of strings

All synonyms for the compound found in PubChem, sorted by popularity.

Tmfloat

Melting temperature [K]

Tbfloat

Boiling temperature [K]

Tcfloat

Critical temperature [K]

Pcfloat

Critical pressure [Pa]

Vcfloat

Critical volume [m^3/mol]

Zcfloat

Critical compressibility [-]

rhocfloat

Critical density [kg/m^3]

rhocmfloat

Critical molar density [mol/m^3]

omegafloat

Acentric factor [-]

StielPolarfloat

Stiel Polar factor, see chemicals.acentric.Stiel_polar_factor for the definition [-]

Ttfloat

Triple temperature, [K]

Ptfloat

Triple pressure, [Pa]

Hfusfloat

Enthalpy of fusion [J/kg]

Hfusmfloat

Molar enthalpy of fusion [J/mol]

Hsubfloat

Enthalpy of sublimation [J/kg]

Hsubmfloat

Molar enthalpy of sublimation [J/mol]

Hfmfloat

Standard state molar enthalpy of formation, [J/mol]

Hffloat

Standard enthalpy of formation in a mass basis, [J/kg]

Hfgmfloat

Ideal-gas molar enthalpy of formation, [J/mol]

Hfgfloat

Ideal-gas enthalpy of formation in a mass basis, [J/kg]

Hcmfloat

Molar higher heat of combustion [J/mol]

Hcfloat

Higher Heat of combustion [J/kg]

Hcm_lowerfloat

Molar lower heat of combustion [J/mol]

Hc_lowerfloat

Lower Heat of combustion [J/kg]

S0mfloat

Standard state absolute molar entropy of the chemical, [J/mol/K]

S0float

Standard state absolute entropy of the chemical, [J/kg/K]

S0gmfloat

Absolute molar entropy in an ideal gas state of the chemical, [J/mol/K]

S0gfloat

Absolute mass entropy in an ideal gas state of the chemical, [J/kg/K]

Gfmfloat

Standard state molar change of Gibbs energy of formation [J/mol]

Gffloat

Standard state change of Gibbs energy of formation [J/kg]

Gfgmfloat

Ideal-gas molar change of Gibbs energy of formation [J/mol]

Gfgfloat

Ideal-gas change of Gibbs energy of formation [J/kg]

Sfmfloat

Standard state molar change of entropy of formation, [J/mol/K]

Sffloat

Standard state change of entropy of formation, [J/kg/K]

Sfgmfloat

Ideal-gas molar change of entropy of formation, [J/mol/K]

Sfgfloat

Ideal-gas change of entropy of formation, [J/kg/K]

Hcgmfloat

Higher molar heat of combustion of the chemical in the ideal gas state, [J/mol]

Hcgfloat

Higher heat of combustion of the chemical in the ideal gas state, [J/kg]

Hcgm_lowerfloat

Lower molar heat of combustion of the chemical in the ideal gas state, [J/mol]

Hcg_lowerfloat

Lower heat of combustion of the chemical in the ideal gas state, [J/kg]

Tflashfloat

Flash point of the chemical, [K]

Tautoignitionfloat

Autoignition point of the chemical, [K]

LFLfloat

Lower flammability limit of the gas in an atmosphere at STP, mole fraction [-]

UFLfloat

Upper flammability limit of the gas in an atmosphere at STP, mole fraction [-]

TWAtuple[quantity, unit]

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

STELtuple[quantity, unit]

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

Ceilingtuple[quantity, unit]

Ceiling limits on worker exposure to dangerous chemicals.

Skinbool

Whether or not a chemical can be absorbed through the skin.

Carcinogenstr or dict

Carcinogen status information.

dipolefloat

Dipole moment in debye, [3.33564095198e-30 ampere*second^2]

Stockmayerfloat

Lennard-Jones depth of potential-energy minimum over k, [K]

molecular_diameterfloat

Lennard-Jones molecular diameter, [angstrom]

GWPfloat

Global warming potential (default 100-year outlook) (impact/mass chemical)/(impact/mass CO2), [-]

ODPfloat

Ozone Depletion potential (impact/mass chemical)/(impact/mass CFC-11), [-]

logPfloat

Octanol-water partition coefficient, [-]

legal_statusstr or dict

Dictionary of legal status indicators for the chemical.

economic_statuslist

Dictionary of economic status indicators for the chemical.

RIfloat

Refractive Index on the Na D line, [-]

RITfloat

Temperature at which refractive index reading was made

conductivityfloat

Electrical conductivity of the fluid, [S/m]

conductivityTfloat

Temperature at which conductivity measurement was made

VaporPressureobject

Instance of thermo.vapor_pressure.VaporPressure, with data and methods loaded for the chemical; performs the actual calculations of vapor pressure of the chemical.

EnthalpyVaporizationobject

Instance of thermo.phase_change.EnthalpyVaporization, with data and methods loaded for the chemical; performs the actual calculations of molar enthalpy of vaporization of the chemical.

VolumeSolidobject

Instance of thermo.volume.VolumeSolid, with data and methods loaded for the chemical; performs the actual calculations of molar volume of the solid phase of the chemical.

VolumeLiquidobject

Instance of thermo.volume.VolumeLiquid, with data and methods loaded for the chemical; performs the actual calculations of molar volume of the liquid phase of the chemical.

VolumeGasobject

Instance of thermo.volume.VolumeGas, with data and methods loaded for the chemical; performs the actual calculations of molar volume of the gas phase of the chemical.

HeatCapacitySolidobject

Instance of thermo.heat_capacity.HeatCapacitySolid, with data and methods loaded for the chemical; performs the actual calculations of molar heat capacity of the solid phase of the chemical.

HeatCapacityLiquidobject

Instance of thermo.heat_capacity.HeatCapacityLiquid, with data and methods loaded for the chemical; performs the actual calculations of molar heat capacity of the liquid phase of the chemical.

HeatCapacityGasobject

Instance of thermo.heat_capacity.HeatCapacityGas, with data and methods loaded for the chemical; performs the actual calculations of molar heat capacity of the gas phase of the chemical.

ViscosityLiquidobject

Instance of thermo.viscosity.ViscosityLiquid, with data and methods loaded for the chemical; performs the actual calculations of viscosity of the liquid phase of the chemical.

ViscosityGasobject

Instance of thermo.viscosity.ViscosityGas, with data and methods loaded for the chemical; performs the actual calculations of viscosity of the gas phase of the chemical.

ThermalConductivityLiquidobject

Instance of thermo.thermal_conductivity.ThermalConductivityLiquid, with data and methods loaded for the chemical; performs the actual calculations of thermal conductivity of the liquid phase of the chemical.

ThermalConductivityGasobject

Instance of thermo.thermal_conductivity.ThermalConductivityGas, with data and methods loaded for the chemical; performs the actual calculations of thermal conductivity of the gas phase of the chemical.

SurfaceTensionobject

Instance of thermo.interface.SurfaceTension, with data and methods loaded for the chemical; performs the actual calculations of surface tension of the chemical.

Permittivityobject

Instance of thermo.permittivity.PermittivityLiquid, with data and methods loaded for the chemical; performs the actual calculations of permittivity of the chemical.

Psat_298float

Vapor pressure of the chemical at 298.15 K, [Pa]

phase_STPstr

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

Vml_Tbfloat

Molar volume of liquid phase at the normal boiling point [m^3/mol]

Vml_Tmfloat

Molar volume of liquid phase at the melting point [m^3/mol]

Vml_STPfloat

Molar volume of liquid phase at 298.15 K and 101325 Pa [m^3/mol]

rhoml_STPfloat

Molar density of liquid phase at 298.15 K and 101325 Pa [mol/m^3]

Vmg_STPfloat

Molar volume of gas phase at 298.15 K and 101325 Pa according to the ideal gas law, [m^3/mol]

Vms_Tmfloat

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

rhos_Tmfloat

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

Hvap_Tbmfloat

Molar enthalpy of vaporization at the normal boiling point [J/mol]

Hvap_Tbfloat

Mass enthalpy of vaporization at the normal boiling point [J/kg]

Hvapm_298float

Molar enthalpy of vaporization at 298.15 K [J/mol]

Hvap_298float

Mass enthalpy of vaporization at 298.15 K [J/kg]

alpha

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

alphag

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

alphal

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

API

API gravity of the liquid phase of the chemical, [degrees].

aromatic_rings

Number of aromatic rings in a chemical, computed with RDKit from a chemical’s SMILES.

atom_fractions

Dictionary of atom:fractional occurence of the elements in a chemical.

Bvirial

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

charge

Charge of a chemical, computed with RDKit from a chemical’s SMILES.

Cp

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

Cpg

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

Cpgm

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

Cpl

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

Cplm

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

Cpm

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

Cps

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

Cpsm

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

Cvg

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

Cvgm

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

eos

Equation of state object held by the chemical; used to calculate excess thermodynamic quantities, and also provides a vapor pressure curve, enthalpy of vaporization curve, fugacity, thermodynamic partial derivatives, and more; see thermo.eos for a full listing.

Hill

Hill formula of a compound.

Hvap

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

Hvapm

Enthalpy of vaporization of the chemical at its current temperature, in units of [J/mol].

isentropic_exponent

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

isobaric_expansion

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

isobaric_expansion_g

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

isobaric_expansion_l

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

JT

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

JTg

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

JTl

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

k

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

kg

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

kl

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

mass_fractions

Dictionary of atom:mass-weighted fractional occurence of elements.

mu

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

mug

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

mul

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

nu

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

nug

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

nul

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

Parachor

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

permittivity

Relative permittivity (dielectric constant) of the chemical at its current temperature, [dimensionless].

Poynting

Poynting correction factor [dimensionless] for use in phase equilibria methods based on activity coefficients or other reference states.

Pr

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

Prg

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

Prl

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

Psat

Vapor pressure of the chemical at its current temperature, in units of [Pa].

PSRK_groups

Dictionary of PSRK subgroup: count groups for the PSRK subgroups, as determined by DDBST’s online service.

rdkitmol

RDKit object of the chemical, without hydrogen.

rdkitmol_Hs

RDKit object of the chemical, with hydrogen.

rho

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

rhog

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

rhogm

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

rhol

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

rholm

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

rhom

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

rhos

Solid-phase mass density of the chemical at its current temperature, in units of [kg/m^3].

rhosm

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

rings

Number of rings in a chemical, computed with RDKit from a chemical’s SMILES.

SG

Specific gravity of the chemical, [dimensionless].

SGg

Specific gravity of the gas phase of the chemical, [dimensionless].

SGl

Specific gravity of the liquid phase of the chemical at the specified temperature and pressure, [dimensionless].

SGs

Specific gravity of the solid phase of the chemical at the specified temperature and pressure, [dimensionless].

sigma

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

solubility_parameter

Solubility parameter of the chemical at its current temperature and pressure, in units of [Pa^0.5].

UNIFAC_Dortmund_groups

Dictionary of Dortmund UNIFAC subgroup: count groups for the Dortmund UNIFAC subgroups, as determined by DDBST’s online service.

UNIFAC_groups

Dictionary of UNIFAC subgroup: count groups for the original UNIFAC subgroups, as determined by DDBST’s online service.

UNIFAC_R

UNIFAC R (normalized Van der Waals volume), dimensionless.

UNIFAC_Q

UNIFAC Q (normalized Van der Waals area), dimensionless.

Van_der_Waals_area

Unnormalized Van der Waals area, in units of [m^2/mol].

Van_der_Waals_volume

Unnormalized Van der Waals volume, in units of [m^3/mol].

Vm

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

Vmg

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

Vml

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

Vms

Solid-phase molar volume of the chemical at its current temperature, in units of [m^3/mol].

Z

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

Zg

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

Zl

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

Zs

Compressibility factor of the chemical in the solid phase at the current temperature and pressure, [dimensionless].

Methods

Reynolds([V, D])

draw_2d([width, height, Hs])

Interface for drawing a 2D image of the molecule.

draw_3d([width, height, style, Hs, atom_labels])

Interface for drawing an interactive 3D view of the molecule.

Bond

Capillary

Grashof

Jakob

Peclet_heat

Tsat

Weber

calc_H

calc_H_excess

calc_S

calc_S_excess

calculate

calculate_PH

calculate_PS

calculate_TH

calculate_TS

set_TP_sources

set_constant_sources

set_constants

set_eos

set_ref

set_thermo

property A

Helmholtz energy of the chemical at its current temperature and pressure, in units of [J/kg].

This property requires that thermo.chemical.set_thermo ran successfully to be accurate. It also depends on the molar volume of the chemical at its current conditions.

property API

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

Examples

>>> Chemical('water').API
10.
property Am

Helmholtz energy of the chemical at its current temperature and pressure, in units of [J/mol].

This property requires that thermo.chemical.set_thermo ran successfully to be accurate. It also depends on the molar volume of the chemical at its current conditions.

Bond(L=None)[source]
property Bvirial

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

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

Examples

>>> Chemical('water', T=500, P=1e5).Bvirial
-0.00017
Capillary(V=None)[source]
property Cp

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

Utilizes the object oriented interfaces thermo.heat_capacity.HeatCapacitySolid, thermo.heat_capacity.HeatCapacityLiquid, and thermo.heat_capacity.HeatCapacityGas to perform the actual calculation of each property. Note that those interfaces provide output in molar units (J/mol/K).

Examples

>>> w = Chemical('water')
>>> w.Cp, w.phase
(4180., 'l')
>>> Chemical('palladium').Cp
234.
property Cpg

Gas-phase heat capacity of the chemical at its current temperature, in units of [J/kg/K]. 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.heat_capacity.HeatCapacityGas; each Chemical instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> w = Chemical('water', T=520)
>>> w.Cpg
1967.
property Cpgm

Gas-phase ideal gas heat capacity of the chemical at its current temperature, in units of [J/mol/K]. 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.heat_capacity.HeatCapacityGas; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('water').Cpgm
33.58
>>> Chemical('water').HeatCapacityGas.T_dependent_property(320)
33.67
>>> Chemical('water').HeatCapacityGas.T_dependent_property_integral(300, 320)
672.
property Cpl

Liquid-phase heat capacity of the chemical at its current temperature, in units of [J/kg/K]. 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.heat_capacity.HeatCapacityLiquid; each Chemical instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> Chemical('water', T=320).Cpl
4177.

Ideal entropy change of water from 280 K to 340 K, output converted back to mass-based units of J/kg/K.

>>> dSm = Chemical('water').HeatCapacityLiquid.T_dependent_property_integral_over_T(280, 340)
>>> property_molar_to_mass(dSm, Chemical('water').MW)
812.
property Cplm

Liquid-phase heat capacity of the chemical at its current temperature, in units of [J/mol/K]. 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.heat_capacity.HeatCapacityLiquid; each Chemical instance creates one to actually perform the calculations.

Notes

Some methods give heat capacity along the saturation line, some at 1 atm but only up to the normal boiling point, and some give heat capacity at 1 atm up to the normal boiling point and then along the saturation line. Real-liquid heat capacity is pressure dependent, but this interface is not.

Examples

>>> Chemical('water').Cplm
75.3
>>> Chemical('water').HeatCapacityLiquid.T_dependent_property(320)
75.2
>>> Chemical('water').HeatCapacityLiquid.T_dependent_property_integral(300, 320)
1505.
property Cpm

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

Utilizes the object oriented interfaces thermo.heat_capacity.HeatCapacitySolid, thermo.heat_capacity.HeatCapacityLiquid, and thermo.heat_capacity.HeatCapacityGas to perform the actual calculation of each property.

Examples

>>> Chemical('cubane').Cpm
137.
>>> Chemical('ethylbenzene', T=550, P=3E6).Cpm
294.
property Cps

Solid-phase heat capacity of the chemical at its current temperature, in units of [J/kg/K]. 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.heat_capacity.HeatCapacitySolid; each Chemical instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> Chemical('palladium', T=400).Cps
241.
>>> Pd = Chemical('palladium', T=400)
>>> Cpsms = [Pd.HeatCapacitySolid.T_dependent_property(T) for T in np.linspace(300, 500, 2)]
>>> [property_molar_to_mass(Cps, Pd.MW) for Cps in Cpsms]
[234., 248.]
property Cpsm

Solid-phase heat capacity of the chemical at its current temperature, in units of [J/mol/K]. 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.heat_capacity.HeatCapacitySolid; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('palladium').Cpsm
24.9
>>> Chemical('palladium').HeatCapacitySolid.T_dependent_property(320)
25.
>>> sorted(Chemical('palladium').HeatCapacitySolid.all_methods)
['CRCSTD', 'LASTOVKA_S', 'PERRY151']
property Cvg

Gas-phase ideal-gas contant-volume heat capacity of the chemical 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

>>> w = Chemical('water', T=520)
>>> w.Cvg
1506.
property Cvgm

Gas-phase ideal-gas contant-volume heat capacity of the chemical at its current temperature, 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

>>> w = Chemical('water', T=520)
>>> w.Cvgm
27.1
Grashof(Tw=None, L=None)[source]
property Hill

Hill formula of a compound. For a description of the Hill system, see chemicals.elements.atoms_to_Hill.

Examples

>>> Chemical('furfuryl alcohol').Hill
'C5H6O2'
property Hvap

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

This property uses the object-oriented interface thermo.phase_change.EnthalpyVaporization, but converts its results from molar to mass units.

Examples

>>> Chemical('water', T=320).Hvap
2389540.
property Hvapm

Enthalpy of vaporization of the chemical at its current temperature, in units of [J/mol]. 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.phase_change.EnthalpyVaporization; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('water', T=320).Hvapm
43048.
>>> Chemical('water').EnthalpyVaporization.T_dependent_property(320)
43048.
>>> sorted(Chemical('water').EnthalpyVaporization.all_methods)
['ALIBAKHSHI', 'CHEN', 'CLAPEYRON', 'COOLPROP', 'CRC_HVAP_298', 'CRC_HVAP_TB', 'DIPPR_PERRY_8E', 'LIU', 'MORGAN_KOBAYASHI', 'PITZER', 'RIEDEL', 'SIVARAMAN_MAGEE_KOBAYASHI', 'VDI_PPDS', 'VELASCO', 'VETERE']
property JT

Joule Thomson coefficient of the chemical at its current phase and temperature, 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

>>> Chemical('water').JT
-2.21e-07
property JTg

Joule Thomson coefficient of the chemical 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)

Utilizes the temperature-derivative method of thermo.volume.VolumeGas and the temperature-dependent heat capacity method thermo.heat_capacity.HeatCapacityGas to obtain the properties required for the actual calculation.

Examples

>>> Chemical('dodecane', T=400, P=1000).JTg
5.4e-05
property JTl

Joule Thomson coefficient of the chemical 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)

Utilizes the temperature-derivative method of thermo.volume.VolumeLiquid and the temperature-dependent heat capacity method thermo.heat_capacity.HeatCapacityLiquid to obtain the properties required for the actual calculation.

Examples

>>> Chemical('dodecane', T=400).JTl
-3.1e-07
Jakob(Tw=None)[source]
property PSRK_groups

Dictionary of PSRK subgroup: count groups for the PSRK subgroups, as determined by DDBST’s online service.

Examples

>>> Chemical('Cumene').PSRK_groups
{1: 2, 9: 5, 13: 1}
property Parachor

Parachor of the chemical 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 phase, and molecular weight. For uses of this property, see thermo.utils.Parachor.

The gas density is calculated using the ideal-gas law.

Examples

>>> Chemical('octane').Parachor
6.2e-05
Peclet_heat(V=None, D=None)[source]
property Poynting

Poynting correction factor [dimensionless] for use in phase equilibria methods based on activity coefficients or other reference states. Performs the shortcut calculation assuming molar volume is independent of pressure.

Poy=exp[Vl(PPsat)RT]\text{Poy} = \exp\left[\frac{V_l (P-P^{sat})}{RT}\right]

The full calculation normally returns values very close to the approximate ones. This property is defined in terms of pure components only.

Notes

The full equation shown below can be used as follows:

Poy=exp[PisatPVildPRT]\text{Poy} = \exp\left[\frac{\int_{P_i^{sat}}^P V_i^l dP}{RT}\right]
>>> from scipy.integrate import quad
>>> c = Chemical('pentane', T=300, P=2e6)
>>> exp(quad(lambda P : c.VolumeLiquid(c.T, P), c.Psat, c.P)[0]/R/c.T)
1.093

Examples

>>> Chemical('pentane', T=300, P=2e6).Poynting
1.09
property Pr

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

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

Examples

>>> Chemical('acetone').Pr
4.1
property Prg

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

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

Utilizes the temperature and pressure dependent object oriented interfaces thermo.viscosity.ViscosityGas, thermo.thermal_conductivity.ThermalConductivityGas, and thermo.heat_capacity.HeatCapacityGas to calculate the actual properties.

Examples

>>> Chemical('NH3').Prg
0.84
property Prl

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

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

Utilizes the temperature and pressure dependent object oriented interfaces thermo.viscosity.ViscosityLiquid, thermo.thermal_conductivity.ThermalConductivityLiquid, and thermo.heat_capacity.HeatCapacityLiquid to calculate the actual properties.

Examples

>>> Chemical('nitrogen', T=70).Prl
2.78
property Psat

Vapor pressure of the chemical at its current temperature, in units of [Pa]. 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.vapor_pressure.VaporPressure; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('water', T=320).Psat
10545.
>>> Chemical('water').VaporPressure.T_dependent_property(320)
10545.
>>> sorted(Chemical('water').VaporPressure.all_methods)
['AMBROSE_WALTON', 'ANTOINE_POLING', 'ANTOINE_WEBBOOK', 'BOILING_CRITICAL', 'COOLPROP', 'DIPPR_PERRY_8E', 'EDALAT', 'EOS', 'IAPWS', 'LEE_KESLER_PSAT', 'SANJARI', 'VDI_PPDS', 'WAGNER_MCGARRY']
property R_specific

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

Examples

>>> Chemical('water').R_specific
461.5
Reynolds(V=None, D=None)[source]
property SG

Specific gravity of the chemical, [dimensionless].

For gas-phase conditions, this is calculated at 15.6 °C (60 °F) and 1 atm for the chemical 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 chemical’s density at the currently specified conditions.

Examples

>>> Chemical('MTBE').SG
0.73
property SGg

Specific gravity of the gas phase of the chemical, [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 chemical both at the reference conditions, not the conditions of the chemical.

Examples

>>> Chemical('argon').SGg
1.3
property SGl

Specific gravity of the liquid phase of the chemical 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

>>> Chemical('water', T=365).SGl
0.965
property SGs

Specific gravity of the solid phase of the chemical 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.

Examples

>>> Chemical('iron').SGs
7.8
Tsat(P)[source]
property U

Internal energy of the chemical at its current temperature and pressure, in units of [J/kg].

This property requires that thermo.chemical.set_thermo ran successfully to be accurate. It also depends on the molar volume of the chemical at its current conditions.

property UNIFAC_Dortmund_groups

Dictionary of Dortmund UNIFAC subgroup: count groups for the Dortmund UNIFAC subgroups, as determined by DDBST’s online service.

Examples

>>> Chemical('Cumene').UNIFAC_Dortmund_groups
{1: 2, 9: 5, 13: 1}
property UNIFAC_Q

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

Examples

>>> Chemical('decane').UNIFAC_Q
6.016
property UNIFAC_R

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

Examples

>>> Chemical('benzene').UNIFAC_R
3.1878
property UNIFAC_groups

Dictionary of UNIFAC subgroup: count groups for the original UNIFAC subgroups, as determined by DDBST’s online service.

Examples

>>> Chemical('Cumene').UNIFAC_groups
{1: 2, 9: 5, 13: 1}
property Um

Internal energy of the chemical at its current temperature and pressure, in units of [J/mol].

This property requires that thermo.chemical.set_thermo ran successfully to be accurate. It also depends on the molar volume of the chemical at its current conditions.

property Van_der_Waals_area

Unnormalized Van der Waals area, in units of [m^2/mol].

Examples

>>> Chemical('hexane').Van_der_Waals_area
964000.0
property Van_der_Waals_volume

Unnormalized Van der Waals volume, in units of [m^3/mol].

Examples

>>> Chemical('hexane').Van_der_Waals_volume
6.8261966e-05
property Vm

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

Utilizes the object oriented interfaces thermo.volume.VolumeSolid, thermo.volume.VolumeLiquid, and thermo.volume.VolumeGas to perform the actual calculation of each property.

Examples

>>> Chemical('ethylbenzene', T=550, P=3E6).Vm
0.00017
property Vmg

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

Examples

Estimate the molar volume of the core of the sun, at 15 million K and 26.5 PetaPascals, assuming pure helium (actually 68% helium):

>>> Chemical('helium', T=15E6, P=26.5E15).Vmg
4.7e-07
property Vmg_ideal

Gas-phase molar volume of the chemical at its current temperature and pressure calculated with the ideal-gas law, in units of [m^3/mol].

Examples

>>> Chemical('helium', T=300.0, P=1e5).Vmg_ideal
0.0249433878544
property Vml

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

Examples

>>> Chemical('cyclobutane', T=225).Vml
7.42395423425395e-05
property Vms

Solid-phase molar volume of the chemical at its current temperature, in units of [m^3/mol]. 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.volume.VolumeSolid; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('iron').Vms
7.09e-06
Weber(V=None, D=None)[source]
property Z

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

Examples

>>> Chemical('MTBE', T=900, P=1E-2).Z
1.
property Zg

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

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

Examples

>>> Chemical('sulfur hexafluoride', T=700, P=1E9).Zg
11.14
property Zl

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

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

Examples

>>> Chemical('water').Zl
0.00073
property Zs

Compressibility factor of the chemical in the solid phase at the current temperature and pressure, [dimensionless].

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

Examples

>>> Chemical('palladium').Z
0.000362
property absolute_permittivity

Absolute permittivity of the chemical at its current temperature, in units of [farad/meter]. Those units are equivalent to ampere^2*second^4/kg/m^3.

Examples

>>> Chemical('water', T=293.15).absolute_permittivity
7.1e-10
property alpha

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

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

Examples

>>> Chemical('furfural').alpha
8.7e-08
property alphag

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

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

Utilizes the temperature and pressure dependent object oriented interfaces thermo.volume.VolumeGas, thermo.thermal_conductivity.ThermalConductivityGas, and thermo.heat_capacity.HeatCapacityGas to calculate the actual properties.

Examples

>>> Chemical('ammonia').alphag
1.69e-05
property alphal

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

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

Utilizes the temperature and pressure dependent object oriented interfaces thermo.volume.VolumeLiquid, thermo.thermal_conductivity.ThermalConductivityLiquid, and thermo.heat_capacity.HeatCapacityLiquid to calculate the actual properties.

Examples

>>> Chemical('nitrogen', T=70).alphal
9.4e-08
property aromatic_rings

Number of aromatic rings in a chemical, computed with RDKit from a chemical’s SMILES. If RDKit is not available, holds None.

Examples

>>> Chemical('Paclitaxel').aromatic_rings
3
property atom_fractions

Dictionary of atom:fractional occurence of the elements in a chemical. Useful when performing element balances. For mass-fraction occurences, see mass_fractions.

Examples

>>> Chemical('Ammonium aluminium sulfate').atom_fractions
{'Al': 0.0625, 'H': 0.25, 'N': 0.0625, 'O': 0.5, 'S': 0.125}
calc_H(T, P)[source]
calc_H_excess(T, P)[source]
calc_S(T, P)[source]
calc_S_excess(T, P)[source]
calculate(T=None, P=None)[source]
calculate_PH(P, H)[source]
calculate_PS(P, S)[source]
calculate_TH(T, H)[source]
calculate_TS(T, S)[source]
property charge

Charge of a chemical, computed with RDKit from a chemical’s SMILES. If RDKit is not available, holds None.

Examples

>>> Chemical('sodium ion').charge
1
draw_2d(width=300, height=300, Hs=False)[source]

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

Parameters
widthint

Number of pixels wide for the view

heightint

Number of pixels tall for the view

Hsbool

Whether or not to show hydrogen

Examples

>>> Chemical('decane').draw_2d() 
<PIL.PngImagePlugin.PngImageFile image mode=RGB size=300x300 at 0x...>
draw_3d(width=300, height=500, style='stick', Hs=True, atom_labels=True)[source]

Interface for drawing an interactive 3D view of the molecule. Requires an HTML5 browser, and the libraries RDKit, pymol3D, and IPython. An exception is raised if all three of these libraries are not installed.

Parameters
widthint

Number of pixels wide for the view, [pixels]

heightint

Number of pixels tall for the view, [pixels]

stylestr

One of ‘stick’, ‘line’, ‘cross’, or ‘sphere’, [-]

Hsbool

Whether or not to show hydrogen, [-]

atom_labelsbool

Whether or not to label the atoms, [-]

Examples

>>> Chemical('cubane').draw_3d()
None
property economic_status

Dictionary of economic status indicators for the chemical.

Examples

>>> Chemical('benzene').economic_status
["US public: {'Manufactured': 6165232.1, 'Imported': 463146.474, 'Exported': 271908.252}", '1,000,000 - 10,000,000 tonnes per annum', 'Intermediate Use Only', 'OECD HPV Chemicals']
property eos

Equation of state object held by the chemical; used to calculate excess thermodynamic quantities, and also provides a vapor pressure curve, enthalpy of vaporization curve, fugacity, thermodynamic partial derivatives, and more; see thermo.eos for a full listing.

Examples

>>> Chemical('methane').eos.V_g
0.0244
property is_acid

Method to return whether or not this chemical is in the category acid, [-]

property is_acyl_halide

Method to return whether or not this chemical is in the category acyl_halide, [-]

property is_alcohol

Method to return whether or not this chemical is in the category alcohol, [-]

property is_aldehyde

Method to return whether or not this chemical is in the category aldehyde, [-]

property is_alkane

Method to return whether or not this chemical is in the category alkane, [-]

property is_alkene

Method to return whether or not this chemical is in the category alkene, [-]

property is_alkylaluminium

Method to return whether or not this chemical is in the category alkylaluminium, [-]

property is_alkyllithium

Method to return whether or not this chemical is in the category alkyllithium, [-]

property is_alkylmagnesium_halide

Method to return whether or not this chemical is in the category alkylmagnesium_halide, [-]

property is_alkyne

Method to return whether or not this chemical is in the category alkyne, [-]

property is_amide

Method to return whether or not this chemical is in the category amide, [-]

property is_amidine

Method to return whether or not this chemical is in the category amidine, [-]

property is_amine

Method to return whether or not this chemical is in the category amine, [-]

property is_anhydride

Method to return whether or not this chemical is in the category anhydride, [-]

property is_aromatic

Method to return whether or not this chemical is in the category aromatic, [-]

property is_azide

Method to return whether or not this chemical is in the category azide, [-]

property is_azo

Method to return whether or not this chemical is in the category azo, [-]

property is_borinic_acid

Method to return whether or not this chemical is in the category borinic_acid, [-]

property is_borinic_ester

Method to return whether or not this chemical is in the category borinic_ester, [-]

property is_boronic_acid

Method to return whether or not this chemical is in the category boronic_acid, [-]

property is_boronic_ester

Method to return whether or not this chemical is in the category boronic_ester, [-]

property is_branched_alkane

Method to return whether or not this chemical is in the category branched_alkane, [-]

property is_bromoalkane

Method to return whether or not this chemical is in the category bromoalkane, [-]

property is_carbamate

Method to return whether or not this chemical is in the category carbamate, [-]

property is_carbodithio

Method to return whether or not this chemical is in the category carbodithio, [-]

property is_carbodithioic_acid

Method to return whether or not this chemical is in the category carbodithioic_acid, [-]

property is_carbonate

Method to return whether or not this chemical is in the category carbonate, [-]

property is_carbothioic_o_acid

Method to return whether or not this chemical is in the category carbothioic_o_acid, [-]

property is_carbothioic_s_acid

Method to return whether or not this chemical is in the category carbothioic_s_acid, [-]

property is_carboxylate

Method to return whether or not this chemical is in the category carboxylate, [-]

property is_carboxylic_acid

Method to return whether or not this chemical is in the category carboxylic_acid, [-]

property is_carboxylic_anhydride

Method to return whether or not this chemical is in the category carboxylic_anhydride, [-]

property is_chloroalkane

Method to return whether or not this chemical is in the category chloroalkane, [-]

property is_cyanate

Method to return whether or not this chemical is in the category cyanate, [-]

property is_cyanide

Method to return whether or not this chemical is in the category cyanide, [-]

property is_cycloalkane

Method to return whether or not this chemical is in the category cycloalkane, [-]

property is_disulfide

Method to return whether or not this chemical is in the category disulfide, [-]

property is_ester

Method to return whether or not this chemical is in the category ester, [-]

property is_ether

Method to return whether or not this chemical is in the category ether, [-]

property is_fluoroalkane

Method to return whether or not this chemical is in the category fluoroalkane, [-]

property is_haloalkane

Method to return whether or not this chemical is in the category haloalkane, [-]

property is_hydroperoxide

Method to return whether or not this chemical is in the category hydroperoxide, [-]

property is_imide

Method to return whether or not this chemical is in the category imide, [-]

property is_imine

Method to return whether or not this chemical is in the category imine, [-]

property is_inorganic

Method to return whether or not this chemical is in the category inorganic, [-]

property is_iodoalkane

Method to return whether or not this chemical is in the category iodoalkane, [-]

property is_is_hydrocarbon

Method to return whether or not this chemical is in the category is_hydrocarbon, [-]

property is_isocyanate

Method to return whether or not this chemical is in the category isocyanate, [-]

property is_isonitrile

Method to return whether or not this chemical is in the category isonitrile, [-]

property is_isothiocyanate

Method to return whether or not this chemical is in the category isothiocyanate, [-]

property is_ketone

Method to return whether or not this chemical is in the category ketone, [-]

property is_mercaptan

Method to return whether or not this chemical is in the category mercaptan, [-]

property is_methylenedioxy

Method to return whether or not this chemical is in the category methylenedioxy, [-]

property is_nitrate

Method to return whether or not this chemical is in the category nitrate, [-]

property is_nitrile

Method to return whether or not this chemical is in the category nitrile, [-]

property is_nitrite

Method to return whether or not this chemical is in the category nitrite, [-]

property is_nitro

Method to return whether or not this chemical is in the category nitro, [-]

property is_nitroso

Method to return whether or not this chemical is in the category nitroso, [-]

property is_organic

Method to return whether or not this chemical is in the category organic, [-]

property is_orthocarbonate_ester

Method to return whether or not this chemical is in the category orthocarbonate_ester, [-]

property is_orthoester

Method to return whether or not this chemical is in the category orthoester, [-]

property is_oxime

Method to return whether or not this chemical is in the category oxime, [-]

property is_peroxide

Method to return whether or not this chemical is in the category peroxide, [-]

property is_phenol

Method to return whether or not this chemical is in the category phenol, [-]

property is_phosphate

Method to return whether or not this chemical is in the category phosphate, [-]

property is_phosphine

Method to return whether or not this chemical is in the category phosphine, [-]

property is_phosphodiester

Method to return whether or not this chemical is in the category phosphodiester, [-]

property is_phosphonic_acid

Method to return whether or not this chemical is in the category phosphonic_acid, [-]

property is_polyol

Method to return whether or not this chemical is in the category polyol, [-]

property is_primary_aldimine

Method to return whether or not this chemical is in the category primary_aldimine, [-]

property is_primary_amine

Method to return whether or not this chemical is in the category primary_amine, [-]

property is_primary_ketimine

Method to return whether or not this chemical is in the category primary_ketimine, [-]

property is_pyridyl

Method to return whether or not this chemical is in the category pyridyl, [-]

property is_quat

Method to return whether or not this chemical is in the category quat, [-]

property is_secondary_aldimine

Method to return whether or not this chemical is in the category secondary_aldimine, [-]

property is_secondary_amine

Method to return whether or not this chemical is in the category secondary_amine, [-]

property is_secondary_ketimine

Method to return whether or not this chemical is in the category secondary_ketimine, [-]

property is_siloxane

Method to return whether or not this chemical is in the category siloxane, [-]

property is_silyl_ether

Method to return whether or not this chemical is in the category silyl_ether, [-]

property is_sulfide

Method to return whether or not this chemical is in the category sulfide, [-]

property is_sulfinic_acid

Method to return whether or not this chemical is in the category sulfinic_acid, [-]

property is_sulfonate_ester

Method to return whether or not this chemical is in the category sulfonate_ester, [-]

property is_sulfone

Method to return whether or not this chemical is in the category sulfone, [-]

property is_sulfonic_acid

Method to return whether or not this chemical is in the category sulfonic_acid, [-]

property is_sulfoxide

Method to return whether or not this chemical is in the category sulfoxide, [-]

property is_tertiary_amine

Method to return whether or not this chemical is in the category tertiary_amine, [-]

property is_thial

Method to return whether or not this chemical is in the category thial, [-]

property is_thiocyanate

Method to return whether or not this chemical is in the category thiocyanate, [-]

property is_thioketone

Method to return whether or not this chemical is in the category thioketone, [-]

property is_thiolester

Method to return whether or not this chemical is in the category thiolester, [-]

property is_thionoester

Method to return whether or not this chemical is in the category thionoester, [-]

property isentropic_exponent

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

Examples

>>> Chemical('hydrogen').isentropic_exponent
1.40
property isobaric_expansion

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

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

Examples

Radical change in value just above and below the critical temperature of water:

>>> Chemical('water', T=647.1, P=22048320.0).isobaric_expansion
0.77
>>> Chemical('water', T=647.2, P=22048320.0).isobaric_expansion
0.39
property isobaric_expansion_g

Isobaric (constant-pressure) expansion of the gas phase of the chemical 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

Utilizes the temperature-derivative method of thermo.VolumeGas to perform the actual calculation. The derivatives are all numerical.

Examples

>>> Chemical('Hexachlorobenzene', T=900).isobaric_expansion_g
0.00115
property isobaric_expansion_l

Isobaric (constant-pressure) expansion of the liquid phase of the chemical 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

Utilizes the temperature-derivative method of thermo.volume.VolumeLiquid to perform the actual calculation. The derivatives are all numerical.

Examples

>>> Chemical('dodecane', T=400).isobaric_expansion_l
0.00116
property k

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

Utilizes the object oriented interfaces thermo.thermal_conductivity.ThermalConductivityLiquid and thermo.thermal_conductivity.ThermalConductivityGas to perform the actual calculation of each property.

Examples

>>> Chemical('ethanol', T=300).kl
0.16
>>> Chemical('ethanol', T=400).kg
0.026
property kg

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

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.ThermalConductivityGas; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('water', T=320, P=100).kg
0.02
property kl

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

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.ThermalConductivityLiquid; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('water', T=320).kl
0.63
property ks

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

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.thermal_conductivity.ThermalConductivitySolid; each Chemical instance creates one to actually perform the calculations.

property legal_status

Dictionary of legal status indicators for the chemical.

Examples

>>> Chemical('benzene').legal_status
{'DSL': 'LISTED', 'TSCA': 'LISTED', 'EINECS': 'LISTED', 'NLP': 'UNLISTED', 'SPIN': 'LISTED'}
property mass_fractions

Dictionary of atom:mass-weighted fractional occurence of elements. Useful when performing mass balances. For atom-fraction occurences, see atom_fractions.

Examples

>>> Chemical('water').mass_fractions
{'H': 0.11189834407236524, 'O': 0.8881016559276347}
property mu

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

Utilizes the object oriented interfaces thermo.viscosity.ViscosityLiquid and thermo.viscosity.ViscosityGas to perform the actual calculation of each property.

Examples

>>> Chemical('ethanol', T=300).mu
0.0010
>>> Chemical('ethanol', T=400).mu
1.18e-05
property mug

Viscosity of the chemical in the gas phase at its current temperature and pressure, 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.ViscosityGas; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('water', T=320, P=100).mug
1.04e-05
property mul

Viscosity of the chemical in the liquid phase at its current temperature and pressure, 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.ViscosityLiquid; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('water', T=320).mul
0.000576
property nu

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

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

Examples

>>> Chemical('argon').nu
1.38e-05
property nug

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

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

Utilizes the temperature and pressure dependent object oriented interfaces thermo.volume.VolumeGas, thermo.viscosity.ViscosityGas to calculate the actual properties.

Examples

>>> Chemical('methane', T=115).nug
2.5e-06
property nul

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

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

Utilizes the temperature and pressure dependent object oriented interfaces thermo.volume.VolumeLiquid, thermo.viscosity.ViscosityLiquid to calculate the actual properties.

Examples

>>> Chemical('methane', T=110).nul
2.858e-07
property permittivity

Relative permittivity (dielectric constant) of the chemical at its current temperature, [dimensionless].

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.permittivity.PermittivityLiquid; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('toluene', T=250).permittivity
2.497
property rdkitmol

RDKit object of the chemical, without hydrogen. If RDKit is not available, holds None.

For examples of what can be done with RDKit, see their website.

property rdkitmol_Hs

RDKit object of the chemical, with hydrogen. If RDKit is not available, holds None.

For examples of what can be done with RDKit, see their website.

property rho

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

Utilizes the object oriented interfaces thermo.volume.VolumeSolid, thermo.volume.VolumeLiquid, and thermo.volume.VolumeGas to perform the actual calculation of each property. Note that those interfaces provide output in units of m^3/mol.

Examples

>>> Chemical('decane', T=550, P=2E6).rho
498.
property rhog

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

Examples

Estimate the density of the core of the sun, at 15 million K and 26.5 PetaPascals, assuming pure helium (actually 68% helium):

>>> Chemical('helium', T=15E6, P=26.5E15).rhog
8519.

Compared to a result on Wikipedia of 150000 kg/m^3, the fundamental equation of state performs poorly.

>>> He = Chemical('helium', T=15E6, P=26.5E15)
>>> He.VolumeGas.method_P = 'IDEAL'
>>> He.rhog
850477.

The ideal-gas law performs somewhat better, but vastly overshoots the density prediction.

property rhogm

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

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

Examples

>>> Chemical('tungsten hexafluoride').rhogm
42.
property rhol

Liquid-phase mass density of the chemical at its current temperature and pressure, in units of [kg/m^3]. 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.volume.VolumeLiquid; each Chemical instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> Chemical('o-xylene', T=297).rhol
876.9946785618097
property rholm

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

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

Examples

>>> Chemical('nitrogen', T=70).rholm
29937.
property rhom

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

Utilizes the object oriented interfaces thermo.volume.VolumeSolid, thermo.volume.VolumeLiquid, and thermo.volume.VolumeGas to perform the actual calculation of each property. Note that those interfaces provide output in units of m^3/mol.

Examples

>>> Chemical('1-hexanol').rhom
7986.
property rhos

Solid-phase mass density of the chemical at its current temperature, in units of [kg/m^3]. 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.volume.VolumeSolid; each Chemical instance creates one to actually perform the calculations. Note that that interface provides output in molar units.

Examples

>>> Chemical('iron').rhos
7870.
property rhosm

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

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

Examples

>>> Chemical('palladium').rhosm
112760.
property rings

Number of rings in a chemical, computed with RDKit from a chemical’s SMILES. If RDKit is not available, holds None.

Examples

>>> Chemical('Paclitaxel').rings
7
set_TP_sources()[source]
set_constant_sources()[source]
set_constants()[source]
set_eos(T, P, eos=<class 'thermo.eos.PR'>)[source]
set_ref(T_ref=298.15, P_ref=101325, phase_ref='calc', H_ref=0, S_ref=0)[source]
set_thermo()[source]
property sigma

Surface tension of the chemical at its current temperature, 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.SurfaceTension; each Chemical instance creates one to actually perform the calculations.

Examples

>>> Chemical('water', T=320).sigma
0.068
>>> Chemical('water', T=320).SurfaceTension.solve_property(0.05)
417.2
property solubility_parameter

Solubility parameter of the chemical at its current temperature and pressure, in units of [Pa^0.5].

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

Calculated based on enthalpy of vaporization and molar volume. Normally calculated at STP. For uses of this property, see thermo.solubility.solubility_parameter.

Examples

>>> Chemical('NH3', T=200).solubility_parameter
31712.