6#ifndef CT_REACTOR_SURFACE_H
7#define CT_REACTOR_SURFACE_H
37 const vector<shared_ptr<ReactorBase>>& reactors,
39 const string&
name=
"(none)");
42 string type()
const override {
43 return "ReactorSurface";
47 double area()
const override;
50 void setArea(
double a)
override;
59 return m_kinetics.get();
64 "Inlets are undefined for reactors of type '{}'.",
type());
69 "Outlets are undefined for reactors of type '{}'.",
type());
83 return m_reactors.size();
89 return m_reactors.at(n);
110 void eval(
double t,
double* LHS,
double* RHS)
override;
111 void evalSteady(
double t,
double* LHS,
double* RHS)
override;
129 shared_ptr<SurfPhase> m_surf;
130 shared_ptr<InterfaceKinetics> m_kinetics;
131 vector<shared_ptr<ReactorBase>> m_reactors;
146 string type()
const override {
return "MoleReactorSurface"; }
150 void eval(
double t,
double* LHS,
double* RHS)
override;
151 void evalSteady(
double t,
double* LHS,
double* RHS)
override;
169 vector<Eigen::SparseMatrix<double>::StorageIndex>
m_kin2net;
186 const vector<shared_ptr<ReactorBase>>& reactors,
188 const string&
name=
"(none)");
191 return "FlowReactorSurface";
196 void evalDae(
double t,
double* y,
double* ydot,
double* residual)
override;
197 void getStateDae(
double* y,
double* ydot)
override;
203 double area()
const override;
210 void setArea(
double A)
override { m_area = A; }
Base class for 'flow devices' (valves, pressure regulators, etc.) connecting reactors.
A surface in contact with a FlowReactor.
double area() const override
Surface area per unit length [m].
double initialAtol() const
Get 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 setInitialMaxErrorFailures(int max_fails)
Set 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 m_ss_rtol
steady-state relative tolerance, used to determine initial surface coverages
int m_max_ss_error_fails
maximum number of steady-state integrator error test failures
void getConstraints(double *constraints) override
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
bool timeIsIndependent() const override
Indicates whether the governing equations for this reactor are functions of time or a spatial variabl...
string type() const override
String indicating the reactor model implemented.
double initialMaxErrorFailures() const
Get the steady state tolerances used to determine the initial state for surface coverages.
void getStateDae(double *y, double *ydot) override
Get the current state and derivative vector of the reactor for a DAE solver.
void setInitialAtol(double atol)
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.
int m_max_ss_steps
maximum number of steady-state coverage integrator-steps
void evalDae(double t, double *y, double *ydot, double *residual) override
Evaluate the reactor governing equations.
void setInitialRtol(double rtol)
Set the steady state tolerances used to determine the initial state for surface coverages.
double m_ss_atol
steady-state absolute tolerance, used to determine initial surface coverages
void setArea(double A) override
Set the reactor surface area per unit length [m].
A kinetics manager for heterogeneous reaction mechanisms.
A surface where the state variables are the total number of moles of each species.
double upperBound(size_t k) const override
Get the upper bound on the k-th component of the local state vector.
vector< double > m_f_energy
Temporary storage for d(energy)/d(moles) scaling factors.
vector< Eigen::SparseMatrix< double >::StorageIndex > m_kin2net
Mapping from InterfaceKinetics species index to ReactorNet state vector index.
void resetBadValues(double *y) override
Reset physically or mathematically problematic values, such as negative species concentrations.
vector< double > m_f_species
Temporary storage for d(moles)/d(moles) scaling factor.
void eval(double t, double *LHS, double *RHS) override
Evaluate the reactor governing equations.
string type() const override
String indicating the reactor model implemented.
vector< ReactorBase * > m_kin2reactor
Mapping from InterfaceKinetics species index to the corresponding reactor.
void evalSteady(double t, double *LHS, double *RHS) override
Evaluate the governing equations with modifications for the steady-state solver.
void getState(double *y) override
Get the current state of the reactor.
void getJacobianElements(vector< Eigen::Triplet< double > > &trips) override
Get Jacobian elements for this reactor within the full reactor network.
double lowerBound(size_t k) const override
Get the lower bound on the k-th component of the local state vector.
void updateState(double *y) override
Set the state of the reactor to correspond to the state vector y.
void initialize(double t0=0.0) override
Initialize the reactor.
vector< double > m_cov_tmp
Temporary storage for coverages.
An error indicating that an unimplemented function has been called.
Base class for reactor objects.
FlowDevice & outlet(size_t n=0)
Return a reference to the n-th outlet FlowDevice connected to this reactor.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
string name() const
Return the name of this reactor.
A surface where reactions can occur that is in contact with the bulk fluid of a Reactor.
double upperBound(size_t k) const override
Get the upper bound on the k-th component of the local state vector.
double area() const override
Returns the surface area [m²].
void resetBadValues(double *y) override
Reset physically or mathematically problematic values, such as negative species concentrations.
void addSurface(ReactorSurface *surf) override
Add a ReactorSurface object to a Reactor object.
void addOutlet(FlowDevice &outlet) override
Connect an outlet FlowDevice to this reactor.
ReactorSurface(shared_ptr< Solution > soln, const vector< shared_ptr< ReactorBase > > &reactors, bool clone, const string &name="(none)")
Create a new ReactorSurface.
void setArea(double a) override
Set the surface area [m²].
void resetSensitivity(double *params) override
Reset the reaction rate multipliers.
void eval(double t, double *LHS, double *RHS) override
Evaluate the reactor governing equations.
string type() const override
String indicating the wall model implemented.
size_t componentIndex(const string &nm) const override
Return the index in the solution vector for this reactor of the component named nm.
void applySensitivity(double *params) override
Set reaction rate multipliers based on the sensitivity variables in params.
void getCoverages(double *cov) const
Get the surface coverages.
void evalSteady(double t, double *LHS, double *RHS) override
Evaluate the governing equations with modifications for the steady-state solver.
void getState(double *y) override
Get the current state of the reactor.
vector< size_t > initializeSteady() override
Initialize the reactor before solving a steady-state problem.
size_t nAdjacent() const
Get the number of Reactor and Reservoir objects adjacent to this surface.
void addSensitivityReaction(size_t rxn) override
Add a sensitivity parameter associated with the reaction number rxn
double lowerBound(size_t k) const override
Get the lower bound on the k-th component of the local state vector.
string componentName(size_t k) override
Return the name of the solution component with index i.
shared_ptr< ReactorBase > adjacent(size_t n)
Access an adjacent Reactor or Reservoir.
void updateState(double *y) override
Set the state of the reactor to correspond to the state vector y.
InterfaceKinetics * kinetics()
Accessor for the InterfaceKinetics object.
void initialize(double t0=0.0) override
Initialize the reactor.
void setCoverages(const double *cov)
Set the surface coverages.
SurfPhase * thermo()
Accessor for the SurfPhase object.
void addInlet(FlowDevice &inlet) override
Connect an inlet FlowDevice to this reactor.
void addWall(WallBase &w, int lr) override
Insert a Wall between this reactor and another reactor.
A simple thermodynamic model for a surface phase, assuming an ideal solution model.
Base class for 'walls' (walls, pistons, etc.) connecting reactors.
Namespace for the Cantera kernel.
map< string, double > Composition
Map from string names to doubles.