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