7#include "cantera/numerics/eigen_dense.h"
12BulkKinetics::BulkKinetics() {
18 shared_ptr<ReactionRate> rate = r->rate();
20 rate->setContext(*r, *
this);
29 for (
const auto& [name, stoich] : r->products) {
32 for (
const auto& [name, stoich] : r->reactants) {
38 string rtype = rate->subType();
58 if (r->thirdBody() !=
nullptr) {
66void BulkKinetics::addThirdBody(shared_ptr<Reaction> r)
68 map<size_t, double> efficiencies;
69 for (
const auto& [name, efficiency] : r->thirdBody()->efficiencies) {
72 efficiencies[k] = efficiency;
74 throw CanteraError(
"BulkKinetics::addThirdBody",
"Found third-body"
75 " efficiency for undefined species '{}' while adding reaction '{}'",
82 r->thirdBody()->default_efficiency,
83 r->thirdBody()->mass_action);
91 shared_ptr<ReactionRate> rate = rNew->rate();
92 string rtype = rate->subType();
100 "Evaluator not available for type '{}'.", rtype);
105 rate->setRateIndex(i);
106 rate->setContext(*rNew, *
this);
118 rates->resize(*
this);
130 m_state.resize(
thermo().stateSize());
133 rates->resize(*
this);
146void BulkKinetics::invalidateCache()
148 Kinetics::invalidateCache();
156 copy(
m_rfn.begin(),
m_rfn.end(), kfwd.begin());
167 vector<double>& delta_gibbs0 =
m_rbuf0;
168 fill(delta_gibbs0.begin(), delta_gibbs0.end(), 0.0);
174 double logStandConc = log(
thermo().standardConcentration());
176 kc[i] = exp(-delta_gibbs0[i] * rrt +
m_dn[i] * logStandConc);
186 if (doIrreversible) {
238 for (
size_t k = 0; k <
m_kk; k++) {
251 for (
size_t k = 0; k <
m_kk; k++) {
261 settings[
"skip-falloff"] = m_jac_skip_falloff;
262 settings[
"rtol-delta"] = m_jac_rtol_delta;
267 bool force = settings.
empty();
268 if (force || settings.
hasKey(
"skip-third-bodies")) {
271 if (force || settings.
hasKey(
"skip-falloff")) {
272 m_jac_skip_falloff = settings.
getBool(
"skip-falloff",
false);
274 if (force || settings.
hasKey(
"rtol-delta")) {
275 m_jac_rtol_delta = settings.
getDouble(
"rtol-delta", 1e-8);
302 Eigen::Map<Eigen::VectorXd> dRevRop(drop.data(),
nReactions());
303 Eigen::Map<Eigen::VectorXd> dRevRop2(m_rbuf2.data(),
nReactions());
316 Eigen::Map<Eigen::VectorXd> dNetRop(drop.data(),
nReactions());
317 Eigen::Map<Eigen::VectorXd> dNetRop2(m_rbuf2.data(),
nReactions());
375 Eigen::Map<Eigen::VectorXd> dNetRop(drop.data(),
nReactions());
378 Eigen::Map<Eigen::VectorXd> dNetRop2(m_rbuf2.data(),
nReactions());
387 vector<double>& rop_rates =
m_rbuf0;
397 vector<double>& rop_rates =
m_rbuf0;
408 vector<double>& rop_rates =
m_rbuf0;
422 vector<double>& rop_rates =
m_rbuf0;
432 vector<double>& rop_rates =
m_rbuf0;
443 vector<double>& rop_rates =
m_rbuf0;
452void BulkKinetics::updateROP()
464 double logStandConc = log(
thermo().standardConcentration());
470 for (
size_t i = 0; i <
m_revindex.size(); i++) {
476 for (
size_t i = 0; i !=
m_irrev.size(); ++i) {
481 if (!last.
validate(T, rho, statenum)) {
494 bool changed = rates->update(
thermo(), *
this);
496 rates->getRateConstants(
m_kf0);
531 for (
size_t i = 0; i <
m_rfn.size(); i++) {
533 "m_rfn[{}] is not finite.", i);
535 "m_ropf[{}] is not finite.", i);
537 "m_ropr[{}] is not finite.", i);
572 vector<double>& grt = m_sbuf0;
573 vector<double>& delta_gibbs0 = m_rbuf1;
574 fill(delta_gibbs0.begin(), delta_gibbs0.end(), 0.0);
583 double Tinv = 1. / T;
584 double rrt_dTinv = rrt * Tinv / m_jac_rtol_delta;
585 double rrtt = rrt * Tinv;
586 for (
size_t i = 0; i <
m_revindex.size(); i++) {
591 drkcn[irxn] *= factor;
594 for (
size_t i = 0; i <
m_irrev.size(); ++i) {
604 copy(in.begin(), in.end(), drop.begin());
606 rates->processRateConstants_ddT(drop,
m_rfn, m_jac_rtol_delta);
613 copy(in.begin(), in.end(), drop.begin());
615 rates->processRateConstants_ddP(drop,
m_rfn, m_jac_rtol_delta);
620 span<double> drop,
bool mass_action)
622 Eigen::Map<Eigen::VectorXd> out(drop.data(),
nReactions());
635 Eigen::Map<Eigen::VectorXd> outM(m_rbuf1.data(),
nReactions());
643 if (!m_jac_skip_falloff) {
647 rates->processRateConstants_ddM(
asSpan(outM),
m_rfn, m_jac_rtol_delta);
656 Eigen::SparseMatrix<double> out;
657 vector<double>& scaled = m_rbuf1;
658 vector<double>& outV = m_rbuf2;
661 copy(in.begin(), in.end(), scaled.begin());
670 copy(scaled.begin(), scaled.end(), outV.begin());
678 copy(scaled.begin(), scaled.end(), outV.begin());
682 if (!m_jac_skip_falloff) {
685 rates->processRateConstants_ddM(outV,
m_rfn, m_jac_rtol_delta,
false);
696 if (!
thermo().isIdeal()) {
698 "Not supported for non-ideal ThermoPhase models.");
Header file for class ThermoPhase, the base class for phases with thermodynamic properties,...
A map of string keys to values whose type can vary at runtime.
double getDouble(const string &key, double default_) const
If key exists, return it as a double, otherwise return default_.
bool hasKey(const string &key) const
Returns true if the map contains an item named key.
bool empty() const
Return boolean indicating whether AnyMap is empty.
bool getBool(const string &key, bool default_) const
If key exists, return it as a bool, otherwise return default_.
bool addReaction(shared_ptr< Reaction > r, bool resize=true) override
Add a single reaction to the mechanism.
Eigen::SparseMatrix< double > netRatesOfProgress_ddX() override
Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant te...
void applyEquilibriumConstants(span< double > rop)
Multiply rate with inverse equilibrium constant.
Eigen::SparseMatrix< double > revRatesOfProgress_ddCi() override
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
void getRevRatesOfProgress_ddP(span< double > drop) override
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
void process_ddT(const span< const double > in, span< double > drop)
Process temperature derivative.
Eigen::SparseMatrix< double > calculateCompositionDerivatives(StoichManagerN &stoich, span< const double > in, bool ddX=true)
Process derivatives.
void getDeltaEnthalpy(span< double > deltaH) override
Return the vector of values for the reactions change in enthalpy.
void getEquilibriumConstants(span< double > kc) override
Return a vector of Equilibrium constants.
Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX() override
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
void getFwdRatesOfProgress_ddT(span< double > drop) override
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
vector< double > m_phys_conc
Physical concentrations, as calculated by ThermoPhase::getConcentrations.
void getRevRatesOfProgress_ddC(span< double > drop) override
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
Eigen::SparseMatrix< double > revRatesOfProgress_ddX() override
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
void getNetRatesOfProgress_ddT(span< double > drop) override
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
void getNetRatesOfProgress_ddP(span< double > drop) override
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
void assertDerivativesValid(const string &name)
Helper function ensuring that all rate derivatives can be calculated.
void getFwdRateConstants(span< double > kfwd) override
Return the forward rate constants.
void process_ddC(StoichManagerN &stoich, span< const double > in, span< double > drop, bool mass_action=true)
Process concentration (molar density) derivative.
void applyEquilibriumConstants_ddT(span< double > drkcn)
Multiply rate with scaled temperature derivatives of the inverse equilibrium constant.
void setMultiplier(size_t i, double f) override
Set the multiplier for reaction i to f.
vector< double > m_dn
Difference between the global reactants order and the global products order.
void process_ddP(const span< const double > in, span< double > drop)
Process pressure derivative.
void getDeltaSSGibbs(span< double > deltaG) override
Return the vector of values for the reaction standard state Gibbs free energy change.
void processThirdBodies(span< double > rop)
Multiply rate with third-body collider concentrations.
bool m_jac_skip_third_bodies
Derivative settings.
void resizeSpecies() override
Resize arrays with sizes that depend on the total number of species.
vector< double > m_grt
Standard chemical potentials for each species.
vector< double > m_act_conc
Activity concentrations, as calculated by ThermoPhase::getActivityConcentrations.
void getFwdRateConstants_ddC(span< double > dkfwd) override
Calculate derivatives for forward rate constants with respect to molar concentration at constant temp...
void getDerivativeSettings(AnyMap &settings) const override
Retrieve derivative settings.
void getNetRatesOfProgress_ddC(span< double > drop) override
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi() override
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
void modifyReaction(size_t i, shared_ptr< Reaction > rNew) override
Modify the rate expression associated with a reaction.
void getRevRateConstants(span< double > krev, bool doIrreversible=false) override
Return the reverse rate constants.
vector< double > m_rbuf0
Buffers for partial rop results with length nReactions()
void getDeltaSSEnthalpy(span< double > deltaH) override
Return the vector of values for the change in the standard state enthalpies of reaction.
void getFwdRatesOfProgress_ddP(span< double > drop) override
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
void getDeltaEntropy(span< double > deltaS) override
Return the vector of values for the reactions change in entropy.
ThirdBodyCalc m_multi_concm
used with MultiRate evaluator
void getDeltaSSEntropy(span< double > deltaS) override
Return the vector of values for the change in the standard state entropies for each reaction.
void resizeReactions() override
Finalize Kinetics object and associated objects.
void getFwdRateConstants_ddP(span< double > dkfwd) override
Calculate derivatives for forward rate constants with respect to pressure at constant temperature,...
void getRevRatesOfProgress_ddT(span< double > drop) override
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
void getThirdBodyConcentrations(span< double > concm) override
Return a vector of values of effective concentrations of third-body collision partners of any reactio...
Eigen::SparseMatrix< double > netRatesOfProgress_ddCi() override
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
void getFwdRateConstants_ddT(span< double > dkfwd) override
Calculate derivatives for forward rate constants with respect to temperature at constant pressure,...
void getDeltaGibbs(span< double > deltaG) override
Return the vector of values for the reaction Gibbs free energy change.
vector< double > m_concm
Third body concentrations.
vector< double > m_kf0
Forward rate constants without perturbation.
void setDerivativeSettings(const AnyMap &settings) override
Set/modify derivative settings.
void getFwdRatesOfProgress_ddC(span< double > drop) override
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
Base class for exceptions thrown by Cantera classes.
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
ValueCache m_cache
Cache for saved calculations within each Kinetics object.
vector< double > m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
vector< double > m_ropf
Forward rate-of-progress for each reaction.
vector< unique_ptr< MultiRateBase > > m_rateHandlers
Vector of rate handlers.
vector< double > m_rkcn
Reciprocal of the equilibrium constant in concentration units.
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
vector< size_t > m_revindex
Indices of reversible reactions.
bool m_skipUndeclaredThirdBodies
See skipUndeclaredThirdBodies()
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
map< string, size_t > m_rateTypes
Mapping of rate handlers.
vector< double > m_ropnet
Net rate-of-progress for each reaction.
virtual void getRevReactionDelta(span< const double > g, span< double > dg) const
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the rev...
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
size_t nReactions() const
Number of reactions in the reaction mechanism.
vector< size_t > m_irrev
Indices of irreversible reactions.
bool m_hasUndeclaredThirdBodies
Flag indicating whether reactions include undeclared third bodies.
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
virtual void getReactionDelta(span< const double > property, span< double > deltaProperty) const
Change in species properties.
virtual void setMultiplier(size_t i, double f)
Set the multiplier for reaction i to f.
vector< double > m_rfn
Forward rate constant for each reaction.
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
vector< double > m_delta_gibbs0
Delta G^0 for all reactions.
An error indicating that an unimplemented function has been called.
virtual double molarDensity() const
Molar density (kmol/m^3).
double temperature() const
Temperature (K).
virtual void getConcentrations(span< double > c) const
Get the species concentrations (kmol/m^3).
virtual double density() const
Density (kg/m^3).
int stateMFNumber() const
Return the State Mole Fraction Number.
virtual void restoreState(span< const double > state)
Restore the state of the phase from a previously saved state vector.
virtual double pressure() const
Return the thermodynamic pressure (Pa).
virtual void saveState(span< double > state) const
Write to array 'state' the current internal state.
This class handles operations involving the stoichiometric coefficients on one side of a reaction (re...
Eigen::SparseMatrix< double > derivatives(span< const double > conc, span< const double > rates)
Calculate derivatives with respect to species concentrations.
void scale(span< const double > in, span< double > out, double factor) const
Scale input by reaction order and factor.
virtual void getPartialMolarEnthalpies(span< double > hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
virtual void setState_TP(double t, double p)
Set the temperature (K) and pressure (Pa)
double RT() const
Return the Gas Constant multiplied by the current temperature.
virtual void getPartialMolarEntropies(span< double > sbar) const
Returns an array of partial molar entropies of the species in the solution.
virtual void getStandardChemPotentials(span< double > mu) const
Get the array of chemical potentials at unit activity for the species at their standard states at the...
virtual void getActivityConcentrations(span< double > c) const
This method returns an array of generalized concentrations.
virtual void getEntropy_R(span< double > sr) const
Get the array of nondimensional Entropy functions for the standard state species at the current T and...
virtual void getEnthalpy_RT(span< double > hrt) const
Get the nondimensional Enthalpy functions for the species at their standard states at the current T a...
virtual void getChemPotentials(span< double > mu) const
Get the species chemical potentials. Units: J/kmol.
bool empty() const
Return boolean indicating whether ThirdBodyCalc is empty.
void install(size_t rxnNumber, const map< size_t, double > &efficiencies, double default_efficiency, bool mass_action)
Install reaction that uses third-body effects in ThirdBodyCalc manager.
void update(span< const double > conc, double ctot, span< double > concm) const
Update third-body concentrations in full vector.
Eigen::SparseMatrix< double > derivatives(span< const double > product)
Calculate derivatives with respect to species concentrations.
void multiply(span< double > output, span< const double > concm)
Multiply output with effective third-body concentration.
void scale(span< const double > in, span< double > out, double factor) const
Scale entries involving third-body collider in law of mass action by factor.
void resizeCoeffs(size_t nSpc, size_t nRxn)
Resize the sparse coefficient matrix.
void scaleM(span< const double > in, span< double > out, span< const double > concm, double factor) const
Scale entries involving third-body collider in rate expression by third-body concentration and factor...
CachedScalar getScalar(int id)
Get a reference to a CachedValue object representing a scalar (double) with the given id.
int getId()
Get a unique id for a cached value.
#define AssertFinite(expr, procedure,...)
Throw an exception if the specified exception is not a finite number.
bool legacy_rate_constants_used()
Returns true if legacy rate constant definition is used.
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"
span< double > asSpan(Eigen::DenseBase< Derived > &v)
Convenience wrapper for accessing Eigen vector/array/map data as a span.
const double BigNumber
largest number to compare to inf.
void checkArraySize(const char *procedure, size_t available, size_t required)
Wrapper for throwing ArraySizeError.
A cached property value and the state at which it was evaluated.
double state2
Value of the second state variable for the state at which value was evaluated, for example density or...
bool validate(double state1New)
Check whether the currently cached value is valid based on a single state variable.
double state1
Value of the first state variable for the state at which value was evaluated, for example temperature...