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