Cantera  2.1.2
ConstPressureReactor.cpp
Go to the documentation of this file.
1 /**
2  * @file ConstPressureReactor.cpp A constant pressure zero-dimensional
3  * reactor
4  */
5 
6 // Copyright 2001 California Institute of Technology
7 
10 #include "cantera/zeroD/Wall.h"
13 
14 using namespace std;
15 
16 namespace Cantera
17 {
18 
19 ConstPressureReactor::ConstPressureReactor() : Reactor() {}
20 
22 getInitialConditions(double t0, size_t leny, double* y)
23 {
24  m_init = true;
25  if (m_thermo == 0) {
26  throw CanteraError("getInitialConditions",
27  "Error: reactor is empty.");
28  }
29  m_thermo->restoreState(m_state);
30 
31  // set the first component to the total mass
32  y[0] = m_thermo->density() * m_vol;
33 
34  // set the second component to the total enthalpy
35  y[1] = m_thermo->enthalpy_mass() * m_thermo->density() * m_vol;
36 
37  // set components y+2 ... y+K+1 to the mass fractions Y_k of each species
38  m_thermo->getMassFractions(y+2);
39 
40  // set the remaining components to the surface species
41  // coverages on the walls
42  size_t loc = m_nsp + 2;
43  SurfPhase* surf;
44  for (size_t m = 0; m < m_nwalls; m++) {
45  surf = m_wall[m]->surface(m_lr[m]);
46  if (surf) {
47  m_wall[m]->getCoverages(m_lr[m], y + loc);
48  loc += surf->nSpecies();
49  }
50  }
51 }
52 
54 {
55  m_thermo->restoreState(m_state);
56  m_sdot.resize(m_nsp, 0.0);
57  m_wdot.resize(m_nsp, 0.0);
58  m_nv = m_nsp + 2;
59  for (size_t w = 0; w < m_nwalls; w++)
60  if (m_wall[w]->surface(m_lr[w])) {
61  m_nv += m_wall[w]->surface(m_lr[w])->nSpecies();
62  }
63 
64  m_enthalpy = m_thermo->enthalpy_mass();
65  m_pressure = m_thermo->pressure();
66  m_intEnergy = m_thermo->intEnergy_mass();
67 
68  size_t nt = 0, maxnt = 0;
69  for (size_t m = 0; m < m_nwalls; m++) {
70  if (m_wall[m]->kinetics(m_lr[m])) {
71  nt = m_wall[m]->kinetics(m_lr[m])->nTotalSpecies();
72  if (nt > maxnt) {
73  maxnt = nt;
74  }
75  if (m_wall[m]->kinetics(m_lr[m])) {
76  if (&m_kin->thermo(0) !=
77  &m_wall[m]->kinetics(m_lr[m])->thermo(0)) {
78  throw CanteraError("ConstPressureReactor::initialize",
79  "First phase of all kinetics managers must be"
80  " the gas.");
81  }
82  }
83  }
84  }
85  m_work.resize(maxnt);
86  std::sort(m_pnum.begin(), m_pnum.end());
87  m_init = true;
88 }
89 
91 {
92  // The components of y are [0] the total mass, [1] the total enthalpy,
93  // [2...K+2) are the mass fractions of each species, and [K+2...] are the
94  // coverages of surface species on each wall.
95  m_mass = y[0];
96  m_thermo->setMassFractions_NoNorm(y+2);
97  if (m_energy) {
98  m_thermo->setState_HP(y[1]/m_mass, m_pressure, 1.0e-4);
99  } else {
100  m_thermo->setPressure(m_pressure);
101  }
102  m_vol = m_mass / m_thermo->density();
103 
104  size_t loc = m_nsp + 2;
105  SurfPhase* surf;
106  for (size_t m = 0; m < m_nwalls; m++) {
107  surf = m_wall[m]->surface(m_lr[m]);
108  if (surf) {
109  m_wall[m]->setCoverages(m_lr[m], y+loc);
110  loc += surf->nSpecies();
111  }
112  }
113 
114  // save parameters needed by other connected reactors
115  m_enthalpy = m_thermo->enthalpy_mass();
116  m_intEnergy = m_thermo->intEnergy_mass();
117  m_thermo->saveState(m_state);
118 }
119 
120 void ConstPressureReactor::evalEqs(doublereal time, doublereal* y,
121  doublereal* ydot, doublereal* params)
122 {
123  size_t nk;
124  m_thermo->restoreState(m_state);
125 
126  Kinetics* kin;
127  size_t npar, ploc;
128  double mult;
129 
130  // process sensitivity parameters
131  if (params) {
132 
133  npar = m_pnum.size();
134  for (size_t n = 0; n < npar; n++) {
135  mult = m_kin->multiplier(m_pnum[n]);
136  m_kin->setMultiplier(m_pnum[n], mult*params[n]);
137  }
138  ploc = npar;
139  for (size_t m = 0; m < m_nwalls; m++) {
140  if (m_nsens_wall[m] > 0) {
141  m_wall[m]->setSensitivityParameters(m_lr[m], params + ploc);
142  ploc += m_nsens_wall[m];
143  }
144  }
145  }
146 
147  m_Q = 0.0;
148 
149  // compute wall terms
150  doublereal rs0, sum, wallarea;
151  double dmdt = 0.0; // dm/dt (gas phase)
152  double* dYdt = ydot + 2;
153 
154  SurfPhase* surf;
155  size_t lr, ns, loc = m_nsp+2, surfloc;
156  fill(m_sdot.begin(), m_sdot.end(), 0.0);
157  for (size_t i = 0; i < m_nwalls; i++) {
158  lr = 1 - 2*m_lr[i];
159  m_Q += lr*m_wall[i]->Q(time);
160  kin = m_wall[i]->kinetics(m_lr[i]);
161  surf = m_wall[i]->surface(m_lr[i]);
162  if (surf && kin) {
163  rs0 = 1.0/surf->siteDensity();
164  nk = surf->nSpecies();
165  sum = 0.0;
166  surf->setTemperature(m_state[0]);
167  m_wall[i]->syncCoverages(m_lr[i]);
168  kin->getNetProductionRates(DATA_PTR(m_work));
169  ns = kin->surfacePhaseIndex();
170  surfloc = kin->kineticsSpeciesIndex(0,ns);
171  for (size_t k = 1; k < nk; k++) {
172  ydot[loc + k] = m_work[surfloc+k]*rs0*surf->size(k);
173  sum -= ydot[loc + k];
174  }
175  ydot[loc] = sum;
176  loc += nk;
177 
178  wallarea = m_wall[i]->area();
179  for (size_t k = 0; k < m_nsp; k++) {
180  m_sdot[k] += m_work[k]*wallarea;
181  }
182  }
183  }
184 
185  const vector_fp& mw = m_thermo->molecularWeights();
186  const doublereal* Y = m_thermo->massFractions();
187 
188  if (m_chem) {
189  m_kin->getNetProductionRates(&m_wdot[0]); // "omega dot"
190  }
191 
192  double mdot_surf = 0.0; // net mass flux from surface
193  for (size_t k = 0; k < m_nsp; k++) {
194  // production in gas phase and from surfaces
195  dYdt[k] = (m_wdot[k] * m_vol + m_sdot[k]) * mw[k] / m_mass;
196  mdot_surf += m_sdot[k] * mw[k];
197  }
198  dmdt += mdot_surf;
199 
200  // external heat transfer
201  double dHdt = - m_Q;
202 
203  for (size_t n = 0; n < m_nsp; n++) {
204  // dilution by net surface mass flux
205  dYdt[n] -= Y[n] * mdot_surf / m_mass;
206  }
207 
208  // add terms for open system
209  if (m_open) {
210  double enthalpy = m_thermo->enthalpy_mass();
211  // outlets
212  for (size_t i = 0; i < m_nOutlets; i++) {
213  double mdot_out = m_outlet[i]->massFlowRate(time); // mass flow out of system
214  dmdt -= mdot_out;
215  dHdt -= mdot_out * enthalpy;
216  }
217 
218  // inlets
219  for (size_t i = 0; i < m_nInlets; i++) {
220  double mdot_in = m_inlet[i]->massFlowRate(time);
221  dmdt += mdot_in; // mass flow into system
222  for (size_t n = 0; n < m_nsp; n++) {
223  double mdot_spec = m_inlet[i]->outletSpeciesMassFlowRate(n);
224  // flow of species into system and dilution by other species
225  dYdt[n] += (mdot_spec - mdot_in * Y[n]) / m_mass;
226  }
227  dHdt += mdot_in * m_inlet[i]->enthalpy_mass();
228  }
229  }
230 
231  ydot[0] = dmdt;
232  if (m_energy) {
233  ydot[1] = dHdt;
234  } else {
235  ydot[1] = 0.0;
236  }
237 
238  // reset sensitivity parameters
239  if (params) {
240  npar = m_pnum.size();
241  for (size_t n = 0; n < npar; n++) {
242  mult = m_kin->multiplier(m_pnum[n]);
243  m_kin->setMultiplier(m_pnum[n], mult/params[n]);
244  }
245  ploc = npar;
246  for (size_t m = 0; m < m_nwalls; m++) {
247  if (m_nsens_wall[m] > 0) {
248  m_wall[m]->resetSensitivityParameters(m_lr[m]);
249  ploc += m_nsens_wall[m];
250  }
251  }
252  }
253 }
254 
255 size_t ConstPressureReactor::componentIndex(const string& nm) const
256 {
257  if (nm == "m") {
258  return 0;
259  }
260  if (nm == "H") {
261  return 1;
262  }
263  // check for a gas species name
264  size_t k = m_thermo->speciesIndex(nm);
265  if (k != npos) {
266  return k + 2;
267  }
268 
269  // check for a wall species
270  size_t walloffset = 0, kp = 0;
271  thermo_t* th;
272  for (size_t m = 0; m < m_nwalls; m++) {
273  if (m_wall[m]->kinetics(m_lr[m])) {
274  kp = m_wall[m]->kinetics(m_lr[m])->reactionPhaseIndex();
275  th = &m_wall[m]->kinetics(m_lr[m])->thermo(kp);
276  k = th->speciesIndex(nm);
277  if (k != npos) {
278  return k + 2 + m_nsp + walloffset;
279  } else {
280  walloffset += th->nSpecies();
281  }
282  }
283  }
284  return npos;
285 }
286 
287 }
virtual void setState_HP(doublereal h, doublereal p, doublereal tol=1.e-4)
Set the internally stored specific enthalpy (J/kg) and pressure (Pa) of the phase.
virtual doublereal density() const
Density (kg/m^3).
Definition: Phase.h:534
void restoreState(const vector_fp &state)
Restore a state saved on a previous call to saveState.
Definition: Phase.cpp:289
virtual size_t componentIndex(const std::string &nm) const
Return the index in the solution vector for this reactor of the component named nm.
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).
Header file for class Wall.
vector_fp m_wdot
Species net molar production rates.
Definition: Reactor.h:156
virtual void getInitialConditions(doublereal t0, size_t leny, doublereal *y)
Called by ReactorNet to get the initial conditions.
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
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
virtual void evalEqs(doublereal t, doublereal *y, doublereal *ydot, doublereal *params)
void setMultiplier(size_t i, doublereal f)
Set the multiplier for reaction i to f.
Definition: Kinetics.h:849
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.
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
virtual void initialize(doublereal t0=0.0)
Initialize the reactor.
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
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
virtual void setMassFractions_NoNorm(const doublereal *const y)
Set the mass fractions to the specified values without normalizing.
Definition: Phase.cpp:388
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
virtual void setPressure(doublereal p)
Set the internally stored pressure (Pa) at constant temperature and composition.
Definition: ThermoPhase.h:331
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
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
Definition: ct_defs.h:36
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
doublereal siteDensity()
Returns the site density.
Definition: SurfPhase.h:415