IdealGasConstPressureMoleReactor is a class for ideal gas constant-pressure reactors which use a state of moles. More...
#include <IdealGasConstPressureMoleReactor.h>
IdealGasConstPressureMoleReactor is a class for ideal gas constant-pressure reactors which use a state of moles.
Definition at line 20 of file IdealGasConstPressureMoleReactor.h.
Public Member Functions | |
| string | type () const override |
| String indicating the reactor model implemented. | |
| void | getState (double *y) override |
| Get the current state of the reactor. | |
| void | initialize (double t0=0.0) override |
| Initialize the reactor. | |
| void | eval (double t, double *LHS, double *RHS) override |
| Evaluate the reactor governing equations. | |
| void | updateState (double *y) override |
| Set the state of the reactor to correspond to the state vector y. | |
| void | getJacobianScalingFactors (double &f_species, double *f_energy) override |
| Get scaling factors for the Jacobian matrix terms proportional to \( d\dot{n}_k/dC_j \). | |
| void | getJacobianElements (vector< Eigen::Triplet< double > > &trips) override |
| Calculate an approximate Jacobian to accelerate preconditioned solvers. | |
| bool | preconditionerSupported () const override |
| Return a false if preconditioning is not supported or true otherwise. | |
| size_t | componentIndex (const string &nm) const override |
| Return the index in the solution vector for this reactor of the component named nm. | |
| string | componentName (size_t k) override |
| Return the name of the solution component with index i. | |
| double | upperBound (size_t k) const override |
| Get the upper bound on the k-th component of the local state vector. | |
| double | lowerBound (size_t k) const override |
| Get the lower bound on the k-th component of the local state vector. | |
| ConstPressureMoleReactor (shared_ptr< Solution > sol, const string &name="(none)") | |
| ConstPressureMoleReactor (shared_ptr< Solution > sol, bool clone, const string &name="(none)") | |
Public Member Functions inherited from ConstPressureMoleReactor | |
| ConstPressureMoleReactor (shared_ptr< Solution > sol, const string &name="(none)") | |
| ConstPressureMoleReactor (shared_ptr< Solution > sol, bool clone, const string &name="(none)") | |
| string | type () const override |
| String indicating the reactor model implemented. | |
| void | getState (double *y) override |
| Get the current state of the reactor. | |
| void | eval (double t, double *LHS, double *RHS) override |
| Evaluate the reactor governing equations. | |
| vector< size_t > | steadyConstraints () const override |
| Get the indices of equations that are algebraic constraints when solving the steady-state problem. | |
| void | updateState (double *y) override |
| Set the state of the reactor to correspond to the state vector y. | |
| size_t | componentIndex (const string &nm) const override |
| Return the index in the solution vector for this reactor of the component named nm. | |
| string | componentName (size_t k) override |
| Return the name of the solution component with index i. | |
| double | upperBound (size_t k) const override |
| Get the upper bound on the k-th component of the local state vector. | |
| double | lowerBound (size_t k) const override |
| Get the lower bound on the k-th component of the local state vector. | |
| void | resetBadValues (double *y) override |
| Reset physically or mathematically problematic values, such as negative species concentrations. | |
Public Member Functions inherited from MoleReactor | |
| MoleReactor (shared_ptr< Solution > sol, const string &name="(none)") | |
| MoleReactor (shared_ptr< Solution > sol, bool clone, const string &name="(none)") | |
| string | type () const override |
| String indicating the reactor model implemented. | |
| void | getState (double *y) override |
| Get the current state of the reactor. | |
| void | updateState (double *y) override |
| Set the state of the reactor to correspond to the state vector y. | |
| void | eval (double t, double *LHS, double *RHS) override |
| Evaluate the reactor governing equations. | |
| size_t | componentIndex (const string &nm) const override |
| Return the index in the solution vector for this reactor of the component named nm. | |
| string | componentName (size_t k) override |
| Return the name of the solution component with index i. | |
| double | upperBound (size_t k) const override |
| Get the upper bound on the k-th component of the local state vector. | |
| double | lowerBound (size_t k) const override |
| Get the lower bound on the k-th component of the local state vector. | |
| void | resetBadValues (double *y) override |
| Reset physically or mathematically problematic values, such as negative species concentrations. | |
Public Member Functions inherited from Reactor | |
| Reactor (shared_ptr< Solution > sol, const string &name="(none)") | |
| Reactor (shared_ptr< Solution > sol, bool clone, const string &name="(none)") | |
| string | type () const override |
| String indicating the reactor model implemented. | |
| virtual bool | isOde () const |
| Indicate whether the governing equations for this reactor type are a system of ODEs or DAEs. | |
| void | setInitialVolume (double vol) override |
| Set the initial reactor volume. | |
| void | setChemistryEnabled (bool cflag=true) override |
| Enable or disable changes in reactor composition due to chemical reactions. | |
| bool | chemistryEnabled () const override |
Returns true if changes in the reactor composition due to chemical reactions are enabled. | |
| void | setEnergyEnabled (bool eflag=true) override |
| Set the energy equation on or off. | |
| bool | energyEnabled () const override |
Returns true if solution of the energy equation is enabled. | |
| void | getState (double *y) override |
| Get the current state of the reactor. | |
| void | initialize (double t0=0.0) override |
| Initialize the reactor. | |
| void | eval (double t, double *LHS, double *RHS) override |
| Evaluate the reactor governing equations. | |
| vector< size_t > | steadyConstraints () const override |
| Get the indices of equations that are algebraic constraints when solving the steady-state problem. | |
| void | updateState (double *y) override |
| Set the state of the reactor to correspond to the state vector y. | |
| void | addSensitivityReaction (size_t rxn) override |
| Add a sensitivity parameter associated with the reaction number rxn | |
| void | addSensitivitySpeciesEnthalpy (size_t k) override |
| Add a sensitivity parameter associated with the enthalpy formation of species k. | |
| size_t | componentIndex (const string &nm) const override |
| Return the index in the solution vector for this reactor of the component named nm. | |
| string | componentName (size_t k) override |
| Return the name of the solution component with index i. | |
| double | upperBound (size_t k) const override |
| Get the upper bound on the k-th component of the local state vector. | |
| double | lowerBound (size_t k) const override |
| Get the lower bound on the k-th component of the local state vector. | |
| void | resetBadValues (double *y) override |
| Reset physically or mathematically problematic values, such as negative species concentrations. | |
| void | setAdvanceLimits (const double *limits) |
| Set absolute step size limits during advance. | |
| bool | hasAdvanceLimits () const |
| Check whether Reactor object uses advance limits. | |
| bool | getAdvanceLimits (double *limits) const |
| Retrieve absolute step size limits during advance. | |
| void | setAdvanceLimit (const string &nm, const double limit) |
| Set individual step size limit for component name nm | |
| Eigen::SparseMatrix< double > | finiteDifferenceJacobian () |
| Calculate the reactor-specific Jacobian using a finite difference method. | |
| void | setDerivativeSettings (AnyMap &settings) override |
| Use this to set the kinetics objects derivative settings. | |
| void | applySensitivity (double *params) override |
| Set reaction rate multipliers based on the sensitivity variables in params. | |
| void | resetSensitivity (double *params) override |
| Reset the reaction rate multipliers. | |
Public Member Functions inherited from ReactorBase | |
| 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 | |
| ReactorBase & | operator= (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< Solution > | phase () |
| Access the Solution object used to represent the contents of this reactor. | |
| shared_ptr< const Solution > | phase () 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. | |
| ReactorNet & | network () |
| 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, 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. | |
| virtual double | area () const |
| Returns an area associated with a reactor [m²]. | |
| virtual void | setArea (double a) |
| Set an area associated with a reactor [m²]. | |
| virtual void | addInlet (FlowDevice &inlet) |
| Connect an inlet FlowDevice to this reactor. | |
| virtual void | addOutlet (FlowDevice &outlet) |
| Connect an outlet FlowDevice to this reactor. | |
| FlowDevice & | inlet (size_t n=0) |
| Return a reference to the n-th inlet FlowDevice connected to this reactor. | |
| FlowDevice & | outlet (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. | |
| WallBase & | wall (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. | |
| ReactorSurface * | surface (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 | getStateDae (double *y, double *ydot) |
| Get the current state and derivative vector of the reactor for a DAE solver. | |
| virtual void | evalDae (double t, double *y, double *ydot, double *residual) |
| Evaluate the reactor governing equations. | |
| virtual void | getConstraints (double *constraints) |
| Given a vector of length neq(), mark which variables should be considered algebraic constraints. | |
| virtual Eigen::SparseMatrix< double > | jacobian () |
| Calculate the Jacobian of a specific reactor specialization. | |
| 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. | |
| 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 Attributes | |
| vector< double > | m_hk |
| Species molar enthalpies. | |
| double | m_TotalCp |
| Total heat capacity ( \( m c_p \)) [J/K]. | |
Protected Attributes inherited from ConstPressureMoleReactor | |
| const size_t | m_sidx = 1 |
Protected Attributes inherited from MoleReactor | |
| const size_t | m_sidx = 2 |
| const value for the species start index | |
Protected Attributes inherited from Reactor | |
| Kinetics * | m_kin = nullptr |
| Pointer to the homogeneous Kinetics object that handles the reactions. | |
| double | m_vdot = 0.0 |
| net rate of volume change from moving walls [m^3/s] | |
| double | m_Qdot = 0.0 |
| net heat transfer into the reactor, through walls [W] | |
| vector< double > | m_wdot |
| Species net molar production rates. | |
| vector< double > | m_uk |
| Species molar internal energies. | |
| vector< double > | m_sdot |
| Total production rate of bulk phase species on surfaces [kmol/s]. | |
| bool | m_chem = false |
| bool | m_energy = true |
| vector< double > | m_advancelimits |
| Advance step limit. | |
Protected Attributes inherited from ReactorBase | |
| size_t | m_nsp = 0 |
| Number of homogeneous species in the mixture. | |
| ThermoPhase * | m_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. | |
| ReactorNet * | m_net = nullptr |
| The ReactorNet that this reactor is part of. | |
| size_t | m_offset = 0 |
| Offset into global ReactorNet state vector. | |
| shared_ptr< Solution > | m_solution |
| Composite thermo/kinetics/transport handler. | |
| vector< SensitivityParameter > | m_sensParams |
Additional Inherited Members | |
Protected Member Functions inherited from MoleReactor | |
| void | getMoles (double *y) |
| Get moles of the system from mass fractions stored by thermo object. | |
| void | setMassFromMoles (double *y) |
| Set internal mass variable based on moles given. | |
Protected Member Functions inherited from Reactor | |
| void | updateSurfaceProductionRates () |
| Update m_sdot to reflect current production rates of bulk phase species due to reactions on adjacent surfaces. | |
| void | evalWalls (double t) override |
| Evaluate terms related to Walls. | |
Protected Member Functions inherited from ReactorBase | |
| ReactorBase (const string &name="(none)") | |
|
inlineoverridevirtual |
String indicating the reactor model implemented.
Usually corresponds to the name of the derived class.
Reimplemented from ConstPressureMoleReactor.
Definition at line 25 of file IdealGasConstPressureMoleReactor.h.
|
overridevirtual |
Get the current state of the reactor.
| [out] | y | state vector representing the initial state of the reactor |
Reimplemented from ConstPressureMoleReactor.
Definition at line 20 of file IdealGasConstPressureMoleReactor.cpp.
|
overridevirtual |
Initialize the reactor.
Called automatically by ReactorNet::initialize.
Reimplemented from ReactorBase.
Definition at line 30 of file IdealGasConstPressureMoleReactor.cpp.
|
overridevirtual |
Evaluate the reactor governing equations.
Called by ReactorNet::eval.
| [in] | t | time. |
| [out] | LHS | pointer to start of vector of left-hand side coefficients for governing equations, length m_nv, default values 1 |
| [out] | RHS | pointer to start of vector of right-hand side coefficients for governing equations, length m_nv, default values 0 |
Reimplemented from ConstPressureMoleReactor.
Definition at line 54 of file IdealGasConstPressureMoleReactor.cpp.
|
overridevirtual |
Set the state of the reactor to correspond to the state vector y.
Reimplemented from ConstPressureMoleReactor.
Definition at line 40 of file IdealGasConstPressureMoleReactor.cpp.
|
overridevirtual |
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.
| f_species | Scaling factor for derivatives appearing in the species equations. Equal to $1/V$. |
| f_energy | Scaling factor for each species term appearing in the energy equation. |
Reimplemented from ReactorBase.
Definition at line 179 of file IdealGasConstPressureMoleReactor.cpp.
|
overridevirtual |
Calculate an approximate Jacobian to accelerate preconditioned solvers.
Neglects derivatives with respect to mole fractions that would generate a fully-dense Jacobian. Currently also neglects terms related to interactions between reactors, for example via inlets and outlets.
Reimplemented from ReactorBase.
Definition at line 109 of file IdealGasConstPressureMoleReactor.cpp.
|
inlineoverridevirtual |
Return a false if preconditioning is not supported or true otherwise.
Reimplemented from ReactorBase.
Definition at line 45 of file IdealGasConstPressureMoleReactor.h.
|
overridevirtual |
Return the index in the solution vector for this reactor of the component named nm.
Reimplemented from ConstPressureMoleReactor.
Definition at line 188 of file IdealGasConstPressureMoleReactor.cpp.
|
overridevirtual |
Return the name of the solution component with index i.
Reimplemented from ConstPressureMoleReactor.
Definition at line 201 of file IdealGasConstPressureMoleReactor.cpp.
|
overridevirtual |
Get the upper bound on the k-th component of the local state vector.
Reimplemented from ConstPressureMoleReactor.
Definition at line 211 of file IdealGasConstPressureMoleReactor.cpp.
|
overridevirtual |
Get the lower bound on the k-th component of the local state vector.
Reimplemented from ConstPressureMoleReactor.
Definition at line 220 of file IdealGasConstPressureMoleReactor.cpp.
| ConstPressureMoleReactor | ( | shared_ptr< Solution > | sol, |
| const string & | name = "(none)" |
||
| ) |
Definition at line 23 of file ConstPressureMoleReactor.cpp.
| ConstPressureMoleReactor | ( | shared_ptr< Solution > | sol, |
| bool | clone, | ||
| const string & | name = "(none)" |
||
| ) |
Definition at line 24 of file ConstPressureMoleReactor.cpp.
|
protected |
Species molar enthalpies.
Definition at line 53 of file IdealGasConstPressureMoleReactor.h.
|
protected |
Total heat capacity ( \( m c_p \)) [J/K].
Definition at line 54 of file IdealGasConstPressureMoleReactor.h.