Core classes#

Combustion Toolbox relies on several core classes to perform its calculations. These classes are the backbone of the software and are responsible for the calculations and data handling. The following sections describe the core classes and their functionalities.

Elements#

class Elements#

Bases: handle

The Elements class is used to store the list of elements that may appear in the database.

The Elements object can be initialized as follows:

elements = Elements()

This creates an instance of the Elements class and initializes it with a predefined list of elements and their indices.

See also: Species(), ChemicalSystem(), Database(), NasaDatabase()

static setElements()#

Set cell with elements name

Returns:

Tuple containing

  • elements (cell): Elements

  • NE (struct): Number of elements

setIndexStableElements()#

The only elements that are stable as diatomic gases are elements 1 (H), 9 (N), 10 (O), 11 (F), and 19 (Cl). The remaining elements that are stable as (monoatomic) gases are the noble gases He (4), Ne (12), Ar (20), Kr (38), Xe (56), and Rn (88), which do not form any compound.

Species#

class Species#

Bases: handle

The Species() class is used to store the properties of a chemical species.

See also: Database(), NasaDatabase()

T = None#

Temperature [K]

Texponents = None#

Exponents polynomials

Tintervals = None#

Number of temperature intervals

Trange = None#

Temperature range intervals

Tref = None#

Temperature reference [K]

W = None#

Molecular weight [kg/mol]

a = None#

Coefficients a polynomials

b = None#

Coefficients b polynomials

comments = None#

Additional comments from database

cpcurve = None#

Gridded interpolant object with specific heat at constant pressure of the individual species

ef = None#

Internal energy of formation [J/mol]

formula = None#

Chemical formula

fullname = None#

Fullname chemical species

g0curve = None#

Gridded interpolant object with Gibbs free energy of the individual species

h0curve = None#

Gridded interpolant object with enthalpy of the individual species

hf = None#

Enthalpy of formation at Tref from its reference species in their standard state [J/mol]

hftoh0 = None#

Enthalpy of formation at Tref relative to molar enthalpy at 0 K for standard state [J/mol]

name = None#

Name chemical species

phase = None#

Phase

refCode = None#

Reference date code

s0curve = None#

Gridded interpolant object with entropy of the individual species

getAdiabaticIndex(obj, T)#

Compute adiabatic index of the species [-] at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation

Parameters:
  • obj (Species) – Species object

  • T (float) – Temperature [K]

Returns:

gamma (float) – Adiabatic index [-]

Example

gamma = getAdiabaticIndex(obj, 300)

getEnthalpy(obj, T)#

Compute enthalpy [J/mol] of the species at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation

Parameters:
  • obj (Species) – Species object

  • T (float) – Temperature [K]

Returns:

h0 (float) – Enthalpy in molar basis [J/mol]

Example

h0 = getEnthalpy(obj, 300)

getEntropy(obj, T)#

Compute entropy [J/(mol-K)] of the species at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation

Parameters:
  • obj (Species) – Species object

  • T (float) – Temperature [K]

Returns:

s0 (float) – Entropy in molar basis [J/(mol-K)]

Example

s0 = getEntropy(obj, 300)

getGibbsEnergy(obj, T)#

Compute Gibbs energy [J/mol] of the species at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation

Parameters:
  • obj (Species) – Species object

  • T (float) – Temperature [K]

Returns:

g0 (float) – Gibbs energy in molar basis [J/mol]

Example

g0 = getGibbsEnergy(obj, 300)

getHeatCapacityPressure(obj, T)#

Compute specific heat at constant pressure [J/(mol-K)] of the species at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation

Parameters:
  • obj (Species) – Species object

  • T (float) – Temperature [K]

Returns:

cp (float) – Specific heat at constant pressure in molar basis [J/(mol-K)]

Example

cp = getHeatCapacityPressure(obj, 300)

getHeatCapacityVolume(obj, T)#

Compute specific heat at constant volume [J/(mol-K)] of the species at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation

Parameters:
  • obj (Species) – Species object

  • T (float) – Temperature [K]

Returns:

cv (float) – Specific heat at constant volume in molar basis [J/(mol-K)]

Example

cv = getHeatCapacityVolume(obj, 300)

getInternalEnergy(obj, T)#

Compute internal energy [J/mol] of the species at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation

Parameters:
  • obj (Species) – Species object

  • T (float) – Temperature [K]

Returns:

e0 (float) – Internal energy in molar basis [J/mol]

Example

e0 = getInternalEnergy(obj, 300)

getThermalEnthalpy(obj, T)#

Compute thermal enthalpy [J/mol] of the species at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation

Parameters:
  • obj (Species) – Species object

  • T (float) – Temperature [K]

Returns:

DhT (float) – Thermal enthalpy in molar basis [J/mol]

Example

DhT = getThermalEnthalpy(obj, 300)

getThermalInternalEnergy(obj, T)#

Compute thermal internal energy [J/mol] of the species at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation

Parameters:
  • obj (Species) – Species object

  • T (float) – Temperature [K]

Returns:

DeT (float) – Thermal internal energy in molar basis [J/mol]

Example

DeT = getThermalInternalEnergy(obj, 300)

ChemicalSystem#

class ChemicalSystem(database, varargin)#

Bases: handle, matlab.mixin.Copyable

The ChemicalSystem() class is used to define a chemical system with a list of species and their thermodynamic properties.

The ChemicalSystem() object can be initialized as follows:

system = ChemicalSystem(DB)
system = ChemicalSystem(DB, {'CO2', 'CO', 'H2O', 'H2', 'O2', 'N2', 'Ar'})
system = ChemicalSystem(DB, 'soot formation extended')

Here DB represents an instance of the NasaDatabase() or BurcatDatabase() class. The first case initializes the chemical system with all the species in the database that may appear as products given the initial mixture. The second case initializes the chemical system with the list of species provided. The third case initialize the chemical system with one of the predefined list of species (see setListSpecies()) that typically appear in hydrocarbon combustion systems.

See also: Database(), NasaDatabase(), findProducts(), setListSpecies()

ChemicalSystem(database, varargin)#

Constructor

FLAG_BURCAT = 'false'#

Find all the combinations of species from the database (without BURCAT’s DB) that can appear as products for the given list of reactants

FLAG_COMPLETE = 'false'#

Flag indicating to compute chemical equilibrium considering a complete combustion

FLAG_CONDENSED = 'true'#

Flag indicating to include condensed species

FLAG_ION = 'false'#

Flag indicating to include ionized species in the automatic finder of species

checkCompleteReaction(equivalenceRatio, equivalenceRatioSoot)#

Check if the list of species corresponds to “complete_reaction” If FLAG_COMPLETE is true, establish the list of species based on the given equivalence ratio (phi)

checkSpecies(species)#
clean()#

Set temperature-dependent matrix properties to zero

cleanMoles()#

Set temperature-dependent matrix properties to zero

static fillPropertiesMatrix(propertiesMatrix, species, moles, T)#

Fill the properties matrix with the data of the mixture

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • propertiesMatrix (float) – Properties matrix

  • species (cell) – Species contained in the system

  • moles (float) – Moles of the species in the mixture [mol]

  • T (float) – Temperature [K]

Returns:

propertiesMatrix (float) – Properties matrix filled

static fillPropertiesMatrixComposition(propertiesMatrix, species, moles)#

Fill the properties matrix with the data of the mixture

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • propertiesMatrix (float) – Properties matrix

  • species (cell) – Species contained in the system

  • moles (float) – Moles of the species in the mixture [mol]

Returns:

propertiesMatrix (float) – Properties matrix filled

static fillPropertiesMatrixCompositionFast(propertiesMatrix, species, moles, index)#

Fill properties matrix with the data of the mixture

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • propertiesMatrix (float) – Properties matrix

  • species (cell) – Species contained in the system

  • moles (float) – Moles of the species in the mixture [mol]

  • index (float) – Vector with the indexes of the species to fill the properties matrix

Returns:

propertiesMatrix (float) – Properties matrix filled

static fillPropertiesMatrixFast(propertiesMatrix, species, moles, T, index)#

Fill properties matrix with the data of the mixture

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • propertiesMatrix (float) – Properties matrix

  • species (cell) – Species contained in the system

  • moles (float) – Moles of the species in the mixture [mol]

  • T (float) – Temperature [K]

  • index (float) – Vector with the indexes of the species to fill the properties matrix

Returns:

propertiesMatrix (float) – Properties matrix filled

static fillPropertiesMatrixFastH0(propertiesMatrix, species, moles, T, index, h0)#

Fill properties matrix with the data of the mixture

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • propertiesMatrix (float) – Properties matrix

  • species (cell) – Species contained in the system

  • moles (float) – Moles of the species in the mixture [mol]

  • T (float) – Temperature [K]

  • index (float) – Vector with the indexes of the species to fill the properties matrix

  • h0 (float) – Enthalpy of formation vector [J/mol]

Returns:

propertiesMatrix (float) – Properties matrix filled

getIndexIons(species)#

Get index of ions for the given list of chemical species

Parameters:

species (cell) – List of chemical species

Returns:

index (float) – Index of ions

getSystemProducts()#

Set List of Species to List of Products

ind_W = '3'#

Index molecular weight

ind_cpi = '7'#

Index specific heat at constant pressure

ind_efi = '2'#

Index internal energy of formation

ind_hfi = '1'#

Index enthalpy of formation

ind_hi = '6'#

Index enthalpy

ind_ni = '5'#

Index number of moles

ind_phase = '4'#

Index phase

ind_si = '8'#

Index entropy

indexCondensed = None#

Indeces condensed species

indexCryogenic = None#

Indeces cryogenic liquified species

indexElements = None#

Index of elements

indexFrozen = None#

Indeces inert/frozen species

indexGas = None#

Indeces gaseous species

indexIons = None#

Indeces ionized species in species

indexReact = None#

Indeces react species

indexSpecies = None#

Index of species

initialization()#

Initialize chemical system

listElements = None#

List of elements

listSpecies = None#

List of species

listSpeciesLean = "{'CO2', 'H2O', 'N2', 'Ar', 'O2'}"#

List of species for a lean complete combustion (equivalence ratio < 1)

listSpeciesRich = "{'CO2', 'H2O', 'N2', 'Ar', 'CO', 'H2'}"#

List of species for a rich complete combustion (equivalence ratio > 1)

listSpeciesSoot = "{'N2', 'Ar', 'CO', 'H2', 'Cbgrb'}"#

List of species for a roch complete combustion with soot formation (equivalence ratio > equivalence ratio soot)

numElements = None#

Number of elements

numProperties = '8'#

Number of properties in propertiesMatrix

numSpecies = None#

Number of species

numSpeciesGas = None#

Number of gaseous species

plus(obj2)#

Overload the plus operator to add propertiesMatrix of two ChemicalSystem objects

propertiesMatrix = None#

Properties matrix

propertyVector = None#

Property vector

setContainedElements()#

Obtain containted elements from the given set of species (reactants and products)

Parameters:

obj (ChemicalSystem) – ChemicalSystem object

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the list of elements contained in the system

setIndexPhaseSpecies()#

Get index of gaseous, condensed and cryogenic species

Parameters:

obj (ChemicalSystem) – ChemicalSystem object

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the index of gaseous, condensed and cryogenic species

setOxidizerReference(listOxidizer)#

Set oxidizer of reference for computations with the equivalence ratio

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • listOxidizer (cell) – List of oxidizers

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the reference oxidizer

setPropertiesMatrix(species, moles, T, varargin)#

Fill the properties matrix with the data of the mixture

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • species (cell) – Species contained in the system

  • moles (float) – Moles of the species in the mixture [mol]

  • T (float) – Temperature [K]

Optional Args:

index (float): Vector with the indexes of the species to fill the properties matrix

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the properties matrix filled

Examples

  • setPropertiesMatrix(obj, {‘N2’, ‘O2’}, [3.76, 1], 300)

  • setPropertiesMatrix(obj, {‘N2’, ‘O2’}, [3.76, 1], 300, [1, 2])

setPropertiesMatrixComposition(species, moles, varargin)#

Fill the properties matrix with the data of the mixture

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • species (cell) – Species contained in the system

  • moles (float) – Moles of the species in the mixture [mol]

  • T (float) – Temperature [K]

Optional Args:

index (float): Vector with the indexes of the species to fill the properties matrix

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the properties matrix filled

Examples

  • setPropertiesMatrixComposition(obj, {‘N2’, ‘O2’}, [3.76, 1])

  • setPropertiesMatrixComposition(obj, {‘N2’, ‘O2’}, [3.76, 1], [1, 2])

setPropertiesMatrixCompositionInitialIndex(species, moles, index)#

Fill the properties matrix with the data of the mixture

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • species (cell) – Species contained in the system

  • moles (float) – Moles of the species in the mixture [mol]

  • T (float) – Temperature [K]

  • index (float) – Vector with the indexes of the species to fill the properties matrix

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the properties matrix filled

Example

setPropertiesMatrixCompositionInitialIndex(obj, {‘N2’, ‘O2’}, [3.76, 1], [1, 2])

setPropertiesMatrixInitialIndex(species, moles, T, index, varargin)#

Fill the properties matrix with the data of the mixture

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • species (cell) – Species contained in the system

  • moles (float) – Moles of the species in the mixture [mol]

  • T (float) – Temperature [K]

  • index (float) – Vector with the indexes of the species to fill the properties matrix

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the properties matrix filled

Examples

  • setPropertiesMatrix(obj, {‘N2’, ‘O2’}, [3.76, 1], 300)

  • setPropertiesMatrix(obj, {‘N2’, ‘O2’}, [3.76, 1], 300, [1, 2])

setPropertiesMatrixInitialize()#

Initialize properties matrix

Parameters:

obj (ChemicalSystem) – ChemicalSystem object

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the properties matrix initialized

setReactIndex(speciesFrozen)#

Set index of react (non-frozen) and frozen species

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • speciesFrozen (char) – Frozen chemical species

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the index of react and frozen species

setStoichiometricMatrix()#

Set stoichiometric matrix of the chemical system

Parameters:

obj (ChemicalSystem) – ChemicalSystem object

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the stoichiometric matrix set

sortIndexPhaseSpecies()#

Reorginize index of gaseous, condensed and cryogenic species

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • LS (cell) – Name list species / list of species

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the index of gaseous, condensed and cryogenic species sorted

sortListSpecies()#

Establish cataloged list of species according to the state of the phase (gaseous or condensed). It also obtains the indices of cryogenic liquid species, i.e., liquified gases.

Parameters:

obj (ChemicalSystem) – ChemicalSystem object

Returns:

obj (ChemicalSystem) – ChemicalSystem object with the list of species sorted

species = None#

Struct with Species objects

stoichiometricMatrix = None#

Stoichiometric matrix

findProducts(obj, listReactants, varargin)#

Find all the combinations of species from DB that can appear as products for the given list of reactants

Parameters:
  • obj (ChemicalSystem) – ChemicalSystem object

  • listReactants (cell) – List of reactants

Optional Name-Value Pairs Args:
  • indexElements_DB (float): Matrix NS_DB x MAX_ELEMENTS with element indeces of the species contained in the database

  • FLAG_BURCAT (bool): Flag indicating to look for species also in Burcat’s database

  • FLAG_ION (bool): Flag indicating to include ionized species

  • FLAG_CONDENSED (bool): Flag indicating to include condensed species

Returns:

Tuple containing

  • listSpecies (cell): List of products

  • indexElements_DB (float): Matrix NS_DB x MAX_ELEMENTS with element indeces of the species contained in the database

Examples

  • [listSpecies, indexElements_DB] = findProducts(ChemicalSystem(DB), {‘O2’, ‘N’, ‘eminus’})

  • [listSpecies, indexElements_DB] = findProducts(ChemicalSystem(DB), {‘O2’, ‘CO’, ‘N’}, ‘flag_burcat’, true)

  • [listSpecies, indexElements_DB] = findProducts(ChemicalSystem(DB), {‘O2’, ‘CO’, ‘N’}, ‘flag_burcat’, true, ‘flag_ion’, true)

  • [listSpecies, indexElements_DB] = findProducts(ChemicalSystem(DB), {‘O2’, ‘CO’, ‘N’}, ‘flag_burcat’, true, ‘flag_ion’, true, ‘flag_condensed’, true, ‘ind’, indexElements_DB)

  • [listSpecies, indexElements_DB] = findProducts(ChemicalSystem(DB), {‘O2’, ‘CO’, ‘N’}, ‘flag_burcat’, true, ‘flag_ion’, true, ‘ind’, indexElements_DB)

setListSpecies(obj, varargin)#

Set list of species in the mixture (products)

Predefined list of species:
  • SOOT FORMATION

  • COMPLETE

  • HC/O2/N2 EXTENDED

  • SOOT FORMATION EXTENDED

  • NASA ALL

  • NASA ALL CONDENSED

  • NASA ALL IONS

  • AIR, DISSOCIATED AIR

  • AIR IONS, AIR_IONS

  • IDEAL_AIR, AIR_IDEAL

  • HYDROGEN

  • HYDROGEN_L, HYDROGEN (L)

  • HC/O2/N2 PROPELLANTS

  • SI/HC/O2/N2 PROPELLANTS

Parameters:

obj (ChemicalSystem) – Chemical system object

Optional Args:
  • listSpecies (cell): Name list species / list of species

  • phi (float): Equivalence ratio

  • phi_c (float): Equivalence ratio in which theoretically appears soot

Returns:

Tuple containing

  • obj (ChemicalSystem): Chemical system object

  • listSpecies (cell): List of species

  • listSpeciesFormula (cell): List with chemical formula of species

Examples

  • [obj, listSpecies, listSpeciesFormula] = setListSpecies(chemicalSystem, ‘soot formation’);

  • [obj, listSpecies, listSpeciesFormula] = setListSpecies(chemicalSystem, ‘soot formation’);

  • [obj, listSpecies, listSpeciesFormula] = setListSpecies(chemicalSystem, ‘complete’, 1.5, 2.5);

Mixture#

class Mixture(chemicalSystem, varargin)#

Bases: handle, matlab.mixin.Copyable

The Mixture class is used to store the properties of a chemical mixture.

The Mixture object can be initialized as follows:

mix = Mixture(chemicalSystem)

This creates an instance of the Mixture class and initializes it with a predefined chemical system.

See also: ChemicalSystem(), Database(), NasaDatabase()

DeT = None#

Thermal internal energy [J]

DhT = None#

Thermal enthalpy [J]

Ds = None#

Entropy of mixing [J/K]

FLAG_TSPECIES = 'false'#

Flag to indicate species-specific initial temperatures are defined (initial mixture)

FLAG_VOLUME = 'false'#

Flag to indicate specific volume is defined (initial mixture)

I_sp = None#

Specific impulse [s]

I_vac = None#

Vacuum impulse [s]

MW = None#

Mean molecular weight [kg/mol]

Mixture(chemicalSystem, varargin)#

Mixture constructor

N = None#

Total number of moles [mol]

T = None#

Temperature [K]

Tspecies = None#

Species-specific initial temperatures [K] (initial mixture)

W = None#

Molecular weight [kg/mol]

Xi = None#

Molar fractions [-]

Yi = None#

Mass fractions [-]

areaRatio = None#

Area ratio = area_i / areaThroat

areaRatioChamber = None#

Area ratio = areaChamber / areaThroat

assignAtomElementsFuel(natomElementsFuel)#

Assign atomic counts of the fuel

Parameters:
  • obj (Mixture) – Mixture object

  • natomElementsFuel (float) – vector containing atomic counts for the fuel

Returns:

obj (Mixture) – Mixture object with updated atomic counts for the fuel

beta = None#

Wave angle [deg]

betaMax = None#

Maximum wave angle [deg]

betaMin = None#

Minimum wave angle [deg]

betaSonic = None#

Wave angle at the sonic point [deg]

cf = None#

Thrust coefficient [-]

chemicalSystem = None#

Chemical system object

cjSpeed = None#

Chapman-Jouguet speed

computeComposition()#

Compute the composition of the mixture

computeEntropyMixing(Ni, N_gas, R0, FLAG_NONZERO)#

Compute entropy of mixing [J/K]

Parameters:
  • obj (Mixture) – Mixture object

  • Ni (float) – Vector of moles of each species [mol]

  • N_gas (float) – Total moles of gas species [mol]

  • R0 (float) – Universal gas constant [J/(mol-K)]

  • FLAG_NONZERO (bool) – Vector of nonzero species

Returns:

DS (float) – Entropy of mixing [J/K]

Note

only nonzero for gaseous species

Example

Ds = computeEntropyMixing(mix, Ni, N_gas, R0, FLAG_NONZERO)

computeEquilibriumTemperature(speciesTemperatures)#

Compute the equilibrium temperature [K] for a mixture of n species, each initially at a specified temperature. A Newton-Raphson method is used to iteratively adjust the temperature until equilibrium is reached

Parameters:

speciesTemperatures (float) – Array containing the initial temperatures [K] for each species

Returns:

T (float) – Equilibrium temperature [K] of the mixture

Example

T = obj.computeEquilibriumTemperature([300, 400, 350])

computeEquivalenceRatio()#

Compute equivalence ratio [-]

Parameters:

obj (Mixture) – Mixture object

Returns:

obj (Mixture) – Mixture object with updated equivalence ratio [-]

computeEquivalenceRatioSoot()#

Compute guess of equivalence ratio in which soot appears considering complete combustion

Parameters:

obj (Mixture) – Mixture object

Returns:

obj (Mixture) – Mixture object with theoretical equivalence ratio at which soot appears [-]

computeMeanMolecularWeight(moles, index)#

Compute Mean Molecular Weight [kg/mol]

computeProperties()#

Compute composition and thermodynamic properties of the mixture

Parameters:

obj (Mixture) – Mixture object

Returns:

obj (Mixture) – Mixture object with the computed properties

Example

mix = computeProperties(obj)

computeRatiosFuelOxidizer(propertiesMatrixFuel, propertiesMatrixOxidizer)#

Compute percentage Fuel, Oxidizer/Fuel ratio and equivalence ratio

Parameters:
  • obj (Mixture) – Mixture object

  • propertiesMatrixFuel (float) – Properties matrix of the fuel

  • propertiesMatrixOxidizer (float) – Properties matrix of the oxidizer

Returns:

obj (Mixture) – Mixture object with updated properties

computeThermodynamics()#

Compute thermodynamic properties of the mixture

Parameters:

obj (Mixture) – Mixture object

Returns:

obj (Mixture) – Mixture object with the computed properties

Example

mix = computeThermodynamics(obj)

config = None#

Mixture configuration object

cp = None#

Specific heat at constant pressure [J/K]

cstar = None#

Characteristic velocity [m/s]

cv = None#

Specific heat at constant volume [J/K]

dVdT_p = None#

Derivative of volume with respect to temperature at constant pressure [-]

dVdp_T = None#

Derivative of volume with respect to pressure at constant temperature [-]

defineF()#

Set Fuel of the mixture

Parameters:

obj (Mixture) – Mixture object

Returns:

obj (Mixture) – Mixture object

defineO()#

Set Oxidizer of the mixture

Parameters:

obj (Mixture) – Mixture object

Returns:

obj (Mixture) – Mixture object

driveFactor = None#

Overdriven/Underdriven factor (detonations)

e = None#

Internal energy [J]

ef = None#

Internal energy of formation [J]

equationState = None#

Equation of State object

equivalenceRatio = None#

Equivalence ratio [-]

equivalenceRatioSoot = None#

Theoretical equivalence ratio at which soot may appear [-]

fuelOxidizerMassRatio = None#

Mass ratio of oxidizer to fuel [-]

g = None#

Gibbs energy [J]

gamma = None#

Adiabatic index [-]

gamma_s = None#

Adiabatic index [-]

static getMass(system, propertiesMatrix)#

Compute mass mixture [kg]

getPropertyListSpecies(fun, temperatures)#

This helper method evaluates a given property function (e.g., specific heat, enthalpy) for each species in the mixture at the specified temperature(s).

Parameters:
  • fun (function handle) – Function to compute the property (e.g., @getHeatCapacityPressure, @getHeatCapacityVolume, @getEnthalpy, or @getInternalEnergy).

  • temperatures (float) – Temperature(s) [K] at which to evaluate the property. If provided as a vector, each species is evaluated at its corresponding temperature.

Returns:

value (vector) – Vector containing the computed property for each species.

Example

cp_values = obj.getPropertyListSpecies(@getHeatCapacityPressure, [300, 400, 350]);

getTypeSpecies()#

Create cell array with the type of species in the mixture

Parameters:

obj (Mixture) – Mixture class

Returns:

typeSpecies (cell) – Cell array with the type of species in the mixture

h = None#

Enthalpy [J]

hf = None#

Enthalpy of formation [J]

indexGas = None#

Index of the gas species (initial mixture)

indexMax = None#

Index of the maximum deflection angle

indexMin = None#

Index of the minimum deflection angle

indexSonic = None#

Index of the sonic point

indexSpecies = None#

Index of the species (initial mixture)

mach = None#

Mach number [-]

mi = None#

Mass mixture [kg]

natomElements = None#

Vector atoms of each element [-]

natomElementsReact = None#

Vector atoms of each element without frozen species [-]

oxidizerFuelMassRatio = None#

Mass ratio of fuel to oxidizer [-]

p = None#

Pressure [bar]

percentageFuel = None#

Percentage of fuel in the mixture [%]

phase = None#

Phase vector [-]

polar = None#

Properties of the polar solution

rho = None#

Density [kg/m3]

s = None#

Entropy [J/K]

s0 = None#

Entropy (frozen) [J/K]

setEquivalenceRatio(equivalenceRatio)#

Set equivalence ratio and compute thermodynamic properties

Parameters:
  • obj (Mixture) – Mixture object

  • equivalenceRatio (float) – Equivalence ratio [-]

Returns:

obj (Mixture) – Mixture object with updated properties

Example

setEquivalenceRatio(obj, 1)

setPressure(p, varargin)#

Set pressure [bar] and compute thermodynamic properties

Parameters:
  • obj (Mixture) – Mixture object

  • p (float) – Pressure [bar]

Returns:

obj (Mixture) – Mixture object with updated properties

Example

setPressure(obj, 1)

setProperties(property, value, varargin)#

Obtain properties at equilibrium for the given thermochemical transformation

Parameters:
  • obj (Mixture) – Mixture object

  • property (char) – Property to be set

  • value (float) – Value of the property

Optional Args (key-value pairs):
  • property (char): Property to be set

  • value (float): Value of the property

Returns:

objArray (Mixture) – Array of Mixture objects with the computed properties

Note

  • Use this method after setting the initial composition of the mixture

Examples

  • mixArray = setProperties(mix, ‘equivalenceRatio’, value)

  • mixArray = setProperties(mix, ‘equivalenceRatio’, value, ‘temperature’, value)

  • mixArray = setProperties(mix, ‘equivalenceRatio’, value, ‘temperature’, value, ‘pressure’, value)

  • mixArray = setProperties(mix, ‘phi’, value, ‘T’, value, ‘p’, value)

setTemperature(T, varargin)#

Set temperature [K] and compute thermodynamic properties

Parameters:
  • obj (Mixture) – Mixture object

  • T (float) – Temperature [K]

Returns:

obj (Mixture) – Mixture object with updated properties

Example

setTemperature(obj, 300)

setTemperatureSpecies(speciesTemperatures)#

Set species-specific temperatures and update equilibrium temperature and properties

Parameters:

speciesTemperatures (float) – Array with the temperatures [K] for each species in the mixture

Returns:

obj (Mixture) – Mixture object with updated equilibrium temperature and properties

Note: The temperature vector is assigned to the species in the same order as they were added to the mixture object

Example

setTemperatureSpecies(obj, [300, 400, 350])

setVolume(vSpecific, varargin)#

Set specific volume [m3/kg] and compute thermodynamic properties

Parameters:
  • obj (Mixture) – Mixture object

  • vSpecific (float) – Specific volume [m3/kg]

Returns:

obj (Mixture) – Mixture object with updated properties

Example

setVolume(obj, 1)

sound = None#

Speed of sound [m/s]

stoichiometricMoles = None#

Theoretical moles of the oxidizer of reference for a stoichiometric combustion

theta = None#

Deflection angle [deg]

thetaMax = None#

Maximum deflection angle [deg]

thetaMin = None#

Minimum eflection angle [deg]

thetaSonic = None#

Deflection angle at the sonic point [deg]

u = None#

Velocity module relative to the shock front [m/s]

uNormal = None#

Normal component of u [m/s]

uShock = None#

Velocity module in the shock tube [m/s]

updateComposition()#

Update the composition of the mixture

Parameters:

obj (Mixture) – Mixture object

Returns:

obj (Mixture) – Mixture object with updated properties

updateIndexSpecies()#

Update index species in the mixture

Parameters:

obj (Mixture) – Mixture object

Returns:

obj (Mixture) – Mixture object with updated indexSpecies

updateThermodynamics()#

Update the thermodynamic state of the mixture

Parameters:

obj (Mixture) – Mixture object

Returns:

obj (Mixture) – Mixture object with updated properties

v = None#

Volume [m3]

vSpecific = None#

Specific volume [m3/kg]

vSpecific2vMolar(vSpecific, moles, molesGas, varargin)#

Compute molar volume [m3/mol] from specific volume [m3/kg]

print(mix, varargin)#

Print properties and composition of the given mixtures in the command window

Parameters:

mix1 (Mixture) – Mixture object with the properties of the mixture

Optional Args:
  • mix2 (Mixture): Mixture object with the properties of the mixture

  • mixi (Mixture): Mixture object with the properties of the mixture

  • mixN (Mixture): Mixture object with the properties of the mixture

Examples

  • print(mix1)

  • print(mix1, mix2)

  • print(mix1, mix2, mix3)

  • print(mix1, mix2, mix3, mix4)

setStagnation(mix, varargin)#

Update mixture properties to stagnation conditions

Parameters:

mix (Mixture) – Mixture object

Returns:

mix (Mixture) – Mixture object at stagnation conditions

Example

mix = setStagnation(mix)

MixtureConfig#

class MixtureConfig(varargin)#

Bases: handle

The MixtureConfig class is used to store configuration settings for a Mixture object.

The MixtureConfig object can be initialized as follows:

config = MixtureConfig('mintolDisplay', 1e-6, ...
                       'compositionUnits', 'molar fraction', ...
                       'FLAG_COMPACT', true)

This creates an instance of the MixtureConfig class with the specified configuration.

See also: Mixture

FLAG_COMPACT = 'true'#

FLAG to print composition in a compact format

MixtureConfig(varargin)#

MixtureConfig contructor

This function accepts name-value pair arguments to override default property values.

Examples

  • config = MixtureConfig();

  • config = MixtureConfig(‘mintolDisplay’, 1e-6);

  • config = MixtureConfig(‘compositionUnits’, ‘molar fraction’);

  • config = MixtureConfig(‘FLAG_COMPACT’, true);

  • config = MixtureConfig(‘mintolDisplay’, 1e-6, ‘compositionUnits’, ‘molar fraction’, ‘FLAG_COMPACT’, true);

compositionUnits = "'molar fraction'"#

Options: ‘mol’, ‘molar fraction’, or ‘mass fraction’

mintolDisplay = '1e-6'#

Minimum tolerance to display species

EquationState#

class EquationState#

Bases: handle

The EquationState() class is an abstract class that defines the interface for computing the pressure and molar volume of a mixture using a specified equation of state.

Subclasses must implement the following abstract methods:
  • getPressure(temperature, molarVolume, varargin)

  • getVolume(temperature, pressure, varargin)

See also: EquationStateIdealGas(), Mixture()

EquationStateIdealGas#

class EquationStateIdealGas#

Bases: combustiontoolbox.core.EquationState

The EquationStateIdealGas() class implements the ideal gas law.

The EquationStateIdealGas() object can be used to compute the pressure and molar volume of a mixture assuming ideal gas behavior.

Example

eos = EquationStateIdealGas();

See also: EquationState() and Mixture()

static getPressure(temperature, molarVolume, varargin)#

Compute pressure [Pa] using the ideal gas law.

Parameters:
  • temperature (float) – Temperature of the mixture [K]

  • molarVolume (float) – Molar volume of the mixture [m3/mol]

Returns:

pressure (float) – Pressure of the mixture [Pa]

Example

P = getPressure(obj, 300, 0.024)

static getVolume(temperature, pressure, varargin)#

Compute molar volume [m3/mol] using the ideal gas law.

Parameters:
  • temperature (float) – Temperature of the mixture [K]

  • pressure (float) – Pressure of the mixture [Pa]

Returns:

molarVolume (float) – Molar volume of the mixture [m3/mol]

Example

V = getVolume(obj, 300, 1e5)