Warning

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

One-Dimensional Reacting Flows

Domain1D

class Domain1D(a, b, c)

Domain1D class constructor.

Parameters:
  • a

    Integer type of domain. Possible values are

    • 1 - Stagnation Flow
    • 2 - Inlet1D
    • 3 - Surf1D
    • 4 - Symm1D
    • 5 - Outlet1D
    • 6 - Reacting Surface
    • 8 - Sim1D
    • -2 - OutletRes
  • b – Instance of class Solution() (for a == 1) or Interface() (for a == 6). Not used for all other valid values of a.
  • c – Integer, either 1 or 2, indicating whether an axisymmetric stagnation flow or a free flame should be created. If not specified, defaults to 1. Ignored if a != 1.
AxiStagnFlow(gas)

Get an axisymmetric stagnation flow domain.

Parameters:gas – Instance of class Solution()
Returns:Domain1D instance representing an axisymmetric stagnation flow.
AxisymmetricFlow(gas, id)

Create an axisymmetric flow domain.

Parameters:
  • gas – Instance of class Solution()
  • id – String, ID of the flow
Returns:

Domain1D instance representing an axisymmetric flow.

Inlet(id)

Create an inlet domain.

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

Parameters:id – String name of the inlet.
Returns:Instance of class Domain1D() representing an inlet.
Outlet(id)

Create an outlet domain.

Parameters:id – String ID of the outlet.
Returns:Instance of Domain1D() representing an outlet.
OutletRes(id)

Create an outlet reservoir domain.

Returns:Instance of Domain1D() representing an outlet reservoir.
Surface(id, surface_mech)

Create a surface domain.

Parameters:
  • id – String ID of surface
  • surface_mech – Instance of class Interface() defining the surface reaction mechanism to be used. Optional.
Returns:

Instance of class Domain1D() representing a non-reacting or reacting surface.

SymmPlane(id)

Create a symmetry plane domain.

Parameters:id – String ID of the symmetry plane.
Returns:Instance of class Domain1D() representing a symmetry plane.
componentIndex(d, name)

Get the index of a component given its name.

Parameters:
  • d – Instance of class Domain1D()
  • name – String name of the component to look up. If a numeric value is passed, it will be returned.
Returns:

Index of the component, or input numeric value.

componentName(d, n)

Get the name of a component given its index.

Parameters:
  • d – Instance of class Domain1D()
  • n – Integer or vector of integers of components’ names to get.
Returns:

Cell array of component names.

disableEnergy(d)

Disable the energy equation.

Parameters:d – Instance of class Domain1D()
domainIndex(d)

Get the domain index.

Parameters:d – Instance of class Domain1D()
Returns:This function returns an integer flag denoting the location of the domain, beginning with 1 at the left.
domainType(d)

Get the type of domain.

Parameters:d – Instance of class Domain1D()
Returns:This function returns an integer flag denoting the domain type.
enableEnergy(d)

Enable the energy equation.

Parameters:d – Instance of class Domain1D()
gridPoints(d, n)

Get grid points from a domain.

Parameters:
  • d – Instance of class Domain1D()
  • n – Optional, vector of grid points to be retrieved.
Returns:

Vector of grid points. Length of n or nPoints().

isFlow(d)

Determine whether a domain is a flow.

Parameters:d – Instance of class Domain1D()
Returns:1 if the domain is a flow domain, and 0 otherwise.
isInlet(d)

Determine whether a domain is an inlet.

Parameters:d – Instance of class Domain1D()
Returns:1 if the domain is an inlet, and 0 otherwise.
isSurface(d)

Determine if a domain is a surface.

Parameters:d – Instance of class Domain1D()
Returns:1 if the domain is a surface, and 0 otherwise.
massFlux(d)

Get the mass flux.

Parameters:d – Instance of class Domain1D()
Returns:The mass flux in the domain.
massFraction(d, k)

Get the mass fraction of a species given its integer index.

This method returns the mass fraction of species k, where k is the integer index of the species in the flow domain to which the boundary domain is attached.

Parameters:
  • d – Instance of class Domain1D()
  • k – Integer species index
Returns:

Mass fraction of species

nComponents(d)

Get the number of components.

Parameters:d – Instance of class Domain1D()
Returns:Number of variables at each grid point
nPoints(d)

Get the number of grid points.

Parameters:d – Instance of class Domain1D()
Returns:Integer number of grid points.
set(a, varargin)

Set properties of a Domain1D.

The properties that may be set are

  • Temperature (T)
  • Pressure (P)
  • Mole Fractions (X)
  • Mass Flux (mdot)
  • tol
  • tol-time
  • grid
  • bounds
  • T_fixed
  • ID

Either the full property name or the symbol may be specified. Mole and mass fractions must be input as vectors (either row or column) with length equal to the number of species.

Examples:

>> set(gas,'Temperature',600.0);
>> set(gas,'T',600.0);
>> set(gas,'T',600.0,'P',2*oneatm,'Y',massfracs);
>> set(gas,'X',ones(nSpecies(gas),1));

Alternatively, individual methods to set properties may be called (setTemperature, setMoleFractions, etc.)

See also: setBounds(), setFixedTempProfile() setID(), setMdot(), setMoleFractions(), setPressure(), setProfile(), setSteadyTolerances(), setTemperature(), setTransientTolerances(), setupGrid()

Parameters:
  • a – Instance of class Domain1D()
  • varargin – Comma separated list of property, value pairs to be set
setBounds(d, component, lower, upper)

Set bounds on the solution components.

Parameters:
  • d – Instance of class Domain1D()
  • component – String, component to set the bounds on
  • lower – Lower bound
  • upper – Upper bound
setCoverageEqs(d, onoff)

Enable or disable solving the coverage equations.

Parameters:
  • d – Instance of class Domain1D()
  • onoff – String, one of 'on' or 'yes' to turn solving the coverage equations on. One of 'off' or 'no' to turn off the coverage equations.
setFixedTempProfile(d, profile)

Set a fixed temperature profile.

Set the temperature profile to use when the energy equation is not being solved. The profile must be entered as an array of positions / temperatures, which may be in rows or columns.

Parameters:
  • d – Instance of class Domain1D()
  • profile – n x 2 or 2 x n array of n points at which the temperature is specified.
setID(d, id)

Set the ID tag for a domain.

Parameters:
  • d – Instance of class Domain1D()
  • id – String ID to assign
setMdot(d, mdot)

Set the mass flow rate.

Parameters:
  • d – Instance of class Domain1D()
  • mdot – Mass flow rate
setMoleFractions(d, x)

Set the mole fractions.

Parameters:
  • d – Instance of class Domain1D()
  • x – String specifying the species and mole fractions in the format 'SPEC:X,SPEC2:X2'.
setPressure(d, p)

Set the pressure.

Parameters:
  • d – Instance of class Domain1D()
  • p – Pressure to be set. Units: Pa
setProfile(d, n, p)

Set the profile of a component.

Convenience function to allow an instance of Domain1D() to have a profile of its components set when it is part of a Stack().

Parameters:
  • d – Instance of class Domain1D()
  • n – Integer index of component, vector of component indices, string of component name, or cell array of strings of component names.
  • p – n x 2 array, whose columns are the relative (normalized) positions and the component values at those points. The number of positions n is arbitrary.
setSteadyTolerances(d, component, rtol, atol)

Set the steady-state tolerances.

Parameters:
  • d – Instance of class Domain1D()
  • component – String or cell array of strings of component values whose tolerances should be set. If 'default' is specified, the tolerance of all components will be set.
  • rtol – Relative tolerance
  • atol – Absolute tolerance
setTemperature(d, t)

Set the temperature.

Parameters:
  • d – Instance of class Domain1D()
  • t – Temperature to be set. Units: K
setTransientTolerances(d, component, rtol, atol)

Set the transient tolerances.

Parameters:
  • d – Instance of class Domain1D()
  • component – String or cell array of strings of component values whose tolerances should be set. If 'default' is specified, the tolerance of all components will be set.
  • rtol – Relative tolerance
  • atol – Absolute tolerance
setupGrid(d, grid)

Set up the solution grid.

Parameters:
temperature(d)

Get the boundary temperature.

Parameters:d – Instance of class Domain1D()
Returns:Temperature. Units: K

Stack

class Stack(domains)

Stack class constructor.

A stack object is a container for one-dimensional domains, which are instances of class Domain1D. The domains are of two types - extended domains, and connector domains.

See also: Domain1D()

Parameters:domains – Vector of domain instances
Returns:Instance of class Stack()
CounterFlowDiffusionFlame(left, flow, right, tp_f, tp_o, oxidizer)

Create a counter flow diffusion flame stack.

Parameters:
  • left – Object representing the left inlet, which must be created using function Inlet().
  • flow – Object representing the flow, created with function AxisymmetricFlow().
  • right – Object representing the right inlet, which must be created using function Inlet().
  • tp_f – Object representing the fuel inlet gas, instance of class Solution(), and an ideal gas.
  • tp_o – Object representing the oxidizer inlet gas, instance of class Solution(), and an ideal gas.
  • oxidizer – String representing the oxidizer species. Most commonly O2.
Returns:

Instance of Stack() object representing the left inlet, flow, and right inlet.

FreeFlame(gas, id)

Create a freely-propagating flat flame.

Parameters:
  • gas – Instance of class Solution()
  • id – String, ID of the flow
Returns:

Domain1D instance representing a freely propagating, adiabatic flame

npflame_init(gas, left, flow, right, fuel, oxidizer, nuox)

Create a non-premixed flame stack.

This function is deprecated in favor of CounterFlowDiffusionFlame() and will be removed after Cantera 2.4.

Parameters:
  • gas – Object representing the gas, instance of class Solution(), and an ideal gas. This object will be used to compute all required thermodynamic, kinetic, and transport properties. The state of this object should be set to an estimate of the gas state emerging from the burner before calling StagnationFlame.
  • left – Object representing the left inlet, which must be created using function Inlet().
  • flow – Object representing the flow, created with function AxisymmetricFlow().
  • right – Object representing the right inlet, which must be created using function Inlet().
  • fuel – String representing the fuel species
  • ox – String representing the oxidizer species
  • nuox – Number of oxidizer molecules required to completely combust one fuel molecule.
Returns:

Instance of Stack() object representing the left inlet, flow, and right inlet.

domainIndex(s, name)

Get the index of a domain in a stack given its name.

Parameters:
  • s – Instance of class Stack()
  • name – If double, the value is returned. Otherwise, the name is looked up and its index is returned.
Returns:

Index of domain

grid(s, name)

Get the grid in one domain.

Parameters:
  • s – Instance of class Stack()
  • name – Name of the domain for which the grid should be retrieved.
Returns:

The grid in domain name

plotSolution(s, domain, component)

Plot a specified solution component.

Parameters:
  • s – Instance of class Stack()
  • domain – Name of domain from which the component should be retrieved
  • component – Name of the component to be plotted
resid(s, domain, rdt, count)

Get the residuals.

Parameters:
  • s – Instance of class Stack()
  • domain – Name of the domain
  • rdt
  • count
Returns:

restore(s, fname, id)

Restore a previously-saved solution.

This method can be used to provide an initial guess for the solution.

See also: save()

Parameters:
  • s – Instance of class Stack()
  • fname – File name of an XML file containing solution information
  • id – ID of the element that should be restored
save(s, fname, id, desc)

Save a solution to a file.

The output file is in a format that can be used by restore()

Parameters:
  • s – Instance of class Stack()
  • fname – File name where XML file should be written
  • id – ID to be assigned to the XML element when it is written
  • desc – Description to be written to the output file
saveSoln(s, fname, id, desc)

Save a solution to a file.

The output file is in a format that can be used by restore()

Parameters:
  • s – Instance of class Stack()
  • fname – File name where XML file should be written
  • id – ID to be assigned to the XML element when it is written
  • desc – Description to be written to the output file
setFlatProfile(s, domain, comp, v)

Set a component to a value across the entire domain.

Parameters:
  • s – Instance of class Stack()
  • domain – Integer ID of the domain
  • comp – Component to be set
  • v – Double, value to be set
setMaxJacAge(s, ss_age, ts_age)

Set the number of times the Jacobian will be used before it is recomputed.

Parameters:
  • s – Instance of class Stack()
  • ss_age – Maximum age of the Jacobian for steady state analysis
  • ts_age – Maximum age of the Jacobian for transient analysis. If not specified, defaults to ss_age.
setProfile(s, name, comp, p)

Specify a profile for one component.

The solution vector values for this component will be linearly interpolated from the discrete function defined by p(:,1) vs. p(:,2). Note that p(1,1) = 0.0 corresponds to the leftmost grid point in the specified domain, and p(1,n) = 1.0 corresponds to the rightmost grid point. This method can be called at any time, but is usually used to set the initial guess for the solution.

Example (assuming s is an instance of Stack()):

>> zr = [0 0.1 0.2 0.4 0.8 1];
>> v  = [500 650 700 730 800 900];
>> setProfile(s, 1, 2, [zr, v]);
Parameters:
  • s – Instance of class Stack()
  • name – Domain name
  • comp – component number
  • p – n x 2 array, whose columns are the relative (normalized) positions and the component values at those points. The number of positions n is arbitrary.
setRefineCriteria(s, n, ratio, slope, curve, prune)

Set the criteria used to refine the grid.

Parameters:
  • s – Instance of class Stack()
  • n – Domain number
  • ratio – Maximum size ratio between adjacent cells
  • slope – Maximum relative difference in value between adjacent points
  • curve – Maximum relative difference in slope between adjacent cells
  • prune – Minimum value for slope or curve for which points will be retained in the grid. If the computed slope or curve value is below prune for all components, it will be deleted, unless either neighboring point is already marked for deletion.
setTimeStep(s, stepsize, steps)

Specify a sequence of time steps.

Parameters:
  • stepsize – Initial step size (s)
  • steps – Vector of number of steps to take before re-attempting solution of steady-state problem. For example, steps = [1, 2, 5, 10] would cause one time step to be taken first the the steady-state solution attempted. If this failed, two time steps would be taken, etc.
setValue(s, n, comp, localPoint, v)

Set the value of a single entry in the solution vector.

Example (assuming s is an instance of Stack()):

setValue(s, 3, 5, 1, 5.6)

This sets component 5 at the leftmost point (local point 1) in domain 3 to the value 5.6. Note that the local index always begins at 1 at the left of each domain, independent of the global index of the point, which depends on the location of this domain in the stack.

Parameters:
  • s – Instance of class Stack()
  • n – Domain number
  • comp – Component number
  • localPoint – Local index of the grid point in the domain
  • v – Value
solution(s, domain, component)

Get a solution component in one domain.

Parameters:
  • s – Instance of class Stack()
  • domain – String, name of the domain from which the solution is desired
  • component – String, component for which the solution is desired. If omitted, solutions for all of the components will be returned in an nPoints() x nComponents() array.
Returns:

Either an nPoints() x 1 vector, or nPoints() x nComponents() array.

solve(s, loglevel, refine_grid)

Solve the problem.

Parameters:
  • s – Instance of class Stack()
  • loglevel – Integer flag controlling the amount of diagnostic output. Zero supresses all output, and 5 produces very verbose output.
  • refine_grid – Integer, 1 to allow grid refinement, 0 to disallow.
writeStats(s)

Print statistics for the current solution.

Prints a summary of the number of function and Jacobian evaluations for each grid, and the CPU time spent on each one.

Parameters:s – Instance of class Stack()