Cantera  3.1.0a2
Loading...
Searching...
No Matches
BulkKinetics.h
Go to the documentation of this file.
1/**
2 * @file BulkKinetics.h
3 */
4
5// This file is part of Cantera. See License.txt in the top-level directory or
6// at https://cantera.org/license.txt for license and copyright information.
7
8#ifndef CT_BULKKINETICS_H
9#define CT_BULKKINETICS_H
10
11#include "Kinetics.h"
12#include "MultiRate.h"
13#include "ThirdBodyCalc.h"
14
15namespace Cantera
16{
17
18//! Specialization of Kinetics for chemistry in a single bulk phase
19//! @ingroup kineticsmgr
20class BulkKinetics : public Kinetics
21{
22public:
23 //! @name Constructors and General Information
24 //! @{
26
27 string kineticsType() const override {
28 return "bulk";
29 }
30
31 bool isReversible(size_t i) override;
32 //! @}
33
34 //! @name Reaction Mechanism Setup Routines
35 //! @{
36 bool addReaction(shared_ptr<Reaction> r, bool resize=true) override;
37 void addThirdBody(shared_ptr<Reaction> r);
38 void modifyReaction(size_t i, shared_ptr<Reaction> rNew) override;
39 void resizeSpecies() override;
40 void resizeReactions() override;
41 void setMultiplier(size_t i, double f) override;
42 void invalidateCache() override;
43 //! @}
44
45 //! @name Reaction rate constants, rates of progress, and thermodynamic properties
46 //! @{
47 void getFwdRateConstants(double* kfwd) override;
48 void getEquilibriumConstants(double* kc) override;
49 void getRevRateConstants(double* krev, bool doIrreversible=false) override;
50
51 void getDeltaGibbs(double* deltaG) override;
52 void getDeltaEnthalpy(double* deltaH) override;
53 void getDeltaEntropy(double* deltaS) override;
54
55 void getDeltaSSGibbs(double* deltaG) override;
56 void getDeltaSSEnthalpy(double* deltaH) override;
57 void getDeltaSSEntropy(double* deltaS) override;
58 //! @}
59
60 //! @name Derivatives of rate constants and rates of progress
61 //! @{
62 void getDerivativeSettings(AnyMap& settings) const override;
63 void setDerivativeSettings(const AnyMap& settings) override;
64 void getFwdRateConstants_ddT(double* dkfwd) override;
65 void getFwdRatesOfProgress_ddT(double* drop) override;
66 void getRevRatesOfProgress_ddT(double* drop) override;
67 void getNetRatesOfProgress_ddT(double* drop) override;
68 void getFwdRateConstants_ddP(double* dkfwd) override;
69 void getFwdRatesOfProgress_ddP(double* drop) override;
70 void getRevRatesOfProgress_ddP(double* drop) override;
71 void getNetRatesOfProgress_ddP(double* drop) override;
72 void getFwdRateConstants_ddC(double* dkfwd) override;
73 void getFwdRatesOfProgress_ddC(double* drop) override;
74 void getRevRatesOfProgress_ddC(double* drop) override;
75 void getNetRatesOfProgress_ddC(double* drop) override;
76 Eigen::SparseMatrix<double> fwdRatesOfProgress_ddX() override;
77 Eigen::SparseMatrix<double> revRatesOfProgress_ddX() override;
78 Eigen::SparseMatrix<double> netRatesOfProgress_ddX() override;
79 Eigen::SparseMatrix<double> fwdRatesOfProgress_ddCi() override;
80 Eigen::SparseMatrix<double> revRatesOfProgress_ddCi() override;
81 Eigen::SparseMatrix<double> netRatesOfProgress_ddCi() override;
82 //! @}
83
84 //! @name Rate calculation intermediate methods
85 //! @{
86
87 void updateROP() override;
88
89 void getThirdBodyConcentrations(double* concm) override;
90 const vector<double>& thirdBodyConcentrations() const override {
91 return m_concm;
92 }
93
94 //! @}
95
96protected:
97 //! @name Internal service methods
98 //!
99 //! @note These methods are for internal use, and seek to avoid code duplication
100 //! while evaluating terms used for rate constants, rates of progress, and
101 //! their derivatives.
102 //! @{
103
104 //! Multiply rate with third-body collider concentrations
105 void processThirdBodies(double* rop);
106
107 //! Multiply rate with inverse equilibrium constant
108 void applyEquilibriumConstants(double* rop);
109
110 //! Multiply rate with scaled temperature derivatives of the inverse
111 //! equilibrium constant
112 /*!
113 * This (scaled) derivative is handled by a finite difference.
114 */
115 void applyEquilibriumConstants_ddT(double* drkcn);
116
117 //! Process temperature derivative
118 //! @param in rate expression used for the derivative calculation
119 //! @param drop pointer to output buffer
120 void process_ddT(const vector<double>& in, double* drop);
121
122 //! Process pressure derivative
123 //! @param in rate expression used for the derivative calculation
124 //! @param drop pointer to output buffer
125 void process_ddP(const vector<double>& in, double* drop);
126
127 //! Process concentration (molar density) derivative
128 //! @param stoich stoichiometry manager
129 //! @param in rate expression used for the derivative calculation
130 //! @param drop pointer to output buffer
131 //! @param mass_action boolean indicating whether law of mass action applies
132 void process_ddC(StoichManagerN& stoich, const vector<double>& in,
133 double* drop, bool mass_action=true);
134
135 //! Process derivatives
136 //! @param stoich stoichiometry manager
137 //! @param in rate expression used for the derivative calculation
138 //! @param ddX true: w.r.t mole fractions false: w.r.t species concentrations
139 //! @return a sparse matrix of derivative contributions for each reaction of
140 //! dimensions nTotalReactions by nTotalSpecies
141 Eigen::SparseMatrix<double> calculateCompositionDerivatives(
142 StoichManagerN& stoich, const vector<double>& in, bool ddX=true);
143
144 //! Helper function ensuring that all rate derivatives can be calculated
145 //! @param name method name used for error output
146 //! @throw CanteraError if ideal gas assumption does not hold
147 void assertDerivativesValid(const string& name);
148
149 //! @}
150
151 //! Vector of rate handlers
152 vector<unique_ptr<MultiRateBase>> m_bulk_rates;
153 map<string, size_t> m_bulk_types; //!< Mapping of rate handlers
154
155 vector<size_t> m_revindex; //!< Indices of reversible reactions
156 vector<size_t> m_irrev; //!< Indices of irreversible reactions
157
158 //! Difference between the global reactants order and the global products
159 //! order. Of type "double" to account for the fact that we can have real-
160 //! valued stoichiometries.
161 vector<double> m_dn;
162
163 ThirdBodyCalc m_multi_concm; //!< used with MultiRate evaluator
164
165 //! Third body concentrations
166 vector<double> m_concm;
167
168 //! Activity concentrations, as calculated by ThermoPhase::getActivityConcentrations
169 vector<double> m_act_conc;
170
171 //! Physical concentrations, as calculated by ThermoPhase::getConcentrations
172 vector<double> m_phys_conc;
173
174 //! Derivative settings
176 bool m_jac_skip_falloff;
177 double m_jac_rtol_delta;
178
179 bool m_ROP_ok = false;
180
181 //! Buffers for partial rop results with length nReactions()
182 vector<double> m_rbuf0;
183 vector<double> m_rbuf1;
184 vector<double> m_rbuf2;
185 vector<double> m_kf0; //!< Forward rate constants without perturbation
186 vector<double> m_sbuf0;
187 vector<double> m_state;
188 vector<double> m_grt; //!< Standard chemical potentials for each species
189};
190
191}
192
193#endif
Base class for kinetics managers and also contains the kineticsmgr module documentation (see Kinetics...
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:427
Specialization of Kinetics for chemistry in a single bulk phase.
bool addReaction(shared_ptr< Reaction > r, bool resize=true) override
Add a single reaction to the mechanism.
void getNetRatesOfProgress_ddC(double *drop) override
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
Eigen::SparseMatrix< double > netRatesOfProgress_ddX() override
Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant te...
void getNetRatesOfProgress_ddP(double *drop) override
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
string kineticsType() const override
Identifies the Kinetics manager type.
Eigen::SparseMatrix< double > revRatesOfProgress_ddCi() override
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
void getFwdRateConstants(double *kfwd) override
Return the forward rate constants.
void getDeltaSSGibbs(double *deltaG) override
Return the vector of values for the reaction standard state Gibbs free energy change.
map< string, size_t > m_bulk_types
Mapping of rate handlers.
Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX() override
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
const vector< double > & thirdBodyConcentrations() const override
Provide direct access to current third-body concentration values.
void getDeltaGibbs(double *deltaG) override
Return the vector of values for the reaction Gibbs free energy change.
vector< double > m_phys_conc
Physical concentrations, as calculated by ThermoPhase::getConcentrations.
Eigen::SparseMatrix< double > revRatesOfProgress_ddX() override
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
bool isReversible(size_t i) override
True if reaction i has been declared to be reversible.
void getFwdRatesOfProgress_ddT(double *drop) override
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
void assertDerivativesValid(const string &name)
Helper function ensuring that all rate derivatives can be calculated.
void getFwdRatesOfProgress_ddP(double *drop) override
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
void getFwdRateConstants_ddT(double *dkfwd) override
Calculate derivatives for forward rate constants with respect to temperature at constant pressure,...
void setMultiplier(size_t i, double f) override
Set the multiplier for reaction i to f.
void getDeltaSSEntropy(double *deltaS) override
Return the vector of values for the change in the standard state entropies for each reaction.
vector< double > m_dn
Difference between the global reactants order and the global products order.
bool m_jac_skip_third_bodies
Derivative settings.
void getDeltaSSEnthalpy(double *deltaH) override
Return the vector of values for the change in the standard state enthalpies of reaction.
Eigen::SparseMatrix< double > calculateCompositionDerivatives(StoichManagerN &stoich, const vector< double > &in, bool ddX=true)
Process derivatives.
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 getRevRateConstants(double *krev, bool doIrreversible=false) override
Return the reverse rate constants.
void getEquilibriumConstants(double *kc) override
Return a vector of Equilibrium constants.
void getNetRatesOfProgress_ddT(double *drop) override
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
vector< size_t > m_revindex
Indices of reversible reactions.
void processThirdBodies(double *rop)
Multiply rate with third-body collider concentrations.
void applyEquilibriumConstants(double *rop)
Multiply rate with inverse equilibrium constant.
void getDerivativeSettings(AnyMap &settings) const override
Retrieve derivative settings.
void process_ddT(const vector< double > &in, double *drop)
Process temperature derivative.
void getDeltaEnthalpy(double *deltaH) override
Return the vector of values for the reactions change in enthalpy.
vector< unique_ptr< MultiRateBase > > m_bulk_rates
Vector of rate handlers.
void getRevRatesOfProgress_ddT(double *drop) override
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
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 getFwdRatesOfProgress_ddC(double *drop) override
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
void process_ddC(StoichManagerN &stoich, const vector< double > &in, double *drop, bool mass_action=true)
Process concentration (molar density) derivative.
void process_ddP(const vector< double > &in, double *drop)
Process pressure derivative.
void getFwdRateConstants_ddP(double *dkfwd) override
Calculate derivatives for forward rate constants with respect to pressure at constant temperature,...
void getRevRatesOfProgress_ddP(double *drop) override
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
void getRevRatesOfProgress_ddC(double *drop) override
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
vector< double > m_rbuf0
Buffers for partial rop results with length nReactions()
vector< size_t > m_irrev
Indices of irreversible reactions.
ThirdBodyCalc m_multi_concm
used with MultiRate evaluator
void applyEquilibriumConstants_ddT(double *drkcn)
Multiply rate with scaled temperature derivatives of the inverse equilibrium constant.
void resizeReactions() override
Finalize Kinetics object and associated objects.
Eigen::SparseMatrix< double > netRatesOfProgress_ddCi() override
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
void getThirdBodyConcentrations(double *concm) override
Return a vector of values of effective concentrations of third-body collision partners of any reactio...
void getDeltaEntropy(double *deltaS) override
Return the vector of values for the reactions change in entropy.
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 getFwdRateConstants_ddC(double *dkfwd) override
Calculate derivatives for forward rate constants with respect to molar concentration at constant temp...
Public interface for kinetics managers.
Definition Kinetics.h:125
This class handles operations involving the stoichiometric coefficients on one side of a reaction (re...
Calculate and apply third-body effects on reaction rates, including non- unity third-body efficiencie...
Namespace for the Cantera kernel.
Definition AnyMap.cpp:564