17#include <boost/math/tools/roots.hpp>
20namespace bmt = boost::math::tools;
25Reactor::Reactor(shared_ptr<Solution> sol,
const string& name)
26 : Reactor(sol, true, name)
30Reactor::Reactor(shared_ptr<Solution> sol,
bool clone,
const string& name)
31 : ReactorBase(sol, clone, name)
33 m_kin = m_solution->kinetics().get();
34 setChemistryEnabled(m_kin->nReactions() > 0);
36 m_sdot.resize(m_nsp, 0.0);
37 m_wdot.resize(m_nsp, 0.0);
41void Reactor::setDerivativeSettings(
AnyMap& settings)
43 m_kin->setDerivativeSettings(settings);
46void Reactor::getState(
double* y)
49 m_mass = m_thermo->density() * m_vol;
56 y[2] = m_thermo->intEnergy_mass() * m_mass;
59 m_thermo->getMassFractions(y+3);
62void Reactor::initialize(
double t0)
64 if (!m_thermo || (m_chem && !m_kin)) {
65 throw CanteraError(
"Reactor::initialize",
"Reactor contents not set"
66 " for reactor '" + m_name +
"'.");
68 updateConnected(
true);
70 for (
size_t n = 0; n < m_wall.size(); n++) {
76void Reactor::updateState(
double* y)
83 m_thermo->setMassFractions_NoNorm(y+3);
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;
93 double T = m_thermo->temperature();
94 boost::uintmax_t maxiter = 100;
95 pair<double, double> TT;
97 TT = bmt::bracket_and_solve_root(
98 u_err, T, 1.2,
true, bmt::eps_tolerance<double>(48), maxiter);
99 }
catch (std::exception&) {
103 TT = bmt::bisect(u_err, m_thermo->minTemp(), m_thermo->maxTemp(),
104 bmt::eps_tolerance<double>(48), maxiter);
105 }
catch (std::exception& err2) {
107 m_thermo->setState_TD(T, m_mass / m_vol);
109 "{}\nat U = {}, rho = {}", err2.what(), U, m_mass / m_vol);
112 if (fabs(TT.first - TT.second) > 1e-7*TT.first) {
113 throw CanteraError(
"Reactor::updateState",
"root finding failed");
115 m_thermo->setState_TD(TT.second, m_mass / m_vol);
117 m_thermo->setDensity(m_mass/m_vol);
120 updateConnected(
true);
123void Reactor::eval(
double time,
double* LHS,
double* RHS)
125 double& dmdt = RHS[0];
126 double* mdYdt = RHS + 3;
129 updateSurfaceProductionRates();
130 const vector<double>& mw = m_thermo->molecularWeights();
131 const double* Y = m_thermo->massFractions();
134 double mdot_surf =
dot(m_sdot.begin(), m_sdot.end(), mw.begin());
141 m_kin->getNetProductionRates(&m_wdot[0]);
144 for (
size_t k = 0; k < m_nsp; k++) {
146 mdYdt[k] = (m_wdot[k] * m_vol + m_sdot[k]) * mw[k];
148 mdYdt[k] -= Y[k] * mdot_surf;
157 RHS[2] = - m_thermo->pressure() * m_vdot + m_Qdot;
158 RHS[2] += m_thermo->intrinsicHeating() * m_vol;
164 for (
auto outlet : m_outlet) {
165 double mdot = outlet->massFlowRate();
168 RHS[2] -= mdot * m_enthalpy;
173 for (
auto inlet : m_inlet) {
174 double mdot = inlet->massFlowRate();
176 for (
size_t n = 0; n < m_nsp; n++) {
177 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
179 mdYdt[n] += (mdot_spec - mdot * Y[n]);
182 RHS[2] += mdot * inlet->enthalpy_mass();
187void Reactor::evalSteady(
double time,
double* LHS,
double* RHS)
189 eval(time, LHS, RHS);
190 RHS[1] = m_vol - m_initialVolume;
193void Reactor::evalWalls(
double t)
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();
205vector<size_t> Reactor::initializeSteady()
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());
213 m_initialVolume = m_vol;
217Eigen::SparseMatrix<double> Reactor::finiteDifferenceJacobian()
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;
224 Eigen::ArrayXd yPerturbed = yCurrent;
225 Eigen::ArrayXd lhsPerturbed(m_nv), lhsCurrent(m_nv);
226 Eigen::ArrayXd rhsPerturbed(m_nv), rhsCurrent(m_nv);
229 updateState(yCurrent.data());
230 eval(time, lhsCurrent.data(), rhsCurrent.data());
232 double rel_perturb = std::sqrt(std::numeric_limits<double>::epsilon());
233 double atol = (m_net !=
nullptr) ? m_net->atol() : 1e-15;
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;
240 updateState(yPerturbed.data());
243 eval(time, lhsPerturbed.data(), rhsPerturbed.data());
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);
255 updateState(yCurrent.data());
257 Eigen::SparseMatrix<double> jac(m_nv, m_nv);
258 jac.setFromTriplets(trips.begin(), trips.end());
262void Reactor::addSensitivityReaction(
size_t rxn)
264 if (!m_chem || rxn >= m_kin->nReactions()) {
266 "Reaction number out of range ({})", rxn);
269 size_t p = network().registerSensitivityParameter(
270 name()+
": "+m_kin->reaction(rxn)->equation(), 1.0, 1.0);
271 m_sensParams.emplace_back(
275void Reactor::addSensitivitySpeciesEnthalpy(
size_t k)
277 if (k >= m_thermo->nSpecies()) {
278 throw CanteraError(
"Reactor::addSensitivitySpeciesEnthalpy",
279 "Species index out of range ({})", k);
282 size_t p = network().registerSensitivityParameter(
283 name() +
": " + m_thermo->speciesName(k) +
" enthalpy",
285 m_sensParams.emplace_back(
287 SensParameterType::enthalpy});
290void Reactor::updateSurfaceProductionRates()
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();
302size_t Reactor::componentIndex(
const string& nm)
const
307 if (nm ==
"volume") {
310 if (nm ==
"int_energy") {
314 return m_thermo->speciesIndex(nm) + 3;
317 "Component '{}' not found", nm);
321string Reactor::componentName(
size_t k) {
328 }
else if (k >= 3 && k < neq()) {
329 return m_thermo->speciesName(k - 3);
331 throw IndexError(
"Reactor::componentName",
"component", k, m_nv);
334double Reactor::upperBound(
size_t k)
const {
341 }
else if (k >= 3 && k < m_nv) {
344 throw CanteraError(
"Reactor::upperBound",
"Index {} is out of bounds.", k);
348double Reactor::lowerBound(
size_t k)
const {
355 }
else if (k >= 3 && k < m_nv) {
358 throw CanteraError(
"Reactor::lowerBound",
"Index {} is out of bounds.", k);
362void Reactor::resetBadValues(
double* y) {
363 for (
size_t k = 3; k < m_nv; k++) {
364 y[k] = std::max(y[k], 0.0);
368void Reactor::applySensitivity(
double* params)
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]);
381 m_thermo->invalidateCache();
383 m_kin->invalidateCache();
387void Reactor::resetSensitivity(
double* params)
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);
399 m_thermo->invalidateCache();
401 m_kin->invalidateCache();
405void Reactor::setAdvanceLimits(
const double *limits)
407 m_advancelimits.assign(limits, limits + m_nv);
410 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
411 [](
double val){return val>0;})) {
412 m_advancelimits.resize(0);
416bool Reactor::getAdvanceLimits(
double *limits)
const
418 bool has_limit = hasAdvanceLimits();
420 std::copy(m_advancelimits.begin(), m_advancelimits.end(), limits);
422 std::fill(limits, limits + m_nv, -1.0);
427void Reactor::setAdvanceLimit(
const string& nm,
const double limit)
429 size_t k = componentIndex(nm);
430 m_advancelimits.resize(m_nv, -1.0);
431 m_advancelimits[k] = limit;
434 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
435 [](
double val){return val>0;})) {
436 m_advancelimits.resize(0);
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.
Base class for exceptions thrown by Cantera classes.
An array index is out of range.
Base class for 'walls' (walls, pistons, etc.) connecting reactors.
virtual void initialize()
Called just before the start of integration.
double dot(InputIter x_begin, InputIter x_end, InputIter2 y_begin)
Function that calculates a templated inner product.
const double GasConstant
Universal Gas Constant [J/kmol/K].
Namespace for the Cantera kernel.
const double Tiny
Small number to compare differences of mole fractions against.
offset
Offsets of solution components in the 1D solution array.
const double BigNumber
largest number to compare to inf.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...