Python Tutorial#

Getting Started#

See also

If you haven’t already installed Cantera, we have instructions for many platforms in our Installation section.

Start by opening an interactive Python session, for example by running IPython. Import the Cantera Python module and NumPy by running:

import cantera as ct
import numpy as np

When using Cantera, the first thing you usually need is an object representing: some phase of matter. Here, we’ll create a gas mixture

gas1 = ct.Solution('gri30.yaml')

To view the state of the mixture, call the gas1 object as if it were a function:

gas1()
  gri30:

       temperature   300 K
          pressure   1.0133e+05 Pa
           density   0.081894 kg/m^3
  mean mol. weight   2.016 kg/kmol
   phase of matter   gas

                          1 kg             1 kmol     
                     ---------------   ---------------
          enthalpy             26469             53361  J
   internal energy       -1.2108e+06        -2.441e+06  J
           entropy             64910        1.3086e+05  J/K
    Gibbs function       -1.9447e+07       -3.9204e+07  J
 heat capacity c_p             14311             28851  J/K
 heat capacity c_v             10187             20536  J/K

                      mass frac. Y      mole frac. X     chem. pot. / RT
                     ---------------   ---------------   ---------------
                H2                 1                 1           -15.717
     [  +52 minor]                 0                 0  

What you have just done is to create an object, gas1 that implements GRI-Mech 3.0, the 53-species, 325-reaction natural gas combustion mechanism developed by Smith et al. [1999].

The gas1 object has properties you would expect for a gas mixture - it has a temperature, a pressure, species mole and mass fractions, etc. As we’ll soon see, it has many more properties.

The summary of the state of gas1 printed above shows that new objects created from the gri30.yaml input file start out with a temperature of 300 K, a pressure of 1 atm, and have a composition that consists of only one species, in this case hydrogen. There is nothing special about H2 - it just happens to be the first species listed in the input file defining GRI-Mech 3.0. In general, whichever species is listed first will initially have a mole fraction of 1.0, and all of the others will be zero.

Setting the State#

The state of the object can easily be changed. For example:

gas1.TP = 1200, 101325

sets the temperature to 1200 K and the pressure to 101325 Pa (Cantera always uses SI units). After this statement, calling gas1() results in:

gas1()
  gri30:

       temperature   1200 K
          pressure   1.0133e+05 Pa
           density   0.020473 kg/m^3
  mean mol. weight   2.016 kg/kmol
   phase of matter   gas

                          1 kg             1 kmol     
                     ---------------   ---------------
          enthalpy        1.3295e+07        2.6802e+07  J
   internal energy        8.3457e+06        1.6825e+07  J
           entropy             85222        1.7181e+05  J/K
    Gibbs function       -8.8972e+07       -1.7937e+08  J
 heat capacity c_p             15377             31000  J/K
 heat capacity c_v             11253             22686  J/K

                      mass frac. Y      mole frac. X     chem. pot. / RT
                     ---------------   ---------------   ---------------
                H2                 1                 1           -17.978
     [  +52 minor]                 0                 0  

Thermodynamics generally requires that two properties in addition to composition information be specified to fix the intensive state of a substance (or mixture). The state of the mixture can be set using several combinations of two properties. The following are all equivalent:

gas1.TP = 1200, 101325          # temperature, pressure
gas1.TD = 1200, 0.020473        # temperature, density
gas1.HP = 1.3295e7, 101325      # specific enthalpy, pressure
gas1.UV = 8.3457e6, 1/0.020473  # specific internal energy, specific volume
gas1.SP = 85222, 101325         # specific entropy, pressure
gas1.SV = 85222, 1/0.020473     # specific entropy, specific volume

In each case, the values of the extensive properties must be entered per unit mass.

Properties may be read independently or together:

gas1.T
1199.9365830924924
gas1.h
13293801.951980166
gas1.UV
(8344978.606604642, 48.84482000683827)

The composition can be set in terms of either mole fractions (X) or mass fractions (Y):

gas1.X = 'CH4:1, O2:2, N2:7.52'

Mass and mole fractions can also be set using dict objects, which is convenient in cases where the composition is stored in a variable or being computed:

phi = 0.8
gas1.X = {'CH4':1, 'O2':2/phi, 'N2':2*3.76/phi}

When the composition alone is changed, the temperature and density are held constant. This means that the pressure and other intensive properties will change. The composition can also be set in conjunction with the intensive properties of the mixture:

gas1.TPX = 1200, 101325, 'CH4:1, O2:2, N2:7.52'
gas1()
  gri30:

       temperature   1200 K
          pressure   1.0133e+05 Pa
           density   0.28063 kg/m^3
  mean mol. weight   27.633 kg/kmol
   phase of matter   gas

                          1 kg             1 kmol     
                     ---------------   ---------------
          enthalpy        8.6193e+05        2.3818e+07  J
   internal energy        5.0087e+05        1.3841e+07  J
           entropy            8914.2        2.4633e+05  J/K
    Gibbs function       -9.8351e+06       -2.7178e+08  J
 heat capacity c_p            1397.3             38611  J/K
 heat capacity c_v            1096.4             30296  J/K

                      mass frac. Y      mole frac. X     chem. pot. / RT
                     ---------------   ---------------   ---------------
                O2           0.22014           0.19011           -28.747
               CH4          0.055187          0.095057           -35.961
                N2           0.72467           0.71483           -25.679
     [  +50 minor]                 0                 0  

The composition above was specified using a string. The format is a comma-separated list of <species name>:<relative mole numbers> pairs. The mole numbers will be normalized to produce the mole fractions, and therefore they are “relative” mole numbers. Mass fractions can be set in this way too by changing X to Y in the above statements.

The composition can also be set using an array, which must have the same size as the number of species. For example, to set all 53 mole fractions to the same value, do this:

gas1.X = np.ones(53)  # NumPy array of 53 ones

Or, to set all the mass fractions to equal values:

gas1.Y = np.ones(53)

When setting the state, you can control what properties are held constant by passing the special value None to the property setter. For example, to change the specific volume to 2.1 m\(^3\)/kg while holding entropy constant:

gas1.SV = None, 2.1

Or to set the mass fractions while holding temperature and pressure constant:

gas1.TPX = None, None, 'CH4:1.0, O2:0.5'

Working with a Subset of Species#

Many properties of a Solution provide values for each species present in the phase. If you want to get values only for a subset of these species, you can use Python’s “slicing” syntax to select data for just the species of interest. To get the mole fractions of just the major species in gas1, in the order specified, you can write:

Xmajor = gas1['CH4','O2','CO2','H2O','N2'].X

If you want to use the same set of species repeatedly, you can keep a reference to the sliced phase object:

major = gas1['CH4','O2','CO2','H2O','N2']
cp_major = major.partial_molar_cp
wdot_major = major.net_production_rates

The slice object and the original object share the same internal state, so modifications to one will affect the other.

Working With Mechanism Files#

In previous example, we created an object that models an ideal gas mixture with the species and reactions of GRI-Mech 3.0, using the gri30.yaml input file included with Cantera. Several other reaction mechanism files are included with Cantera, including ones that model high- temperature air, a hydrogen/oxygen reaction mechanism, and a few surface reaction mechanisms. These files are usually located in the data subdirectory of the Cantera installation directory, for example C:\Program Files\Cantera\data on Windows or /usr/local/cantera/data/ on Unix/Linux/Mac OS X machines, depending on how you installed Cantera and the options you specified.

If for some reason Cantera has difficulty finding where these files are on your system, set environment variable CANTERA_DATA to the directory or directories (separated using ; on Windows or : on other operating systems) where they are located. Alternatively, you can call function add_directory() to add a directory to the Cantera search path:

ct.add_directory('~/cantera/my_data_files')

Cantera input files are plain text files, and can be created with any text editor. See the page Introduction to YAML Input Files for more information.

A Cantera input file may contain more than one phase specification, or may contain specifications of interfaces (surfaces). Here we import definitions of two bulk phases and the interface between them from file diamond.yaml:

gas2 = ct.Solution('diamond.yaml', 'gas')
diamond = ct.Solution('diamond.yaml', 'diamond')
diamond_surf = ct.Interface('diamond.yaml' , 'diamond_100',
                            [gas2, diamond])

Note that the bulk (3D) phases that participate in the surface reactions must also be passed as arguments to Interface.

Converting CK-format files#

See the page Converting Chemkin Format Files for information on how to convert from CK-format to Cantera’s YAML format.

Getting Help#

In addition to the Sphinx-generated Python documentation, documentation of the Python classes and their methods can be accessed from within the Python interpreter as well.

Suppose you have created a Cantera object and want to know what methods are available for it, and get help on using the methods:

g = ct.Solution('gri30.yaml')

To get help on the Python class that this object is an instance of:

help(g)
Hide code cell output
Help on Solution in module cantera.composite object:

class Solution(cantera.transport.Transport, cantera.kinetics.Kinetics, cantera.thermo.ThermoPhase)
 |  A class for chemically-reacting solutions. Instances can be created to
 |  represent any type of solution -- a mixture of gases, a liquid solution, or
 |  a solid solution, for example.
 |  
 |  Class `Solution` derives from classes `ThermoPhase`, `Kinetics`, and
 |  `Transport`.  It defines no methods of its own, and is provided so that a
 |  single object can be used to compute thermodynamic, kinetic, and transport
 |  properties of a solution.
 |  
 |  To skip initialization of the Transport object, pass the keyword argument
 |  ``transport_model=None`` to the `Solution` constructor.
 |  
 |  The most common way to instantiate `Solution` objects is by using a phase
 |  definition, species and reactions defined in an input file::
 |  
 |      gas = ct.Solution('gri30.yaml')
 |  
 |  If an input file defines multiple phases, the corresponding key in the
 |  ``phases`` map can be used to specify the desired phase via the ``name`` keyword
 |  argument of the constructor::
 |  
 |      gas = ct.Solution('diamond.yaml', name='gas')
 |      diamond = ct.Solution('diamond.yaml', name='diamond')
 |  
 |  The name of the `Solution` object defaults to the *phase* identifier
 |  specified in the input file. Upon initialization of a `Solution` object,
 |  a custom name can assigned via::
 |  
 |      gas.name = 'my_custom_name'
 |  
 |  `Solution` objects can also be constructed using `Species` and `Reaction`
 |  objects which can themselves either be imported from input files or defined
 |  directly in Python::
 |  
 |      spec = ct.Species.list_from_file("gri30.yaml")
 |      spec_gas = ct.Solution(thermo='ideal-gas', species=spec)
 |      rxns = ct.Reaction.list_from_file("gri30.yaml", spec_gas)
 |      gas = ct.Solution(thermo='ideal-tas', kinetics='gas',
 |                        species=spec, reactions=rxns, name='my_custom_name')
 |  
 |  where the ``thermo`` and ``kinetics`` keyword arguments are strings
 |  specifying the thermodynamic and kinetics model, respectively, and
 |  ``species`` and ``reactions`` keyword arguments are lists of `Species` and
 |  `Reaction` objects, respectively. Note that importing the reactions from a
 |  YAML input file requires a `Kinetics` object containing the species, as
 |  shown.
 |  
 |  Types of underlying models that form the composite `Solution` object are
 |  queried using the ``thermo_model``, ``kinetics_model`` and
 |  ``transport_model`` attributes; further, the ``composite`` attribute is a
 |  shorthand returning a tuple containing the types of the three constitutive
 |  models.
 |  
 |  For non-trivial uses cases of this functionality, see the examples
 |  `extract_submechanism.py <https://cantera.org/examples/python/kinetics/extract_submechanism.py.html>`_
 |  and `mechanism_reduction.py <https://cantera.org/examples/python/kinetics/mechanism_reduction.py.html>`_.
 |  
 |  In addition, `Solution` objects can be constructed by passing the text of
 |  the YAML phase definition in directly, using the ``yaml`` keyword
 |  argument::
 |  
 |      yaml_def = '''
 |      phases:
 |      - name: gas
 |        thermo: ideal-gas
 |        kinetics: gas
 |        elements: [O, H, Ar]
 |        species:
 |        - gri30.yaml/species: all
 |        reactions:
 |        - gri30.yaml/reactions: declared-species
 |        skip-undeclared-elements: true
 |        skip-undeclared-third-bodies: true
 |        state: {T: 300, P: 1 atm}
 |      '''
 |      gas = ct.Solution(yaml=yaml_def)
 |  
 |  Method resolution order:
 |      Solution
 |      cantera.transport.Transport
 |      cantera.kinetics.Kinetics
 |      cantera.thermo.ThermoPhase
 |      cantera.solutionbase._SolutionBase
 |      builtins.object
 |  
 |  Methods inherited from cantera.transport.Transport:
 |  
 |  __init__(self, /, *args, **kwargs)
 |      Initialize self.  See help(type(self)) for accurate signature.
 |  
 |  __reduce_cython__(...)
 |  
 |  __setstate_cython__(...)
 |  
 |  get_binary_diff_coeffs_polynomial(self, i, j)
 |      Get the polynomial fit to the logarithm of temperature for
 |      the binary diffusion coefficient of species ``i`` and ``j``.
 |  
 |  get_collision_integral_polynomials(self, i, j)
 |      Get the polynomial fit to the logarithm of temperature for
 |      the collision integral of species ``i`` and ``j``.
 |  
 |  get_thermal_conductivity_polynomial(self, i)
 |      Get the polynomial fit to the logarithm of temperature for
 |      the thermal conductivity of species ``i``.
 |  
 |  get_viscosity_polynomial(self, i)
 |      Get the polynomial fit to the logarithm of temperature for
 |      the viscosity of species ``i``.
 |  
 |  set_binary_diff_coeffs_polynomial(self, i, j, values)
 |      Set the polynomial fit to the logarithm of temperature for
 |      the binary diffusion coefficient of species ``i`` and ``j``.
 |  
 |  set_collision_integral_polynomial(self, i, j, avalues, bvalues, cvalues, actualT=True)
 |      Get the polynomial fit to the logarithm of temperature for
 |      the collision integral of species ``i`` and ``j``.
 |  
 |  set_thermal_conductivity_polynomial(self, i, values)
 |      Set the polynomial fit to the logarithm of temperature for
 |      the thermal conductivity of species ``i``.
 |  
 |  set_viscosity_polynomial(self, i, values)
 |      Set the polynomial fit to the logarithm of temperature for
 |      the viscosity of species ``i``.
 |  
 |  ----------------------------------------------------------------------
 |  Static methods inherited from cantera.transport.Transport:
 |  
 |  __new__(*args, **kwargs) from builtins.type
 |      Create and return a new object.  See help(type) for accurate signature.
 |  
 |  ----------------------------------------------------------------------
 |  Data descriptors inherited from cantera.transport.Transport:
 |  
 |  CK_mode
 |      Boolean to indicate if the chemkin interpretation is used.
 |  
 |  binary_diff_coeffs
 |      Binary diffusion coefficients [m^2/s].
 |  
 |  electrical_conductivity
 |      Electrical conductivity. [S/m].
 |  
 |  mix_diff_coeffs
 |      Mixture-averaged diffusion coefficients [m^2/s] relating the
 |      mass-averaged diffusive fluxes (with respect to the mass averaged
 |      velocity) to gradients in the species mole fractions.
 |  
 |  mix_diff_coeffs_mass
 |      Mixture-averaged diffusion coefficients [m^2/s] relating the
 |      diffusive mass fluxes to gradients in the species mass fractions.
 |  
 |  mix_diff_coeffs_mole
 |      Mixture-averaged diffusion coefficients [m^2/s] relating the
 |      molar diffusive fluxes to gradients in the species mole fractions.
 |  
 |  mobilities
 |      Electrical mobilities of charged species [m^2/s-V]
 |  
 |  multi_diff_coeffs
 |      Multicomponent diffusion coefficients, D[i,j], the diffusion
 |      coefficient for species i due to concentration gradients in
 |      species j [m**2/s].
 |  
 |  species_viscosities
 |      Pure species viscosities [Pa-s]
 |  
 |  thermal_conductivity
 |      Thermal conductivity. [W/m/K]
 |  
 |  thermal_diff_coeffs
 |      Return a one-dimensional array of the species thermal diffusion
 |      coefficients [kg/m/s].
 |  
 |  transport_model
 |      Get/Set the string specifying the transport model, such as `multicomponent`,
 |      `mixture-averaged`, or `unity-Lewis-number`.
 |      
 |      Setting a new transport model deletes the underlying C++ Transport
 |      object and replaces it with a new one implementing the specified model.
 |  
 |  viscosity
 |      Viscosity [Pa-s].
 |  
 |  ----------------------------------------------------------------------
 |  Methods inherited from cantera.kinetics.Kinetics:
 |  
 |  add_reaction(self, rxn)
 |      Add a new reaction to this phase.
 |  
 |  kinetics_species_index(self, species, phase=0)
 |      The index of species ``species`` of phase ``phase`` within arrays returned
 |      by methods of class `Kinetics`. If ``species`` is a string, the ``phase``
 |      argument is unused.
 |  
 |  kinetics_species_name(self, k)
 |      Name of the species with index ``k`` in the arrays returned by methods
 |      of class `Kinetics`.
 |  
 |  modify_reaction(self, irxn, rxn)
 |      Modify the `Reaction` with index ``irxn`` to have the same rate parameters as
 |      ``rxn``. ``rxn`` must have the same reactants and products and use the same rate
 |      parameterization (for example, `ArrheniusRate`, `FalloffRate`, `PlogRate`, etc.)
 |      as the existing reaction. This method does not modify the third-body
 |      efficiencies, reaction orders, or reversibility of the reaction.
 |  
 |  multiplier(self, i_reaction)
 |      A scaling factor applied to the rate coefficient for reaction
 |      ``i_reaction``. Can be used to carry out sensitivity analysis or to
 |      selectively disable a particular reaction. See `set_multiplier`.
 |  
 |  product_stoich_coeff(self, k_spec, i_reaction)
 |      The stoichiometric coefficient of species ``k_spec`` as a product in
 |      reaction ``i_reaction``.
 |  
 |  reactant_stoich_coeff(self, k_spec, i_reaction)
 |      The stoichiometric coefficient of species ``k_spec`` as a reactant in
 |      reaction ``i_reaction``.
 |  
 |  reaction(self, i_reaction)
 |      Return a `Reaction` object representing the reaction with index
 |      ``i_reaction``. Changes to this object do not affect the `Kinetics` or
 |      `Solution` object until the `modify_reaction` function is called.
 |  
 |  reaction_equations(self, indices=None)
 |      Returns a list containing the reaction equation for all reactions in the
 |      mechanism if ``indices`` is unspecified, or the equations for each
 |      reaction in the sequence ``indices``. For example::
 |      
 |          >>> gas.reaction_equations()
 |          ['2 O + M <=> O2 + M', 'O + H + M <=> OH + M', 'O + H2 <=> H + OH', ...]
 |          >>> gas.reaction_equations([2,3])
 |          ['O + H + M <=> OH + M', 'O + H2 <=> H + OH']
 |      
 |      See also `reaction_equation`.
 |  
 |  reactions(self)
 |      Return a list of all `Reaction` objects. Changes to these objects do not
 |      affect the `Kinetics` or `Solution` object until the `modify_reaction`
 |      function is called.
 |  
 |  set_multiplier(self, value, i_reaction=-1)
 |      Set the multiplier for for reaction ``i_reaction`` to ``value``.
 |      If ``i_reaction`` is not specified, then the multiplier for all reactions
 |      is set to ``value``. See `multiplier`.
 |  
 |  ----------------------------------------------------------------------
 |  Data descriptors inherited from cantera.kinetics.Kinetics:
 |  
 |  creation_rates
 |      Creation rates for each species. [kmol/m^3/s] for bulk phases or
 |      [kmol/m^2/s] for surface phases.
 |  
 |  creation_rates_ddC
 |      Calculate derivatives of species creation rates with respect to molar
 |      concentration at constant temperature, pressure and mole fractions.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  creation_rates_ddCi
 |      Calculate derivatives for species creation rates with respect to species
 |      concentration at constant temperature, pressure, and concentration of all other
 |      species. For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      The method returns a matrix with `n_total_species` rows and `n_total_species`
 |      columns.
 |      
 |      For a derivative with respect to :math: `c_i`, all other :math: `c_i` are
 |      held constant.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |      
 |      .. versionadded:: 3.0
 |  
 |  creation_rates_ddP
 |      Calculate derivatives of species creation rates with respect to pressure
 |      at constant temperature, molar concentration and mole fractions.
 |  
 |  creation_rates_ddT
 |      Calculate derivatives of species creation rates with respect to temperature
 |      at constant pressure, molar concentration and mole fractions.
 |  
 |  creation_rates_ddX
 |      Calculate derivatives for species creation rates with respect to species
 |      concentrations at constant temperature, pressure and molar concentration.
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      Note that for derivatives with respect to :math:`X_i`, all other :math:`X_j`
 |      are held constant, rather than enforcing :math:`\sum X_j = 1`.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  delta_enthalpy
 |      Change in enthalpy for each reaction [J/kmol].
 |  
 |  delta_entropy
 |      Change in entropy for each reaction [J/kmol/K].
 |  
 |  delta_gibbs
 |      Change in Gibbs free energy for each reaction [J/kmol].
 |  
 |  delta_standard_enthalpy
 |      Change in standard-state enthalpy (independent of composition) for
 |      each reaction [J/kmol].
 |  
 |  delta_standard_entropy
 |      Change in standard-state entropy (independent of composition) for
 |      each reaction [J/kmol/K].
 |  
 |  delta_standard_gibbs
 |      Change in standard-state Gibbs free energy (independent of composition)
 |      for each reaction [J/kmol].
 |  
 |  derivative_settings
 |      Property setting behavior of derivative evaluation.
 |      
 |      For :ct:`BulkKinetics`, the following keyword/value pairs are supported:
 |      
 |      -  ``skip-third-bodies`` (boolean) ... if `False` (default), third body
 |         concentrations are considered for the evaluation of derivatives
 |      
 |      -  ``skip-falloff`` (boolean) ... if `True` (default), third-body effects
 |         on reaction rates are not considered.
 |      
 |      -  ``rtol-delta`` (double) ... relative tolerance used to perturb properties
 |         when calculating numerical derivatives. The default value is 1e-8.
 |      
 |      Derivative settings are updated using a dictionary::
 |      
 |          >>> gas.derivative_settings = {"skip-falloff": True}
 |      
 |      Passing an empty dictionary will reset all values to their defaults.
 |  
 |  destruction_rates
 |      Destruction rates for each species. [kmol/m^3/s] for bulk phases or
 |      [kmol/m^2/s] for surface phases.
 |  
 |  destruction_rates_ddC
 |      Calculate derivatives of species destruction rates with respect to molar
 |      concentration at constant temperature, pressure and mole fractions.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  destruction_rates_ddCi
 |      Calculate derivatives for species destruction rates with respect to species
 |      concentration at constant temperature, pressure, and concentration of all other
 |      species. For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      The method returns a matrix with `n_total_species` rows and `n_total_species` columns.
 |      For a derivative with respect to :math: `c_i`, all other :math: `c_i` are
 |      held constant.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |      
 |      .. versionadded:: 3.0
 |  
 |  destruction_rates_ddP
 |      Calculate derivatives of species destruction rates with respect to pressure
 |      at constant temperature, molar concentration and mole fractions.
 |  
 |  destruction_rates_ddT
 |      Calculate derivatives of species destruction rates with respect to temperature
 |      at constant pressure, molar concentration and mole fractions.
 |  
 |  destruction_rates_ddX
 |      Calculate derivatives for species destruction rates with respect to species
 |      concentrations at constant temperature, pressure and molar concentration.
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      Note that for derivatives with respect to :math:`X_i`, all other :math:`X_j`
 |      are held constant, rather than enforcing :math:`\sum X_j = 1`.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  equilibrium_constants
 |      Equilibrium constants in concentration units for all reactions.
 |  
 |  forward_rate_constants
 |      Forward rate constants for all reactions.
 |      
 |      The computed values include all temperature-dependent and pressure-dependent
 |      contributions. By default, third-body concentrations are only considered if
 |      they are part of the reaction rate definition; for a legacy implementation that
 |      includes third-body concentrations, see `use_legacy_rate_constants`.
 |  
 |  forward_rate_constants_ddC
 |      Calculate derivatives for forward rate constants with respect to molar
 |      concentration at constant temperature, pressure and mole fractions.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  forward_rate_constants_ddP
 |      Calculate derivatives for forward rate constants with respect to pressure
 |      at constant temperature, molar concentration and mole fractions.
 |  
 |  forward_rate_constants_ddT
 |      Calculate derivatives for forward rate constants with respect to temperature
 |      at constant pressure, molar concentration and mole fractions.
 |  
 |  forward_rates_of_progress
 |      Forward rates of progress for the reactions. [kmol/m^3/s] for bulk
 |      phases or [kmol/m^2/s] for surface phases.
 |  
 |  forward_rates_of_progress_ddC
 |      Calculate derivatives for forward rates-of-progress with respect to molar
 |      concentration at constant temperature, pressure and mole fractions.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  forward_rates_of_progress_ddCi
 |      Calculate derivatives for forward rates-of-progress with respect to species
 |      concentrations at constant temperature, pressure and remaining species
 |      concentrations.
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      Note that for derivatives with respect to :math:`c_i`, all other :math:`c_j`
 |      are held constant.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |      
 |      .. versionadded:: 3.0
 |  
 |  forward_rates_of_progress_ddP
 |      Calculate derivatives for forward rates-of-progress with respect to pressure
 |      at constant temperature, molar concentration and mole fractions.
 |  
 |  forward_rates_of_progress_ddT
 |      Calculate derivatives for forward rates-of-progress with respect to temperature
 |      at constant pressure, molar concentration and mole fractions.
 |  
 |  forward_rates_of_progress_ddX
 |      Calculate derivatives for forward rates-of-progress with respect to species
 |      concentrations at constant temperature, pressure and molar concentration.
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      Note that for derivatives with respect to :math:`X_i`, all other :math:`X_j`
 |      are held constant, rather than enforcing :math:`\sum X_j = 1`.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  heat_production_rates
 |      Get the volumetric heat production rates [W/m^3] on a per-reaction
 |      basis. The sum over all reactions results in the total volumetric heat
 |      release rate.
 |      Example: C. K. Law: Combustion Physics (2006), Fig. 7.8.6
 |      
 |      >>> gas.heat_production_rates[1]  # heat production rate of the 2nd reaction
 |  
 |  heat_release_rate
 |      Get the total volumetric heat release rate [W/m^3].
 |  
 |  kinetics_model
 |      Return type of kinetics.
 |  
 |  kinetics_species_names
 |      A list of all species names, corresponding to the arrays returned by
 |      methods of class `Kinetics`.
 |  
 |  n_phases
 |      Number of phases in the reaction mechanism.
 |  
 |  n_reactions
 |      Number of reactions in the reaction mechanism.
 |  
 |  n_total_species
 |      Total number of species in all phases participating in the kinetics
 |      mechanism.
 |  
 |  net_production_rates
 |      Net production rates for each species. [kmol/m^3/s] for bulk phases or
 |      [kmol/m^2/s] for surface phases.
 |  
 |  net_production_rates_ddC
 |      Calculate derivatives of species net production rates with respect to molar
 |      density at constant temperature, pressure and mole fractions.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  net_production_rates_ddCi
 |      Calculate derivatives for species net production rates with respect to species
 |      concentration at constant temperature, pressure, and concentration of all other
 |      species. For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      The method returns a matrix with `n_total_species` rows and `n_total_species` columns.
 |      For a derivative with respect to :math: `c_i`, all other :math: `c_i` are
 |      held constant.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |      
 |      .. versionadded:: 3.0
 |  
 |  net_production_rates_ddP
 |      Calculate derivatives of species net production rates with respect to pressure
 |      at constant temperature, molar concentration and mole fractions.
 |  
 |  net_production_rates_ddT
 |      Calculate derivatives of species net production rates with respect to
 |      temperature at constant pressure, molar concentration and mole fractions.
 |  
 |  net_production_rates_ddX
 |      Calculate derivatives for species net production rates with respect to species
 |      concentrations at constant temperature, pressure and molar concentration.
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      Note that for derivatives with respect to :math:`X_i`, all other :math:`X_j`
 |      are held constant, rather than enforcing :math:`\sum X_j = 1`.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  net_rates_of_progress
 |      Net rates of progress for the reactions. [kmol/m^3/s] for bulk phases
 |      or [kmol/m^2/s] for surface phases.
 |  
 |  net_rates_of_progress_ddC
 |      Calculate derivatives for net rates-of-progress with respect to molar
 |      concentration at constant temperature, pressure and mole fractions.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  net_rates_of_progress_ddCi
 |      Calculate derivatives for net rates-of-progress with respect to species
 |      concentrations at constant temperature, pressure and remaining species
 |      concentrations. For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      Note that for derivatives with respect to :math:`c_i`, all other :math:`c_j`
 |      are held constant.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |      
 |      .. versionadded:: 3.0
 |  
 |  net_rates_of_progress_ddP
 |      Calculate derivatives for net rates-of-progress with respect to pressure
 |      at constant temperature, molar concentration and mole fractions.
 |  
 |  net_rates_of_progress_ddT
 |      Calculate derivatives for net rates-of-progress with respect to temperature
 |      at constant pressure, molar concentration and mole fractions.
 |  
 |  net_rates_of_progress_ddX
 |      Calculate derivatives for net rates-of-progress with respect to species
 |      concentrations at constant temperature, pressure and molar concentration.
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      Note that for derivatives with respect to :math:`X_i`, all other :math:`X_j`
 |      are held constant, rather than enforcing :math:`\sum X_j = 1`.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  product_stoich_coeffs
 |      The array of product stoichiometric coefficients. Element ``[k,i]`` of
 |      this array is the product stoichiometric coefficient of species ``k`` in
 |      reaction ``i``.
 |      
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      .. versionchanged:: 3.0
 |      
 |          Method was changed to a property in Cantera 3.0.
 |  
 |  product_stoich_coeffs_reversible
 |      The array of product stoichiometric coefficients of reversible reactions.
 |      Element ``[k,i]`` of this array is the product stoichiometric coefficient
 |      of species ``k`` in reaction ``i``.
 |      
 |      For sparse output, set ``ct.use_sparse(True)``.
 |  
 |  reactant_stoich_coeffs
 |      The array of reactant stoichiometric coefficients. Element ``[k,i]`` of
 |      this array is the reactant stoichiometric coefficient of species ``k`` in
 |      reaction ``i``.
 |      
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      .. versionchanged:: 3.0
 |      
 |          Method was changed to a property in Cantera 3.0.
 |  
 |  reaction_phase_index
 |      The index of the phase where the reactions occur.
 |      
 |      .. deprecated:: 3.0
 |      
 |          After Cantera 3.0, the reacting phase is always the first phase associated
 |          with the Kinetics object. This method will be removed after Cantera 3.1.
 |  
 |  reverse_rate_constants
 |      Reverse rate constants for all reactions.
 |      
 |      The computed values include all temperature-dependent and pressure-dependent
 |      contributions. By default, third-body concentrations are only considered if
 |      they are part of the reaction rate definition; for a legacy implementation that
 |      includes third-body concentrations, see `use_legacy_rate_constants`.
 |  
 |  reverse_rates_of_progress
 |      Reverse rates of progress for the reactions. [kmol/m^3/s] for bulk
 |      phases or [kmol/m^2/s] for surface phases.
 |  
 |  reverse_rates_of_progress_ddC
 |      Calculate derivatives for reverse rates-of-progress with respect to molar
 |      concentration at constant temperature, pressure and mole fractions.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  reverse_rates_of_progress_ddCi
 |      Calculate derivatives for reverse rates-of-progress with respect to species
 |      concentrations at constant temperature, pressure and remaining species
 |      concentrations.
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      Note that for derivatives with respect to :math:`c_i`, all other :math:`c_j`
 |      are held constant.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |      
 |      .. versionadded:: 3.0
 |  
 |  reverse_rates_of_progress_ddP
 |      Calculate derivatives for reverse rates-of-progress with respect to pressure
 |      at constant temperature, molar concentration and mole fractions.
 |  
 |  reverse_rates_of_progress_ddT
 |      Calculate derivatives for reverse rates-of-progress with respect to temperature
 |      at constant pressure, molar concentration and mole fractions.
 |  
 |  reverse_rates_of_progress_ddX
 |      Calculate derivatives for reverse rates-of-progress with respect to species
 |      concentrations at constant temperature, pressure and molar concentration.
 |      For sparse output, set ``ct.use_sparse(True)``.
 |      
 |      Note that for derivatives with respect to :math:`X_i`, all other :math:`X_j`
 |      are held constant, rather than enforcing :math:`\sum X_j = 1`.
 |      
 |      .. warning::
 |      
 |          This property is an experimental part of the Cantera API and
 |          may be changed or removed without notice.
 |  
 |  third_body_concentrations
 |      Effective third-body concentrations used by individual reactions; values
 |      are only defined for reactions involving third-bodies and are set to
 |      not-a-number otherwise.
 |  
 |  ----------------------------------------------------------------------
 |  Methods inherited from cantera.thermo.ThermoPhase:
 |  
 |  __call__(self, /, *args, **kwargs)
 |      Call self as a function.
 |  
 |  add_species(self, species)
 |      Add a new species to this phase. Missing elements will be added
 |      automatically.
 |  
 |  add_species_alias(self, name, alias)
 |      Add the alternate species name ``alias`` for an original species ``name``.
 |  
 |  atomic_weight(self, m)
 |      Atomic weight [kg/kmol] of element ``m``
 |  
 |  element_index(self, element)
 |      The index of element ``element``, which may be specified as a string or
 |      an integer. In the latter case, the index is checked for validity and
 |      returned. If no such element is present, an exception is thrown.
 |  
 |  element_name(self, m)
 |      Name of the element with index ``m``.
 |  
 |  elemental_mass_fraction(self, m)
 |      Get the elemental mass fraction :math:`Z_{\mathrm{mass},m}` of element
 |      :math:`m` as defined by:
 |      
 |      .. math:: Z_{\mathrm{mass},m} = \sum_k \frac{a_{m,k} M_m}{M_k} Y_k
 |      
 |      with :math:`a_{m,k}` being the number of atoms of element :math:`m` in
 |      species :math:`k`, :math:`M_m` the atomic weight of element :math:`m`,
 |      :math:`M_k` the molecular weight of species :math:`k`, and :math:`Y_k`
 |      the mass fraction of species :math:`k`::
 |      
 |          >>> phase.elemental_mass_fraction('H')
 |          1.0
 |      
 |      :param m:
 |          Base element, may be specified by name or by index.
 |  
 |  elemental_mole_fraction(self, m)
 |      Get the elemental mole fraction :math:`Z_{\mathrm{mole},m}` of element
 |      :math:`m` (the number of atoms of element m divided by the total number
 |      of atoms) as defined by:
 |      
 |      .. math:: Z_{\mathrm{mole},m} = \frac{\sum_k a_{m,k} X_k}
 |                                           {\sum_k \sum_j a_{j,k} X_k}
 |      
 |      with :math:`a_{m,k}` being the number of atoms of element :math:`m` in
 |      species :math:`k`, :math:`\sum_j` being a sum over all elements, and
 |      :math:`X_k` being the mole fraction of species :math:`k`::
 |      
 |          >>> phase.elemental_mole_fraction('H')
 |          1.0
 |      
 |      :param m:
 |          Base element, may be specified by name or by index.
 |  
 |  equilibrate(self, XY, solver='auto', rtol=1e-09, max_steps=1000, max_iter=100, estimate_equil=0, log_level=0)
 |      Set to a state of chemical equilibrium holding property pair
 |      ``XY`` constant.
 |      
 |      :param XY:
 |          A two-letter string, which must be one of the set::
 |      
 |              ['TP','TV','HP','SP','SV','UV']
 |      
 |      :param solver:
 |          Specifies the equilibrium solver to use. May be one of the following:
 |      
 |          * ``'element_potential'`` - a fast solver using the element potential
 |            method
 |          * ``'gibbs'`` - a slower but more robust Gibbs minimization solver
 |          * ``'vcs'`` - the VCS non-ideal equilibrium solver
 |          * ``'auto'`` - The element potential solver will be tried first, then
 |            if it fails the Gibbs solver will be tried.
 |      :param rtol:
 |          The relative error tolerance.
 |      :param max_steps:
 |          The maximum number of steps in composition to take to find a converged
 |          solution.
 |      :param max_iter:
 |          For the Gibbs minimization solver, this specifies the number of
 |          outer iterations on T or P when some property pair other
 |          than TP is specified.
 |      :param estimate_equil:
 |          Integer indicating whether the solver should estimate its own
 |          initial condition. If 0, the initial mole fraction vector in the
 |          `ThermoPhase` object is used as the initial condition. If 1, the
 |          initial mole fraction vector is used if the element abundances are
 |          satisfied. If -1, the initial mole fraction vector is thrown out,
 |          and an estimate is formulated.
 |      :param log_level:
 |          Set to a value greater than 0 to write diagnostic output.
 |  
 |  equivalence_ratio(self, fuel=None, oxidizer=None, basis='mole', include_species=None)
 |      Get the equivalence ratio :math:`\phi` of the current mixture, which is a
 |      conserved quantity. Considers the oxidation of C to CO2, H to H2O
 |      and S to SO2. Other elements are assumed not to participate in oxidation
 |      (that is, N ends up as N2). If fuel and oxidizer are not specified, the
 |      equivalence ratio is computed from the available oxygen and the
 |      required oxygen for complete oxidation:
 |      
 |      .. math:: \phi = \frac{Z_{\mathrm{mole},C} + Z_{\mathrm{mole},S}
 |                + \frac{1}{4}Z_{\mathrm{mole},H}} {\frac{1}{2}Z_{\mathrm{mole},O}}
 |      
 |      where :math:`Z_{\mathrm{mole},e}` is the elemental mole fraction of element
 |      :math:`e`. If the fuel and oxidizer compositions are specified, :math:`\phi` is
 |      computed from:
 |      
 |      .. math:: \phi = \frac{Z}{1-Z}\frac{1-Z_{\mathrm{st}}}{Z_{\mathrm{st}}}
 |      
 |      where :math:`Z` is the Bilger mixture fraction and :math:`Z_{\mathrm{st}}`
 |      the Bilger mixture fraction at stoichiometric conditions.
 |      The ``basis`` determines the composition of fuel and oxidizer:
 |      ``basis='mole'`` (default) means mole fractions, ``basis='mass'`` means
 |      mass fractions. Note that this definition takes all species into account.
 |      In case certain species like inert diluents should be ignored, a
 |      list of species can be provided with ``include_species``. This means that
 |      only these species are considered for the computation of the equivalence
 |      ratio. For more information, see `Python example
 |      <https://cantera.org/examples/python/thermo/equivalenceRatio.py.html>`_ ::
 |      
 |          >>> gas.set_equivalence_ratio(0.5, fuel='CH3:0.5, CH3OH:.5, N2:0.125', oxidizer='O2:0.21, N2:0.79, NO:0.01')
 |          >>> gas.equivalence_ratio(fuel='CH3:0.5, CH3OH:.5, N2:0.125', oxidizer='O2:0.21, N2:0.79, NO:0.01')
 |          0.5
 |      
 |      :param fuel:
 |          Fuel species name or mole/mass fractions as string, array, or dict.
 |      :param oxidizer:
 |          Oxidizer species name or mole/mass fractions as a string, array, or dict.
 |      :param basis:
 |          Determines if ``fuel`` and ``oxidizer`` are given in mole fractions
 |          (``basis="mole"``) or mass fractions (``basis="mass"``)
 |      :param include_species:
 |          List of species names (optional). Only these species are considered for the
 |          computation of the equivalence ratio. By default, all species are considered
 |  
 |  find_isomers(self, comp)
 |      Find species/isomers matching a composition specified by ``comp``.
 |  
 |  mass_fraction_dict(self, threshold=0.0)
 |      Return a dictionary giving the mass fraction for each species by name where the
 |      mass fraction is greater than ``threshold``.
 |  
 |  mixture_fraction(self, fuel, oxidizer, basis='mole', element='Bilger')
 |      Get the mixture fraction of the current mixture in
 |      (kg fuel / (kg oxidizer + kg fuel)). This is a quantity that is conserved after
 |      oxidation. Considers the oxidation of C to CO2, H to H2O and S to SO2. Other
 |      elements are assumed not to participate in oxidation (that is, N ends up as N2).
 |      The ``basis`` determines the composition of fuel and oxidizer:
 |      ``basis="mole"`` (default) means mole fractions, ``basis="mass"`` means mass
 |      fractions. The mixture fraction can be computed from a single element (for
 |      example, carbon with ``element="C"``)
 |      
 |      .. math:: Z_m = \frac{Z_{\mathrm{mass},m}-Z_{\mathrm{mass},m,\mathrm{ox}}}
 |          {Z_{\mathrm{mass},\mathrm{fuel}}-Z_{\mathrm{mass},m,\mathrm{ox}}}
 |      
 |      where :math:`Z_{\mathrm{mass},m}` is the elemental mass fraction of
 |      element :math:`m` in the mixture, and :math:`Z_{\mathrm{mass},m,\mathrm{ox}}`
 |      and :math:`Z_{\mathrm{mass},\mathrm{fuel}}` are the elemental mass fractions of
 |      the oxidizer and fuel, or from the Bilger mixture fraction
 |      (``element="Bilger"``), which considers the elements C, S, H and O
 |      :cite:p:`bilger1979`. The Bilger mixture fraction is computed by default:
 |      
 |      .. math:: Z_m = Z_{\mathrm{Bilger}} = \frac{\beta-\beta_{\mathrm{ox}}}
 |          {\beta_{\mathrm{fuel}}-\beta_{\mathrm{ox}}}
 |      
 |      with
 |      
 |      .. math:: \beta = 2\frac{Z_C}{M_C}+2\frac{Z_S}{M_S}+\frac{1}{2}\frac{Z_H}{M_H}
 |          - \frac{Z_O}{M_O}
 |      
 |      and :math:`M_m` the atomic weight of element :math:`m`.
 |      For more information, see `Python example
 |      <https://cantera.org/examples/python/thermo/equivalenceRatio.py.html>`_.::
 |      
 |          >>> gas.set_mixture_fraction(0.5, 'CH3:0.5, CH3OH:0.5, N2:0.125', 'O2:0.21, N2:0.79, NO:0.01')
 |          >>> gas.mixture_fraction('CH3:0.5, CH3OH:0.5, N2:0.125', 'O2:0.21, N2:0.79, NO:.01')
 |          0.5
 |      
 |      :param fuel:
 |          Fuel species name or mole/mass fractions as string, array, or dict.
 |      :param oxidizer:
 |          Oxidizer species name or mole/mass fractions as a string, array, or
 |          dict.
 |      :param basis:
 |          Determines if ``fuel`` and ``oxidizer`` are given in mole
 |          fractions (``basis='mole'``) or mass fractions (``basis='mass'``)
 |      :param element:
 |          Computes the mixture fraction from the specified elemental
 |          mass fraction (given by element name or element index) or as
 |          the Bilger mixture fraction (default)
 |  
 |  modify_species(self, k, species)
 |      Modify the thermodynamic data associated with a species. The species name,
 |      elemental composition, and type of thermo parameterization must be unchanged.
 |  
 |  mole_fraction_dict(self, threshold=0.0)
 |      Return a dictionary giving the mole fraction for each species by name where the
 |      mole fraction is greater than ``threshold``.
 |  
 |  n_atoms(self, species, element)
 |      Number of atoms of element ``element`` in species ``species``. The element
 |      and species may be specified by name or by index.
 |      
 |      >>> phase.n_atoms('CH4','H')
 |      4
 |  
 |  report(self, show_thermo=True, threshold=1e-14)
 |      Generate a report describing the thermodynamic state of this phase. To
 |      print the report to the terminal, simply call the phase object. The
 |      following two statements are equivalent::
 |      
 |      >>> phase()
 |      >>> print(phase.report())
 |      
 |      :param show_thermo:
 |          A Boolean argument specifying whether to show phase thermodynamic
 |          information in the ouptut.
 |      :param threshold:
 |          The threshold used to clip data in the output. Values below the threshold
 |          are not displayed.
 |  
 |  set_discretized_electron_energy_distribution(self, levels, distribution)
 |      Set electron energy distribution. When this method is used, electron
 |      temperature is calculated from the distribution.
 |      
 |      :param levels:
 |          vector of electron energy levels [eV]
 |      :param distribution:
 |          vector of distribution
 |  
 |  set_equivalence_ratio(self, phi, fuel, oxidizer, basis='mole', *, diluent=None, fraction=None)
 |      Set the composition to a mixture of ``fuel`` and ``oxidizer`` at the
 |      specified equivalence ratio ``phi``, holding temperature and pressure
 |      constant. Considers the oxidation of C to CO2, H to H2O and S to SO2.
 |      Other elements are assumed not to participate in oxidation (that is,
 |      N ends up as N2). The ``basis`` determines the fuel and oxidizer
 |      compositions: ``basis='mole'`` means mole fractions (default),
 |      ``basis='mass'`` means mass fractions. The fuel/oxidizer mixture can be
 |      be diluted by a ``diluent`` based on a mixing ``fraction``. The amount of
 |      diluent is quantified as a fraction of fuel, oxidizer or the fuel/oxidizer
 |      mixture. For more information, see `Python example
 |      <https://cantera.org/examples/python/thermo/equivalenceRatio.py.html>`_ ::
 |      
 |          >>> gas.set_equivalence_ratio(0.5, 'CH4', 'O2:1.0, N2:3.76', basis='mole')
 |          >>> gas.mass_fraction_dict()
 |          {'CH4': 0.02837633052851, 'N2': 0.7452356312613, 'O2': 0.22638803821018}
 |          >>> gas.set_equivalence_ratio(1.2, 'NH3:0.8,CO:0.2', 'O2:1', basis='mole')
 |          >>> gas.mass_fraction_dict()
 |          {'CO': 0.14784006249290, 'NH3': 0.35956645545401, 'O2': 0.49259348205308}
 |      
 |      :param phi:
 |          Equivalence ratio
 |      :param fuel:
 |          Fuel species name or mole/mass fractions as string, array, or dict.
 |      :param oxidizer:
 |          Oxidizer species name or mole/mass fractions as a string, array, or dict.
 |      :param basis:
 |          Determines if ``fuel`` and ``oxidizer`` are given in mole
 |          fractions (``basis='mole'``) or mass fractions (``basis='mass'``).
 |      :param diluent:
 |          Optional parameter. Required if dilution is used. Specifies the composition
 |          of the diluent in mole/mass fractions as a string, array or dict.
 |      :param fraction:
 |          Optional parameter. Dilutes the fuel/oxidizer mixture with the diluent
 |          according to ``fraction``. Fraction can refer to the fraction of diluent in
 |          the  mixture (for example ``fraction="diluent:0.7"`` will create a mixture
 |          with 30 % fuel/oxidizer and 70 % diluent), the fraction of fuel in the
 |          mixture (for example ``fraction="fuel:0.1"`` means that the mixture contains
 |          10 % fuel. The amount of oxidizer is determined from the equivalence ratio
 |          and the remaining mixture is the diluent) or fraction of oxidizer in the
 |          mixture (for example ``fraction="oxidizer:0.1"``). The fraction itself is
 |          interpreted as mole or mass fraction based on ``basis``. The diluent is not
 |          considered in the computation of the equivalence ratio. Default is no
 |          dilution or ``fraction=None``. May be given as string or dictionary (for
 |          example ``fraction={"fuel":0.7}``).
 |  
 |  set_mixture_fraction(self, mixture_fraction, fuel, oxidizer, basis='mole')
 |      Set the composition to a mixture of ``fuel`` and ``oxidizer`` at the
 |      specified mixture fraction ``mixture_fraction`` (kg fuel / kg mixture), holding
 |      temperature and pressure constant. Considers the oxidation of C to CO2,
 |      H to H2O and S to SO2. Other elements are assumed not to participate in
 |      oxidation (that is, N ends up as N2). The ``basis`` determines the composition
 |      of fuel and oxidizer: ``basis='mole'`` (default) means mole fractions,
 |      ``basis='mass'`` means mass fractions. For more information, see `Python
 |      example
 |      <https://cantera.org/examples/python/thermo/equivalenceRatio.py.html>`_ ::
 |      
 |          >>> gas.set_mixture_fraction(0.5, 'CH4', 'O2:1.0, N2:3.76')
 |          >>> gas.mass_fraction_dict()
 |          {'CH4': 0.5, 'N2': 0.38350014242997776, 'O2': 0.11649985757002226}
 |          >>> gas.set_mixture_fraction(0.5, {'NH3':0.8, 'CO':0.2}, 'O2:1.0')
 |          >>> gas.mass_fraction_dict()
 |          {'CO': 0.145682068778996, 'NH3': 0.354317931221004, 'O2': 0.5}
 |      
 |      :param mixture_fraction:
 |          Mixture fraction (kg fuel / kg mixture)
 |      :param fuel:
 |          Fuel species name or mole/mass fractions as string, array, or dict.
 |      :param oxidizer:
 |          Oxidizer species name or mole/mass fractions as a string, array, or
 |          dict.
 |      :param basis: determines if ``fuel`` and ``oxidizer`` are given in mole
 |          fractions (``basis='mole'``) or mass fractions (``basis='mass'``)
 |  
 |  set_unnormalized_mass_fractions(self, Y)
 |      Set the mass fractions without normalizing to force ``sum(Y) == 1.0``.
 |      Useful primarily when calculating derivatives with respect to ``Y[k]`` by
 |      finite difference.
 |  
 |  set_unnormalized_mole_fractions(self, X)
 |      Set the mole fractions without normalizing to force ``sum(X) == 1.0``.
 |      Useful primarily when calculating derivatives with respect to ``X[k]``
 |      by finite difference.
 |  
 |  species(self, k=None)
 |      Return the `Species` object for species ``k``, where ``k`` is either the
 |      species index or the species name. If ``k`` is not specified, a list of
 |      all species objects is returned. Changes to this object do not affect
 |      the `ThermoPhase` or `Solution` object until the `modify_species`
 |      function is called.
 |  
 |  species_index(self, species)
 |      The index of species ``species``, which may be specified as a string or
 |      an integer. In the latter case, the index is checked for validity and
 |      returned. If no such species is present, an exception is thrown.
 |  
 |  species_name(self, k)
 |      Name of the species with index ``k``.
 |  
 |  stoich_air_fuel_ratio(self, fuel, oxidizer, basis='mole')
 |      Get the stoichiometric air to fuel ratio (kg oxidizer / kg fuel). Considers the
 |      oxidation of C to CO2, H to H2O and S to SO2. Other elements are assumed
 |      not to participate in oxidation (that is, N ends up as N2).
 |      The ``basis`` determines the composition of fuel and oxidizer: ``basis='mole'`` (default)
 |      means mole fractions, ``basis='mass'`` means mass fractions::
 |      
 |          >>> gas.set_mixture_fraction(0.5, 'CH3:0.5, CH3OH:.5, N2:0.125', 'O2:0.21, N2:0.79, NO:0.01')
 |          >>> gas.stoich_air_fuel_ratio('CH3:0.5, CH3OH:.5, N2:0.125', 'O2:0.21, N2:0.79, NO:0.01')
 |          8.148040722239438
 |      
 |      :param fuel:
 |          Fuel species name or mole/mass fractions as string, array, or dict.
 |      :param oxidizer:
 |          Oxidizer species name or mole/mass fractions as a string, array, or
 |          dict.
 |      :param basis:
 |          Determines if ``fuel`` and ``oxidizer`` are given in mole
 |          fractions (``basis='mole'``) or mass fractions (``basis='mass'``)
 |  
 |  ----------------------------------------------------------------------
 |  Data descriptors inherited from cantera.thermo.ThermoPhase:
 |  
 |  DP
 |      Get/Set density [kg/m^3] and pressure [Pa].
 |  
 |  DPX
 |      Get/Set density [kg/m^3], pressure [Pa], and mole fractions.
 |  
 |  DPY
 |      Get/Set density [kg/m^3], pressure [Pa], and mass fractions.
 |  
 |  HP
 |      Get/Set enthalpy [J/kg or J/kmol] and pressure [Pa].
 |  
 |  HPX
 |      Get/Set enthalpy [J/kg or J/kmol], pressure [Pa] and mole fractions.
 |  
 |  HPY
 |      Get/Set enthalpy [J/kg or J/kmol], pressure [Pa] and mass fractions.
 |  
 |  P
 |      Pressure [Pa].
 |  
 |  P_sat
 |      Saturation pressure [Pa] at the current temperature.
 |  
 |  Pe
 |      Get electron Pressure [Pa].
 |  
 |  SP
 |      Get/Set entropy [J/kg/K or J/kmol/K] and pressure [Pa].
 |  
 |  SPX
 |      Get/Set entropy [J/kg/K or J/kmol/K], pressure [Pa], and mole fractions.
 |  
 |  SPY
 |      Get/Set entropy [J/kg/K or J/kmol/K], pressure [Pa], and mass fractions.
 |  
 |  SV
 |      Get/Set entropy [J/kg/K or J/kmol/K] and specific volume [m^3/kg or
 |      m^3/kmol].
 |  
 |  SVX
 |      Get/Set entropy [J/kg/K or J/kmol/K], specific volume [m^3/kg or
 |      m^3/kmol], and mole fractions.
 |  
 |  SVY
 |      Get/Set entropy [J/kg/K or J/kmol/K], specific volume [m^3/kg or
 |      m^3/kmol], and mass fractions.
 |  
 |  T
 |      Temperature [K].
 |  
 |  TD
 |      Get/Set temperature [K] and density [kg/m^3 or kmol/m^3].
 |  
 |  TDX
 |      Get/Set temperature [K], density [kg/m^3 or kmol/m^3], and mole
 |      fractions.
 |  
 |  TDY
 |      Get/Set temperature [K] and density [kg/m^3 or kmol/m^3], and mass
 |      fractions.
 |  
 |  TP
 |      Get/Set temperature [K] and pressure [Pa].
 |  
 |  TPX
 |      Get/Set temperature [K], pressure [Pa], and mole fractions.
 |  
 |  TPY
 |      Get/Set temperature [K], pressure [Pa], and mass fractions.
 |  
 |  T_sat
 |      Saturation temperature [K] at the current pressure.
 |  
 |  Te
 |      Get/Set electron Temperature [K].
 |  
 |  UV
 |      Get/Set internal energy [J/kg or J/kmol] and specific volume
 |      [m^3/kg or m^3/kmol].
 |  
 |  UVX
 |      Get/Set internal energy [J/kg or J/kmol], specific volume
 |      [m^3/kg or m^3/kmol], and mole fractions.
 |  
 |  UVY
 |      Get/Set internal energy [J/kg or J/kmol], specific volume
 |      [m^3/kg or m^3/kmol], and mass fractions.
 |  
 |  X
 |      Get/Set the species mole fractions. Can be set as an array, as a dictionary,
 |      or as a string. Always returns an array::
 |      
 |          >>> phase.X = [0.1, 0, 0, 0.4, 0, 0, 0, 0, 0.5]
 |          >>> phase.X = {'H2':0.1, 'O2':0.4, 'AR':0.5}
 |          >>> phase.X = 'H2:0.1, O2:0.4, AR:0.5'
 |          >>> phase.X
 |          array([0.1, 0, 0, 0.4, 0, 0, 0, 0, 0.5])
 |  
 |  Y
 |      Get/Set the species mass fractions. Can be set as an array, as a dictionary,
 |      or as a string. Always returns an array::
 |      
 |          >>> phase.Y = [0.1, 0, 0, 0.4, 0, 0, 0, 0, 0.5]
 |          >>> phase.Y = {'H2':0.1, 'O2':0.4, 'AR':0.5}
 |          >>> phase.Y = 'H2:0.1, O2:0.4, AR:0.5'
 |          >>> phase.Y
 |          array([0.1, 0, 0, 0.4, 0, 0, 0, 0, 0.5])
 |  
 |  activities
 |      Array of nondimensional activities. Returns either molar or molal
 |      activities depending on the convention of the thermodynamic model.
 |  
 |  activity_coefficients
 |      Array of nondimensional, molar activity coefficients.
 |  
 |  atomic_weights
 |      Array of atomic weight [kg/kmol] for each element in the mixture.
 |  
 |  basis
 |      Determines whether intensive thermodynamic properties are treated on a
 |      ``mass`` (per kg) or ``molar`` (per kmol) basis. This affects the values
 |      returned by the properties `h`, `u`, `s`, `g`, `v`, `density`, `cv`,
 |      and `cp`, as well as the values used with the state-setting properties
 |      such as `HPX` and `UV`.
 |  
 |  case_sensitive_species_names
 |      Enforce case-sensitivity for look up of species names
 |  
 |  charges
 |      Array of species charges [elem. charge].
 |  
 |  chemical_potentials
 |      Array of species chemical potentials [J/kmol].
 |  
 |  concentrations
 |      Get/Set the species concentrations. Units are kmol/m^3 for bulk phases, kmol/m^2
 |      for surface phases, and kmol/m for edge phases.
 |  
 |  cp
 |      Heat capacity at constant pressure [J/kg/K or J/kmol/K] depending
 |      on `basis`.
 |  
 |  cp_mass
 |      Specific heat capacity at constant pressure [J/kg/K].
 |  
 |  cp_mole
 |      Molar heat capacity at constant pressure [J/kmol/K].
 |  
 |  critical_density
 |      Critical density [kg/m^3 or kmol/m^3] depending on `basis`.
 |  
 |  critical_pressure
 |      Critical pressure [Pa].
 |  
 |  critical_temperature
 |      Critical temperature [K].
 |  
 |  cv
 |      Heat capacity at constant volume [J/kg/K or J/kmol/K] depending on
 |      `basis`.
 |  
 |  cv_mass
 |      Specific heat capacity at constant volume [J/kg/K].
 |  
 |  cv_mole
 |      Molar heat capacity at constant volume [J/kmol/K].
 |  
 |  density
 |      Density [kg/m^3 or kmol/m^3] depending on `basis`.
 |  
 |  density_mass
 |      (Mass) density [kg/m^3].
 |  
 |  density_mole
 |      Molar density [kmol/m^3].
 |  
 |  electric_potential
 |      Get/Set the electric potential [V] for this phase.
 |  
 |  electrochemical_potentials
 |      Array of species electrochemical potentials [J/kmol].
 |  
 |  electron_energy_distribution
 |      Electron energy distribution
 |  
 |  electron_energy_distribution_type
 |      Electron energy distribution type
 |  
 |  electron_energy_levels
 |      Electron energy levels [eV]
 |  
 |  element_names
 |      A list of all the element names.
 |  
 |  enthalpy_mass
 |      Specific enthalpy [J/kg].
 |  
 |  enthalpy_mole
 |      Molar enthalpy [J/kmol].
 |  
 |  entropy_mass
 |      Specific entropy [J/kg/K].
 |  
 |  entropy_mole
 |      Molar entropy [J/kmol/K].
 |  
 |  g
 |      Gibbs free energy [J/kg or J/kmol] depending on `basis`.
 |  
 |  gibbs_mass
 |      Specific Gibbs free energy [J/kg].
 |  
 |  gibbs_mole
 |      Molar Gibbs free energy [J/kmol].
 |  
 |  h
 |      Enthalpy [J/kg or J/kmol] depending on `basis`.
 |  
 |  has_phase_transition
 |      Returns true if the phase represents a substance with phase transitions
 |  
 |  int_energy_mass
 |      Specific internal energy [J/kg].
 |  
 |  int_energy_mole
 |      Molar internal energy [J/kmol].
 |  
 |  is_compressible
 |      Returns true if the density of the phase is an independent variable defining
 |      the thermodynamic state of a substance
 |  
 |  is_pure
 |      Returns true if the phase represents a pure (fixed composition) substance
 |  
 |  isothermal_compressibility
 |      Isothermal compressibility [1/Pa].
 |  
 |  isotropic_shape_factor
 |      Shape factor of isotropic-velocity distribution for electron energy
 |  
 |  max_temp
 |      Maximum temperature for which the thermodynamic data for the phase are
 |      valid.
 |  
 |  mean_electron_energy
 |      Mean electron energy [eV]
 |  
 |  mean_molecular_weight
 |      The mean molecular weight (molar mass) [kg/kmol].
 |  
 |  min_temp
 |      Minimum temperature for which the thermodynamic data for the phase are
 |      valid.
 |  
 |  molecular_weights
 |      Array of species molecular weights (molar masses) [kg/kmol].
 |  
 |  n_electron_energy_levels
 |      Number of electron energy levels
 |  
 |  n_elements
 |      Number of elements.
 |  
 |  n_selected_species
 |      Number of species selected for output (by slicing of Solution object)
 |  
 |  n_species
 |      Number of species.
 |  
 |  normalize_electron_energy_distribution_enabled
 |      Automatically normalize electron energy distribution
 |  
 |  partial_molar_cp
 |      Array of species partial molar specific heat capacities at constant
 |      pressure [J/kmol/K].
 |  
 |  partial_molar_enthalpies
 |      Array of species partial molar enthalpies [J/kmol].
 |  
 |  partial_molar_entropies
 |      Array of species partial molar entropies [J/kmol/K].
 |  
 |  partial_molar_int_energies
 |      Array of species partial molar internal energies [J/kmol].
 |  
 |  partial_molar_volumes
 |      Array of species partial molar volumes [m^3/kmol].
 |  
 |  phase_of_matter
 |      Get the thermodynamic phase (gas, liquid, etc.) at the current conditions.
 |  
 |  quadrature_method
 |      Quadrature method
 |  
 |  reference_pressure
 |      Reference state pressure [Pa].
 |  
 |  s
 |      Entropy [J/kg/K or J/kmol/K] depending on `basis`.
 |  
 |  sound_speed
 |      Speed of sound [m/s].
 |  
 |  species_names
 |      A list of all the species names.
 |  
 |  standard_concentration_units
 |      Get standard concentration units for this phase.
 |  
 |  standard_cp_R
 |      Array of nondimensional species standard-state specific heat capacities
 |      at constant pressure at the current temperature and pressure.
 |  
 |  standard_enthalpies_RT
 |      Array of nondimensional species standard-state enthalpies at the
 |      current temperature and pressure.
 |  
 |  standard_entropies_R
 |      Array of nondimensional species standard-state entropies at the
 |      current temperature and pressure.
 |  
 |  standard_gibbs_RT
 |      Array of nondimensional species standard-state Gibbs free energies at
 |      the current temperature and pressure.
 |  
 |  standard_int_energies_RT
 |      Array of nondimensional species standard-state internal energies at the
 |      current temperature and pressure.
 |  
 |  state
 |      Get/Set the full thermodynamic state as a single array, arranged as
 |      [temperature, density, mass fractions] for most phases. Useful mainly
 |      in cases where it is desired to store many states in a multidimensional
 |      array.
 |  
 |  state_size
 |      Return size of vector defining internal state of the phase.
 |  
 |  thermal_expansion_coeff
 |      Thermal expansion coefficient [1/K].
 |  
 |  thermo_model
 |      Return thermodynamic model describing phase.
 |  
 |  u
 |      Internal energy in [J/kg or J/kmol].
 |  
 |  v
 |      Specific volume [m^3/kg or m^3/kmol] depending on `basis`.
 |  
 |  volume_mass
 |      Specific volume [m^3/kg].
 |  
 |  volume_mole
 |      Molar volume [m^3/kmol].
 |  
 |  ----------------------------------------------------------------------
 |  Data and other attributes inherited from cantera.thermo.ThermoPhase:
 |  
 |  __pyx_vtable__ = <capsule object NULL>
 |  
 |  ----------------------------------------------------------------------
 |  Methods inherited from cantera.solutionbase._SolutionBase:
 |  
 |  __copy__(self)
 |  
 |  __del__(...)
 |  
 |  __getitem__(self, key, /)
 |      Return self[key].
 |  
 |  __getstate__(self)
 |      Save complete information of the current phase for pickling.
 |  
 |  __setstate__(self, pkl)
 |      Restore Solution from pickled information.
 |  
 |  clear_user_data(self)
 |      Clear all saved input data, so that the data given by `input_data` or
 |      `write_yaml` will only include values generated by Cantera based on the
 |      current object state.
 |  
 |  clear_user_header(self)
 |      Clear all saved header data, so that the data given by `input_header` or
 |      `write_yaml` will only include values generated by Cantera based on the
 |      current object state.
 |  
 |  update_user_data(self, data)
 |      Add the contents of the provided `dict` as additional fields when generating
 |      YAML phase definition files with `write_yaml` or in the data returned by
 |      `input_data`. Existing keys with matching names are overwritten.
 |  
 |  update_user_header(self, data)
 |      Add the contents of the provided `dict` as additional top-level YAML fields
 |      when generating files with `write_yaml` or in the data returned by
 |      `input_header`. Existing keys with matching names are overwritten.
 |  
 |  write_chemkin(self, mechanism_path=None, thermo_path=None, transport_path=None, sort_species=None, sort_elements=None, overwrite=False, quiet=False)
 |      Write this `~cantera.Solution` instance to one or more Chemkin-format files.
 |      See the documentation for `cantera.yaml2ck.convert` for information about the
 |      arguments to this function.
 |  
 |  write_yaml(self, filename=None, phases=None, units=None, precision=None, skip_user_defined=None, header=True)
 |      Write the definition for this phase, any additional phases specified,
 |      and their species and reactions to the specified file.
 |      
 |      :param filename:
 |          The name of the output file; if ``None``, a YAML string is returned
 |      :param phases:
 |          Additional ThermoPhase / Solution objects to be included in the
 |          output file
 |      :param units:
 |          A `UnitSystem` object or dictionary of the units to be used for
 |          each dimension.
 |          See `YamlWriter.output_units <cantera.YamlWriter.output_units>`.
 |      :param precision:
 |          For output floating point values, the maximum number of digits to
 |          the right of the decimal point. The default is 15 digits.
 |      :param skip_user_defined:
 |          If `True`, user-defined fields which are not used by Cantera will
 |          be stripped from the output. These additional contents can also be
 |          controlled using the `update_user_data` and `clear_user_data` functions.
 |      :param header:
 |          If `True`, fields of the `input_header` will be added to the YAML header;
 |          note that fields name ``generator``, ``cantera-version``, ``git-commit``
 |          and ``date`` are reserved, which means that any existing data are
 |          replaced by automatically generated content when the file is written.
 |  
 |  ----------------------------------------------------------------------
 |  Data descriptors inherited from cantera.solutionbase._SolutionBase:
 |  
 |  composite
 |      Returns tuple of thermo/kinetics/transport models associated with
 |      this SolutionBase object.
 |  
 |  input_data
 |      Get input data corresponding to the current state of this Solution,
 |      along with any user-specified data provided with its input (YAML)
 |      definition.
 |  
 |  input_header
 |      Retrieve input header data not associated with the current state of this
 |      Solution, which corresponds to fields at the root level of the YAML input
 |      that are not required for the instantiation of Cantera objects.
 |  
 |  name
 |      The name assigned to this object. The default value corresponds
 |      to the YAML input file phase entry.
 |  
 |  selected_species
 |      Get/set the set of species that are included when returning results that have
 |      a value for each species, such as `species_names <cantera.ThermoPhase.species_names>`,
 |      `partial_molar_enthalpies <cantera.ThermoPhase.partial_molar_enthalpies>`, or
 |      `net_production_rates <cantera.Kinetics.net_production_rates>`. The list of
 |      selected species can be set by name or index. This property returns the
 |      species by index.::
 |      
 |         >>> gas.selected_species = ["H2", "O2"]
 |         >>> print(gas.molecular_weights)
 |         [ 2.016 31.998]
 |      
 |      This method is often used implicitly by using an indexing expression on a
 |      `Solution` object::
 |      
 |         >>> print(gas["H2", "O2"].molecular_weights)
 |         [ 2.016 31.998]
 |  
 |  source
 |      The source of this object (such as a file name).

For a simple list of the properties and methods of this object:

dir(g)
Hide code cell output
['CK_mode',
 'DP',
 'DPX',
 'DPY',
 'HP',
 'HPX',
 'HPY',
 'P',
 'P_sat',
 'Pe',
 'SP',
 'SPX',
 'SPY',
 'SV',
 'SVX',
 'SVY',
 'T',
 'TD',
 'TDX',
 'TDY',
 'TP',
 'TPX',
 'TPY',
 'T_sat',
 'Te',
 'UV',
 'UVX',
 'UVY',
 'X',
 'Y',
 '_ThermoPhase__composition_to_array',
 '__call__',
 '__class__',
 '__copy__',
 '__del__',
 '__delattr__',
 '__dir__',
 '__doc__',
 '__eq__',
 '__format__',
 '__ge__',
 '__getattribute__',
 '__getitem__',
 '__getstate__',
 '__gt__',
 '__hash__',
 '__init__',
 '__init_subclass__',
 '__le__',
 '__lt__',
 '__module__',
 '__ne__',
 '__new__',
 '__pyx_vtable__',
 '__reduce__',
 '__reduce_cython__',
 '__reduce_ex__',
 '__repr__',
 '__setattr__',
 '__setstate__',
 '__setstate_cython__',
 '__sizeof__',
 '__slots__',
 '__str__',
 '__subclasshook__',
 '_check_kinetics_species_index',
 '_check_phase_index',
 '_check_reaction_index',
 '_cinit',
 '_custom_rates',
 '_enable_plasma',
 '_full_states',
 '_init_parts',
 '_init_yaml',
 '_native_mode',
 '_native_state',
 '_partial_states',
 '_references',
 'activities',
 'activity_coefficients',
 'add_reaction',
 'add_species',
 'add_species_alias',
 'atomic_weight',
 'atomic_weights',
 'basis',
 'binary_diff_coeffs',
 'case_sensitive_species_names',
 'charges',
 'chemical_potentials',
 'clear_user_data',
 'clear_user_header',
 'composite',
 'concentrations',
 'cp',
 'cp_mass',
 'cp_mole',
 'creation_rates',
 'creation_rates_ddC',
 'creation_rates_ddCi',
 'creation_rates_ddP',
 'creation_rates_ddT',
 'creation_rates_ddX',
 'critical_density',
 'critical_pressure',
 'critical_temperature',
 'cv',
 'cv_mass',
 'cv_mole',
 'delta_enthalpy',
 'delta_entropy',
 'delta_gibbs',
 'delta_standard_enthalpy',
 'delta_standard_entropy',
 'delta_standard_gibbs',
 'density',
 'density_mass',
 'density_mole',
 'derivative_settings',
 'destruction_rates',
 'destruction_rates_ddC',
 'destruction_rates_ddCi',
 'destruction_rates_ddP',
 'destruction_rates_ddT',
 'destruction_rates_ddX',
 'electric_potential',
 'electrical_conductivity',
 'electrochemical_potentials',
 'electron_energy_distribution',
 'electron_energy_distribution_type',
 'electron_energy_levels',
 'element_index',
 'element_name',
 'element_names',
 'elemental_mass_fraction',
 'elemental_mole_fraction',
 'enthalpy_mass',
 'enthalpy_mole',
 'entropy_mass',
 'entropy_mole',
 'equilibrate',
 'equilibrium_constants',
 'equivalence_ratio',
 'find_isomers',
 'forward_rate_constants',
 'forward_rate_constants_ddC',
 'forward_rate_constants_ddP',
 'forward_rate_constants_ddT',
 'forward_rates_of_progress',
 'forward_rates_of_progress_ddC',
 'forward_rates_of_progress_ddCi',
 'forward_rates_of_progress_ddP',
 'forward_rates_of_progress_ddT',
 'forward_rates_of_progress_ddX',
 'g',
 'get_binary_diff_coeffs_polynomial',
 'get_collision_integral_polynomials',
 'get_thermal_conductivity_polynomial',
 'get_viscosity_polynomial',
 'gibbs_mass',
 'gibbs_mole',
 'h',
 'has_phase_transition',
 'heat_production_rates',
 'heat_release_rate',
 'input_data',
 'input_header',
 'int_energy_mass',
 'int_energy_mole',
 'is_compressible',
 'is_pure',
 'isothermal_compressibility',
 'isotropic_shape_factor',
 'kinetics_model',
 'kinetics_species_index',
 'kinetics_species_name',
 'kinetics_species_names',
 'mass_fraction_dict',
 'max_temp',
 'mean_electron_energy',
 'mean_molecular_weight',
 'min_temp',
 'mix_diff_coeffs',
 'mix_diff_coeffs_mass',
 'mix_diff_coeffs_mole',
 'mixture_fraction',
 'mobilities',
 'modify_reaction',
 'modify_species',
 'mole_fraction_dict',
 'molecular_weights',
 'multi_diff_coeffs',
 'multiplier',
 'n_atoms',
 'n_electron_energy_levels',
 'n_elements',
 'n_phases',
 'n_reactions',
 'n_selected_species',
 'n_species',
 'n_total_species',
 'name',
 'net_production_rates',
 'net_production_rates_ddC',
 'net_production_rates_ddCi',
 'net_production_rates_ddP',
 'net_production_rates_ddT',
 'net_production_rates_ddX',
 'net_rates_of_progress',
 'net_rates_of_progress_ddC',
 'net_rates_of_progress_ddCi',
 'net_rates_of_progress_ddP',
 'net_rates_of_progress_ddT',
 'net_rates_of_progress_ddX',
 'normalize_electron_energy_distribution_enabled',
 'partial_molar_cp',
 'partial_molar_enthalpies',
 'partial_molar_entropies',
 'partial_molar_int_energies',
 'partial_molar_volumes',
 'phase_of_matter',
 'product_stoich_coeff',
 'product_stoich_coeffs',
 'product_stoich_coeffs_reversible',
 'quadrature_method',
 'reactant_stoich_coeff',
 'reactant_stoich_coeffs',
 'reaction',
 'reaction_equations',
 'reaction_phase_index',
 'reactions',
 'reference_pressure',
 'report',
 'reverse_rate_constants',
 'reverse_rates_of_progress',
 'reverse_rates_of_progress_ddC',
 'reverse_rates_of_progress_ddCi',
 'reverse_rates_of_progress_ddP',
 'reverse_rates_of_progress_ddT',
 'reverse_rates_of_progress_ddX',
 's',
 'selected_species',
 'set_binary_diff_coeffs_polynomial',
 'set_collision_integral_polynomial',
 'set_discretized_electron_energy_distribution',
 'set_equivalence_ratio',
 'set_mixture_fraction',
 'set_multiplier',
 'set_thermal_conductivity_polynomial',
 'set_unnormalized_mass_fractions',
 'set_unnormalized_mole_fractions',
 'set_viscosity_polynomial',
 'sound_speed',
 'source',
 'species',
 'species_index',
 'species_name',
 'species_names',
 'species_viscosities',
 'standard_concentration_units',
 'standard_cp_R',
 'standard_enthalpies_RT',
 'standard_entropies_R',
 'standard_gibbs_RT',
 'standard_int_energies_RT',
 'state',
 'state_size',
 'stoich_air_fuel_ratio',
 'thermal_conductivity',
 'thermal_diff_coeffs',
 'thermal_expansion_coeff',
 'thermo_model',
 'third_body_concentrations',
 'transport_model',
 'u',
 'update_user_data',
 'update_user_header',
 'v',
 'viscosity',
 'volume_mass',
 'volume_mole',
 'write_chemkin',
 'write_yaml']

To get help on a specific method, such as the species_index method:

help(g.species_index)
Help on method species_index in module cantera.thermo:

species_index(species) method of cantera.composite.Solution instance
    The index of species ``species``, which may be specified as a string or
    an integer. In the latter case, the index is checked for validity and
    returned. If no such species is present, an exception is thrown.

For properties, getting the documentation is slightly trickier, as the usual method will give you the help for the result. For example:

help(g.T)
Hide code cell output
Help on float object:

class float(object)
 |  float(x=0, /)
 |  
 |  Convert a string or number to a floating point number, if possible.
 |  
 |  Methods defined here:
 |  
 |  __abs__(self, /)
 |      abs(self)
 |  
 |  __add__(self, value, /)
 |      Return self+value.
 |  
 |  __bool__(self, /)
 |      True if self else False
 |  
 |  __ceil__(self, /)
 |      Return the ceiling as an Integral.
 |  
 |  __divmod__(self, value, /)
 |      Return divmod(self, value).
 |  
 |  __eq__(self, value, /)
 |      Return self==value.
 |  
 |  __float__(self, /)
 |      float(self)
 |  
 |  __floor__(self, /)
 |      Return the floor as an Integral.
 |  
 |  __floordiv__(self, value, /)
 |      Return self//value.
 |  
 |  __format__(self, format_spec, /)
 |      Formats the float according to format_spec.
 |  
 |  __ge__(self, value, /)
 |      Return self>=value.
 |  
 |  __getattribute__(self, name, /)
 |      Return getattr(self, name).
 |  
 |  __getnewargs__(self, /)
 |  
 |  __gt__(self, value, /)
 |      Return self>value.
 |  
 |  __hash__(self, /)
 |      Return hash(self).
 |  
 |  __int__(self, /)
 |      int(self)
 |  
 |  __le__(self, value, /)
 |      Return self<=value.
 |  
 |  __lt__(self, value, /)
 |      Return self<value.
 |  
 |  __mod__(self, value, /)
 |      Return self%value.
 |  
 |  __mul__(self, value, /)
 |      Return self*value.
 |  
 |  __ne__(self, value, /)
 |      Return self!=value.
 |  
 |  __neg__(self, /)
 |      -self
 |  
 |  __pos__(self, /)
 |      +self
 |  
 |  __pow__(self, value, mod=None, /)
 |      Return pow(self, value, mod).
 |  
 |  __radd__(self, value, /)
 |      Return value+self.
 |  
 |  __rdivmod__(self, value, /)
 |      Return divmod(value, self).
 |  
 |  __repr__(self, /)
 |      Return repr(self).
 |  
 |  __rfloordiv__(self, value, /)
 |      Return value//self.
 |  
 |  __rmod__(self, value, /)
 |      Return value%self.
 |  
 |  __rmul__(self, value, /)
 |      Return value*self.
 |  
 |  __round__(self, ndigits=None, /)
 |      Return the Integral closest to x, rounding half toward even.
 |      
 |      When an argument is passed, work like built-in round(x, ndigits).
 |  
 |  __rpow__(self, value, mod=None, /)
 |      Return pow(value, self, mod).
 |  
 |  __rsub__(self, value, /)
 |      Return value-self.
 |  
 |  __rtruediv__(self, value, /)
 |      Return value/self.
 |  
 |  __sub__(self, value, /)
 |      Return self-value.
 |  
 |  __truediv__(self, value, /)
 |      Return self/value.
 |  
 |  __trunc__(self, /)
 |      Return the Integral closest to x between 0 and x.
 |  
 |  as_integer_ratio(self, /)
 |      Return integer ratio.
 |      
 |      Return a pair of integers, whose ratio is exactly equal to the original float
 |      and with a positive denominator.
 |      
 |      Raise OverflowError on infinities and a ValueError on NaNs.
 |      
 |      >>> (10.0).as_integer_ratio()
 |      (10, 1)
 |      >>> (0.0).as_integer_ratio()
 |      (0, 1)
 |      >>> (-.25).as_integer_ratio()
 |      (-1, 4)
 |  
 |  conjugate(self, /)
 |      Return self, the complex conjugate of any float.
 |  
 |  hex(self, /)
 |      Return a hexadecimal representation of a floating-point number.
 |      
 |      >>> (-0.1).hex()
 |      '-0x1.999999999999ap-4'
 |      >>> 3.14159.hex()
 |      '0x1.921f9f01b866ep+1'
 |  
 |  is_integer(self, /)
 |      Return True if the float is an integer.
 |  
 |  ----------------------------------------------------------------------
 |  Class methods defined here:
 |  
 |  __getformat__(typestr, /) from builtins.type
 |      You probably don't want to use this function.
 |      
 |        typestr
 |          Must be 'double' or 'float'.
 |      
 |      It exists mainly to be used in Python's test suite.
 |      
 |      This function returns whichever of 'unknown', 'IEEE, big-endian' or 'IEEE,
 |      little-endian' best describes the format of floating point numbers used by the
 |      C type named by typestr.
 |  
 |  fromhex(string, /) from builtins.type
 |      Create a floating-point number from a hexadecimal string.
 |      
 |      >>> float.fromhex('0x1.ffffp10')
 |      2047.984375
 |      >>> float.fromhex('-0x1p-1074')
 |      -5e-324
 |  
 |  ----------------------------------------------------------------------
 |  Static methods defined here:
 |  
 |  __new__(*args, **kwargs) from builtins.type
 |      Create and return a new object.  See help(type) for accurate signature.
 |  
 |  ----------------------------------------------------------------------
 |  Data descriptors defined here:
 |  
 |  imag
 |      the imaginary part of a complex number
 |  
 |  real
 |      the real part of a complex number

will provide help on Python’s float class. To get the help for the temperature property, ask for the attribute of the class object itself:

help(g.__class__.T)
Help on getset descriptor cantera.thermo.ThermoPhase.T:

T
    Temperature [K].

If you are using the IPython shell, help can also be obtained using the ? syntax:

In[1]: g.species_index?

Chemical Equilibrium#

To set a gas mixture to a state of chemical equilibrium, use the ThermoPhase.equilibrate() method:

import cantera as ct
g = ct.Solution('gri30.yaml')
g.TPX = 300.0, ct.one_atm, 'CH4:0.95,O2:2,N2:7.52'
g.equilibrate('TP')

The above statement sets the state of object g to the state of chemical equilibrium holding temperature and pressure fixed. Alternatively, the specific enthalpy and pressure can be held fixed:

g.TPX = 300.0, ct.one_atm, 'CH4:0.95,O2:2,N2:7.52'
g.equilibrate('HP')

Other options are:

  • UV: fixed specific internal energy and specific volume

  • SV: fixed specific entropy and specific volume

  • SP: fixed specific entropy and pressure

How can you tell if equilibrate has correctly found the chemical equilibrium state? One way is verify that the net rates of progress of all reversible reactions are zero. Here is the code to do this:

g.TPX = 300.0, ct.one_atm, 'CH4:0.95,O2:2,N2:7.52'
g.equilibrate('HP')

rf = g.forward_rates_of_progress
rr = g.reverse_rates_of_progress
for i in range(g.n_reactions):
    if g.reaction(i).reversible and rf[i] != 0.0:
        print(' %4i  %10.4g  ' % (i, (rf[i] - rr[i])/rf[i]))
Hide code cell output
    0  -8.808e-16  
    1   1.111e-15  
    2   3.739e-15  
    3  -5.209e-15  
    4    8.33e-15  
    5   5.387e-15  
    6   5.236e-15  
    7  -5.431e-15  
    8   5.028e-15  
    9   5.077e-15  
   10   8.513e-15  
   11    7.97e-15  
   12   -1.72e-15  
   13   1.505e-15  
   14   2.791e-15  
   15   5.203e-15  
   16   5.405e-15  
   17  -3.531e-15  
   18  -1.305e-15  
   19   -2.18e-15  
   20  -1.006e-15  
   21   6.253e-15  
   22   1.064e-14  
   23   1.233e-14  
   24  -1.781e-15  
   25   5.316e-15  
   26  -1.762e-14  
   27   1.296e-14  
   28  -6.816e-15  
   29   7.145e-15  
   30   9.262e-15  
   31   1.956e-15  
   32   9.953e-15  
   33   1.015e-14  
   34   9.966e-15  
   35   9.966e-15  
   37   7.012e-15  
   38  -2.277e-15  
   39  -2.401e-15  
   40  -2.423e-15  
   41  -2.581e-15  
   42   3.011e-15  
   43  -1.752e-15  
   44  -5.306e-15  
   45  -4.992e-15  
   46   1.022e-15  
   47   2.279e-16  
   48   1.981e-15  
   49  -4.373e-15  
   50  -8.975e-15  
   51  -4.161e-15  
   52   1.347e-15  
   53  -4.325e-15  
   54  -7.039e-15  
   55  -7.608e-15  
   56    -5.9e-15  
   57   2.753e-15  
   58   1.161e-15  
   59   5.358e-15  
   60   6.685e-15  
   61   1.377e-14  
   62  -5.757e-16  
   63  -1.747e-15  
   64   3.591e-15  
   65   4.656e-15  
   66   1.165e-14  
   67  -3.484e-15  
   68  -1.702e-15  
   69  -7.812e-15  
   70  -3.591e-15  
   71     1.9e-14  
   72   3.902e-15  
   73  -1.927e-14  
   74  -2.238e-14  
   75   1.501e-14  
   76   1.593e-14  
   77   -1.85e-14  
   78   1.341e-14  
   79  -7.199e-15  
   80   1.979e-15  
   81   7.139e-16  
   82   2.875e-15  
   83   9.593e-15  
   84   1.757e-15  
   85   2.505e-15  
   86   -1.77e-15  
   87   1.095e-14  
   88   1.095e-14  
   89   3.434e-15  
   90   1.414e-14  
   91   3.323e-15  
   92   1.007e-15  
   93   1.528e-15  
   94  -6.289e-15  
   95   5.257e-15  
   96   7.135e-15  
   97   1.759e-14  
   98   1.953e-15  
   99   6.951e-15  
  100    1.25e-14  
  101   1.425e-14  
  102   1.409e-14  
  103  -8.699e-16  
  104   1.167e-15  
  105  -7.184e-15  
  106   6.157e-15  
  107   6.038e-15  
  108   1.221e-14  
  109   8.676e-15  
  110    1.44e-14  
  111  -1.246e-14  
  112  -8.757e-15  
  113  -4.051e-15  
  114  -1.271e-14  
  115  -1.273e-14  
  116  -7.259e-15  
  117  -7.377e-15  
  118  -9.833e-15  
  119  -1.762e-15  
  120   1.502e-14  
  121   1.702e-15  
  122  -1.526e-15  
  123  -7.325e-15  
  124   1.583e-14  
  125   1.025e-14  
  126   1.575e-15  
  127  -3.489e-15  
  128  -7.316e-15  
  129   1.777e-14  
  130  -4.969e-15  
  131           0  
  132  -2.988e-16  
  133           0  
  135  -2.423e-15  
  136  -1.405e-14  
  137   7.162e-15  
  138  -1.285e-15  
  139  -7.158e-15  
  140  -7.135e-15  
  141  -2.061e-15  
  143   4.048e-15  
  144   1.063e-14  
  145  -5.855e-16  
  146  -1.326e-14  
  147  -1.958e-15  
  148   8.771e-15  
  149  -2.705e-15  
  150  -1.948e-15  
  151  -8.177e-15  
  152  -1.435e-14  
  153  -2.211e-14  
  154   1.713e-15  
  155   5.103e-15  
  156   6.333e-15  
  157   7.849e-15  
  158   -7.38e-15  
  159  -1.615e-15  
  160   8.253e-15  
  161  -4.898e-15  
  162  -3.129e-15  
  163  -2.362e-14  
  164  -2.581e-14  
  165  -5.618e-15  
  166  -5.607e-15  
  167  -1.852e-15  
  168   4.266e-15  
  169    2.33e-15  
  170   3.585e-15  
  171  -6.001e-15  
  172   2.122e-14  
  173  -1.551e-14  
  174   1.612e-14  
  175   1.148e-14  
  176   7.937e-15  
  177           0  
  178   3.777e-15  
  179           0  
  180  -1.804e-15  
  181  -1.651e-16  
  182  -3.395e-15  
  183   2.556e-15  
  184  -4.498e-16  
  185  -6.941e-15  
  186   2.693e-15  
  187  -5.172e-15  
  188  -5.105e-15  
  189   1.751e-15  
  190   -1.87e-15  
  191   7.998e-15  
  192   8.029e-15  
  193   8.106e-15  
  194   5.656e-15  
  195  -3.527e-15  
  196  -8.369e-15  
  197    3.73e-15  
  198    4.98e-15  
  199  -6.266e-15  
  200   1.876e-15  
  201  -3.632e-15  
  202    2.74e-15  
  203  -3.539e-15  
  204  -3.564e-15  
  205   1.487e-16  
  206   1.975e-15  
  207  -1.737e-15  
  208   1.739e-16  
  209   1.089e-14  
  210   1.831e-16  
  211   1.622e-15  
  212  -3.438e-15  
  213  -1.053e-14  
  214   1.224e-16  
  215   8.503e-15  
  216   1.483e-15  
  217   -4.93e-16  
  218  -2.089e-15  
  219   6.373e-15  
  220   4.205e-15  
  221   1.719e-15  
  222   3.595e-15  
  223   3.627e-15  
  224           0  
  225   3.678e-15  
  226   1.065e-15  
  227   9.008e-15  
  228   1.051e-14  
  229  -2.538e-15  
  230   -8.03e-16  
  231   3.416e-15  
  232  -2.607e-16  
  233   -2.64e-15  
  234   2.575e-16  
  235   2.563e-15  
  236  -1.526e-15  
  237  -1.391e-16  
  238  -1.957e-16  
  239   8.211e-15  
  240   3.917e-15  
  241   1.687e-15  
  242  -6.453e-16  
  243  -9.985e-16  
  244   3.177e-15  
  245   5.331e-15  
  246    5.19e-15  
  247   1.492e-14  
  248   3.319e-15  
  249   3.565e-15  
  250  -3.859e-15  
  251    3.33e-15  
  252   3.555e-15  
  253  -5.533e-15  
  254  -3.849e-15  
  255  -2.854e-15  
  256   2.685e-15  
  257   5.415e-15  
  258   2.055e-14  
  259   1.496e-14  
  260   1.802e-16  
  261   5.062e-15  
  262   3.634e-15  
  263   -1.21e-15  
  264   2.761e-15  
  265  -1.863e-15  
  266   8.085e-15  
  267   4.211e-15  
  268  -1.151e-15  
  269   8.965e-15  
  270   7.137e-15  
  271   1.064e-14  
  272    2.55e-15  
  273    1.39e-14  
  274   2.607e-16  
  275           0  
  276  -7.676e-15  
  277   1.749e-15  
  278  -9.985e-16  
  279  -8.825e-16  
  280  -5.619e-15  
  281           0  
  282  -3.391e-15  
  284   2.556e-15  
  285   3.219e-15  
  286  -1.732e-15  
  288     6.4e-15  
  290    3.42e-15  
  293   1.807e-14  
  294   3.866e-15  
  295   1.703e-14  
  298    1.67e-14  
  303   1.131e-14  
  307  -3.482e-15  
  308   -7.28e-15  
  309  -1.068e-14  
  310   -4.05e-15  
  311   1.084e-14  
  312   -1.15e-14  
  313  -5.482e-15  
  314  -2.877e-15  
  315   1.271e-14  
  316  -6.703e-15  
  317  -1.861e-14  
  318   5.405e-15  
  319   1.007e-14  
  320   -9.18e-16  
  321  -1.179e-15  
  322           0  
  324  -1.507e-14  

If the magnitudes of the numbers in this list are all very small, then each reversible reaction is very nearly equilibrated, which only occurs if the gas is in chemical equilibrium.

You might be wondering how equilibrate works. (Then again, you might not). Method equilibrate invokes Cantera’s chemical equilibrium solver, which uses an element potential method. The element potential method is one of a class of equivalent nonstoichiometric methods that all have the characteristic that the problem reduces to solving a set of \(M\) nonlinear algebraic equations, where \(M\) is the number of elements (not species). The so-called stoichiometric methods, on the other hand, (including Gibbs minimization), require solving \(K\) nonlinear equations, where \(K\) is the number of species (usually \(K >> M\)). See Smith and Missen [1982] for more information on the various algorithms and their characteristics.

Cantera uses a damped Newton method to solve these equations, and does a few other things to generate a good starting guess and to produce a reasonably robust algorithm. If you want to know more about the details, look at the C++ class ChemEquil.

Chemical Kinetics#

Solution objects are also Kinetics objects, and provide all of the methods necessary to compute the thermodynamic quantities associated with each reaction, reaction rates, and species creation and destruction rates. They also provide methods to inspect the quantities that define each reaction such as the rate constants and the stoichiometric coefficients. The rate calculation functions are used extensively within Cantera’s reactor network model and 1D flame model.

Information about individual reactions that is independent of the thermodynamic state can be obtained by accessing Reaction objects with the Kinetics.reaction() method:

g = ct.Solution('gri30.yaml')
r = g.reaction(2) # get a Reaction object
r
H2 + O <=> H + OH    <Arrhenius>
r.reactants
{'H2': 1.0, 'O': 1.0}
r.products
{'H': 1.0, 'OH': 1.0}

Information about specific reaction rate parameterizations should be queried using the input_data property, which returns a YAML-compatible dictionary that represents input data needed to create the corresponding rate object:

r.rate
<ArrheniusRate at 110195d30>
r.rate.input_data
{'type': 'Arrhenius', 'rate-constant': {'A': 38.7, 'b': 2.7, 'Ea': 26191840.0}}

If we are interested in only certain types of reactions, we can use this information to filter the full list of reactions to find the just the ones of interest. For example, here we find the indices of just those reactions which convert CO into CO2:

II = [i for i,r in enumerate(g.reactions())
      if 'CO' in r.reactants and 'CO2' in r.products]
for i in II:
    print(g.reaction(i).equation)
CO + O (+M) <=> CO2 (+M)
CO + O2 <=> CO2 + O
CO + OH <=> CO2 + H
CO + HO2 <=> CO2 + OH

(Actually, we should also include reactions where the reaction is written such that CO2 is a reactant and CO is a product, but for this example, we’ll just stick to this smaller set of reactions.) Now, let’s set the composition to an interesting equilibrium state:

g.TPX = 300, 101325, {'CH4':0.6, 'O2':1.0, 'N2':3.76}
g.equilibrate('HP')

We can verify that this is an equilibrium state by seeing that the net reaction rates are essentially zero:

g.net_rates_of_progress[II]
array([-1.35525272e-20, -5.29395592e-21, -1.66533454e-15, -8.47032947e-21])

Now, let’s see what happens if we decrease the temperature of the mixture:

g.TP = g.T-100, None
g.net_rates_of_progress[II]
array([3.18644907e-05, 5.00489883e-08, 1.05964910e-01, 2.89502678e-06])

All of the reaction rates are positive, favoring the formation of CO2 from CO, with the third reaction, CO + OH <=> CO2 + H proceeding the fastest. If we look at the enthalpy change associated with each of these reactions:

g.delta_enthalpy[II]
array([-5.33034690e+08, -2.23248529e+07, -8.76650141e+07, -2.49169644e+08])

we see that the change is negative in each case, indicating a net release of thermal energy. The total heat release rate can be computed either from the reaction rates:

np.dot(g.net_rates_of_progress, g.delta_enthalpy)
-58013369.59816061

or from the species production rates:

np.dot(g.net_production_rates, g.partial_molar_enthalpies)
-58013369.598160654

The contribution from just the selected reactions is:

np.dot(g.net_rates_of_progress[II], g.delta_enthalpy[II])
-9307122.692794647

Or about 16% of the total heat release rate.

Next Steps#

Congratulations! You have finished the Cantera Python tutorial. You should now be ready to begin using Cantera on your own. Please see additional sections of the Cantera User Guide for assistance with intermediate and advanced Cantera functionality. Good luck!