Cantera  4.0.0a1
Loading...
Searching...
No Matches
ReactorBase Class Reference

Base class for reactor objects. More...

#include <ReactorBase.h>

Inheritance diagram for ReactorBase:
[legend]

Detailed Description

Base class for reactor objects.

Allows using any substance model, with arbitrary inflow, outflow, heat loss/gain, surface chemistry, and volume change, whenever defined.

Definition at line 50 of file ReactorBase.h.

Public Member Functions

 ReactorBase (shared_ptr< Solution > sol, const string &name="(none)")
 Instantiate a ReactorBase object with Solution contents.
 
 ReactorBase (shared_ptr< Solution > sol, bool clone, const string &name="(none)")
 Instantiate a ReactorBase object with Solution contents.
 
 ReactorBase (const ReactorBase &)=delete
 
ReactorBaseoperator= (const ReactorBase &)=delete
 
virtual string type () const
 String indicating the reactor model implemented.
 
string name () const
 Return the name of this reactor.
 
void setName (const string &name)
 Set the name of this reactor.
 
bool setDefaultName (map< string, int > &counts)
 Set the default name of a reactor. Returns false if it was previously set.
 
shared_ptr< Solutionphase ()
 Access the Solution object used to represent the contents of this reactor.
 
shared_ptr< const Solutionphase () const
 Access the Solution object used to represent the contents of this reactor.
 
virtual bool timeIsIndependent () const
 Indicates whether the governing equations for this reactor are functions of time or a spatial variable.
 
size_t neq ()
 Number of equations (state variables) for this reactor.
 
virtual void syncState ()
 Set the state of the reactor to the associated ThermoPhase object.
 
virtual void updateConnected (bool updatePressure)
 Update state information needed by connected reactors, flow devices, and walls.
 
double residenceTime ()
 Return the residence time (s) of the contents of this reactor, based on the outlet mass flow rates and the mass of the reactor contents.
 
ReactorNetnetwork ()
 The ReactorNet that this reactor belongs to.
 
void setNetwork (ReactorNet *net)
 Set the ReactorNet that this reactor belongs to.
 
size_t offset () const
 Get the starting offset for this reactor's state variables within the global state vector of the ReactorNet.
 
void setOffset (size_t offset)
 Set the starting offset for this reactor's state variables within the global state vector of the ReactorNet.
 
size_t speciesOffset () const
 Offset of the first species in the local state vector.
 
virtual void getJacobianScalingFactors (double &f_species, span< double > f_energy)
 Get scaling factors for the Jacobian matrix terms proportional to \( d\dot{n}_k/dC_j \).
 
virtual void addSensitivityReaction (size_t rxn)
 Add a sensitivity parameter associated with the reaction number rxn
 
virtual size_t nSensParams () const
 Number of sensitivity parameters associated with this reactor.
 
Methods to set up a simulation
virtual void setInitialVolume (double vol)
 Set the initial reactor volume.
 
virtual double area () const
 Returns an area associated with a reactor [m²].
 
virtual void evalWalls (double t)
 Evaluate contributions from walls connected to this reactor.
 
virtual void setArea (double a)
 Set an area associated with a reactor [m²].
 
virtual bool chemistryEnabled () const
 Returns true if changes in the reactor composition due to chemical reactions are enabled.
 
virtual void setChemistryEnabled (bool cflag=true)
 Enable or disable changes in reactor composition due to chemical reactions.
 
virtual bool energyEnabled () const
 Returns true if solution of the energy equation is enabled.
 
virtual void setEnergyEnabled (bool eflag=true)
 Set the energy equation on or off.
 
virtual void addInlet (FlowDevice &inlet)
 Connect an inlet FlowDevice to this reactor.
 
virtual void addOutlet (FlowDevice &outlet)
 Connect an outlet FlowDevice to this reactor.
 
FlowDeviceinlet (size_t n=0)
 Return a reference to the n-th inlet FlowDevice connected to this reactor.
 
FlowDeviceoutlet (size_t n=0)
 Return a reference to the n-th outlet FlowDevice connected to this reactor.
 
size_t nInlets ()
 Return the number of inlet FlowDevice objects connected to this reactor.
 
size_t nOutlets ()
 Return the number of outlet FlowDevice objects connected to this reactor.
 
size_t nWalls ()
 Return the number of Wall objects connected to this reactor.
 
virtual void addWall (WallBase &w, int lr)
 Insert a Wall between this reactor and another reactor.
 
WallBasewall (size_t n)
 Return a reference to the n-th Wall connected to this reactor.
 
virtual void addSurface (ReactorSurface *surf)
 Add a ReactorSurface object to a Reactor object.
 
ReactorSurfacesurface (size_t n)
 Return a reference to the n-th ReactorSurface connected to this reactor.
 
virtual size_t nSurfs () const
 Return the number of surfaces in a reactor.
 
span< const double > surfaceProductionRates () const
 Production rates on surfaces.
 
virtual void initialize (double t0=0.0)
 Initialize the reactor.
 
virtual void getState (span< double > y)
 Get the current state of the reactor.
 
void setAbsoluteTolerances (span< const double > atol)
 Set absolute tolerances for this reactor's state variables.
 
bool getAbsoluteTolerances (span< double > atol) const
 Get the absolute tolerances set for this reactor's state variables.
 
void clearAbsoluteTolerances ()
 Clear user-specified absolute tolerances for this reactor.
 
bool hasUserTolerances () const
 Returns true if user-specified absolute tolerances are set for this reactor.
 
virtual void updateDefaultTolerances (span< double > atol, double baseAtol)
 Update default absolute tolerances based on this reactor's state variables.
 
virtual void getStateDae (span< double > y, span< double > ydot)
 Get the current state and derivative vector of the reactor for a DAE solver.
 
virtual void eval (double t, span< double > LHS, span< double > RHS)
 Evaluate the reactor governing equations.
 
virtual void evalDae (double t, span< const double > y, span< const double > ydot, span< double > residual)
 Evaluate the reactor governing equations.
 
virtual void evalSteady (double t, span< double > LHS, span< double > RHS)
 Evaluate the governing equations with modifications for the steady-state solver.
 
virtual void getConstraints (span< double > constraints)
 Given a vector of length neq(), mark which variables should be considered algebraic constraints.
 
virtual vector< size_t > initializeSteady ()
 Initialize the reactor before solving a steady-state problem.
 
virtual void updateState (span< const double > y)
 Set the state of the reactor to correspond to the state vector y.
 
virtual void addSensitivitySpeciesEnthalpy (size_t k)
 Add a sensitivity parameter associated with the enthalpy formation of species k.
 
virtual size_t componentIndex (const string &nm) const
 Return the index in the solution vector for this reactor of the component named nm.
 
virtual string componentName (size_t k)
 Return the name of the solution component with index i.
 
virtual double upperBound (size_t k) const
 Get the upper bound on the k-th component of the local state vector.
 
virtual double lowerBound (size_t k) const
 Get the lower bound on the k-th component of the local state vector.
 
virtual void resetBadValues (span< double > y)
 Reset physically or mathematically problematic values, such as negative species concentrations.
 
virtual void getJacobianElements (SparseTriplets &trips)
 Get Jacobian elements for this reactor within the full reactor network.
 
virtual void addPressureJacobian (SparseTriplets &trips, size_t row, double coeff, bool includeSpecies=true) const
 Add terms proportional to derivatives with respect to this reactor's pressure.
 
virtual void addTemperatureJacobian (SparseTriplets &trips, size_t row, double coeff) const
 Add terms proportional to derivatives with respect to this reactor's temperature.
 
virtual void addSpeciesMassFractionJacobian (SparseTriplets &trips, size_t row, size_t k, double coeff) const
 Add terms proportional to derivatives with respect to a species mass fraction.
 
virtual void addEnthalpyJacobian (SparseTriplets &trips, size_t row, double coeff, bool includeComposition=true) const
 Add terms proportional to derivatives of this reactor's specific enthalpy.
 
virtual Eigen::SparseMatrix< double > jacobian ()
 Calculate the Jacobian of a specific reactor specialization.
 
virtual void setDerivativeSettings (AnyMap &settings)
 Use this to set the kinetics objects derivative settings.
 
virtual void applySensitivity (span< const double > params)
 Set reaction rate multipliers based on the sensitivity variables in params.
 
virtual void resetSensitivity (span< const double > params)
 Reset the reaction rate multipliers.
 
virtual bool preconditionerSupported () const
 Return a false if preconditioning is not supported or true otherwise.
 
Solution components

The values returned are those after the last call to ReactorNet::advance or ReactorNet::step.

double volume () const
 Returns the current volume (m^3) of the reactor.
 
double density () const
 Returns the current density (kg/m^3) of the reactor's contents.
 
double temperature () const
 Returns the current temperature (K) of the reactor's contents.
 
double enthalpy_mass () const
 Returns the current enthalpy (J/kg) of the reactor's contents.
 
double pressure () const
 Returns the current pressure (Pa) of the reactor.
 
double mass () const
 Returns the mass (kg) of the reactor's contents.
 
span< const double > massFractions () const
 Return the vector of species mass fractions.
 
double massFraction (size_t k) const
 Return the mass fraction of the k-th species.
 

Protected Member Functions

 ReactorBase (const string &name="(none)")
 
bool isJacobianLocalRow (size_t row) const
 Check whether a global Jacobian row belongs to this reactor.
 

Protected Attributes

size_t m_nsp = 0
 Number of homogeneous species in the mixture.
 
ThermoPhasem_thermo = nullptr
 
double m_vol = 0.0
 Current volume of the reactor [m^3].
 
double m_mass = 0.0
 Current mass of the reactor [kg].
 
double m_enthalpy = 0.0
 Current specific enthalpy of the reactor [J/kg].
 
double m_pressure = 0.0
 Current pressure in the reactor [Pa].
 
vector< double > m_state
 
vector< FlowDevice * > m_inlet
 
vector< FlowDevice * > m_outlet
 
vector< WallBase * > m_wall
 
vector< ReactorSurface * > m_surfaces
 
vector< double > m_sdot
 species production rates on surfaces
 
vector< int > m_lr
 Vector of length nWalls(), indicating whether this reactor is on the left (0) or right (1) of each wall.
 
string m_name
 Reactor name.
 
bool m_defaultNameSet = false
 true if default name has been previously set.
 
size_t m_nv = 0
 Number of state variables for this reactor.
 
ReactorNetm_net = nullptr
 The ReactorNet that this reactor is part of.
 
size_t m_offset = 0
 Offset into global ReactorNet state vector.
 
vector< double > m_userAtol
 User-specified absolute tolerances for this reactor's state variables.
 
shared_ptr< Solutionm_solution
 Composite thermo/kinetics/transport handler.
 
vector< SensitivityParameterm_sensParams
 

Constructor & Destructor Documentation

◆ ReactorBase() [1/3]

ReactorBase ( shared_ptr< Solution sol,
const string &  name = "(none)" 
)

Instantiate a ReactorBase object with Solution contents.

Parameters
solSolution object to be set.
nameName of the reactor.
Since
New in Cantera 3.1.

Definition at line 23 of file ReactorBase.cpp.

◆ ReactorBase() [2/3]

ReactorBase ( shared_ptr< Solution sol,
bool  clone,
const string &  name = "(none)" 
)

Instantiate a ReactorBase object with Solution contents.

Parameters
solSolution object representing the contents of this reactor
cloneDetermines whether to clone sol so that the internal state of this reactor is independent of the original Solution object and any Solution objects used by other reactors in the network.
nameName of the reactor.
Since
Added the clone argument in Cantera 3.2. If not specified, the default behavior in Cantera 3.2 is not to clone the Solution object. This will change after Cantera 3.2 to default to true.

Definition at line 28 of file ReactorBase.cpp.

◆ ~ReactorBase()

~ReactorBase ( )
virtual

Definition at line 47 of file ReactorBase.cpp.

◆ ReactorBase() [3/3]

ReactorBase ( const string &  name = "(none)")
explicitprotected

Definition at line 19 of file ReactorBase.cpp.

Member Function Documentation

◆ type()

virtual string type ( ) const
inlinevirtual

String indicating the reactor model implemented.

Usually corresponds to the name of the derived class.

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, MoleReactor, Reactor, ReactorSurface, MoleReactorSurface, FlowReactorSurface, and Reservoir.

Definition at line 76 of file ReactorBase.h.

◆ name()

string name ( ) const
inline

Return the name of this reactor.

Definition at line 81 of file ReactorBase.h.

◆ setName()

void setName ( const string &  name)
inline

Set the name of this reactor.

Definition at line 86 of file ReactorBase.h.

◆ setDefaultName()

bool setDefaultName ( map< string, int > &  counts)

Set the default name of a reactor. Returns false if it was previously set.

Definition at line 54 of file ReactorBase.cpp.

◆ phase() [1/2]

shared_ptr< Solution > phase ( )
inline

Access the Solution object used to represent the contents of this reactor.

Since
New in Cantera 3.2

Definition at line 95 of file ReactorBase.h.

◆ phase() [2/2]

shared_ptr< const Solution > phase ( ) const
inline

Access the Solution object used to represent the contents of this reactor.

Since
New in Cantera 3.2

Definition at line 99 of file ReactorBase.h.

◆ timeIsIndependent()

virtual bool timeIsIndependent ( ) const
inlinevirtual

Indicates whether the governing equations for this reactor are functions of time or a spatial variable.

All reactors in a network must have the same value.

Reimplemented in FlowReactor, and FlowReactorSurface.

Definition at line 103 of file ReactorBase.h.

◆ neq()

size_t neq ( )
inline

Number of equations (state variables) for this reactor.

Definition at line 108 of file ReactorBase.h.

◆ setInitialVolume()

virtual void setInitialVolume ( double  vol)
inlinevirtual

Set the initial reactor volume.

Reimplemented in Reactor.

Definition at line 116 of file ReactorBase.h.

◆ area()

virtual double area ( ) const
inlinevirtual

Returns an area associated with a reactor [m²].

Examples: surface area of ReactorSurface or cross section area of FlowReactor.

Reimplemented in FlowReactor, ReactorSurface, and FlowReactorSurface.

Definition at line 123 of file ReactorBase.h.

◆ evalWalls()

virtual void evalWalls ( double  t)
inlinevirtual

Evaluate contributions from walls connected to this reactor.

Reimplemented in Reactor.

Definition at line 129 of file ReactorBase.h.

◆ setArea()

virtual void setArea ( double  a)
inlinevirtual

Set an area associated with a reactor [m²].

Examples: surface area of ReactorSurface or cross section area of FlowReactor.

Reimplemented in ReactorSurface, FlowReactorSurface, and FlowReactor.

Definition at line 136 of file ReactorBase.h.

◆ chemistryEnabled()

virtual bool chemistryEnabled ( ) const
inlinevirtual

Returns true if changes in the reactor composition due to chemical reactions are enabled.

Since
New in Cantera 3.2.

Reimplemented in Reactor.

Definition at line 144 of file ReactorBase.h.

◆ setChemistryEnabled()

virtual void setChemistryEnabled ( bool  cflag = true)
inlinevirtual

Enable or disable changes in reactor composition due to chemical reactions.

Since
New in Cantera 3.2.

Reimplemented in Reactor.

Definition at line 151 of file ReactorBase.h.

◆ energyEnabled()

virtual bool energyEnabled ( ) const
inlinevirtual

Returns true if solution of the energy equation is enabled.

Since
New in Cantera 3.2.

Reimplemented in Reactor.

Definition at line 158 of file ReactorBase.h.

◆ setEnergyEnabled()

virtual void setEnergyEnabled ( bool  eflag = true)
inlinevirtual

Set the energy equation on or off.

Since
New in Cantera 3.2.

Reimplemented in Reactor.

Definition at line 165 of file ReactorBase.h.

◆ addInlet()

void addInlet ( FlowDevice inlet)
virtual

Connect an inlet FlowDevice to this reactor.

Reimplemented in ReactorSurface.

Definition at line 67 of file ReactorBase.cpp.

◆ addOutlet()

void addOutlet ( FlowDevice outlet)
virtual

Connect an outlet FlowDevice to this reactor.

Reimplemented in ReactorSurface.

Definition at line 77 of file ReactorBase.cpp.

◆ inlet()

FlowDevice & inlet ( size_t  n = 0)

Return a reference to the n-th inlet FlowDevice connected to this reactor.

Definition at line 252 of file ReactorBase.cpp.

◆ outlet()

FlowDevice & outlet ( size_t  n = 0)

Return a reference to the n-th outlet FlowDevice connected to this reactor.

Definition at line 256 of file ReactorBase.cpp.

◆ nInlets()

size_t nInlets ( )
inline

Return the number of inlet FlowDevice objects connected to this reactor.

Definition at line 183 of file ReactorBase.h.

◆ nOutlets()

size_t nOutlets ( )
inline

Return the number of outlet FlowDevice objects connected to this reactor.

Definition at line 188 of file ReactorBase.h.

◆ nWalls()

size_t nWalls ( )
inline

Return the number of Wall objects connected to this reactor.

Definition at line 193 of file ReactorBase.h.

◆ addWall()

void addWall ( WallBase w,
int  lr 
)
virtual

Insert a Wall between this reactor and another reactor.

lr = 0 if this reactor is to the left of the wall and lr = 1 if this reactor is to the right of the wall. This method is called automatically for both the left and right reactors by WallBase::install.

Reimplemented in ReactorSurface.

Definition at line 87 of file ReactorBase.cpp.

◆ wall()

WallBase & wall ( size_t  n)

Return a reference to the n-th Wall connected to this reactor.

Definition at line 102 of file ReactorBase.cpp.

◆ addSurface()

void addSurface ( ReactorSurface surf)
virtual

Add a ReactorSurface object to a Reactor object.

Attention
This method should generally not be called directly by users. Reactor and ReactorSurface objects should be connected by providing adjacent reactors to the newReactorSurface factory function.

Reimplemented in ReactorSurface.

Definition at line 107 of file ReactorBase.cpp.

◆ surface()

ReactorSurface * surface ( size_t  n)

Return a reference to the n-th ReactorSurface connected to this reactor.

Definition at line 119 of file ReactorBase.cpp.

◆ nSurfs()

virtual size_t nSurfs ( ) const
inlinevirtual

Return the number of surfaces in a reactor.

Definition at line 218 of file ReactorBase.h.

◆ surfaceProductionRates()

span< const double > surfaceProductionRates ( ) const
inline

Production rates on surfaces.

For bulk reactors, this contains the total production rates [kmol/s] of bulk phase species due to reactions on all adjacent species.

For surfaces, this contains the production rates [kmol/m²/s] of species on the surface and all adjacent phases, in the order defined by the InterfaceKinetics object.

Definition at line 230 of file ReactorBase.h.

◆ initialize()

virtual void initialize ( double  t0 = 0.0)
inlinevirtual

◆ getState()

virtual void getState ( span< double >  y)
inlinevirtual

Get the current state of the reactor.

Parameters
[out]ystate vector representing the initial state of the reactor

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, MoleReactor, Reactor, ReactorSurface, and MoleReactorSurface.

Definition at line 245 of file ReactorBase.h.

◆ setAbsoluteTolerances()

void setAbsoluteTolerances ( span< const double >  atol)

Set absolute tolerances for this reactor's state variables.

Parameters
atolUser-specified absolute tolerances. The length must be equal to neq(), and entries are ordered according to componentName().

ReactorNet uses these values for this reactor's portion of the state vector instead of the network-level scalar absolute tolerance or reactor-specific default tolerances.

Since
New in Cantera 4.0.

Definition at line 146 of file ReactorBase.cpp.

◆ getAbsoluteTolerances()

bool getAbsoluteTolerances ( span< double >  atol) const

Get the absolute tolerances set for this reactor's state variables.

Parameters
[out]atolUser-specified absolute tolerances. The length must be equal to neq().
Returns
true if user-specified absolute tolerances are set for this reactor; false otherwise.
Since
New in Cantera 4.0.

Definition at line 158 of file ReactorBase.cpp.

◆ clearAbsoluteTolerances()

void clearAbsoluteTolerances ( )

Clear user-specified absolute tolerances for this reactor.

Since
New in Cantera 4.0.

Definition at line 171 of file ReactorBase.cpp.

◆ hasUserTolerances()

bool hasUserTolerances ( ) const
inline

Returns true if user-specified absolute tolerances are set for this reactor.

Since
New in Cantera 4.0.

Definition at line 276 of file ReactorBase.h.

◆ updateDefaultTolerances()

virtual void updateDefaultTolerances ( span< double >  atol,
double  baseAtol 
)
inlinevirtual

Update default absolute tolerances based on this reactor's state variables.

Parameters
[in,out]atolMutable slice of ReactorNet's absolute tolerance vector for this reactor. The length is neq(), and entries are ordered according to componentName(). On entry, all values are set to baseAtol.
baseAtolThe network-level scalar default absolute tolerance.

ReactorNet calls this hook only when neither this reactor nor the network has user-specified absolute tolerances.

Since
New in Cantera 4.0.

Reimplemented in MoleReactor, and MoleReactorSurface.

Definition at line 290 of file ReactorBase.h.

◆ getStateDae()

virtual void getStateDae ( span< double >  y,
span< double >  ydot 
)
inlinevirtual

Get the current state and derivative vector of the reactor for a DAE solver.

Parameters
[out]ystate vector representing the initial state of the reactor
[out]ydotstate vector representing the initial derivatives of the reactor

Reimplemented in FlowReactor, and FlowReactorSurface.

Definition at line 298 of file ReactorBase.h.

◆ eval()

virtual void eval ( double  t,
span< double >  LHS,
span< double >  RHS 
)
inlinevirtual

Evaluate the reactor governing equations.

Called by ReactorNet::eval.

Parameters
[in]ttime.
[out]LHSpointer to start of vector of left-hand side coefficients for governing equations, length m_nv, default values 1
[out]RHSpointer to start of vector of right-hand side coefficients for governing equations, length m_nv, default values 0

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, MoleReactor, Reactor, ReactorSurface, and MoleReactorSurface.

Definition at line 308 of file ReactorBase.h.

◆ evalDae()

virtual void evalDae ( double  t,
span< const double >  y,
span< const double >  ydot,
span< double >  residual 
)
inlinevirtual

Evaluate the reactor governing equations.

Called by ReactorNet::eval.

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

Reimplemented in FlowReactor, and FlowReactorSurface.

Definition at line 319 of file ReactorBase.h.

◆ evalSteady()

virtual void evalSteady ( double  t,
span< double >  LHS,
span< double >  RHS 
)
inlinevirtual

Evaluate the governing equations with modifications for the steady-state solver.

This method calls the standard eval() method then modifies elements of RHS that correspond to algebraic constraints.

Since
New in Cantera 4.0.

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, Reactor, ReactorSurface, and MoleReactorSurface.

Definition at line 330 of file ReactorBase.h.

◆ getConstraints()

virtual void getConstraints ( span< double >  constraints)
inlinevirtual

Given a vector of length neq(), mark which variables should be considered algebraic constraints.

Reimplemented in FlowReactor, and FlowReactorSurface.

Definition at line 337 of file ReactorBase.h.

◆ initializeSteady()

virtual vector< size_t > initializeSteady ( )
inlinevirtual

Initialize the reactor before solving a steady-state problem.

This method is responsible for storing the initial value for any algebraic constraints and returning the indices of those constraints.

Returns
Indices of equations that are algebraic constraints when solving the steady-state problem.
Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.
Since
New in Cantera 3.2. Renamed from steadyConstraints in Cantera 4.0.

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, Reactor, and ReactorSurface.

Definition at line 352 of file ReactorBase.h.

◆ updateState()

virtual void updateState ( span< const double >  y)
inlinevirtual

◆ addSensitivitySpeciesEnthalpy()

virtual void addSensitivitySpeciesEnthalpy ( size_t  k)
inlinevirtual

Add a sensitivity parameter associated with the enthalpy formation of species k.

Reimplemented in Reactor.

Definition at line 363 of file ReactorBase.h.

◆ componentIndex()

virtual size_t componentIndex ( const string &  nm) const
inlinevirtual

Return the index in the solution vector for this reactor of the component named nm.

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, FlowReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, MoleReactor, Reactor, and ReactorSurface.

Definition at line 369 of file ReactorBase.h.

◆ componentName()

virtual string componentName ( size_t  k)
inlinevirtual

◆ upperBound()

virtual double upperBound ( size_t  k) const
inlinevirtual

Get the upper bound on the k-th component of the local state vector.

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, MoleReactor, Reactor, ReactorSurface, and MoleReactorSurface.

Definition at line 380 of file ReactorBase.h.

◆ lowerBound()

virtual double lowerBound ( size_t  k) const
inlinevirtual

Get the lower bound on the k-th component of the local state vector.

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, IdealGasConstPressureMoleReactor, IdealGasConstPressureReactor, IdealGasMoleReactor, IdealGasReactor, MoleReactor, Reactor, ReactorSurface, and MoleReactorSurface.

Definition at line 385 of file ReactorBase.h.

◆ resetBadValues()

virtual void resetBadValues ( span< double >  y)
inlinevirtual

Reset physically or mathematically problematic values, such as negative species concentrations.

Parameters
[in,out]ycurrent state vector, to be updated; length neq()

Reimplemented in ConstPressureMoleReactor, ConstPressureReactor, MoleReactor, Reactor, ReactorSurface, and MoleReactorSurface.

Definition at line 392 of file ReactorBase.h.

◆ getJacobianElements()

virtual void getJacobianElements ( SparseTriplets &  trips)
inlinevirtual

Get Jacobian elements for this reactor within the full reactor network.

Indices within trips are global indices within the full reactor network. The reactor is responsible for providing all elements of the Jacobian in the rows corresponding to its state variables, that is, all derivatives of its state variables with respect to all state variables in the network.

Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.

Reimplemented in IdealGasConstPressureMoleReactor, IdealGasMoleReactor, and MoleReactorSurface.

Definition at line 405 of file ReactorBase.h.

◆ addPressureJacobian()

virtual void addPressureJacobian ( SparseTriplets &  trips,
size_t  row,
double  coeff,
bool  includeSpecies = true 
) const
inlinevirtual

Add terms proportional to derivatives with respect to this reactor's pressure.

This method is used by connectors to apply the chain rule. For example, a valve may know a row contribution as coeff * d(mdot)/d(P_in - P_out), while the reactor knows how its pressure depends on its own state variables. Calling this method appends all entries coeff * dP/dy_j for this reactor's state variables y_j.

Parameters
[in,out]tripsSparse Jacobian entries. Implementations append entries using global row and column indices in the reactor network.
rowGlobal row index receiving these chain-rule terms. This may be a row belonging to another reactor when a connector creates cross-reactor coupling.
coeffMultiplicative factor applied to dP/dy_j.
includeSpeciesInclude pressure derivatives with respect to species state variables. These terms may be dense for mole reactors and are controlled by preconditioner sparsity settings.
Since
New in Cantera 4.0.

Reimplemented in IdealGasMoleReactor.

Definition at line 425 of file ReactorBase.h.

◆ addTemperatureJacobian()

virtual void addTemperatureJacobian ( SparseTriplets &  trips,
size_t  row,
double  coeff 
) const
inlinevirtual

Add terms proportional to derivatives with respect to this reactor's temperature.

This method is used by connectors to apply the chain rule. Calling this method appends all entries coeff * dT/dy_j for this reactor's state variables y_j.

Parameters
[in,out]tripsSparse Jacobian entries. Implementations append entries using global row and column indices in the reactor network.
rowGlobal row index receiving these chain-rule terms. This may be a row belonging to another reactor when a connector creates cross-reactor coupling.
coeffMultiplicative factor applied to dT/dy_j.
Since
New in Cantera 4.0.

Reimplemented in IdealGasConstPressureMoleReactor, and IdealGasMoleReactor.

Definition at line 441 of file ReactorBase.h.

◆ addSpeciesMassFractionJacobian()

virtual void addSpeciesMassFractionJacobian ( SparseTriplets &  trips,
size_t  row,
size_t  k,
double  coeff 
) const
inlinevirtual

Add terms proportional to derivatives with respect to a species mass fraction.

Flow devices transport species according to upstream mass fractions. This method lets a connector request chain-rule entries for coeff * dY_k/dy_j while the reactor implementation handles how Y_k depends on its native state variables, such as species moles.

Parameters
[in,out]tripsSparse Jacobian entries. Implementations append entries using global row and column indices in the reactor network.
rowGlobal row index receiving these chain-rule terms. This may be a row belonging to another reactor when a connector creates cross-reactor coupling.
kSpecies index in this reactor's phase.
coeffMultiplicative factor applied to dY_k/dy_j.
Since
New in Cantera 4.0.

Reimplemented in IdealGasConstPressureMoleReactor, and IdealGasMoleReactor.

Definition at line 459 of file ReactorBase.h.

◆ addEnthalpyJacobian()

virtual void addEnthalpyJacobian ( SparseTriplets &  trips,
size_t  row,
double  coeff,
bool  includeComposition = true 
) const
inlinevirtual

Add terms proportional to derivatives of this reactor's specific enthalpy.

This method is used by flow device connectors to apply the chain rule for inlet enthalpy contributions. It appends entries coeff * dh/dy_j for this reactor's state variables y_j. The temperature contribution is coeff * cp_mass; the optional composition contribution adds one entry per species.

Parameters
[in,out]tripsSparse Jacobian entries. Implementations append entries using global row and column indices in the reactor network.
rowGlobal row index receiving these chain-rule terms. This may be a row belonging to another reactor when the inlet energy source is cross-reactor.
coeffMultiplicative factor applied to dh/dy_j.
includeCompositionInclude derivatives of specific enthalpy with respect to species state variables. These terms may be dense and are controlled by preconditioner sparsity settings.
Since
New in Cantera 4.0.

Reimplemented in IdealGasConstPressureMoleReactor, and IdealGasMoleReactor.

Definition at line 478 of file ReactorBase.h.

◆ setDerivativeSettings()

virtual void setDerivativeSettings ( AnyMap settings)
inlinevirtual

Use this to set the kinetics objects derivative settings.

Reimplemented in Reactor.

Definition at line 492 of file ReactorBase.h.

◆ applySensitivity()

virtual void applySensitivity ( span< const double >  params)
inlinevirtual

Set reaction rate multipliers based on the sensitivity variables in params.

Reimplemented in Reactor, and ReactorSurface.

Definition at line 498 of file ReactorBase.h.

◆ resetSensitivity()

virtual void resetSensitivity ( span< const double >  params)
inlinevirtual

Reset the reaction rate multipliers.

Reimplemented in Reactor, and ReactorSurface.

Definition at line 503 of file ReactorBase.h.

◆ preconditionerSupported()

virtual bool preconditionerSupported ( ) const
inlinevirtual

Return a false if preconditioning is not supported or true otherwise.

Warning
This method is an experimental part of the Cantera API and may be changed or removed without notice.
Since
New in Cantera 3.0

Reimplemented in IdealGasConstPressureMoleReactor, and IdealGasMoleReactor.

Definition at line 514 of file ReactorBase.h.

◆ syncState()

void syncState ( )
virtual

Set the state of the reactor to the associated ThermoPhase object.

This method will trigger integrator reinitialization.

Deprecated:
To be removed after Cantera 4.0. Use ReactorNet::reinitialize to indicate a change in state that requires integrator reinitialization.

Definition at line 136 of file ReactorBase.cpp.

◆ updateConnected()

void updateConnected ( bool  updatePressure)
virtual

Update state information needed by connected reactors, flow devices, and walls.

Called from updateState() for normal reactor types, and from ReactorNet::updateState for Reservoir.

Parameters
updatePressureIndicates whether to update m_pressure. Should true for reactors where the pressure is a dependent property, calculated from the state, and false when the pressure is constant or an independent variable.

Definition at line 179 of file ReactorBase.cpp.

◆ residenceTime()

double residenceTime ( )

Return the residence time (s) of the contents of this reactor, based on the outlet mass flow rates and the mass of the reactor contents.

Definition at line 223 of file ReactorBase.cpp.

◆ volume()

double volume ( ) const
inline

Returns the current volume (m^3) of the reactor.

Definition at line 546 of file ReactorBase.h.

◆ density()

double density ( ) const

Returns the current density (kg/m^3) of the reactor's contents.

Definition at line 232 of file ReactorBase.cpp.

◆ temperature()

double temperature ( ) const

Returns the current temperature (K) of the reactor's contents.

Definition at line 237 of file ReactorBase.cpp.

◆ enthalpy_mass()

double enthalpy_mass ( ) const
inline

Returns the current enthalpy (J/kg) of the reactor's contents.

Definition at line 557 of file ReactorBase.h.

◆ pressure()

double pressure ( ) const
inline

Returns the current pressure (Pa) of the reactor.

Definition at line 562 of file ReactorBase.h.

◆ mass()

double mass ( ) const
inline

Returns the mass (kg) of the reactor's contents.

Definition at line 567 of file ReactorBase.h.

◆ massFractions()

span< const double > massFractions ( ) const

Return the vector of species mass fractions.

Definition at line 242 of file ReactorBase.cpp.

◆ massFraction()

double massFraction ( size_t  k) const

Return the mass fraction of the k-th species.

Definition at line 247 of file ReactorBase.cpp.

◆ network()

ReactorNet & network ( )

The ReactorNet that this reactor belongs to.

Definition at line 204 of file ReactorBase.cpp.

◆ setNetwork()

void setNetwork ( ReactorNet net)

Set the ReactorNet that this reactor belongs to.

Definition at line 214 of file ReactorBase.cpp.

◆ offset()

size_t offset ( ) const
inline

Get the starting offset for this reactor's state variables within the global state vector of the ReactorNet.

Definition at line 587 of file ReactorBase.h.

◆ setOffset()

void setOffset ( size_t  offset)
inline

Set the starting offset for this reactor's state variables within the global state vector of the ReactorNet.

Definition at line 591 of file ReactorBase.h.

◆ speciesOffset()

size_t speciesOffset ( ) const
inline

Offset of the first species in the local state vector.

Definition at line 594 of file ReactorBase.h.

◆ getJacobianScalingFactors()

virtual void getJacobianScalingFactors ( double &  f_species,
span< double >  f_energy 
)
inlinevirtual

Get scaling factors for the Jacobian matrix terms proportional to \( d\dot{n}_k/dC_j \).

Used to determine contribution of surface phases to the Jacobian.

Parameters
f_speciesScaling factor for derivatives appearing in the species equations. Equal to $1/V$.
f_energyScaling factor for each species term appearing in the energy equation.

Reimplemented in IdealGasConstPressureMoleReactor, and IdealGasMoleReactor.

Definition at line 607 of file ReactorBase.h.

◆ addSensitivityReaction()

virtual void addSensitivityReaction ( size_t  rxn)
inlinevirtual

Add a sensitivity parameter associated with the reaction number rxn

Reimplemented in Reactor, and ReactorSurface.

Definition at line 612 of file ReactorBase.h.

◆ nSensParams()

virtual size_t nSensParams ( ) const
inlinevirtual

Number of sensitivity parameters associated with this reactor.

Definition at line 617 of file ReactorBase.h.

◆ isJacobianLocalRow()

bool isJacobianLocalRow ( size_t  row) const
inlineprotected

Check whether a global Jacobian row belongs to this reactor.

Used by Jacobian helper methods when local rows are already covered by the reactor's own block approximation, while cross-reactor rows need explicit connector entries.

Definition at line 629 of file ReactorBase.h.

Member Data Documentation

◆ m_nsp

size_t m_nsp = 0
protected

Number of homogeneous species in the mixture.

Definition at line 634 of file ReactorBase.h.

◆ m_thermo

ThermoPhase* m_thermo = nullptr
protected

Definition at line 636 of file ReactorBase.h.

◆ m_vol

double m_vol = 0.0
protected

Current volume of the reactor [m^3].

Definition at line 637 of file ReactorBase.h.

◆ m_mass

double m_mass = 0.0
protected

Current mass of the reactor [kg].

Definition at line 638 of file ReactorBase.h.

◆ m_enthalpy

double m_enthalpy = 0.0
protected

Current specific enthalpy of the reactor [J/kg].

Definition at line 639 of file ReactorBase.h.

◆ m_pressure

double m_pressure = 0.0
protected

Current pressure in the reactor [Pa].

Definition at line 640 of file ReactorBase.h.

◆ m_state

vector<double> m_state
protected

Definition at line 641 of file ReactorBase.h.

◆ m_inlet

vector<FlowDevice*> m_inlet
protected

Definition at line 642 of file ReactorBase.h.

◆ m_outlet

vector<FlowDevice*> m_outlet
protected

Definition at line 642 of file ReactorBase.h.

◆ m_wall

vector<WallBase*> m_wall
protected

Definition at line 644 of file ReactorBase.h.

◆ m_surfaces

vector<ReactorSurface*> m_surfaces
protected

Definition at line 645 of file ReactorBase.h.

◆ m_sdot

vector<double> m_sdot
protected

species production rates on surfaces

Definition at line 646 of file ReactorBase.h.

◆ m_lr

vector<int> m_lr
protected

Vector of length nWalls(), indicating whether this reactor is on the left (0) or right (1) of each wall.

Definition at line 650 of file ReactorBase.h.

◆ m_name

string m_name
protected

Reactor name.

Definition at line 651 of file ReactorBase.h.

◆ m_defaultNameSet

bool m_defaultNameSet = false
protected

true if default name has been previously set.

Definition at line 652 of file ReactorBase.h.

◆ m_nv

size_t m_nv = 0
protected

Number of state variables for this reactor.

Definition at line 653 of file ReactorBase.h.

◆ m_net

ReactorNet* m_net = nullptr
protected

The ReactorNet that this reactor is part of.

Definition at line 655 of file ReactorBase.h.

◆ m_offset

size_t m_offset = 0
protected

Offset into global ReactorNet state vector.

Definition at line 656 of file ReactorBase.h.

◆ m_userAtol

vector<double> m_userAtol
protected

User-specified absolute tolerances for this reactor's state variables.

Definition at line 659 of file ReactorBase.h.

◆ m_solution

shared_ptr<Solution> m_solution
protected

Composite thermo/kinetics/transport handler.

Definition at line 662 of file ReactorBase.h.

◆ m_sensParams

vector<SensitivityParameter> m_sensParams
protected

Definition at line 665 of file ReactorBase.h.


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