Cantera  4.0.0a1
Loading...
Searching...
No Matches
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
16
17#include <boost/math/tools/roots.hpp>
18
19using namespace std;
20namespace bmt = boost::math::tools;
21
22namespace Cantera
23{
24
25Reactor::Reactor(shared_ptr<Solution> sol, const string& name)
26 : Reactor(sol, true, name)
27{
28}
29
30Reactor::Reactor(shared_ptr<Solution> sol, bool clone, const string& name)
31 : ReactorBase(sol, clone, name)
32{
33 m_kin = m_solution->kinetics().get();
34 setChemistryEnabled(m_kin->nReactions() > 0);
35 m_vol = 1.0; // By default, the volume is set to 1.0 m^3.
36 m_sdot.resize(m_nsp, 0.0);
37 m_wdot.resize(m_nsp, 0.0);
38 m_nv = 3 + m_nsp; // mass, volume, internal energy, and species mass fractions
39}
40
41void Reactor::setDerivativeSettings(AnyMap& settings)
42{
43 m_kin->setDerivativeSettings(settings);
44}
45
46void Reactor::getState(double* y)
47{
48 // set the first component to the total mass
49 m_mass = m_thermo->density() * m_vol;
50 y[0] = m_mass;
51
52 // set the second component to the total volume
53 y[1] = m_vol;
54
55 // set the third component to the total internal energy
56 y[2] = m_thermo->intEnergy_mass() * m_mass;
57
58 // set components y+3 ... y+K+2 to the mass fractions of each species
59 m_thermo->getMassFractions(y+3);
60}
61
62void Reactor::initialize(double t0)
63{
64 if (!m_thermo || (m_chem && !m_kin)) {
65 throw CanteraError("Reactor::initialize", "Reactor contents not set"
66 " for reactor '" + m_name + "'.");
67 }
68 updateConnected(true);
69
70 for (size_t n = 0; n < m_wall.size(); n++) {
71 WallBase* W = m_wall[n];
72 W->initialize();
73 }
74}
75
76void Reactor::updateState(double* y)
77{
78 // The components of y are [0] the total mass, [1] the total volume,
79 // [2] the total internal energy, [3...K+3] are the mass fractions of each
80 // species, and [K+3...] are the coverages of surface species on each wall.
81 m_mass = y[0];
82 m_vol = y[1];
83 m_thermo->setMassFractions_NoNorm(y+3);
84
85 if (m_energy) {
86 double U = y[2];
87 // Residual function: error in internal energy as a function of T
88 auto u_err = [this, U](double T) {
89 m_thermo->setState_TD(T, m_mass / m_vol);
90 return m_thermo->intEnergy_mass() * m_mass - U;
91 };
92
93 double T = m_thermo->temperature();
94 boost::uintmax_t maxiter = 100;
95 pair<double, double> TT;
96 try {
97 TT = bmt::bracket_and_solve_root(
98 u_err, T, 1.2, true, bmt::eps_tolerance<double>(48), maxiter);
99 } catch (std::exception&) {
100 // Try full-range bisection if bracketing fails (for example, near
101 // temperature limits for the phase's equation of state)
102 try {
103 TT = bmt::bisect(u_err, m_thermo->minTemp(), m_thermo->maxTemp(),
104 bmt::eps_tolerance<double>(48), maxiter);
105 } catch (std::exception& err2) {
106 // Set m_thermo back to a reasonable state if root finding fails
107 m_thermo->setState_TD(T, m_mass / m_vol);
108 throw CanteraError("Reactor::updateState",
109 "{}\nat U = {}, rho = {}", err2.what(), U, m_mass / m_vol);
110 }
111 }
112 if (fabs(TT.first - TT.second) > 1e-7*TT.first) {
113 throw CanteraError("Reactor::updateState", "root finding failed");
114 }
115 m_thermo->setState_TD(TT.second, m_mass / m_vol);
116 } else {
117 m_thermo->setDensity(m_mass/m_vol);
118 }
119
120 updateConnected(true);
121}
122
123void Reactor::eval(double time, double* LHS, double* RHS)
124{
125 double& dmdt = RHS[0];
126 double* mdYdt = RHS + 3; // mass * dY/dt
127
128 evalWalls(time);
129 updateSurfaceProductionRates();
130 const vector<double>& mw = m_thermo->molecularWeights();
131 const double* Y = m_thermo->massFractions();
132
133 // mass added to gas phase from surface reactions
134 double mdot_surf = dot(m_sdot.begin(), m_sdot.end(), mw.begin());
135 dmdt = mdot_surf;
136
137 // volume equation
138 RHS[1] = m_vdot;
139
140 if (m_chem) {
141 m_kin->getNetProductionRates(&m_wdot[0]); // "omega dot"
142 }
143
144 for (size_t k = 0; k < m_nsp; k++) {
145 // production in gas phase and from surfaces
146 mdYdt[k] = (m_wdot[k] * m_vol + m_sdot[k]) * mw[k];
147 // dilution by net surface mass flux
148 mdYdt[k] -= Y[k] * mdot_surf;
149 LHS[k+3] = m_mass;
150 }
151
152 // Energy equation.
153 // @f[
154 // \dot U = -P\dot V + A \dot q + \dot m_{in} h_{in} - \dot m_{out} h.
155 // @f]
156 if (m_energy) {
157 RHS[2] = - m_thermo->pressure() * m_vdot + m_Qdot;
158 RHS[2] += m_thermo->intrinsicHeating() * m_vol;
159 } else {
160 RHS[2] = 0.0;
161 }
162
163 // add terms for outlets
164 for (auto outlet : m_outlet) {
165 double mdot = outlet->massFlowRate();
166 dmdt -= mdot; // mass flow out of system
167 if (m_energy) {
168 RHS[2] -= mdot * m_enthalpy;
169 }
170 }
171
172 // add terms for inlets
173 for (auto inlet : m_inlet) {
174 double mdot = inlet->massFlowRate();
175 dmdt += mdot; // mass flow into system
176 for (size_t n = 0; n < m_nsp; n++) {
177 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
178 // flow of species into system and dilution by other species
179 mdYdt[n] += (mdot_spec - mdot * Y[n]);
180 }
181 if (m_energy) {
182 RHS[2] += mdot * inlet->enthalpy_mass();
183 }
184 }
185}
186
187void Reactor::evalSteady(double time, double* LHS, double* RHS)
188{
189 eval(time, LHS, RHS);
190 RHS[1] = m_vol - m_initialVolume;
191}
192
193void Reactor::evalWalls(double t)
194{
195 // time is currently unused
196 m_vdot = 0.0;
197 m_Qdot = 0.0;
198 for (size_t i = 0; i < m_wall.size(); i++) {
199 int f = 2 * m_lr[i] - 1;
200 m_vdot -= f * m_wall[i]->expansionRate();
201 m_Qdot += f * m_wall[i]->heatRate();
202 }
203}
204
205vector<size_t> Reactor::initializeSteady()
206{
207 if (!energyEnabled()) {
208 throw CanteraError("Reactor::initializeSteady", "Steady state solver cannot"
209 " be used with {0} when energy equation is disabled."
210 "\nConsider using IdealGas{0} instead.\n"
211 "See https://github.com/Cantera/enhancements/issues/234", type());
212 }
213 m_initialVolume = m_vol;
214 return {1}; // volume
215}
216
217Eigen::SparseMatrix<double> Reactor::finiteDifferenceJacobian()
218{
219 vector<Eigen::Triplet<double>> trips;
220 Eigen::ArrayXd yCurrent(m_nv);
221 getState(yCurrent.data());
222 double time = (m_net != nullptr) ? m_net->time() : 0.0;
223
224 Eigen::ArrayXd yPerturbed = yCurrent;
225 Eigen::ArrayXd lhsPerturbed(m_nv), lhsCurrent(m_nv);
226 Eigen::ArrayXd rhsPerturbed(m_nv), rhsCurrent(m_nv);
227 lhsCurrent = 1.0;
228 rhsCurrent = 0.0;
229 updateState(yCurrent.data());
230 eval(time, lhsCurrent.data(), rhsCurrent.data());
231
232 double rel_perturb = std::sqrt(std::numeric_limits<double>::epsilon());
233 double atol = (m_net != nullptr) ? m_net->atol() : 1e-15;
234
235 for (size_t j = 0; j < m_nv; j++) {
236 yPerturbed = yCurrent;
237 double delta_y = std::max(std::abs(yCurrent[j]), 1000 * atol) * rel_perturb;
238 yPerturbed[j] += delta_y;
239
240 updateState(yPerturbed.data());
241 lhsPerturbed = 1.0;
242 rhsPerturbed = 0.0;
243 eval(time, lhsPerturbed.data(), rhsPerturbed.data());
244
245 // d ydot_i/dy_j
246 for (size_t i = 0; i < m_nv; i++) {
247 double ydotPerturbed = rhsPerturbed[i] / lhsPerturbed[i];
248 double ydotCurrent = rhsCurrent[i] / lhsCurrent[i];
249 if (ydotCurrent != ydotPerturbed) {
250 trips.emplace_back(static_cast<int>(i), static_cast<int>(j),
251 (ydotPerturbed - ydotCurrent) / delta_y);
252 }
253 }
254 }
255 updateState(yCurrent.data());
256
257 Eigen::SparseMatrix<double> jac(m_nv, m_nv);
258 jac.setFromTriplets(trips.begin(), trips.end());
259 return jac;
260}
261
262void Reactor::addSensitivityReaction(size_t rxn)
263{
264 if (!m_chem || rxn >= m_kin->nReactions()) {
265 throw CanteraError("Reactor::addSensitivityReaction",
266 "Reaction number out of range ({})", rxn);
267 }
268
269 size_t p = network().registerSensitivityParameter(
270 name()+": "+m_kin->reaction(rxn)->equation(), 1.0, 1.0);
271 m_sensParams.emplace_back(
272 SensitivityParameter{rxn, p, 1.0, SensParameterType::reaction});
273}
274
275void Reactor::addSensitivitySpeciesEnthalpy(size_t k)
276{
277 if (k >= m_thermo->nSpecies()) {
278 throw CanteraError("Reactor::addSensitivitySpeciesEnthalpy",
279 "Species index out of range ({})", k);
280 }
281
282 size_t p = network().registerSensitivityParameter(
283 name() + ": " + m_thermo->speciesName(k) + " enthalpy",
284 0.0, GasConstant * 298.15);
285 m_sensParams.emplace_back(
286 SensitivityParameter{k, p, m_thermo->Hf298SS(k),
287 SensParameterType::enthalpy});
288}
289
290void Reactor::updateSurfaceProductionRates()
291{
292 m_sdot.assign(m_nsp, 0.0);
293 for (auto& S : m_surfaces) {
294 const auto& sdot = S->surfaceProductionRates();
295 size_t offset = S->kinetics()->kineticsSpeciesIndex(m_thermo->speciesName(0));
296 for (size_t k = 0; k < m_nsp; k++) {
297 m_sdot[k] += sdot[offset + k] * S->area();
298 }
299 }
300}
301
302size_t Reactor::componentIndex(const string& nm) const
303{
304 if (nm == "mass") {
305 return 0;
306 }
307 if (nm == "volume") {
308 return 1;
309 }
310 if (nm == "int_energy") {
311 return 2;
312 }
313 try {
314 return m_thermo->speciesIndex(nm) + 3;
315 } catch (const CanteraError&) {
316 throw CanteraError("Reactor::componentIndex",
317 "Component '{}' not found", nm);
318 }
319}
320
321string Reactor::componentName(size_t k) {
322 if (k == 0) {
323 return "mass";
324 } else if (k == 1) {
325 return "volume";
326 } else if (k == 2) {
327 return "int_energy";
328 } else if (k >= 3 && k < neq()) {
329 return m_thermo->speciesName(k - 3);
330 }
331 throw IndexError("Reactor::componentName", "component", k, m_nv);
332}
333
334double Reactor::upperBound(size_t k) const {
335 if (k == 0) {
336 return BigNumber; // mass
337 } else if (k == 1) {
338 return BigNumber; // volume
339 } else if (k == 2) {
340 return BigNumber; // internal energy
341 } else if (k >= 3 && k < m_nv) {
342 return 1.0; // species mass fraction or surface coverage
343 } else {
344 throw CanteraError("Reactor::upperBound", "Index {} is out of bounds.", k);
345 }
346}
347
348double Reactor::lowerBound(size_t k) const {
349 if (k == 0) {
350 return 0; // mass
351 } else if (k == 1) {
352 return 0; // volume
353 } else if (k == 2) {
354 return -BigNumber; // internal energy
355 } else if (k >= 3 && k < m_nv) {
356 return -Tiny; // species mass fraction or surface coverage
357 } else {
358 throw CanteraError("Reactor::lowerBound", "Index {} is out of bounds.", k);
359 }
360}
361
362void Reactor::resetBadValues(double* y) {
363 for (size_t k = 3; k < m_nv; k++) {
364 y[k] = std::max(y[k], 0.0);
365 }
366}
367
368void Reactor::applySensitivity(double* params)
369{
370 if (!params) {
371 return;
372 }
373 for (auto& p : m_sensParams) {
374 if (p.type == SensParameterType::reaction) {
375 p.value = m_kin->multiplier(p.local);
376 m_kin->setMultiplier(p.local, p.value*params[p.global]);
377 } else if (p.type == SensParameterType::enthalpy) {
378 m_thermo->modifyOneHf298SS(p.local, p.value + params[p.global]);
379 }
380 }
381 m_thermo->invalidateCache();
382 if (m_kin) {
383 m_kin->invalidateCache();
384 }
385}
386
387void Reactor::resetSensitivity(double* params)
388{
389 if (!params) {
390 return;
391 }
392 for (auto& p : m_sensParams) {
393 if (p.type == SensParameterType::reaction) {
394 m_kin->setMultiplier(p.local, p.value);
395 } else if (p.type == SensParameterType::enthalpy) {
396 m_thermo->resetHf298(p.local);
397 }
398 }
399 m_thermo->invalidateCache();
400 if (m_kin) {
401 m_kin->invalidateCache();
402 }
403}
404
405void Reactor::setAdvanceLimits(const double *limits)
406{
407 m_advancelimits.assign(limits, limits + m_nv);
408
409 // resize to zero length if no limits are set
410 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
411 [](double val){return val>0;})) {
412 m_advancelimits.resize(0);
413 }
414}
415
416bool Reactor::getAdvanceLimits(double *limits) const
417{
418 bool has_limit = hasAdvanceLimits();
419 if (has_limit) {
420 std::copy(m_advancelimits.begin(), m_advancelimits.end(), limits);
421 } else {
422 std::fill(limits, limits + m_nv, -1.0);
423 }
424 return has_limit;
425}
426
427void Reactor::setAdvanceLimit(const string& nm, const double limit)
428{
429 size_t k = componentIndex(nm);
430 m_advancelimits.resize(m_nv, -1.0);
431 m_advancelimits[k] = limit;
432
433 // resize to zero length if no limits are set
434 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
435 [](double val){return val>0;})) {
436 m_advancelimits.resize(0);
437 }
438}
439
440}
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:431
Base class for exceptions thrown by Cantera classes.
An array index is out of range.
Base class for 'walls' (walls, pistons, etc.) connecting reactors.
Definition Wall.h:23
virtual void initialize()
Called just before the start of integration.
Definition Wall.h:63
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:122
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const double Tiny
Small number to compare differences of mole fractions against.
Definition ct_defs.h:175
offset
Offsets of solution components in the 1D solution array.
Definition Flow1D.h:25
const double BigNumber
largest number to compare to inf.
Definition ct_defs.h:162
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...