Cantera  3.0.0
Loading...
Searching...
No Matches
ThirdBodyCalc.h
Go to the documentation of this file.
1/**
2 * @file ThirdBodyCalc.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_THIRDBODYCALC_H
9#define CT_THIRDBODYCALC_H
10
12
13namespace Cantera
14{
15
16//! Calculate and apply third-body effects on reaction rates, including non-
17//! unity third-body efficiencies.
18//! @ingroup rateEvaluators
20{
21public:
22 //! Install reaction that uses third-body effects in ThirdBodyCalc manager
23 void install(size_t rxnNumber, const map<size_t, double>& efficiencies,
24 double default_efficiency, bool mass_action) {
25 m_reaction_index.push_back(rxnNumber);
26 m_default.push_back(default_efficiency);
27
28 if (mass_action) {
29 m_mass_action_index.push_back(m_reaction_index.size() - 1);
30 } else {
31 m_no_mass_action_index.push_back(m_reaction_index.size() - 1);
32 }
33
34 m_species.emplace_back();
35 m_eff.emplace_back();
36 for (const auto& [k, efficiency] : efficiencies) {
37 AssertTrace(k != npos);
38 m_species.back().push_back(k);
39 m_eff.back().push_back(efficiency - default_efficiency);
40 m_efficiencyList.emplace_back(
41 static_cast<int>(rxnNumber),
42 static_cast<int>(k), efficiency - default_efficiency);
43 }
44 }
45
46 //! Resize the sparse coefficient matrix
47 void resizeCoeffs(size_t nSpc, size_t nRxn) {
48 // Sparse Efficiency coefficient matrix
49 Eigen::SparseMatrix<double> efficiencies;
50 efficiencies.setZero();
51 efficiencies.resize(nRxn, nSpc);
52 efficiencies.reserve(m_efficiencyList.size());
53 efficiencies.setFromTriplets(
54 m_efficiencyList.begin(), m_efficiencyList.end());
55
56 // derivative matrix multipliers
57 vector<Eigen::Triplet<double>> triplets;
58 triplets.reserve(m_reaction_index.size() * nSpc);
59 for (size_t i = 0; i < m_default.size(); i++) {
60 if (m_default[i] != 0) {
61 for (size_t j = 0; j < nSpc; j++) {
62 triplets.emplace_back(
63 static_cast<int>(m_reaction_index[i]),
64 static_cast<int>(j), m_default[i]);
65 }
66 }
67 }
68 Eigen::SparseMatrix<double> defaults(nRxn, nSpc);
69 defaults.reserve(triplets.size());
70 defaults.setFromTriplets(triplets.begin(), triplets.end());
71 m_multipliers = efficiencies + defaults;
72 }
73
74 //! Update third-body concentrations in full vector
75 void update(const vector<double>& conc, double ctot, double* concm) const {
76 for (size_t i = 0; i < m_reaction_index.size(); i++) {
77 double sum = 0.0;
78 for (size_t j = 0; j < m_species[i].size(); j++) {
79 sum += m_eff[i][j] * conc[m_species[i][j]];
80 }
81 concm[m_reaction_index[i]] = m_default[i] * ctot + sum;
82 }
83 }
84
85 //! Multiply output with effective third-body concentration
86 void multiply(double* output, const double* concm) {
87 for (size_t i = 0; i < m_mass_action_index.size(); i++) {
89 output[ix] *= concm[ix];
90 }
91 }
92
93 //! Calculate derivatives with respect to species concentrations.
94 /*!
95 * @param product Product of law of mass action and rate terms.
96 */
97 Eigen::SparseMatrix<double> derivatives(const double* product) {
98 Eigen::Map<const Eigen::VectorXd> mapped(product, m_multipliers.rows());
99 return mapped.asDiagonal() * m_multipliers;
100 }
101
102 //! Scale entries involving third-body collider in law of mass action by factor
103 void scale(const double* in, double* out, double factor) const {
104 for (size_t i = 0; i < m_mass_action_index.size(); i++) {
106 out[ix] = factor * in[ix];
107 }
108 }
109
110 //! Scale entries involving third-body collider in rate expression
111 //! by third-body concentration and factor
112 void scaleM(const double* in, double* out,
113 const double* concm, double factor) const
114 {
115 for (size_t i = 0; i < m_no_mass_action_index.size(); i++) {
117 out[ix] = factor * concm[ix] * in[ix];
118 }
119 }
120
121 //! Return boolean indicating whether ThirdBodyCalc is empty
122 bool empty() const {
123 return m_reaction_index.empty();
124 }
125
126protected:
127 //! Indices of reactions that use third-bodies within vector of concentrations
128 vector<size_t> m_reaction_index;
129
130 //! Indices within m_reaction_index of reactions that consider third-body effects
131 //! in the law of mass action
132 vector<size_t> m_mass_action_index;
133
134 //! Indices within m_reaction_index of reactions that consider third-body effects
135 //! in the rate expression
137
138 //! m_species[i][j] is the index of the j-th species in reaction i.
139 vector<vector<size_t>> m_species;
140
141 //! m_eff[i][j] is the efficiency of the j-th species in reaction i.
142 vector<vector<double>> m_eff;
143
144 //! The default efficiency for each reaction
145 vector<double> m_default;
146
147 //! Sparse efficiency matrix (compensated for defaults)
148 //! Each triplet corresponds to (reaction index, species index, efficiency)
149 vector<Eigen::Triplet<double>> m_efficiencyList;
150
151 //! Sparse derivative multiplier matrix
152 Eigen::SparseMatrix<double> m_multipliers;
153};
154
155using ThirdBodyCalc3 = ThirdBodyCalc; // @todo: remove after Cantera 3.0
156
157}
158
159#endif
Calculate and apply third-body effects on reaction rates, including non- unity third-body efficiencie...
void scaleM(const double *in, double *out, const double *concm, double factor) const
Scale entries involving third-body collider in rate expression by third-body concentration and factor...
Eigen::SparseMatrix< double > m_multipliers
Sparse derivative multiplier matrix.
vector< Eigen::Triplet< double > > m_efficiencyList
Sparse efficiency matrix (compensated for defaults) Each triplet corresponds to (reaction index,...
vector< vector< size_t > > m_species
m_species[i][j] is the index of the j-th species in reaction i.
vector< vector< double > > m_eff
m_eff[i][j] is the efficiency of the j-th species in reaction i.
vector< double > m_default
The default efficiency for each reaction.
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.
vector< size_t > m_no_mass_action_index
Indices within m_reaction_index of reactions that consider third-body effects in the rate expression.
void multiply(double *output, const double *concm)
Multiply output with effective third-body concentration.
void update(const vector< double > &conc, double ctot, double *concm) const
Update third-body concentrations in full vector.
Eigen::SparseMatrix< double > derivatives(const double *product)
Calculate derivatives with respect to species concentrations.
void scale(const double *in, double *out, double factor) const
Scale entries involving third-body collider in law of mass action by factor.
vector< size_t > m_mass_action_index
Indices within m_reaction_index of reactions that consider third-body effects in the law of mass acti...
void resizeCoeffs(size_t nSpc, size_t nRxn)
Resize the sparse coefficient matrix.
vector< size_t > m_reaction_index
Indices of reactions that use third-bodies within vector of concentrations.
This file contains definitions of constants, types and terms that are used in internal routines and a...
#define AssertTrace(expr)
Assertion must be true or an error is thrown.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:564
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:195