Ember
|
Class which manages the main integration loop. More...
#include <flameSolver.h>
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< SourceSystem > | sourceTerms |
vector< bool > | useCVODE |
boost::ptr_vector< DiffusionSystem > | diffusionTerms |
boost::ptr_vector< TridiagonalIntegrator > | diffusionSolvers |
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 |
ScalarFunction * | strainfunc |
Function which describes strain rate a(t) and its derivative. | |
ScalarFunction * | rateMultiplierFunction |
Function which describes a multiplier for the chemical reaction rates. | |
LoggerCallback * | stateWriter |
LoggerCallback * | timeseriesWriter |
IntegratorCallback * | heatLossFunction |
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< CanteraGas > | gases |
std::mutex | gasInitMutex |
std::unique_ptr< tbb::global_control > | tbbTaskSched |
std::shared_ptr< BilinearInterpolator > | vzInterp |
Data for solving quasi-2D method-of-lines problems. | |
std::shared_ptr< BilinearInterpolator > | vrInterp |
std::shared_ptr< BilinearInterpolator > | TInterp |
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 | |
dvec & | x |
The coordinates of the grid points [m]. | |
dvec & | r |
"radius" at x[j]. | |
dvec & | rphalf |
"radius" at x[j+1/2]. | |
dvec & | hh |
Grid spacing between x[j] and x[j+1]. | |
dvec & | dlj |
Average of left and right grid spacing. | |
dvec & | cfm |
Coefficient for y[j-1] in first centered difference. | |
dvec & | cf |
Coefficient for y[j] in first centered difference. | |
dvec & | cfp |
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 ) | |
Class which manages the main integration loop.
Contains the split solvers and is responsible for the large-scale time integration.
FlameSolver::FlameSolver | ( | ) |
|
virtual |
void FlameSolver::setOptions | ( | const ConfigOptions & | options | ) |
Set options read from the configuration file.
void FlameSolver::initialize | ( | void | ) |
call to generate profiles and perform one-time setup
void FlameSolver::finalize | ( | ) |
|
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.
bool FlameSolver::checkTerminationCondition | ( | void | ) |
Determine whether to terminate integration based on reaching steady-state solution.
|
virtual |
Create prof.h5 file.
fileName | The name of the output file. |
errorFile | Setting this to 'true' will generate a more verbose output file |
updateDerivatives | Setting 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.
void FlameSolver::saveTimeSeriesData | ( | const std::string & | filename, |
bool | write | ||
) |
Collect info for out.h5 file.
double FlameSolver::getHeatReleaseRate | ( | void | ) |
Compute integral heat release rate [W/m^2].
Assumes that qDot has already been updated.
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.
double FlameSolver::getFlamePosition | ( | void | ) |
Compute position of flame (heat release centroid) [m].
|
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.
|
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.
|
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.
void FlameSolver::resizeAuxiliary | ( | ) |
Handle resizing of data structures as grid size changes.
void FlameSolver::resizeMappedArrays | ( | ) |
update data that shadows SplitSolver arrays
void FlameSolver::updateCrossTerms | ( | ) |
calculates values of cross-component terms: jSoret, sumcpj, and jCorr
void FlameSolver::updateChemicalProperties | ( | ) |
Update thermodynamic, transport, and kinetic properties.
void FlameSolver::updateChemicalProperties | ( | size_t | j1, |
size_t | j2 | ||
) |
Update thermodynamic, transport, and kinetic properties.
void FlameSolver::updateBC | ( | ) |
Set boundary condition for left edge of domain.
void FlameSolver::calculateQdot | ( | ) |
Compute heat release rate using the current temperature and mass fractions.
void FlameSolver::correctMassFractions | ( | ) |
Correct the drift of the total mass fractions and reset any negative mass fractions.
void FlameSolver::setDiffusionSolverState | ( | double | tInitial | ) |
void FlameSolver::setConvectionSolverState | ( | double | tInitial | ) |
void FlameSolver::setProductionSolverState | ( | double | tInitial | ) |
|
virtual |
Implements SplitSolver.
|
virtual |
Implements SplitSolver.
|
virtual |
Implements SplitSolver.
void FlameSolver::integrateProductionTerms | ( | size_t | j1, |
size_t | j2 | ||
) |
void FlameSolver::integrateDiffusionTerms | ( | size_t | k1, |
size_t | k2 | ||
) |
void FlameSolver::update_xStag | ( | const double | t, |
const bool | updateIntError | ||
) |
double FlameSolver::targetFlamePosition | ( | double | t | ) |
[m]
void FlameSolver::printPerformanceStats | ( | void | ) |
void FlameSolver::printPerfString | ( | std::ostream & | stats, |
const std::string & | label, | ||
const PerfTimer & | T | ||
) |
dvector FlameSolver::timeVector |
Time [s] at the end of every ConfigOptions::outputStepInterval timesteps.
dvector FlameSolver::heatReleaseRate |
Integral heat release rate [W/m^2] at the times in timeVector.
long int FlameSolver::nTotal |
total number of timesteps taken
int FlameSolver::nRegrid |
number of time steps since regridding/adaptation
int FlameSolver::nOutput |
number of time steps since storing integral flame parameters
int FlameSolver::nProfile |
number of time steps since saving flame profiles
int FlameSolver::nTerminate |
number of steps since last checking termination condition
int FlameSolver::nCurrentState |
number of time steps since profNow.h5 and out.h5 were written
double FlameSolver::terminationCondition |
Current value to be compared with terminationTolerance.
double FlameSolver::tOutput |
time of next integral flame parameters output (this step)
double FlameSolver::tRegrid |
time of next regridding
double FlameSolver::tProfile |
time of next profile output
boost::ptr_vector<SourceSystem> FlameSolver::sourceTerms |
vector<bool> FlameSolver::useCVODE |
boost::ptr_vector<DiffusionSystem> FlameSolver::diffusionTerms |
boost::ptr_vector<TridiagonalIntegrator> FlameSolver::diffusionSolvers |
ConvectionSystemSplit FlameSolver::convectionSystem |
System representing the convection terms and containing their integrators.
double FlameSolver::rhou |
double FlameSolver::rhob |
double FlameSolver::rhoLeft |
double FlameSolver::rhoRight |
double FlameSolver::Tleft |
double FlameSolver::Tright |
dvec FlameSolver::Yleft |
dvec FlameSolver::Yright |
size_t FlameSolver::nSpec |
Number of chemical species.
size_t FlameSolver::nVars |
Number of state variables at each grid point (nSpec + 2)
VecMap FlameSolver::U |
normalized tangential velocity (u*a/u_inf) [1/s]
VecMap FlameSolver::T |
temperature [K]
dvec FlameSolver::rho |
density [kg/m^3]
dvec FlameSolver::drhodt |
time derivative of density [kg/m^3*s]
dvec FlameSolver::jCorr |
Correction to ensure sum of mass fractions = 1.
dvec FlameSolver::sumcpj |
part of the enthalpy flux term
dvec FlameSolver::qDot |
Heat release rate [W/m^3].
dmatrix FlameSolver::wDot |
species production rates [kmol/m^3*s]
dvec FlameSolver::Wmx |
mixture molecular weight [kg/kmol]
dvec FlameSolver::W |
species molecular weights [kg/kmol]
dvec FlameSolver::mu |
dynamic viscosity [Pa*s]
dvec FlameSolver::lambda |
thermal conductivity [W/m*K]
dvec FlameSolver::cp |
mixture heat capacity [J/kg*K]
dmatrix FlameSolver::cpSpec |
species molar heat capacities [J/kmol*K]
dmatrix FlameSolver::rhoD |
density * diffusivity [kg/m*s]
dmatrix FlameSolver::Dkt |
thermal diffusivity
dmatrix FlameSolver::hk |
species molar enthalpies [J/kmol]
dmatrix FlameSolver::jFick |
Fickian mass flux [kg/m^2*s].
dmatrix FlameSolver::jSoret |
Soret mass flux [kg/m^2*s].
DiffusionSystem FlameSolver::jCorrSystem |
TridiagonalIntegrator FlameSolver::jCorrSolver |
ScalarFunction* FlameSolver::strainfunc |
Function which describes strain rate a(t) and its derivative.
ScalarFunction* FlameSolver::rateMultiplierFunction |
Function which describes a multiplier for the chemical reaction rates.
LoggerCallback* FlameSolver::stateWriter |
LoggerCallback* FlameSolver::timeseriesWriter |
IntegratorCallback* FlameSolver::heatLossFunction |
double FlameSolver::rVcenter |
mass flux at centerline [kg/m^2 or kg/m*rad*s]
double FlameSolver::rVzero |
mass flux at j=0
double FlameSolver::tFlamePrev |
double FlameSolver::tFlameNext |
double FlameSolver::xFlameTarget |
double FlameSolver::xFlameActual |
double FlameSolver::flamePosIntegralError |
CanteraGas FlameSolver::gas |
Cantera data.
tbb::enumerable_thread_specific<CanteraGas> FlameSolver::gases |
std::mutex FlameSolver::gasInitMutex |
std::unique_ptr<tbb::global_control> FlameSolver::tbbTaskSched |
std::shared_ptr<BilinearInterpolator> FlameSolver::vzInterp |
Data for solving quasi-2D method-of-lines problems.
std::shared_ptr<BilinearInterpolator> FlameSolver::vrInterp |
std::shared_ptr<BilinearInterpolator> FlameSolver::TInterp |
PerfTimer FlameSolver::totalTimer |
PerfTimer FlameSolver::setupTimer |
PerfTimer FlameSolver::reactionTimer |
PerfTimer FlameSolver::diffusionTimer |
PerfTimer FlameSolver::convectionTimer |
PerfTimer FlameSolver::regridTimer |
PerfTimer FlameSolver::reactionRatesTimer |
PerfTimer FlameSolver::transportTimer |
PerfTimer FlameSolver::thermoTimer |
PerfTimer FlameSolver::jacobianTimer |
PerfTimer FlameSolver::conductivityTimer |
PerfTimer FlameSolver::viscosityTimer |
PerfTimer FlameSolver::diffusivityTimer |