17 const vector<shared_ptr<ReactorBase>>& reactors,
22 if (reactors.empty()) {
24 "Surface requires at least one adjacent reactor.");
26 vector<shared_ptr<Solution>>
adjacent;
27 for (
auto R : reactors) {
29 m_reactors.push_back(R);
38 m_surf = std::dynamic_pointer_cast<SurfPhase>(
m_solution->thermo());
41 "Solution object must have a SurfPhase object as the thermo manager.");
44 if (!soln->kinetics() ) {
45 throw CanteraError(
"ReactorSurface::ReactorSurface",
46 "Solution object must have kinetics manager.");
47 }
else if (!std::dynamic_pointer_cast<InterfaceKinetics>(soln->kinetics())) {
48 throw CanteraError(
"ReactorSurface::ReactorSurface",
49 "Kinetics manager must be an InterfaceKinetics object.");
51 m_kinetics = std::dynamic_pointer_cast<InterfaceKinetics>(
m_solution->kinetics());
52 m_thermo = m_surf.get();
54 m_sdot.resize(m_kinetics->nTotalSpecies(), 0.0);
69 m_surf->setCoveragesNoNorm(cov);
74 m_surf->setCoveragesByName(cov);
79 m_surf->setCoveragesByName(cov);
84 m_surf->getCoverages(cov);
89 m_surf->getCoverages(y);
105 m_surf->setCoveragesNoNorm(y);
107 m_kinetics->getNetProductionRates(
m_sdot.data());
112 size_t nsp = m_surf->nSpecies();
113 double rs0 = 1.0 / m_surf->siteDensity();
115 for (
size_t k = 1; k < nsp; k++) {
116 RHS[k] =
m_sdot[k] * rs0 * m_surf->size(k);
125 vector<double> cov(
m_nsp);
126 m_surf->getCoverages(cov.data());
128 for (
size_t k = 0; k <
m_nsp; k++) {
139 for (
auto& p : m_sensParams) {
140 if (p.type == SensParameterType::reaction) {
141 p.value = m_kinetics->multiplier(p.local);
142 m_kinetics->setMultiplier(p.local, p.value*params[p.global]);
143 }
else if (p.type == SensParameterType::enthalpy) {
148 m_kinetics->invalidateCache();
156 for (
auto& p : m_sensParams) {
157 if (p.type == SensParameterType::reaction) {
158 m_kinetics->setMultiplier(p.local, p.value);
159 }
else if (p.type == SensParameterType::enthalpy) {
164 m_kinetics->invalidateCache();
169 if (rxn >= m_kinetics->nReactions()) {
170 throw CanteraError(
"ReactorSurface::addSensitivityReaction",
171 "Reaction number out of range ({})", rxn);
173 size_t p = m_reactors[0]->network().registerSensitivityParameter(
174 m_kinetics->reaction(rxn)->equation(), 1.0, 1.0);
175 m_sensParams.emplace_back(
181 return m_surf->speciesIndex(nm);
186 return m_surf->speciesName(k);
201 for (
size_t k = 0; k <
m_nsp; k++) {
202 y[k] = std::max(y[k], 0.0);
204 m_surf->setCoverages(y);
205 m_surf->getCoverages(y);
214 m_f_energy.resize(m_kinetics->nTotalSpecies(), 0.0);
215 m_f_species.resize(m_kinetics->nTotalSpecies(), 0.0);
216 m_kin2net.resize(m_kinetics->nTotalSpecies(), -1);
219 for (
size_t k = 0; k <
m_nsp; k++) {
222 for (
auto R : m_reactors) {
223 size_t nsp = R->phase()->thermo()->nSpecies();
224 size_t k0 = m_kinetics->speciesOffset(*R->phase()->thermo());
225 int offset =
static_cast<int>(R->offset() + R->speciesOffset());
226 for (
size_t k = 0; k < nsp; k++) {
235 m_surf->getCoverages(y);
236 double totalSites = m_surf->siteDensity() * m_area;
237 for (
size_t k = 0; k <
m_nsp; k++) {
238 y[k] *= totalSites / m_surf->size(k);
245 double totalSites = m_surf->siteDensity() * m_area;
246 for (
size_t k = 0; k <
m_nsp; k++) {
247 m_cov_tmp[k] *= m_surf->size(k) / totalSites;
249 m_surf->setCoveragesNoNorm(
m_cov_tmp.data());
251 m_kinetics->getNetProductionRates(
m_sdot.data());
256 for (
size_t k = 0; k <
m_nsp; k++) {
257 RHS[k] =
m_sdot[k] * m_area / m_surf->size(k);
264 double cov_sum = 0.0;
265 for (
size_t k = 0; k <
m_nsp; k++) {
268 RHS[0] = 1.0 - cov_sum;
283 for (
size_t k = 0; k <
m_nsp; k++) {
284 y[k] = std::max(y[k], 0.0);
290 auto sdot_ddC = m_kinetics->netProductionRates_ddCi();
291 for (
auto R : m_reactors) {
293 size_t nsp_R = R->phase()->thermo()->nSpecies();
294 size_t k0 = m_kinetics->speciesOffset(*R->phase()->thermo());
295 R->getJacobianScalingFactors(f_species, &
m_f_energy[k0]);
299 for (
int k = 0; k < sdot_ddC.outerSize(); k++) {
304 for (Eigen::SparseMatrix<double>::InnerIterator it(sdot_ddC, k); it; ++it) {
306 if (row == -1 || it.value() == 0.0) {
310 trips.emplace_back(row, col, it.value() *
m_f_species[k]);
312 trips.emplace_back(R->
offset(), col,
322 const vector<shared_ptr<ReactorBase>>& reactors,
337 return 2.0 * sqrt(
Pi * m_reactors[0]->
area());
342 size_t nsp = m_surf->nSpecies();
344 for (
size_t k = 1; k < nsp; k++) {
348 residual[0] = sum - 1.0;
364 std::fill(constraints, constraints +
m_nsp, 1.0);
Base class for kinetics managers and also contains the kineticsmgr module documentation (see Kinetics...
Header file for class ReactorSurface.
Header for a simple thermodynamics model of a surface phase derived from ThermoPhase,...
Base class for exceptions thrown by Cantera classes.
double area() const override
Surface area per unit length [m].
FlowReactorSurface(shared_ptr< Solution > soln, const vector< shared_ptr< ReactorBase > > &reactors, bool clone, const string &name="(none)")
Create a new ReactorSurface.
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.
void getStateDae(double *y, double *ydot) override
Get the current state and derivative vector of the reactor for a DAE solver.
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.
double m_ss_atol
steady-state absolute tolerance, used to determine initial surface coverages
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.
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.
Base class for reactor objects.
shared_ptr< Solution > m_solution
Composite thermo/kinetics/transport handler.
double pressure() const
Returns the current pressure (Pa) of the reactor.
size_t m_nv
Number of state variables for this reactor.
double temperature() const
Returns the current temperature (K) of the reactor's contents.
size_t m_offset
Offset into global ReactorNet state vector.
vector< double > m_sdot
species production rates on surfaces
size_t m_nsp
Number of homogeneous species in the mixture.
size_t offset() const
Get the starting offset for this reactor's state variables within the global state vector of the Reac...
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.
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.
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.
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.
void initialize(double t0=0.0) override
Initialize the reactor.
void setCoverages(const double *cov)
Set the surface coverages.
virtual void setState_TP(double t, double p)
Set the temperature (K) and pressure (Pa)
virtual void modifyOneHf298SS(const size_t k, const double Hf298New)
Modify the value of the 298 K Heat of Formation of one species in the phase (J kmol-1)
void invalidateCache() override
Invalidate any cached values which are normally updated only when a change in state is detected.
virtual void resetHf298(const size_t k=npos)
Restore the original heat of formation of one or more species.
Namespace for the Cantera kernel.
const double Tiny
Small number to compare differences of mole fractions against.
offset
Offsets of solution components in the 1D solution array.
const double BigNumber
largest number to compare to inf.
map< string, double > Composition
Map from string names to doubles.