Ignis
Loading...
Searching...
No Matches
ignis Class Reference

#include <ignis.h>

Collaboration diagram for ignis:
[legend]

Public Member Functions

void setIC (const std::string icType, const std::string fname="")
 
void storeState ()
 
void setFluxesUnity ()
 
void setFluxes ()
 
void setGrid (double _L)
 
void writeFileHdf5 (const std::string gname, const std::string timeType)
 
void writeFile (const std::string fname)
 
void solveSS ()
 
void setChi (const double _chi0)
 
void solveUnsteady (const double nTauRun, const int nsteps, const bool doWriteTime=true, const double Tmin=0, const double Tmax=0)
 
int Func (const double *vars, double *F)
 
int rhsf (const double *vars, double *dvarsdt)
 
int rhsf_flamelet (const double *vars, double *dvarsdt)
 
void setQrad (std::vector< double > &Q)
 
void setTprof (const std::vector< double > &_Tprof_h, const std::vector< double > &_Tprof_T)
 
void setDerivative2 (const double vL, const double vR, const std::vector< double > &v, std::vector< double > &d2vdx2)
 
void setDerivative (const double vL, const double vR, const std::vector< double > &v, std::vector< double > &dvdx)
 
size_t I (size_t i, size_t k)
 
size_t Ia (size_t i, size_t k)
 
 ignis (const bool _isPremixed, const bool _doEnergyEqn, const bool _isFlamelet, const bool _doSoot, const size_t _ngrd, const double _L, const double _P, std::shared_ptr< Cantera::Solution > csol, std::string _radType, const std::vector< double > &_yLbc, const std::vector< double > &_yRbc, const double _TLbc, const double _TRbc, std::shared_ptr< soot::sootModel > _SM, std::shared_ptr< soot::state > _SMstate)
 

Public Attributes

size_t ngrd
 number of interior grid points
 
size_t nsp
 number of gas species
 
size_t nsoot
 number of soot variables
 
size_t nvar
 number of transported variables at each grid point
 
size_t nvarA
 nvar * ngrd
 
double P
 system pressure, uniform (Pa)
 
std::vector< std::vector< double > > y
 mass fractions: y[igrid][isp]
 
std::vector< double > T
 temperature (K)
 
std::vector< std::vector< double > > sootvars
 soot moments or sections; sootvars[igrid][isoot]
 
double Pstore
 stored system pressure (for initializing from stored state)
 
std::vector< std::vector< double > > ystore
 stored mass fractions
 
std::vector< double > Tstore
 stored temperature
 
std::vector< std::vector< double > > sootstore
 stored soot variables
 
std::vector< double > yLbc
 
std::vector< double > yRbc
 y boundary values: left and right (as needed)
 
double TLbc
 
double TRbc
 T boundary values: left and right (as needed)
 
double hLbc
 
double hRbc
 h boundary values: left and right (as needed)
 
double cpLbc
 
double cpRbc
 cp boundary values: left and right (as needed)
 
std::vector< double > hspLbc
 species enthalpies on left boundary
 
std::vector< double > hspRbc
 species enthalpies on right boundary
 
double L
 domain size (m)
 
std::vector< double > x
 grid position values (m)
 
std::vector< double > dx
 grid spacing (m), nonuniform is fine
 
std::vector< double > fl
 fractions for interpolation
 
std::vector< double > fr
 fractions for interpolation
 
double Tscale
 scaling value for temperature (for solvers)
 
double hscale
 scaling value for enthalpy (for solvers)
 
std::vector< double > sootScales
 scaling value for soot variables (for solvers)
 
std::vector< double > vars0
 for homotopy approaches
 
std::vector< double > F0
 for homotopy approaches
 
double s
 homotopy variable
 
std::shared_ptr< Cantera::ThermoPhase > gas
 Cantera thermo object.
 
std::shared_ptr< Cantera::Kinetics > kin
 Cantera kinetics object.
 
std::shared_ptr< Cantera::Transport > trn
 Cantera transport object.
 
std::shared_ptr< streamsstrm
 streams object (mixture fraction, etc.)
 
std::shared_ptr< rad > radProps
 radiation object
 
std::string radType
 radiation model name
 
bool doRadiation
 radiation flag
 
std::vector< double > kabs_sur
 kabs for surroundings
 
std::vector< double > awts_sur
 awts for surroundings
 
bool doLe1 = false
 true if doing unity Lewis numbers (default false)

 
double Ttarget
 for unsteady cases, run until this max T instead of for a given time
 
double dT
 delta T increment for unsteady cases
 
int isave
 file counter for save during unsteady cases
 
std::vector< std::vector< double > > flux_y
 species fluxes: [I(igrid, ksp)] I(igrid,ksp) maps 2D onto 1D
 
std::vector< std::vector< double > > flux_soot
 species fluxes: [I(igrid, ksoot)]
 
std::vector< double > flux_h
 species fluxes: [igrid]
 
bool isFlamelet = false
 true for laminar flamelet (mixture fraction coordinate)
 
std::vector< double > chi
 dissipation rate profile
 
double chi0
 
bool isPremixed = false
 true of the case is a premixed flame, (only left boundary condition, constant mass flux through domain)
 
double mflux = 0.0
 premixed flame mass flux (kg/m2*s)
 
bool doEnergyEqn = true
 for premixed flames: can solve energy equation or set T profile
 
std::shared_ptr< linearInterpLI
 interpolator for specified temperature profiles
 
std::vector< double > Tprof_h
 temperature profile position (h is height above burner (m))
 
std::vector< double > Tprof_T
 temperature profile T values
 
bool doSoot = false
 soot flag
 
std::shared_ptr< soot::sootModel > SM
 soot model
 
std::shared_ptr< soot::state > SMstate
 holds state variables (gas and soot) for soot model
 
std::shared_ptr< HighFive::File > fh5
 hdf5 file pointer
 

Detailed Description

One-dimensional flame solver: burner stabilized premixed and diffusion flames Premixed flames can use a fixed temperature profile, or solve the energy equation. Diffusion flames are diffusion only, with Dirichlet boundary conditions. Flame stain occurs through the domain length.

Definition at line 29 of file ignis.h.

Constructor & Destructor Documentation

◆ ignis()

ignis::ignis ( const bool  _isPremixed,
const bool  _doEnergyEqn,
const bool  _isFlamelet,
const bool  _doSoot,
const size_t  _ngrd,
const double  _L,
const double  _P,
std::shared_ptr< Cantera::Solution >  csol,
std::string  _radType,
const std::vector< double > &  _yLbc,
const std::vector< double > &  _yRbc,
const double  _TLbc,
const double  _TRbc,
std::shared_ptr< soot::sootModel >  _SM,
std::shared_ptr< soot::state >  _SMstate 
)

Constructor function

Parameters
_isPremixedinput: flag indicating if running a premixed or diffusion flame
_doEnergyEqninput: flag indicating if solving the energy equation or using a given T profile (for premixed flames)
_doSootinput: flag indicating if solving with soot
_ngrdinput: number of interior grid points
_Linput: domain size (m)
_Pinput: system pressure (Pa)
csolinput: Cantera Solution object (used to set gas, kin, trn)
_yLbcinput: mass fractions on left boundary
_yRbcinput: mass fractions on right boundary (not used for premixed flames)
_TLbcinput: temperature (K) on left boundary
_TRbcinput: temperature (K) on right boundary (not used for premixed flames)
_SMinput: soot model object
_SMstateinput: soot model state object (gas and soot properties for input to SM)

Definition at line 45 of file ignis.cc.

Here is the call graph for this function:

Member Function Documentation

◆ setIC()

void ignis::setIC ( const std::string  icType,
const std::string  fname = "" 
)

Set initial condition for solution (or initial guess for steady solvers).

Parameters
icTypeinput: string indicating the type of the initial condition to use
fnameinput: read input from this file name (defaults to empty string)

Definition at line 485 of file ignis.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ storeState()

void ignis::storeState ( )

Store system state so that we can use this state later as an initial condition. See setIC

Definition at line 469 of file ignis.cc.

Here is the caller graph for this function:

◆ setFluxesUnity()

void ignis::setFluxesUnity ( )

Compute fluxes for all transported variables assuming unity Lewis numbers for energy and species. Soot uses a thermophoretic flux. Only called if doLe1 is false (the default).

Definition at line 590 of file ignis.cc.

Here is the caller graph for this function:

◆ setFluxes()

void ignis::setFluxes ( )

Compute fluxes for all transported variables assuming mixture-averaged transport for species and energy. Soot uses a thermophoretic flux. Only called if doLe1 is false (the default).

Definition at line 738 of file ignis.cc.

Here is the caller graph for this function:

◆ setGrid()

void ignis::setGrid ( double  _L)

Set spatial grid (x and dx)

Parameters
_Linput: domain length (resets class data member L, so cases can be run with different L)

| * | * | * | * | * | dx dx dx dx dx

Definition at line 202 of file ignis.cc.

Here is the caller graph for this function:

◆ writeFileHdf5()

void ignis::writeFileHdf5 ( const std::string  gname,
const std::string  timeType 
)

There is a single hdf5 file, with multiple groups. Whenever you would normally write a file with x, mixf, T, h, rho, y, etc., instead, you make a new group inside the one hdf5 file.

Parameters
gnameinput: write output to this group name

Definition at line 264 of file ignis.cc.

Here is the caller graph for this function:

◆ writeFile()

void ignis::writeFile ( const std::string  fname)

Write output file

Parameters
fnameinput: write output to this file name (includes path)

Definition at line 345 of file ignis.cc.

Here is the caller graph for this function:

◆ solveSS()

void ignis::solveSS ( )

Solve steady state problem. Uses Sundials Kinsol. It is more robust to solve the unsteady problem to steady state. Solves F(vars) = 0, where vars are the vector of variables at all grid points and F is the equation for each of them. See ignis::Func.

Definition at line 942 of file ignis.cc.

Here is the call graph for this function:

◆ setChi()

void ignis::setChi ( const double  _chi0)

Definition at line 1607 of file ignis.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ solveUnsteady()

void ignis::solveUnsteady ( const double  nTauRun,
const int  nSteps,
const bool  doWriteTime = true,
const double  Tmin = 0,
const double  Tmax = 0 
)

Solve unsteady ignis problems. Assumes y, T are initialized Two modes: write on every time step of size dt, or write on temperature steps of size dT. Default is in time --> Tmin, Tmax are zero --> dT = 0

Parameters
nTauRuninput: number of characteristic times to solve for.
nStepsinput: number of steps to take during the solution.
doWriteTimeinput: if true (default) then write solution in time
Tmininput: minimum temperature (default is 0, see code) (for stepping to desired T)
Tmaxinput: maximum temperature (default is 0, see code) (for stepping to desired T)

Definition at line 1122 of file ignis.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Func()

int ignis::Func ( const double *  vars,
double *  F 
)

System of equations solved by the steady solver, F(vars) = 0

Parameters
varsinput: pointer to array of variables at every grid point.
Foutput: pointer to array of function values for each variable at each grid point 1D vectors are accessed by convenience function Ia(igrid, kvar), as in vars[Ia(i,k)].

Definition at line 1008 of file ignis.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ rhsf()

int ignis::rhsf ( const double *  vars,
double *  dvarsdt 
)

Right hand side function for unsteady solver, as in dvar/dt = rhsf(var)

Parameters
varsinput: current value of variables at each grid point
dvarsdtoutput: rate of each variable at each grid point

Definition at line 1215 of file ignis.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ rhsf_flamelet()

int ignis::rhsf_flamelet ( const double *  vars,
double *  dvarsdt 
)

Right hand side function for unsteady solver for flamelet equations, as in dvar/dt = rhsf(var)

Parameters
varsinput: current value of variables at each grid point
dvarsdtoutput: rate of each variable at each grid point

Definition at line 1317 of file ignis.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setQrad()

void ignis::setQrad ( std::vector< double > &  Q)

Set radiative source term for the energy equation.

Parameters
Qoutput: source term at each grid point (W/m3)

Definition at line 1082 of file ignis.cc.

Here is the caller graph for this function:

◆ setTprof()

void ignis::setTprof ( const std::vector< double > &  _Tprof_h,
const std::vector< double > &  _Tprof_T 
)
inline

Definition at line 130 of file ignis.h.

Here is the caller graph for this function:

◆ setDerivative2()

void ignis::setDerivative2 ( const double  vL,
const double  vR,
const std::vector< double > &  v,
std::vector< double > &  d2vdx2 
)

Compute second derivative of a profile at cell centers

Parameters
vLinput: variable on left face of domain
vRinput: variable on right face of domain
vinput: variable grid points
d2vdx2output: derivative on grid points

Get derivative on faces, then take difference between faces | * | * | * | * |

Definition at line 1568 of file ignis.cc.

Here is the caller graph for this function:

◆ setDerivative()

void ignis::setDerivative ( const double  vL,
const double  vR,
const std::vector< double > &  v,
std::vector< double > &  dvdx 
)

Compute derivative of a profile at cell centers

Parameters
vLinput: variable on left face of domain
vRinput: variable on right face of domain
vinput: variable grid points
dvdxoutput: derivative on grid points

Interpolate variable to faces, then take the derivative between faces | * | * | * | * |

Definition at line 1538 of file ignis.cc.

Here is the caller graph for this function:

◆ I()

size_t ignis::I ( size_t  i,
size_t  k 
)
inline

Definition at line 142 of file ignis.h.

◆ Ia()

size_t ignis::Ia ( size_t  i,
size_t  k 
)
inline

Definition at line 143 of file ignis.h.

Here is the caller graph for this function:

Member Data Documentation

◆ ngrd

size_t ignis::ngrd

number of interior grid points

Definition at line 33 of file ignis.h.

◆ nsp

size_t ignis::nsp

number of gas species

Definition at line 34 of file ignis.h.

◆ nsoot

size_t ignis::nsoot

number of soot variables

Definition at line 35 of file ignis.h.

◆ nvar

size_t ignis::nvar

number of transported variables at each grid point

Definition at line 36 of file ignis.h.

◆ nvarA

size_t ignis::nvarA

nvar * ngrd

Definition at line 37 of file ignis.h.

◆ P

double ignis::P

system pressure, uniform (Pa)

Definition at line 39 of file ignis.h.

◆ y

std::vector<std::vector<double> > ignis::y

mass fractions: y[igrid][isp]

Definition at line 40 of file ignis.h.

◆ T

std::vector<double> ignis::T

temperature (K)

Definition at line 41 of file ignis.h.

◆ sootvars

std::vector<std::vector<double> > ignis::sootvars

soot moments or sections; sootvars[igrid][isoot]

Definition at line 42 of file ignis.h.

◆ Pstore

double ignis::Pstore

stored system pressure (for initializing from stored state)

Definition at line 44 of file ignis.h.

◆ ystore

std::vector<std::vector<double> > ignis::ystore

stored mass fractions

Definition at line 45 of file ignis.h.

◆ Tstore

std::vector<double> ignis::Tstore

stored temperature

Definition at line 46 of file ignis.h.

◆ sootstore

std::vector<std::vector<double> > ignis::sootstore

stored soot variables

Definition at line 47 of file ignis.h.

◆ yLbc

std::vector<double> ignis::yLbc

Definition at line 49 of file ignis.h.

◆ yRbc

std::vector<double> ignis::yRbc

y boundary values: left and right (as needed)

Definition at line 49 of file ignis.h.

◆ TLbc

double ignis::TLbc

Definition at line 50 of file ignis.h.

◆ TRbc

double ignis::TRbc

T boundary values: left and right (as needed)

Definition at line 50 of file ignis.h.

◆ hLbc

double ignis::hLbc

Definition at line 51 of file ignis.h.

◆ hRbc

double ignis::hRbc

h boundary values: left and right (as needed)

Definition at line 51 of file ignis.h.

◆ cpLbc

double ignis::cpLbc

Definition at line 52 of file ignis.h.

◆ cpRbc

double ignis::cpRbc

cp boundary values: left and right (as needed)

Definition at line 52 of file ignis.h.

◆ hspLbc

std::vector<double> ignis::hspLbc

species enthalpies on left boundary

Definition at line 53 of file ignis.h.

◆ hspRbc

std::vector<double> ignis::hspRbc

species enthalpies on right boundary

Definition at line 54 of file ignis.h.

◆ L

double ignis::L

domain size (m)

Definition at line 56 of file ignis.h.

◆ x

std::vector<double> ignis::x

grid position values (m)

Definition at line 57 of file ignis.h.

◆ dx

std::vector<double> ignis::dx

grid spacing (m), nonuniform is fine

Definition at line 58 of file ignis.h.

◆ fl

std::vector<double> ignis::fl

fractions for interpolation

Definition at line 59 of file ignis.h.

◆ fr

std::vector<double> ignis::fr

fractions for interpolation

Definition at line 60 of file ignis.h.

◆ Tscale

double ignis::Tscale

scaling value for temperature (for solvers)

Definition at line 62 of file ignis.h.

◆ hscale

double ignis::hscale

scaling value for enthalpy (for solvers)

Definition at line 63 of file ignis.h.

◆ sootScales

std::vector<double> ignis::sootScales

scaling value for soot variables (for solvers)

Definition at line 64 of file ignis.h.

◆ vars0

std::vector<double> ignis::vars0

for homotopy approaches

Definition at line 66 of file ignis.h.

◆ F0

std::vector<double> ignis::F0

for homotopy approaches

Definition at line 67 of file ignis.h.

◆ s

double ignis::s

homotopy variable

Definition at line 68 of file ignis.h.

◆ gas

std::shared_ptr<Cantera::ThermoPhase> ignis::gas

Cantera thermo object.

Definition at line 70 of file ignis.h.

◆ kin

std::shared_ptr<Cantera::Kinetics> ignis::kin

Cantera kinetics object.

Definition at line 71 of file ignis.h.

◆ trn

std::shared_ptr<Cantera::Transport> ignis::trn

Cantera transport object.

Definition at line 72 of file ignis.h.

◆ strm

std::shared_ptr<streams> ignis::strm

streams object (mixture fraction, etc.)

Definition at line 74 of file ignis.h.

◆ radProps

std::shared_ptr<rad> ignis::radProps

radiation object

Definition at line 75 of file ignis.h.

◆ radType

std::string ignis::radType

radiation model name

Definition at line 76 of file ignis.h.

◆ doRadiation

bool ignis::doRadiation

radiation flag

Definition at line 77 of file ignis.h.

◆ kabs_sur

std::vector<double> ignis::kabs_sur

kabs for surroundings

Definition at line 78 of file ignis.h.

◆ awts_sur

std::vector<double> ignis::awts_sur

awts for surroundings

Definition at line 79 of file ignis.h.

◆ doLe1

bool ignis::doLe1 = false

true if doing unity Lewis numbers (default false)

Definition at line 81 of file ignis.h.

◆ Ttarget

double ignis::Ttarget

for unsteady cases, run until this max T instead of for a given time

Definition at line 83 of file ignis.h.

◆ dT

double ignis::dT

delta T increment for unsteady cases

Definition at line 84 of file ignis.h.

◆ isave

int ignis::isave

file counter for save during unsteady cases

Definition at line 85 of file ignis.h.

◆ flux_y

std::vector<std::vector<double> > ignis::flux_y

species fluxes: [I(igrid, ksp)] I(igrid,ksp) maps 2D onto 1D

Definition at line 87 of file ignis.h.

◆ flux_soot

std::vector<std::vector<double> > ignis::flux_soot

species fluxes: [I(igrid, ksoot)]

Definition at line 88 of file ignis.h.

◆ flux_h

std::vector<double> ignis::flux_h

species fluxes: [igrid]

Definition at line 89 of file ignis.h.

◆ isFlamelet

bool ignis::isFlamelet = false

true for laminar flamelet (mixture fraction coordinate)

Definition at line 91 of file ignis.h.

◆ chi

std::vector<double> ignis::chi

dissipation rate profile

Definition at line 92 of file ignis.h.

◆ chi0

double ignis::chi0

Definition at line 93 of file ignis.h.

◆ isPremixed

bool ignis::isPremixed = false

true of the case is a premixed flame, (only left boundary condition, constant mass flux through domain)

Definition at line 95 of file ignis.h.

◆ mflux

double ignis::mflux = 0.0

premixed flame mass flux (kg/m2*s)

Definition at line 96 of file ignis.h.

◆ doEnergyEqn

bool ignis::doEnergyEqn = true

for premixed flames: can solve energy equation or set T profile

Definition at line 98 of file ignis.h.

◆ LI

std::shared_ptr<linearInterp> ignis::LI

interpolator for specified temperature profiles

Definition at line 99 of file ignis.h.

◆ Tprof_h

std::vector<double> ignis::Tprof_h

temperature profile position (h is height above burner (m))

Definition at line 100 of file ignis.h.

◆ Tprof_T

std::vector<double> ignis::Tprof_T

temperature profile T values

Definition at line 101 of file ignis.h.

◆ doSoot

bool ignis::doSoot = false

soot flag

Definition at line 105 of file ignis.h.

◆ SM

std::shared_ptr<soot::sootModel> ignis::SM

soot model

Definition at line 106 of file ignis.h.

◆ SMstate

std::shared_ptr<soot::state> ignis::SMstate

holds state variables (gas and soot) for soot model

Definition at line 107 of file ignis.h.

◆ fh5

std::shared_ptr<HighFive::File> ignis::fh5

hdf5 file pointer

Definition at line 111 of file ignis.h.


The documentation for this class was generated from the following files: