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

Represents a system of equations used to integrate the (chemical) source term at a single point using the CVODE integrator. More...

#include <sourceSystem.h>

Inheritance diagram for SourceSystemCVODE:
[legend]

Public Member Functions

 SourceSystemCVODE ()
 
int f (const realtype t, const sdVector &y, sdVector &ydot)
 The ODE function: ydot = f(t,y)
 
int denseJacobian (const realtype t, const sdVector &y, const sdVector &ydot, sdMatrix &J)
 Calculate the Jacobian matrix: J = df/dy.
 
int fdJacobian (const realtype t, const sdVector &y, const sdVector &ydot, sdMatrix &J)
 A simpler finite difference based Jacobian.
 
void setState (double tInitial, double uu, double tt, const dvec &yy)
 Set the initial condition for integrator.
 
int integrateToTime (double tf)
 Take as many steps as needed to reach tf.
 
int integrateOneStep (double tf)
 Take one step toward tf without stepping past it.
 
double time () const
 Current integrator time, relative to tInitial.
 
void unroll_y (const sdVector &y, double t)
 fill in the current state variables from y
 
void unroll_y ()
 fill in the current state variables from the integrator state
 
void roll_y (sdVector &y) const
 fill in y with the values of the current state variables
 
void roll_ydot (sdVector &ydot) const
 fill in ydot with the values of the current time derivatives
 
std::string getStats ()
 Return a string indicating the number of internal timesteps taken.
 
void initialize (size_t nSpec)
 Resize internal arrays for a problem of the specified size (nSpec+2)
 
void setOptions (ConfigOptions &options)
 Set integrator tolerances and other parameters.
 
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.
 
void writeJacobian (std::ostream &out)
 Print the current Jacobian matrix ot the specified stream.
 
- Public Member Functions inherited from SourceSystem
 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 dUdt
 time derivative of the tangential velocity
 
double dTdt
 time derivative of the temperature
 
dvec dYdt
 time derivative of the species mass fractions
 
dvec wDot
 species net production rates [kmol/m^3*s]
 
- Public Attributes inherited from SourceSystem
double U
 tangential velocity
 
double T
 temperature
 
dvec Y
 species mass fraction
 
dvec splitConst
 Extra constant term introduced by splitting.
 

Private Attributes

std::unique_ptr< SundialsCvodeintegrator
 
SundialsContext sunContext
 

Additional Inherited Members

- Protected Attributes inherited from SourceSystem
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.
 
- Private Member Functions inherited from sdODE
virtual int f (const realtype t, const sdVector &y, sdVector &ydot)=0
 Evaluate the ODE function.
 
virtual int g (realtype t, sdVector &y, realtype *gOut)
 Optional function for CVODE to find roots of.
 
virtual int denseJacobian (const realtype t, const sdVector &y, const sdVector &ydot, sdMatrix &J)
 Optional function used to compute the dense Jacobian.
 
virtual int bandedJacobian (const realtype t, const sdVector &y, const sdVector &ydot, sdBandMatrix &J)
 Optional function used to compute a banded Jacobian.
 
virtual ~sdODE ()
 

Detailed Description

Represents a system of equations used to integrate the (chemical) source term at a single point using the CVODE integrator.

Constructor & Destructor Documentation

◆ SourceSystemCVODE()

SourceSystemCVODE::SourceSystemCVODE ( )
inline

Member Function Documentation

◆ f()

int SourceSystemCVODE::f ( const realtype  t,
const sdVector y,
sdVector ydot 
)
virtual

The ODE function: ydot = f(t,y)

Implements sdODE.

◆ denseJacobian()

int SourceSystemCVODE::denseJacobian ( const realtype  t,
const sdVector y,
const sdVector ydot,
sdMatrix J 
)
virtual

Calculate the Jacobian matrix: J = df/dy.

Reimplemented from sdODE.

◆ fdJacobian()

int SourceSystemCVODE::fdJacobian ( const realtype  t,
const sdVector y,
const sdVector ydot,
sdMatrix J 
)

A simpler finite difference based Jacobian.

◆ setState()

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

Set the initial condition for integrator.

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

Implements SourceSystem.

◆ integrateToTime()

int SourceSystemCVODE::integrateToTime ( double  tf)
virtual

Take as many steps as needed to reach tf.

tf is relative to tInitial.

Implements SourceSystem.

◆ integrateOneStep()

int SourceSystemCVODE::integrateOneStep ( double  tf)
virtual

Take one step toward tf without stepping past it.

tf is relative to tInitial.

Implements SourceSystem.

◆ time()

double SourceSystemCVODE::time ( ) const
virtual

Current integrator time, relative to tInitial.

Implements SourceSystem.

◆ unroll_y() [1/2]

void SourceSystemCVODE::unroll_y ( const sdVector y,
double  t 
)

fill in the current state variables from y

◆ unroll_y() [2/2]

void SourceSystemCVODE::unroll_y ( )
inlinevirtual

fill in the current state variables from the integrator state

Implements SourceSystem.

◆ roll_y()

void SourceSystemCVODE::roll_y ( sdVector y) const

fill in y with the values of the current state variables

◆ roll_ydot()

void SourceSystemCVODE::roll_ydot ( sdVector ydot) const

fill in ydot with the values of the current time derivatives

◆ getStats()

std::string SourceSystemCVODE::getStats ( )
virtual

Return a string indicating the number of internal timesteps taken.

Implements SourceSystem.

◆ initialize()

void SourceSystemCVODE::initialize ( size_t  nSpec)
virtual

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

Reimplemented from SourceSystem.

◆ setOptions()

void SourceSystemCVODE::setOptions ( ConfigOptions options_)
virtual

Set integrator tolerances and other parameters.

Reimplemented from SourceSystem.

◆ writeState()

void SourceSystemCVODE::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 from SourceSystem.

◆ writeJacobian()

void SourceSystemCVODE::writeJacobian ( std::ostream &  out)
virtual

Print the current Jacobian matrix ot the specified stream.

Reimplemented from SourceSystem.

Member Data Documentation

◆ dUdt

double SourceSystemCVODE::dUdt

time derivative of the tangential velocity

◆ dTdt

double SourceSystemCVODE::dTdt

time derivative of the temperature

◆ dYdt

dvec SourceSystemCVODE::dYdt

time derivative of the species mass fractions

◆ wDot

dvec SourceSystemCVODE::wDot

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

◆ integrator

std::unique_ptr<SundialsCvode> SourceSystemCVODE::integrator
private

◆ sunContext

SundialsContext SourceSystemCVODE::sunContext
private

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