A surface in contact with a FlowReactor. More...
#include <ReactorSurface.h>
A surface in contact with a FlowReactor.
May represent the reactor wall or a catalytic surface within the reactor.
Definition at line 164 of file ReactorSurface.h.
Public Member Functions | |
| FlowReactorSurface (shared_ptr< Solution > soln, const vector< shared_ptr< ReactorBase > > &reactors, bool clone, const string &name="(none)") | |
| Create a new ReactorSurface. | |
| string | type () const override |
| String indicating the reactor model implemented. | |
| void | evalDae (double t, double *y, double *ydot, double *residual) override |
| Evaluate the reactor governing equations. | |
| void | getStateDae (double *y, double *ydot) override |
| Get the current state and derivative vector of the reactor for a DAE solver. | |
| void | getConstraints (double *constraints) override |
| Given a vector of length neq(), mark which variables should be considered algebraic constraints. | |
| double | area () const override |
| Surface area per unit length [m]. | |
| void | setArea (double A) override |
| Set the reactor surface area per unit length [m]. | |
| double | initialAtol () const |
| Get the steady state tolerances used to determine the initial state for surface coverages. | |
| void | setInitialAtol (double atol) |
| Set the steady state tolerances used to determine the initial state for surface coverages. | |
| double | initialRtol () const |
| Get the steady state tolerances used to determine the initial state for surface coverages. | |
| void | setInitialRtol (double rtol) |
| Set the steady state tolerances used to determine the initial state for surface coverages. | |
| double | initialMaxSteps () const |
| Get the steady state tolerances used to determine the initial state for surface coverages. | |
| void | setInitialMaxSteps (int max_steps) |
| Set the steady state tolerances used to determine the initial state for surface coverages. | |
| double | initialMaxErrorFailures () const |
| Get the steady state tolerances used to determine the initial state for surface coverages. | |
| void | setInitialMaxErrorFailures (int max_fails) |
| Set the steady state tolerances used to determine the initial state for surface coverages. | |
Public Member Functions inherited from ReactorSurface | |
| ReactorSurface (shared_ptr< Solution > soln, const vector< shared_ptr< ReactorBase > > &reactors, bool clone, const string &name="(none)") | |
| Create a new ReactorSurface. | |
| string | type () const override |
| String indicating the wall model implemented. | |
| double | area () const override |
| Returns the surface area [m²]. | |
| void | setArea (double a) override |
| Set the surface area [m²]. | |
| SurfPhase * | thermo () |
| Accessor for the SurfPhase object. | |
| InterfaceKinetics * | kinetics () |
| Accessor for the InterfaceKinetics object. | |
| void | addInlet (FlowDevice &inlet) override |
| Connect an inlet FlowDevice to this reactor. | |
| void | addOutlet (FlowDevice &outlet) override |
| Connect an outlet FlowDevice to this reactor. | |
| void | addWall (WallBase &w, int lr) override |
| Insert a Wall between this reactor and another reactor. | |
| void | addSurface (ReactorSurface *surf) override |
| Add a ReactorSurface object to a Reactor object. | |
| void | setCoverages (const double *cov) |
| Set the surface coverages. | |
| void | setCoverages (const Composition &cov) |
| Set the surface coverages by name. | |
| void | setCoverages (const string &cov) |
| Set the surface coverages by name. | |
| void | getCoverages (double *cov) const |
| Get the surface coverages. | |
| void | getState (double *y) override |
| Get the current state of the reactor. | |
| void | initialize (double t0=0.0) override |
| Initialize 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. | |
| void | addSensitivityReaction (size_t rxn) override |
| Add a sensitivity parameter associated with the reaction number rxn | |
| 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. | |
| 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. | |
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 void | setInitialVolume (double vol) |
| Set the initial reactor volume. | |
| virtual void | evalWalls (double t) |
| Evaluate contributions from walls connected to this reactor. | |
| 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. | |
| 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. | |
| WallBase & | wall (size_t n) |
| Return a reference to the n-th Wall connected to this reactor. | |
| 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 vector< size_t > | steadyConstraints () const |
| Get the indices of equations that are algebraic constraints when solving the steady-state problem. | |
| virtual void | addSensitivitySpeciesEnthalpy (size_t k) |
| Add a sensitivity parameter associated with the enthalpy formation of species k. | |
| 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 (double *y) |
| Reset physically or mathematically problematic values, such as negative species concentrations. | |
| virtual void | getJacobianElements (vector< Eigen::Triplet< double > > &trips) |
| Get Jacobian elements for this reactor within the full reactor network. | |
| 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 bool | preconditionerSupported () const |
| Return a false if preconditioning is not supported or true otherwise. | |
| 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 | |
| double | m_ss_rtol = 1e-7 |
| steady-state relative tolerance, used to determine initial surface coverages | |
| double | m_ss_atol = 1e-14 |
| steady-state absolute tolerance, used to determine initial surface coverages | |
| int | m_max_ss_steps = 20000 |
| maximum number of steady-state coverage integrator-steps | |
| int | m_max_ss_error_fails = 10 |
| maximum number of steady-state integrator error test failures | |
Protected Attributes inherited from ReactorSurface | |
| double | m_area = 1.0 |
| shared_ptr< SurfPhase > | m_surf |
| shared_ptr< InterfaceKinetics > | m_kinetics |
| vector< ReactorBase * > | m_reactors |
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 ReactorBase | |
| ReactorBase (const string &name="(none)") | |
| FlowReactorSurface | ( | shared_ptr< Solution > | soln, |
| const vector< shared_ptr< ReactorBase > > & | reactors, | ||
| bool | clone, | ||
| const string & | name = "(none)" |
||
| ) |
Create a new ReactorSurface.
| soln | Thermodynamic and kinetic model representing species and reactions on the surface |
| clone | Determines whether to clone soln so that the internal state of this reactor surface is independent of the original Solution (Interface) object and any Solution objects used by other reactors in the network. |
| reactors | One or more reactors whose phases participate in reactions occurring on the surface. For the purpose of rate evaluation, the temperature of the surface is set equal to the temperature of the first reactor specified. |
| name | Name used to identify the surface |
reactors and clone arguments introduced in Cantera 3.2. Definition at line 258 of file ReactorSurface.cpp.
|
inlineoverridevirtual |
String indicating the reactor model implemented.
Usually corresponds to the name of the derived class.
Reimplemented from ReactorBase.
Definition at line 173 of file ReactorSurface.h.
|
overridevirtual |
Evaluate the reactor governing equations.
Called by ReactorNet::eval.
| [in] | t | time. |
| [in] | y | solution vector, length neq() |
| [in] | ydot | rate of change of solution vector, length neq() |
| [out] | residual | residuals vector, length neq() |
Reimplemented from ReactorBase.
Definition at line 277 of file ReactorSurface.cpp.
|
overridevirtual |
Get the current state and derivative vector of the reactor for a DAE solver.
| [out] | y | state vector representing the initial state of the reactor |
| [out] | ydot | state vector representing the initial derivatives of the reactor |
Reimplemented from ReactorBase.
Definition at line 288 of file ReactorSurface.cpp.
|
overridevirtual |
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
Reimplemented from ReactorBase.
Definition at line 298 of file ReactorSurface.cpp.
|
overridevirtual |
Surface area per unit length [m].
Reimplemented from ReactorBase.
Definition at line 267 of file ReactorSurface.cpp.
|
inlineoverridevirtual |
Set the reactor surface area per unit length [m].
If the surface is the wall of the reactor, then this is the perimeter of the reactor. If the surface represents something like a catalyst monolith, this is the inverse of the surface area to volume ratio.
Reimplemented from ReactorBase.
Definition at line 191 of file ReactorSurface.h.
|
inline |
Get the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 195 of file ReactorSurface.h.
|
inline |
Set the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 201 of file ReactorSurface.h.
|
inline |
Get the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 207 of file ReactorSurface.h.
|
inline |
Set the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 213 of file ReactorSurface.h.
|
inline |
Get the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 219 of file ReactorSurface.h.
|
inline |
Set the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 225 of file ReactorSurface.h.
|
inline |
Get the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 231 of file ReactorSurface.h.
|
inline |
Set the steady state tolerances used to determine the initial state for surface coverages.
Definition at line 237 of file ReactorSurface.h.
|
protected |
steady-state relative tolerance, used to determine initial surface coverages
Definition at line 243 of file ReactorSurface.h.
|
protected |
steady-state absolute tolerance, used to determine initial surface coverages
Definition at line 245 of file ReactorSurface.h.
|
protected |
maximum number of steady-state coverage integrator-steps
Definition at line 247 of file ReactorSurface.h.
|
protected |
maximum number of steady-state integrator error test failures
Definition at line 249 of file ReactorSurface.h.