16#include "cantera/numerics/eigen_dense.h"
18#include <boost/math/tools/roots.hpp>
21namespace bmt = boost::math::tools;
26Reactor::Reactor(shared_ptr<Solution> sol,
const string& name)
27 : Reactor(sol, true, name)
31Reactor::Reactor(shared_ptr<Solution> sol,
bool clone,
const string& name)
32 : ReactorBase(sol, clone, name)
34 m_kin = m_solution->kinetics().get();
35 setChemistryEnabled(m_kin->nReactions() > 0);
37 m_sdot.resize(m_nsp, 0.0);
38 m_wdot.resize(m_nsp, 0.0);
42void Reactor::setDerivativeSettings(
AnyMap& settings)
44 m_kin->setDerivativeSettings(settings);
47void Reactor::getState(span<double> y)
50 m_mass = m_thermo->density() * m_vol;
57 y[2] = m_thermo->intEnergy_mass() * m_mass;
60 m_thermo->getMassFractions(y.subspan(3, m_nsp));
63void Reactor::initialize(
double t0)
65 if (!m_thermo || (m_chem && !m_kin)) {
66 throw CanteraError(
"Reactor::initialize",
"Reactor contents not set"
67 " for reactor '" + m_name +
"'.");
69 updateConnected(
true);
71 for (
size_t n = 0; n < m_wall.size(); n++) {
77void Reactor::updateState(span<const double> y)
84 m_thermo->setMassFractions_NoNorm(y.subspan(3, m_nsp));
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;
94 double T = m_thermo->temperature();
95 boost::uintmax_t maxiter = 100;
96 pair<double, double> TT;
98 TT = bmt::bracket_and_solve_root(
99 u_err, T, 1.2,
true, bmt::eps_tolerance<double>(48), maxiter);
100 }
catch (std::exception&) {
104 TT = bmt::bisect(u_err, m_thermo->minTemp(), m_thermo->maxTemp(),
105 bmt::eps_tolerance<double>(48), maxiter);
106 }
catch (std::exception& err2) {
108 m_thermo->setState_TD(T, m_mass / m_vol);
110 "{}\nat U = {}, rho = {}", err2.what(), U, m_mass / m_vol);
113 if (fabs(TT.first - TT.second) > 1e-7*TT.first) {
114 throw CanteraError(
"Reactor::updateState",
"root finding failed");
116 m_thermo->setState_TD(TT.second, m_mass / m_vol);
118 m_thermo->setDensity(m_mass/m_vol);
121 updateConnected(
true);
124void Reactor::eval(
double time, span<double> LHS, span<double> RHS)
126 double& dmdt = RHS[0];
127 auto mdYdt = RHS.subspan(3);
130 updateSurfaceProductionRates();
131 auto mw = m_thermo->molecularWeights();
132 auto Y = m_thermo->massFractions();
135 double mdot_surf =
dot(m_sdot.begin(), m_sdot.end(), mw.begin());
142 m_kin->getNetProductionRates(m_wdot);
145 for (
size_t k = 0; k < m_nsp; k++) {
147 mdYdt[k] = (m_wdot[k] * m_vol + m_sdot[k]) * mw[k];
149 mdYdt[k] -= Y[k] * mdot_surf;
158 RHS[2] = - m_thermo->pressure() * m_vdot + m_Qdot;
159 RHS[2] += m_thermo->intrinsicHeating() * m_vol;
165 for (
auto outlet : m_outlet) {
166 double mdot = outlet->massFlowRate();
169 RHS[2] -= mdot * m_enthalpy;
174 for (
auto inlet : m_inlet) {
175 double mdot = inlet->massFlowRate();
177 for (
size_t n = 0; n < m_nsp; n++) {
178 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
180 mdYdt[n] += (mdot_spec - mdot * Y[n]);
183 RHS[2] += mdot * inlet->enthalpy_mass();
188void Reactor::evalSteady(
double time, span<double> LHS, span<double> RHS)
190 eval(time, LHS, RHS);
191 RHS[1] = m_vol - m_initialVolume;
194void Reactor::evalWalls(
double t)
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();
206vector<size_t> Reactor::initializeSteady()
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());
214 m_initialVolume = m_vol;
218Eigen::SparseMatrix<double> Reactor::finiteDifferenceJacobian()
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;
225 Eigen::ArrayXd yPerturbed = yCurrent;
226 Eigen::ArrayXd lhsPerturbed(m_nv), lhsCurrent(m_nv);
227 Eigen::ArrayXd rhsPerturbed(m_nv), rhsCurrent(m_nv);
230 updateState(
asSpan(yCurrent));
233 double rel_perturb = std::sqrt(std::numeric_limits<double>::epsilon());
234 double atol = (m_net !=
nullptr) ? m_net->atol() : 1e-15;
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;
241 updateState(
asSpan(yPerturbed));
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);
256 updateState(
asSpan(yCurrent));
258 Eigen::SparseMatrix<double> jac(m_nv, m_nv);
259 jac.setFromTriplets(trips.begin(), trips.end());
263void Reactor::addSensitivityReaction(
size_t rxn)
265 if (!m_chem || rxn >= m_kin->nReactions()) {
267 "Reaction number out of range ({})", rxn);
270 size_t p = network().registerSensitivityParameter(
271 name()+
": "+m_kin->reaction(rxn)->equation(), 1.0, 1.0);
272 m_sensParams.emplace_back(
276void Reactor::addSensitivitySpeciesEnthalpy(
size_t k)
278 if (k >= m_thermo->nSpecies()) {
279 throw CanteraError(
"Reactor::addSensitivitySpeciesEnthalpy",
280 "Species index out of range ({})", k);
283 size_t p = network().registerSensitivityParameter(
284 name() +
": " + m_thermo->speciesName(k) +
" enthalpy",
286 m_sensParams.emplace_back(
288 SensParameterType::enthalpy});
291void Reactor::updateSurfaceProductionRates()
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();
303size_t Reactor::componentIndex(
const string& nm)
const
308 if (nm ==
"volume") {
311 if (nm ==
"int_energy") {
315 return m_thermo->speciesIndex(nm) + 3;
318 "Component '{}' not found", nm);
322string Reactor::componentName(
size_t k) {
329 }
else if (k >= 3 && k < neq()) {
330 return m_thermo->speciesName(k - 3);
332 throw IndexError(
"Reactor::componentName",
"component", k, m_nv);
335double Reactor::upperBound(
size_t k)
const {
342 }
else if (k >= 3 && k < m_nv) {
345 throw CanteraError(
"Reactor::upperBound",
"Index {} is out of bounds.", k);
349double Reactor::lowerBound(
size_t k)
const {
356 }
else if (k >= 3 && k < m_nv) {
359 throw CanteraError(
"Reactor::lowerBound",
"Index {} is out of bounds.", k);
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);
369void Reactor::applySensitivity(span<const double> params)
371 if (params.empty()) {
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]);
382 m_thermo->invalidateCache();
384 m_kin->invalidateCache();
388void Reactor::resetSensitivity(span<const double> params)
390 if (params.empty()) {
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);
400 m_thermo->invalidateCache();
402 m_kin->invalidateCache();
406void Reactor::setAdvanceLimits(span<const double> limits)
408 m_advancelimits.assign(limits.begin(), limits.end());
411 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
412 [](
double val){return val>0;})) {
413 m_advancelimits.resize(0);
417bool Reactor::getAdvanceLimits(span<double> limits)
const
419 bool has_limit = hasAdvanceLimits();
421 std::copy(m_advancelimits.begin(), m_advancelimits.end(), limits.begin());
423 std::fill(limits.begin(), limits.end(), -1.0);
428void Reactor::setAdvanceLimit(
const string& nm,
const double limit)
430 size_t k = componentIndex(nm);
431 m_advancelimits.resize(m_nv, -1.0);
432 m_advancelimits[k] = limit;
435 if (std::none_of(m_advancelimits.begin(), m_advancelimits.end(),
436 [](
double val){return val>0;})) {
437 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.
span< double > asSpan(Eigen::DenseBase< Derived > &v)
Convenience wrapper for accessing Eigen vector/array/map data as a span.
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...