Warning
This documentation is for an old version of Cantera. You can find docs for newer versions here.
Objects Representing Phases¶
Composite Phase Objects¶
These classes are composite representations of a substance which has thermodynamic, chemical kinetic, and (optionally) transport properties.
-
class
cantera.
Solution
(infile='', phaseid='', source=None, thermo=None, species=(), kinetics=None, reactions=(), **kwargs)¶ Bases:
cantera._cantera.ThermoPhase
,cantera._cantera.Kinetics
,cantera._cantera.Transport
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 classesThermoPhase
,Kinetics
, andTransport
. 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 theSolution
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.cti')
If an input file defines multiple phases, the phase name (in CTI) or id (in XML) can be used to specify the desired phase:
gas = ct.Solution('diamond.cti', 'gas') diamond = ct.Solution('diamond.cti', 'diamond')
Solution
objects can also be constructed usingSpecies
andReaction
objects which can themselves either be imported from input files or defined directly in Python:spec = ct.Species.listFromFile('gri30.cti') rxns = ct.Reaction.listFromFile('gri30.cti') gas = ct.Solution(thermo='IdealGas', kinetics='GasKinetics', species=spec, reactions=rxns)
where the
thermo
andkinetics
keyword arguments are strings specifying the thermodynamic and kinetics model, respectively, andspecies
andreactions
keyword arguments are lists ofSpecies
andReaction
objects, respectively.For non-trivial uses cases of this functionality, see the examples extract_submechanism.py and mechanism_reduction.py.
In addition,
Solution
objects can be constructed by passing the text of the CTI or XML phase definition in directly, using thesource
keyword argument:cti_def = ''' ideal_gas(name='gas', elements='O H Ar', species='gri30: all', reactions='gri30: all', options=['skip_undeclared_elements', 'skip_undeclared_species', 'skip_undeclared_third_bodies'], initial_state=state(temperature=300, pressure=101325))''' gas = ct.Solution(source=cti_def)
-
class
cantera.
Interface
(infile='', phaseid='', phases=(), thermo=None, species=(), kinetics=None, reactions=())¶ Bases:
cantera._cantera.InterfacePhase
,cantera._cantera.InterfaceKinetics
Two-dimensional interfaces.
Instances of class
Interface
represent reacting 2D interfaces between bulk 3D phases. ClassInterface
defines no methods of its own. All of its methods derive from eitherInterfacePhase
orInterfaceKinetics
.To construct an
Interface
object, adjacent bulk phases which participate in reactions need to be created and then passed in as a list in thephases
argument to the constructor:gas = ct.Solution('diamond.cti', 'gas') diamond = ct.Solution('diamond.cti', 'diamond') diamond_surf = ct.Interface('diamond.cti', 'diamond_100', [gas, diamond])
-
class
cantera.
DustyGas
(infile, phaseid='')¶ Bases:
cantera._cantera.ThermoPhase
,cantera._cantera.Kinetics
,cantera._cantera.DustyGasTransport
A composite class which models a gas in a stationary, solid, porous medium.
The only transport properties computed are the multicomponent diffusion coefficients. The model does not compute viscosity or thermal conductivity.
Pure Fluid Phases¶
The following convenience functions can be used to create PureFluid
objects
with the indicated equation of state:
Representing Quantities of Phases¶
-
class
cantera.
Quantity
(phase, mass=None, moles=None, constant='UV')¶ Bases:
object
A class representing a specific quantity of a
Solution
. In addition to the properties which can be computed for classSolution
, classQuantity
provides several additional capabilities. AQuantity
object is created from aSolution
with either the mass or number of moles specified:>>> gas = ct.Solution('gri30.xml') >>> gas.TPX = 300, 5e5, 'O2:1.0, N2:3.76' >>> q1 = ct.Quantity(gas, mass=5) # 5 kg of air
The state of a
Quantity
can be changed in the same way as aSolution
:>>> q1.TP = 500, 101325
Quantities have properties which provide access to extensive properties:
>>> q1.volume 7.1105094 >>> q1.enthalpy 1032237.84
The size of a
Quantity
can be changed by setting the mass or number of moles:>>> q1.moles = 3 >>> q1.mass 86.552196 >>> q1.volume 123.086
or by multiplication:
>>> q1 *= 2 >>> q1.moles 6.0
Finally, Quantities can be added, providing an easy way of calculating the state resulting from mixing two substances:
>>> q1.mass = 5 >>> q2 = ct.Quantity(gas) >>> q2.TPX = 300, 101325, 'CH4:1.0' >>> q2.mass = 1 >>> q3 = q1 + q2 # combine at constant UV >>> q3.T 432.31234 >>> q3.P 97974.9871 >>> q3.mole_fraction_dict() {'CH4': 0.26452900448117395, 'N2': 0.5809602821745349, 'O2': 0.1545107133442912}
If a different property pair should be held constant when combining, this can be specified as follows:
>>> q1.constant = q2.constant = 'HP' >>> q3 = q1 + q2 # combine at constant HP >>> q3.T 436.03320 >>> q3.P 101325.0
-
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.
-
ID
¶ The ID of the phase. The default is taken from the CTI/XML input file.
-
P
¶ Pressure [Pa].
-
P_sat
¶ Saturation pressure [Pa] at the current temperature.
-
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.
-
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])
-
add_reaction
¶ Add a new reaction to this phase.
-
add_species
¶ Add a new species to this phase. Missing elements will be added automatically.
-
atomic_weight
¶ Atomic weight [kg/kmol] of element m
-
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) ormolar
(per kmol) basis. This affects the values returned by the propertiesh
,u
,s
,g
,v
,density
,cv
, andcp
, as well as the values used with the state-setting properties such asHPX
andUV
.
-
binary_diff_coeffs
¶ Binary diffusion coefficients [m^2/s].
-
chemical_potentials
¶ Array of species chemical potentials [J/kmol].
-
concentrations
¶ Get/Set the species concentrations [kmol/m^3].
-
cp_mass
¶ Specific heat capacity at constant pressure [J/kg/K].
-
cp_mole
¶ Molar heat capacity at constant pressure [J/kmol/K].
-
creation_rates
¶ Creation rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.
-
critical_pressure
¶ Critical pressure [Pa].
-
critical_temperature
¶ Critical temperature [K].
-
cv_mass
¶ Specific heat capacity at constant volume [J/kg/K].
-
cv_mole
¶ Molar heat capacity at constant volume [J/kmol/K].
-
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].
-
density_mass
¶ (Mass) density [kg/m^3].
-
density_mole
¶ Molar density [kmol/m^3].
-
destruction_rates
¶ Destruction rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.
-
electric_potential
¶ Get/Set the electric potential [V] for this phase.
-
electrical_conductivity
¶ Electrical conductivity. [S/m].
-
electrochemical_potentials
¶ Array of species electrochemical potentials [J/kmol].
-
element_index
¶ 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
¶ Name of the element with index m.
-
element_names
¶ A list of all the element names.
-
element_potentials
¶ Get the array of element potentials. The element potentials are only defined for equilibrium states. This method first sets the composition to a state of equilibrium at constant T and P, then computes the element potentials for this equilibrium state.
-
elemental_mass_fraction
¶ Get the elemental mass fraction \(Z_{\mathrm{mass},m}\) of element \(m\) as defined by:
\[Z_{\mathrm{mass},m} = \sum_k \frac{a_{m,k} M_m}{M_k} Y_k\]with \(a_{m,k}\) being the number of atoms of element \(m\) in species \(k\), \(M_m\) the atomic weight of element \(m\), \(M_k\) the molecular weight of species \(k\), and \(Y_k\) the mass fraction of species \(k\).
Parameters: m – Base element, may be specified by name or by index. >>> phase.elemental_mass_fraction('H') 1.0
-
elemental_mole_fraction
¶ Get the elemental mole fraction \(Z_{\mathrm{mole},m}\) of element \(m\) (the number of atoms of element m divided by the total number of atoms) as defined by:
\[Z_{\mathrm{mole},m} = \frac{\sum_k a_{m,k} X_k} {\sum_k \sum_j a_{j,k} X_k}\]with \(a_{m,k}\) being the number of atoms of element \(m\) in species \(k\), \(\sum_j\) being a sum over all elements, and \(X_k\) being the mole fraction of species \(k\).
Parameters: m – Base element, may be specified by name or by index. >>> phase.elemental_mole_fraction('H') 1.0
-
enthalpy_mass
¶ Specific enthalpy [J/kg].
-
enthalpy_mole
¶ Molar enthalpy [J/kmol].
-
entropy_mass
¶ Specific entropy [J/kg].
-
entropy_mole
¶ Molar entropy [J/kmol/K].
-
equilibrate
(XY=None, *args, **kwargs)¶ Set the state to equilibrium. By default, the property pair
self.constant
is held constant. SeeThermoPhase.equilibrate
.
-
equilibrium_constants
¶ Equilibrium constants in concentration units for all reactions.
-
forward_rate_constants
¶ Forward rate constants for all reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.
-
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.
-
gibbs_mass
¶ Specific Gibbs free energy [J/kg].
-
gibbs_mole
¶ Molar Gibbs free energy [J/kmol].
-
int_energy_mass
¶ Specific internal energy [J/kg].
-
int_energy_mole
¶ Molar internal energy [J/kmol].
-
is_reversible
¶ True if reaction
i_reaction
is reversible.
-
isothermal_compressibility
¶ Isothermal compressibility [1/Pa].
-
kinetics_species_index
¶ 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.
-
mass_fraction_dict
¶
-
max_temp
¶ Maximum temperature for which the thermodynamic data for the phase are valid.
-
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.
-
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.
-
modify_reaction
¶ Modify the
Reaction
with indexirxn
to have the same rate parameters asrxn
.rxn
must have the same reactants and products and be of the same type (i.e.ElementaryReaction
,FalloffReaction
,PlogReaction
, etc.) as the existing reaction. This method does not modify the third-body efficiencies, reaction orders, or reversibility of the reaction.
-
modify_species
¶
-
mole_fraction_dict
¶
-
molecular_weights
¶ Array of species molecular weights (molar masses) [kg/kmol].
-
multi_diff_coeffs
¶ Multicomponent diffusion coefficients [m^2/s].
-
multiplier
¶ 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
.
-
n_atoms
¶ 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
-
n_elements
¶ Number of elements.
-
n_phases
¶ Number of phases in the reaction mechanism.
-
n_reactions
¶ Number of reactions in the reaction mechanism.
-
n_selected_species
¶ Number of species selected for output (by slicing of Solution object)
-
n_species
¶ Number of species.
-
n_total_species
¶ Total number of species in all phases participating in the kinetics mechanism.
-
name
¶ The name assigned to this phase. The default is taken from the CTI/XML input file.
-
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_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.
-
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
¶ Get the underlying
Solution
object, with the state set to match the wrappingQuantity
object.
-
product_stoich_coeff
¶ The stoichiometric coefficient of species k_spec as a product in reaction i_reaction.
-
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.
-
products
¶ The products portion of the reaction equation
-
reactant_stoich_coeff
¶ The stoichiometric coefficient of species k_spec as a reactant in reaction i_reaction.
-
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.
-
reactants
¶ The reactants portion of the reaction equation
-
reaction_equation
¶ The equation for the specified reaction. See also
reaction_equations
.
-
reaction_equations
¶ 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
.
-
reaction_phase_index
¶ The index of the phase where the reactions occur.
-
reaction_type
¶ Type of reaction i_reaction.
-
reference_pressure
¶ Reference state pressure [Pa].
-
report
¶ 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())
-
reverse_rate_constants
¶ Reverse rate constants for all reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.
-
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.
-
selected_species
¶
-
set_equivalence_ratio
¶ 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 and H to CO2 and H2O. Other elements are assumed not to participate in oxidation (i.e. N ends up as N2):
>>> gas.set_equivalence_ratio(0.5, 'CH4', 'O2:1.0, N2:3.76') >>> gas.mole_fraction_dict() {'CH4': 0.049900199, 'N2': 0.750499001, 'O2': 0.199600798} >>> gas.set_equivalence_ratio(1.2, {'NH3;:0.8, 'CO':0.2}, 'O2:1.0') >>> gas.mole_fraction_dict() {'CO': 0.1263157894, 'NH3': 0.505263157, 'O2': 0.36842105}
Parameters: - phi – Equivalence ratio
- fuel – Fuel species name or molar composition as string, array, or dict.
- oxidizer – Oxidizer species name or molar composition as a string, array, or dict.
-
set_multiplier
¶ 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
.
-
set_unnormalized_mass_fractions
¶ 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
¶ 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
¶ 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.
-
species_index
¶ 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
¶ Name of the species with index k.
-
species_names
¶ A list of all the species names.
-
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.
-
thermal_conductivity
¶ Thermal conductivity. [W/m/K].
-
thermal_diff_coeffs
¶ Return a one-dimensional array of the species thermal diffusion coefficients [kg/m/s].
-
thermal_expansion_coeff
¶ Thermal expansion coefficient [1/K].
-
transport_model
¶ Get/Set the transport model associated with this transport model.
Setting a new transport model deletes the underlying C++ Transport object and replaces it with a new one implementing the specified model.
-
u
¶ Internal energy in [J/kg or J/kmol].
-
viscosity
¶ Viscosity [Pa-s].
-
volume_mass
¶ Specific volume [m^3/kg].
-
volume_mole
¶ Molar volume [m^3/kmol].
-
Representing Multiple States¶
-
class
cantera.
SolutionArray
(phase, shape=(0, ), states=None, extra=None)¶ Bases:
object
A class providing a convenient interface for representing many thermodynamic states using the same
Solution
object and computing properties for that array of states.SolutionArray can represent both 1D and multi-dimensional arrays of states, with shapes described in the same way as Numpy arrays. All of the states can be set in a single call.
>>> gas = ct.Solution('gri30.cti') >>> states = ct.SolutionArray(gas, (6, 10)) >>> T = np.linspace(300, 1000, 10) # row vector >>> P = ct.one_atm * np.linspace(0.1, 5.0, 6)[:,np.newaxis] # column vector >>> X = 'CH4:1.0, O2:1.0, N2:3.76' >>> states.TPX = T, P, X
Similar to Numpy arrays, input with fewer non-singleton dimensions than the SolutionArray is ‘broadcast’ to generate input of the appropriate shape. In the above example, the single value for the mole fraction input is applied to each input, while each row has a constant temperature and each column has a constant pressure.
Computed properties are returned as Numpy arrays with the same shape as the array of states, with additional dimensions appended as necessary for non- scalar output (e.g. per-species or per-reaction properties):
>>> h = states.enthalpy_mass >>> h[i,j] # -> enthalpy at P[i] and T[j] >>> sk = states.partial_molar_entropies >>> sk[i,:,k] # -> entropy of species k at P[i] and each temperature >>> ropnet = states.net_rates_of_progress >>> ropnet[i,j,n] # -> net reaction rate for reaction n at P[i] and T[j]
In the case of 1D arrays, additional states can be appended one at a time:
>>> states = ct.SolutionArray(gas) # creates an empty SolutionArray >>> for phi in np.linspace(0.5, 2.0, 20): ... states.append(T=300, P=ct.one_atm, ... X={'CH4': phi, 'O2': 2, 'N2': 2*3.76}) >>> # 'states' now contains 20 elements >>> states.equilibrate('HP') >>> states.T # -> adiabatic flame temperature at various equivalence ratios
SolutionArray objects can also be ‘sliced’ like Numpy arrays, which can be used both for accessing and setting properties:
>>> states = ct.SolutionArray(gas, (6, 10)) >>> states[0].TP = 400, None # set the temperature of the first row to 400 K >>> cp = states[:,1].cp_mass # heat capacity of the second column
If many slices or elements of a property are going to be accessed (i.e. within a loop), it is generally more efficient to compute the property array once and access this directly, rather than repeatedly slicing the SolutionArray object, e.g.:
>>> mu = states.viscosity >>> for i,j in np.ndindex(mu.shape): ... # do something with mu[i,j]
Information about a subset of species may also be accessed, using parentheses to specify the species:
>>> states('CH4').Y # -> mass fraction of CH4 in each state >>> states('H2','O2').partial_molar_cp # -> cp for H2 and O2
Properties and functions which are not dependent on the thermodynamic state act as pass-throughs to the underlying
Solution
object, and are not converted to arrays:>>> states.element_names ['O', 'H', 'C', 'N', 'Ar'] >>> s.reaction_equation(10) 'CH4 + O <=> CH3 + OH'
Data represnted by a SolutionArray can be extracted and saved to a CSV file using the
write_csv
method:>>> states.write_csv('somefile.csv', cols=('T','P','X','net_rates_of_progress'))
Parameters: - phase – The
Solution
object used to compute the thermodynamic, kinetic, and transport properties - shape – A tuple or integer indicating the dimensions of the
SolutionArray. If the shape is 1D, the array may be extended using the
append
method. Otherwise, the shape is fixed. - states – The initial array of states. Used internally to provide slicing support.
-
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.
-
ID
¶ The ID of the phase. The default is taken from the CTI/XML input file.
-
P
¶ Pressure [Pa].
-
PV
¶ Get/Set the pressure [Pa] and specific volume [m^3/kg] of a PureFluid.
-
PX
¶ Get/Set the pressure [Pa] and vapor fraction of a two-phase state.
-
P_sat
¶ Saturation pressure [Pa] at the current temperature.
-
SH
¶ Get/Set the specific entropy [J/kg/K] and the specific enthalpy [J/kg] of a PureFluid.
-
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.
-
ST
¶ Get/Set the entropy [J/kg/K] and temperature [K] of a PureFluid.
-
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.
-
TH
¶ Get/Set the temperature [K] and the specific enthalpy [J/kg] of a PureFluid.
-
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.
-
TV
¶ Get/Set the temperature [K] and specific volume [m^3/kg] of a PureFluid.
-
TX
¶ Get/Set the temperature [K] and vapor fraction of a two-phase state.
-
T_sat
¶ Saturation temperature [K] at the current pressure.
-
UP
¶ Get/Set the specific internal energy [J/kg] and the pressure [Pa] of a PureFluid.
-
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.
-
VH
¶ Get/Set the specfic volume [m^3/kg] and the specific enthalpy [J/kg] of a PureFluid.
-
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])
-
append
(state=None, **kwargs)¶ Append an element to the array with the specified state. Elements can only be appended in cases where the array of states is one-dimensional.
The state may be specified in one of three ways:
as the array of [temperature, density, mass fractions] which is returned by
Solution.state
:mystates.append(gas.state)
as a tuple of three elements that corresponds to any of the full-state setters of
Solution
, e.g.TPY
orHPX
:mystates.append(TPX=(300, 101325, 'O2:1.0, N2:3.76'))
as separate keywords for each of the elements corresponding to one of the full-state setters:
mystates.append(T=300, P=101325, X={'O2':1.0, 'N2':3.76})
-
atomic_weight
¶ Atomic weight [kg/kmol] of element m
-
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) ormolar
(per kmol) basis. This affects the values returned by the propertiesh
,u
,s
,g
,v
,density
,cv
, andcp
, as well as the values used with the state-setting properties such asHPX
andUV
.
-
binary_diff_coeffs
¶ Binary diffusion coefficients [m^2/s].
-
chemical_potentials
¶ Array of species chemical potentials [J/kmol].
-
collect_data
(cols=('extra', 'T', 'density', 'Y'), threshold=0, species='Y')¶ Returns the data specified by cols in a single 2D Numpy array, along with a list of column labels.
Parameters: - cols – A list of any properties of the solution that are scalars or which have a value for each species or reaction. If species names are specified, then either the mass or mole fraction of that species will be taken, depending on the value of species. cols may also include any arrays which were specified as ‘extra’ variables when defining the SolutionArray object. The special value ‘extra’ can be used to include all ‘extra’ variables.
- threshold – Relative tolerance for including a particular column. The tolerance is applied by comparing the maximum absolute value for a particular column to the maximum absolute value in all columns for the same variable (e.g. mass fraction).
- species – Specifies whether to use mass (‘Y’) or mole (‘X’) fractions for individual species specified in ‘cols’
-
concentrations
¶ Get/Set the species concentrations [kmol/m^3].
-
coverages
¶ Get/Set the fraction of sites covered by each species.
-
cp_mass
¶ Specific heat capacity at constant pressure [J/kg/K].
-
cp_mole
¶ Molar heat capacity at constant pressure [J/kmol/K].
-
creation_rates
¶ Creation rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.
-
critical_pressure
¶ Critical pressure [Pa].
-
critical_temperature
¶ Critical temperature [K].
-
cv_mass
¶ Specific heat capacity at constant volume [J/kg/K].
-
cv_mole
¶ Molar heat capacity at constant volume [J/kmol/K].
-
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].
-
density_mass
¶ (Mass) density [kg/m^3].
-
density_mole
¶ Molar density [kmol/m^3].
-
destruction_rates
¶ Destruction rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.
-
electric_potential
¶ Get/Set the electric potential [V] for this phase.
-
electrical_conductivity
¶ Electrical conductivity. [S/m].
-
electrochemical_potentials
¶ Array of species electrochemical potentials [J/kmol].
-
element_index
¶ 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
¶ Name of the element with index m.
-
element_names
¶ A list of all the element names.
-
elemental_mass_fraction
(*args, **kwargs)¶
-
elemental_mole_fraction
(*args, **kwargs)¶
-
enthalpy_mass
¶ Specific enthalpy [J/kg].
-
enthalpy_mole
¶ Molar enthalpy [J/kmol].
-
entropy_mass
¶ Specific entropy [J/kg].
-
entropy_mole
¶ Molar entropy [J/kmol/K].
-
equilibrate
(*args, **kwargs)¶
-
equilibrium_constants
¶ Equilibrium constants in concentration units for all reactions.
-
forward_rate_constants
¶ Forward rate constants for all reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.
-
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.
-
gibbs_mass
¶ Specific Gibbs free energy [J/kg].
-
gibbs_mole
¶ Molar Gibbs free energy [J/kmol].
-
int_energy_mass
¶ Specific internal energy [J/kg].
-
int_energy_mole
¶ Molar internal energy [J/kmol].
-
is_reversible
¶ True if reaction
i_reaction
is reversible.
-
isothermal_compressibility
¶ Isothermal compressibility [1/Pa].
-
kinetics_species_index
¶ 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.
-
max_temp
¶ Maximum temperature for which the thermodynamic data for the phase are valid.
-
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.
-
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.
-
modify_reaction
¶ Modify the
Reaction
with indexirxn
to have the same rate parameters asrxn
.rxn
must have the same reactants and products and be of the same type (i.e.ElementaryReaction
,FalloffReaction
,PlogReaction
, etc.) as the existing reaction. This method does not modify the third-body efficiencies, reaction orders, or reversibility of the reaction.
-
molecular_weights
¶ Array of species molecular weights (molar masses) [kg/kmol].
-
multi_diff_coeffs
¶ Multicomponent diffusion coefficients [m^2/s].
-
multiplier
¶ 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
.
-
n_atoms
¶ 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
-
n_elements
¶ Number of elements.
-
n_phases
¶ Number of phases in the reaction mechanism.
-
n_reactions
¶ Number of reactions in the reaction mechanism.
-
n_species
¶ Number of species.
-
n_total_species
¶ Total number of species in all phases participating in the kinetics mechanism.
-
name
¶ The name assigned to this phase. The default is taken from the CTI/XML input file.
-
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_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.
-
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].
-
product_stoich_coeff
¶ The stoichiometric coefficient of species k_spec as a product in reaction i_reaction.
-
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.
-
products
¶ The products portion of the reaction equation
-
reactant_stoich_coeff
¶ The stoichiometric coefficient of species k_spec as a reactant in reaction i_reaction.
-
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.
-
reactants
¶ The reactants portion of the reaction equation
-
reaction_equation
¶ The equation for the specified reaction. See also
reaction_equations
.
-
reaction_equations
¶ 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
.
-
reaction_phase_index
¶ The index of the phase where the reactions occur.
-
reaction_type
¶ Type of reaction i_reaction.
-
reference_pressure
¶ Reference state pressure [Pa].
-
reverse_rate_constants
¶ Reverse rate constants for all reactions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.
-
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.
-
set_multiplier
¶ 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
.
-
site_density
¶ Get/Set the site density. [kmol/m^2] for surface phases; [kmol/m] for edge phases.
-
species
¶ 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.
-
species_index
¶ 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
¶ Name of the species with index k.
-
species_names
¶ A list of all the species names.
-
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.
-
thermal_conductivity
¶ Thermal conductivity. [W/m/K].
-
thermal_diff_coeffs
¶ Return a one-dimensional array of the species thermal diffusion coefficients [kg/m/s].
-
thermal_expansion_coeff
¶ Thermal expansion coefficient [1/K].
-
transport_model
¶ Get/Set the transport model associated with this transport model.
Setting a new transport model deletes the underlying C++ Transport object and replaces it with a new one implementing the specified model.
-
u
¶ Internal energy in [J/kg or J/kmol].
-
viscosity
¶ Viscosity [Pa-s].
-
volume_mass
¶ Specific volume [m^3/kg].
-
volume_mole
¶ Molar volume [m^3/kmol].
-
write_csv
(filename, cols=('extra', 'T', 'density', 'Y'), *args, **kwargs)¶ Write a CSV file named filename containing the data specified by cols. The first row of the CSV file will contain column labels.
Additional arguments are passed on to
collect_data
. This method works only with 1D SolutionArray objects.
- phase – The
Utility Functions¶
-
cantera.
add_directory
(directory)¶ Add a directory to search for Cantera data files.
-
cantera.
get_data_directories
()¶ Get a list of the directories Cantera searches for data files.