Cantera 2.6.0
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
15
16#include <boost/math/tools/roots.hpp>
17
18using namespace std;
19namespace bmt = boost::math::tools;
20
21namespace Cantera
22{
23Reactor::Reactor() :
24 m_kin(0),
25 m_vdot(0.0),
26 m_Q(0.0),
27 m_Qdot(0.0),
28 m_mass(0.0),
29 m_chem(false),
30 m_energy(true),
31 m_nv(0)
32{}
33
34void Reactor::insert(shared_ptr<Solution> sol) {
35 setThermoMgr(*sol->thermo());
36 setKineticsMgr(*sol->kinetics());
37}
38
39void 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
49void 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
75void 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
84void Reactor::initialize(doublereal 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
112size_t Reactor::nSensParams()
113{
114 size_t ns = m_sensParams.size();
115 for (auto& S : m_surfaces) {
116 ns += S->nSensParams();
117 }
118 return ns;
119}
120
121void Reactor::syncState()
122{
123 ReactorBase::syncState();
124 m_mass = m_thermo->density() * m_vol;
125}
126
127void Reactor::updateState(doublereal* 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_TR(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 std::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_TR(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_TR(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
175void 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
184void 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 = (m_net != nullptr) ? m_net->time() : 0.0;
195 for (size_t i = 0; i < m_outlet.size(); i++) {
196 m_outlet[i]->updateMassFlowRate(time);
197 }
198 for (size_t i = 0; i < m_inlet.size(); i++) {
199 m_inlet[i]->updateMassFlowRate(time);
200 }
201}
202
203void Reactor::eval(double time, double* LHS, double* RHS)
204{
205 double& dmdt = RHS[0];
206 double* mdYdt = RHS + 3; // mass * dY/dt
207
208 evalWalls(time);
209 m_thermo->restoreState(m_state);
210 const vector_fp& mw = m_thermo->molecularWeights();
211 const doublereal* Y = m_thermo->massFractions();
212
213 evalSurfaces(LHS + m_nsp + 3, RHS + m_nsp + 3, m_sdot.data());
214 // mass added to gas phase from surface reactions
215 double mdot_surf = dot(m_sdot.begin(), m_sdot.end(), mw.begin());
216 dmdt = mdot_surf;
217
218 // volume equation
219 RHS[1] = m_vdot;
220
221 if (m_chem) {
222 m_kin->getNetProductionRates(&m_wdot[0]); // "omega dot"
223 }
224
225 for (size_t k = 0; k < m_nsp; k++) {
226 // production in gas phase and from surfaces
227 mdYdt[k] = (m_wdot[k] * m_vol + m_sdot[k]) * mw[k];
228 // dilution by net surface mass flux
229 mdYdt[k] -= Y[k] * mdot_surf;
230 LHS[k+3] = m_mass;
231 }
232
233 // Energy equation.
234 // \f[
235 // \dot U = -P\dot V + A \dot q + \dot m_{in} h_{in} - \dot m_{out} h.
236 // \f]
237 if (m_energy) {
238 RHS[2] = - m_thermo->pressure() * m_vdot + m_Qdot;
239 } else {
240 RHS[2] = 0.0;
241 }
242
243 // add terms for outlets
244 for (auto outlet : m_outlet) {
245 double mdot = outlet->massFlowRate();
246 dmdt -= mdot; // mass flow out of system
247 if (m_energy) {
248 RHS[2] -= mdot * m_enthalpy;
249 }
250 }
251
252 // add terms for inlets
253 for (auto inlet : m_inlet) {
254 double mdot = inlet->massFlowRate();
255 dmdt += mdot; // mass flow into system
256 for (size_t n = 0; n < m_nsp; n++) {
257 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
258 // flow of species into system and dilution by other species
259 mdYdt[n] += (mdot_spec - mdot * Y[n]);
260 }
261 if (m_energy) {
262 RHS[2] += mdot * inlet->enthalpy_mass();
263 }
264 }
265}
266
267void Reactor::evalWalls(double t)
268{
269 m_vdot = 0.0;
270 m_Qdot = 0.0;
271 for (size_t i = 0; i < m_wall.size(); i++) {
272 int f = 2 * m_lr[i] - 1;
273 m_vdot -= f * m_wall[i]->vdot(t);
274 m_Qdot += f * m_wall[i]->Q(t);
275 }
276 m_Q = -m_Qdot;
277}
278
279void Reactor::evalSurfaces(double* LHS, double* RHS, double* sdot)
280{
281 fill(sdot, sdot + m_nsp, 0.0);
282 size_t loc = 0; // offset into ydot
283
284 for (auto S : m_surfaces) {
285 Kinetics* kin = S->kinetics();
286 SurfPhase* surf = S->thermo();
287
288 double rs0 = 1.0/surf->siteDensity();
289 size_t nk = surf->nSpecies();
290 double sum = 0.0;
291 S->syncState();
292 kin->getNetProductionRates(&m_work[0]);
293 size_t ns = kin->surfacePhaseIndex();
294 size_t surfloc = kin->kineticsSpeciesIndex(0,ns);
295 for (size_t k = 1; k < nk; k++) {
296 LHS[loc] = 1.0;
297 RHS[loc + k] = m_work[surfloc + k] * rs0 * surf->size(k);
298 sum -= RHS[loc + k];
299 }
300 LHS[loc] = 1.0;
301 RHS[loc] = sum;
302 loc += nk;
303
304 size_t bulkloc = kin->kineticsSpeciesIndex(m_thermo->speciesName(0));
305 double wallarea = S->area();
306 for (size_t k = 0; k < m_nsp; k++) {
307 sdot[k] += m_work[bulkloc + k] * wallarea;
308 }
309 }
310}
311
312void Reactor::addSensitivityReaction(size_t rxn)
313{
314 if (!m_chem || rxn >= m_kin->nReactions()) {
315 throw CanteraError("Reactor::addSensitivityReaction",
316 "Reaction number out of range ({})", rxn);
317 }
318
319 size_t p = network().registerSensitivityParameter(
320 name()+": "+m_kin->reactionString(rxn), 1.0, 1.0);
321 m_sensParams.emplace_back(
322 SensitivityParameter{rxn, p, 1.0, SensParameterType::reaction});
323}
324
325void Reactor::addSensitivitySpeciesEnthalpy(size_t k)
326{
327 if (k >= m_thermo->nSpecies()) {
328 throw CanteraError("Reactor::addSensitivitySpeciesEnthalpy",
329 "Species index out of range ({})", k);
330 }
331
332 size_t p = network().registerSensitivityParameter(
333 name() + ": " + m_thermo->speciesName(k) + " enthalpy",
334 0.0, GasConstant * 298.15);
335 m_sensParams.emplace_back(
336 SensitivityParameter{k, p, m_thermo->Hf298SS(k),
337 SensParameterType::enthalpy});
338}
339
340size_t Reactor::speciesIndex(const string& nm) const
341{
342 // check for a gas species name
343 size_t k = m_thermo->speciesIndex(nm);
344 if (k != npos) {
345 return k;
346 }
347
348 // check for a wall species
349 size_t offset = m_nsp;
350 for (auto& S : m_surfaces) {
351 ThermoPhase* th = S->thermo();
352 k = th->speciesIndex(nm);
353 if (k != npos) {
354 return k + offset;
355 } else {
356 offset += th->nSpecies();
357 }
358 }
359 return npos;
360}
361
362size_t Reactor::componentIndex(const string& nm) const
363{
364 size_t k = speciesIndex(nm);
365 if (k != npos) {
366 return k + 3;
367 } else if (nm == "mass") {
368 return 0;
369 } else if (nm == "volume") {
370 return 1;
371 } else if (nm == "int_energy") {
372 return 2;
373 } else {
374 return npos;
375 }
376}
377
378std::string Reactor::componentName(size_t k) {
379 if (k == 0) {
380 return "mass";
381 } else if (k == 1) {
382 return "volume";
383 } else if (k == 2) {
384 return "int_energy";
385 } else if (k >= 3 && k < neq()) {
386 k -= 3;
387 if (k < m_thermo->nSpecies()) {
388 return m_thermo->speciesName(k);
389 } else {
390 k -= m_thermo->nSpecies();
391 }
392 for (auto& S : m_surfaces) {
393 ThermoPhase* th = S->thermo();
394 if (k < th->nSpecies()) {
395 return th->speciesName(k);
396 } else {
397 k -= th->nSpecies();
398 }
399 }
400 }
401 throw CanteraError("Reactor::componentName", "Index is out of bounds.");
402}
403
404void Reactor::applySensitivity(double* params)
405{
406 if (!params) {
407 return;
408 }
409 for (auto& p : m_sensParams) {
410 if (p.type == SensParameterType::reaction) {
411 p.value = m_kin->multiplier(p.local);
412 m_kin->setMultiplier(p.local, p.value*params[p.global]);
413 } else if (p.type == SensParameterType::enthalpy) {
414 m_thermo->modifyOneHf298SS(p.local, p.value + params[p.global]);
415 }
416 }
417 for (auto& S : m_surfaces) {
418 S->setSensitivityParameters(params);
419 }
420 m_thermo->invalidateCache();
421 if (m_kin) {
422 m_kin->invalidateCache();
423 }
424}
425
426void Reactor::resetSensitivity(double* params)
427{
428 if (!params) {
429 return;
430 }
431 for (auto& p : m_sensParams) {
432 if (p.type == SensParameterType::reaction) {
433 m_kin->setMultiplier(p.local, p.value);
434 } else if (p.type == SensParameterType::enthalpy) {
435 m_thermo->resetHf298(p.local);
436 }
437 }
438 for (auto& S : m_surfaces) {
439 S->resetSensitivityParameters();
440 }
441 m_thermo->invalidateCache();
442 if (m_kin) {
443 m_kin->invalidateCache();
444 }
445}
446
447void Reactor::setAdvanceLimits(const double *limits)
448{
449 if (m_thermo == 0) {
450 throw CanteraError("Reactor::setAdvanceLimits",
451 "Error: reactor is empty.");
452 }
453 m_advancelimits.assign(limits, limits + m_nv);
454
455 // resize to zero length if no limits are set
456 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
457 [](double val){return val>0;})) {
458 m_advancelimits.resize(0);
459 }
460}
461
462bool Reactor::getAdvanceLimits(double *limits)
463{
464 bool has_limit = hasAdvanceLimits();
465 if (has_limit) {
466 std::copy(m_advancelimits.begin(), m_advancelimits.end(), limits);
467 } else {
468 std::fill(limits, limits + m_nv, -1.0);
469 }
470 return has_limit;
471}
472
473void Reactor::setAdvanceLimit(const string& nm, const double limit)
474{
475 size_t k = componentIndex(nm);
476 if (k == npos) {
477 throw CanteraError("Reactor::setAdvanceLimit", "No component named '{}'", nm);
478 }
479
480 if (m_thermo == 0) {
481 throw CanteraError("Reactor::setAdvanceLimit",
482 "Error: reactor is empty.");
483 }
484 if (m_nv == 0) {
485 if (m_net == 0) {
486 throw CanteraError("Reactor::setAdvanceLimit",
487 "Cannot set limit on a reactor that is not "
488 "assigned to a ReactorNet object.");
489 } else {
490 m_net->initialize();
491 }
492 } else if (k > m_nv) {
493 throw CanteraError("Reactor::setAdvanceLimit",
494 "Index out of bounds.");
495 }
496 m_advancelimits.resize(m_nv, -1.0);
497 m_advancelimits[k] = limit;
498
499 // resize to zero length if no limits are set
500 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
501 [](double val){return val>0;})) {
502 m_advancelimits.resize(0);
503 }
504}
505
506}
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.
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
Public interface for kinetics managers.
Definition: Kinetics.h:114
virtual void getNetProductionRates(doublereal *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition: Kinetics.cpp:484
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition: Kinetics.h:265
size_t surfacePhaseIndex() const
This returns the integer index of the phase which has ThermoPhase type cSurf.
Definition: Kinetics.h:206
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:273
std::string speciesName(size_t k) const
Name of the species with index k.
Definition: Phase.cpp:200
size_t speciesIndex(const std::string &name) const
Returns the index of a species named 'name' within the Phase object.
Definition: Phase.cpp:187
A simple thermodynamic model for a surface phase, assuming an ideal solution model.
Definition: SurfPhase.h:125
virtual double size(size_t k) const
Returns the number of sites occupied by one molecule of species k.
Definition: SurfPhase.h:313
double siteDensity() const
Returns the site density.
Definition: SurfPhase.h:308
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:102
Base class for 'walls' (walls, pistons, etc.) connecting reactors.
Definition: Wall.h:25
virtual void initialize()
Called just before the start of integration.
Definition: Wall.h:69
Namespace for the Cantera kernel.
Definition: AnyMap.h:29
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:192
doublereal dot(InputIter x_begin, InputIter x_end, InputIter2 y_begin)
Function that calculates a templated inner product.
Definition: utilities.h:77
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:184
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition: ct_defs.h:113
Various templated functions that carry out common vector operations (see Templated Utility Functions)...