Warning

This documentation is for an old version of Cantera. You can find docs for newer versions here.

One-dimensional Reacting Flows

Composite Domains

FreeFlame

class cantera.FreeFlame(gas, grid=None, width=None)

Bases: cantera.onedim.FlameBase

A freely-propagating flat flame.

A domain of type IdealGasFlow named ‘flame’ will be created to represent the flame and set to free flow. The three domains comprising the stack are stored as self.inlet, self.flame, and self.outlet.

Parameters:
  • grid – A list of points to be used as the initial grid. Not recommended unless solving only on a fixed grid; Use the width parameter instead.
  • width – Defines a grid on the interval [0, width] with internal points determined automatically by the solver.
flame
get_flame_speed_reaction_sensitivities()

Compute the normalized sensitivities of the laminar flame speed \(S_u\) with respect to the reaction rate constants \(k_i\):

\[s_i = \frac{k_i}{S_u} \frac{dS_u}{dk_i}\]
inlet
outlet
set_initial_guess(locs=[0.0, 0.3, 0.5, 1.0])

Set the initial guess for the solution. The adiabatic flame temperature and equilibrium composition are computed for the inlet gas composition.

Parameters:locs – A list of four locations to define the temperature and mass fraction profiles. Profiles rise linearly between the second and third location. Locations are given as a fraction of the entire domain
solve(loglevel=1, refine_grid=True, auto=False)

Solve the problem.

Parameters:
  • loglevel – integer flag controlling the amount of diagnostic output. Zero suppresses all output, and 5 produces very verbose output.
  • refine_grid – if True, enable grid refinement.
  • auto – if True, sequentially execute the different solution stages and attempt to automatically recover from errors. Attempts to first solve on the initial grid with energy enabled. If that does not succeed, a fixed-temperature solution will be tried followed by enabling the energy equation, and then with grid refinement enabled. If non-default tolerances have been specified or multicomponent transport is enabled, an additional solution using these options will be calculated.

BurnerFlame

class cantera.BurnerFlame(gas, grid=None, width=None)

Bases: cantera.onedim.FlameBase

A burner-stabilized flat flame.

Parameters:
  • gasSolution (using the IdealGas thermodynamic model) used to evaluate all gas properties and reaction rates.
  • grid – A list of points to be used as the initial grid. Not recommended unless solving only on a fixed grid; Use the width parameter instead.
  • width – Defines a grid on the interval [0, width] with internal points determined automatically by the solver.

A domain of class IdealGasFlow named flame will be created to represent the flame and set to axisymmetric stagnation flow. The three domains comprising the stack are stored as self.burner, self.flame, and self.outlet.

burner
flame
outlet
set_initial_guess()

Set the initial guess for the solution. The adiabatic flame temperature and equilibrium composition are computed for the burner gas composition. The temperature profile rises linearly in the first 20% of the flame to Tad, then is flat. The mass fraction profiles are set similarly.

solve(loglevel=1, refine_grid=True, auto=False)

Solve the problem.

Parameters:
  • loglevel – integer flag controlling the amount of diagnostic output. Zero suppresses all output, and 5 produces very verbose output.
  • refine_grid – if True, enable grid refinement.
  • auto – if True, sequentially execute the different solution stages and attempt to automatically recover from errors. Attempts to first solve on the initial grid with energy enabled. If that does not succeed, a fixed-temperature solution will be tried followed by enabling the energy equation, and then with grid refinement enabled. If non-default tolerances have been specified or multicomponent transport is enabled, an additional solution using these options will be calculated.

CounterflowDiffusionFlame

class cantera.CounterflowDiffusionFlame(gas, grid=None, width=None)

Bases: cantera.onedim.FlameBase

A counterflow diffusion flame

Parameters:
  • gasSolution (using the IdealGas thermodynamic model) used to evaluate all gas properties and reaction rates.
  • grid – A list of points to be used as the initial grid. Not recommended unless solving only on a fixed grid; Use the width parameter instead.
  • width – Defines a grid on the interval [0, width] with internal points determined automatically by the solver.

A domain of class IdealGasFlow named flame will be created to represent the flame and set to axisymmetric stagnation flow. The three domains comprising the stack are stored as self.fuel_inlet, self.flame, and self.oxidizer_inlet.

extinct(self)

Method overloaded for some flame types to indicate if the flame has been extinguished. Base class method always returns ‘False’

flame
fuel_inlet
mixture_fraction(m)

Compute the mixture fraction based on element m

The mixture fraction is computed from the elemental mass fraction of element m, normalized by its values on the fuel and oxidizer inlets:

\[Z = \frac{Z_{\mathrm{mass},m}(z) - Z_{\mathrm{mass},m}(z_\mathrm{oxidizer})} {Z_{\mathrm{mass},m}(z_\mathrm{fuel}) - Z_{\mathrm{mass},m}(z_\mathrm{oxidizer})} \]
Parameters:m – The element based on which the mixture fraction is computed, may be specified by name or by index
>>> f.mixture_fraction('H')
oxidizer_inlet
set_initial_guess()

Set the initial guess for the solution. The initial guess is generated by assuming infinitely-fast chemistry.

solve(loglevel=1, refine_grid=True, auto=False)

Solve the problem.

Parameters:
  • loglevel – integer flag controlling the amount of diagnostic output. Zero suppresses all output, and 5 produces very verbose output.
  • refine_grid – if True, enable grid refinement.
  • auto – if True, sequentially execute the different solution stages and attempt to automatically recover from errors. Attempts to first solve on the initial grid with energy enabled. If that does not succeed, a fixed-temperature solution will be tried followed by enabling the energy equation, and then with grid refinement enabled. If non-default tolerances have been specified or multicomponent transport is enabled, an additional solution using these options will be calculated.
strain_rate(definition, fuel=None, oxidizer='O2', stoich=None)

Return the axial strain rate of the counterflow diffusion flame in 1/s.

Parameters:
  • definition – The definition of the strain rate to be calculated. Options are: mean, max, stoichiometric, potential_flow_fuel, and potential_flow_oxidizer.
  • fuel – The fuel species. Used only if definition is stoichiometric.
  • oxidizer – The oxidizer species, default O2. Used only if definition is stoichiometric.
  • stoich – The molar stoichiometric oxidizer-to-fuel ratio. Can be omitted if the oxidizer is O2. Used only if definition is stoichiometric.

The parameter definition sets the method to compute the strain rate. Possible options are:

mean:
The mean axial velocity gradient in the entire domain
\[a_{mean} = \left| \frac{\Delta u}{\Delta z} \right| \]
max:

The maximum axial velocity gradient

\[a_{max} = \max \left( \left| \frac{du}{dz} \right| \right) \]
stoichiometric:

The axial velocity gradient at the stoichiometric surface.

\[a_{stoichiometric} = \left| \left. \frac{du}{dz} \right|_{\phi=1} \right|\]

This method uses the additional keyword arguments fuel, oxidizer, and stoich.

>>> f.strain_rate('stoichiometric', fuel='H2', oxidizer='O2',
                  stoich=0.5)
potential_flow_fuel:

The corresponding axial strain rate for a potential flow boundary condition at the fuel inlet.

\[a_{f} = \sqrt{-\frac{\Lambda}{\rho_{f}}} \]
potential_flow_oxidizer:

The corresponding axial strain rate for a potential flow boundary condition at the oxidizer inlet.

\[a_{o} = \sqrt{-\frac{\Lambda}{\rho_{o}}} \]

CounterflowPremixedFlame

class cantera.CounterflowPremixedFlame(gas, grid=None, width=None)

Bases: cantera.onedim.FlameBase

A premixed counterflow flame

Parameters:
  • gasSolution (using the IdealGas thermodynamic model) used to evaluate all gas properties and reaction rates.
  • grid – Array of initial grid points. Not recommended unless solving only on a fixed grid; Use the width parameter instead.
  • width – Defines a grid on the interval [0, width] with internal points determined automatically by the solver.

A domain of class IdealGasFlow named flame will be created to represent the flame and set to axisymmetric stagnation flow. The three domains comprising the stack are stored as self.reactants, self.flame, and self.products.

flame
products
reactants
set_initial_guess(equilibrate=True)

Set the initial guess for the solution.

If equilibrate is True, then the products composition and temperature will be set to the equilibrium state of the reactants mixture.

ImpingingJet

class cantera.ImpingingJet(gas, grid=None, width=None)

Bases: cantera.onedim.FlameBase

An axisymmetric flow impinging on a surface at normal incidence.

Parameters:
  • gasSolution (using the IdealGas thermodynamic model) used to evaluate all gas properties and reaction rates.
  • grid – A list of points to be used as the initial grid. Not recommended unless solving only on the initial grid; Use the width parameter instead.
  • width – Defines a grid on the interval [0, width] with internal points determined automatically by the solver.
  • surface – A Kinetics object used to compute any surface reactions.

A domain of class IdealGasFlow named flame will be created to represent the flame and set to axisymmetric stagnation flow. The three domains comprising the stack are stored as self.inlet, self.flame, and self.surface.

flame
inlet
set_initial_guess(products='inlet')

Set the initial guess for the solution. If products = ‘equil’, then the equilibrium composition at the adiabatic flame temperature will be used to form the initial guess. Otherwise the inlet composition will be used.

surface

IonFreeFlame

class cantera.IonFreeFlame(gas, grid=None, width=None)

Bases: cantera.onedim.IonFlameBase, cantera.onedim.FreeFlame

A freely-propagating flame with ionized gas.

E

Array containing the electric field strength at each point.

electric_field_enabled

Get/Set whether or not to solve the Poisson’s equation.

solve(self, loglevel=1, refine_grid=True, auto=False)

Solve the problem.

Parameters:
  • loglevel – integer flag controlling the amount of diagnostic output. Zero suppresses all output, and 5 produces very verbose output.
  • refine_grid – if True, enable grid refinement.
  • auto – if True, sequentially execute the different solution stages and attempt to automatically recover from errors. Attempts to first solve on the initial grid with energy enabled. If that does not succeed, a fixed-temperature solution will be tried followed by enabling the energy equation, and then with grid refinement enabled. If non-default tolerances have been specified or multicomponent transport is enabled, an additional solution using these options will be calculated.
flame
inlet
outlet

IonBurnerFlame

class cantera.IonBurnerFlame(gas, grid=None, width=None)

Bases: cantera.onedim.IonFlameBase, cantera.onedim.BurnerFlame

A burner-stabilized flat flame with ionized gas.

E

Array containing the electric field strength at each point.

electric_field_enabled

Get/Set whether or not to solve the Poisson’s equation.

solve(self, loglevel=1, refine_grid=True, auto=False)

Solve the problem.

Parameters:
  • loglevel – integer flag controlling the amount of diagnostic output. Zero suppresses all output, and 5 produces very verbose output.
  • refine_grid – if True, enable grid refinement.
  • auto – if True, sequentially execute the different solution stages and attempt to automatically recover from errors. Attempts to first solve on the initial grid with energy enabled. If that does not succeed, a fixed-temperature solution will be tried followed by enabling the energy equation, and then with grid refinement enabled. If non-default tolerances have been specified or multicomponent transport is enabled, an additional solution using these options will be calculated.
burner
flame
outlet

Flow Domains

IdealGasFlow

class cantera.IdealGasFlow(thermo)

Bases: cantera._cantera._FlowBase

An ideal gas flow domain. Functions set_free_flow and set_axisymmetric_flow can be used to set different type of flow.

For the type of axisymmetric flow, the equations solved are the similarity equations for the flow in a finite-height gap of infinite radial extent. The solution variables are:

u
axial velocity
V
radial velocity divided by radius
T
temperature
lambda
(1/r)(dP/dr)
Y_k
species mass fractions

It may be shown that if the boundary conditions on these variables are independent of radius, then a similarity solution to the exact governing equations exists in which these variables are all independent of radius. This solution holds only in in low-Mach-number limit, in which case (dP/dz) = 0, and lambda is a constant. (Lambda is treated as a spatially- varying solution variable for numerical reasons, but in the final solution it is always independent of z.) As implemented here, the governing equations assume an ideal gas mixture. Arbitrary chemistry is allowed, as well as arbitrary variation of the transport properties.

P

Pressure [Pa]

bounds(self, component)

Return the (lower, upper) bounds for a solution component.

>>> d.bounds('T')
(200.0, 5000.0)
component_index(self, str name)

Index of the component with name ‘name’

component_name(self, int n)

Name of the nth component.

component_names

List of the names of all components of this domain.

energy_enabled

Determines whether or not to solve the energy equation.

grid

The grid for this domain

have_user_tolerances

have_user_tolerances: __builtin__.bool

index

Index of this domain in a stack. Returns -1 if this domain is not part of a stack.

n_components

Number of solution components at each grid point.

n_points

Number of grid points belonging to this domain.

name

The name / id of this domain

radiation_enabled

Determines whether or not to include radiative heat transfer

set_axisymmetric_flow(self)

Set flow configuration for axisymmetric counterflow or burner-stabilized flames, using specified inlet mass fluxes.

set_boundary_emissivities(self, e_left, e_right)
set_bounds(self, *, default=None, Y=None, **kwargs)

Set the lower and upper bounds on the solution.

The argument list should consist of keyword/value pairs, with component names as keywords and (lower_bound, upper_bound) tuples as the values. The keyword default may be used to specify default bounds for all unspecified components. The keyword Y can be used to stand for all species mass fractions in flow domains.

>>> d.set_bounds(default=(0, 1), Y=(-1.0e-5, 2.0))
set_default_tolerances(self)

Set all tolerances to their default values

set_fixed_temp_profile(self, pos, T)

Set the fixed temperature profile. This profile is used whenever the energy equation is disabled.

Parameters:
  • pos – arrray of relative positions from 0 to 1
  • temp – array of temperature values
>>> d.set_fixed_temp_profile(array([0.0, 0.5, 1.0]),
...                          array([500.0, 1500.0, 2000.0])
set_free_flow(self)

Set flow configuration for freely-propagating flames, using an internal point with a fixed temperature as the condition to determine the inlet mass flux.

set_steady_tolerances(self, *, default=None, Y=None, abs=None, rel=None, **kwargs)

Set the error tolerances for the steady-state problem.

The argument list should consist of keyword/value pairs, with component names as keywords and (rtol, atol) tuples as the values. The keyword default may be used to specify default bounds for all unspecified components. The keyword Y can be used to stand for all species mass fractions in flow domains. Alternatively, the keywords abs and rel can be used to specify arrays for the absolute and relative tolerances for each solution component.

set_transient_tolerances(self, *, default=None, Y=None, abs=None, rel=None, **kwargs)

Set the error tolerances for the steady-state problem.

The argument list should consist of keyword/value pairs, with component names as keywords and (rtol, atol) tuples as the values. The keyword default may be used to specify default bounds for all unspecified components. The keyword Y can be used to stand for all species mass fractions in flow domains. Alternatively, the keywords abs and rel can be used to specify arrays for the absolute and relative tolerances for each solution component.

set_transport(self, _SolutionBase phase)

Set the Solution object used for calculating transport properties.

soret_enabled

Determines whether or not to include diffusive mass fluxes due to the Soret effect. Enabling this option works only when using the multicomponent transport model.

steady_abstol(self, component=None)

Return the absolute error tolerance for the steady state problem for a specified solution component, or all components if none is specified.

steady_reltol(self, component=None)

Return the relative error tolerance for the steady state problem for a specified solution component, or all components if none is specified.

tolerances(self, component)

Return the (relative, absolute) error tolerances for a solution component.

>>> rtol, atol = d.tolerances('u')
transient_abstol(self, component=None)

Return the absolute error tolerance for the transient problem for a specified solution component, or all components if none is specified.

transient_reltol(self, component=None)

Return the relative error tolerance for the transient problem for a specified solution component, or all components if none is specified.

IonFlow

class cantera.IonFlow(thermo)

Bases: cantera._cantera._FlowBase

An ion flow domain.

In an ion flow dommain, the electric drift is added to the diffusion flux

electric_field_enabled

Determines whether or not to solve the energy equation.

set_default_tolerances(self)
set_solving_stage(self, stage)

FreeFlow

class cantera.FreeFlow(thermo)

Bases: cantera._cantera.IdealGasFlow

FreeFlow(*args, **kwargs)

Deprecated since version 2.4: To be removed after Cantera 2.4. Use class IdealGasFlow instead and call the set_free_flow() method.

AxisymmetricStagnationFlow

class cantera.AxisymmetricStagnationFlow(thermo)

Bases: cantera._cantera.IdealGasFlow

AxisymmetricStagnationFlow(*args, **kwargs)

Deprecated since version 2.4: To be removed after Cantera 2.4. Use class IdealGasFlow instead and call the set_axisymmetric_flow() method.

Boundaries

Inlet1D

class cantera.Inlet1D(phase, name=None)

Bases: cantera._cantera.Boundary1D

A one-dimensional inlet. Note that an inlet can only be a terminal domain - it must be either the leftmost or rightmost domain in a stack.

spread_rate

Get/set the tangential velocity gradient [1/s] at this boundary.

Outlet1D

class cantera.Outlet1D(phase, name=None)

Bases: cantera._cantera.Boundary1D

A one-dimensional outlet. An outlet imposes a zero-gradient boundary condition on the flow.

OutletReservoir1D

class cantera.OutletReservoir1D(phase, name=None)

Bases: cantera._cantera.Boundary1D

A one-dimensional outlet into a reservoir.

SymmetryPlane1D

class cantera.SymmetryPlane1D(phase, name=None)

Bases: cantera._cantera.Boundary1D

A symmetry plane.

Surface1D

class cantera.Surface1D(phase, name=None)

Bases: cantera._cantera.Boundary1D

A solid surface.

ReactingSurface1D

class cantera.ReactingSurface1D(phase, name=None)

Bases: cantera._cantera.Boundary1D

A reacting solid surface.

coverage_enabled

Controls whether or not to solve the surface coverage equations.

set_kinetics(self, Kinetics kin)

Set the kinetics manager (surface reaction mechanism object).

Base Classes

Domain1D

class cantera.Domain1D(name=None)

Bases: object

Domain1D(_SolutionBase phase, name=None, *args, **kwargs)

bounds(self, component)

Return the (lower, upper) bounds for a solution component.

>>> d.bounds('T')
(200.0, 5000.0)
component_index(self, str name)

Index of the component with name ‘name’

component_name(self, int n)

Name of the nth component.

component_names

List of the names of all components of this domain.

grid

The grid for this domain

have_user_tolerances

have_user_tolerances: __builtin__.bool

index

Index of this domain in a stack. Returns -1 if this domain is not part of a stack.

n_components

Number of solution components at each grid point.

n_points

Number of grid points belonging to this domain.

name

The name / id of this domain

set_bounds(self, *, default=None, Y=None, **kwargs)

Set the lower and upper bounds on the solution.

The argument list should consist of keyword/value pairs, with component names as keywords and (lower_bound, upper_bound) tuples as the values. The keyword default may be used to specify default bounds for all unspecified components. The keyword Y can be used to stand for all species mass fractions in flow domains.

>>> d.set_bounds(default=(0, 1), Y=(-1.0e-5, 2.0))
set_default_tolerances(self)

Set all tolerances to their default values

set_steady_tolerances(self, *, default=None, Y=None, abs=None, rel=None, **kwargs)

Set the error tolerances for the steady-state problem.

The argument list should consist of keyword/value pairs, with component names as keywords and (rtol, atol) tuples as the values. The keyword default may be used to specify default bounds for all unspecified components. The keyword Y can be used to stand for all species mass fractions in flow domains. Alternatively, the keywords abs and rel can be used to specify arrays for the absolute and relative tolerances for each solution component.

set_transient_tolerances(self, *, default=None, Y=None, abs=None, rel=None, **kwargs)

Set the error tolerances for the steady-state problem.

The argument list should consist of keyword/value pairs, with component names as keywords and (rtol, atol) tuples as the values. The keyword default may be used to specify default bounds for all unspecified components. The keyword Y can be used to stand for all species mass fractions in flow domains. Alternatively, the keywords abs and rel can be used to specify arrays for the absolute and relative tolerances for each solution component.

steady_abstol(self, component=None)

Return the absolute error tolerance for the steady state problem for a specified solution component, or all components if none is specified.

steady_reltol(self, component=None)

Return the relative error tolerance for the steady state problem for a specified solution component, or all components if none is specified.

tolerances(self, component)

Return the (relative, absolute) error tolerances for a solution component.

>>> rtol, atol = d.tolerances('u')
transient_abstol(self, component=None)

Return the absolute error tolerance for the transient problem for a specified solution component, or all components if none is specified.

transient_reltol(self, component=None)

Return the relative error tolerance for the transient problem for a specified solution component, or all components if none is specified.

Boundary1D

class cantera.Boundary1D(phase, name=None)

Bases: cantera._cantera.Domain1D

Boundary1D(*args, **kwargs)

Base class for boundary domains.

Parameters:phase – The (gas) phase corresponding to the adjacent flow domain
T

The temperature [K] at this boundary.

X

Species mole fractions at this boundary. May be set as either a string or as an array.

Y

Species mass fractions at this boundary. May be set as either a string or as an array.

mdot

The mass flow rate per unit area [kg/m^2]

Sim1D

class cantera.Sim1D(domains)

Bases: object

Sim1D(domains, *args, **kwargs)

Class Sim1D is a container for one-dimensional domains. It also holds the multi-domain solution vector, and controls the process of finding the solution.

Domains are ordered left-to-right, with domain number 0 at the left.

clear_stats(self)

Clear solver statistics.

domain_index(self, dom)

Get the index of a domain, specified either by name or as a Domain1D object.

domains
eval(self, rdt=0.0)

Evaluate the governing equations using the current solution estimate, storing the residual in the array which is accessible with the work_value function.

Parameters:rdt – Reciprocal of the time-step
eval_count_stats

Return number of non-Jacobian function evaluations made in each call to solve()

eval_time_stats

Return CPU time spent on non-Jacobian function evaluations in each call to solve()

extinct(self)

Method overloaded for some flame types to indicate if the flame has been extinguished. Base class method always returns ‘False’

get_max_grid_points(self, domain)

Get the maximum number of grid points in the specified domain.

get_refine_criteria(self, domain)

Get a dictionary of the criteria used to refine one domain. The items in the dictionary are the ratio, slope, curve, and prune, as defined in set_refine_criteria.

Parameters:domain – domain object, index, or name
>>> s.set_refine_criteria(d, ratio=5.0, slope=0.2, curve=0.3, prune=0.03)
>>> s.get_refine_criteria(d)
{'ratio': 5.0, 'slope': 0.2, 'curve': 0.3, 'prune': 0.03}
grid_size_stats

Return total grid size in each call to solve()

jacobian_count_stats

Return number of Jacobian evaluations made in each call to solve()

jacobian_time_stats

Return CPU time spent evaluating Jacobians in each call to solve()

max_time_step_count

Get/Set the maximum number of time steps allowed before reaching the steady-state solution

profile(self, domain, component)

Spatial profile of one component in one domain.

Parameters:
  • domain – Domain1D object, name, or index
  • component – component name or index
>>> T = s.profile(flow, 'T')
refine(self, loglevel=1)

Refine the grid, adding points where solution is not adequately resolved.

restore(self, filename='soln.xml', name='solution', loglevel=2)

Set the solution vector to a previously-saved solution.

Parameters:
  • filename – solution file
  • name – solution name within the file
  • loglevel – Amount of logging information to display while restoring, from 0 (disabled) to 2 (most verbose).
>>> s.restore(filename='save.xml', name='energy_off')
restore_steady_solution(self)

Set the current solution vector to the last successful steady-state solution. This can be used to examine the solver progress after a failure during grid refinement.

restore_time_stepping_solution(self)

Set the current solution vector to the last successful time-stepping solution. This can be used to examine the solver progress after a failed integration.

save(self, filename='soln.xml', name='solution', description='none', loglevel=1)

Save the solution in XML format.

Parameters:
  • filename – solution file
  • name – solution name within the file
  • description – custom description text
>>> s.save(filename='save.xml', name='energy_off',
...        description='solution with energy eqn. disabled')
set_fixed_temperature(self, T)

Set the temperature used to fix the spatial location of a freely propagating flame.

set_flat_profile(self, domain, component, value)

Set a flat profile for one component in one domain.

Parameters:
  • domain – Domain1D object, name, or index
  • component – component name or index
  • v – value
>>> s.set_flat_profile(d, 'u', -3.0)
set_grid_min(self, dz, domain=None)

Set the minimum grid spacing on domain. If domain is None, then set the grid spacing for all domains.

set_initial_guess(self, *args, **kwargs)

Set the initial guess for the solution. Derived classes extend this function to set approximations for the temperature and composition profiles.

set_interrupt(self, f)

Set an interrupt function to be called each time that OneDim::eval is called. The signature of f is float f(float). The default interrupt function is used to trap KeyboardInterrupt exceptions so that ctrl-c can be used to break out of the C++ solver loop.

set_max_grid_points(self, domain, npmax)

Set the maximum number of grid points in the specified domain.

set_max_jac_age(self, ss_age, ts_age)

Set the maximum number of times the Jacobian will be used before it must be re-evaluated.

Parameters:
  • ss_age – age criterion during steady-state mode
  • ts_age – age criterion during time-stepping mode
set_max_time_step(self, tsmax)

Set the maximum time step.

set_min_time_step(self, tsmin)

Set the minimum time step.

set_profile(self, domain, component, positions, values)

Set an initial estimate for a profile of one component in one domain.

Parameters:
  • domain – Domain1D object, name, or index
  • component – component name or index
  • positions – sequence of relative positions, from 0 on the left to 1 on the right
  • values – sequence of values at the relative positions specified in positions
>>> s.set_profile(d, 'T', [0.0, 0.2, 1.0], [400.0, 800.0, 1500.0])
set_refine_criteria(self, domain, ratio=10.0, slope=0.8, curve=0.8, prune=0.05)

Set the criteria used to refine one domain.

Parameters:
  • domain – domain object, index, or name
  • ratio – additional points will be added if the ratio of the spacing on either side of a grid point exceeds this value
  • slope – maximum difference in value between two adjacent points, scaled by the maximum difference in the profile (0.0 < slope < 1.0). Adds points in regions of high slope.
  • curve – maximum difference in slope between two adjacent intervals, scaled by the maximum difference in the profile (0.0 < curve < 1.0). Adds points in regions of high curvature.
  • prune – if the slope or curve criteria are satisfied to the level of ‘prune’, the grid point is assumed not to be needed and is removed. Set prune significantly smaller than ‘slope’ and ‘curve’. Set to zero to disable pruning the grid.
>>> s.set_refine_criteria(d, ratio=5.0, slope=0.2, curve=0.3, prune=0.03)
set_steady_callback(self, f)

Set a callback function to be called after each successful steady-state solve, before regridding. The signature of f is float f(float). The argument passed to f is “0” and the output is ignored.

set_time_step(self, stepsize, n_steps)

Set the sequence of time steps to try when Newton fails.

Parameters:
  • stepsize – initial time step size [s]
  • n_steps – sequence of integer step numbers
>>> s.set_time_step(1.0e-5, [1, 2, 5, 10])
set_time_step_callback(self, f)

Set a callback function to be called after each successful timestep. The signature of f is float f(float). The argument passed to f is the size of the timestep. The output is ignored.

set_time_step_factor(self, tfactor)

Set the factor by which the time step will be increased after a successful step, or decreased after an unsuccessful one.

set_value(self, domain, component, point, value)

Set the value of one component in one domain at one point to ‘value’.

Parameters:
  • domain – Domain1D object, name, or index
  • component – component name or index
  • point – grid point number within domain starting with 0 on the left
  • value – numerical value
>>> s.set(d, 3, 5, 6.7)
>>> s.set(1, 0, 5, 6.7)
>>> s.set('flow', 'T', 5, 500)
show_solution(self)

print the current solution.

show_stats(self, print_time=True)

Show the statistics for the last solution.

If invoked with no arguments or with a non-zero argument, the timing statistics will be printed. Otherwise, the timing will not be printed.

solve(self, loglevel=1, refine_grid=True, auto=False)

Solve the problem.

Parameters:
  • loglevel – integer flag controlling the amount of diagnostic output. Zero suppresses all output, and 5 produces very verbose output.
  • refine_grid – if True, enable grid refinement.
  • auto – if True, sequentially execute the different solution stages and attempt to automatically recover from errors. Attempts to first solve on the initial grid with energy enabled. If that does not succeed, a fixed-temperature solution will be tried followed by enabling the energy equation, and then with grid refinement enabled. If non-default tolerances have been specified or multicomponent transport is enabled, an additional solution using these options will be calculated.
solve_adjoint(self, perturb, n_params, dgdx, g=None, dp=1e-5)

Find the sensitivities of an objective function using an adjoint method.

For an objective function \(g(x, p)\) where \(x\) is the state vector of the system and \(p\) is a vector of parameters, this computes the vector of sensitivities \(dg/dp\). This assumes that the system of equations has already been solved to find \(x\).

Parameters:
  • perturb

    A function with the signature perturb(sim, i, dp) which perturbs parameter i by a relative factor of dp. To perturb a reaction rate constant, this function could be defined as:

    def perturb(sim, i, dp):
        sim.gas.set_multiplier(1+dp, i)
    

    Calling perturb(sim, i, 0) should restore that parameter to its default value.

  • n_params – The length of the vector of sensitivity parameters
  • dgdx – The vector of partial derivatives of the function \(g(x, p)\) with respect to the system state \(x\).
  • g – A function with the signature value = g(sim) which computes the value of \(g(x,p)\) at the current system state. This is used to compute \(\partial g/\partial p\). If this is identically zero (i.e. \(g\) is independent of \(p\)) then this argument may be omitted.
  • dp – A relative value by which to perturb each parameter
time_step_stats

Return number of time steps taken in each call to solve()

value(self, domain, component, point)

Solution value at one point

Parameters:
  • domain – Domain1D object, name, or index
  • component – component name or index
  • point – grid point number within domain starting with 0 on the left
>>> t = s.value('flow', 'T', 6)
work_value(self, domain, component, point)

Internal work array value at one point. After calling eval, this array contains the values of the residual function.

Parameters:
  • domain – Domain1D object, name, or index
  • component – component name or index
  • point – grid point number in the domain, starting with zero at the left
>>> t = s.value(flow, 'T', 6)

FlameBase

class cantera.FlameBase(domains, gas, grid=None)

Bases: cantera._cantera.Sim1D

Base class for flames with a single flow domain

Parameters:
  • gas – object to use to evaluate all gas properties and reaction rates
  • grid – array of initial grid points
L

Array containing the radial pressure gradient (1/r)(dP/dr) [N/m^4] at each point. Note: This value is named ‘lambda’ in the C++ code.

P

Get/Set the pressure of the flame [Pa]

T

Array containing the temperature [K] at each grid point.

V

Array containing the tangential velocity gradient [1/s] at each point.

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])

Returns an array of size n_species x n_points.

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])

Returns an array of size n_species x n_points.

chemical_potentials

Array of species chemical potentials [J/kmol].

Returns an array of size n_species x n_points.

concentrations

Get/Set the species concentrations [kmol/m^3].

Returns an array of size n_species x n_points.

cp

Heat capacity at constant pressure [J/kg/K or J/kmol/K] depending on basis.

Returns an array of length n_points.

cp_mass

Specific heat capacity at constant pressure [J/kg/K].

Returns an array of length n_points.

cp_mole

Molar heat capacity at constant pressure [J/kmol/K].

Returns an array of length n_points.

creation_rates

Creation rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.

Returns an array of size n_species x n_points.

cv

Heat capacity at constant volume [J/kg/K or J/kmol/K] depending on basis.

Returns an array of length n_points.

cv_mass

Specific heat capacity at constant volume [J/kg/K].

Returns an array of length n_points.

cv_mole

Molar heat capacity at constant volume [J/kmol/K].

Returns an array of length n_points.

delta_enthalpy

Change in enthalpy for each reaction [J/kmol].

Returns an array of size n_reactions x n_points.

delta_entropy

Change in entropy for each reaction [J/kmol/K].

Returns an array of size n_reactions x n_points.

delta_gibbs

Change in Gibbs free energy for each reaction [J/kmol].

Returns an array of size n_reactions x n_points.

delta_standard_enthalpy

Change in standard-state enthalpy (independent of composition) for each reaction [J/kmol].

Returns an array of size n_reactions x n_points.

delta_standard_entropy

Change in standard-state entropy (independent of composition) for each reaction [J/kmol/K].

Returns an array of size n_reactions x n_points.

delta_standard_gibbs

Change in standard-state Gibbs free energy (independent of composition) for each reaction [J/kmol].

Returns an array of size n_reactions x n_points.

density

Density [kg/m^3 or kmol/m^3] depending on basis.

Returns an array of length n_points.

density_mass

(Mass) density [kg/m^3].

Returns an array of length n_points.

density_mole

Molar density [kmol/m^3].

Returns an array of length n_points.

destruction_rates

Destruction rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.

Returns an array of size n_species x n_points.

electrochemical_potentials

Array of species electrochemical potentials [J/kmol].

Returns an array of size n_species x n_points.

elemental_mass_fraction(m)

Get the elemental mass fraction \(Z_{\mathrm{mass},m}\) of element \(m\) at each grid point, which is defined as:

\[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, ..., 0.0]
elemental_mole_fraction(m)

Get the elemental mole fraction \(Z_{\mathrm{mole},m}\) of element \(m\) at each grid point, which is defined as:

\[Z_{\mathrm{mole},m} = \sum_k \frac{a_{m,k}}{\sum_j a_{j,k}} X_k \]

with \(a_{m,k}\) being the number of atoms of element \(m\) in species \(k\) and \(X_k\) 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, ..., 0.0]
energy_enabled

Get/Set whether or not to solve the energy equation.

enthalpy_mass

Specific enthalpy [J/kg].

Returns an array of length n_points.

enthalpy_mole

Molar enthalpy [J/kmol].

Returns an array of length n_points.

entropy_mass

Specific entropy [J/kg].

Returns an array of length n_points.

entropy_mole

Molar entropy [J/kmol/K].

Returns an array of length n_points.

equilibrium_constants

Equilibrium constants in concentration units for all reactions.

Returns an array of size n_reactions x n_points.

forward_rate_constants

Forward rate constants for all reactions. The computed values include all temperature-dependent, pressure-dependent, and third body contributions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.

Returns an array of size n_reactions x n_points.

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.

Returns an array of size n_reactions x n_points.

g

Gibbs free energy [J/kg or J/kmol] depending on basis.

Returns an array of length n_points.

gas
get_refine_criteria()

Get a dictionary of the criteria used for grid refinement. The items in the dictionary are the ratio, slope, curve, and prune, as defined in set_refine_criteria.

>>> f.set_refine_criteria(ratio=3.0, slope=0.1, curve=0.2, prune=0)
>>> f.get_refine_criteria()
{'ratio': 3.0, 'slope': 0.1, 'curve': 0.2, 'prune': 0.0}
gibbs_mass

Specific Gibbs free energy [J/kg].

Returns an array of length n_points.

gibbs_mole

Molar Gibbs free energy [J/kmol].

Returns an array of length n_points.

grid

Array of grid point positions along the flame.

h

Enthalpy [J/kg or J/kmol] depending on basis.

Returns an array of length n_points.

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

>>> f.heat_production_rates[2]  # heat production rate of the 2nd reaction
heat_release_rate

Get the total volumetric heat release rate [W/m^3].

int_energy

Internal energy in [J/kg or J/kmol].

Returns an array of length n_points.

int_energy_mass

Specific internal energy [J/kg].

Returns an array of length n_points.

int_energy_mole

Molar internal energy [J/kmol].

Returns an array of length n_points.

isothermal_compressibility

Isothermal compressibility [1/Pa].

Returns an array of length n_points.

max_grid_points

Get/Set the maximum number of grid points used in the solution of this flame.

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.

Returns an array of size n_species x n_points.

mix_diff_coeffs_mass

Mixture-averaged diffusion coefficients [m^2/s] relating the diffusive mass fluxes to gradients in the species mass fractions.

Returns an array of size n_species x n_points.

mix_diff_coeffs_mole

Mixture-averaged diffusion coefficients [m^2/s] relating the molar diffusive fluxes to gradients in the species mole fractions.

Returns an array of size n_species x n_points.

net_production_rates

Net production rates for each species. [kmol/m^3/s] for bulk phases or [kmol/m^2/s] for surface phases.

Returns an array of size n_species x n_points.

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.

Returns an array of size n_reactions x n_points.

partial_molar_cp

Array of species partial molar specific heat capacities at constant pressure [J/kmol/K].

Returns an array of size n_species x n_points.

partial_molar_enthalpies

Array of species partial molar enthalpies [J/kmol].

Returns an array of size n_species x n_points.

partial_molar_entropies

Array of species partial molar entropies [J/kmol/K].

Returns an array of size n_species x n_points.

partial_molar_int_energies

Array of species partial molar internal energies [J/kmol].

Returns an array of size n_species x n_points.

partial_molar_volumes

Array of species partial molar volumes [m^3/kmol].

Returns an array of size n_species x n_points.

radiation_enabled

Get/Set whether or not to include radiative heat transfer

reverse_rate_constants

Reverse rate constants for all reactions. The computed values include all temperature-dependent, pressure-dependent, and third body contributions. Units are a combination of kmol, m^3 and s, that depend on the rate expression for the reaction.

Returns an array of size n_reactions x n_points.

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.

Returns an array of size n_reactions x n_points.

s

Entropy [J/kg/K or J/kmol/K] depending on basis.

Returns an array of length n_points.

set_boundary_emissivities(e_left, e_right)
set_gas_state(point)

Set the state of the the Solution object used for calculations, self.gas, to the temperature and composition at the point with index point.

set_profile(component, locations, values)

Set an initial estimate for a profile of one component.

Parameters:
  • component – component name or index
  • positions – sequence of relative positions, from 0 on the left to 1 on the right
  • values – sequence of values at the relative positions specified in positions
>>> f.set_profile('T', [0.0, 0.2, 1.0], [400.0, 800.0, 1500.0])
set_refine_criteria(ratio=10.0, slope=0.8, curve=0.8, prune=0.0)

Set the criteria used for grid refinement.

Parameters:
  • ratio – additional points will be added if the ratio of the spacing on either side of a grid point exceeds this value
  • slope – maximum difference in value between two adjacent points, scaled by the maximum difference in the profile (0.0 < slope < 1.0). Adds points in regions of high slope.
  • curve – maximum difference in slope between two adjacent intervals, scaled by the maximum difference in the profile (0.0 < curve < 1.0). Adds points in regions of high curvature.
  • prune – if the slope or curve criteria are satisfied to the level of ‘prune’, the grid point is assumed not to be needed and is removed. Set prune significantly smaller than ‘slope’ and ‘curve’. Set to zero to disable pruning the grid.
>>> f.set_refine_criteria(ratio=3.0, slope=0.1, curve=0.2, prune=0)
solution(component, point=None)

Get the solution at one point or for the full flame domain (if point=None) for the specified component. The component can be specified by name or index.

soret_enabled

Get/Set whether or not to include diffusive mass fluxes due to the Soret effect. Enabling this option works only when using the multicomponent transport model.

standard_cp_R

Array of nondimensional species standard-state specific heat capacities at constant pressure at the current temperature and pressure.

Returns an array of size n_species x n_points.

standard_enthalpies_RT

Array of nondimensional species standard-state enthalpies at the current temperature and pressure.

Returns an array of size n_species x n_points.

standard_entropies_R

Array of nondimensional species standard-state entropies at the current temperature and pressure.

Returns an array of size n_species x n_points.

standard_gibbs_RT

Array of nondimensional species standard-state Gibbs free energies at the current temperature and pressure.

Returns an array of size n_species x n_points.

standard_int_energies_RT

Array of nondimensional species standard-state internal energies at the current temperature and pressure.

Returns an array of size n_species x n_points.

thermal_conductivity

Thermal conductivity. [W/m/K].

Returns an array of length n_points.

thermal_diff_coeffs

Return a one-dimensional array of the species thermal diffusion coefficients [kg/m/s].

Returns an array of size n_species x n_points.

thermal_expansion_coeff

Thermal expansion coefficient [1/K].

Returns an array of length n_points.

transport_model

Get/Set the transport model used by the Solution object used for this simulation.

u

Array containing the velocity [m/s] normal to the flame at each point.

viscosity

Viscosity [Pa-s].

Returns an array of length n_points.

volume

Specific volume [m^3/kg or m^3/kmol] depending on basis.

Returns an array of length n_points.

volume_mass

Specific volume [m^3/kg].

Returns an array of length n_points.

volume_mole

Molar volume [m^3/kmol].

Returns an array of length n_points.

write_csv(filename, species='X', quiet=True)

Write the velocity, temperature, density, and species profiles to a CSV file.

Parameters:
  • filename – Output file name
  • species – Attribute to use obtaining species profiles, e.g. X for mole fractions or Y for mass fractions.