Cantera  2.1.2
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"
13 
14 #include <cfloat>
15 
16 using namespace std;
17 
18 namespace Cantera
19 {
20 Reactor::Reactor() : ReactorBase(),
21  m_kin(0),
22  m_vdot(0.0),
23  m_Q(0.0),
24  m_chem(false),
25  m_energy(true),
26  m_nsens(npos)
27 {}
28 
29 void Reactor::getInitialConditions(double t0, size_t leny, double* y)
30 {
31  m_init = true;
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
53  size_t loc = m_nsp + 3;
54  SurfPhase* surf;
55  for (size_t m = 0; m < m_nwalls; m++) {
56  surf = m_wall[m]->surface(m_lr[m]);
57  if (surf) {
58  m_wall[m]->getCoverages(m_lr[m], y + loc);
59  loc += surf->nSpecies();
60  }
61  }
62 }
63 
64 void Reactor::initialize(doublereal t0)
65 {
66  m_thermo->restoreState(m_state);
67  m_sdot.resize(m_nsp, 0.0);
68  m_wdot.resize(m_nsp, 0.0);
69  m_nv = m_nsp + 3;
70  for (size_t w = 0; w < m_nwalls; w++)
71  if (m_wall[w]->surface(m_lr[w])) {
72  m_nv += m_wall[w]->surface(m_lr[w])->nSpecies();
73  }
74 
75  m_enthalpy = m_thermo->enthalpy_mass();
76  m_pressure = m_thermo->pressure();
77  m_intEnergy = m_thermo->intEnergy_mass();
78 
79  size_t nt = 0, maxnt = 0;
80  for (size_t m = 0; m < m_nwalls; m++) {
81  m_wall[m]->initialize();
82  if (m_wall[m]->kinetics(m_lr[m])) {
83  nt = m_wall[m]->kinetics(m_lr[m])->nTotalSpecies();
84  if (nt > maxnt) {
85  maxnt = nt;
86  }
87  if (m_wall[m]->kinetics(m_lr[m])) {
88  if (&m_kin->thermo(0) !=
89  &m_wall[m]->kinetics(m_lr[m])->thermo(0)) {
90  throw CanteraError("Reactor::initialize",
91  "First phase of all kinetics managers must be"
92  " the gas.");
93  }
94  }
95  }
96  }
97  m_work.resize(maxnt);
98  std::sort(m_pnum.begin(), m_pnum.end());
99  m_init = true;
100 }
101 
103 {
104  if (m_nsens == npos) {
105  // determine the number of sensitivity parameters
106  size_t m, ns;
107  m_nsens = m_pnum.size();
108  for (m = 0; m < m_nwalls; m++) {
109  ns = m_wall[m]->nSensParams(m_lr[m]);
110  m_nsens_wall.push_back(ns);
111  m_nsens += ns;
112  }
113  }
114  return m_nsens;
115 }
116 
117 void Reactor::updateState(doublereal* y)
118 {
119  for (size_t i = 0; i < m_nv; i++) {
120  AssertFinite(y[i], "Reactor::updateState",
121  "y[" + int2str(i) + "] is not finite");
122  }
123 
124  // The components of y are [0] the total mass, [1] the total volume,
125  // [2] the total internal energy, [3...K+3] are the mass fractions of each
126  // species, and [K+3...] are the coverages of surface species on each wall.
127  m_mass = y[0];
128  m_vol = y[1];
129 
130  m_thermo->setMassFractions_NoNorm(y+3);
131 
132  if (m_energy) {
133  // Use Newton's method to determine the mixture temperature. Tight
134  // tolerances are required both for Jacobian evaluation and for
135  // sensitivity analysis to work correctly.
136 
137  doublereal U = y[2];
138  doublereal T = temperature();
139  double dT = 100;
140 
141  int i = 0;
142  while (abs(dT / T) > 10 * DBL_EPSILON) {
143  m_thermo->setState_TR(T, m_mass / m_vol);
144  double dUdT = m_thermo->cv_mass() * m_mass;
145  dT = (m_thermo->intEnergy_mass() * m_mass - U) / dUdT;
146  dT = std::min(dT, 0.5 * T);
147  T -= dT;
148  i++;
149  if (i > 100) {
150  std::string message = "no convergence";
151  message += "\nU/m = " + fp2str(U / m_mass);
152  message += "\nT = " + fp2str(T);
153  message += "\nrho = " + fp2str(m_mass / m_vol);
154  message += "\n";
155  throw CanteraError("Reactor::updateState", message);
156  }
157  }
158  } else {
159  m_thermo->setDensity(m_mass/m_vol);
160  }
161 
162  size_t loc = m_nsp + 3;
163  SurfPhase* surf;
164  for (size_t m = 0; m < m_nwalls; m++) {
165  surf = m_wall[m]->surface(m_lr[m]);
166  if (surf) {
167  m_wall[m]->setCoverages(m_lr[m], y+loc);
168  loc += surf->nSpecies();
169  }
170  }
171 
172  // save parameters needed by other connected reactors
173  m_enthalpy = m_thermo->enthalpy_mass();
174  m_pressure = m_thermo->pressure();
175  m_intEnergy = m_thermo->intEnergy_mass();
176  m_thermo->saveState(m_state);
177 }
178 
179 void Reactor::evalEqs(doublereal time, doublereal* y,
180  doublereal* ydot, doublereal* params)
181 {
182  m_thermo->restoreState(m_state);
183 
184  // process sensitivity parameters
185  if (params) {
186  size_t npar = m_pnum.size();
187  for (size_t n = 0; n < npar; n++) {
188  double mult = m_kin->multiplier(m_pnum[n]);
189  m_kin->setMultiplier(m_pnum[n], mult*params[n]);
190  }
191  size_t ploc = npar;
192  for (size_t m = 0; m < m_nwalls; m++) {
193  if (m_nsens_wall[m] > 0) {
194  m_wall[m]->setSensitivityParameters(m_lr[m], params + ploc);
195  ploc += m_nsens_wall[m];
196  }
197  }
198  }
199 
200  m_vdot = 0.0;
201  m_Q = 0.0;
202  double dmdt = 0.0; // dm/dt (gas phase)
203  double* dYdt = ydot + 3;
204 
205  // compute wall terms
206  size_t loc = m_nsp+3;
207  fill(m_sdot.begin(), m_sdot.end(), 0.0);
208  for (size_t i = 0; i < m_nwalls; i++) {
209  int lr = 1 - 2*m_lr[i];
210  double vdot = lr*m_wall[i]->vdot(time);
211  m_vdot += vdot;
212  m_Q += lr*m_wall[i]->Q(time);
213  Kinetics* kin = m_wall[i]->kinetics(m_lr[i]);
214  SurfPhase* surf = m_wall[i]->surface(m_lr[i]);
215  if (surf && kin) {
216  double rs0 = 1.0/surf->siteDensity();
217  size_t nk = surf->nSpecies();
218  double sum = 0.0;
219  surf->setTemperature(m_state[0]);
220  m_wall[i]->syncCoverages(m_lr[i]);
221  kin->getNetProductionRates(DATA_PTR(m_work));
222  size_t ns = kin->surfacePhaseIndex();
223  size_t surfloc = kin->kineticsSpeciesIndex(0,ns);
224  for (size_t k = 1; k < nk; k++) {
225  ydot[loc + k] = m_work[surfloc+k]*rs0*surf->size(k);
226  sum -= ydot[loc + k];
227  }
228  ydot[loc] = sum;
229  loc += nk;
230 
231  double wallarea = m_wall[i]->area();
232  for (size_t k = 0; k < m_nsp; k++) {
233  m_sdot[k] += m_work[k]*wallarea;
234  }
235  }
236  }
237 
238  // volume equation
239  ydot[1] = m_vdot;
240 
241  const vector_fp& mw = m_thermo->molecularWeights();
242  const doublereal* Y = m_thermo->massFractions();
243 
244  if (m_chem) {
245  m_kin->getNetProductionRates(&m_wdot[0]); // "omega dot"
246  }
247 
248  double mdot_surf = 0.0; // net mass flux from surfaces
249  for (size_t k = 0; k < m_nsp; k++) {
250  // production in gas phase and from surfaces
251  dYdt[k] = (m_wdot[k] * m_vol + m_sdot[k]) * mw[k] / m_mass;
252  mdot_surf += m_sdot[k] * mw[k];
253  }
254  dmdt += mdot_surf; // mass added to gas phase from surface reations
255 
256  for (size_t k = 0; k < m_nsp; k++) {
257  // dilution by net surface mass flux
258  dYdt[k] -= Y[k] * mdot_surf / m_mass;
259  }
260 
261  /*
262  * Energy equation.
263  * \f[
264  * \dot U = -P\dot V + A \dot q + \dot m_{in} h_{in}
265  * - \dot m_{out} h.
266  * \f]
267  */
268  if (m_energy) {
269  ydot[2] = - m_thermo->pressure() * m_vdot - m_Q;
270  } else {
271  ydot[2] = 0.0;
272  }
273 
274  // add terms for open system
275  if (m_open) {
276  doublereal enthalpy = m_thermo->enthalpy_mass();
277 
278  // outlets
279  for (size_t i = 0; i < m_nOutlets; i++) {
280  double mdot_out = m_outlet[i]->massFlowRate(time);
281  dmdt -= mdot_out; // mass flow out of system
282  if (m_energy) {
283  ydot[2] -= mdot_out * enthalpy;
284  }
285  }
286 
287  // inlets
288  for (size_t i = 0; i < m_nInlets; i++) {
289  double mdot_in = m_inlet[i]->massFlowRate(time);
290  dmdt += mdot_in; // mass flow into system
291  for (size_t n = 0; n < m_nsp; n++) {
292  double mdot_spec = m_inlet[i]->outletSpeciesMassFlowRate(n);
293  // flow of species into system and dilution by other species
294  dYdt[n] += (mdot_spec - mdot_in * Y[n]) / m_mass;
295  }
296  if (m_energy) {
297  ydot[2] += mdot_in * m_inlet[i]->enthalpy_mass();
298  }
299  }
300  }
301 
302  ydot[0] = dmdt;
303 
304  for (size_t i = 0; i < m_nv; i++) {
305  AssertFinite(ydot[i], "Reactor::evalEqs",
306  "ydot[" + int2str(i) + "] is not finite");
307  }
308 
309  // reset sensitivity parameters
310  if (params) {
311  size_t npar = m_pnum.size();
312  for (size_t n = 0; n < npar; n++) {
313  double mult = m_kin->multiplier(m_pnum[n]);
314  m_kin->setMultiplier(m_pnum[n], mult/params[n]);
315  }
316  size_t ploc = npar;
317  for (size_t m = 0; m < m_nwalls; m++) {
318  if (m_nsens_wall[m] > 0) {
319  m_wall[m]->resetSensitivityParameters(m_lr[m]);
320  ploc += m_nsens_wall[m];
321  }
322  }
323  }
324 }
325 
327 {
328  if (rxn >= m_kin->nReactions())
329  throw CanteraError("Reactor::addSensitivityReaction",
330  "Reaction number out of range ("+int2str(rxn)+")");
331 
333  name()+": "+m_kin->reactionString(rxn));
334  m_pnum.push_back(rxn);
335  m_mult_save.push_back(1.0);
336 }
337 
338 std::vector<std::pair<void*, int> > Reactor::getSensitivityOrder() const
339 {
340  std::vector<std::pair<void*, int> > order;
341  order.push_back(std::make_pair(const_cast<Reactor*>(this), 0));
342  for (size_t n = 0; n < m_nwalls; n++) {
343  if (m_nsens_wall[n]) {
344  order.push_back(std::make_pair(m_wall[n], m_lr[n]));
345  }
346  }
347  return order;
348 }
349 
350 size_t Reactor::componentIndex(const string& nm) const
351 {
352  if (nm == "m") {
353  return 0;
354  }
355  if (nm == "V") {
356  return 1;
357  }
358  if (nm == "U") {
359  return 2;
360  }
361 
362  // check for a gas species name
363  size_t k = m_thermo->speciesIndex(nm);
364  if (k != npos) {
365  return k + 3;
366  }
367 
368  // check for a wall species
369  size_t walloffset = 0, kp = 0;
370  thermo_t* th;
371  for (size_t m = 0; m < m_nwalls; m++) {
372  if (m_wall[m]->kinetics(m_lr[m])) {
373  kp = m_wall[m]->kinetics(m_lr[m])->reactionPhaseIndex();
374  th = &m_wall[m]->kinetics(m_lr[m])->thermo(kp);
375  k = th->speciesIndex(nm);
376  if (k != npos) {
377  return k + 3 + m_nsp + walloffset;
378  } else {
379  walloffset += th->nSpecies();
380  }
381  }
382  }
383  return npos;
384 }
385 
386 }
virtual void getInitialConditions(doublereal t0, size_t leny, doublereal *y)
Called by ReactorNet to get the initial conditions.
Definition: Reactor.cpp:29
virtual doublereal density() const
Density (kg/m^3).
Definition: Phase.h:534
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:40
void restoreState(const vector_fp &state)
Restore a state saved on a previous call to saveState.
Definition: Phase.cpp:289
doublereal m_mass
total mass
Definition: Reactor.h:153
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:102
Header file for class Wall.
vector_fp m_wdot
Species net molar production rates.
Definition: Reactor.h:156
void getMassFractions(doublereal *const y) const
Get the species mass fractions.
Definition: Phase.cpp:561
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:288
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:173
virtual void getNetProductionRates(doublereal *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition: Kinetics.h:600
virtual void addSensitivityReaction(size_t rxn)
Add a sensitivity parameter associated with the reaction number rxn (in the homogeneous phase)...
Definition: Reactor.cpp:326
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:398
doublereal multiplier(size_t i) const
The current value of the multiplier for reaction i.
Definition: Kinetics.h:840
void setMultiplier(size_t i, doublereal f)
Set the multiplier for reaction i to f.
Definition: Kinetics.h:849
#define AssertFinite(expr, procedure, message)
Throw an exception if the specified exception is not a finite number.
Definition: ctexceptions.h:254
virtual void initialize(doublereal t0=0.0)
Initialize the reactor.
Definition: Reactor.cpp:64
doublereal intEnergy_mass() const
Specific internal energy.
Definition: ThermoPhase.h:940
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:101
virtual void updateState(doublereal *y)
Set the state of the reactor to correspond to the state vector y.
Definition: Reactor.cpp:117
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
Definition: Reactor.h:149
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:229
size_t surfacePhaseIndex()
This returns the integer index of the phase which has ThermoPhase type cSurf.
Definition: Kinetics.h:263
std::string fp2str(const double x, const std::string &fmt)
Convert a double into a c++ string.
Definition: stringUtils.cpp:29
const doublereal * massFractions() const
Return a const pointer to the mass fraction array.
Definition: Phase.h:454
Public interface for kinetics managers.
Definition: Kinetics.h:131
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:68
doublereal cv_mass() const
Specific heat at constant volume.
Definition: ThermoPhase.h:968
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition: Kinetics.h:331
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:388
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:270
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:252
virtual doublereal pressure() const
Return the thermodynamic pressure (Pa).
Definition: ThermoPhase.h:314
const vector_fp & molecularWeights() const
Return a const reference to the internal vector of molecular weights.
Definition: Phase.cpp:505
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:165
virtual void setTemperature(const doublereal temp)
Set the internally stored temperature of the phase (K).
Definition: Phase.h:562
doublereal enthalpy_mass() const
Specific enthalpy.
Definition: ThermoPhase.h:933
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition: Kinetics.h:196
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
Definition: ct_defs.h:36
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:338
size_t m_nsp
Number of homogeneous species in the mixture.
Definition: ReactorBase.h:213
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:277
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:350
virtual std::string reactionString(size_t i) const
Return a string representing the reaction.
Definition: Kinetics.h:716
virtual void evalEqs(doublereal t, doublereal *y, doublereal *ydot, doublereal *params)
Definition: Reactor.cpp:179
doublereal m_vdot
Tolerance on the temperature.
Definition: Reactor.h:152
void setState_TR(doublereal t, doublereal rho)
Set the internally stored temperature (K) and density (kg/m^3)
Definition: Phase.cpp:454
doublereal temperature() const
Returns the current temperature (K) of the reactor's contents.
Definition: ReactorBase.h:164
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:549
doublereal siteDensity()
Returns the site density.
Definition: SurfPhase.h:415