Cantera 2.6.0
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ReactorNet Class Reference

A class representing a network of connected reactors. More...

#include <ReactorNet.h>

Inheritance diagram for ReactorNet:
[legend]
Collaboration diagram for ReactorNet:
[legend]

Public Member Functions

 ReactorNet (const ReactorNet &)=delete
 
ReactorNetoperator= (const ReactorNet &)=delete
 
doublereal time ()
 Current value of the simulation time. More...
 
doublereal rtol ()
 Relative tolerance. More...
 
doublereal atol ()
 Absolute integration tolerance. More...
 
doublereal rtolSensitivity () const
 Relative sensitivity tolerance. More...
 
doublereal atolSensitivity () const
 Absolute sensitivity tolerance. More...
 
void advance (doublereal time)
 Advance the state of all reactors in time. More...
 
double advance (double time, bool applylimit)
 Advance the state of all reactors in time. More...
 
double step ()
 Advance the state of all reactors in time. More...
 
void addReactor (Reactor &r)
 Add the reactor r to this reactor network. More...
 
Reactorreactor (int n)
 Return a reference to the n-th reactor in this network. More...
 
bool verbose () const
 Returns true if verbose logging output is enabled. More...
 
void setVerbose (bool v=true)
 Enable or disable verbose logging while setting up and integrating the reactor network. More...
 
Integratorintegrator ()
 Return a reference to the integrator. More...
 
void updateState (doublereal *y)
 Update the state of all the reactors in the network to correspond to the values in the solution vector y. More...
 
double sensitivity (size_t k, size_t p)
 Return the sensitivity of the k-th solution component with respect to the p-th sensitivity parameter. More...
 
double sensitivity (const std::string &component, size_t p, int reactor=0)
 Return the sensitivity of the component named component with respect to the p-th sensitivity parameter. More...
 
void evalJacobian (doublereal t, doublereal *y, doublereal *ydot, doublereal *p, Array2D *j)
 Evaluate the Jacobian matrix for the reactor network. More...
 
virtual size_t neq ()
 Number of equations. More...
 
virtual void eval (doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
 Evaluate the right-hand-side function. More...
 
virtual void getState (doublereal *y)
 Fill in the vector y with the current state of the system. More...
 
virtual void getDerivative (int k, double *dky)
 Return k-th derivative at the current time. More...
 
virtual size_t nparams ()
 Number of sensitivity parameters. More...
 
size_t globalComponentIndex (const std::string &component, size_t reactor=0)
 Return the index corresponding to the component named component in the reactor with index reactor in the global state vector for the reactor network. More...
 
std::string componentName (size_t i) const
 Return the name of the i-th component of the global state vector. More...
 
size_t registerSensitivityParameter (const std::string &name, double value, double scale)
 Used by Reactor and Wall objects to register the addition of sensitivity parameters so that the ReactorNet can keep track of the order in which sensitivity parameters are added. More...
 
const std::string & sensitivityParameterName (size_t p)
 The name of the p-th sensitivity parameter added to this ReactorNet. More...
 
void initialize ()
 Initialize the reactor network. More...
 
void reinitialize ()
 Reinitialize the integrator. More...
 
void setNeedsReinit ()
 Called to trigger integrator reinitialization before further integration. More...
 
virtual void setMaxSteps (int nmax)
 Set the maximum number of internal integration time-steps the integrator will take before reaching the next output time. More...
 
virtual int maxSteps ()
 Returns the maximum number of internal integration time-steps the integrator will take before reaching the next output time. More...
 
void setAdvanceLimits (const double *limits)
 Set absolute step size limits during advance. More...
 
bool hasAdvanceLimits ()
 Check whether ReactorNet object uses advance limits. More...
 
bool getAdvanceLimits (double *limits)
 Retrieve absolute step size limits during advance. More...
 
Methods to set up a simulation.
void setInitialTime (double time)
 Set initial time. More...
 
double maxTimeStep ()
 Get the maximum time step. More...
 
void setMaxTimeStep (double maxstep)
 Set the maximum time step. More...
 
void setMaxErrTestFails (int nmax)
 Set the maximum number of error test failures permitted by the CVODES integrator in a single time step. More...
 
void setTolerances (double rtol, double atol)
 Set the relative and absolute tolerances for the integrator. More...
 
void setSensitivityTolerances (double rtol, double atol)
 Set the relative and absolute tolerances for integrating the sensitivity equations. More...
 
- Public Member Functions inherited from FuncEval
int eval_nothrow (double t, double *y, double *ydot)
 Evaluate the right-hand side using return code to indicate status. More...
 
void suppressErrors (bool suppress)
 Enable or disable suppression of errors when calling eval() More...
 
bool suppressErrors () const
 Get current state of error suppression. More...
 
std::string getErrors () const
 Return a string containing the text of any suppressed errors. More...
 
void clearErrors ()
 Clear any previously-stored suppressed errors. More...
 

Protected Member Functions

virtual void getEstimate (double time, int k, double *yest)
 Estimate a future state based on current derivatives. More...
 
virtual int lastOrder ()
 Returns the order used for last solution step of the ODE integrator The function is intended for internal use by ReactorNet::advance and deliberately not exposed in external interfaces. More...
 

Protected Attributes

std::vector< Reactor * > m_reactors
 
std::unique_ptr< Integratorm_integ
 
doublereal m_time
 
bool m_init
 
bool m_integrator_init
 True if integrator initialization is current. More...
 
size_t m_nv
 
std::vector< size_t > m_start
 m_start[n] is the starting point in the state vector for reactor n More...
 
vector_fp m_atol
 
doublereal m_rtol
 
doublereal m_rtolsens
 
doublereal m_atols
 
doublereal m_atolsens
 
doublereal m_maxstep
 Maximum integrator internal timestep. Default of 0.0 means infinity. More...
 
int m_maxErrTestFails
 
bool m_verbose
 
std::vector< std::string > m_paramNames
 Names corresponding to each sensitivity parameter. More...
 
vector_fp m_ydot
 
vector_fp m_yest
 
vector_fp m_advancelimits
 
vector_fp m_LHS
 m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation More...
 
vector_fp m_RHS
 
bool m_checked_eval_deprecation
 
std::vector< bool > m_have_deprecated_eval
 
- Protected Attributes inherited from FuncEval
bool m_suppress_errors
 
std::vector< std::string > m_errors
 Errors occurring during function evaluations. More...
 

Additional Inherited Members

- Public Attributes inherited from FuncEval
vector_fp m_sens_params
 Values for the problem parameters for which sensitivities are computed This is the array which is perturbed and passed back as the fourth argument to eval(). More...
 
vector_fp m_paramScales
 Scaling factors for each sensitivity parameter. More...
 

Detailed Description

A class representing a network of connected reactors.

This class is used to integrate the time-dependent governing equations for a network of reactors (Reactor, ConstPressureReactor) connected by various means, for example Wall, MassFlowController, Valve, or PressureController.

Definition at line 26 of file ReactorNet.h.

Constructor & Destructor Documentation

◆ ReactorNet()

Definition at line 18 of file ReactorNet.cpp.

◆ ~ReactorNet()

~ReactorNet ( )
virtual

Definition at line 34 of file ReactorNet.cpp.

Member Function Documentation

◆ setInitialTime()

void setInitialTime ( double  time)

Set initial time.

Default = 0.0 s. Restarts integration from this time using the current mixture state as the initial condition.

Definition at line 38 of file ReactorNet.cpp.

◆ maxTimeStep()

double maxTimeStep ( )
inline

Get the maximum time step.

Definition at line 42 of file ReactorNet.h.

References ReactorNet::m_maxstep.

◆ setMaxTimeStep()

void setMaxTimeStep ( double  maxstep)

Set the maximum time step.

Definition at line 44 of file ReactorNet.cpp.

◆ setMaxErrTestFails()

void setMaxErrTestFails ( int  nmax)

Set the maximum number of error test failures permitted by the CVODES integrator in a single time step.

Definition at line 50 of file ReactorNet.cpp.

◆ setTolerances()

void setTolerances ( double  rtol,
double  atol 
)

Set the relative and absolute tolerances for the integrator.

Definition at line 56 of file ReactorNet.cpp.

◆ setSensitivityTolerances()

void setSensitivityTolerances ( double  rtol,
double  atol 
)

Set the relative and absolute tolerances for integrating the sensitivity equations.

Definition at line 67 of file ReactorNet.cpp.

◆ time()

doublereal time ( )
inline

Current value of the simulation time.

Definition at line 63 of file ReactorNet.h.

◆ rtol()

doublereal rtol ( )
inline

Relative tolerance.

Definition at line 68 of file ReactorNet.h.

◆ atol()

doublereal atol ( )
inline

Absolute integration tolerance.

Definition at line 73 of file ReactorNet.h.

◆ rtolSensitivity()

doublereal rtolSensitivity ( ) const
inline

Relative sensitivity tolerance.

Definition at line 78 of file ReactorNet.h.

◆ atolSensitivity()

doublereal atolSensitivity ( ) const
inline

Absolute sensitivity tolerance.

Definition at line 83 of file ReactorNet.h.

◆ advance() [1/2]

void advance ( doublereal  time)

Advance the state of all reactors in time.

Take as many internal timesteps as necessary to reach time.

Parameters
timeTime to advance to (s).

Definition at line 143 of file ReactorNet.cpp.

◆ advance() [2/2]

double advance ( double  time,
bool  applylimit 
)

Advance the state of all reactors in time.

Take as many internal timesteps as necessary towards time. If applylimit is true, the advance step will be automatically reduced if needed to stay within limits (set by setAdvanceLimit). Returns the time at the end of integration.

Parameters
timeTime to advance to (s).
applylimitLimit advance step (boolean).

Definition at line 155 of file ReactorNet.cpp.

◆ step()

double step ( )

Advance the state of all reactors in time.

Definition at line 210 of file ReactorNet.cpp.

◆ addReactor()

void addReactor ( Reactor r)

Add the reactor r to this reactor network.

Definition at line 247 of file ReactorNet.cpp.

References ReactorBase::setNetwork().

◆ reactor()

Reactor & reactor ( int  n)
inline

Return a reference to the n-th reactor in this network.

The reactor indices are determined by the order in which the reactors were added to the reactor network.

Definition at line 114 of file ReactorNet.h.

Referenced by ReactorNet::sensitivity().

◆ verbose()

bool verbose ( ) const
inline

Returns true if verbose logging output is enabled.

Definition at line 119 of file ReactorNet.h.

◆ setVerbose()

void setVerbose ( bool  v = true)
inline

Enable or disable verbose logging while setting up and integrating the reactor network.

Definition at line 125 of file ReactorNet.h.

References FuncEval::suppressErrors().

◆ integrator()

Integrator & integrator ( )
inline

Return a reference to the integrator.

Definition at line 131 of file ReactorNet.h.

◆ updateState()

void updateState ( doublereal *  y)

Update the state of all the reactors in the network to correspond to the values in the solution vector y.

Definition at line 348 of file ReactorNet.cpp.

References Cantera::checkFinite().

◆ sensitivity() [1/2]

double sensitivity ( size_t  k,
size_t  p 
)

Return the sensitivity of the k-th solution component with respect to the p-th sensitivity parameter.

The sensitivity coefficient \( S_{ki} \) of solution variable \( y_k \) with respect to sensitivity parameter \( p_i \) is defined as:

\[ S_{ki} = \frac{1}{y_k} \frac{\partial y_k}{\partial p_i} \]

For reaction sensitivities, the parameter is a multiplier on the forward rate constant (and implicitly on the reverse rate constant for reversible reactions) which has a nominal value of 1.0, and the sensitivity is nondimensional.

For species enthalpy sensitivities, the parameter is a perturbation to the molar enthalpy of formation, such that the dimensions of the sensitivity are kmol/J.

Definition at line 309 of file ReactorNet.cpp.

References Cantera::SmallNumber.

Referenced by ReactorNet::sensitivity().

◆ sensitivity() [2/2]

double sensitivity ( const std::string &  component,
size_t  p,
int  reactor = 0 
)
inline

Return the sensitivity of the component named component with respect to the p-th sensitivity parameter.

The sensitivity coefficient \( S_{ki} \) of solution variable \( y_k \) with respect to sensitivity parameter \( p_i \) is defined as:

\[ S_{ki} = \frac{1}{y_k} \frac{\partial y_k}{\partial p_i} \]

For reaction sensitivities, the parameter is a multiplier on the forward rate constant (and implicitly on the reverse rate constant for reversible reactions) which has a nominal value of 1.0, and the sensitivity is nondimensional.

For species enthalpy sensitivities, the parameter is a perturbation to the molar enthalpy of formation, such that the dimensions of the sensitivity are kmol/J.

Definition at line 161 of file ReactorNet.h.

References ReactorNet::globalComponentIndex(), ReactorNet::reactor(), and ReactorNet::sensitivity().

◆ evalJacobian()

void evalJacobian ( doublereal  t,
doublereal *  y,
doublereal *  ydot,
doublereal *  p,
Array2D j 
)

Evaluate the Jacobian matrix for the reactor network.

Parameters
[in]tTime at which to evaluate the Jacobian
[in]yGlobal state vector at time t
[out]ydotTime derivative of the state vector evaluated at t.
[in]psensitivity parameter vector (unused?)
[out]jJacobian matrix, size neq() by neq().

Definition at line 325 of file ReactorNet.cpp.

◆ neq()

virtual size_t neq ( )
inlinevirtual

Number of equations.

Implements FuncEval.

Definition at line 178 of file ReactorNet.h.

◆ eval()

void eval ( doublereal  t,
doublereal *  y,
doublereal *  ydot,
doublereal *  p 
)
virtual

Evaluate the right-hand-side function.

Called by the integrator.

Parameters
[in]ttime.
[in]ysolution vector, length neq()
[out]ydotrate of change of solution vector, length neq()
[in]psensitivity parameter vector, length nparams()

Implements FuncEval.

Definition at line 253 of file ReactorNet.cpp.

◆ getState()

void getState ( doublereal *  y)
virtual

Fill in the vector y with the current state of the system.

Reimplemented from FuncEval.

Definition at line 356 of file ReactorNet.cpp.

◆ getDerivative()

void getDerivative ( int  k,
double *  dky 
)
virtual

Return k-th derivative at the current time.

Definition at line 363 of file ReactorNet.cpp.

◆ nparams()

virtual size_t nparams ( )
inlinevirtual

Number of sensitivity parameters.

Reimplemented from FuncEval.

Definition at line 189 of file ReactorNet.h.

References FuncEval::m_sens_params.

◆ globalComponentIndex()

size_t globalComponentIndex ( const std::string &  component,
size_t  reactor = 0 
)

Return the index corresponding to the component named component in the reactor with index reactor in the global state vector for the reactor network.

Definition at line 399 of file ReactorNet.cpp.

Referenced by ReactorNet::sensitivity().

◆ componentName()

std::string componentName ( size_t  i) const

Return the name of the i-th component of the global state vector.

The name returned includes both the name of the reactor and the specific component, for example ‘'reactor1: CH4’`.

Definition at line 407 of file ReactorNet.cpp.

◆ registerSensitivityParameter()

size_t registerSensitivityParameter ( const std::string &  name,
double  value,
double  scale 
)

Used by Reactor and Wall objects to register the addition of sensitivity parameters so that the ReactorNet can keep track of the order in which sensitivity parameters are added.

Parameters
nameA name describing the parameter, for example the reaction string
valueThe nominal value of the parameter
scaleA scaling factor to be applied to the sensitivity coefficient
Returns
the index of this parameter in the vector of sensitivity parameters (global across all reactors)

Definition at line 419 of file ReactorNet.cpp.

References Cantera::scale().

◆ sensitivityParameterName()

const std::string & sensitivityParameterName ( size_t  p)
inline

The name of the p-th sensitivity parameter added to this ReactorNet.

Definition at line 216 of file ReactorNet.h.

References ReactorNet::m_paramNames.

◆ initialize()

void initialize ( )

Initialize the reactor network.

Called automatically the first time advance or step is called.

Definition at line 78 of file ReactorNet.cpp.

◆ reinitialize()

void reinitialize ( )

Reinitialize the integrator.

Used to solve a new problem (different initial conditions) but with the same configuration of the reactor network. Can be called manually, or automatically after calling setInitialTime or modifying a reactor's contents.

Definition at line 122 of file ReactorNet.cpp.

References Cantera::debuglog().

◆ setNeedsReinit()

void setNeedsReinit ( )
inline

Called to trigger integrator reinitialization before further integration.

Definition at line 232 of file ReactorNet.h.

References ReactorNet::m_integrator_init.

◆ setMaxSteps()

void setMaxSteps ( int  nmax)
virtual

Set the maximum number of internal integration time-steps the integrator will take before reaching the next output time.

Parameters
nmaxThe maximum number of steps, setting this value to zero disables this option.

Definition at line 133 of file ReactorNet.cpp.

◆ maxSteps()

int maxSteps ( )
virtual

Returns the maximum number of internal integration time-steps the integrator will take before reaching the next output time.

Definition at line 138 of file ReactorNet.cpp.

◆ setAdvanceLimits()

void setAdvanceLimits ( const double *  limits)

Set absolute step size limits during advance.

Definition at line 371 of file ReactorNet.cpp.

◆ hasAdvanceLimits()

bool hasAdvanceLimits ( )

Check whether ReactorNet object uses advance limits.

Definition at line 381 of file ReactorNet.cpp.

◆ getAdvanceLimits()

bool getAdvanceLimits ( double *  limits)

Retrieve absolute step size limits during advance.

Definition at line 390 of file ReactorNet.cpp.

◆ getEstimate()

void getEstimate ( double  time,
int  k,
double *  yest 
)
protectedvirtual

Estimate a future state based on current derivatives.

The function is intended for internal use by ReactorNet::advance and deliberately not exposed in external interfaces.

Definition at line 222 of file ReactorNet.cpp.

◆ lastOrder()

int lastOrder ( )
protectedvirtual

Returns the order used for last solution step of the ODE integrator The function is intended for internal use by ReactorNet::advance and deliberately not exposed in external interfaces.

Definition at line 242 of file ReactorNet.cpp.

Member Data Documentation

◆ m_reactors

std::vector<Reactor*> m_reactors
protected

Definition at line 268 of file ReactorNet.h.

◆ m_integ

std::unique_ptr<Integrator> m_integ
protected

Definition at line 269 of file ReactorNet.h.

◆ m_time

doublereal m_time
protected

Definition at line 270 of file ReactorNet.h.

◆ m_init

bool m_init
protected

Definition at line 271 of file ReactorNet.h.

◆ m_integrator_init

bool m_integrator_init
protected

True if integrator initialization is current.

Definition at line 272 of file ReactorNet.h.

Referenced by ReactorNet::setNeedsReinit().

◆ m_nv

size_t m_nv
protected

Definition at line 273 of file ReactorNet.h.

◆ m_start

std::vector<size_t> m_start
protected

m_start[n] is the starting point in the state vector for reactor n

Definition at line 276 of file ReactorNet.h.

◆ m_atol

vector_fp m_atol
protected

Definition at line 278 of file ReactorNet.h.

◆ m_rtol

doublereal m_rtol
protected

Definition at line 279 of file ReactorNet.h.

◆ m_rtolsens

doublereal m_rtolsens
protected

Definition at line 279 of file ReactorNet.h.

◆ m_atols

doublereal m_atols
protected

Definition at line 280 of file ReactorNet.h.

◆ m_atolsens

doublereal m_atolsens
protected

Definition at line 280 of file ReactorNet.h.

◆ m_maxstep

doublereal m_maxstep
protected

Maximum integrator internal timestep. Default of 0.0 means infinity.

Definition at line 283 of file ReactorNet.h.

Referenced by ReactorNet::maxTimeStep().

◆ m_maxErrTestFails

int m_maxErrTestFails
protected

Definition at line 285 of file ReactorNet.h.

◆ m_verbose

bool m_verbose
protected

Definition at line 286 of file ReactorNet.h.

◆ m_paramNames

std::vector<std::string> m_paramNames
protected

Names corresponding to each sensitivity parameter.

Definition at line 289 of file ReactorNet.h.

Referenced by ReactorNet::sensitivityParameterName().

◆ m_ydot

vector_fp m_ydot
protected

Definition at line 291 of file ReactorNet.h.

◆ m_yest

vector_fp m_yest
protected

Definition at line 292 of file ReactorNet.h.

◆ m_advancelimits

vector_fp m_advancelimits
protected

Definition at line 293 of file ReactorNet.h.

◆ m_LHS

vector_fp m_LHS
protected

m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation

Definition at line 296 of file ReactorNet.h.

◆ m_RHS

vector_fp m_RHS
protected

Definition at line 297 of file ReactorNet.h.

◆ m_checked_eval_deprecation

bool m_checked_eval_deprecation
protected
Todo:
Remove after Cantera 2.6

Definition at line 298 of file ReactorNet.h.

◆ m_have_deprecated_eval

std::vector<bool> m_have_deprecated_eval
protected
Todo:
Remove after Cantera 2.6

Definition at line 299 of file ReactorNet.h.


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