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

A surface where reactions can occur that is in contact with the bulk fluid of a Reactor. More...

#include <ReactorSurface.h>

Inheritance diagram for ReactorSurface:
[legend]

Detailed Description

A surface where reactions can occur that is in contact with the bulk fluid of a Reactor.

Definition at line 20 of file ReactorSurface.h.

Public Member Functions

 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²].
 
SurfPhasethermo ()
 Accessor for the SurfPhase object.
 
InterfaceKineticskinetics ()
 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.
 
size_t nAdjacent () const
 Get the number of Reactor and Reservoir objects adjacent to this surface.
 
shared_ptr< ReactorBaseadjacent (size_t n)
 Access an adjacent Reactor or Reservoir.
 
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.
 
vector< size_t > initializeSteady () override
 Initialize the reactor before solving a steady-state problem.
 
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 evalSteady (double t, double *LHS, double *RHS) override
 Evaluate the governing equations with modifications for the steady-state solver.
 
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.
 
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 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
 
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, 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.
 
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.
 
WallBasewall (size_t n)
 Return a reference to the n-th Wall connected to this reactor.
 
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 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 void addSensitivitySpeciesEnthalpy (size_t k)
 Add a sensitivity parameter associated with the enthalpy formation of species k.
 
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_area = 1.0
 
shared_ptr< SurfPhasem_surf
 
shared_ptr< InterfaceKineticsm_kinetics
 
vector< shared_ptr< ReactorBase > > m_reactors
 
- Protected Attributes inherited from ReactorBase
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.
 
shared_ptr< Solutionm_solution
 Composite thermo/kinetics/transport handler.
 
vector< SensitivityParameterm_sensParams
 

Additional Inherited Members

- Protected Member Functions inherited from ReactorBase
 ReactorBase (const string &name="(none)")
 

Constructor & Destructor Documentation

◆ ReactorSurface()

ReactorSurface ( shared_ptr< Solution soln,
const vector< shared_ptr< ReactorBase > > &  reactors,
bool  clone,
const string &  name = "(none)" 
)

Create a new ReactorSurface.

Parameters
solnThermodynamic and kinetic model representing species and reactions on the surface
cloneDetermines 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.
reactorsOne 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.
nameName used to identify the surface
Since
Constructor signature including reactors and clone arguments introduced in Cantera 3.2.

Definition at line 16 of file ReactorSurface.cpp.

Member Function Documentation

◆ type()

string type ( ) const
inlineoverridevirtual

String indicating the wall model implemented.

Reimplemented from ReactorBase.

Definition at line 42 of file ReactorSurface.h.

◆ area()

double area ( ) const
overridevirtual

Returns the surface area [m²].

Reimplemented from ReactorBase.

Definition at line 57 of file ReactorSurface.cpp.

◆ setArea()

void setArea ( double  a)
overridevirtual

Set the surface area [m²].

Reimplemented from ReactorBase.

Definition at line 62 of file ReactorSurface.cpp.

◆ thermo()

SurfPhase * thermo ( )
inline

Accessor for the SurfPhase object.

Definition at line 53 of file ReactorSurface.h.

◆ kinetics()

InterfaceKinetics * kinetics ( )
inline

Accessor for the InterfaceKinetics object.

Definition at line 58 of file ReactorSurface.h.

◆ addInlet()

void addInlet ( FlowDevice inlet)
inlineoverridevirtual

Connect an inlet FlowDevice to this reactor.

Reimplemented from ReactorBase.

Definition at line 62 of file ReactorSurface.h.

◆ addOutlet()

void addOutlet ( FlowDevice outlet)
inlineoverridevirtual

Connect an outlet FlowDevice to this reactor.

Reimplemented from ReactorBase.

Definition at line 67 of file ReactorSurface.h.

◆ addWall()

void addWall ( WallBase w,
int  lr 
)
inlineoverridevirtual

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 from ReactorBase.

Definition at line 72 of file ReactorSurface.h.

◆ addSurface()

void addSurface ( ReactorSurface surf)
inlineoverridevirtual

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 from ReactorBase.

Definition at line 76 of file ReactorSurface.h.

◆ nAdjacent()

size_t nAdjacent ( ) const
inline

Get the number of Reactor and Reservoir objects adjacent to this surface.

Since
New in Cantera 4.0.

Definition at line 82 of file ReactorSurface.h.

◆ adjacent()

shared_ptr< ReactorBase > adjacent ( size_t  n)
inline

Access an adjacent Reactor or Reservoir.

Since
New in Cantera 4.0.

Definition at line 88 of file ReactorSurface.h.

◆ setCoverages() [1/3]

void setCoverages ( const double *  cov)

Set the surface coverages.

Array cov has length equal to the number of surface species.

Definition at line 67 of file ReactorSurface.cpp.

◆ setCoverages() [2/3]

void setCoverages ( const Composition cov)

Set the surface coverages by name.

Definition at line 72 of file ReactorSurface.cpp.

◆ setCoverages() [3/3]

void setCoverages ( const string &  cov)

Set the surface coverages by name.

Definition at line 77 of file ReactorSurface.cpp.

◆ getCoverages()

void getCoverages ( double *  cov) const

Get the surface coverages.

Array cov should have length equal to the number of surface species.

Definition at line 82 of file ReactorSurface.cpp.

◆ getState()

void getState ( double *  y)
overridevirtual

Get the current state of the reactor.

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

Reimplemented from ReactorBase.

Definition at line 87 of file ReactorSurface.cpp.

◆ initialize()

void initialize ( double  t0 = 0.0)
overridevirtual

Initialize the reactor.

Called automatically by ReactorNet::initialize.

Reimplemented from ReactorBase.

Definition at line 92 of file ReactorSurface.cpp.

◆ initializeSteady()

vector< size_t > initializeSteady ( )
overridevirtual

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 from ReactorBase.

Definition at line 98 of file ReactorSurface.cpp.

◆ updateState()

void updateState ( double *  y)
overridevirtual

Set the state of the reactor to correspond to the state vector y.

Reimplemented from ReactorBase.

Definition at line 103 of file ReactorSurface.cpp.

◆ eval()

void eval ( double  t,
double *  LHS,
double *  RHS 
)
overridevirtual

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 from ReactorBase.

Definition at line 110 of file ReactorSurface.cpp.

◆ evalSteady()

void evalSteady ( double  t,
double *  LHS,
double *  RHS 
)
overridevirtual

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 from ReactorBase.

Definition at line 122 of file ReactorSurface.cpp.

◆ addSensitivityReaction()

void addSensitivityReaction ( size_t  rxn)
overridevirtual

Add a sensitivity parameter associated with the reaction number rxn

Reimplemented from ReactorBase.

Definition at line 167 of file ReactorSurface.cpp.

◆ applySensitivity()

void applySensitivity ( double *  params)
overridevirtual

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

Reimplemented from ReactorBase.

Definition at line 134 of file ReactorSurface.cpp.

◆ resetSensitivity()

void resetSensitivity ( double *  params)
overridevirtual

Reset the reaction rate multipliers.

Reimplemented from ReactorBase.

Definition at line 151 of file ReactorSurface.cpp.

◆ componentIndex()

size_t componentIndex ( const string &  nm) const
overridevirtual

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

Reimplemented from ReactorBase.

Definition at line 179 of file ReactorSurface.cpp.

◆ componentName()

string componentName ( size_t  k)
overridevirtual

Return the name of the solution component with index i.

See also
componentIndex()

Reimplemented from ReactorBase.

Definition at line 184 of file ReactorSurface.cpp.

◆ upperBound()

double upperBound ( size_t  k) const
overridevirtual

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

Reimplemented from ReactorBase.

Definition at line 189 of file ReactorSurface.cpp.

◆ lowerBound()

double lowerBound ( size_t  k) const
overridevirtual

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

Reimplemented from ReactorBase.

Definition at line 194 of file ReactorSurface.cpp.

◆ resetBadValues()

void resetBadValues ( double *  y)
overridevirtual

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

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

Reimplemented from ReactorBase.

Definition at line 199 of file ReactorSurface.cpp.

Member Data Documentation

◆ m_area

double m_area = 1.0
protected

Definition at line 127 of file ReactorSurface.h.

◆ m_surf

shared_ptr<SurfPhase> m_surf
protected

Definition at line 129 of file ReactorSurface.h.

◆ m_kinetics

shared_ptr<InterfaceKinetics> m_kinetics
protected

Definition at line 130 of file ReactorSurface.h.

◆ m_reactors

vector<shared_ptr<ReactorBase> > m_reactors
protected

Definition at line 131 of file ReactorSurface.h.


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