The Flamelet Table Tool

The Flamelet Table Tool has evolved to have more options that previous versions, and so this section serves to document the various options that can be enabled via the input control file to the Flamelet Table Tool.

For compressible tables (FPVC, FPIC), the critical elements are the following:

  1. A Cantera formatted yaml mechanism file for the mechanism that was used to generate the flamelet solutions.

  2. A simple text file that maps the names of the species in the Cantera mechanism file to those in the flamelet solution files if there are any. Only the species that have different names in the flamelet solution files need to be included in this file.

  3. A version of Cantera installed in the local Python environment (A Python virtual environment is the easiest way to set this up).

Python Virtual Environment with Cantera

Cantera can be compiled from scratch and installed to a virtual environment if a user has a custom version of Cantera. Cantera also maintains Python pip installation option where a user can install Cantera by simply using:

pip install cantera

A tutorial for installing from the source code is here

To create a python virtual environment, do the following:

  1. python3 -m venv cantera-env

  2. Enter environment: source cantera-env/bin/activate

  3. Leave environment: deactivate

The build method requires a separate source and build directory. Place the downloaded Cantera source into a directory such as cantera/``and create an empty directory called ``cantera_build.

Enter the Python Cantera environment before compiling. Use the following Scons options when configuring the Cantera installation:

Scons options for building Cantera. Note put all of these options on one line.
scons build prefix = '<path>/cantera_build/' python_package = 'full'
python_prefix = '<path_to_python_venv>/lib/python3.4/site-packages'
python3_package = 'full'
python3_prefix = '<path_to_python_venv>/'

Then run scons install, which should install the Cantera package into the Python virtual environment.

Compressible Flamelet Tables

The compressible flamelet tables require a few additional inputs from non-compressible tables. Namely the following options:

  1. MECHANISM - Name of the Cantera formatted mechanism

  2. SPECIESNAMEREMAPFILE - Map between flamelet species names & Cantera mechanism species names (only ones that are different need to be included)

  3. RETAIN_CANTERA_NAMES - If TRUE, final species variable names in the outputs use the Cantera names (via the mapping file). If not set, species keep the flamelet names (except optional MDL renaming).

  4. EOS - The type of equation of state to use for encoding EOS data into the compressible flamelet table

The file contents below show a sample input file for generating a compressible diffusion flamelet table.

Sample input control file for generating a compressible diffusion flamelet table using the Flamelet Table Tool.
# Author name to appear in flamelet table
AUTHOR Christopher Robbin

# Path to the flamelet calculations
PREFIX <path_to_flamelets>

# Describes the folders and name of files of the flamelets. Multiple folders
# can be added
FLAMELETPATHS H2*

# Type of flamelet table (FPV, FPVT, FPVC, FPIC)
TABLETYPE FPVC
CLOSURETYPE ThickenedFlame
MECHANISM burke.cti

#Map between flamelet species name & Cantera mechanism species names
SPECIESNAMEREMAPFILE cantera_burke_species_names.txt
RETAIN_CANTERA_NAMES TRUE
EOS ideal #ideal or pengrobinson

# Type of flamelet file to use for creating table (FlameMaster, Cantera)
FLAMELETTYPE Flamemaster

# Characteristics of the flamelet table
NZMEAN 401
NCMEAN 401

ZST 0.22858
ZSPACING zst
CSPACING homogeneous

# Desired outputs
OUTPUTNAME burke_fpvc_ideal
OUTPUTVARIABLES W_OH W_H2O2
OUTPUTTYPE hdf5 vtk
DEFINEPROGVAR Y_H2O
TABLE_VERIFICATION TRUE

Progress Variable Definition and Optimization

Most flamelet tables generated by this tool are parameterized by mixture fraction, Z, and a reaction progress variable, C. The progress variable is defined from the species listed in DEFINEPROGVAR. For the usual unity-weight or manually weighted definitions, the tool uses the direct weighted sum

\[C = \sum_i w_i Y_i\]

where \(Y_i\) is the species mass fraction and \(w_i\) is the corresponding progress-variable weight. The source term stored in the table as W_C is the same weighted sum applied to the species source terms,

\[W_C = \sum_i w_i W_i\]

where \(W_i\) is the source term for species \(i\).

Optimized Progress Variables: Math

The PROGVARWEIGHTMETHOD optimal option is based on the progress-variable regularization idea described by Ihme, Shunn, and Zhang, “Regularization of reaction progress variable for application to flamelet-based combustion models,” Journal of Computational Physics, 2012, DOI: 10.1016/j.jcp.2012.06.029.

The goal is to choose weights that reduce overlap between neighboring flamelets when the manifold is viewed in (Z, C) space. At each optimizer evaluation, the tool computes a candidate raw progress variable,

\[C_\mathrm{raw} = \sum_i w_i Y_i\]

and normalizes it onto the table coordinate interval,

\[C = \left(C_\mathrm{raw} - C_\mathrm{min}\right)\frac{1}{C_\mathrm{max} - C_\mathrm{min}}\]

The normalization is stored in the log as offset and scale:

\[\mathrm{offset} = C_\mathrm{min}, \qquad \mathrm{scale} = \frac{1}{C_\mathrm{max} - C_\mathrm{min}}\]

This affine normalization only affects optimized weights. Unity and manual weights use offset=0 and scale=1. For an optimized progress variable, the table uses

\[C = \left(\sum_i w_i Y_i - \mathrm{offset}\right)\mathrm{scale}\]

and

\[W_C = \mathrm{scale}\sum_i w_i W_i\]

The offset does not appear in W_C because it is constant.

The optimization loss penalizes regions where the progress variable is not monotone with respect to the stoichiometric progress value, \(C_\mathrm{st}\), across the flamelet family. In implementation form, the monotonicity part of the loss is

\[J_\mathrm{mono} = \int_{C_\mathrm{st}}\int_Z C(Z, C_\mathrm{st})\, H\left(\epsilon - \frac{\partial C}{\partial C_\mathrm{st}}\right)\,dZ\,dC_\mathrm{st}\]

where \(H\) is a Heaviside function and \(\epsilon\) is a small slope tolerance. A smaller value indicates less overlap/non-monotonicity in the candidate progress coordinate.

The tool uses SciPy’s differential-evolution optimizer with each weight bounded in [-1, 1]. Differential evolution is a stochastic population method, so individual loss evaluations are not expected to decrease monotonically. Watch the best loss in the log rather than the current loss.

Progress-variable optimization is mathematical rather than chemical. It can find a coordinate that separates flamelets well but is dominated by an intermediate species. Such a result can be useful in some manifolds, but it changes the meaning of the transported progress variable and its source term. Product-heavy definitions such as Y_H2O, Y_CO2, Y_CO, and sometimes Y_H2 usually remain easier to interpret as a burnedness coordinate. If the optimizer reports that one species dominates the weight norm, inspect the resulting PROG and SRC_PROG fields before using the table in production.

Using Progress-variable Weights in the Control File

If no weighting method is specified, the tool defaults to unity weights:

DEFINEPROGVAR Y_H2O Y_H2 Y_CO Y_CO2
PROGVARWEIGHTMETHOD unity

This is equivalent to omitting PROGVARWEIGHTMETHOD. The resulting progress variable is Y_H2O + Y_H2 + Y_CO + Y_CO2.

Manual weights are supplied in the same order as DEFINEPROGVAR:

DEFINEPROGVAR Y_H2O Y_H2 Y_CO Y_CO2
PROGVARWEIGHTMETHOD manual
PROGVARWEIGHTS 1.0 0.5 1.0 1.0

Optimized weights are enabled with:

DEFINEPROGVAR Y_H2O Y_H2 Y_CO Y_CO2
PROGVARWEIGHTMETHOD optimal

The optimized run logs the final species weights, the optimized loss compared with the unity-weight loss, the normalization offset and scale, and a warning if the result is dominated by one species. A typical optimized trace contains entries like:

optimization loss evaluation 250: current=7.164e-03 best=1.466e-03 monotonicity=7.164e-03
Optimized progress variable weights: Y_H2O=-0.025, Y_H2=0.796, Y_CO=0.027, Y_CO2=0.017
Optimized progress variable normalization: offset=-0.000277 scale=20.7171

The optimizer also supports an optional reference-weight regularization term:

PROGVAROPTIMIZER_WEIGHT_REGULARIZATION 0.03

This term biases the solution toward the unity-weight direction while still minimizing flamelet overlap:

\[J = J_\mathrm{mono} + \lambda\left(1 - \frac{\mathbf{w}\cdot\mathbf{1}}{\|\mathbf{w}\|\|\mathbf{1}\|}\right)^2\]

where \(\lambda\) is PROGVAROPTIMIZER_WEIGHT_REGULARIZATION. The default is 0.0, so the term is inactive unless explicitly requested. This option is useful when the unconstrained optimizer finds a low-overlap but chemically awkward definition, for example one dominated by a single intermediate species.

Progress Source-term Scaling

Compressible tables can optionally carry two extra coefficients for correcting the tabulated progress-variable source term outside the baseline flamelet state:

\[\dot{\omega}_C = \dot{\omega}_{C,0} \left(\frac{\rho}{\rho_0}\right)^\alpha \exp\left[-T_a\left(\frac{1}{T} - \frac{1}{T_0}\right)\right]\]

The table stores \(\alpha\) as SRC_PROG_ALPHA and \(T_a\) as SRC_PROG_TA. The baseline source term is still written as SRC_PROG. The CFD solver applies the correction using its local rho and T together with the table values rho0 and T0.

The table tool treats these coefficients as a generic source-scaling surface. ZND fitting for detonation work is one way to generate the surface, but the same interface can be used for coefficients fitted from pressure sweeps, constant-volume reactors, or other calibration data.

Generating ZND Profiles with SDToolbox

The vendored SDToolbox copy can be installed as an editable Python package from the repository root:

python -m pip install -e third-party/SDToolbox

It is Python-only and uses Cantera YAML mechanisms. Mechanisms bundled with SDToolbox can be loaded by filename, for example Li2015.yaml or Mevel2017.yaml. A minimal ZND profile generation script is:

make_znd_profile.py
import cantera as ct

from sdtoolbox.mechanism import Solution
from sdtoolbox.postshock import CJspeed, PostShock_fr
from sdtoolbox.utilities import znd_profile_to_csv
from sdtoolbox.znd import zndsolve


pressure = 100000.0
temperature = 300.0
composition = "H2:2 O2:1 N2:3.76"
mechanism = "Li2015.yaml"

cj_speed = CJspeed(pressure, temperature, composition, mechanism)

gas1 = Solution(mechanism)
gas1.TPX = temperature, pressure, composition

gas = PostShock_fr(cj_speed, pressure, temperature, composition, mechanism)
znd_output = zndsolve(
    gas,
    gas1,
    cj_speed,
    t_end=1.0e-6,
    max_step=1.0e-8,
    advanced_output=True,
)

znd_profile_to_csv("znd_profile.csv", znd_output)

The exported CSV contains distance, time, pressure, temperature, density, velocity, thermicity, heat-release rate, specific heats, enthalpy, internal energy, species mass fractions (Y_*), species mole fractions (X_*), and net production rates (WDOT_*).

Wiring ZND Fits into a Flamelet Table

Raw ZND profiles are not passed directly to the flamelet table generator. They are calibration data. A preprocessing step should reduce one or more ZND calculations into the coefficient surface used by the CFD-side correction:

\[\log\left(\frac{\dot{\omega}_C}{\dot{\omega}_{C,0}}\right) = \alpha \log\left(\frac{\rho}{\rho_0}\right) - T_a\left(\frac{1}{T} - \frac{1}{T_0}\right)\]

For each table coordinate where you want a correction, collect a baseline source term \(\dot{\omega}_{C,0}\) and one or more perturbed source terms \(\dot{\omega}_C\), then fit SRC_PROG_ALPHA and SRC_PROG_TA with a linear least-squares solve. The current table-generator interface expects the reduced coefficient file, not the full calibration profile:

Z,C,SRC_PROG_ALPHA,SRC_PROG_TA
0.02,0.10,0.70,11500.0
0.02,0.50,0.82,12800.0
0.02,0.90,0.76,12100.0
0.10,0.10,0.68,10900.0
0.10,0.50,0.79,12400.0
0.10,0.90,0.73,11800.0

Use coordinate columns that the table tool can evaluate from each flamelet. Common choices are Z C for FPVC-like tables, Phi C for FPIC-like premixed tables, or P C for pressure-dependent calibration. The coordinate names are not hard-coded beyond the built-in aliases; if a column matches a flamelet variable or property name, it can be used.

For a constant scaling coefficient, use:

PROG_SOURCE_SCALING constant
PROG_SOURCE_SCALING_ALPHA 0.7
PROG_SOURCE_SCALING_TA 12000.0

For a coefficient surface, provide a CSV file with named coordinate columns and the two coefficient columns:

Z,C,SRC_PROG_ALPHA,SRC_PROG_TA
0.0,0.0,0.0,1000.0
1.0,0.0,0.8,1100.0
0.0,1.0,0.1,1010.0
1.0,1.0,0.9,1110.0

and enable it with:

PROG_SOURCE_SCALING surface
PROG_SOURCE_SCALING_FILE source_scaling.csv
PROG_SOURCE_SCALING_COORDINATES Z C

PROG_SOURCE_SCALING znd is accepted as an alias for a surface generated from ZND preprocessing. The table code does not otherwise care where the coefficients came from.

ZND-native Tables

The longer-term, more detonation-aligned path is to use ZND profiles as the table manifold itself rather than using premixed flamelets with an external source-term correction. This is exposed as a first implementation slice through:

TABLETYPE FPZND
CLOSURETYPE ThickenedFlame
FLAMELETTYPE ZND

In this mode each input CSV file is one ZND trajectory at a fixed inlet mixture state. The profile direction is mapped to the progress-variable coordinate C exactly like a premixed flamelet, while the family coordinate is stored as z_inlet. The table is therefore a two-dimensional (Z, C) manifold built directly from ZND profiles.

Touched/new classes for the first ZND-native path are:

sdtoolbox.utilities.znd_profile_to_csv

Exports SDToolbox ZND solutions to a single CSV and now accepts optional comment metadata such as z_inlet or equivalence_ratio.

flamelet_tool.flamelet.ZNDProfileFlamelet

A ZND profile flamelet object. It reuses the premixed C-grid behavior because the table marches along the profile by progress variable.

flamelet_tool.flamelet_reader.ZNDProfileFlameletReader

Reads SDToolbox-style ZND CSV files, converts Y_* columns to species mass fractions, converts WDOT_* columns from kmol/m^3/s to PRODRATE_* in kg/m^3/s, and fills transport properties from Cantera.

flamelet_tool.flamelet_table.ZNDTableCoordinateCalculator

Builds the table family-coordinate grid over the available ZND z_inlet range by default. Use ZSPACING full if a 0..1 grid is desired.

flamelet_tool.thickened_table.FPZNDThickenedFlameFlameletTable

Builds the ZND-native thickened table. It intentionally skips FPIC’s synthetic pure-oxidizer/pure-fuel flamelets so the manifold is only what the ZND data supports.

flamelet_tool.table_writer.FPZNDThickenedFlameTableWriter

Writes the table with Combustion Model = FPZND and Table Type = ZND Native Thickened Flame.

flamelet_tool.main.Main and FlameletTableFactory

Dispatch TABLETYPE FPZND and FLAMELETTYPE ZND to the reader/table classes above.

The reader needs each profile to define its family coordinate. The simplest option is to write metadata comments before the CSV header:

# z_inlet: 0.1111111111
# equivalence_ratio: 1.0
distance,time,pressure,temperature,density,...

Equivalently, provide a manifest:

path,z_inlet,equivalence_ratio
znd_phi_080.csv,0.0909,0.8
znd_phi_100.csv,0.1111,1.0
znd_phi_120.csv,0.1304,1.2

and set:

ZND_PROFILE_MANIFEST znd_manifest.csv

A minimal FPZND table input looks like:

AUTHOR znd native example
PREFIX znd_profiles
FLAMELETPATHS znd_phi_*.csv

TABLETYPE FPZND
CLOSURETYPE ThickenedFlame
FLAMELETTYPE ZND
MECHANISM H2.yaml
EOS ideal

ZST 0.1111111111111111
FUEL H2
NZMEAN 21
NCMEAN 101

OUTPUTNAME fpznd_h2_air
OUTPUTTYPE hdf5 vtk
OUTPUTVARIABLES Z HEATRELEASE W_C

DEFINEPROGVAR Y_H2O
PROGVARWEIGHTMETHOD unity

The repository also includes a methane/oxygen playground that generates the ZND profiles and then builds the FPZND table:

cd tests/integration_tests/fpznd_thickened_test/generate_znd_profiles
python generate_znd_profiles.py

cd ../table
flamelet_table_tool input_file

That case uses gri30.yaml, initial methane/oxygen mixtures at 1 MPa, and a generated ZND_PROFILE_MANIFEST so the profile metadata is easy to inspect.

Current limitations are deliberate: this first path is thickened-flame/HDF5/VTK focused, it expects pre-generated ZND CSV profiles, and it does not yet include a high-level batch ZND generator or beta-closure FPZND table. Those are the next natural pieces once we settle the profile family coordinates and the progress-variable definition for the RDE cases.

Useful Concepts in the Flamelet Table Tool

The following are some helpful general concepts and hints to aid users in using the Flamelet Table Tool.

  1. For all flamelet table types, the value of the OUTPUTVARIABLES option does not need to be specified as all the required output variables that are needed for that flamelet table type are automatically included in the output. However, if the user wants to include additional variables from the flamelet solutions, they can be specified in the OUTPUTVARIABLES option.

  2. The TABLE_VERIFICATION option is used to trigger a few sanity checks on the data within the flamelet table to ensure that the data is consistent and that the table is usable. This option is useful for debugging purposes. It’s value can be TRUE or FALSE.

  3. The GENERATE_DIAGNOSTICS option can be used to generate additional diagnostic information about the flamelet table generation process. This option can generate a lot of output, so it is recommended to use this option only when debugging the flamelet table generation process. The outputs that are generated however are very useful in visualizing the set of flamelets that have been provided to the flamelet table tool and seeing how well the table approximates the data from the flamelets.

  4. There are some cases where the GENERATE_DIAGNOSTICS option generates a great deal of data and slows down the table generation process. If a user still wants to get a subset of the diagnostic information, they can use the DIAGNOSTIC_OUTPUTS option to specify the type of diagnostic information that they want to generate. Valid values for this variable are special and variables. The special option generates various plots that give a user an idea of the state of the flamelet table. The variables option generates plots of the flamelet data that is used to generate the table, and is useful for getting a feel for what data is being provided to the flamelet table tool.

Inert Flamelet Type (Synthetic Boundary Mixing)

The FLAMELETTYPE inert option constructs a synthetic non-reacting diffusion flamelet directly from user-specified boundary states at Z=0 and Z=1. This workflow does not read flamelet files from disk, so PREFIX and FLAMELETPATHS are not required.

This workflow is provided as a reproducible example in this repository at tests/integration_tests/inert_fpvc_pure_mixing_table.

Use the following options to generate an inert FPVC pure-mixing table:

Example control file (inert_fpvc_table_input.dat) for an inert FPVC pure-mixing table.
AUTHOR Wanda Daradar
TABLETYPE FPVC_PURE_MIXING
CLOSURETYPE ThickenedFlame
MECHANISM he_n2.yaml
EOS ideal
FLAMELETTYPE inert

# Required by current C_st lookup path
ZST 0.5
FUEL N2

# Boundary states used to synthesize the inert flamelet
BOUNDARY_MIXING_P 6.0e6
BOUNDARY_MIXING_Z0_T 900.0
BOUNDARY_MIXING_Z1_T 363.0
BOUNDARY_MIXING_Z0_Y HE:1.0
BOUNDARY_MIXING_Z1_Y N2:1.0
BOUNDARY_MIXING_NPOINTS 101

NZMEAN 101
NCMEAN 101
DEFINEPROGVAR Y_N2
OUTPUTNAME he_n2_inert_fpvc_pure_mixing
OUTPUTTYPE hdf5 vtk
TABLE_VERIFICATION TRUE
GENERATE_DIAGNOSTICS TRUE
LOGGING_LEVEL INFO

Run the case with:

cd tests/integration_tests/inert_fpvc_pure_mixing_table
flamelet_table_tool inert_fpvc_table_input.dat

Notes:

  1. Z=0 is treated as the oxidizer-side boundary and Z=1 as the fuel-side boundary in the generated flamelet metadata.

  2. The synthetic inert flamelet has no chemistry; species production rates and heat release are set to zero.

  3. For a 3-D inert table with dimensions (Z, C, Y_inert), use TABLETYPE FPVC_INERT with inert-level flamelet groups (see the inert section in the FlameMaster guide).