Ember
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Protected Attributes | List of all members
SourceSystem Class Referenceabstract

Base class used to integrate the chemical source term at a single point. More...

#include <sourceSystem.h>

Inheritance diagram for SourceSystem:
[legend]

Public Member Functions

 SourceSystem ()
 
virtual ~SourceSystem ()
 
virtual void setState (double tInitial, double uu, double tt, const dvec &yy)=0
 Set the initial condition for integrator.
 
virtual int integrateToTime (double tf)=0
 Take as many steps as needed to reach tf.
 
virtual int integrateOneStep (double tf)=0
 Take one step toward tf without stepping past it.
 
virtual double time () const =0
 Current integrator time, relative to tInitial.
 
virtual void unroll_y ()=0
 Extract current internal integrator state into U, T, and Y.
 
virtual void setDebug (bool debug_)
 
virtual std::string getStats ()=0
 Return a string indicating the number of internal timesteps taken.
 
void updateThermo ()
 Compute thermodynamic properties (density, enthalpies, heat capacities) from the current state variables.
 
double getQdotIgniter (double t)
 Calculate the heat release rate associated with a possible external ignition source.
 
void setGas (CanteraGas *_gas)
 Set the CanteraGas object to use for thermodynamic and kinetic property calculations.
 
virtual void initialize (size_t nSpec)
 Resize internal arrays for a problem of the specified size (nSpec+2)
 
virtual void setOptions (ConfigOptions &options_)
 Set integrator tolerances and other parameters.
 
void setTimers (PerfTimer *reactionRates, PerfTimer *thermo, PerfTimer *jacobian)
 
void setPosition (size_t j, double x)
 Set the index j and position x that are represented by this system.
 
void setStrainFunction (ScalarFunction *f)
 Set the function used to compute the strain rate as a function of time.
 
void setRateMultiplierFunction (ScalarFunction *f)
 Set the function used to compute the reaction rate multiplier.
 
void setHeatLossFunction (IntegratorCallback *f)
 Set the function used to compute the heat loss rate to the environment.
 
void setRhou (double _rhou)
 Set the density of the unburned mixture.
 
void resetSplitConstants ()
 Set all the balanced splitting constants to zero.
 
void setupQuasi2d (std::shared_ptr< BilinearInterpolator > vzInterp, std::shared_ptr< BilinearInterpolator > TInterp)
 Assign the interpolators used for solving quasi-2D problems.
 
virtual void writeState (std::ostream &out, bool init)
 Write the current values of the state variables, formatted to be read by Python, to the specified stream.
 
virtual void writeJacobian (std::ostream &out)
 

Public Attributes

double U
 tangential velocity
 
double T
 temperature
 
dvec Y
 species mass fraction
 
dvec splitConst
 Extra constant term introduced by splitting.
 

Protected Attributes

bool debug
 
ConfigOptionsoptions
 
CanteraGasgas
 Cantera data.
 
PerfTimerreactionRatesTimer
 Timer for time spent evaluating reaction rates.
 
PerfTimerthermoTimer
 Timer for time spent evaluating thermodynamic properties.
 
PerfTimerjacobianTimer
 Timer for time spent evaluating and factorizing the Jacobian.
 
ScalarFunctionstrainFunction
 A class that provides the strain rate and its time derivative.
 
ScalarFunctionrateMultiplierFunction
 Provides a multiplier (optional) for the production terms.
 
IntegratorCallbackheatLoss
 Heat loss rate.
 
size_t nSpec
 number of species
 
int j
 grid index for this system
 
double x
 grid position for this system
 
double rhou
 density of the unburned gas
 
double qDot
 heat release rate per unit volume [W/m^3]
 
double qLoss
 heat loss to the environment [W/m^3]
 
double rho
 density [kg/m^3]
 
double cp
 specific heat capacity (average) [J/kg*K]
 
dvec cpSpec
 species specific heat capacity [J/mol*K]
 
double Wmx
 mixture molecular weight [kg/mol]
 
dvec W
 species molecular weights [kg/kmol]
 
dvec hk
 species enthalpies [J/kmol]
 
bool quasi2d
 Flag set to 'true' when solving a quasi-2D problem with prescribed velocity and temperature fields.
 
std::shared_ptr< BilinearInterpolatorvzInterp
 An interpolator for computing the axial (z) velocity when solving a quasi-2D problem.
 
std::shared_ptr< BilinearInterpolatorTInterp
 An interpolator for computing the temperature when solving a quasi-2D problem.
 

Detailed Description

Base class used to integrate the chemical source term at a single point.

Constructor & Destructor Documentation

◆ SourceSystem()

SourceSystem::SourceSystem ( )

◆ ~SourceSystem()

virtual SourceSystem::~SourceSystem ( )
inlinevirtual

Member Function Documentation

◆ setState()

virtual void SourceSystem::setState ( double  tInitial,
double  uu,
double  tt,
const dvec yy 
)
pure virtual

Set the initial condition for integrator.

Parameters
tInitialintegrator start time
uutangential velocity
tttemperature
yyvector of species mass fractions

Implemented in SourceSystemCVODE, and SourceSystemQSS.

◆ integrateToTime()

virtual int SourceSystem::integrateToTime ( double  tf)
pure virtual

Take as many steps as needed to reach tf.

tf is relative to tInitial.

Implemented in SourceSystemCVODE, and SourceSystemQSS.

◆ integrateOneStep()

virtual int SourceSystem::integrateOneStep ( double  tf)
pure virtual

Take one step toward tf without stepping past it.

tf is relative to tInitial.

Implemented in SourceSystemCVODE, and SourceSystemQSS.

◆ time()

virtual double SourceSystem::time ( ) const
pure virtual

Current integrator time, relative to tInitial.

Implemented in SourceSystemCVODE, and SourceSystemQSS.

◆ unroll_y()

virtual void SourceSystem::unroll_y ( )
pure virtual

Extract current internal integrator state into U, T, and Y.

Implemented in SourceSystemCVODE, and SourceSystemQSS.

◆ setDebug()

virtual void SourceSystem::setDebug ( bool  debug_)
inlinevirtual

◆ getStats()

virtual std::string SourceSystem::getStats ( )
pure virtual

Return a string indicating the number of internal timesteps taken.

Implemented in SourceSystemCVODE, and SourceSystemQSS.

◆ updateThermo()

void SourceSystem::updateThermo ( )

Compute thermodynamic properties (density, enthalpies, heat capacities) from the current state variables.

◆ getQdotIgniter()

double SourceSystem::getQdotIgniter ( double  t)

Calculate the heat release rate associated with a possible external ignition source.

◆ setGas()

void SourceSystem::setGas ( CanteraGas _gas)
inline

Set the CanteraGas object to use for thermodynamic and kinetic property calculations.

◆ initialize()

void SourceSystem::initialize ( size_t  nSpec)
virtual

Resize internal arrays for a problem of the specified size (nSpec+2)

Reimplemented in SourceSystemCVODE, and SourceSystemQSS.

◆ setOptions()

void SourceSystem::setOptions ( ConfigOptions options_)
virtual

Set integrator tolerances and other parameters.

Reimplemented in SourceSystemCVODE, and SourceSystemQSS.

◆ setTimers()

void SourceSystem::setTimers ( PerfTimer reactionRates,
PerfTimer thermo,
PerfTimer jacobian 
)

◆ setPosition()

void SourceSystem::setPosition ( size_t  j,
double  x 
)

Set the index j and position x that are represented by this system.

◆ setStrainFunction()

void SourceSystem::setStrainFunction ( ScalarFunction f)
inline

Set the function used to compute the strain rate as a function of time.

◆ setRateMultiplierFunction()

void SourceSystem::setRateMultiplierFunction ( ScalarFunction f)
inline

Set the function used to compute the reaction rate multiplier.

◆ setHeatLossFunction()

void SourceSystem::setHeatLossFunction ( IntegratorCallback f)
inline

Set the function used to compute the heat loss rate to the environment.

◆ setRhou()

void SourceSystem::setRhou ( double  _rhou)
inline

Set the density of the unburned mixture.

This value appears in the source term of the momentum equation.

◆ resetSplitConstants()

void SourceSystem::resetSplitConstants ( )
inline

Set all the balanced splitting constants to zero.

◆ setupQuasi2d()

void SourceSystem::setupQuasi2d ( std::shared_ptr< BilinearInterpolator vzInterp,
std::shared_ptr< BilinearInterpolator TInterp 
)

Assign the interpolators used for solving quasi-2D problems.

◆ writeState()

void SourceSystem::writeState ( std::ostream &  out,
bool  init 
)
virtual

Write the current values of the state variables, formatted to be read by Python, to the specified stream.

Call with init=true when first called to include initializers for the variables.

Reimplemented in SourceSystemCVODE.

◆ writeJacobian()

virtual void SourceSystem::writeJacobian ( std::ostream &  out)
inlinevirtual

Reimplemented in SourceSystemCVODE.

Member Data Documentation

◆ U

double SourceSystem::U

tangential velocity

◆ T

double SourceSystem::T

temperature

◆ Y

dvec SourceSystem::Y

species mass fraction

◆ splitConst

dvec SourceSystem::splitConst

Extra constant term introduced by splitting.

◆ debug

bool SourceSystem::debug
protected

◆ options

ConfigOptions* SourceSystem::options
protected

◆ gas

CanteraGas* SourceSystem::gas
protected

Cantera data.

◆ reactionRatesTimer

PerfTimer* SourceSystem::reactionRatesTimer
protected

Timer for time spent evaluating reaction rates.

◆ thermoTimer

PerfTimer* SourceSystem::thermoTimer
protected

Timer for time spent evaluating thermodynamic properties.

◆ jacobianTimer

PerfTimer* SourceSystem::jacobianTimer
protected

Timer for time spent evaluating and factorizing the Jacobian.

◆ strainFunction

ScalarFunction* SourceSystem::strainFunction
protected

A class that provides the strain rate and its time derivative.

◆ rateMultiplierFunction

ScalarFunction* SourceSystem::rateMultiplierFunction
protected

Provides a multiplier (optional) for the production terms.

◆ heatLoss

IntegratorCallback* SourceSystem::heatLoss
protected

Heat loss rate.

◆ nSpec

size_t SourceSystem::nSpec
protected

number of species

◆ j

int SourceSystem::j
protected

grid index for this system

◆ x

double SourceSystem::x
protected

grid position for this system

◆ rhou

double SourceSystem::rhou
protected

density of the unburned gas

◆ qDot

double SourceSystem::qDot
protected

heat release rate per unit volume [W/m^3]

◆ qLoss

double SourceSystem::qLoss
protected

heat loss to the environment [W/m^3]

◆ rho

double SourceSystem::rho
protected

density [kg/m^3]

◆ cp

double SourceSystem::cp
protected

specific heat capacity (average) [J/kg*K]

◆ cpSpec

dvec SourceSystem::cpSpec
protected

species specific heat capacity [J/mol*K]

◆ Wmx

double SourceSystem::Wmx
protected

mixture molecular weight [kg/mol]

◆ W

dvec SourceSystem::W
protected

species molecular weights [kg/kmol]

◆ hk

dvec SourceSystem::hk
protected

species enthalpies [J/kmol]

◆ quasi2d

bool SourceSystem::quasi2d
protected

Flag set to 'true' when solving a quasi-2D problem with prescribed velocity and temperature fields.

◆ vzInterp

std::shared_ptr<BilinearInterpolator> SourceSystem::vzInterp
protected

An interpolator for computing the axial (z) velocity when solving a quasi-2D problem.

◆ TInterp

std::shared_ptr<BilinearInterpolator> SourceSystem::TInterp
protected

An interpolator for computing the temperature when solving a quasi-2D problem.


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