Database classes#
Combustion Toolbox generates its own databases using an up-to-date version of NASA’s 9-coefficient polynomial fits [McBride, 2002]. This information is stored in the Database
abstract class. Three main databases are used in the code:
NasaDatabase
: contains data extracted from NASA’s database [McBride, 2002].BurcatDatabase
: extracts data from the Third Millennium database [Burcat and Ruscic, 2005].SolarAbundances
: stores the solar abundances of the elements [Asplund et al., 2009].
Database#
- class Database(varargin)#
Bases:
handle
The
Database
abstract class contains the common methods between database objects. This class is used as a base class for theNasaDatabase
andBurcatDatabase
classes.See also:
NasaDatabase()
,BurcatDatabase()
- Database(varargin)#
Constructor
- addSpecies(speciesName, DB_master)#
Add species to the database
- extrapolationMethod = None#
Extrapolation method for the griddedInterpolant objects
- filename = None#
Database filename
- getProperty(listSpecies, property)#
Gets the vector of the defined property for the given set of species
- Parameters:
obj (
Database
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fitslistSpecies (
cell
) – List of speciesproperty (
char
) – Property name
- Returns:
value (float) – Property vector
Example
value = getProperty(obj, {‘H2O’, ‘CO2’}, ‘hf’)
- interpolationMethod = None#
Interpolation method for the griddedInterpolant objects
- load(varargin)#
Load database from file
- Parameters:
obj (
Database
) – Database object
- Optional Args:
filename (char): Filename of the database
- Returns:
obj (Database) – Database object
- name = None#
Database name
- pointsTemperature = None#
Number of points to use in the interpolation of the thermodynamic data
- species = None#
Struct with Species objects
- temperatureReference = None#
Default temperature of reference
- thermoFile = None#
Name with path for the thermodynamic database file (raw data)
- units = None#
Units thermodynamic data [molar or mass]
NasaDatabase#
- class NasaDatabase(varargin)#
Bases:
combustiontoolbox.databases.Database
,handle
The
NasaDatabase()
class is used to store thermodynamic data from NASA’s database using NASA’s 9 coefficient polynomial fits.The
NasaDatabase()
object can be initialized as follows:database = NasaDatabase()
This creates an instance of the
NasaDatabase()
class and initializes it with the chemical species contained in NASA’s database.See also:
BurcatDatabase()
,Database()
- NasaDatabase(varargin)#
Constructor
- Optional Args:
varargin (optional): key-value pairs to initialize the database
- Returns:
obj (NasaDatabase) – Object with NASA’s database
Examples
db = combustiontoolbox.databases.NasaDatabase();
db = combustiontoolbox.databases.NasaDatabase(‘filename’, ‘DB.mat’);
- static fullname2name(species)#
Get full name of the given species
- Parameters:
species (
char
) – Chemical species- Returns:
name (char) – Full name of the given species
- generateDatabase(varargin)#
Generate Database with thermochemical interpolation curves from the data extracted from the thermoFile
- Parameters:
obj (
NasaDatabase
) – NasaDatabase object
- Optional Args:
listSpecies (cell): List of species to generate the database
- Returns:
obj (NasaDatabase) – NasaDatabase object with thermochemical interpolation curves
Examples
DB = generateDatabase(NasaDatabase());
DB = generateDatabase(NasaDatabase(), {‘N2’, ‘O2’, ‘NO’, ‘O’, ‘N’});
- generateDatabaseMaster()#
Generate Master Database (DB_master) with the thermodynamic data of the chemical species
- Parameters:
obj (
NasaDatabase
) – NasaDatabase object- Returns:
DB_master (struct) – Database with the thermodynamic data of the chemical species
Example
DB_master = generateDatabaseMaster(‘thermo_CT.inp’)
- static getChangeMolesGasReaction(elements, elementMatrix, phase)#
In order to compute the internal energy of formation from the enthalpy of formation of a given species, we must determine the change in moles of gases during the formation reaction of a mole of that species starting from the elements in their reference state.
Notes
The only elements that are stable as diatomic gases are elements 1 (H), 8 (N), 9 (O), 10 (F), and 18 (Cl). The remaining elements that are stable as (monoatomic) gases are the noble gases He (3), Ne (11), Ar (19), Kr (37), Xe (55), and Rn (87), which do not form any compound.
- Parameters:
elements (
Elements
) – Elements objectelementMatrix (
float
) – Element matrix of the speciesphase (
float
) – 0 or 1 indicating gas or condensed species
- Optional Args:
elements (Elements): Elements object
- Returns:
Delta_n (float) – Change in moles of gases during the formation reaction of a mole of that species starting from the elements in their reference state
Example
Delta_n = getChangeMolesGasReaction(elements, element_matrix, phase)
- static getCoefficients(species, DB)#
Unpack NASA’s polynomials coefficients from database
- Parameters:
species (
char
) – Chemical speciesDB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
Tuple containing
a (cell): Temperature coefficients
b (cell): Integration constants
Trange (cell): Ranges of temperatures [K]
Texponents (cell): Exponent coefficients
Tintervals (float): Number of intervals of temperatures
phase (float): 0 or 1 indicating gas or condensed phase, respectively
hf0 (float): Enthalpy of formation [J/mol]
W (float): Molecular weight [kg/mol]
FLAG_REFERENCE (bool): Flag indicating species is a reference element/species
Example
[a, b, Trange, Texponents, Tintervals, phase, hf0, W, FLAG_REFERENCE] = getCoefficients(‘H2O’, DB)
- getDatabaseMaster(thermoFile)#
Generate Master Database (DB_master) with the thermodynamic data of the chemical species
- Parameters:
thermoFile (
char
) – path of thermoFile- Returns:
DB_master (struct) – Database with the thermodynamic data of the chemical species
- static getIndexTempereratureInterval(species, T, DB)#
Get interval of the NASA’s polynomials from the Database (DB) for the given species and temperature [K].
- Parameters:
species (
char
) – Chemical speciesT (
float
) – Temperature [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
Tinterval (float) – Index of the interval of temperatures
- getSpeciesThermo(DB, species, temperature, units)#
Calculates the thermodynamic properties of any species included in the NASA database
- Parameters:
DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fitsspecies (
char
) – Chemical speciesT (
float
) – Temperature [K]units (
char
) – Label indicating mass [kg] or molar [mol] units
- Returns:
Tuple containing
cP0 (float): Specific heat at constant pressure [J/(mol-k)]
hf0 (float): Enthalpy of formation [J/mol]
h0 (float): Enthalpy [J/mol]
ef0 (float): Internal energy of formation [J/mol]
s0 (float): Entropy [J/(mol-k)]
Dg0 (float): Gibbs energy [J/mol]
Example
[cp0, hf0, h0, ef0, s0, g0] = getSpeciesThermo(obj, DB, ‘CO’, 1000, ‘molar’)
cp0 = 33.1788
hf0 = -1.1054e+05
h0 = -8.8848e+04
ef0 = -1.1177e+05
s0 = 234.5409
- getSpeciesThermoFull(DB, species, temperature, units)#
Compute thermodynamic function using NASA’s 9 polynomials
- Parameters:
species (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fitsunits (
char
) – Label indicating mass [kg] or molar [mol] units
- Returns:
Tuple containing
cp (float): Specific heat at constant pressure in units basis [J/(units-K)]
cv (float): Specific heat at constant volume in units basis [J/(units-K)]
h0 (float): Enthalpy in units basis [J/units]
DhT (float): Thermal enthalpy in units basis [J/units]
e0 (float): Internal energy in units basis [J/units]
DeT (float): Thermal internal energy in units basis [J/units]
s0 (float): Entropy in units basis [J/(units-K)]
g0 (float): Gibbs energy in units basis [J/units]
Example
[cp, cv, h0, DhT, e0, DeT, s0, g0] = speciesThermoFull(‘H2O’, 300:100:6000, DB)
- species_DeT_NASA(species, temperature, DB)#
Compute thermal internal energy [J/mol] of the species at the given temperature [K] using NASA’s 9 polynomials
- Parameters:
species (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
DeT (float) – Thermal internal energy in molar basis [J/mol]
Example
DeT = species_DeT_NASA(‘H2O’, 300:100:6000, DB)
- species_DhT_NASA(species, temperature, DB)#
Compute thermal enthalpy [J/mol] of the species at the given temperature [K] using NASA’s 9 polynomials
- Parameters:
species (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
DhT (float) – Thermal enthalpy in molar basis [J/mol]
Example
DhT = species_DhT_NASA(‘H2O’, 300:100:6000, DB)
- species_cP_NASA(obj, species, temperature)#
Compute specific heats at constant pressure and at constant volume [J/(mol-K)] of the species at the given temperature [K] using NASA’s 9 polynomials
- Parameters:
obj (
NasaDatabase
) – NasaDatabase objectspecies (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]
- Returns:
Tuple containing
cp (float): Specific heat at constant pressure in molar basis [J/(mol-K)]
cv (float): Specific heat at constant volume in molar basis [J/(mol-K)]
Example
[cp, cv] = species_cP_NASA(‘H2O’, 300:100:6000, DB)
- species_cV_NASA(species, temperature, DB)#
Compute specific heat at constant volume [J/(mol-K)] of the species at the given temperature [K] using NASA’s 9 polynomials
- Parameters:
species (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
cV (float) – Specific heat at constant volume in molar basis [J/(mol-K)]
Example
cV = species_cV_NASA(‘H2O’, 300:100:6000, DB)
- species_e0_NASA(species, temperature, DB)#
Compute internal energy and the thermal internal energy [J/mol] of the species at the given temperature [K] using NASA’s 9 polynomials
- Parameters:
species (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
Tuple containing
e0 (float): Internal energy in molar basis [J/mol]
DeT (float): Thermal internal energy in molar basis [J/mol]
Example
[e0, DeT] = species_e0_NASA(‘H2O’, 300:100:6000, DB)
- species_g0_NASA(obj, species, temperature)#
Compute Compute Gibbs energy [J/mol] of the species at the given temperature [K] using NASA’s 9 polynomials
- Parameters:
obj (
NasaDatabase
) – NasaDatabase objectspecies (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]
- Returns:
g0 (float) – Gibbs energy in molar basis [J/mol]
Example
g0 = species_g0_NASA(NasaDatabase(), ‘H2O’, 300:100:6000)
- species_gamma_NASA(species, T, DB)#
Compute adiabatic index of the species [-] at the given temperature [K] using piecewise cubic Hermite interpolating polynomials and linear extrapolation
- Parameters:
species (
char
) – Chemical speciesT (
float
) – Range of temperatures to evaluate [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
gamma (float) – Adiabatic index [-]
Example
gamma = species_gamma_NASA(‘H2O’, 300:100:6000, DB)
- species_h0_NASA(species, temperature, DB)#
Compute enthalpy and thermal enthalpy [kJ/mol] of the species at the given temperature [K] using NASA’s 9 polynomials
- Parameters:
species (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
Tuple containing
h0 (float): Enthalpy in molar basis [J/mol]
DhT (float): Thermal enthalpy in molar basis [J/mol]
Example
[h0, DhT] = species_h0_NASA(‘H2O’, 300:100:6000, DB)
- species_s0_NASA(species, temperature, DB)#
Compute entropy [J/(mol-K)] of the species at the given temperature [K] using NASA’s 9 polynomials
- Parameters:
species (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
s0 (float) – Entropy in molar basis [J/(mol-K)]
Example
s0 = species_s0_NASA(‘H2O’, 300:100:6000, DB)
- species_thermo_NASA(species, temperature, DB)#
Compute thermodynamic function using NASA’s 9 polynomials
- Parameters:
species (
char
) – Chemical speciestemperature (
float
) – Range of temperatures to evaluate [K]DB (
struct
) – Database with custom thermodynamic polynomials functions generated from NASAs 9 polynomials fits
- Returns:
Tuple containing
cP (float): Specific heat at constant pressure in molar basis [J/(mol-K)]
cV (float): Specific heat at constant volume in molar basis [J/(mol-K)]
h0 (float): Enthalpy in molar basis [J/mol]
DhT (float): Thermal enthalpy in molar basis [J/mol]
e0 (float): Internal energy in molar basis [J/mol]
DeT (float): Thermal internal energy in molar basis [J/mol]
s0 (float): Entropy in molar basis [J/(mol-K)]
g0 (float): Gibbs energy in molar basis [J/mol]
Example
[cP, cV, h0, DhT, e0, DeT, s0, g0] = species_thermo_NASA(‘H2O’, 300:100:6000, DB)
BurcatDatabase#
- class BurcatDatabase(varargin)#
Bases:
combustiontoolbox.databases.Database
,handle
The
BurcatDatabase()
class is used to store thermodynamic data from Burcat’s database using NASA’s 9 coefficient polynomial fits.The
BurcatDatabase()
object can be initialized as follows:database = BurcatDatabase()
This creates an instance of the
BurcatDatabase()
class and initializes it with the chemical species contained in Burcat’s database.See also:
NasaDatabase()
,Database()
- BurcatDatabase(varargin)#
Constructor
- static thermoMillennium_2_thermoNASA9(filenameInput, varargin)#
Read Extended Third Millennium Thermodynamic Database of New NASA Polynomials with Active Thermochemical Tables update and write a new file compatible with thermo NASA 9 format
- Parameters:
filenameInput (
char
) – Filename of the thermoMillennium datafilenameOutput (
char
) – Filename of the thermoMillennium data as NASA9
SolarAbundances#
- class SolarAbundances(varargin)#
The
SolarAbundances
class is used to read solar abundances and compute the initial molar composition of the mixtureSee also:
abundances2moles()
- SolarAbundances(varargin)#
Class constructor
- abundances2moles(elements, varargin)#
Read solar abundances in log 10 scale and compute the initial molar fractions in the mixture [-]
- Parameters:
obj (
SolarAbundances
) – SolarAbundances objectelements (
cell
) – List with the given elementsfilename (
file
) – Filename with the data
- Optional Args:
metallicity (float): Metallicity
- Returns:
moles (float) – moles relative to H of the remaining elements in the mixture
Examples
moles = SolarAbundances.abundances2moles({‘H’, ‘He’, ‘C’, ‘N’, ‘O’, ‘Ne’, ‘Ar’, ‘S’, ‘Cl’, ‘Fe’})
moles = SolarAbundances.abundances2moles({‘H’, ‘He’, ‘C’, ‘N’, ‘O’, ‘Ne’, ‘Ar’, ‘S’, ‘Cl’, ‘Fe’}, 10)
- static read_abundances(filename)#
Read solar abundances file
Format: [number element, element, abundance, name, molar mass (g/mol)]
- Parameters:
filename (
file
) – Filename with the data- Returns:
Tuple containing
abundances (float): Vector with the logarithmic base 10 solar abundances
elements (cell): List with the given elements