Cantera  4.0.0a1
Loading...
Searching...
No Matches
ReactorSurface.cpp
Go to the documentation of this file.
1//! @file ReactorSurface.cpp
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
12
13namespace Cantera
14{
15
16ReactorSurface::ReactorSurface(shared_ptr<Solution> soln,
17 const vector<shared_ptr<ReactorBase>>& reactors,
18 bool clone,
19 const string& name)
20 : ReactorBase(name)
21{
22 if (reactors.empty()) {
23 throw CanteraError("ReactorSurface::ReactorSurface",
24 "Surface requires at least one adjacent reactor.");
25 }
26 vector<shared_ptr<Solution>> adjacent;
27 for (auto R : reactors) {
28 adjacent.push_back(R->phase());
29 m_reactors.push_back(R);
30 R->addSurface(this);
31 }
32 if (clone) {
33 m_solution = soln->clone(adjacent, true, false);
34 } else {
35 m_solution = soln;
36 }
37 m_solution->thermo()->addSpeciesLock();
38 m_surf = std::dynamic_pointer_cast<SurfPhase>(m_solution->thermo());
39 if (!m_surf) {
40 throw CanteraError("ReactorSurface::ReactorSurface",
41 "Solution object must have a SurfPhase object as the thermo manager.");
42 }
43
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.");
50 }
51 m_kinetics = std::dynamic_pointer_cast<InterfaceKinetics>(m_solution->kinetics());
52 m_thermo = m_surf.get();
53 m_nsp = m_nv = m_surf->nSpecies();
54 m_sdot.resize(m_kinetics->nTotalSpecies(), 0.0);
55}
56
58{
59 return m_area;
60}
61
63{
64 m_area = a;
65}
66
67void ReactorSurface::setCoverages(const double* cov)
68{
69 m_surf->setCoveragesNoNorm(cov);
70}
71
73{
74 m_surf->setCoveragesByName(cov);
75}
76
77void ReactorSurface::setCoverages(const string& cov)
78{
79 m_surf->setCoveragesByName(cov);
80}
81
82void ReactorSurface::getCoverages(double* cov) const
83{
84 m_surf->getCoverages(cov);
85}
86
88{
89 m_surf->getCoverages(y);
90}
91
93{
94 // Sync the surface temperature and pressure to that of the first adjacent reactor
95 m_thermo->setState_TP(m_reactors[0]->temperature(), m_reactors[0]->pressure());
96}
97
99{
100 return {0}; // sum of coverages constraint
101}
102
104{
105 m_surf->setCoveragesNoNorm(y);
106 m_thermo->setState_TP(m_reactors[0]->temperature(), m_reactors[0]->pressure());
107 m_kinetics->getNetProductionRates(m_sdot.data());
108}
109
110void ReactorSurface::eval(double t, double* LHS, double* RHS)
111{
112 size_t nsp = m_surf->nSpecies();
113 double rs0 = 1.0 / m_surf->siteDensity();
114 double sum = 0.0;
115 for (size_t k = 1; k < nsp; k++) {
116 RHS[k] = m_sdot[k] * rs0 * m_surf->size(k);
117 sum -= RHS[k];
118 }
119 RHS[0] = sum;
120}
121
122void ReactorSurface::evalSteady(double t, double* LHS, double* RHS)
123{
124 eval(t, LHS, RHS);
125 vector<double> cov(m_nsp);
126 m_surf->getCoverages(cov.data());
127 double sum = 0.0;
128 for (size_t k = 0; k < m_nsp; k++) {
129 sum += cov[k];
130 }
131 RHS[0] = 1.0 - sum;
132}
133
135{
136 if (!params) {
137 return;
138 }
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) {
144 m_thermo->modifyOneHf298SS(p.local, p.value + params[p.global]);
145 }
146 }
147 m_thermo->invalidateCache();
148 m_kinetics->invalidateCache();
149}
150
152{
153 if (!params) {
154 return;
155 }
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) {
160 m_thermo->resetHf298(p.local);
161 }
162 }
163 m_thermo->invalidateCache();
164 m_kinetics->invalidateCache();
165}
166
168{
169 if (rxn >= m_kinetics->nReactions()) {
170 throw CanteraError("ReactorSurface::addSensitivityReaction",
171 "Reaction number out of range ({})", rxn);
172 }
173 size_t p = m_reactors[0]->network().registerSensitivityParameter(
174 m_kinetics->reaction(rxn)->equation(), 1.0, 1.0);
175 m_sensParams.emplace_back(
176 SensitivityParameter{rxn, p, 1.0, SensParameterType::reaction});
177}
178
179size_t ReactorSurface::componentIndex(const string& nm) const
180{
181 return m_surf->speciesIndex(nm);
182}
183
185{
186 return m_surf->speciesName(k);
187}
188
189double ReactorSurface::upperBound(size_t k) const
190{
191 return 1.0;
192}
193
194double ReactorSurface::lowerBound(size_t k) const
195{
196 return -Tiny;
197}
198
200{
201 for (size_t k = 0; k < m_nsp; k++) {
202 y[k] = std::max(y[k], 0.0);
203 }
204 m_surf->setCoverages(y);
205 m_surf->getCoverages(y);
206}
207
208// ------ MoleReactorSurface methods ------
209
211{
213 m_cov_tmp.resize(m_nsp);
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);
217 m_kin2reactor.resize(m_kinetics->nTotalSpecies(), nullptr);
218
219 for (size_t k = 0; k < m_nsp; k++) {
220 m_kin2net[k] = static_cast<int>(m_offset + k);
221 }
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++) {
227 m_kin2net[k + k0] = offset + k;
228 m_kin2reactor[k + k0] = R.get();
229 }
230 }
231}
232
234{
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);
239 }
240}
241
243{
244 std::copy(y, y + m_nsp, m_cov_tmp.data());
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;
248 }
249 m_surf->setCoveragesNoNorm(m_cov_tmp.data());
250 m_thermo->setState_TP(m_reactors[0]->temperature(), m_reactors[0]->pressure());
251 m_kinetics->getNetProductionRates(m_sdot.data());
252}
253
254void MoleReactorSurface::eval(double t, double* LHS, double* RHS)
255{
256 for (size_t k = 0; k < m_nsp; k++) {
257 RHS[k] = m_sdot[k] * m_area / m_surf->size(k);
258 }
259}
260
261void MoleReactorSurface::evalSteady(double t, double* LHS, double* RHS)
262{
263 eval(t, LHS, RHS);
264 double cov_sum = 0.0;
265 for (size_t k = 0; k < m_nsp; k++) {
266 cov_sum += m_cov_tmp[k];
267 }
268 RHS[0] = 1.0 - cov_sum;
269}
270
271double MoleReactorSurface::upperBound(size_t k) const
272{
273 return BigNumber;
274}
275
276double MoleReactorSurface::lowerBound(size_t k) const
277{
278 return -Tiny;
279}
280
282{
283 for (size_t k = 0; k < m_nsp; k++) {
284 y[k] = std::max(y[k], 0.0);
285 }
286}
287
288void MoleReactorSurface::getJacobianElements(vector<Eigen::Triplet<double>>& trips)
289{
290 auto sdot_ddC = m_kinetics->netProductionRates_ddCi();
291 for (auto R : m_reactors) {
292 double f_species;
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]);
296 std::fill(&m_f_species[k0], &m_f_species[k0 + nsp_R], m_area * f_species);
297 }
298 std::fill(&m_f_species[0], &m_f_species[m_nsp], 1.0); // surface species
299 for (int k = 0; k < sdot_ddC.outerSize(); k++) {
300 int col = m_kin2net[k];
301 if (col == -1) {
302 continue;
303 }
304 for (Eigen::SparseMatrix<double>::InnerIterator it(sdot_ddC, k); it; ++it) {
305 int row = m_kin2net[it.row()];
306 if (row == -1 || it.value() == 0.0) {
307 continue;
308 }
309 ReactorBase* R = m_kin2reactor[it.row()];
310 trips.emplace_back(row, col, it.value() * m_f_species[k]);
311 if (m_f_energy[it.row()] != 0.0) {
312 trips.emplace_back(R->offset(), col,
313 it.value() * m_f_energy[it.row()] * m_f_species[k]);
314 }
315 }
316 }
317}
318
319// ------ FlowReactorSurface methods ------
320
322 const vector<shared_ptr<ReactorBase>>& reactors,
323 bool clone,
324 const string& name)
325 : ReactorSurface(soln, reactors, clone, name)
326{
327 m_area = -1.0; // default to perimeter of cylindrical reactor
328}
329
331 if (m_area > 0) {
332 return m_area;
333 }
334
335 // Assuming a cylindrical cross section, P = 2 * pi * r and A = pi * r^2, so
336 // P(A) = 2 * sqrt(pi * A)
337 return 2.0 * sqrt(Pi * m_reactors[0]->area());
338}
339
340void FlowReactorSurface::evalDae(double t, double* y, double* ydot, double* residual)
341{
342 size_t nsp = m_surf->nSpecies();
343 double sum = y[0];
344 for (size_t k = 1; k < nsp; k++) {
345 residual[k] = m_sdot[k];
346 sum += y[k];
347 }
348 residual[0] = sum - 1.0;
349}
350
351void FlowReactorSurface::getStateDae(double* y, double* ydot)
352{
353 // Advance the surface to steady state to get consistent initial coverages
354 m_kinetics->advanceCoverages(100.0, m_ss_rtol, m_ss_atol, 0, m_max_ss_steps,
356 getCoverages(y);
357 // Update the values in m_sdot that are needed by the adjacent FlowReactor
358 updateState(y);
359}
360
361void FlowReactorSurface::getConstraints(double* constraints)
362{
363 // The species coverages are algebraic constraints
364 std::fill(constraints, constraints + m_nsp, 1.0);
365}
366
367}
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.
Definition ReactorBase.h:51
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.
const double Pi
Pi.
Definition ct_defs.h:70
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const double Tiny
Small number to compare differences of mole fractions against.
Definition ct_defs.h:175
offset
Offsets of solution components in the 1D solution array.
Definition Flow1D.h:25
const double BigNumber
largest number to compare to inf.
Definition ct_defs.h:162
map< string, double > Composition
Map from string names to doubles.
Definition ct_defs.h:179