17#include <boost/math/tools/roots.hpp>
20namespace bmt = boost::math::tools;
27Reactor::Reactor(shared_ptr<Solution> sol,
const string& name)
30 m_kin = m_solution->kinetics().get();
31 setChemistryEnabled(m_kin->nReactions() > 0);
35Reactor::Reactor(shared_ptr<Solution> sol,
bool clone,
const string& name)
36 : ReactorBase(sol, clone, name)
38 m_kin = m_solution->kinetics().get();
39 setChemistryEnabled(m_kin->nReactions() > 0);
45 m_kin->setDerivativeSettings(settings);
47 for (
auto S : m_surfaces) {
48 S->kinetics()->setDerivativeSettings(settings);
55 "After Cantera 3.2, a change of reactor contents after instantiation "
65 "Error: reactor is empty.");
67 m_thermo->restoreState(m_state);
77 y[2] = m_thermo->intEnergy_mass() *
m_mass;
80 m_thermo->getMassFractions(y+3);
90 for (
auto& S : m_surfaces) {
91 S->getCoverages(y + loc);
92 loc += S->thermo()->nSpecies();
98 if (!m_thermo || (m_chem && !
m_kin)) {
99 throw CanteraError(
"Reactor::initialize",
"Reactor contents not set"
100 " for reactor '" +
m_name +
"'.");
102 m_thermo->restoreState(m_state);
107 for (
size_t n = 0; n < m_wall.size(); n++) {
115 for (
auto& S : m_surfaces) {
116 m_nv_surf += S->thermo()->nSpecies();
117 size_t nt = S->kinetics()->nTotalSpecies();
118 maxnt = std::max(maxnt, nt);
121 m_work.resize(maxnt);
126 size_t ns = m_sensParams.size();
127 for (
auto& S : m_surfaces) {
128 ns += S->nSensParams();
140 m_thermo->setMassFractions_NoNorm(y+3);
145 auto u_err = [
this, U](
double T) {
147 return m_thermo->intEnergy_mass() *
m_mass - U;
150 double T = m_thermo->temperature();
151 boost::uintmax_t maxiter = 100;
152 pair<double, double> TT;
154 TT = bmt::bracket_and_solve_root(
155 u_err, T, 1.2,
true, bmt::eps_tolerance<double>(48), maxiter);
156 }
catch (std::exception&) {
160 TT = bmt::bisect(u_err, m_thermo->minTemp(), m_thermo->maxTemp(),
161 bmt::eps_tolerance<double>(48), maxiter);
162 }
catch (std::exception& err2) {
166 "{}\nat U = {}, rho = {}", err2.what(), U,
m_mass /
m_vol);
169 if (fabs(TT.first - TT.second) > 1e-7*TT.first) {
170 throw CanteraError(
"Reactor::updateState",
"root finding failed");
184 for (
auto& S : m_surfaces) {
185 S->setCoverages(y+loc);
186 loc += S->thermo()->nSpecies();
193 if (updatePressure) {
201 m_thermo->saveState(m_state);
205 if (
m_net !=
nullptr) {
208 for (
size_t i = 0; i < m_outlet.size(); i++) {
209 m_outlet[i]->setSimTime(time);
210 m_outlet[i]->updateMassFlowRate(time);
212 for (
size_t i = 0; i < m_inlet.size(); i++) {
213 m_inlet[i]->setSimTime(time);
214 m_inlet[i]->updateMassFlowRate(time);
216 for (
size_t i = 0; i < m_wall.size(); i++) {
217 m_wall[i]->setSimTime(time);
223 double& dmdt = RHS[0];
224 double* mdYdt = RHS + 3;
227 m_thermo->restoreState(m_state);
228 const vector<double>& mw = m_thermo->molecularWeights();
229 const double* Y = m_thermo->massFractions();
243 for (
size_t k = 0; k <
m_nsp; k++) {
247 mdYdt[k] -= Y[k] * mdot_surf;
262 for (
auto outlet : m_outlet) {
263 double mdot =
outlet->massFlowRate();
271 for (
auto inlet : m_inlet) {
272 double mdot =
inlet->massFlowRate();
274 for (
size_t n = 0; n <
m_nsp; n++) {
275 double mdot_spec =
inlet->outletSpeciesMassFlowRate(n);
277 mdYdt[n] += (mdot_spec - mdot * Y[n]);
280 RHS[2] += mdot *
inlet->enthalpy_mass();
290 for (
size_t i = 0; i < m_wall.size(); i++) {
291 int f = 2 *
m_lr[i] - 1;
292 m_vdot -= f * m_wall[i]->expansionRate();
293 m_Qdot += f * m_wall[i]->heatRate();
299 fill(sdot, sdot +
m_nsp, 0.0);
302 for (
auto S : m_surfaces) {
311 for (
size_t k = 1; k < nk; k++) {
312 RHS[loc + k] = m_work[k] * rs0 * surf->
size(k);
319 double wallarea = S->area();
320 for (
size_t k = 0; k <
m_nsp; k++) {
321 sdot[k] += m_work[bulkloc + k] * wallarea;
329 throw CanteraError(
"Reactor::steadyConstraints",
"Steady state solver cannot"
330 " be used with {0} when energy equation is disabled."
331 "\nConsider using IdealGas{0} instead.\n"
332 "See https://github.com/Cantera/enhancements/issues/234",
type());
335 throw CanteraError(
"Reactor::steadyConstraints",
"Steady state solver cannot"
336 " currently be used when reactor surfaces are present.\n"
337 "See https://github.com/Cantera/enhancements/issues/234.");
346 "Reactor must be initialized first.");
351 Eigen::ArrayXd yCurrent(m_nv);
353 double time = (
m_net !=
nullptr) ?
m_net->time() : 0.0;
355 Eigen::ArrayXd yPerturbed = yCurrent;
356 Eigen::ArrayXd lhsPerturbed(m_nv), lhsCurrent(m_nv);
357 Eigen::ArrayXd rhsPerturbed(m_nv), rhsCurrent(m_nv);
361 eval(time, lhsCurrent.data(), rhsCurrent.data());
363 double rel_perturb = std::sqrt(std::numeric_limits<double>::epsilon());
364 double atol = (
m_net !=
nullptr) ?
m_net->atol() : 1e-15;
366 for (
size_t j = 0; j < m_nv; j++) {
367 yPerturbed = yCurrent;
368 double delta_y = std::max(std::abs(yCurrent[j]), 1000 * atol) * rel_perturb;
369 yPerturbed[j] += delta_y;
374 eval(time, lhsPerturbed.data(), rhsPerturbed.data());
377 for (
size_t i = 0; i < m_nv; i++) {
378 double ydotPerturbed = rhsPerturbed[i] / lhsPerturbed[i];
379 double ydotCurrent = rhsCurrent[i] / lhsCurrent[i];
380 if (ydotCurrent != ydotPerturbed) {
382 static_cast<int>(i),
static_cast<int>(j),
383 (ydotPerturbed - ydotCurrent) / delta_y);
389 Eigen::SparseMatrix<double> jac(m_nv, m_nv);
397 fill(sdot, sdot + m_nsp, 0.0);
400 for (
auto S : m_surfaces) {
409 for (
size_t k = 1; k < nk; k++) {
410 RHS[loc + k] = m_work[k] * rs0 * surf->
size(k);
417 double wallarea = S->area();
418 for (
size_t k = 0; k < m_nsp; k++) {
419 sdot[k] += m_work[bulkloc + k] * wallarea;
426 if (!m_chem || rxn >=
m_kin->nReactions()) {
428 "Reaction number out of range ({})", rxn);
431 size_t p =
network().registerSensitivityParameter(
432 name()+
": "+
m_kin->reaction(rxn)->equation(), 1.0, 1.0);
433 m_sensParams.emplace_back(
439 if (k >= m_thermo->nSpecies()) {
440 throw CanteraError(
"Reactor::addSensitivitySpeciesEnthalpy",
441 "Species index out of range ({})", k);
444 size_t p =
network().registerSensitivityParameter(
445 name() +
": " + m_thermo->speciesName(k) +
" enthalpy",
447 m_sensParams.emplace_back(
449 SensParameterType::enthalpy});
455 size_t k = m_thermo->speciesIndex(nm,
false);
462 for (
auto& S : m_surfaces) {
472 "Species '{}' not found", nm);
480 if (nm ==
"volume") {
483 if (nm ==
"int_energy") {
490 "Component '{}' not found", nm);
501 }
else if (k >= 3 && k <
neq()) {
503 if (k < m_thermo->nSpecies()) {
504 return m_thermo->speciesName(k);
506 k -= m_thermo->nSpecies();
508 for (
auto& S : m_surfaces) {
510 if (k < th->nSpecies()) {
517 throw IndexError(
"Reactor::componentName",
"component", k, m_nv);
527 }
else if (k >= 3 && k < m_nv) {
530 throw CanteraError(
"Reactor::upperBound",
"Index {} is out of bounds.", k);
541 }
else if (k >= 3 && k < m_nv) {
544 throw CanteraError(
"Reactor::lowerBound",
"Index {} is out of bounds.", k);
549 for (
size_t k = 3; k < m_nv; k++) {
550 y[k] = std::max(y[k], 0.0);
559 for (
auto& p : m_sensParams) {
560 if (p.type == SensParameterType::reaction) {
561 p.value =
m_kin->multiplier(p.local);
562 m_kin->setMultiplier(p.local, p.value*params[p.global]);
563 }
else if (p.type == SensParameterType::enthalpy) {
564 m_thermo->modifyOneHf298SS(p.local, p.value + params[p.global]);
567 for (
auto& S : m_surfaces) {
568 S->setSensitivityParameters(params);
570 m_thermo->invalidateCache();
572 m_kin->invalidateCache();
581 for (
auto& p : m_sensParams) {
582 if (p.type == SensParameterType::reaction) {
583 m_kin->setMultiplier(p.local, p.value);
584 }
else if (p.type == SensParameterType::enthalpy) {
585 m_thermo->resetHf298(p.local);
588 for (
auto& S : m_surfaces) {
589 S->resetSensitivityParameters();
591 m_thermo->invalidateCache();
593 m_kin->invalidateCache();
601 "Error: reactor is empty.");
607 [](
double val){return val>0;})) {
618 std::fill(limits, limits + m_nv, -1.0);
629 "Error: reactor is empty.");
634 "Cannot set limit on a reactor that is not "
635 "assigned to a ReactorNet object.");
639 }
else if (k > m_nv) {
641 "Index out of bounds.");
648 [](
double val){return val>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.
Public interface for kinetics managers.
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
virtual void getNetProductionRates(double *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
An error indicating that an unimplemented function has been called.
size_t nSpecies() const
Returns the number of species in the phase.
string speciesName(size_t k) const
Name of the species with index k.
size_t speciesIndex(const string &name) const
Returns the index of a species named 'name' within the Phase object.
Base class for reactor objects.
FlowDevice & outlet(size_t n=0)
Return a reference to the n-th outlet FlowDevice connected to this reactor.
double m_pressure
Current pressure in the reactor [Pa].
ReactorNet * m_net
The ReactorNet that this reactor is part of.
virtual size_t nSurfs() const
Return the number of surfaces in a reactor.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
vector< int > m_lr
Vector of length nWalls(), indicating whether this reactor is on the left (0) or right (1) of each wa...
double m_vol
Current volume of the reactor [m^3].
double m_intEnergy
Current internal energy of the reactor [J/kg].
double m_mass
Current mass of the reactor [kg].
size_t m_nsp
Number of homogeneous species in the mixture.
string m_name
Reactor name.
ReactorNet & network()
The ReactorNet that this reactor belongs to.
double m_enthalpy
Current specific enthalpy of the reactor [J/kg].
string name() const
Return the name of this reactor.
virtual double lowerBound(size_t k) const
Get the lower bound on the k-th component of the local state vector.
virtual string componentName(size_t k)
Return the name of the solution component with index i.
void setChemistryEnabled(bool cflag=true) override
Enable or disable changes in reactor composition due to chemical reactions.
virtual void evalSurfaces(double *LHS, double *RHS, double *sdot)
Evaluate terms related to surface reactions.
virtual size_t componentIndex(const string &nm) const
Return the index in the solution vector for this reactor of the component named nm.
virtual void updateSurfaceState(double *y)
Update the state of SurfPhase objects attached to this reactor.
virtual void applySensitivity(double *params)
Set reaction rate multipliers based on the sensitivity variables in params.
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
vector< double > m_wdot
Species net molar production rates.
virtual void eval(double t, double *LHS, double *RHS)
Evaluate the reactor governing equations.
virtual void evalWalls(double t)
Evaluate terms related to Walls.
Eigen::SparseMatrix< double > finiteDifferenceJacobian()
Calculate the reactor-specific Jacobian using a finite difference method.
bool energyEnabled() const override
Returns true if solution of the energy equation is enabled.
size_t neq()
Number of equations (state variables) for this reactor.
virtual void addSensitivitySpeciesEnthalpy(size_t k)
Add a sensitivity parameter associated with the enthalpy formation of species k (in the homogeneous p...
size_t nSensParams() const override
Number of sensitivity parameters associated with this reactor.
string type() const override
String indicating the reactor model implemented.
double m_Qdot
net heat transfer into the reactor, through walls [W]
void setKinetics(Kinetics &kin) override
vector< double > m_advancelimits
!< Number of variables associated with reactor surfaces
virtual vector< size_t > steadyConstraints() const
Get the indices of equations that are algebraic constraints when solving the steady-state problem.
virtual double upperBound(size_t k) const
Get the upper bound on the k-th component of the local state vector.
virtual void setDerivativeSettings(AnyMap &settings)
Use this to set the kinetics objects derivative settings.
virtual void updateState(double *y)
Set the state of the reactor to correspond to the state vector y.
void setAdvanceLimits(const double *limits)
Set absolute step size limits during advance.
vector< Eigen::Triplet< double > > m_jac_trips
Vector of triplets representing the jacobian.
void setAdvanceLimit(const string &nm, const double limit)
Set individual step size limit for component name nm
virtual void resetSensitivity(double *params)
Reset the reaction rate multipliers.
void addSensitivityReaction(size_t rxn) override
Add a sensitivity parameter associated with the reaction number rxn
vector< double > m_sdot
Production rates of gas phase species on surfaces [kmol/s].
virtual void getState(double *y)
Get the the current state of the reactor.
bool hasAdvanceLimits() const
Check whether Reactor object uses advance limits.
double m_vdot
net rate of volume change from moving walls [m^3/s]
virtual void getSurfaceInitialConditions(double *y)
Get initial conditions for SurfPhase objects attached to this reactor.
void initialize(double t0=0.0) override
Initialize the reactor.
bool getAdvanceLimits(double *limits) const
Retrieve absolute step size limits during advance.
virtual size_t speciesIndex(const string &nm) const
Return the index in the solution vector for this reactor of the species named nm, in either the homog...
virtual bool timeIsIndependent() const
Indicates whether the governing equations for this reactor are functions of time or a spatial variabl...
virtual void resetBadValues(double *y)
Reset physically or mathematically problematic values, such as negative species concentrations.
virtual void updateConnected(bool updatePressure)
Update the state information needed by connected reactors, flow devices, and reactor walls.
A simple thermodynamic model for a surface phase, assuming an ideal solution model.
double size(size_t k) const
Returns the number of sites occupied by one molecule of species k.
double siteDensity() const
Returns the site density.
Base class for a phase with thermodynamic properties.
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 size_t npos
index returned by functions to indicate "no position"
const double Tiny
Small number to compare differences of mole fractions against.
offset
Offsets of solution components in the 1D solution array.
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.
const double BigNumber
largest number to compare to inf.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...