Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Reactor.cpp
Go to the documentation of this file.
1 /**
2  * @file Reactor.cpp A zero-dimensional reactor
3  */
4 
5 // Copyright 2001 California Institute of Technology
6 
9 #include "cantera/zeroD/Wall.h"
12 
13 #include <cfloat>
14 
15 using namespace std;
16 
17 namespace Cantera
18 {
19 Reactor::Reactor() :
20  m_kin(0),
21  m_vdot(0.0),
22  m_Q(0.0),
23  m_mass(0.0),
24  m_chem(false),
25  m_energy(true),
26  m_nv(0),
27  m_nsens(npos)
28 {}
29 
30 void Reactor::getInitialConditions(double t0, size_t leny, double* y)
31 {
32  if (m_thermo == 0) {
33  cout << "Error: reactor is empty." << endl;
34  return;
35  }
36  m_thermo->restoreState(m_state);
37 
38  // set the first component to the total mass
39  m_mass = m_thermo->density() * m_vol;
40  y[0] = m_mass;
41 
42  // set the second component to the total volume
43  y[1] = m_vol;
44 
45  // set the third component to the total internal energy
46  y[2] = m_thermo->intEnergy_mass() * m_mass;
47 
48  // set components y+3 ... y+K+2 to the mass fractions of each species
49  m_thermo->getMassFractions(y+3);
50 
51  // set the remaining components to the surface species
52  // coverages on the walls
54 }
55 
57 {
58  size_t loc = 0;
59  for (size_t m = 0; m < m_wall.size(); m++) {
60  SurfPhase* surf = m_wall[m]->surface(m_lr[m]);
61  if (surf) {
62  m_wall[m]->getCoverages(m_lr[m], y + loc);
63  loc += surf->nSpecies();
64  }
65  }
66 }
67 
68 void Reactor::initialize(doublereal t0)
69 {
70  if (!m_thermo || !m_kin) {
71  throw CanteraError("Reactor::initialize", "Reactor contents not set"
72  " for reactor '" + m_name + "'.");
73  }
74  m_thermo->restoreState(m_state);
75  m_sdot.resize(m_nsp, 0.0);
76  m_wdot.resize(m_nsp, 0.0);
77  m_nv = m_nsp + 3;
78  for (size_t w = 0; w < m_wall.size(); w++)
79  if (m_wall[w]->surface(m_lr[w])) {
80  m_nv += m_wall[w]->surface(m_lr[w])->nSpecies();
81  }
82 
83  m_enthalpy = m_thermo->enthalpy_mass();
84  m_pressure = m_thermo->pressure();
85  m_intEnergy = m_thermo->intEnergy_mass();
86 
87  size_t nt = 0, maxnt = 0;
88  for (size_t m = 0; m < m_wall.size(); m++) {
89  m_wall[m]->initialize();
90  if (m_wall[m]->kinetics(m_lr[m])) {
91  nt = m_wall[m]->kinetics(m_lr[m])->nTotalSpecies();
92  maxnt = std::max(maxnt, nt);
93  if (m_wall[m]->kinetics(m_lr[m])) {
94  if (&m_kin->thermo(0) !=
95  &m_wall[m]->kinetics(m_lr[m])->thermo(0)) {
96  throw CanteraError("Reactor::initialize",
97  "First phase of all kinetics managers must be"
98  " the gas.");
99  }
100  }
101  }
102  }
103  m_work.resize(maxnt);
104  std::sort(m_pnum.begin(), m_pnum.end());
105 }
106 
108 {
109  if (m_nsens == npos) {
110  // determine the number of sensitivity parameters
111  size_t m, ns;
112  m_nsens = m_pnum.size();
113  for (m = 0; m < m_wall.size(); m++) {
114  ns = m_wall[m]->nSensParams(m_lr[m]);
115  m_nsens_wall.push_back(ns);
116  m_nsens += ns;
117  }
118  }
119  return m_nsens;
120 }
121 
123 {
125  m_mass = m_thermo->density() * m_vol;
126 }
127 
128 void Reactor::updateState(doublereal* y)
129 {
130  for (size_t i = 0; i < m_nv; i++) {
131  AssertFinite(y[i], "Reactor::updateState",
132  "y[" + int2str(i) + "] is not finite");
133  }
134 
135  // The components of y are [0] the total mass, [1] the total volume,
136  // [2] the total internal energy, [3...K+3] are the mass fractions of each
137  // species, and [K+3...] are the coverages of surface species on each wall.
138  m_mass = y[0];
139  m_vol = y[1];
140 
141  m_thermo->setMassFractions_NoNorm(y+3);
142 
143  if (m_energy) {
144  // Use a damped Newton's method to determine the mixture temperature.
145  // Tight tolerances are required both for Jacobian evaluation and for
146  // sensitivity analysis to work correctly.
147 
148  doublereal U = y[2];
149  doublereal T = temperature();
150  double dT = 100;
151  double dUprev = 1e10;
152  double dU = 1e10;
153 
154  int i = 0;
155  double damp = 1.0;
156  while (abs(dT / T) > 10 * DBL_EPSILON) {
157  dUprev = dU;
158  m_thermo->setState_TR(T, m_mass / m_vol);
159  double dUdT = m_thermo->cv_mass() * m_mass;
160  dU = m_thermo->intEnergy_mass() * m_mass - U;
161  dT = dU / dUdT;
162  // Reduce the damping coefficient if the magnitude of the error
163  // isn't decreasing
164  if (std::abs(dU) < std::abs(dUprev)) {
165  damp = 1.0;
166  } else {
167  damp *= 0.8;
168  }
169  dT = std::min(dT, 0.5 * T) * damp;
170  T -= dT;
171  i++;
172  if (i > 100) {
173  std::string message = "no convergence";
174  message += "\nU/m = " + fp2str(U / m_mass);
175  message += "\nT = " + fp2str(T);
176  message += "\nrho = " + fp2str(m_mass / m_vol);
177  message += "\n";
178  throw CanteraError("Reactor::updateState", message);
179  }
180  }
181  } else {
182  m_thermo->setDensity(m_mass/m_vol);
183  }
184 
185  updateSurfaceState(y + m_nsp + 3);
186 
187  // save parameters needed by other connected reactors
188  m_enthalpy = m_thermo->enthalpy_mass();
189  m_pressure = m_thermo->pressure();
190  m_intEnergy = m_thermo->intEnergy_mass();
191  m_thermo->saveState(m_state);
192 }
193 
195 {
196  size_t loc = 0;
197  for (size_t m = 0; m < m_wall.size(); m++) {
198  SurfPhase* surf = m_wall[m]->surface(m_lr[m]);
199  if (surf) {
200  m_wall[m]->setCoverages(m_lr[m], y+loc);
201  loc += surf->nSpecies();
202  }
203  }
204 }
205 
206 void Reactor::evalEqs(doublereal time, doublereal* y,
207  doublereal* ydot, doublereal* params)
208 {
209  double dmdt = 0.0; // dm/dt (gas phase)
210  double* dYdt = ydot + 3;
211 
212  m_thermo->restoreState(m_state);
213  applySensitivity(params);
214  evalWalls(time);
215  double mdot_surf = evalSurfaces(time, ydot + m_nsp + 3);
216  dmdt += mdot_surf; // mass added to gas phase from surface reations
217 
218  // volume equation
219  ydot[1] = m_vdot;
220 
221  const vector_fp& mw = m_thermo->molecularWeights();
222  const doublereal* Y = m_thermo->massFractions();
223 
224  if (m_chem) {
225  m_kin->getNetProductionRates(&m_wdot[0]); // "omega dot"
226  }
227 
228  for (size_t k = 0; k < m_nsp; k++) {
229  // production in gas phase and from surfaces
230  dYdt[k] = (m_wdot[k] * m_vol + m_sdot[k]) * mw[k] / m_mass;
231  // dilution by net surface mass flux
232  dYdt[k] -= Y[k] * mdot_surf / m_mass;
233  }
234 
235  /*
236  * Energy equation.
237  * \f[
238  * \dot U = -P\dot V + A \dot q + \dot m_{in} h_{in}
239  * - \dot m_{out} h.
240  * \f]
241  */
242  if (m_energy) {
243  ydot[2] = - m_thermo->pressure() * m_vdot - m_Q;
244  } else {
245  ydot[2] = 0.0;
246  }
247 
248  // add terms for outlets
249  for (size_t i = 0; i < m_outlet.size(); i++) {
250  double mdot_out = m_outlet[i]->massFlowRate(time);
251  dmdt -= mdot_out; // mass flow out of system
252  if (m_energy) {
253  ydot[2] -= mdot_out * m_enthalpy;
254  }
255  }
256 
257  // add terms for inlets
258  for (size_t i = 0; i < m_inlet.size(); i++) {
259  double mdot_in = m_inlet[i]->massFlowRate(time);
260  dmdt += mdot_in; // mass flow into system
261  for (size_t n = 0; n < m_nsp; n++) {
262  double mdot_spec = m_inlet[i]->outletSpeciesMassFlowRate(n);
263  // flow of species into system and dilution by other species
264  dYdt[n] += (mdot_spec - mdot_in * Y[n]) / m_mass;
265  }
266  if (m_energy) {
267  ydot[2] += mdot_in * m_inlet[i]->enthalpy_mass();
268  }
269  }
270 
271  ydot[0] = dmdt;
272 
273  for (size_t i = 0; i < m_nv; i++) {
274  AssertFinite(ydot[i], "Reactor::evalEqs",
275  "ydot[" + int2str(i) + "] is not finite");
276  }
277 
278  resetSensitivity(params);
279 }
280 
281 void Reactor::evalWalls(double t)
282 {
283  m_vdot = 0.0;
284  m_Q = 0.0;
285  for (size_t i = 0; i < m_wall.size(); i++) {
286  int lr = 1 - 2*m_lr[i];
287  m_vdot += lr*m_wall[i]->vdot(t);
288  m_Q += lr*m_wall[i]->Q(t);
289  }
290 }
291 
292 double Reactor::evalSurfaces(double t, double* ydot)
293 {
294  const vector_fp& mw = m_thermo->molecularWeights();
295 
296  fill(m_sdot.begin(), m_sdot.end(), 0.0);
297  size_t loc = 0; // offset into ydot
298  double mdot_surf = 0.0; // net mass flux from surface
299 
300  for (size_t i = 0; i < m_wall.size(); i++) {
301  Kinetics* kin = m_wall[i]->kinetics(m_lr[i]);
302  SurfPhase* surf = m_wall[i]->surface(m_lr[i]);
303  if (surf && kin) {
304  double rs0 = 1.0/surf->siteDensity();
305  size_t nk = surf->nSpecies();
306  double sum = 0.0;
307  surf->setTemperature(m_state[0]);
308  m_wall[i]->syncCoverages(m_lr[i]);
309  kin->getNetProductionRates(&m_work[0]);
310  size_t ns = kin->surfacePhaseIndex();
311  size_t surfloc = kin->kineticsSpeciesIndex(0,ns);
312  for (size_t k = 1; k < nk; k++) {
313  ydot[loc + k] = m_work[surfloc+k]*rs0*surf->size(k);
314  sum -= ydot[loc + k];
315  }
316  ydot[loc] = sum;
317  loc += nk;
318 
319  double wallarea = m_wall[i]->area();
320  for (size_t k = 0; k < m_nsp; k++) {
321  m_sdot[k] += m_work[k]*wallarea;
322  mdot_surf += m_sdot[k] * mw[k];
323  }
324  }
325  }
326  return mdot_surf;
327 }
328 
330 {
331  if (rxn >= m_kin->nReactions())
332  throw CanteraError("Reactor::addSensitivityReaction",
333  "Reaction number out of range ("+int2str(rxn)+")");
334 
336  name()+": "+m_kin->reactionString(rxn));
337  m_pnum.push_back(rxn);
338  m_mult_save.push_back(1.0);
339 }
340 
341 std::vector<std::pair<void*, int> > Reactor::getSensitivityOrder() const
342 {
343  std::vector<std::pair<void*, int> > order;
344  order.push_back(std::make_pair(const_cast<Reactor*>(this), 0));
345  for (size_t n = 0; n < m_wall.size(); n++) {
346  if (m_nsens_wall[n]) {
347  order.push_back(std::make_pair(m_wall[n], m_lr[n]));
348  }
349  }
350  return order;
351 }
352 
353 size_t Reactor::speciesIndex(const string& nm) const
354 {
355  // check for a gas species name
356  size_t k = m_thermo->speciesIndex(nm);
357  if (k != npos) {
358  return k;
359  }
360 
361  // check for a wall species
362  size_t walloffset = 0, kp = 0;
363  thermo_t* th;
364  for (size_t m = 0; m < m_wall.size(); m++) {
365  if (m_wall[m]->kinetics(m_lr[m])) {
366  kp = m_wall[m]->kinetics(m_lr[m])->reactionPhaseIndex();
367  th = &m_wall[m]->kinetics(m_lr[m])->thermo(kp);
368  k = th->speciesIndex(nm);
369  if (k != npos) {
370  return k + m_nsp + walloffset;
371  } else {
372  walloffset += th->nSpecies();
373  }
374  }
375  }
376  return npos;
377 }
378 
379 size_t Reactor::componentIndex(const string& nm) const
380 {
381  size_t k = speciesIndex(nm);
382  if (k != npos) {
383  return k + 3;
384  } else if (nm == "m" || nm == "mass") {
385  return 0;
386  } else if (nm == "V" || nm == "volume") {
387  return 1;
388  } else if (nm == "U" || nm == "int_energy") {
389  return 2;
390  } else {
391  return npos;
392  }
393 }
394 
395 void Reactor::applySensitivity(double* params)
396 {
397  if (!params) {
398  return;
399  }
400  size_t npar = m_pnum.size();
401  for (size_t n = 0; n < npar; n++) {
402  double mult = m_kin->multiplier(m_pnum[n]);
403  m_kin->setMultiplier(m_pnum[n], mult*params[n]);
404  }
405  size_t ploc = npar;
406  for (size_t m = 0; m < m_wall.size(); m++) {
407  if (m_nsens_wall[m] > 0) {
408  m_wall[m]->setSensitivityParameters(m_lr[m], params + ploc);
409  ploc += m_nsens_wall[m];
410  }
411  }
412 
413 }
414 
415 void Reactor::resetSensitivity(double* params)
416 {
417  if (!params) {
418  return;
419  }
420  size_t npar = m_pnum.size();
421  for (size_t n = 0; n < npar; n++) {
422  double mult = m_kin->multiplier(m_pnum[n]);
423  m_kin->setMultiplier(m_pnum[n], mult/params[n]);
424  }
425  size_t ploc = npar;
426  for (size_t m = 0; m < m_wall.size(); m++) {
427  if (m_nsens_wall[m] > 0) {
428  m_wall[m]->resetSensitivityParameters(m_lr[m]);
429  ploc += m_nsens_wall[m];
430  }
431  }
432 }
433 
434 }
const std::string & reactionString(size_t i) const
Return a string representing the reaction.
Definition: Kinetics.h:706
virtual void getInitialConditions(doublereal t0, size_t leny, doublereal *y)
Called by ReactorNet to get the initial conditions.
Definition: Reactor.cpp:30
virtual doublereal density() const
Density (kg/m^3).
Definition: Phase.h:608
vector_fp m_sdot
Production rates of gas phase species on surfaces [kmol/s].
Definition: Reactor.h:189
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:39
void restoreState(const vector_fp &state)
Restore a state saved on a previous call to saveState.
Definition: Phase.cpp:314
doublereal m_mass
total mass
Definition: Reactor.h:185
Header for a simple thermodynamics model of a surface phase derived from ThermoPhase, assuming an ideal solution model (see Thermodynamic Properties and class SurfPhase).
virtual size_t nSensParams()
Number of sensitivity parameters associated with this reactor (including walls)
Definition: Reactor.cpp:107
Header file for class Wall.
vector_fp m_wdot
Species net molar production rates.
Definition: Reactor.h:191
void getMassFractions(doublereal *const y) const
Get the species mass fractions.
Definition: Phase.cpp:598
thermo_t & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism...
Definition: Kinetics.h:285
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:165
virtual void addSensitivityReaction(size_t rxn)
Add a sensitivity parameter associated with the reaction number rxn (in the homogeneous phase)...
Definition: Reactor.cpp:329
ReactorNet & network()
The ReactorNet that this reactor belongs to.
Definition: ReactorBase.cpp:72
doublereal size(size_t k) const
This routine returns the size of species k.
Definition: Phase.h:411
virtual void evalWalls(double t)
Evaluate terms related to Walls Calculates m_vdot and m_Q based on wall movement and heat transfer...
Definition: Reactor.cpp:281
virtual void resetSensitivity(double *params)
Reset the reaction rate multipliers.
Definition: Reactor.cpp:415
doublereal multiplier(size_t i) const
The current value of the multiplier for reaction i.
Definition: Kinetics.h:875
#define AssertFinite(expr, procedure, message)
Throw an exception if the specified exception is not a finite number.
Definition: ctexceptions.h:308
virtual void initialize(doublereal t0=0.0)
Initialize the reactor.
Definition: Reactor.cpp:68
doublereal intEnergy_mass() const
Specific internal energy.
Definition: ThermoPhase.h:899
virtual void getNetProductionRates(doublereal *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition: Kinetics.cpp:501
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:97
virtual void updateState(doublereal *y)
Set the state of the reactor to correspond to the state vector y.
Definition: Reactor.cpp:128
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
Definition: Reactor.h:181
A simple thermodynamic model for a surface phase, assuming an ideal solution model.
Definition: SurfPhase.h:143
size_t speciesIndex(const std::string &name) const
Returns the index of a species named 'name' within the Phase object.
Definition: Phase.cpp:257
size_t surfacePhaseIndex()
This returns the integer index of the phase which has ThermoPhase type cSurf.
Definition: Kinetics.h:260
virtual void applySensitivity(double *params)
Set reaction rate multipliers based on the sensitivity variables in params.
Definition: Reactor.cpp:395
std::string fp2str(const double x, const std::string &fmt)
Convert a double into a c++ string.
Definition: stringUtils.cpp:28
doublereal m_Q
net heat transfer through walls [W]
Definition: Reactor.h:184
virtual void syncState()
Set the state of the reactor to correspond to the state of the associated ThermoPhase object...
Definition: ReactorBase.cpp:36
const doublereal * massFractions() const
Return a const pointer to the mass fraction array.
Definition: Phase.h:482
virtual double evalSurfaces(double t, double *ydot)
Evaluate terms related to surface reactions Calculates m_sdot and rate of change in surface species c...
Definition: Reactor.cpp:292
Public interface for kinetics managers.
Definition: Kinetics.h:128
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:99
virtual size_t speciesIndex(const std::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:353
doublereal cv_mass() const
Specific heat at constant volume.
Definition: ThermoPhase.h:927
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition: Kinetics.h:323
std::string name() const
Return the name of this reactor.
Definition: ReactorBase.h:42
virtual void setMassFractions_NoNorm(const doublereal *const y)
Set the mass fractions to the specified values without normalizing.
Definition: Phase.cpp:404
void registerSensitivityReaction(void *reactor, size_t reactionIndex, const std::string &name, int leftright=0)
Used by Reactor and Wall objects to register the addition of sensitivity reactions so that the Reacto...
Definition: ReactorNet.cpp:234
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:265
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
Definition: ThermoPhase.h:278
const vector_fp & molecularWeights() const
Return a const reference to the internal vector of molecular weights.
Definition: Phase.cpp:515
virtual void syncState()
Set the state of the reactor to correspond to the state of the associated ThermoPhase object...
Definition: Reactor.cpp:122
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:157
virtual void setTemperature(const doublereal temp)
Set the internally stored temperature of the phase (K).
Definition: Phase.h:638
doublereal enthalpy_mass() const
Specific enthalpy.
Definition: ThermoPhase.h:892
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition: Kinetics.h:193
std::vector< std::pair< void *, int > > getSensitivityOrder() const
Return a vector specifying the ordering of objects to use when determining sensitivity parameter indi...
Definition: Reactor.cpp:341
size_t m_nsp
Number of homogeneous species in the mixture.
Definition: ReactorBase.h:217
void saveState(vector_fp &state) const
Save the current internal state of the phase Write to vector 'state' the current internal state...
Definition: Phase.cpp:302
virtual size_t componentIndex(const std::string &nm) const
Return the index in the solution vector for this reactor of the component named nm.
Definition: Reactor.cpp:379
virtual void setMultiplier(size_t i, doublereal f)
Set the multiplier for reaction i to f.
Definition: Kinetics.h:884
virtual void evalEqs(doublereal t, doublereal *y, doublereal *ydot, doublereal *params)
Definition: Reactor.cpp:206
doublereal m_vdot
net rate of volume change from moving walls [m^3/s]
Definition: Reactor.h:183
void setState_TR(doublereal t, doublereal rho)
Set the internally stored temperature (K) and density (kg/m^3)
Definition: Phase.cpp:464
virtual void updateSurfaceState(double *y)
Update the state of SurfPhase objects attached to this reactor.
Definition: Reactor.cpp:194
doublereal temperature() const
Returns the current temperature (K) of the reactor's contents.
Definition: ReactorBase.h:173
virtual void setDensity(const doublereal density_)
Set the internally stored density (kg/m^3) of the phase Note the density of a phase is an independent...
Definition: Phase.h:623
virtual void getSurfaceInitialConditions(double *y)
Get initial conditions for SurfPhase objects attached to this reactor.
Definition: Reactor.cpp:56
doublereal siteDensity()
Returns the site density.
Definition: SurfPhase.h:410