Cantera  3.1.0a1
Loading...
Searching...
No Matches
Reactor.cpp
Go to the documentation of this file.
1//! @file Reactor.cpp A zero-dimensional reactor
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
16
17#include <boost/math/tools/roots.hpp>
18
19using namespace std;
20namespace bmt = boost::math::tools;
21
22namespace Cantera
23{
24
25Reactor::Reactor(shared_ptr<Solution> sol, const string& name)
26{
27 m_solution = sol;
28 m_name = name;
29 setThermoMgr(*sol->thermo());
30 setKineticsMgr(*sol->kinetics());
31}
32
34{
36 // translate settings to surfaces
37 for (auto S : m_surfaces) {
38 S->kinetics()->setDerivativeSettings(settings);
39 }
40}
41
43{
44 m_kin = &kin;
45 if (m_kin->nReactions() == 0) {
46 setChemistry(false);
47 } else {
48 setChemistry(true);
49 }
50}
51
52void Reactor::getState(double* y)
53{
54 if (m_thermo == 0) {
55 throw CanteraError("Reactor::getState",
56 "Error: reactor is empty.");
57 }
58 m_thermo->restoreState(m_state);
59
60 // set the first component to the total mass
61 m_mass = m_thermo->density() * m_vol;
62 y[0] = m_mass;
63
64 // set the second component to the total volume
65 y[1] = m_vol;
66
67 // set the third component to the total internal energy
68 y[2] = m_thermo->intEnergy_mass() * m_mass;
69
70 // set components y+3 ... y+K+2 to the mass fractions of each species
71 m_thermo->getMassFractions(y+3);
72
73 // set the remaining components to the surface species
74 // coverages on the walls
76}
77
79{
80 size_t loc = 0;
81 for (auto& S : m_surfaces) {
82 S->getCoverages(y + loc);
83 loc += S->thermo()->nSpecies();
84 }
85}
86
87void Reactor::initialize(double t0)
88{
89 if (!m_thermo || (m_chem && !m_kin)) {
90 throw CanteraError("Reactor::initialize", "Reactor contents not set"
91 " for reactor '" + m_name + "'.");
92 }
93 m_thermo->restoreState(m_state);
94 m_sdot.resize(m_nsp, 0.0);
95 m_wdot.resize(m_nsp, 0.0);
96 updateConnected(true);
97
98 for (size_t n = 0; n < m_wall.size(); n++) {
99 WallBase* W = m_wall[n];
100 W->initialize();
101 }
102
103 m_nv = m_nsp + 3;
104 m_nv_surf = 0;
105 size_t maxnt = 0;
106 for (auto& S : m_surfaces) {
107 m_nv_surf += S->thermo()->nSpecies();
108 size_t nt = S->kinetics()->nTotalSpecies();
109 maxnt = std::max(maxnt, nt);
110 }
111 m_nv += m_nv_surf;
112 m_work.resize(maxnt);
113}
114
116{
117 size_t ns = m_sensParams.size();
118 for (auto& S : m_surfaces) {
119 ns += S->nSensParams();
120 }
121 return ns;
122}
123
125{
127 m_mass = m_thermo->density() * m_vol;
128}
129
130void Reactor::updateState(double* y)
131{
132 // The components of y are [0] the total mass, [1] the total volume,
133 // [2] the total internal energy, [3...K+3] are the mass fractions of each
134 // species, and [K+3...] are the coverages of surface species on each wall.
135 m_mass = y[0];
136 m_vol = y[1];
137 m_thermo->setMassFractions_NoNorm(y+3);
138
139 if (m_energy) {
140 double U = y[2];
141 // Residual function: error in internal energy as a function of T
142 auto u_err = [this, U](double T) {
143 m_thermo->setState_TD(T, m_mass / m_vol);
144 return m_thermo->intEnergy_mass() * m_mass - U;
145 };
146
147 double T = m_thermo->temperature();
148 boost::uintmax_t maxiter = 100;
149 pair<double, double> TT;
150 try {
151 TT = bmt::bracket_and_solve_root(
152 u_err, T, 1.2, true, bmt::eps_tolerance<double>(48), maxiter);
153 } catch (std::exception&) {
154 // Try full-range bisection if bracketing fails (for example, near
155 // temperature limits for the phase's equation of state)
156 try {
157 TT = bmt::bisect(u_err, m_thermo->minTemp(), m_thermo->maxTemp(),
158 bmt::eps_tolerance<double>(48), maxiter);
159 } catch (std::exception& err2) {
160 // Set m_thermo back to a reasonable state if root finding fails
161 m_thermo->setState_TD(T, m_mass / m_vol);
162 throw CanteraError("Reactor::updateState",
163 "{}\nat U = {}, rho = {}", err2.what(), U, m_mass / m_vol);
164 }
165 }
166 if (fabs(TT.first - TT.second) > 1e-7*TT.first) {
167 throw CanteraError("Reactor::updateState", "root finding failed");
168 }
169 m_thermo->setState_TD(TT.second, m_mass / m_vol);
170 } else {
171 m_thermo->setDensity(m_mass/m_vol);
172 }
173
174 updateConnected(true);
175 updateSurfaceState(y + m_nsp + 3);
176}
177
179{
180 size_t loc = 0;
181 for (auto& S : m_surfaces) {
182 S->setCoverages(y+loc);
183 loc += S->thermo()->nSpecies();
184 }
185}
186
187void Reactor::updateConnected(bool updatePressure) {
188 // save parameters needed by other connected reactors
189 m_enthalpy = m_thermo->enthalpy_mass();
190 if (updatePressure) {
191 m_pressure = m_thermo->pressure();
192 }
193 m_intEnergy = m_thermo->intEnergy_mass();
194 m_thermo->saveState(m_state);
195
196 // Update the mass flow rate of connected flow devices
197 double time = 0.0;
198 if (m_net != nullptr) {
199 time = (timeIsIndependent()) ? m_net->time() : m_net->distance();
200 }
201 for (size_t i = 0; i < m_outlet.size(); i++) {
202 m_outlet[i]->setSimTime(time);
203 m_outlet[i]->updateMassFlowRate(time);
204 }
205 for (size_t i = 0; i < m_inlet.size(); i++) {
206 m_inlet[i]->setSimTime(time);
207 m_inlet[i]->updateMassFlowRate(time);
208 }
209 for (size_t i = 0; i < m_wall.size(); i++) {
210 m_wall[i]->setSimTime(time);
211 }
212}
213
214void Reactor::eval(double time, double* LHS, double* RHS)
215{
216 double& dmdt = RHS[0];
217 double* mdYdt = RHS + 3; // mass * dY/dt
218
219 evalWalls(time);
220 m_thermo->restoreState(m_state);
221 const vector<double>& mw = m_thermo->molecularWeights();
222 const double* Y = m_thermo->massFractions();
223
224 evalSurfaces(LHS + m_nsp + 3, RHS + m_nsp + 3, m_sdot.data());
225 // mass added to gas phase from surface reactions
226 double mdot_surf = dot(m_sdot.begin(), m_sdot.end(), mw.begin());
227 dmdt = mdot_surf;
228
229 // volume equation
230 RHS[1] = m_vdot;
231
232 if (m_chem) {
233 m_kin->getNetProductionRates(&m_wdot[0]); // "omega dot"
234 }
235
236 for (size_t k = 0; k < m_nsp; k++) {
237 // production in gas phase and from surfaces
238 mdYdt[k] = (m_wdot[k] * m_vol + m_sdot[k]) * mw[k];
239 // dilution by net surface mass flux
240 mdYdt[k] -= Y[k] * mdot_surf;
241 LHS[k+3] = m_mass;
242 }
243
244 // Energy equation.
245 // @f[
246 // \dot U = -P\dot V + A \dot q + \dot m_{in} h_{in} - \dot m_{out} h.
247 // @f]
248 if (m_energy) {
249 RHS[2] = - m_thermo->pressure() * m_vdot + m_Qdot;
250 } else {
251 RHS[2] = 0.0;
252 }
253
254 // add terms for outlets
255 for (auto outlet : m_outlet) {
256 double mdot = outlet->massFlowRate();
257 dmdt -= mdot; // mass flow out of system
258 if (m_energy) {
259 RHS[2] -= mdot * m_enthalpy;
260 }
261 }
262
263 // add terms for inlets
264 for (auto inlet : m_inlet) {
265 double mdot = inlet->massFlowRate();
266 dmdt += mdot; // mass flow into system
267 for (size_t n = 0; n < m_nsp; n++) {
268 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
269 // flow of species into system and dilution by other species
270 mdYdt[n] += (mdot_spec - mdot * Y[n]);
271 }
272 if (m_energy) {
273 RHS[2] += mdot * inlet->enthalpy_mass();
274 }
275 }
276}
277
278void Reactor::evalWalls(double t)
279{
280 // time is currently unused
281 m_vdot = 0.0;
282 m_Qdot = 0.0;
283 for (size_t i = 0; i < m_wall.size(); i++) {
284 int f = 2 * m_lr[i] - 1;
285 m_vdot -= f * m_wall[i]->expansionRate();
286 m_Qdot += f * m_wall[i]->heatRate();
287 }
288}
289
290void Reactor::evalSurfaces(double* LHS, double* RHS, double* sdot)
291{
292 fill(sdot, sdot + m_nsp, 0.0);
293 size_t loc = 0; // offset into ydot
294
295 for (auto S : m_surfaces) {
296 Kinetics* kin = S->kinetics();
297 SurfPhase* surf = S->thermo();
298
299 double rs0 = 1.0/surf->siteDensity();
300 size_t nk = surf->nSpecies();
301 double sum = 0.0;
302 S->syncState();
303 kin->getNetProductionRates(&m_work[0]);
304 for (size_t k = 1; k < nk; k++) {
305 RHS[loc + k] = m_work[k] * rs0 * surf->size(k);
306 sum -= RHS[loc + k];
307 }
308 RHS[loc] = sum;
309 loc += nk;
310
311 size_t bulkloc = kin->kineticsSpeciesIndex(m_thermo->speciesName(0));
312 double wallarea = S->area();
313 for (size_t k = 0; k < m_nsp; k++) {
314 sdot[k] += m_work[bulkloc + k] * wallarea;
315 }
316 }
317}
318
319Eigen::SparseMatrix<double> Reactor::finiteDifferenceJacobian()
320{
321 if (m_nv == 0) {
322 throw CanteraError("Reactor::finiteDifferenceJacobian",
323 "Reactor must be initialized first.");
324 }
325 // clear former jacobian elements
326 m_jac_trips.clear();
327
328 Eigen::ArrayXd yCurrent(m_nv);
329 getState(yCurrent.data());
330 double time = (m_net != nullptr) ? m_net->time() : 0.0;
331
332 Eigen::ArrayXd yPerturbed = yCurrent;
333 Eigen::ArrayXd lhsPerturbed(m_nv), lhsCurrent(m_nv);
334 Eigen::ArrayXd rhsPerturbed(m_nv), rhsCurrent(m_nv);
335 lhsCurrent = 1.0;
336 rhsCurrent = 0.0;
337 updateState(yCurrent.data());
338 eval(time, lhsCurrent.data(), rhsCurrent.data());
339
340 double rel_perturb = std::sqrt(std::numeric_limits<double>::epsilon());
341 double atol = (m_net != nullptr) ? m_net->atol() : 1e-15;
342
343 for (size_t j = 0; j < m_nv; j++) {
344 yPerturbed = yCurrent;
345 double delta_y = std::max(std::abs(yCurrent[j]), 1000 * atol) * rel_perturb;
346 yPerturbed[j] += delta_y;
347
348 updateState(yPerturbed.data());
349 lhsPerturbed = 1.0;
350 rhsPerturbed = 0.0;
351 eval(time, lhsPerturbed.data(), rhsPerturbed.data());
352
353 // d ydot_i/dy_j
354 for (size_t i = 0; i < m_nv; i++) {
355 double ydotPerturbed = rhsPerturbed[i] / lhsPerturbed[i];
356 double ydotCurrent = rhsCurrent[i] / lhsCurrent[i];
357 if (ydotCurrent != ydotPerturbed) {
358 m_jac_trips.emplace_back(
359 static_cast<int>(i), static_cast<int>(j),
360 (ydotPerturbed - ydotCurrent) / delta_y);
361 }
362 }
363 }
364 updateState(yCurrent.data());
365
366 Eigen::SparseMatrix<double> jac(m_nv, m_nv);
367 jac.setFromTriplets(m_jac_trips.begin(), m_jac_trips.end());
368 return jac;
369}
370
371
372void Reactor::evalSurfaces(double* RHS, double* sdot)
373{
374 fill(sdot, sdot + m_nsp, 0.0);
375 size_t loc = 0; // offset into ydot
376
377 for (auto S : m_surfaces) {
378 Kinetics* kin = S->kinetics();
379 SurfPhase* surf = S->thermo();
380
381 double rs0 = 1.0/surf->siteDensity();
382 size_t nk = surf->nSpecies();
383 double sum = 0.0;
384 S->syncState();
385 kin->getNetProductionRates(&m_work[0]);
386 for (size_t k = 1; k < nk; k++) {
387 RHS[loc + k] = m_work[k] * rs0 * surf->size(k);
388 sum -= RHS[loc + k];
389 }
390 RHS[loc] = sum;
391 loc += nk;
392
393 size_t bulkloc = kin->kineticsSpeciesIndex(m_thermo->speciesName(0));
394 double wallarea = S->area();
395 for (size_t k = 0; k < m_nsp; k++) {
396 sdot[k] += m_work[bulkloc + k] * wallarea;
397 }
398 }
399}
400
402{
403 if (!m_chem || rxn >= m_kin->nReactions()) {
404 throw CanteraError("Reactor::addSensitivityReaction",
405 "Reaction number out of range ({})", rxn);
406 }
407
409 name()+": "+m_kin->reaction(rxn)->equation(), 1.0, 1.0);
410 m_sensParams.emplace_back(
411 SensitivityParameter{rxn, p, 1.0, SensParameterType::reaction});
412}
413
415{
416 if (k >= m_thermo->nSpecies()) {
417 throw CanteraError("Reactor::addSensitivitySpeciesEnthalpy",
418 "Species index out of range ({})", k);
419 }
420
422 name() + ": " + m_thermo->speciesName(k) + " enthalpy",
423 0.0, GasConstant * 298.15);
424 m_sensParams.emplace_back(
425 SensitivityParameter{k, p, m_thermo->Hf298SS(k),
426 SensParameterType::enthalpy});
427}
428
429size_t Reactor::speciesIndex(const string& nm) const
430{
431 // check for a gas species name
432 size_t k = m_thermo->speciesIndex(nm);
433 if (k != npos) {
434 return k;
435 }
436
437 // check for a wall species
438 size_t offset = m_nsp;
439 for (auto& S : m_surfaces) {
440 ThermoPhase* th = S->thermo();
441 k = th->speciesIndex(nm);
442 if (k != npos) {
443 return k + offset;
444 } else {
445 offset += th->nSpecies();
446 }
447 }
448 return npos;
449}
450
451size_t Reactor::componentIndex(const string& nm) const
452{
453 size_t k = speciesIndex(nm);
454 if (k != npos) {
455 return k + 3;
456 } else if (nm == "mass") {
457 return 0;
458 } else if (nm == "volume") {
459 return 1;
460 } else if (nm == "int_energy") {
461 return 2;
462 } else {
463 return npos;
464 }
465}
466
467string Reactor::componentName(size_t k) {
468 if (k == 0) {
469 return "mass";
470 } else if (k == 1) {
471 return "volume";
472 } else if (k == 2) {
473 return "int_energy";
474 } else if (k >= 3 && k < neq()) {
475 k -= 3;
476 if (k < m_thermo->nSpecies()) {
477 return m_thermo->speciesName(k);
478 } else {
479 k -= m_thermo->nSpecies();
480 }
481 for (auto& S : m_surfaces) {
482 ThermoPhase* th = S->thermo();
483 if (k < th->nSpecies()) {
484 return th->speciesName(k);
485 } else {
486 k -= th->nSpecies();
487 }
488 }
489 }
490 throw CanteraError("Reactor::componentName", "Index is out of bounds.");
491}
492
493void Reactor::applySensitivity(double* params)
494{
495 if (!params) {
496 return;
497 }
498 for (auto& p : m_sensParams) {
499 if (p.type == SensParameterType::reaction) {
500 p.value = m_kin->multiplier(p.local);
501 m_kin->setMultiplier(p.local, p.value*params[p.global]);
502 } else if (p.type == SensParameterType::enthalpy) {
503 m_thermo->modifyOneHf298SS(p.local, p.value + params[p.global]);
504 }
505 }
506 for (auto& S : m_surfaces) {
507 S->setSensitivityParameters(params);
508 }
509 m_thermo->invalidateCache();
510 if (m_kin) {
511 m_kin->invalidateCache();
512 }
513}
514
515void Reactor::resetSensitivity(double* params)
516{
517 if (!params) {
518 return;
519 }
520 for (auto& p : m_sensParams) {
521 if (p.type == SensParameterType::reaction) {
522 m_kin->setMultiplier(p.local, p.value);
523 } else if (p.type == SensParameterType::enthalpy) {
524 m_thermo->resetHf298(p.local);
525 }
526 }
527 for (auto& S : m_surfaces) {
528 S->resetSensitivityParameters();
529 }
530 m_thermo->invalidateCache();
531 if (m_kin) {
532 m_kin->invalidateCache();
533 }
534}
535
536void Reactor::setAdvanceLimits(const double *limits)
537{
538 if (m_thermo == 0) {
539 throw CanteraError("Reactor::setAdvanceLimits",
540 "Error: reactor is empty.");
541 }
542 m_advancelimits.assign(limits, limits + m_nv);
543
544 // resize to zero length if no limits are set
545 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
546 [](double val){return val>0;})) {
547 m_advancelimits.resize(0);
548 }
549}
550
551bool Reactor::getAdvanceLimits(double *limits) const
552{
553 bool has_limit = hasAdvanceLimits();
554 if (has_limit) {
555 std::copy(m_advancelimits.begin(), m_advancelimits.end(), limits);
556 } else {
557 std::fill(limits, limits + m_nv, -1.0);
558 }
559 return has_limit;
560}
561
562void Reactor::setAdvanceLimit(const string& nm, const double limit)
563{
564 size_t k = componentIndex(nm);
565 if (k == npos) {
566 throw CanteraError("Reactor::setAdvanceLimit", "No component named '{}'", nm);
567 }
568
569 if (m_thermo == 0) {
570 throw CanteraError("Reactor::setAdvanceLimit",
571 "Error: reactor is empty.");
572 }
573 if (m_nv == 0) {
574 if (m_net == 0) {
575 throw CanteraError("Reactor::setAdvanceLimit",
576 "Cannot set limit on a reactor that is not "
577 "assigned to a ReactorNet object.");
578 } else {
579 m_net->initialize();
580 }
581 } else if (k > m_nv) {
582 throw CanteraError("Reactor::setAdvanceLimit",
583 "Index out of bounds.");
584 }
585 m_advancelimits.resize(m_nv, -1.0);
586 m_advancelimits[k] = limit;
587
588 // resize to zero length if no limits are set
589 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
590 [](double val){return val>0;})) {
591 m_advancelimits.resize(0);
592 }
593}
594
595}
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,...
Header file for base class WallBase.
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:427
Base class for exceptions thrown by Cantera classes.
double outletSpeciesMassFlowRate(size_t k)
Mass flow rate (kg/s) of outlet species k.
double enthalpy_mass()
specific enthalpy
double massFlowRate()
Mass flow rate (kg/s).
Definition FlowDevice.h:39
Public interface for kinetics managers.
Definition Kinetics.h:125
double multiplier(size_t i) const
The current value of the multiplier for reaction i.
Definition Kinetics.h:1360
shared_ptr< Reaction > reaction(size_t i)
Return the Reaction object for reaction i.
Definition Kinetics.cpp:684
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition Kinetics.h:152
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition Kinetics.h:276
virtual void setMultiplier(size_t i, double f)
Set the multiplier for reaction i to f.
Definition Kinetics.h:1369
virtual void getNetProductionRates(double *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:363
void restoreState(const vector< double > &state)
Restore a state saved on a previous call to saveState.
Definition Phase.cpp:260
size_t nSpecies() const
Returns the number of species in the phase.
Definition Phase.h:231
virtual void setMassFractions_NoNorm(const double *const y)
Set the mass fractions to the specified values without normalizing.
Definition Phase.cpp:355
void setState_TD(double t, double rho)
Set the internally stored temperature (K) and density (kg/m^3)
Definition Phase.cpp:377
double temperature() const
Temperature (K).
Definition Phase.h:562
void saveState(vector< double > &state) const
Save the current internal state of the phase.
Definition Phase.cpp:236
string speciesName(size_t k) const
Name of the species with index k.
Definition Phase.cpp:142
virtual void setDensity(const double density_)
Set the internally stored density (kg/m^3) of the phase.
Definition Phase.cpp:586
const double * massFractions() const
Return a const pointer to the mass fraction array.
Definition Phase.h:442
const vector< double > & molecularWeights() const
Return a const reference to the internal vector of molecular weights.
Definition Phase.cpp:395
size_t speciesIndex(const string &name) const
Returns the index of a species named 'name' within the Phase object.
Definition Phase.cpp:129
virtual double density() const
Density (kg/m^3).
Definition Phase.h:587
void getMassFractions(double *const y) const
Get the species mass fractions.
Definition Phase.cpp:471
virtual double pressure() const
Return the thermodynamic pressure (Pa).
Definition Phase.h:580
FlowDevice & outlet(size_t n=0)
Return a reference to the n-th outlet FlowDevice connected to this reactor.
shared_ptr< Solution > m_solution
Composite thermo/kinetics/transport handler.
double m_pressure
Current pressure in the reactor [Pa].
ReactorNet * m_net
The ReactorNet that this reactor is part of.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
vector< int > m_lr
Vector of length nWalls(), indicating whether this reactor is on the left (0) or right (1) of each wa...
double m_vol
Current volume of the reactor [m^3].
virtual void syncState()
Set the state of the reactor to correspond to the state of the associated ThermoPhase object.
double m_intEnergy
Current internal energy of the reactor [J/kg].
size_t m_nsp
Number of homogeneous species in the mixture.
virtual void setThermoMgr(ThermoPhase &thermo)
Specify the mixture contained in the reactor.
ReactorNet & network()
The ReactorNet that this reactor belongs to.
double m_enthalpy
Current specific enthalpy of the reactor [J/kg].
string name() const
Return the name of this reactor.
Definition ReactorBase.h:70
void initialize()
Initialize the reactor network.
double time()
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
size_t registerSensitivityParameter(const string &name, double value, double scale)
Used by Reactor and Wall objects to register the addition of sensitivity parameters so that the React...
double distance()
Current position [m] along the length of the reactor network, for reactors that are solved as a funct...
double atol()
Absolute integration tolerance.
Definition ReactorNet.h:95
virtual string componentName(size_t k)
Return the name of the solution component with index i.
Definition Reactor.cpp:467
void setKineticsMgr(Kinetics &kin) override
Definition Reactor.cpp:42
virtual void evalSurfaces(double *LHS, double *RHS, double *sdot)
Evaluate terms related to surface reactions.
Definition Reactor.cpp:290
virtual size_t componentIndex(const string &nm) const
Return the index in the solution vector for this reactor of the component named nm.
Definition Reactor.cpp:451
virtual void updateSurfaceState(double *y)
Update the state of SurfPhase objects attached to this reactor.
Definition Reactor.cpp:178
virtual void applySensitivity(double *params)
Set reaction rate multipliers based on the sensitivity variables in params.
Definition Reactor.cpp:493
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
Definition Reactor.h:283
vector< double > m_wdot
Species net molar production rates.
Definition Reactor.h:295
virtual void eval(double t, double *LHS, double *RHS)
Evaluate the reactor governing equations.
Definition Reactor.cpp:214
virtual void evalWalls(double t)
Evaluate terms related to Walls.
Definition Reactor.cpp:278
Eigen::SparseMatrix< double > finiteDifferenceJacobian()
Calculate the reactor-specific Jacobian using a finite difference method.
Definition Reactor.cpp:319
size_t neq()
Number of equations (state variables) for this reactor.
Definition Reactor.h:109
virtual void addSensitivitySpeciesEnthalpy(size_t k)
Add a sensitivity parameter associated with the enthalpy formation of species k (in the homogeneous p...
Definition Reactor.cpp:414
double m_Qdot
net heat transfer into the reactor, through walls [W]
Definition Reactor.h:287
vector< double > m_advancelimits
!< Number of variables associated with reactor surfaces
Definition Reactor.h:302
virtual void addSensitivityReaction(size_t rxn)
Add a sensitivity parameter associated with the reaction number rxn (in the homogeneous phase).
Definition Reactor.cpp:401
virtual size_t nSensParams() const
Number of sensitivity parameters associated with this reactor (including walls)
Definition Reactor.cpp:115
virtual void setDerivativeSettings(AnyMap &settings)
Use this to set the kinetics objects derivative settings.
Definition Reactor.cpp:33
virtual void updateState(double *y)
Set the state of the reactor to correspond to the state vector y.
Definition Reactor.cpp:130
void setChemistry(bool cflag=true) override
Enable or disable changes in reactor composition due to chemical reactions.
Definition Reactor.h:85
void setAdvanceLimits(const double *limits)
Set absolute step size limits during advance.
Definition Reactor.cpp:536
double m_mass
total mass
Definition Reactor.h:289
vector< Eigen::Triplet< double > > m_jac_trips
Vector of triplets representing the jacobian.
Definition Reactor.h:308
void setAdvanceLimit(const string &nm, const double limit)
Set individual step size limit for component name nm
Definition Reactor.cpp:562
virtual void resetSensitivity(double *params)
Reset the reaction rate multipliers.
Definition Reactor.cpp:515
vector< double > m_sdot
Production rates of gas phase species on surfaces [kmol/s].
Definition Reactor.h:293
virtual void getState(double *y)
Get the the current state of the reactor.
Definition Reactor.cpp:52
bool hasAdvanceLimits() const
Check whether Reactor object uses advance limits.
Definition Reactor.h:192
double m_vdot
net rate of volume change from moving walls [m^3/s]
Definition Reactor.h:285
void syncState() override
Set the state of the reactor to correspond to the state of the associated ThermoPhase object.
Definition Reactor.cpp:124
virtual void getSurfaceInitialConditions(double *y)
Get initial conditions for SurfPhase objects attached to this reactor.
Definition Reactor.cpp:78
void initialize(double t0=0.0) override
Initialize the reactor.
Definition Reactor.cpp:87
bool getAdvanceLimits(double *limits) const
Retrieve absolute step size limits during advance.
Definition Reactor.cpp:551
virtual size_t speciesIndex(const string &nm) const
Return the index in the solution vector for this reactor of the species named nm, in either the homog...
Definition Reactor.cpp:429
virtual bool timeIsIndependent() const
Indicates whether the governing equations for this reactor are functions of time or a spatial variabl...
Definition Reactor.h:65
virtual void updateConnected(bool updatePressure)
Update the state information needed by connected reactors, flow devices, and reactor walls.
Definition Reactor.cpp:187
A simple thermodynamic model for a surface phase, assuming an ideal solution model.
Definition SurfPhase.h:98
double size(size_t k) const
Returns the number of sites occupied by one molecule of species k.
Definition SurfPhase.h:221
double siteDensity() const
Returns the site density.
Definition SurfPhase.h:216
Base class for a phase with thermodynamic properties.
virtual double minTemp(size_t k=npos) const
Minimum temperature for which the thermodynamic data for the species or phase are valid.
virtual double maxTemp(size_t k=npos) const
Maximum temperature for which the thermodynamic data for the species are valid.
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.
double intEnergy_mass() const
Specific internal energy. Units: J/kg.
double Hf298SS(const size_t k) const
Report the 298 K Heat of Formation of the standard state of one species (J kmol-1)
double enthalpy_mass() const
Specific enthalpy. Units: J/kg.
Base class for 'walls' (walls, pistons, etc.) connecting reactors.
Definition Wall.h:22
virtual void initialize()
Called just before the start of integration.
Definition Wall.h:93
virtual void setDerivativeSettings(const AnyMap &settings)
Set/modify derivative settings.
Definition Kinetics.h:703
double dot(InputIter x_begin, InputIter x_end, InputIter2 y_begin)
Function that calculates a templated inner product.
Definition utilities.h:82
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:120
Namespace for the Cantera kernel.
Definition AnyMap.cpp:564
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180
offset
Offsets of solution components in the 1D solution array.
Definition StFlow.h:24
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...