From 761f54f24151108b3cb9b3506fe12df8026a6142 Mon Sep 17 00:00:00 2001 From: Ray Speth Date: Thu, 16 Jul 2015 15:37:41 -0400 Subject: [PATCH] Eliminate unnecessary counter variable Kinetics::m_ii --- include/cantera/kinetics/InterfaceKinetics.h | 4 +- include/cantera/kinetics/Kinetics.h | 47 +++++++++----------- src/kinetics/AqueousKinetics.cpp | 6 +-- src/kinetics/BulkKinetics.cpp | 6 +-- src/kinetics/GasKinetics.cpp | 6 +-- src/kinetics/InterfaceKinetics.cpp | 24 +++++----- src/kinetics/Kinetics.cpp | 16 +++---- 7 files changed, 51 insertions(+), 58 deletions(-) diff --git a/include/cantera/kinetics/InterfaceKinetics.h b/include/cantera/kinetics/InterfaceKinetics.h index d48fb0d8dd..15c25fb7ba 100644 --- a/include/cantera/kinetics/InterfaceKinetics.h +++ b/include/cantera/kinetics/InterfaceKinetics.h @@ -114,8 +114,8 @@ class InterfaceKinetics : public Kinetics * * units = J kmol-1 * - * @param deltaG Output vector of deltaG's for reactions Length: m_ii. - * If 0, this updates the internally stored values only. + * @param deltaG Output vector of deltaG's for reactions Length: + * nReactions(). If 0, this updates the internally stored values only. */ virtual void getDeltaGibbs(doublereal* deltaG); diff --git a/include/cantera/kinetics/Kinetics.h b/include/cantera/kinetics/Kinetics.h index 1d2b56833a..015c1fbf2b 100644 --- a/include/cantera/kinetics/Kinetics.h +++ b/include/cantera/kinetics/Kinetics.h @@ -189,7 +189,7 @@ class Kinetics //! Number of reactions in the reaction mechanism. size_t nReactions() const { - return m_ii; + return m_reactions.size(); } //! Check that the specified reaction index is in range @@ -401,7 +401,7 @@ class Kinetics * least as large as the total number of reactions. * * @param fwdROP Output vector containing forward rates - * of progress of the reactions. Length: m_ii. + * of progress of the reactions. Length: nReactions(). */ virtual void getFwdRatesOfProgress(doublereal* fwdROP); @@ -411,7 +411,7 @@ class Kinetics * dimensioned at least as large as the total number of reactions. * * @param revROP Output vector containing reverse rates - * of progress of the reactions. Length: m_ii. + * of progress of the reactions. Length: nReactions(). */ virtual void getRevRatesOfProgress(doublereal* revROP); @@ -420,7 +420,7 @@ class Kinetics * progress in array netROP, which must be dimensioned at least as large * as the total number of reactions. * - * @param netROP Output vector of the net ROP. Length: m_ii. + * @param netROP Output vector of the net ROP. Length: nReactions(). */ virtual void getNetRatesOfProgress(doublereal* netROP); @@ -435,7 +435,7 @@ class Kinetics * \f] * * @param kc Output vector containing the equilibrium constants. - * Length: m_ii. + * Length: nReactions(). */ virtual void getEquilibriumConstants(doublereal* kc) { throw NotImplementedError("Kinetics::getEquilibriumConstants"); @@ -455,7 +455,7 @@ class Kinetics * reaction. * * @param property Input vector of property value. Length: m_kk. - * @param deltaProperty Output vector of deltaRxn. Length: m_ii. + * @param deltaProperty Output vector of deltaRxn. Length: nReactions(). */ virtual void getReactionDelta(const doublereal* property, doublereal* deltaProperty); @@ -479,7 +479,8 @@ class Kinetics * * units = J kmol-1 * - * @param deltaG Output vector of deltaG's for reactions Length: m_ii. + * @param deltaG Output vector of deltaG's for reactions Length: + * nReactions(). */ virtual void getDeltaGibbs(doublereal* deltaG) { throw NotImplementedError("Kinetics::getDeltaGibbs"); @@ -493,7 +494,8 @@ class Kinetics * * units = J kmol-1 * - * @param deltaM Output vector of deltaM's for reactions Length: m_ii. + * @param deltaM Output vector of deltaM's for reactions Length: + * nReactions(). */ virtual void getDeltaElectrochemPotentials(doublereal* deltaM) { throw NotImplementedError("Kinetics::getDeltaElectrochemPotentials"); @@ -505,7 +507,8 @@ class Kinetics * * units = J kmol-1 * - * @param deltaH Output vector of deltaH's for reactions Length: m_ii. + * @param deltaH Output vector of deltaH's for reactions Length: + * nReactions(). */ virtual void getDeltaEnthalpy(doublereal* deltaH) { throw NotImplementedError("Kinetics::getDeltaEnthalpy"); @@ -517,7 +520,8 @@ class Kinetics * * units = J kmol-1 Kelvin-1 * - * @param deltaS Output vector of deltaS's for reactions Length: m_ii. + * @param deltaS Output vector of deltaS's for reactions Length: + * nReactions(). */ virtual void getDeltaEntropy(doublereal* deltaS) { throw NotImplementedError("Kinetics::getDeltaEntropy"); @@ -530,7 +534,8 @@ class Kinetics * * units = J kmol-1 * - * @param deltaG Output vector of ss deltaG's for reactions Length: m_ii. + * @param deltaG Output vector of ss deltaG's for reactions Length: + * nReactions(). */ virtual void getDeltaSSGibbs(doublereal* deltaG) { throw NotImplementedError("Kinetics::getDeltaSSGibbs"); @@ -543,7 +548,8 @@ class Kinetics * * units = J kmol-1 * - * @param deltaH Output vector of ss deltaH's for reactions Length: m_ii. + * @param deltaH Output vector of ss deltaH's for reactions Length: + * nReactions(). */ virtual void getDeltaSSEnthalpy(doublereal* deltaH) { throw NotImplementedError("Kinetics::getDeltaSSEnthalpy"); @@ -556,7 +562,8 @@ class Kinetics * * units = J kmol-1 Kelvin-1 * - * @param deltaS Output vector of ss deltaS's for reactions Length: m_ii. + * @param deltaS Output vector of ss deltaS's for reactions Length: + * nReactions(). */ virtual void getDeltaSSEntropy(doublereal* deltaS) { throw NotImplementedError("Kinetics::getDeltaSSEntropy"); @@ -695,7 +702,7 @@ class Kinetics * length is the number of reactions. units depends on many issues. * * @param kfwd Output vector containing the forward reaction rate - * constants. Length: m_ii. + * constants. Length: nReactions(). */ virtual void getFwdRateConstants(doublereal* kfwd) { throw NotImplementedError("Kinetics::getFwdRateConstants"); @@ -842,15 +849,6 @@ class Kinetics //@} - /** - * Increment the number of reactions in the mechanism by one. - * @todo Should be protected? - */ - void incrementRxnCount() { - m_ii++; - m_perturb.push_back(1.0); - } - /** * Returns true if the kinetics manager has been properly * initialized and finalized. @@ -929,9 +927,6 @@ class Kinetics StoichManagerN m_irrevProductStoich; //@} - //! Number of reactions in the mechanism - size_t m_ii; - //! The number of species in all of the phases //! that participate in this kinetics mechanism. size_t m_kk; diff --git a/src/kinetics/AqueousKinetics.cpp b/src/kinetics/AqueousKinetics.cpp index ddc2df0e9f..a3a6e21a24 100644 --- a/src/kinetics/AqueousKinetics.cpp +++ b/src/kinetics/AqueousKinetics.cpp @@ -87,7 +87,7 @@ void AqueousKinetics::getEquilibriumConstants(doublereal* kc) getReactionDelta(&m_grt[0], &m_rkcn[0]); doublereal rrt = 1.0/(GasConstant * thermo().temperature()); - for (size_t i = 0; i < m_ii; i++) { + for (size_t i = 0; i < nReactions(); i++) { kc[i] = exp(-m_rkcn[i]*rrt); } @@ -124,7 +124,7 @@ void AqueousKinetics::updateROP() // for reversible reactions, multiply ropr by concentration products m_revProductStoich.multiply(&m_conc[0], &m_ropr[0]); - for (size_t j = 0; j != m_ii; ++j) { + for (size_t j = 0; j != nReactions(); ++j) { m_ropnet[j] = m_ropf[j] - m_ropr[j]; } @@ -142,7 +142,7 @@ void AqueousKinetics::getFwdRateConstants(doublereal* kfwd) // multiply by perturbation factor multiply_each(m_ropf.begin(), m_ropf.end(), m_perturb.begin()); - for (size_t i = 0; i < m_ii; i++) { + for (size_t i = 0; i < nReactions(); i++) { kfwd[i] = m_ropf[i]; } } diff --git a/src/kinetics/BulkKinetics.cpp b/src/kinetics/BulkKinetics.cpp index 9097d20374..f7c6dd52f5 100644 --- a/src/kinetics/BulkKinetics.cpp +++ b/src/kinetics/BulkKinetics.cpp @@ -104,12 +104,12 @@ void BulkKinetics::getRevRateConstants(doublereal* krev, bool doIrreversible) if (doIrreversible) { getEquilibriumConstants(&m_ropnet[0]); - for (size_t i = 0; i < m_ii; i++) { + for (size_t i = 0; i < nReactions(); i++) { krev[i] /= m_ropnet[i]; } } else { // m_rkcn[] is zero for irreversible reactions - for (size_t i = 0; i < m_ii; i++) { + for (size_t i = 0; i < nReactions(); i++) { krev[i] *= m_rkcn[i]; } } @@ -168,7 +168,7 @@ void BulkKinetics::finalize() // Guarantee that these arrays can be converted to double* even in the // special case where there are no reactions defined. - if (!m_ii) { + if (!nReactions()) { m_perturb.resize(1, 1.0); m_ropf.resize(1, 0.0); m_ropr.resize(1, 0.0); diff --git a/src/kinetics/GasKinetics.cpp b/src/kinetics/GasKinetics.cpp index b46e6bc927..562da7d94e 100644 --- a/src/kinetics/GasKinetics.cpp +++ b/src/kinetics/GasKinetics.cpp @@ -127,7 +127,7 @@ void GasKinetics::getEquilibriumConstants(doublereal* kc) getReactionDelta(&m_grt[0], &m_rkcn[0]); doublereal rrt = 1.0/(GasConstant * thermo().temperature()); - for (size_t i = 0; i < m_ii; i++) { + for (size_t i = 0; i < nReactions(); i++) { kc[i] = exp(-m_rkcn[i]*rrt + m_dn[i]*m_logStandConc); } @@ -199,7 +199,7 @@ void GasKinetics::updateROP() // for reversible reactions, multiply ropr by concentration products m_revProductStoich.multiply(&m_conc[0], &m_ropr[0]); - for (size_t j = 0; j != m_ii; ++j) { + for (size_t j = 0; j != nReactions(); ++j) { m_ropnet[j] = m_ropf[j] - m_ropr[j]; } @@ -235,7 +235,7 @@ void GasKinetics::getFwdRateConstants(doublereal* kfwd) // multiply by perturbation factor multiply_each(m_ropf.begin(), m_ropf.end(), m_perturb.begin()); - for (size_t i = 0; i < m_ii; i++) { + for (size_t i = 0; i < nReactions(); i++) { kfwd[i] = m_ropf[i]; } } diff --git a/src/kinetics/InterfaceKinetics.cpp b/src/kinetics/InterfaceKinetics.cpp index a227f27a2d..48ed3ce724 100644 --- a/src/kinetics/InterfaceKinetics.cpp +++ b/src/kinetics/InterfaceKinetics.cpp @@ -294,11 +294,11 @@ void InterfaceKinetics::getEquilibriumConstants(doublereal* kc) updateMu0(); doublereal rrt = 1.0 / (GasConstant * thermo(0).temperature()); - std::fill(kc, kc + m_ii, 0.0); + std::fill(kc, kc + nReactions(), 0.0); getReactionDelta(DATA_PTR(m_mu0_Kc), kc); - for (size_t i = 0; i < m_ii; i++) { + for (size_t i = 0; i < nReactions(); i++) { kc[i] = exp(-kc[i]*rrt); } } @@ -333,7 +333,7 @@ void InterfaceKinetics::updateExchangeCurrentQuantities() getReactionDelta(DATA_PTR(m_mu0), DATA_PTR(m_deltaG0)); // Calculate the product of the standard concentrations of the reactants - for (size_t i = 0; i < m_ii; i++) { + for (size_t i = 0; i < nReactions(); i++) { m_ProdStanConcReac[i] = 1.0; } m_reactantStoich.multiply(DATA_PTR(m_StandardConc), DATA_PTR(m_ProdStanConcReac)); @@ -463,7 +463,7 @@ void InterfaceKinetics::getRevRateConstants(doublereal* krev, bool doIrreversibl getFwdRateConstants(krev); if (doIrreversible) { getEquilibriumConstants(&m_ropnet[0]); - for (size_t i = 0; i < m_ii; i++) { + for (size_t i = 0; i < nReactions(); i++) { krev[i] /= m_ropnet[i]; } } else { @@ -506,7 +506,7 @@ void InterfaceKinetics::updateROP() // Fix up these calculations for cases where the above formalism doesn't hold double OCV = 0.0; - for (size_t jrxn = 0; jrxn != m_ii; ++jrxn) { + for (size_t jrxn = 0; jrxn != nReactions(); ++jrxn) { int reactionType = m_rxntype[jrxn]; if (reactionType == BUTLERVOLMER_RXN) { // @@ -543,7 +543,7 @@ void InterfaceKinetics::updateROP() } } - for (size_t j = 0; j != m_ii; ++j) { + for (size_t j = 0; j != nReactions(); ++j) { m_ropnet[j] = m_ropf[j] - m_ropr[j]; } @@ -553,7 +553,7 @@ void InterfaceKinetics::updateROP() * phases that are stoichiometric phases containing one species with a unity activity */ if (m_phaseExistsCheck) { - for (size_t j = 0; j != m_ii; ++j) { + for (size_t j = 0; j != nReactions(); ++j) { if ((m_ropr[j] > m_ropf[j]) && (m_ropr[j] > 0.0)) { for (size_t p = 0; p < nPhases(); p++) { if (m_rxnPhaseIsProduct[j][p]) { @@ -624,7 +624,7 @@ void InterfaceKinetics::getDeltaGibbs(doublereal* deltaG) // Use the stoichiometric manager to find deltaG for each reaction. getReactionDelta(DATA_PTR(m_mu), DATA_PTR(m_deltaG)); if (deltaG != 0 && (DATA_PTR(m_deltaG) != deltaG)) { - for (size_t j = 0; j < m_ii; ++j) { + for (size_t j = 0; j < nReactions(); ++j) { deltaG[j] = m_deltaG[j]; } } @@ -967,7 +967,7 @@ void InterfaceKinetics::init() void InterfaceKinetics::finalize() { Kinetics::finalize(); - size_t safe_reaction_size = std::max(m_ii, 1); + size_t safe_reaction_size = std::max(nReactions(), 1); deltaElectricEnergy_.resize(safe_reaction_size); size_t ks = reactionPhaseIndex(); if (ks == npos) throw CanteraError("InterfaceKinetics::finalize", @@ -992,7 +992,7 @@ void InterfaceKinetics::finalize() // Guarantee that these arrays can be converted to double* even in the // special case where there are no reactions defined. - if (!m_ii) { + if (!nReactions()) { m_perturb.resize(1, 1.0); m_ropf.resize(1, 0.0); m_ropr.resize(1, 0.0); @@ -1170,7 +1170,7 @@ void EdgeKinetics::finalize() // handled by the InterfaceKinetics::finalize() call. Kinetics::finalize(); - size_t safe_reaction_size = std::max(m_ii, 1); + size_t safe_reaction_size = std::max(nReactions(), 1); deltaElectricEnergy_.resize(safe_reaction_size); size_t ks = reactionPhaseIndex(); if (ks == npos) throw CanteraError("EdgeKinetics::finalize", @@ -1195,7 +1195,7 @@ void EdgeKinetics::finalize() // Guarantee that these arrays can be converted to double* even in the // special case where there are no reactions defined. - if (!m_ii) { + if (!nReactions()) { m_perturb.resize(1, 1.0); m_ropf.resize(1, 0.0); m_ropr.resize(1, 0.0); diff --git a/src/kinetics/Kinetics.cpp b/src/kinetics/Kinetics.cpp index d6a5763801..76c85bc308 100644 --- a/src/kinetics/Kinetics.cpp +++ b/src/kinetics/Kinetics.cpp @@ -16,7 +16,6 @@ using namespace std; namespace Cantera { Kinetics::Kinetics() : - m_ii(0), m_kk(0), m_thermo(0), m_surfphase(npos), @@ -49,7 +48,6 @@ Kinetics& Kinetics::operator=(const Kinetics& right) m_reactantStoich = right.m_reactantStoich; m_revProductStoich = right.m_revProductStoich; m_irrevProductStoich = right.m_irrevProductStoich; - m_ii = right.m_ii; m_kk = right.m_kk; m_perturb = right.m_perturb; m_reactions = right.m_reactions; @@ -94,15 +92,15 @@ int Kinetics::type() const void Kinetics::checkReactionIndex(size_t i) const { - if (i >= m_ii) { - throw IndexError("checkReactionIndex", "reactions", i, m_ii-1); + if (i >= nReactions()) { + throw IndexError("checkReactionIndex", "reactions", i, nReactions()-1); } } void Kinetics::checkReactionArraySize(size_t ii) const { - if (m_ii > ii) { - throw ArraySizeError("checkReactionArraySize", ii, m_ii); + if (nReactions() > ii) { + throw ArraySizeError("checkReactionArraySize", ii, nReactions()); } } @@ -452,7 +450,7 @@ void Kinetics::getNetRatesOfProgress(doublereal* netROP) void Kinetics::getReactionDelta(const double* prop, double* deltaProp) { - fill(deltaProp, deltaProp + m_ii, 0.0); + fill(deltaProp, deltaProp + nReactions(), 0.0); // products add m_revProductStoich.incrementReactions(prop, deltaProp); m_irrevProductStoich.incrementReactions(prop, deltaProp); @@ -462,7 +460,7 @@ void Kinetics::getReactionDelta(const double* prop, double* deltaProp) void Kinetics::getRevReactionDelta(const double* prop, double* deltaProp) { - fill(deltaProp, deltaProp + m_ii, 0.0); + fill(deltaProp, deltaProp + nReactions(), 0.0); // products add m_revProductStoich.incrementReactions(prop, deltaProp); // reactants subtract @@ -652,7 +650,6 @@ bool Kinetics::addReaction(shared_ptr r) m_irrevProductStoich.add(irxn, pk, pstoich, pstoich); } - incrementRxnCount(); m_reactions.push_back(r); m_rxneqn.push_back(r->equation()); m_reactantStrings.push_back(r->reactantString()); @@ -663,6 +660,7 @@ bool Kinetics::addReaction(shared_ptr r) m_ropf.push_back(0.0); m_ropr.push_back(0.0); m_ropnet.push_back(0.0); + m_perturb.push_back(1.0); return true; }