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