Ember
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
FlameSolver Class Reference

Class which manages the main integration loop. More...

#include <flameSolver.h>

Inheritance diagram for FlameSolver:
[legend]

Public Member Functions

 FlameSolver ()
 
virtual ~FlameSolver ()
 
void setOptions (const ConfigOptions &options)
 Set options read from the configuration file.
 
void initialize ()
 call to generate profiles and perform one-time setup
 
void finalize ()
 
virtual void loadProfile ()
 Load initial temperature, mass fraction and velocity profiles.
 
bool checkTerminationCondition (void)
 Determine whether to terminate integration based on reaching steady-state solution.
 
void writeStateFile (const std::string &fileName="", bool errorFile=false, bool updateDerivatives=true)
 Create prof.h5 file.
 
void saveTimeSeriesData (const std::string &filename, bool write)
 Collect info for out.h5 file.
 
double getHeatReleaseRate (void)
 Compute integral heat release rate [W/m^2].
 
double getConsumptionSpeed (void)
 Compute flame consumption speed [m/s].
 
double getFlamePosition (void)
 Compute position of flame (heat release centroid) [m].
 
void setupStep ()
 Take actions at the start of the time step.
 
void prepareIntegrators ()
 Prepare split term integrators.
 
int finishStep ()
 Take actions at the end of a time step.
 
void resizeAuxiliary ()
 Handle resizing of data structures as grid size changes.
 
void resizeMappedArrays ()
 update data that shadows SplitSolver arrays
 
void updateCrossTerms ()
 calculates values of cross-component terms: jSoret, sumcpj, and jCorr
 
void updateChemicalProperties ()
 Update thermodynamic, transport, and kinetic properties.
 
void updateChemicalProperties (size_t j1, size_t j2)
 Update thermodynamic, transport, and kinetic properties.
 
void updateBC ()
 Set boundary condition for left edge of domain.
 
void calculateQdot ()
 Compute heat release rate using the current temperature and mass fractions.
 
void correctMassFractions ()
 Correct the drift of the total mass fractions and reset any negative mass fractions.
 
void setDiffusionSolverState (double tInitial)
 
void setConvectionSolverState (double tInitial)
 
void setProductionSolverState (double tInitial)
 
void integrateConvectionTerms ()
 
void integrateProductionTerms ()
 
void integrateDiffusionTerms ()
 
void integrateProductionTerms (size_t j1, size_t j2)
 
void integrateDiffusionTerms (size_t k1, size_t k2)
 
void rollVectorVector (vector< dvector > &vv, const dmatrix &M) const
 
void unrollVectorVector (vector< dvector > &vv, dmatrix &M, size_t i) const
 
void update_xStag (const double t, const bool updateIntError)
 
double targetFlamePosition (double t)
 [m]
 
void printPerformanceStats (void)
 
void printPerfString (std::ostream &stats, const std::string &label, const PerfTimer &T)
 
- Public Member Functions inherited from GridBased
 GridBased ()
 
virtual ~GridBased ()
 
virtual void setGrid (const OneDimGrid &grid)
 Copy the specified grid to this object.
 
- Public Member Functions inherited from SplitSolver
 SplitSolver ()
 
virtual ~SplitSolver ()
 
void setOptions (const ConfigOptions &_options)
 Set configuration options needed by the solver.
 
void resize (index_t nRows, index_t nCols)
 Set the size of all internal arrays when the problem size changes.
 
void calculateTimeDerivatives (double dt)
 Compute estimated time derivatives for each split term based on deltas for each integrator phase.
 
int step ()
 Take one global timestep.
 
virtual void setupStep ()=0
 Take actions at the start of the time step.
 
virtual int finishStep ()=0
 Take actions at the end of a time step.
 
virtual void prepareIntegrators ()=0
 Prepare split term integrators.
 
void integrateSplitTerms (double t, double dt)
 Take one global timestep by integrating the split terms and computing the deltas for each integrator phase.
 
virtual void writeStateFile (const std::string &fileName="", bool errorFile=false, bool updateDerivatives=true)
 Create prof.h5 file.
 

Public Attributes

dvector timeVector
 Time [s] at the end of every ConfigOptions::outputStepInterval timesteps.
 
dvector heatReleaseRate
 Integral heat release rate [W/m^2] at the times in timeVector.
 
long int nTotal
 total number of timesteps taken
 
int nRegrid
 number of time steps since regridding/adaptation
 
int nOutput
 number of time steps since storing integral flame parameters
 
int nProfile
 number of time steps since saving flame profiles
 
int nTerminate
 number of steps since last checking termination condition
 
int nCurrentState
 number of time steps since profNow.h5 and out.h5 were written
 
double terminationCondition
 Current value to be compared with terminationTolerance.
 
double tOutput
 time of next integral flame parameters output (this step)
 
double tRegrid
 time of next regridding
 
double tProfile
 time of next profile output
 
boost::ptr_vector< SourceSystemsourceTerms
 
vector< bool > useCVODE
 
boost::ptr_vector< DiffusionSystemdiffusionTerms
 
boost::ptr_vector< TridiagonalIntegratordiffusionSolvers
 
ConvectionSystemSplit convectionSystem
 System representing the convection terms and containing their integrators.
 
double rhou
 
double rhob
 
double rhoLeft
 
double rhoRight
 
double Tleft
 
double Tright
 
dvec Yleft
 
dvec Yright
 
size_t nSpec
 Number of chemical species.
 
size_t nVars
 Number of state variables at each grid point (nSpec + 2)
 
VecMap U
 normalized tangential velocity (u*a/u_inf) [1/s]
 
VecMap T
 temperature [K]
 
MatrixMap Y
 species mass fractions, Y(k,j) [-]
 
dvec rho
 density [kg/m^3]
 
dvec drhodt
 time derivative of density [kg/m^3*s]
 
dvec jCorr
 Correction to ensure sum of mass fractions = 1.
 
dvec sumcpj
 part of the enthalpy flux term
 
dvec qDot
 Heat release rate [W/m^3].
 
dmatrix wDot
 species production rates [kmol/m^3*s]
 
dvec Wmx
 mixture molecular weight [kg/kmol]
 
dvec W
 species molecular weights [kg/kmol]
 
dvec mu
 dynamic viscosity [Pa*s]
 
dvec lambda
 thermal conductivity [W/m*K]
 
dvec cp
 mixture heat capacity [J/kg*K]
 
dmatrix cpSpec
 species molar heat capacities [J/kmol*K]
 
dmatrix rhoD
 density * diffusivity [kg/m*s]
 
dmatrix Dkt
 thermal diffusivity
 
dmatrix hk
 species molar enthalpies [J/kmol]
 
dmatrix jFick
 Fickian mass flux [kg/m^2*s].
 
dmatrix jSoret
 Soret mass flux [kg/m^2*s].
 
DiffusionSystem jCorrSystem
 
TridiagonalIntegrator jCorrSolver
 
ScalarFunctionstrainfunc
 Function which describes strain rate a(t) and its derivative.
 
ScalarFunctionrateMultiplierFunction
 Function which describes a multiplier for the chemical reaction rates.
 
LoggerCallbackstateWriter
 
LoggerCallbacktimeseriesWriter
 
IntegratorCallbackheatLossFunction
 
double rVcenter
 mass flux at centerline [kg/m^2 or kg/m*rad*s]
 
double rVzero
 mass flux at j=0
 
double tFlamePrev
 
double tFlameNext
 
double xFlameTarget
 
double xFlameActual
 
double flamePosIntegralError
 
CanteraGas gas
 Cantera data.
 
tbb::enumerable_thread_specific< CanteraGasgases
 
std::mutex gasInitMutex
 
std::unique_ptr< tbb::global_controltbbTaskSched
 
std::shared_ptr< BilinearInterpolatorvzInterp
 Data for solving quasi-2D method-of-lines problems.
 
std::shared_ptr< BilinearInterpolatorvrInterp
 
std::shared_ptr< BilinearInterpolatorTInterp
 
PerfTimer totalTimer
 
PerfTimer setupTimer
 
PerfTimer reactionTimer
 
PerfTimer diffusionTimer
 
PerfTimer convectionTimer
 
PerfTimer regridTimer
 
PerfTimer reactionRatesTimer
 
PerfTimer transportTimer
 
PerfTimer thermoTimer
 
PerfTimer jacobianTimer
 
PerfTimer conductivityTimer
 
PerfTimer viscosityTimer
 
PerfTimer diffusivityTimer
 
- Public Attributes inherited from GridBased
OneDimGrid grid
 the actual grid
 
- Public Attributes inherited from SplitSolver
dmatrix state
 Current state.
 
dmatrix startState
 State at the start of the current integrator stage.
 
dmatrix deltaConv
 Change in state due to convection terms during the current time step.
 
dmatrix deltaDiff
 Change in state due to diffusion terms during the current time step.
 
dmatrix deltaProd
 Change in state due to production terms during the current time step.
 
dmatrix ddtConv
 Estimated time derivatives of convection terms, based on the deltaConv of the previous time step.
 
dmatrix ddtDiff
 Estimated time derivatives of the diffusion terms, based on the deltaDiff of the previous time step.
 
dmatrix ddtProd
 Estimated time derivatives of production terms, based on the deltaProd of the previous time step.
 
dmatrix ddtCross
 Time derivative associated with terms taken as constant over a time step.
 
dmatrix splitConstConv
 Splitting constants for the convection terms.
 
dmatrix splitConstDiff
 Splitting constants for the diffusion terms.
 
dmatrix splitConstProd
 Splitting constants for the production terms.
 
PerfTimer splitTimer
 Timer for amount of time used by splitting procedure.
 
ConfigOptions options
 Options read from the input file.
 
double tStart
 Integrator start time.
 
double tEnd
 Integrator termination time (upper bound)
 
double tNow
 Current time reached by the integrator.
 
double tStageStart
 Start time for the active integrator stage.
 
double tStageEnd
 End time for the active integrator stage.
 
double t
 start time of the current global timestep
 
double dt
 Global timestep used by the split solver.
 

Additional Inherited Members

- Protected Attributes inherited from GridBased
dvecx
 The coordinates of the grid points [m].
 
dvecr
 "radius" at x[j].
 
dvecrphalf
 "radius" at x[j+1/2].
 
dvechh
 Grid spacing between x[j] and x[j+1].
 
dvecdlj
 Average of left and right grid spacing.
 
dveccfm
 Coefficient for y[j-1] in first centered difference.
 
dveccf
 Coefficient for y[j] in first centered difference.
 
dveccfp
 Coefficient for y[j+1] in first centered difference.
 
int & alpha
 curved grid exponent.
 
int & beta
 curved grid exponent.
 
size_t & nPoints
 number of grid point
 
size_t & jj
 index of last grid point (== nPoints-1)
 

Detailed Description

Class which manages the main integration loop.

Contains the split solvers and is responsible for the large-scale time integration.

Constructor & Destructor Documentation

◆ FlameSolver()

FlameSolver::FlameSolver ( )

◆ ~FlameSolver()

FlameSolver::~FlameSolver ( )
virtual

Member Function Documentation

◆ setOptions()

void FlameSolver::setOptions ( const ConfigOptions options)

Set options read from the configuration file.

◆ initialize()

void FlameSolver::initialize ( void  )

call to generate profiles and perform one-time setup

◆ finalize()

void FlameSolver::finalize ( )

◆ loadProfile()

void FlameSolver::loadProfile ( void  )
virtual

Load initial temperature, mass fraction and velocity profiles.

Profiles are loaded either from an HDF5 "restart" file or from the the ConfigOptions object, populated from the Python "Config" object.

◆ checkTerminationCondition()

bool FlameSolver::checkTerminationCondition ( void  )

Determine whether to terminate integration based on reaching steady-state solution.

◆ writeStateFile()

void FlameSolver::writeStateFile ( const std::string &  fileName = "",
bool  errorFile = false,
bool  updateDerivatives = true 
)
virtual

Create prof.h5 file.

Parameters
fileNameThe name of the output file.
errorFileSetting this to 'true' will generate a more verbose output file
updateDerivativesSetting this to 'true' will compute time derivatives based on the current state. This flag should only set if all of the solvers and auxillary arrays are sized corresponding to the current state, e.g. after integration but but before regridding.

Reimplemented from SplitSolver.

◆ saveTimeSeriesData()

void FlameSolver::saveTimeSeriesData ( const std::string &  filename,
bool  write 
)

Collect info for out.h5 file.

◆ getHeatReleaseRate()

double FlameSolver::getHeatReleaseRate ( void  )

Compute integral heat release rate [W/m^2].

Assumes that qDot has already been updated.

◆ getConsumptionSpeed()

double FlameSolver::getConsumptionSpeed ( void  )

Compute flame consumption speed [m/s].

Computed by integrating the energy equation:

\[ S_c = \frac{\int_{-\infty}^\infty \dot{q}/c_p}{\rho_u(T_b-T_u)} \]

Assumes that qDot has already been updated.

◆ getFlamePosition()

double FlameSolver::getFlamePosition ( void  )

Compute position of flame (heat release centroid) [m].

◆ setupStep()

void FlameSolver::setupStep ( )
virtual

Take actions at the start of the time step.

Used by derived classes to perform actions that need to be done at the start of a time step, before the splitting constants are computed.

Implements SplitSolver.

◆ prepareIntegrators()

void FlameSolver::prepareIntegrators ( )
virtual

Prepare split term integrators.

Used by derived classes to assign split constants to the integrators for the individual terms and perform any other actions that should take place

Implements SplitSolver.

◆ finishStep()

int FlameSolver::finishStep ( )
virtual

Take actions at the end of a time step.

Used by derived classes to perform any actions that need to be done at the end of the time step, after the time derivatives have been calculated based on the results of integrating the split equations.

Implements SplitSolver.

◆ resizeAuxiliary()

void FlameSolver::resizeAuxiliary ( )

Handle resizing of data structures as grid size changes.

◆ resizeMappedArrays()

void FlameSolver::resizeMappedArrays ( )

update data that shadows SplitSolver arrays

◆ updateCrossTerms()

void FlameSolver::updateCrossTerms ( )

calculates values of cross-component terms: jSoret, sumcpj, and jCorr

◆ updateChemicalProperties() [1/2]

void FlameSolver::updateChemicalProperties ( )

Update thermodynamic, transport, and kinetic properties.

◆ updateChemicalProperties() [2/2]

void FlameSolver::updateChemicalProperties ( size_t  j1,
size_t  j2 
)

Update thermodynamic, transport, and kinetic properties.

◆ updateBC()

void FlameSolver::updateBC ( )

Set boundary condition for left edge of domain.

◆ calculateQdot()

void FlameSolver::calculateQdot ( )

Compute heat release rate using the current temperature and mass fractions.

◆ correctMassFractions()

void FlameSolver::correctMassFractions ( )

Correct the drift of the total mass fractions and reset any negative mass fractions.

◆ setDiffusionSolverState()

void FlameSolver::setDiffusionSolverState ( double  tInitial)

◆ setConvectionSolverState()

void FlameSolver::setConvectionSolverState ( double  tInitial)

◆ setProductionSolverState()

void FlameSolver::setProductionSolverState ( double  tInitial)

◆ integrateConvectionTerms()

void FlameSolver::integrateConvectionTerms ( )
virtual

Implements SplitSolver.

◆ integrateProductionTerms() [1/2]

void FlameSolver::integrateProductionTerms ( )
virtual

Implements SplitSolver.

◆ integrateDiffusionTerms() [1/2]

void FlameSolver::integrateDiffusionTerms ( )
virtual

Implements SplitSolver.

◆ integrateProductionTerms() [2/2]

void FlameSolver::integrateProductionTerms ( size_t  j1,
size_t  j2 
)

◆ integrateDiffusionTerms() [2/2]

void FlameSolver::integrateDiffusionTerms ( size_t  k1,
size_t  k2 
)

◆ rollVectorVector()

void FlameSolver::rollVectorVector ( vector< dvector > &  vv,
const dmatrix M 
) const

◆ unrollVectorVector()

void FlameSolver::unrollVectorVector ( vector< dvector > &  vv,
dmatrix M,
size_t  i 
) const

◆ update_xStag()

void FlameSolver::update_xStag ( const double  t,
const bool  updateIntError 
)

◆ targetFlamePosition()

double FlameSolver::targetFlamePosition ( double  t)

[m]

◆ printPerformanceStats()

void FlameSolver::printPerformanceStats ( void  )

◆ printPerfString()

void FlameSolver::printPerfString ( std::ostream &  stats,
const std::string &  label,
const PerfTimer T 
)

Member Data Documentation

◆ timeVector

dvector FlameSolver::timeVector

Time [s] at the end of every ConfigOptions::outputStepInterval timesteps.

◆ heatReleaseRate

dvector FlameSolver::heatReleaseRate

Integral heat release rate [W/m^2] at the times in timeVector.

◆ nTotal

long int FlameSolver::nTotal

total number of timesteps taken

◆ nRegrid

int FlameSolver::nRegrid

number of time steps since regridding/adaptation

◆ nOutput

int FlameSolver::nOutput

number of time steps since storing integral flame parameters

◆ nProfile

int FlameSolver::nProfile

number of time steps since saving flame profiles

◆ nTerminate

int FlameSolver::nTerminate

number of steps since last checking termination condition

◆ nCurrentState

int FlameSolver::nCurrentState

number of time steps since profNow.h5 and out.h5 were written

◆ terminationCondition

double FlameSolver::terminationCondition

Current value to be compared with terminationTolerance.

◆ tOutput

double FlameSolver::tOutput

time of next integral flame parameters output (this step)

◆ tRegrid

double FlameSolver::tRegrid

time of next regridding

◆ tProfile

double FlameSolver::tProfile

time of next profile output

◆ sourceTerms

boost::ptr_vector<SourceSystem> FlameSolver::sourceTerms

◆ useCVODE

vector<bool> FlameSolver::useCVODE

◆ diffusionTerms

boost::ptr_vector<DiffusionSystem> FlameSolver::diffusionTerms

◆ diffusionSolvers

boost::ptr_vector<TridiagonalIntegrator> FlameSolver::diffusionSolvers

◆ convectionSystem

ConvectionSystemSplit FlameSolver::convectionSystem

System representing the convection terms and containing their integrators.

◆ rhou

double FlameSolver::rhou

◆ rhob

double FlameSolver::rhob

◆ rhoLeft

double FlameSolver::rhoLeft

◆ rhoRight

double FlameSolver::rhoRight

◆ Tleft

double FlameSolver::Tleft

◆ Tright

double FlameSolver::Tright

◆ Yleft

dvec FlameSolver::Yleft

◆ Yright

dvec FlameSolver::Yright

◆ nSpec

size_t FlameSolver::nSpec

Number of chemical species.

◆ nVars

size_t FlameSolver::nVars

Number of state variables at each grid point (nSpec + 2)

◆ U

VecMap FlameSolver::U

normalized tangential velocity (u*a/u_inf) [1/s]

◆ T

VecMap FlameSolver::T

temperature [K]

◆ Y

MatrixMap FlameSolver::Y

species mass fractions, Y(k,j) [-]

◆ rho

dvec FlameSolver::rho

density [kg/m^3]

◆ drhodt

dvec FlameSolver::drhodt

time derivative of density [kg/m^3*s]

◆ jCorr

dvec FlameSolver::jCorr

Correction to ensure sum of mass fractions = 1.

◆ sumcpj

dvec FlameSolver::sumcpj

part of the enthalpy flux term

◆ qDot

dvec FlameSolver::qDot

Heat release rate [W/m^3].

◆ wDot

dmatrix FlameSolver::wDot

species production rates [kmol/m^3*s]

◆ Wmx

dvec FlameSolver::Wmx

mixture molecular weight [kg/kmol]

◆ W

dvec FlameSolver::W

species molecular weights [kg/kmol]

◆ mu

dvec FlameSolver::mu

dynamic viscosity [Pa*s]

◆ lambda

dvec FlameSolver::lambda

thermal conductivity [W/m*K]

◆ cp

dvec FlameSolver::cp

mixture heat capacity [J/kg*K]

◆ cpSpec

dmatrix FlameSolver::cpSpec

species molar heat capacities [J/kmol*K]

◆ rhoD

dmatrix FlameSolver::rhoD

density * diffusivity [kg/m*s]

◆ Dkt

dmatrix FlameSolver::Dkt

thermal diffusivity

◆ hk

dmatrix FlameSolver::hk

species molar enthalpies [J/kmol]

◆ jFick

dmatrix FlameSolver::jFick

Fickian mass flux [kg/m^2*s].

◆ jSoret

dmatrix FlameSolver::jSoret

Soret mass flux [kg/m^2*s].

◆ jCorrSystem

DiffusionSystem FlameSolver::jCorrSystem

◆ jCorrSolver

TridiagonalIntegrator FlameSolver::jCorrSolver

◆ strainfunc

ScalarFunction* FlameSolver::strainfunc

Function which describes strain rate a(t) and its derivative.

◆ rateMultiplierFunction

ScalarFunction* FlameSolver::rateMultiplierFunction

Function which describes a multiplier for the chemical reaction rates.

◆ stateWriter

LoggerCallback* FlameSolver::stateWriter

◆ timeseriesWriter

LoggerCallback* FlameSolver::timeseriesWriter

◆ heatLossFunction

IntegratorCallback* FlameSolver::heatLossFunction

◆ rVcenter

double FlameSolver::rVcenter

mass flux at centerline [kg/m^2 or kg/m*rad*s]

◆ rVzero

double FlameSolver::rVzero

mass flux at j=0

◆ tFlamePrev

double FlameSolver::tFlamePrev

◆ tFlameNext

double FlameSolver::tFlameNext

◆ xFlameTarget

double FlameSolver::xFlameTarget

◆ xFlameActual

double FlameSolver::xFlameActual

◆ flamePosIntegralError

double FlameSolver::flamePosIntegralError

◆ gas

CanteraGas FlameSolver::gas

Cantera data.

◆ gases

◆ gasInitMutex

std::mutex FlameSolver::gasInitMutex

◆ tbbTaskSched

std::unique_ptr<tbb::global_control> FlameSolver::tbbTaskSched

◆ vzInterp

std::shared_ptr<BilinearInterpolator> FlameSolver::vzInterp

Data for solving quasi-2D method-of-lines problems.

◆ vrInterp

std::shared_ptr<BilinearInterpolator> FlameSolver::vrInterp

◆ TInterp

std::shared_ptr<BilinearInterpolator> FlameSolver::TInterp

◆ totalTimer

PerfTimer FlameSolver::totalTimer

◆ setupTimer

PerfTimer FlameSolver::setupTimer

◆ reactionTimer

PerfTimer FlameSolver::reactionTimer

◆ diffusionTimer

PerfTimer FlameSolver::diffusionTimer

◆ convectionTimer

PerfTimer FlameSolver::convectionTimer

◆ regridTimer

PerfTimer FlameSolver::regridTimer

◆ reactionRatesTimer

PerfTimer FlameSolver::reactionRatesTimer

◆ transportTimer

PerfTimer FlameSolver::transportTimer

◆ thermoTimer

PerfTimer FlameSolver::thermoTimer

◆ jacobianTimer

PerfTimer FlameSolver::jacobianTimer

◆ conductivityTimer

PerfTimer FlameSolver::conductivityTimer

◆ viscosityTimer

PerfTimer FlameSolver::viscosityTimer

◆ diffusivityTimer

PerfTimer FlameSolver::diffusivityTimer

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