Cantera  3.2.0a4
Loading...
Searching...
No Matches
Kinetics.h
Go to the documentation of this file.
1/**
2 * @file Kinetics.h
3 * Base class for kinetics managers and also contains the kineticsmgr
4 * module documentation (see @ref kineticsmgr and class
5 * @link Cantera::Kinetics Kinetics@endlink).
6 */
7
8// This file is part of Cantera. See License.txt in the top-level directory or
9// at https://cantera.org/license.txt for license and copyright information.
10
11#ifndef CT_KINETICS_H
12#define CT_KINETICS_H
13
14#include "StoichManager.h"
16#include "MultiRate.h"
17
18namespace Cantera
19{
20
21class ThermoPhase;
22class Reaction;
23class Solution;
24class AnyMap;
25
26//! @defgroup derivGroup Derivative Calculations
27//! @details Methods for calculating analytical and/or numerical derivatives.
28
29/**
30 * @defgroup chemkinetics Chemical Kinetics
31 */
32
33//! @defgroup reactionGroup Reactions and Reaction Rates
34//! Classes for handling reactions and reaction rates.
35//! @ingroup chemkinetics
36
37//! @defgroup kineticsmgr Kinetics Managers
38//! Classes implementing models for chemical kinetics.
39//! @section kinmodman Models and Managers
40//!
41//! A kinetics manager is a C++ class that implements a kinetics model; a
42//! kinetics model is a set of mathematical equation describing how various
43//! kinetic quantities are to be computed -- reaction rates, species production
44//! rates, etc. Many different kinetics models might be defined to handle
45//! different types of kinetic processes. For example, one kinetics model might
46//! use expressions valid for elementary reactions in ideal gas mixtures. It
47//! might, for example, require the reaction orders to be integral and equal to
48//! the forward stoichiometric coefficients, require that each reaction be
49//! reversible with a reverse rate satisfying detailed balance, include
50//! pressure-dependent unimolecular reactions, etc. Another kinetics model might
51//! be designed for heterogeneous chemistry at interfaces, and might allow
52//! empirical reaction orders, coverage-dependent activation energies,
53//! irreversible reactions, and include effects of potential differences across
54//! the interface on reaction rates.
55//!
56//! A kinetics manager implements a kinetics model. Since the model equations
57//! may be complex and expensive to evaluate, a kinetics manager may adopt
58//! various strategies to 'manage' the computation and evaluate the expressions
59//! efficiently. For example, if there are rate coefficients or other quantities
60//! that depend only on temperature, a manager class may choose to store these
61//! quantities internally, and re-evaluate them only when the temperature has
62//! actually changed. Or a manager designed for use with reaction mechanisms
63//! with a few repeated activation energies might precompute the terms @f$
64//! \exp(-E/RT) @f$, instead of evaluating the exponential repeatedly for each
65//! reaction. There are many other possible 'management styles', each of which
66//! might be better suited to some reaction mechanisms than others.
67//!
68//! But however a manager structures the internal computation, the tasks the
69//! manager class must perform are, for the most part, the same. It must be able
70//! to compute reaction rates, species production rates, equilibrium constants,
71//! etc. Therefore, all kinetics manager classes should have a common set of
72//! public methods, but differ in how they implement these methods.
73//!
74//! A kinetics manager computes reaction rates of progress, species production
75//! rates, equilibrium constants, and similar quantities for a reaction
76//! mechanism. All kinetics manager classes derive from class Kinetics, which
77//! defines a common public interface for all kinetics managers. Each derived
78//! class overloads the virtual methods of Kinetics to implement a particular
79//! kinetics model.
80//!
81//! For example, class BulkKinetics implements reaction rate expressions appropriate for
82//! homogeneous reactions, and class InterfaceKinetics implements expressions
83//! appropriate for heterogeneous mechanisms at interfaces, including how to handle
84//! reactions involving charged species of phases with different electric potentials ---
85//! something that class BulkKinetics doesn't deal with at all.
86//!
87//! Many of the methods of class Kinetics write into arrays the values of some
88//! quantity for each species, for example the net production rate. These
89//! methods always write the results into flat arrays, ordered by phase in the
90//! order the phase was added, and within a phase in the order the species were
91//! added to the phase (which is the same ordering as in the input file).
92//! Example: suppose a heterogeneous mechanism involves three phases -- a bulk
93//! phase 'a', another bulk phase 'b', and the surface phase 'a:b' at the a/b
94//! interface. Phase 'a' contains 12 species, phase 'b' contains 3, and at the
95//! interface there are 5 adsorbed species defined in phase 'a:b'. Then methods
96//! like getNetProductionRates(double* net) will write and output array of
97//! length 20, beginning at the location pointed to by 'net'. The first 12
98//! values will be the net production rates for all 12 species of phase 'a'
99//! (even if some do not participate in the reactions), the next 3 will be for
100//! phase 'b', and finally the net production rates for the surface species will
101//! occupy the last 5 locations.
102//! @ingroup chemkinetics
103
104//! @defgroup rateEvaluators Rate Evaluators
105//! These classes are used to evaluate the rates of reactions.
106//! @ingroup chemkinetics
107
108
109//! Public interface for kinetics managers.
110/*!
111 * This class serves as a base class to derive 'kinetics managers', which are
112 * classes that manage homogeneous chemistry within one phase, or heterogeneous
113 * chemistry at one interface. The virtual methods of this class are meant to be
114 * overloaded in subclasses. The non-virtual methods perform generic functions
115 * and are implemented in Kinetics. They should not be overloaded. Only those
116 * methods required by a subclass need to be overloaded; the rest will throw
117 * exceptions if called.
118 *
119 * When the nomenclature "kinetics species index" is used below, this means that
120 * the species index ranges over all species in all phases handled by the
121 * kinetics manager.
122 *
123 * @ingroup kineticsmgr
124 */
126{
127public:
128 //! @name Constructors and General Information about Mechanism
129 //! @{
130
131 //! Default constructor.
132 Kinetics() = default;
133
134 virtual ~Kinetics() = default;
135
136 //! Kinetics objects are not copyable or assignable
137 Kinetics(const Kinetics&) = delete;
138 Kinetics& operator=(const Kinetics&)= delete;
139
140 //! Create a new Kinetics object with the same kinetics model and reactions as
141 //! this one.
142 //! @param phases Phases used to specify the state for the newly cloned Kinetics
143 //! object. These can be created from the phases used by this object using the
144 //! ThermoPhase::clone() method.
145 //! @since New in %Cantera 3.2.
146 shared_ptr<Kinetics> clone(const vector<shared_ptr<ThermoPhase>>& phases) const;
147
148 //! Identifies the Kinetics manager type.
149 //! Each class derived from Kinetics should override this method to return
150 //! a meaningful identifier.
151 //! @since Starting in %Cantera 3.0, the name returned by this method corresponds
152 //! to the canonical name used in the YAML input format.
153 virtual string kineticsType() const {
154 return "none";
155 }
156
157 //! Finalize Kinetics object and associated objects
158 virtual void resizeReactions();
159
160 //! Number of reactions in the reaction mechanism.
161 size_t nReactions() const {
162 return m_reactions.size();
163 }
164
165 //! Check that the specified reaction index is in range
166 //! Throws an exception if i is greater than nReactions()
167 void checkReactionIndex(size_t m) const;
168
169 //! Check that an array size is at least nReactions()
170 //! Throws an exception if ii is less than nReactions(). Used before calls
171 //! which take an array pointer.
172 void checkReactionArraySize(size_t ii) const;
173
174 //! Check that the specified species index is in range
175 //! Throws an exception if k is greater than nSpecies()-1
176 void checkSpeciesIndex(size_t k) const;
177
178 //! Check that an array size is at least nSpecies()
179 //! Throws an exception if kk is less than nSpecies(). Used before calls
180 //! which take an array pointer.
181 void checkSpeciesArraySize(size_t mm) const;
182
183 //! @}
184 //! @name Information/Lookup Functions about Phases and Species
185 //! @{
186
187 /**
188 * The number of phases participating in the reaction mechanism. For a
189 * homogeneous reaction mechanism, this will always return 1, but for a
190 * heterogeneous mechanism it will return the total number of phases in the
191 * mechanism.
192 */
193 size_t nPhases() const {
194 return m_thermo.size();
195 }
196
197 //! Check that the specified phase index is in range
198 //! Throws an exception if m is greater than nPhases()
199 void checkPhaseIndex(size_t m) const;
200
201 //! Check that an array size is at least nPhases()
202 //! Throws an exception if mm is less than nPhases(). Used before calls
203 //! which take an array pointer.
204 void checkPhaseArraySize(size_t mm) const;
205
206 /**
207 * Return the phase index of a phase in the list of phases defined within
208 * the object.
209 *
210 * @param ph string name of the phase
211 *
212 * If a -1 is returned, then the phase is not defined in the Kinetics
213 * object.
214 */
215 size_t phaseIndex(const string& ph) const {
216 if (m_phaseindex.find(ph) == m_phaseindex.end()) {
217 return npos;
218 } else {
219 return m_phaseindex.at(ph) - 1;
220 }
221 }
222
223 /**
224 * Return pointer to phase where the reactions occur.
225 * @since New in %Cantera 3.0
226 */
227 shared_ptr<ThermoPhase> reactionPhase() const;
228
229 /**
230 * Return pointer to phase associated with Kinetics by index.
231 * @param n Index of the ThermoPhase being sought.
232 * @since New in %Cantera 3.2.
233 * @see thermo
234 */
235 shared_ptr<ThermoPhase> phase(size_t n=0) const {
236 return m_thermo[n];
237 }
238
239 /**
240 * This method returns a reference to the nth ThermoPhase object defined
241 * in this kinetics mechanism. It is typically used so that member
242 * functions of the ThermoPhase object may be called. For homogeneous
243 * mechanisms, there is only one object, and this method can be called
244 * without an argument to access it.
245 *
246 * @param n Index of the ThermoPhase being sought.
247 */
248 ThermoPhase& thermo(size_t n=0) {
249 return *m_thermo[n];
250 }
251 const ThermoPhase& thermo(size_t n=0) const {
252 return *m_thermo[n];
253 }
254
255 /**
256 * The total number of species in all phases participating in the kinetics
257 * mechanism. This is useful to dimension arrays for use in calls to
258 * methods that return the species production rates, for example.
259 */
260 size_t nTotalSpecies() const {
261 return m_kk;
262 }
263
264 /**
265 * The location of species k of phase n in species arrays. Kinetics manager
266 * classes return species production rates in flat arrays, with the species
267 * of each phases following one another, in the order the phases were added.
268 * This method is useful to find the value for a particular species of a
269 * particular phase in arrays returned from methods like getCreationRates
270 * that return an array of species-specific quantities.
271 *
272 * Example: suppose a heterogeneous mechanism involves three phases. The
273 * first contains 12 species, the second 26, and the third 3. Then species
274 * arrays must have size at least 41, and positions 0 - 11 are the values
275 * for the species in the first phase, positions 12 - 37 are the values for
276 * the species in the second phase, etc. Then kineticsSpeciesIndex(7, 0) =
277 * 7, kineticsSpeciesIndex(4, 1) = 16, and kineticsSpeciesIndex(2, 2) = 40.
278 *
279 * @param k species index
280 * @param n phase index for the species
281 */
282 size_t kineticsSpeciesIndex(size_t k, size_t n) const {
283 return m_start[n] + k;
284 }
285
286 //! Return the name of the kth species in the kinetics manager.
287 /*!
288 * k is an integer from 0 to ktot - 1, where ktot is the number of
289 * species in the kinetics manager, which is the sum of the number of
290 * species in all phases participating in the kinetics manager. If k is
291 * out of bounds, the string "<unknown>" is returned.
292 *
293 * @param k species index
294 */
295 string kineticsSpeciesName(size_t k) const;
296
297 /**
298 * This routine will look up a species number based on the input
299 * string nm. The lookup of species will occur for all phases
300 * listed in the kinetics object.
301 *
302 * return
303 * - If a match is found, the position in the species list is returned.
304 * - If no match is found, the value -1 is returned.
305 *
306 * @param nm Input string name of the species
307 */
308 size_t kineticsSpeciesIndex(const string& nm) const;
309
310 /**
311 * This function looks up the name of a species and returns a
312 * reference to the ThermoPhase object of the phase where the species
313 * resides. Will throw an error if the species doesn't match.
314 *
315 * @param nm String containing the name of the species.
316 */
317 ThermoPhase& speciesPhase(const string& nm);
318 const ThermoPhase& speciesPhase(const string& nm) const;
319
320 /**
321 * This function takes as an argument the kineticsSpecies index
322 * (that is, the list index in the list of species in the kinetics
323 * manager) and returns the species' owning ThermoPhase object.
324 *
325 * @param k Species index
326 */
328 return thermo(speciesPhaseIndex(k));
329 }
330
331 /**
332 * This function takes as an argument the kineticsSpecies index (that is, the
333 * list index in the list of species in the kinetics manager) and returns
334 * the index of the phase owning the species.
335 *
336 * @param k Species index
337 */
338 size_t speciesPhaseIndex(size_t k) const;
339
340 //! @}
341 //! @name Reaction Rates Of Progress
342 //! @{
343
344 //! Return the forward rates of progress of the reactions
345 /*!
346 * Forward rates of progress. Return the forward rates of
347 * progress in array fwdROP, which must be dimensioned at
348 * least as large as the total number of reactions.
349 *
350 * @param fwdROP Output vector containing forward rates
351 * of progress of the reactions. Length: nReactions().
352 */
353 virtual void getFwdRatesOfProgress(double* fwdROP);
354
355 //! Return the Reverse rates of progress of the reactions
356 /*!
357 * Return the reverse rates of progress in array revROP, which must be
358 * dimensioned at least as large as the total number of reactions.
359 *
360 * @param revROP Output vector containing reverse rates
361 * of progress of the reactions. Length: nReactions().
362 */
363 virtual void getRevRatesOfProgress(double* revROP);
364
365 /**
366 * Net rates of progress. Return the net (forward - reverse) rates of
367 * progress in array netROP, which must be dimensioned at least as large
368 * as the total number of reactions.
369 *
370 * @param netROP Output vector of the net ROP. Length: nReactions().
371 */
372 virtual void getNetRatesOfProgress(double* netROP);
373
374 //! Return a vector of Equilibrium constants.
375 /*!
376 * Return the equilibrium constants of the reactions in concentration
377 * units in array kc, which must be dimensioned at least as large as the
378 * total number of reactions.
379 *
380 * @f[
381 * Kc_i = \exp [ \Delta G_{ss,i} ] \prod(Cs_k) \exp(\sum_k \nu_{k,i} F \phi_n)
382 * @f]
383 *
384 * @param kc Output vector containing the equilibrium constants.
385 * Length: nReactions().
386 */
387 virtual void getEquilibriumConstants(double* kc) {
388 throw NotImplementedError("Kinetics::getEquilibriumConstants");
389 }
390
391 /**
392 * Change in species properties. Given an array of molar species property
393 * values @f$ z_k, k = 1, \dots, K @f$, return the array of reaction values
394 * @f[
395 * \Delta Z_i = \sum_k \nu_{k,i} z_k, i = 1, \dots, I.
396 * @f]
397 * For example, if this method is called with the array of standard-state
398 * molar Gibbs free energies for the species, then the values returned in
399 * array @c deltaProperty would be the standard-state Gibbs free energies of
400 * reaction for each reaction.
401 *
402 * @param property Input vector of property value. Length: #m_kk.
403 * @param deltaProperty Output vector of deltaRxn. Length: nReactions().
404 */
405 virtual void getReactionDelta(const double* property, double* deltaProperty) const;
406
407 /**
408 * Given an array of species properties 'g', return in array 'dg' the
409 * change in this quantity in the reversible reactions. Array 'g' must
410 * have a length at least as great as the number of species, and array
411 * 'dg' must have a length as great as the total number of reactions.
412 * This method only computes 'dg' for the reversible reactions, and the
413 * entries of 'dg' for the irreversible reactions are unaltered. This is
414 * primarily designed for use in calculating reverse rate coefficients
415 * from thermochemistry for reversible reactions.
416 */
417 virtual void getRevReactionDelta(const double* g, double* dg) const;
418
419 //! Return the vector of values for the reaction Gibbs free energy change.
420 /*!
421 * (virtual from Kinetics.h)
422 * These values depend upon the concentration of the solution.
423 *
424 * units = J kmol-1
425 *
426 * @param deltaG Output vector of deltaG's for reactions Length:
427 * nReactions().
428 */
429 virtual void getDeltaGibbs(double* deltaG) {
430 throw NotImplementedError("Kinetics::getDeltaGibbs");
431 }
432
433 //! Return the vector of values for the reaction electrochemical free
434 //! energy change.
435 /*!
436 * These values depend upon the concentration of the solution and the
437 * voltage of the phases
438 *
439 * units = J kmol-1
440 *
441 * @param deltaM Output vector of deltaM's for reactions Length:
442 * nReactions().
443 */
444 virtual void getDeltaElectrochemPotentials(double* deltaM) {
445 throw NotImplementedError("Kinetics::getDeltaElectrochemPotentials");
446 }
447
448 /**
449 * Return the vector of values for the reactions change in enthalpy.
450 * These values depend upon the concentration of the solution.
451 *
452 * units = J kmol-1
453 *
454 * @param deltaH Output vector of deltaH's for reactions Length:
455 * nReactions().
456 */
457 virtual void getDeltaEnthalpy(double* deltaH) {
458 throw NotImplementedError("Kinetics::getDeltaEnthalpy");
459 }
460
461 /**
462 * Return the vector of values for the reactions change in entropy. These
463 * values depend upon the concentration of the solution.
464 *
465 * units = J kmol-1 Kelvin-1
466 *
467 * @param deltaS Output vector of deltaS's for reactions Length:
468 * nReactions().
469 */
470 virtual void getDeltaEntropy(double* deltaS) {
471 throw NotImplementedError("Kinetics::getDeltaEntropy");
472 }
473
474 /**
475 * Return the vector of values for the reaction standard state Gibbs free
476 * energy change. These values don't depend upon the concentration of the
477 * solution.
478 *
479 * units = J kmol-1
480 *
481 * @param deltaG Output vector of ss deltaG's for reactions Length:
482 * nReactions().
483 */
484 virtual void getDeltaSSGibbs(double* deltaG) {
485 throw NotImplementedError("Kinetics::getDeltaSSGibbs");
486 }
487
488 /**
489 * Return the vector of values for the change in the standard state
490 * enthalpies of reaction. These values don't depend upon the concentration
491 * of the solution.
492 *
493 * units = J kmol-1
494 *
495 * @param deltaH Output vector of ss deltaH's for reactions Length:
496 * nReactions().
497 */
498 virtual void getDeltaSSEnthalpy(double* deltaH) {
499 throw NotImplementedError("Kinetics::getDeltaSSEnthalpy");
500 }
501
502 /**
503 * Return the vector of values for the change in the standard state
504 * entropies for each reaction. These values don't depend upon the
505 * concentration of the solution.
506 *
507 * units = J kmol-1 Kelvin-1
508 *
509 * @param deltaS Output vector of ss deltaS's for reactions Length:
510 * nReactions().
511 */
512 virtual void getDeltaSSEntropy(double* deltaS) {
513 throw NotImplementedError("Kinetics::getDeltaSSEntropy");
514 }
515
516 /**
517 * Return a vector of values of effective concentrations of third-body
518 * collision partners of any reaction. Entries for reactions not involving
519 * third-body collision partners are not defined and set to not-a-number.
520 *
521 * @param concm Output vector of effective third-body concentrations.
522 * Length: nReactions().
523 */
524 virtual void getThirdBodyConcentrations(double* concm) {
525 throw NotImplementedError("Kinetics::getThirdBodyConcentrations",
526 "Not applicable/implemented for Kinetics object of type '{}'",
527 kineticsType());
528 }
529
530 /**
531 * Provide direct access to current third-body concentration values.
532 * @see getThirdBodyConcentrations.
533 */
534 virtual const vector<double>& thirdBodyConcentrations() const {
535 throw NotImplementedError("Kinetics::thirdBodyConcentrations",
536 "Not applicable/implemented for Kinetics object of type '{}'",
537 kineticsType());
538 }
539
540 //! @}
541 //! @name Species Production Rates
542 //! @{
543
544 /**
545 * Species creation rates [kmol/m^3/s or kmol/m^2/s]. Return the species
546 * creation rates in array cdot, which must be dimensioned at least as
547 * large as the total number of species in all phases. @see nTotalSpecies.
548 *
549 * @param cdot Output vector of creation rates. Length: #m_kk.
550 */
551 virtual void getCreationRates(double* cdot);
552
553 /**
554 * Species destruction rates [kmol/m^3/s or kmol/m^2/s]. Return the species
555 * destruction rates in array ddot, which must be dimensioned at least as
556 * large as the total number of species. @see nTotalSpecies.
557 *
558 * @param ddot Output vector of destruction rates. Length: #m_kk.
559 */
560 virtual void getDestructionRates(double* ddot);
561
562 /**
563 * Species net production rates [kmol/m^3/s or kmol/m^2/s]. Return the
564 * species net production rates (creation - destruction) in array wdot,
565 * which must be dimensioned at least as large as the total number of
566 * species. @see nTotalSpecies.
567 *
568 * @param wdot Output vector of net production rates. Length: #m_kk.
569 */
570 virtual void getNetProductionRates(double* wdot);
571
572 //! @}
573
574 //! @addtogroup derivGroup
575 //! @{
576
577 /**
578 * @anchor kinDerivs
579 * @par Routines to Calculate Kinetics Derivatives (Jacobians)
580 * @name
581 *
582 * Kinetics derivatives are calculated with respect to temperature, pressure,
583 * molar concentrations and species mole fractions for forward/reverse/net rates
584 * of progress as well as creation/destruction and net production of species.
585 *
586 * The following suffixes are used to indicate derivatives:
587 * - `_ddT`: derivative with respect to temperature (a vector)
588 * - `_ddP`: derivative with respect to pressure (a vector)
589 * - `_ddC`: derivative with respect to molar concentration (a vector)
590 * - `_ddX`: derivative with respect to species mole fractions (a matrix)
591 * - `_ddCi`: derivative with respect to species concentrations (a matrix)
592 *
593 * @since New in Cantera 2.6
594 *
595 * @warning The calculation of kinetics derivatives is an experimental part of the
596 * %Cantera API and may be changed or removed without notice.
597 *
598 * Source term derivatives are based on a generic rate-of-progress expression
599 * for the @f$ i @f$-th reaction @f$ R_i @f$, which is a function of temperature
600 * @f$ T @f$, pressure @f$ P @f$ and molar concentrations @f$ C_j @f$:
601 * @f[
602 * R_i = k_{f,i} C_M^{\nu_{M,i}} \prod_j C_j^{\nu_{ji}^\prime} -
603 * k_{r,i} C_M^{\nu_{M,i}} \prod_j C_j^{\nu_{ji}^{\prime\prime}}
604 * @f]
605 * Forward/reverse rate expressions @f$ k_{f,i} @f$ and @f$ k_{r,i} @f$ are
606 * implemented by ReactionRate specializations; forward/reverse stoichiometric
607 * coefficients are @f$ \nu_{ji}^\prime @f$ and @f$ \nu_{ji}^{\prime\prime} @f$.
608 * Unless the reaction involves third-body colliders, @f$ \nu_{M,i} = 0 @f$.
609 * For three-body reactions, effective ThirdBody collider concentrations @f$ C_M @f$
610 * are considered with @f$ \nu_{M,i} = 1 @f$. For more detailed information on
611 * relevant theory, see, for example, Perini, et al. @cite perini2012 or Niemeyer,
612 * et al. @cite niemeyer2017, although specifics of %Cantera's implementation may
613 * differ.
614 *
615 * Partial derivatives are obtained from the product rule, where resulting terms
616 * consider reaction rate derivatives, derivatives of the concentration product
617 * term, and, if applicable, third-body term derivatives. ReactionRate
618 * specializations may implement exact derivatives (example:
619 * ArrheniusRate::ddTScaledFromStruct) or approximate them numerically (examples:
620 * ReactionData::perturbTemperature, PlogData::perturbPressure,
621 * FalloffData::perturbThirdBodies). Derivatives of concentration and third-body
622 * terms are based on analytic expressions.
623 *
624 * %Species creation and destruction rates are obtained by multiplying
625 * rate-of-progress vectors by stoichiometric coefficient matrices. As this is a
626 * linear operation, it is possible to calculate derivatives the same way.
627 *
628 * All derivatives are calculated for source terms while holding other properties
629 * constant, independent of whether equation of state or @f$ \sum X_k = 1 @f$
630 * constraints are satisfied. Thus, derivatives deviate from Jacobians and
631 * numerical derivatives that implicitly enforce these constraints. Depending
632 * on application and equation of state, derivatives can nevertheless be used to
633 * obtain Jacobians, for example:
634 *
635 * - The Jacobian of net production rates @f$ \dot{\omega}_{k,\mathrm{net}} @f$
636 * with respect to temperature at constant pressure needs to consider changes
637 * of molar density @f$ C @f$ due to temperature
638 * @f[
639 * \left.
640 * \frac{\partial \dot{\omega}_{k,\mathrm{net}}}{\partial T}
641 * \right|_{P=\mathrm{const}} =
642 * \frac{\partial \dot{\omega}_{k,\mathrm{net}}}{\partial T} +
643 * \frac{\partial \dot{\omega}_{k,\mathrm{net}}}{\partial C}
644 * \left. \frac{\partial C}{\partial T} \right|_{P=\mathrm{const}}
645 * @f]
646 * where for an ideal gas @f$ \partial C / \partial T = - C / T @f$. The
647 * remaining partial derivatives are obtained from getNetProductionRates_ddT()
648 * and getNetProductionRates_ddC(), respectively.
649 *
650 * - The Jacobian of @f$ \dot{\omega}_{k,\mathrm{net}} @f$ with respect to
651 * temperature at constant volume needs to consider pressure changes due to
652 * temperature
653 * @f[
654 * \left.
655 * \frac{\partial \dot{\omega}_{k,\mathrm{net}}}{\partial T}
656 * \right|_{V=\mathrm{const}} =
657 * \frac{\partial \dot{\omega}_{k,\mathrm{net}}}{\partial T} +
658 * \frac{\partial \dot{\omega}_{k,\mathrm{net}}}{\partial P}
659 * \left. \frac{\partial P}{\partial T} \right|_{V=\mathrm{const}}
660 * @f]
661 * where for an ideal gas @f$ \partial P / \partial T = P / T @f$. The
662 * remaining partial derivatives are obtained from getNetProductionRates_ddT()
663 * and getNetProductionRates_ddP(), respectively.
664 *
665 * - Similar expressions can be derived for other derivatives and source terms.
666 *
667 * While some applications require exact derivatives, others can tolerate
668 * approximate derivatives that neglect terms to increase computational speed
669 * and/or improve Jacobian sparsity (example: AdaptivePreconditioner).
670 * Derivative evaluations settings are accessible by keyword/value pairs
671 * using the methods getDerivativeSettings() and setDerivativeSettings().
672 *
673 * For BulkKinetics, the following keyword/value pairs are supported:
674 * - `skip-third-bodies` (boolean): if `false` (default), third body
675 * concentrations are considered for the evaluation of Jacobians
676 * - `skip-falloff` (boolean): if `false` (default), third-body effects
677 * on rate constants are considered for the evaluation of derivatives.
678 * - `rtol-delta` (double): relative tolerance used to perturb properties
679 * when calculating numerical derivatives. The default value is 1e-8.
680 *
681 * For InterfaceKinetics, the following keyword/value pairs are supported:
682 * - `skip-coverage-dependence` (boolean): if `false` (default), rate constant
683 * coverage dependence is not considered when evaluating derivatives.
684 * - `skip-electrochemistry` (boolean): if `false` (default), electrical charge
685 * is not considered in evaluating the derivatives and these reactions are
686 * treated as normal surface reactions.
687 * - `rtol-delta` (double): relative tolerance used to perturb properties
688 * when calculating numerical derivatives. The default value is 1e-8.
689 *
690 * @{
691 */
692
693 /**
694 * Retrieve derivative settings.
695 *
696 * @param settings AnyMap containing settings determining derivative evaluation.
697 */
698 virtual void getDerivativeSettings(AnyMap& settings) const
699 {
700 throw NotImplementedError("Kinetics::getDerivativeSettings",
701 "Not implemented for kinetics type '{}'.", kineticsType());
702 }
703
704 /**
705 * Set/modify derivative settings.
706 *
707 * @param settings AnyMap containing settings determining derivative evaluation.
708 */
709 virtual void setDerivativeSettings(const AnyMap& settings)
710 {
711 throw NotImplementedError("Kinetics::setDerivativeSettings",
712 "Not implemented for kinetics type '{}'.", kineticsType());
713 }
714
715 /**
716 * Calculate derivatives for forward rate constants with respect to temperature
717 * at constant pressure, molar concentration and mole fractions
718 * @param[out] dkfwd Output vector of derivatives. Length: nReactions().
719 */
720 virtual void getFwdRateConstants_ddT(double* dkfwd)
721 {
722 throw NotImplementedError("Kinetics::getFwdRateConstants_ddT",
723 "Not implemented for kinetics type '{}'.", kineticsType());
724 }
725
726 /**
727 * Calculate derivatives for forward rate constants with respect to pressure
728 * at constant temperature, molar concentration and mole fractions.
729 * @param[out] dkfwd Output vector of derivatives. Length: nReactions().
730 */
731 virtual void getFwdRateConstants_ddP(double* dkfwd)
732 {
733 throw NotImplementedError("Kinetics::getFwdRateConstants_ddP",
734 "Not implemented for kinetics type '{}'.", kineticsType());
735 }
736
737 /**
738 * Calculate derivatives for forward rate constants with respect to molar
739 * concentration at constant temperature, pressure and mole fractions.
740 * @param[out] dkfwd Output vector of derivatives. Length: nReactions().
741 *
742 * @warning This method is an experimental part of the %Cantera API and
743 * may be changed or removed without notice.
744 */
745 virtual void getFwdRateConstants_ddC(double* dkfwd)
746 {
747 throw NotImplementedError("Kinetics::getFwdRateConstants_ddC",
748 "Not implemented for kinetics type '{}'.", kineticsType());
749 }
750
751 /**
752 * Calculate derivatives for forward rates-of-progress with respect to temperature
753 * at constant pressure, molar concentration and mole fractions.
754 * @param[out] drop Output vector of derivatives. Length: nReactions().
755 */
756 virtual void getFwdRatesOfProgress_ddT(double* drop)
757 {
758 throw NotImplementedError("Kinetics::getFwdRatesOfProgress_ddT",
759 "Not implemented for kinetics type '{}'.", kineticsType());
760 }
761
762 /**
763 * Calculate derivatives for forward rates-of-progress with respect to pressure
764 * at constant temperature, molar concentration and mole fractions.
765 * @param[out] drop Output vector of derivatives. Length: nReactions().
766 */
767 virtual void getFwdRatesOfProgress_ddP(double* drop)
768 {
769 throw NotImplementedError("Kinetics::getFwdRatesOfProgress_ddP",
770 "Not implemented for kinetics type '{}'.", kineticsType());
771 }
772
773 /**
774 * Calculate derivatives for forward rates-of-progress with respect to molar
775 * concentration at constant temperature, pressure and mole fractions.
776 * @param[out] drop Output vector of derivatives. Length: nReactions().
777 *
778 * @warning This method is an experimental part of the %Cantera API and
779 * may be changed or removed without notice.
780 */
781 virtual void getFwdRatesOfProgress_ddC(double* drop)
782 {
783 throw NotImplementedError("Kinetics::getFwdRatesOfProgress_ddC",
784 "Not implemented for kinetics type '{}'.", kineticsType());
785 }
786
787 /**
788 * Calculate derivatives for forward rates-of-progress with respect to species
789 * mole fractions at constant temperature, pressure and molar concentration.
790 *
791 * The method returns a matrix with nReactions() rows and nTotalSpecies() columns.
792 * For a derivative with respect to @f$ X_i @f$, all other @f$ X_j @f$ are held
793 * constant, rather than enforcing @f$ \sum X_j = 1 @f$.
794 *
795 * @warning This method is an experimental part of the %Cantera API and
796 * may be changed or removed without notice.
797 */
798 virtual Eigen::SparseMatrix<double> fwdRatesOfProgress_ddX()
799 {
800 throw NotImplementedError("Kinetics::fwdRatesOfProgress_ddX",
801 "Not implemented for kinetics type '{}'.", kineticsType());
802 }
803
804 /**
805 * Calculate derivatives for forward rates-of-progress with respect to species
806 * concentration at constant temperature, pressure and remaining species
807 * concentrations.
808 *
809 * The method returns a matrix with nReactions() rows and nTotalSpecies() columns.
810 * For a derivative with respect to @f$ c_i @f$, all other @f$ c_j @f$ are held
811 * constant.
812 *
813 * @warning This method is an experimental part of the %Cantera API and
814 * may be changed or removed without notice.
815 *
816 * @since New in %Cantera 3.0.
817 */
818 virtual Eigen::SparseMatrix<double> fwdRatesOfProgress_ddCi()
819 {
820 throw NotImplementedError("Kinetics::fwdRatesOfProgress_ddCi",
821 "Not implemented for kinetics type '{}'.", kineticsType());
822 }
823
824 /**
825 * Calculate derivatives for reverse rates-of-progress with respect to temperature
826 * at constant pressure, molar concentration and mole fractions.
827 * @param[out] drop Output vector of derivatives. Length: nReactions().
828 */
829 virtual void getRevRatesOfProgress_ddT(double* drop)
830 {
831 throw NotImplementedError("Kinetics::getRevRatesOfProgress_ddT",
832 "Not implemented for kinetics type '{}'.", kineticsType());
833 }
834
835 /**
836 * Calculate derivatives for reverse rates-of-progress with respect to pressure
837 * at constant temperature, molar concentration and mole fractions.
838 * @param[out] drop Output vector of derivatives. Length: nReactions().
839 */
840 virtual void getRevRatesOfProgress_ddP(double* drop)
841 {
842 throw NotImplementedError("Kinetics::getRevRatesOfProgress_ddP",
843 "Not implemented for kinetics type '{}'.", kineticsType());
844 }
845
846 /**
847 * Calculate derivatives for reverse rates-of-progress with respect to molar
848 * concentration at constant temperature, pressure and mole fractions.
849 * @param[out] drop Output vector of derivatives. Length: nReactions().
850 *
851 * @warning This method is an experimental part of the %Cantera API and
852 * may be changed or removed without notice.
853 */
854 virtual void getRevRatesOfProgress_ddC(double* drop)
855 {
856 throw NotImplementedError("Kinetics::getRevRatesOfProgress_ddC",
857 "Not implemented for kinetics type '{}'.", kineticsType());
858 }
859
860 /**
861 * Calculate derivatives for reverse rates-of-progress with respect to species
862 * mole fractions at constant temperature, pressure and molar concentration.
863 *
864 * The method returns a matrix with nReactions() rows and nTotalSpecies() columns.
865 * For a derivative with respect to @f$ X_i @f$, all other @f$ X_j @f$ are held
866 * constant, rather than enforcing @f$ \sum X_j = 1 @f$.
867 *
868 * @warning This method is an experimental part of the %Cantera API and
869 * may be changed or removed without notice.
870 */
871 virtual Eigen::SparseMatrix<double> revRatesOfProgress_ddX()
872 {
873 throw NotImplementedError("Kinetics::revRatesOfProgress_ddX",
874 "Not implemented for kinetics type '{}'.", kineticsType());
875 }
876
877 /**
878 * Calculate derivatives for forward rates-of-progress with respect to species
879 * concentration at constant temperature, pressure and remaining species
880 * concentrations.
881 *
882 * The method returns a matrix with nReactions() rows and nTotalSpecies() columns.
883 * For a derivative with respect to @f$ c_i @f$, all other @f$ c_j @f$ are held
884 * constant.
885 *
886 * @warning This method is an experimental part of the %Cantera API and
887 * may be changed or removed without notice.
888 *
889 * @since New in %Cantera 3.0.
890 */
891 virtual Eigen::SparseMatrix<double> revRatesOfProgress_ddCi()
892 {
893 throw NotImplementedError("Kinetics::revRatesOfProgress_ddCi",
894 "Not implemented for kinetics type '{}'.", kineticsType());
895 }
896
897 /**
898 * Calculate derivatives for net rates-of-progress with respect to temperature
899 * at constant pressure, molar concentration and mole fractions.
900 * @param[out] drop Output vector of derivatives. Length: nReactions().
901 */
902 virtual void getNetRatesOfProgress_ddT(double* drop)
903 {
904 throw NotImplementedError("Kinetics::getNetRatesOfProgress_ddT",
905 "Not implemented for kinetics type '{}'.", kineticsType());
906 }
907
908 /**
909 * Calculate derivatives for net rates-of-progress with respect to pressure
910 * at constant temperature, molar concentration and mole fractions.
911 * @param[out] drop Output vector of derivatives. Length: nReactions().
912 */
913 virtual void getNetRatesOfProgress_ddP(double* drop)
914 {
915 throw NotImplementedError("Kinetics::getNetRatesOfProgress_ddP",
916 "Not implemented for kinetics type '{}'.", kineticsType());
917 }
918
919 /**
920 * Calculate derivatives for net rates-of-progress with respect to molar
921 * concentration at constant temperature, pressure and mole fractions.
922 * @param[out] drop Output vector of derivatives. Length: nReactions().
923 *
924 * @warning This method is an experimental part of the %Cantera API and
925 * may be changed or removed without notice.
926 */
927 virtual void getNetRatesOfProgress_ddC(double* drop)
928 {
929 throw NotImplementedError("Kinetics::getNetRatesOfProgress_ddC",
930 "Not implemented for kinetics type '{}'.", kineticsType());
931 }
932
933 /**
934 * Calculate derivatives for net rates-of-progress with respect to species
935 * mole fractions at constant temperature, pressure and molar concentration.
936 *
937 * The method returns a matrix with nReactions() rows and nTotalSpecies() columns.
938 * For a derivative with respect to @f$ X_i @f$, all other @f$ X_j @f$ are held
939 * constant, rather than enforcing @f$ \sum X_j = 1 @f$.
940 *
941 * @warning This method is an experimental part of the %Cantera API and
942 * may be changed or removed without notice.
943 */
944 virtual Eigen::SparseMatrix<double> netRatesOfProgress_ddX()
945 {
946 throw NotImplementedError("Kinetics::netRatesOfProgress_ddX",
947 "Not implemented for kinetics type '{}'.", kineticsType());
948 }
949
950 /**
951 * Calculate derivatives for net rates-of-progress with respect to species
952 * concentration at constant temperature, pressure, and remaining species
953 * concentrations.
954 *
955 * The method returns a matrix with nReactions() rows and nTotalSpecies() columns.
956 * For a derivative with respect to @f$ c_i @f$, all other @f$ c_j @f$ are held
957 * constant.
958 *
959 * @warning This method is an experimental part of the %Cantera API and
960 * may be changed or removed without notice.
961 *
962 * @since New in %Cantera 3.0.
963 */
964 virtual Eigen::SparseMatrix<double> netRatesOfProgress_ddCi()
965 {
966 throw NotImplementedError("Kinetics::netRatesOfProgress_ddCi",
967 "Not implemented for kinetics type '{}'.", kineticsType());
968 }
969
970 /**
971 * Calculate derivatives for species creation rates with respect to temperature
972 * at constant pressure, molar concentration and mole fractions.
973 * @param[out] dwdot Output vector of derivatives. Length: #m_kk.
974 */
975 void getCreationRates_ddT(double* dwdot);
976
977 /**
978 * Calculate derivatives for species creation rates with respect to pressure
979 * at constant temperature, molar concentration and mole fractions.
980 * @param[out] dwdot Output vector of derivatives. Length: #m_kk.
981 */
982 void getCreationRates_ddP(double* dwdot);
983
984 /**
985 * Calculate derivatives for species creation rates with respect to molar
986 * concentration at constant temperature, pressure and mole fractions.
987 * @param[out] dwdot Output vector of derivatives. Length: #m_kk.
988 *
989 * @warning This method is an experimental part of the %Cantera API and
990 * may be changed or removed without notice.
991 */
992 void getCreationRates_ddC(double* dwdot);
993
994 /**
995 * Calculate derivatives for species creation rates with respect to species
996 * mole fractions at constant temperature, pressure and molar concentration.
997 *
998 * The method returns a square matrix with nTotalSpecies() rows and columns.
999 * For a derivative with respect to @f$ X_i @f$, all other @f$ X_j @f$ are held
1000 * constant, rather than enforcing @f$ \sum X_j = 1 @f$.
1001 *
1002 * @warning This method is an experimental part of the %Cantera API and
1003 * may be changed or removed without notice.
1004 */
1005 Eigen::SparseMatrix<double> creationRates_ddX();
1006
1007 /**
1008 * Calculate derivatives for species creation rates with respect to species
1009 * concentration at constant temperature, pressure, and concentration of all other
1010 * species.
1011 *
1012 * The method returns a square matrix with nTotalSpecies() rows and columns.
1013 * For a derivative with respect to @f$ c_i @f$, all other @f$ c_j @f$ are held
1014 * constant.
1015 *
1016 * @warning This method is an experimental part of the %Cantera API and
1017 * may be changed or removed without notice.
1018 *
1019 * @since New in %Cantera 3.0.
1020 */
1021 Eigen::SparseMatrix<double> creationRates_ddCi();
1022
1023 /**
1024 * Calculate derivatives for species destruction rates with respect to temperature
1025 * at constant pressure, molar concentration and mole fractions.
1026 * @param[out] dwdot Output vector of derivatives. Length: #m_kk.
1027 */
1028 void getDestructionRates_ddT(double* dwdot);
1029
1030 /**
1031 * Calculate derivatives for species destruction rates with respect to pressure
1032 * at constant temperature, molar concentration and mole fractions.
1033 * @param[out] dwdot Output vector of derivatives. Length: #m_kk.
1034 */
1035 void getDestructionRates_ddP(double* dwdot);
1036
1037 /**
1038 * Calculate derivatives for species destruction rates with respect to molar
1039 * concentration at constant temperature, pressure and mole fractions.
1040 * @param[out] dwdot Output vector of derivatives. Length: #m_kk.
1041 *
1042 * @warning This method is an experimental part of the %Cantera API and
1043 * may be changed or removed without notice.
1044 */
1045 void getDestructionRates_ddC(double* dwdot);
1046
1047 /**
1048 * Calculate derivatives for species destruction rates with respect to species
1049 * mole fractions at constant temperature, pressure and molar concentration.
1050 *
1051 * The method returns a square matrix with nTotalSpecies() rows and columns.
1052 * For a derivative with respect to @f$ X_i @f$, all other @f$ X_j @f$ are held
1053 * constant, rather than enforcing @f$ \sum X_j = 1 @f$.
1054 *
1055 * @warning This method is an experimental part of the %Cantera API and
1056 * may be changed or removed without notice.
1057 */
1058 Eigen::SparseMatrix<double> destructionRates_ddX();
1059
1060 /**
1061 * Calculate derivatives for species destruction rates with respect to species
1062 * concentration at constant temperature, pressure, and concentration of all other
1063 * species.
1064 *
1065 * The method returns a square matrix with nTotalSpecies() rows and columns.
1066 * For a derivative with respect to @f$ c_i @f$, all other @f$ c_j @f$ are held
1067 * constant.
1068 *
1069 * @warning This method is an experimental part of the %Cantera API and
1070 * may be changed or removed without notice.
1071 *
1072 * @since New in %Cantera 3.0.
1073 */
1074 Eigen::SparseMatrix<double> destructionRates_ddCi();
1075
1076 /**
1077 * Calculate derivatives for species net production rates with respect to
1078 * temperature at constant pressure, molar concentration and mole fractions.
1079 * @param[out] dwdot Output vector of derivatives. Length: #m_kk.
1080 */
1081 void getNetProductionRates_ddT(double* dwdot);
1082
1083 /**
1084 * Calculate derivatives for species net production rates with respect to pressure
1085 * at constant temperature, molar concentration and mole fractions.
1086 * @param[out] dwdot Output vector of derivatives. Length: #m_kk.
1087 */
1088 void getNetProductionRates_ddP(double* dwdot);
1089
1090 /**
1091 * Calculate derivatives for species net production rates with respect to molar
1092 * concentration at constant temperature, pressure and mole fractions.
1093 * @param[out] dwdot Output vector of derivatives. Length: #m_kk.
1094 *
1095 * @warning This method is an experimental part of the %Cantera API and
1096 * may be changed or removed without notice.
1097 */
1098 void getNetProductionRates_ddC(double* dwdot);
1099
1100 /**
1101 * Calculate derivatives for species net production rates with respect to species
1102 * mole fractions at constant temperature, pressure and molar concentration.
1103 *
1104 * The method returns a square matrix with nTotalSpecies() rows and columns.
1105 * For a derivative with respect to @f$ X_i @f$, all other @f$ X_j @f$ are held
1106 * constant, rather than enforcing @f$ \sum X_j = 1 @f$.
1107 *
1108 * @warning This method is an experimental part of the %Cantera API and
1109 * may be changed or removed without notice.
1110 */
1111 Eigen::SparseMatrix<double> netProductionRates_ddX();
1112
1113 /**
1114 * Calculate derivatives for species net production rates with respect to species
1115 * concentration at constant temperature, pressure, and concentration of all other
1116 * species.
1117 *
1118 * The method returns a square matrix with nTotalSpecies() rows and columns.
1119 * For a derivative with respect to @f$ c_i @f$, all other @f$ c_j @f$ are held
1120 * constant.
1121 *
1122 * @warning This method is an experimental part of the %Cantera API and
1123 * may be changed or removed without notice.
1124 *
1125 * @since New in %Cantera 3.0.
1126 */
1127 Eigen::SparseMatrix<double> netProductionRates_ddCi();
1128
1129 /** @} End of Kinetics Derivatives */
1130 //! @} End of addtogroup derivGroup
1131
1132 //! @name Reaction Mechanism Informational Query Routines
1133 //! @{
1134
1135 /**
1136 * Stoichiometric coefficient of species k as a reactant in reaction i.
1137 *
1138 * @param k kinetic species index
1139 * @param i reaction index
1140 */
1141 virtual double reactantStoichCoeff(size_t k, size_t i) const;
1142
1143 /**
1144 * Stoichiometric coefficient matrix for reactants.
1145 */
1146 Eigen::SparseMatrix<double> reactantStoichCoeffs() const {
1148 }
1149
1150 /**
1151 * Stoichiometric coefficient of species k as a product in reaction i.
1152 *
1153 * @param k kinetic species index
1154 * @param i reaction index
1155 */
1156 virtual double productStoichCoeff(size_t k, size_t i) const;
1157
1158 /**
1159 * Stoichiometric coefficient matrix for products.
1160 */
1161 Eigen::SparseMatrix<double> productStoichCoeffs() const {
1163 }
1164
1165 /**
1166 * Stoichiometric coefficient matrix for products of reversible reactions.
1167 */
1168 Eigen::SparseMatrix<double> revProductStoichCoeffs() const {
1170 }
1171
1172 //! Reactant order of species k in reaction i.
1173 /*!
1174 * This is the nominal order of the activity concentration in
1175 * determining the forward rate of progress of the reaction
1176 *
1177 * @param k kinetic species index
1178 * @param i reaction index
1179 */
1180 virtual double reactantOrder(size_t k, size_t i) const {
1181 throw NotImplementedError("Kinetics::reactantOrder");
1182 }
1183
1184 //! product Order of species k in reaction i.
1185 /*!
1186 * This is the nominal order of the activity concentration of species k in
1187 * determining the reverse rate of progress of the reaction i
1188 *
1189 * For irreversible reactions, this will all be zero.
1190 *
1191 * @param k kinetic species index
1192 * @param i reaction index
1193 */
1194 virtual double productOrder(int k, int i) const {
1195 throw NotImplementedError("Kinetics::productOrder");
1196 }
1197
1198 //! Get the vector of activity concentrations used in the kinetics object
1199 /*!
1200 * @param[out] conc Vector of activity concentrations. Length is equal
1201 * to the number of species in the kinetics object
1202 */
1203 virtual void getActivityConcentrations(double* const conc) {
1204 throw NotImplementedError("Kinetics::getActivityConcentrations");
1205 }
1206
1207 /**
1208 * True if reaction i has been declared to be reversible. If isReversible(i)
1209 * is false, then the reverse rate of progress for reaction i is always
1210 * zero.
1211 *
1212 * @param i reaction index
1213 */
1214 bool isReversible(size_t i) const {
1215 return std::find(m_revindex.begin(), m_revindex.end(), i) < m_revindex.end();
1216 }
1217
1218 /**
1219 * Return the forward rate constants
1220 *
1221 * The computed values include all temperature-dependent and pressure-dependent
1222 * contributions. By default, third-body concentrations are only considered if
1223 * they are part of the reaction rate definition; for a legacy implementation that
1224 * includes third-body concentrations see Cantera::use_legacy_rate_constants().
1225 * Length is the number of reactions. Units are a combination of kmol, m^3 and s,
1226 * that depend on the rate expression for the reaction.
1227 *
1228 * @param kfwd Output vector containing the forward reaction rate
1229 * constants. Length: nReactions().
1230 */
1231 virtual void getFwdRateConstants(double* kfwd) {
1232 throw NotImplementedError("Kinetics::getFwdRateConstants");
1233 }
1234
1235 /**
1236 * Return the reverse rate constants.
1237 *
1238 * The computed values include all temperature-dependent and pressure-dependent
1239 * contributions. By default, third-body concentrations are only considered if
1240 * they are part of the reaction rate definition; for a legacy implementation that
1241 * includes third-body concentrations see Cantera::use_legacy_rate_constants().
1242 * Length is the number of reactions. Units are a combination of kmol, m^3 and s,
1243 * that depend on the rate expression for the reaction. Note, this routine will
1244 * return rate constants for irreversible reactions if the default for
1245 * `doIrreversible` is overridden.
1246 *
1247 * @param krev Output vector of reverse rate constants
1248 * @param doIrreversible boolean indicating whether irreversible reactions
1249 * should be included.
1250 */
1251 virtual void getRevRateConstants(double* krev,
1252 bool doIrreversible = false) {
1253 throw NotImplementedError("Kinetics::getRevRateConstants");
1254 }
1255
1256 //! @}
1257 //! @name Reaction Mechanism Construction
1258 //! @{
1259
1260 //! Add a phase to the kinetics manager object.
1261 /*!
1262 * This must be done before the function init() is called or before any
1263 * reactions are input. The following fields are updated:
1264 *
1265 * - #m_start -> vector of integers, containing the starting position of
1266 * the species for each phase in the kinetics mechanism.
1267 * - #m_thermo -> vector of pointers to ThermoPhase phases that
1268 * participate in the kinetics mechanism.
1269 * - #m_phaseindex -> map containing the string id of each
1270 * ThermoPhase phase as a key and the index of the phase within the
1271 * kinetics manager object as the value.
1272 *
1273 * @param thermo Reference to the ThermoPhase to be added.
1274 * @since New in %Cantera 3.0. Replaces addPhase.
1275 */
1276 virtual void addThermo(shared_ptr<ThermoPhase> thermo);
1277
1278 /**
1279 * Prepare the class for the addition of reactions, after all phases have
1280 * been added. This method is called automatically when the first reaction
1281 * is added. It needs to be called directly only in the degenerate case
1282 * where there are no reactions. The base class method does nothing, but
1283 * derived classes may use this to perform any initialization (allocating
1284 * arrays, etc.) that requires knowing the phases.
1285 */
1286 virtual void init() {}
1287
1288 //! Set kinetics-related parameters from an AnyMap phase description.
1289 //! @since New in %Cantera 3.2.
1290 virtual void setParameters(const AnyMap& phaseNode);
1291
1292 //! Return the parameters for a phase definition which are needed to
1293 //! reconstruct an identical object using the newKinetics function. This
1294 //! excludes the reaction definitions, which are handled separately.
1295 AnyMap parameters() const;
1296
1297 /**
1298 * Resize arrays with sizes that depend on the total number of species.
1299 * Automatically called before adding each Reaction and Phase.
1300 */
1301 virtual void resizeSpecies();
1302
1303 /**
1304 * Add a single reaction to the mechanism. Derived classes should call the
1305 * base class method in addition to handling their own specialized behavior.
1306 *
1307 * @param r Pointer to the Reaction object to be added.
1308 * @param resize If `true`, resizeReactions is called after reaction is added.
1309 * @return `true` if the reaction is added or `false` if it was skipped
1310 */
1311 virtual bool addReaction(shared_ptr<Reaction> r, bool resize=true);
1312
1313 /**
1314 * Modify the rate expression associated with a reaction. The
1315 * stoichiometric equation, type of the reaction, reaction orders, third
1316 * body efficiencies, reversibility, etc. must be unchanged.
1317 *
1318 * @param i Index of the reaction to be modified
1319 * @param rNew Reaction with the new rate expressions
1320 */
1321 virtual void modifyReaction(size_t i, shared_ptr<Reaction> rNew);
1322
1323 /**
1324 * Return the Reaction object for reaction *i*. Changes to this object do
1325 * not affect the Kinetics object until the #modifyReaction function is
1326 * called.
1327 */
1328 shared_ptr<Reaction> reaction(size_t i);
1329
1330 shared_ptr<const Reaction> reaction(size_t i) const;
1331
1332 //! Determine behavior when adding a new reaction that contains species not
1333 //! defined in any of the phases associated with this kinetics manager. If
1334 //! set to true, the reaction will silently be ignored. If false, (the
1335 //! default) an exception will be raised.
1336 void skipUndeclaredSpecies(bool skip) {
1338 }
1339 bool skipUndeclaredSpecies() const {
1341 }
1342
1343 //! Determine behavior when adding a new reaction that contains third-body
1344 //! efficiencies for species not defined in any of the phases associated
1345 //! with this kinetics manager. If set to true, the given third-body
1346 //! efficiency will be ignored. If false, (the default) an exception will be
1347 //! raised.
1350 }
1351 bool skipUndeclaredThirdBodies() const {
1353 }
1354
1355 //! Specify how to handle duplicate third body reactions where one reaction
1356 //! has an explicit third body and the other has the generic third body with a
1357 //! non-zero efficiency for the former third body. Options are "warn" (default),
1358 //! "error", "mark-duplicate", and "modify-efficiency".
1359 void setExplicitThirdBodyDuplicateHandling(const string& flag);
1360 string explicitThirdBodyDuplicateHandling() const {
1361 return m_explicit_third_body_duplicates;
1362 }
1363
1364 //! @}
1365 //! @name Altering Reaction Rates
1366 //!
1367 //! These methods alter reaction rates. They are designed primarily for
1368 //! carrying out sensitivity analysis, but may be used for any purpose
1369 //! requiring dynamic alteration of rate constants. For each reaction, a
1370 //! real-valued multiplier may be defined that multiplies the reaction rate
1371 //! coefficient. The multiplier may be set to zero to completely remove a
1372 //! reaction from the mechanism.
1373 //! @{
1374
1375 //! The current value of the multiplier for reaction i.
1376 /*!
1377 * @param i index of the reaction
1378 */
1379 double multiplier(size_t i) const {
1380 return m_perturb[i];
1381 }
1382
1383 //! Set the multiplier for reaction i to f.
1384 /*!
1385 * @param i index of the reaction
1386 * @param f value of the multiplier.
1387 */
1388 virtual void setMultiplier(size_t i, double f) {
1389 m_perturb[i] = f;
1390 }
1391
1392 virtual void invalidateCache() {
1393 m_cache.clear();
1394 };
1395
1396 //! @}
1397 //! Check for unmarked duplicate reactions and unmatched marked duplicates
1398 //!
1399 //! @param throw_err If `true`, raise an exception that identifies any unmarked
1400 //! duplicate reactions and any reactions marked as duplicate that do not
1401 //! actually have a matching reaction.
1402 //! @param fix If `true` (and if `throw_err` is false), update the `duplicate`
1403 //! flag on all reactions to correctly indicate whether or not they are
1404 //! duplicates.
1405 //! @return If `throw_err` and `fix` are `false`, the indices of the first pair
1406 //! of duplicate reactions or the index of an unmatched duplicate as both
1407 //! elements of the `pair`. Otherwise, `(npos, npos)` if no errors were detected
1408 //! or if the errors were fixed.
1409 //! @since The `fix` argument was added in %Cantera 3.2.
1410 virtual pair<size_t, size_t> checkDuplicates(bool throw_err=true, bool fix=false);
1411
1412 //! Set root Solution holding all phase information
1413 virtual void setRoot(shared_ptr<Solution> root) {
1414 m_root = root;
1415 }
1416
1417 //! Get the Solution object containing this Kinetics object and associated
1418 //! ThermoPhase objects
1419 shared_ptr<Solution> root() const {
1420 return m_root.lock();
1421 }
1422
1423 //! Register a function to be called if reaction is added.
1424 //! @param id A unique ID corresponding to the object affected by the callback.
1425 //! Typically, this is a pointer to an object that also holds a reference to the
1426 //! Kinetics object.
1427 //! @param callback The callback function to be called after any reaction is added.
1428 //! When the callback becomes invalid (for example, the corresponding object is
1429 //! being deleted, the removeReactionAddedCallback() method must be invoked.
1430 //! @since New in %Cantera 3.1
1431 void registerReactionAddedCallback(void* id, const function<void()>& callback)
1432 {
1433 m_reactionAddedCallbacks[id] = callback;
1434 }
1435
1436 //! Remove the reaction-changed callback function associated with the specified object.
1437 //! @since New in %Cantera 3.1
1439 {
1440 m_reactionAddedCallbacks.erase(id);
1441 }
1442
1443protected:
1444 //! Cache for saved calculations within each Kinetics object.
1446
1447 // Update internal rate-of-progress variables #m_ropf and #m_ropr.
1448 virtual void updateROP() {
1449 throw NotImplementedError("Kinetics::updateROP");
1450 }
1451
1452 //! Check whether `r1` and `r2` represent duplicate stoichiometries
1453 //! This function returns a ratio if two reactions are duplicates of
1454 //! one another, and 0.0 otherwise.
1455 /*!
1456 * `r1` and `r2` are maps of species key to stoichiometric coefficient, one
1457 * for each reaction, where the species key is `1+k` for reactants and
1458 * `-1-k` for products and `k` is the species index.
1459 *
1460 * @return 0.0 if the stoichiometries are not multiples of one another
1461 * Otherwise, it returns the ratio of the stoichiometric coefficients.
1462 */
1463 double checkDuplicateStoich(map<int, double>& r1, map<int, double>& r2) const;
1464
1465 //! Vector of rate handlers
1466 vector<unique_ptr<MultiRateBase>> m_rateHandlers;
1467 map<string, size_t> m_rateTypes; //!< Mapping of rate handlers
1468
1469 //! @name Stoichiometry management
1470 //!
1471 //! These objects and functions handle turning reaction extents into species
1472 //! production rates and also handle turning thermo properties into reaction
1473 //! thermo properties.
1474 //! @{
1475
1476 //! Stoichiometry manager for the reactants for each reaction
1478
1479 //! Stoichiometry manager for the products for each reaction
1481
1482 //! Stoichiometry manager for the products of reversible reactions
1484
1485 //! Net stoichiometry (products - reactants)
1486 Eigen::SparseMatrix<double> m_stoichMatrix;
1487 //! @}
1488
1489 //! The number of species in all of the phases
1490 //! that participate in this kinetics mechanism.
1491 size_t m_kk = 0;
1492
1493 //! Vector of perturbation factors for each reaction's rate of
1494 //! progress vector. It is initialized to one.
1495 vector<double> m_perturb;
1496
1497 //! Default values for perturbations defined in the phase definition's
1498 //! `rate-multipliers` field. Key -1 contains the default multiplier.
1499 map<size_t, double> m_defaultPerturb = {{-1, 1.0}};
1500
1501 //! Vector of Reaction objects represented by this Kinetics manager
1502 vector<shared_ptr<Reaction>> m_reactions;
1503
1504 //! m_thermo is a vector of pointers to ThermoPhase objects that are
1505 //! involved with this kinetics operator
1506 /*!
1507 * For homogeneous kinetics applications, this vector will only have one
1508 * entry. For interfacial reactions, this vector will consist of multiple
1509 * entries; some of them will be surface phases, and the other ones will be
1510 * bulk phases. The order that the objects are listed determines the order
1511 * in which the species comprising each phase are listed in the source term
1512 * vector, originating from the reaction mechanism.
1513 */
1514 vector<shared_ptr<ThermoPhase>> m_thermo;
1515
1516 /**
1517 * m_start is a vector of integers specifying the beginning position for the
1518 * species vector for the n'th phase in the kinetics class.
1519 */
1520 vector<size_t> m_start;
1521
1522 /**
1523 * Mapping of the phase name to the position of the phase within the
1524 * kinetics object. Positions start with the value of 1. The member
1525 * function, phaseIndex() decrements by one before returning the index
1526 * value, so that missing phases return -1.
1527 */
1528 map<string, size_t> m_phaseindex;
1529
1530 //! number of spatial dimensions of lowest-dimensional phase.
1531 size_t m_mindim = 4;
1532
1533 //! Forward rate constant for each reaction
1534 vector<double> m_rfn;
1535
1536 //! Delta G^0 for all reactions
1537 vector<double> m_delta_gibbs0;
1538
1539 //! Reciprocal of the equilibrium constant in concentration units
1540 vector<double> m_rkcn;
1541
1542 //! Forward rate-of-progress for each reaction
1543 vector<double> m_ropf;
1544
1545 //! Reverse rate-of-progress for each reaction
1546 vector<double> m_ropr;
1547
1548 //! Net rate-of-progress for each reaction
1549 vector<double> m_ropnet;
1550
1551 vector<size_t> m_revindex; //!< Indices of reversible reactions
1552 vector<size_t> m_irrev; //!< Indices of irreversible reactions
1553
1554 //! The enthalpy change for each reaction to calculate Blowers-Masel rates
1555 vector<double> m_dH;
1556
1557 //! Buffer used for storage of intermediate reaction-specific results
1558 vector<double> m_rbuf;
1559
1560 //! See skipUndeclaredSpecies()
1562
1563 //! See skipUndeclaredThirdBodies()
1565
1566 //! Flag indicating whether reactions include undeclared third bodies
1568
1569
1570 string m_explicit_third_body_duplicates = "warn";
1571
1572 //! reference to Solution
1573 std::weak_ptr<Solution> m_root;
1574
1575 //! Callback functions that are invoked when the reaction is added.
1576 map<void*, function<void()>> m_reactionAddedCallbacks;
1577};
1578
1579}
1580
1581#endif
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
Public interface for kinetics managers.
Definition Kinetics.h:126
shared_ptr< ThermoPhase > phase(size_t n=0) const
Return pointer to phase associated with Kinetics by index.
Definition Kinetics.h:235
virtual void resizeReactions()
Finalize Kinetics object and associated objects.
Definition Kinetics.cpp:49
Eigen::SparseMatrix< double > reactantStoichCoeffs() const
Stoichiometric coefficient matrix for reactants.
Definition Kinetics.h:1146
double multiplier(size_t i) const
The current value of the multiplier for reaction i.
Definition Kinetics.h:1379
void checkPhaseIndex(size_t m) const
Check that the specified phase index is in range Throws an exception if m is greater than nPhases()
Definition Kinetics.cpp:74
void checkSpeciesArraySize(size_t mm) const
Check that an array size is at least nSpecies() Throws an exception if kk is less than nSpecies().
Definition Kinetics.cpp:100
virtual void getFwdRatesOfProgress(double *fwdROP)
Return the forward rates of progress of the reactions.
Definition Kinetics.cpp:367
virtual void setParameters(const AnyMap &phaseNode)
Set kinetics-related parameters from an AnyMap phase description.
Definition Kinetics.cpp:602
virtual void getEquilibriumConstants(double *kc)
Return a vector of Equilibrium constants.
Definition Kinetics.h:387
Kinetics(const Kinetics &)=delete
Kinetics objects are not copyable or assignable.
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition Kinetics.h:248
vector< shared_ptr< Reaction > > m_reactions
Vector of Reaction objects represented by this Kinetics manager.
Definition Kinetics.h:1502
void checkSpeciesIndex(size_t k) const
Check that the specified species index is in range Throws an exception if k is greater than nSpecies(...
Definition Kinetics.cpp:93
double checkDuplicateStoich(map< int, double > &r1, map< int, double > &r2) const
Check whether r1 and r2 represent duplicate stoichiometries This function returns a ratio if two reac...
Definition Kinetics.cpp:260
ValueCache m_cache
Cache for saved calculations within each Kinetics object.
Definition Kinetics.h:1445
vector< double > m_perturb
Vector of perturbation factors for each reaction's rate of progress vector.
Definition Kinetics.h:1495
vector< double > m_ropf
Forward rate-of-progress for each reaction.
Definition Kinetics.h:1543
virtual const vector< double > & thirdBodyConcentrations() const
Provide direct access to current third-body concentration values.
Definition Kinetics.h:534
virtual void getRevReactionDelta(const double *g, double *dg) const
Given an array of species properties 'g', return in array 'dg' the change in this quantity in the rev...
Definition Kinetics.cpp:394
void setExplicitThirdBodyDuplicateHandling(const string &flag)
Specify how to handle duplicate third body reactions where one reaction has an explicit third body an...
Definition Kinetics.cpp:107
virtual void getNetRatesOfProgress(double *netROP)
Net rates of progress.
Definition Kinetics.cpp:379
size_t phaseIndex(const string &ph) const
Return the phase index of a phase in the list of phases defined within the object.
Definition Kinetics.h:215
virtual void getDeltaSSEntropy(double *deltaS)
Return the vector of values for the change in the standard state entropies for each reaction.
Definition Kinetics.h:512
vector< unique_ptr< MultiRateBase > > m_rateHandlers
Vector of rate handlers.
Definition Kinetics.h:1466
vector< size_t > m_start
m_start is a vector of integers specifying the beginning position for the species vector for the n'th...
Definition Kinetics.h:1520
virtual string kineticsType() const
Identifies the Kinetics manager type.
Definition Kinetics.h:153
virtual void getDeltaElectrochemPotentials(double *deltaM)
Return the vector of values for the reaction electrochemical free energy change.
Definition Kinetics.h:444
vector< double > m_rkcn
Reciprocal of the equilibrium constant in concentration units.
Definition Kinetics.h:1540
shared_ptr< ThermoPhase > reactionPhase() const
Return pointer to phase where the reactions occur.
Definition Kinetics.cpp:88
size_t m_kk
The number of species in all of the phases that participate in this kinetics mechanism.
Definition Kinetics.h:1491
virtual pair< size_t, size_t > checkDuplicates(bool throw_err=true, bool fix=false)
Check for unmarked duplicate reactions and unmatched marked duplicates.
Definition Kinetics.cpp:119
virtual void getReactionDelta(const double *property, double *deltaProperty) const
Change in species properties.
Definition Kinetics.cpp:385
void checkPhaseArraySize(size_t mm) const
Check that an array size is at least nPhases() Throws an exception if mm is less than nPhases().
Definition Kinetics.cpp:81
vector< double > m_dH
The enthalpy change for each reaction to calculate Blowers-Masel rates.
Definition Kinetics.h:1555
vector< double > m_ropr
Reverse rate-of-progress for each reaction.
Definition Kinetics.h:1546
size_t nPhases() const
The number of phases participating in the reaction mechanism.
Definition Kinetics.h:193
vector< shared_ptr< ThermoPhase > > m_thermo
m_thermo is a vector of pointers to ThermoPhase objects that are involved with this kinetics operator
Definition Kinetics.h:1514
map< void *, function< void()> > m_reactionAddedCallbacks
Callback functions that are invoked when the reaction is added.
Definition Kinetics.h:1576
AnyMap parameters() const
Return the parameters for a phase definition which are needed to reconstruct an identical object usin...
Definition Kinetics.cpp:619
virtual bool addReaction(shared_ptr< Reaction > r, bool resize=true)
Add a single reaction to the mechanism.
Definition Kinetics.cpp:675
virtual void getThirdBodyConcentrations(double *concm)
Return a vector of values of effective concentrations of third-body collision partners of any reactio...
Definition Kinetics.h:524
vector< size_t > m_revindex
Indices of reversible reactions.
Definition Kinetics.h:1551
void skipUndeclaredSpecies(bool skip)
Determine behavior when adding a new reaction that contains species not defined in any of the phases ...
Definition Kinetics.h:1336
string kineticsSpeciesName(size_t k) const
Return the name of the kth species in the kinetics manager.
Definition Kinetics.cpp:303
virtual void getDestructionRates(double *ddot)
Species destruction rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:417
bool isReversible(size_t i) const
True if reaction i has been declared to be reversible.
Definition Kinetics.h:1214
Kinetics()=default
Default constructor.
virtual void getDeltaEntropy(double *deltaS)
Return the vector of values for the reactions change in entropy.
Definition Kinetics.h:470
virtual void init()
Prepare the class for the addition of reactions, after all phases have been added.
Definition Kinetics.h:1286
virtual void getFwdRateConstants(double *kfwd)
Return the forward rate constants.
Definition Kinetics.h:1231
shared_ptr< Solution > root() const
Get the Solution object containing this Kinetics object and associated ThermoPhase objects.
Definition Kinetics.h:1419
bool m_skipUndeclaredThirdBodies
See skipUndeclaredThirdBodies()
Definition Kinetics.h:1564
virtual double reactantOrder(size_t k, size_t i) const
Reactant order of species k in reaction i.
Definition Kinetics.h:1180
virtual void modifyReaction(size_t i, shared_ptr< Reaction > rNew)
Modify the rate expression associated with a reaction.
Definition Kinetics.cpp:768
map< string, size_t > m_rateTypes
Mapping of rate handlers.
Definition Kinetics.h:1467
vector< double > m_ropnet
Net rate-of-progress for each reaction.
Definition Kinetics.h:1549
void skipUndeclaredThirdBodies(bool skip)
Determine behavior when adding a new reaction that contains third-body efficiencies for species not d...
Definition Kinetics.h:1348
virtual double productStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a product in reaction i.
Definition Kinetics.cpp:362
virtual void addThermo(shared_ptr< ThermoPhase > thermo)
Add a phase to the kinetics manager object.
Definition Kinetics.cpp:585
Eigen::SparseMatrix< double > productStoichCoeffs() const
Stoichiometric coefficient matrix for products.
Definition Kinetics.h:1161
vector< double > m_rbuf
Buffer used for storage of intermediate reaction-specific results.
Definition Kinetics.h:1558
Eigen::SparseMatrix< double > m_stoichMatrix
Net stoichiometry (products - reactants)
Definition Kinetics.h:1486
virtual void getActivityConcentrations(double *const conc)
Get the vector of activity concentrations used in the kinetics object.
Definition Kinetics.h:1203
bool m_skipUndeclaredSpecies
See skipUndeclaredSpecies()
Definition Kinetics.h:1561
map< string, size_t > m_phaseindex
Mapping of the phase name to the position of the phase within the kinetics object.
Definition Kinetics.h:1528
size_t m_mindim
number of spatial dimensions of lowest-dimensional phase.
Definition Kinetics.h:1531
StoichManagerN m_productStoich
Stoichiometry manager for the products for each reaction.
Definition Kinetics.h:1480
virtual void getDeltaSSGibbs(double *deltaG)
Return the vector of values for the reaction standard state Gibbs free energy change.
Definition Kinetics.h:484
std::weak_ptr< Solution > m_root
reference to Solution
Definition Kinetics.h:1573
virtual void getDeltaSSEnthalpy(double *deltaH)
Return the vector of values for the change in the standard state enthalpies of reaction.
Definition Kinetics.h:498
StoichManagerN m_revProductStoich
Stoichiometry manager for the products of reversible reactions.
Definition Kinetics.h:1483
shared_ptr< Kinetics > clone(const vector< shared_ptr< ThermoPhase > > &phases) const
Create a new Kinetics object with the same kinetics model and reactions as this one.
Definition Kinetics.cpp:27
Eigen::SparseMatrix< double > revProductStoichCoeffs() const
Stoichiometric coefficient matrix for products of reversible reactions.
Definition Kinetics.h:1168
size_t nReactions() const
Number of reactions in the reaction mechanism.
Definition Kinetics.h:161
map< size_t, double > m_defaultPerturb
Default values for perturbations defined in the phase definition's rate-multipliers field.
Definition Kinetics.h:1499
void removeReactionAddedCallback(void *id)
Remove the reaction-changed callback function associated with the specified object.
Definition Kinetics.h:1438
virtual void getRevRatesOfProgress(double *revROP)
Return the Reverse rates of progress of the reactions.
Definition Kinetics.cpp:373
vector< size_t > m_irrev
Indices of irreversible reactions.
Definition Kinetics.h:1552
virtual void setRoot(shared_ptr< Solution > root)
Set root Solution holding all phase information.
Definition Kinetics.h:1413
void registerReactionAddedCallback(void *id, const function< void()> &callback)
Register a function to be called if reaction is added.
Definition Kinetics.h:1431
bool m_hasUndeclaredThirdBodies
Flag indicating whether reactions include undeclared third bodies.
Definition Kinetics.h:1567
size_t kineticsSpeciesIndex(size_t k, size_t n) const
The location of species k of phase n in species arrays.
Definition Kinetics.h:282
virtual double productOrder(int k, int i) const
product Order of species k in reaction i.
Definition Kinetics.h:1194
virtual void setMultiplier(size_t i, double f)
Set the multiplier for reaction i to f.
Definition Kinetics.h:1388
vector< double > m_rfn
Forward rate constant for each reaction.
Definition Kinetics.h:1534
virtual void getDeltaEnthalpy(double *deltaH)
Return the vector of values for the reactions change in enthalpy.
Definition Kinetics.h:457
StoichManagerN m_reactantStoich
Stoichiometry manager for the reactants for each reaction.
Definition Kinetics.h:1477
virtual void getRevRateConstants(double *krev, bool doIrreversible=false)
Return the reverse rate constants.
Definition Kinetics.h:1251
virtual void resizeSpecies()
Resize arrays with sizes that depend on the total number of species.
Definition Kinetics.cpp:663
void checkReactionArraySize(size_t ii) const
Check that an array size is at least nReactions() Throws an exception if ii is less than nReactions()...
Definition Kinetics.cpp:66
size_t nTotalSpecies() const
The total number of species in all phases participating in the kinetics mechanism.
Definition Kinetics.h:260
virtual double reactantStoichCoeff(size_t k, size_t i) const
Stoichiometric coefficient of species k as a reactant in reaction i.
Definition Kinetics.cpp:357
void checkReactionIndex(size_t m) const
Check that the specified reaction index is in range Throws an exception if i is greater than nReactio...
Definition Kinetics.cpp:42
virtual void getDeltaGibbs(double *deltaG)
Return the vector of values for the reaction Gibbs free energy change.
Definition Kinetics.h:429
virtual void getNetProductionRates(double *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:428
ThermoPhase & speciesPhase(const string &nm)
This function looks up the name of a species and returns a reference to the ThermoPhase object of the...
Definition Kinetics.cpp:325
vector< double > m_delta_gibbs0
Delta G^0 for all reactions.
Definition Kinetics.h:1537
virtual void getCreationRates(double *cdot)
Species creation rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:403
ThermoPhase & speciesPhase(size_t k)
This function takes as an argument the kineticsSpecies index (that is, the list index in the list of ...
Definition Kinetics.h:327
size_t speciesPhaseIndex(size_t k) const
This function takes as an argument the kineticsSpecies index (that is, the list index in the list of ...
Definition Kinetics.cpp:346
An error indicating that an unimplemented function has been called.
This class handles operations involving the stoichiometric coefficients on one side of a reaction (re...
const Eigen::SparseMatrix< double > & stoichCoeffs() const
Return matrix containing stoichiometric coefficients.
Base class for a phase with thermodynamic properties.
Storage for cached values.
Definition ValueCache.h:153
void clear()
Clear all cached values.
virtual void setDerivativeSettings(const AnyMap &settings)
Set/modify derivative settings.
Definition Kinetics.h:709
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi()
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
Definition Kinetics.h:818
Eigen::SparseMatrix< double > creationRates_ddCi()
Calculate derivatives for species creation rates with respect to species concentration at constant te...
Definition Kinetics.cpp:485
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddCi()
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
Definition Kinetics.h:964
void getCreationRates_ddT(double *dwdot)
Calculate derivatives for species creation rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:439
virtual Eigen::SparseMatrix< double > netRatesOfProgress_ddX()
Calculate derivatives for net rates-of-progress with respect to species mole fractions at constant te...
Definition Kinetics.h:944
Eigen::SparseMatrix< double > destructionRates_ddX()
Calculate derivatives for species destruction rates with respect to species mole fractions at constan...
Definition Kinetics.cpp:531
void getCreationRates_ddC(double *dwdot)
Calculate derivatives for species creation rates with respect to molar concentration at constant temp...
Definition Kinetics.cpp:463
void getDestructionRates_ddP(double *dwdot)
Calculate derivatives for species destruction rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:507
void getNetProductionRates_ddC(double *dwdot)
Calculate derivatives for species net production rates with respect to molar concentration at constan...
Definition Kinetics.cpp:567
void getDestructionRates_ddT(double *dwdot)
Calculate derivatives for species destruction rates with respect to temperature at constant pressure,...
Definition Kinetics.cpp:495
virtual void getFwdRatesOfProgress_ddP(double *drop)
Calculate derivatives for forward rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:767
Eigen::SparseMatrix< double > netProductionRates_ddX()
Calculate derivatives for species net production rates with respect to species mole fractions at cons...
Definition Kinetics.cpp:575
Eigen::SparseMatrix< double > creationRates_ddX()
Calculate derivatives for species creation rates with respect to species mole fractions at constant t...
Definition Kinetics.cpp:475
Eigen::SparseMatrix< double > destructionRates_ddCi()
Calculate derivatives for species destruction rates with respect to species concentration at constant...
Definition Kinetics.cpp:541
void getNetProductionRates_ddT(double *dwdot)
Calculate derivatives for species net production rates with respect to temperature at constant pressu...
Definition Kinetics.cpp:551
virtual Eigen::SparseMatrix< double > fwdRatesOfProgress_ddX()
Calculate derivatives for forward rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:798
virtual void getRevRatesOfProgress_ddT(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:829
virtual void getFwdRateConstants_ddP(double *dkfwd)
Calculate derivatives for forward rate constants with respect to pressure at constant temperature,...
Definition Kinetics.h:731
virtual void getFwdRateConstants_ddT(double *dkfwd)
Calculate derivatives for forward rate constants with respect to temperature at constant pressure,...
Definition Kinetics.h:720
void getCreationRates_ddP(double *dwdot)
Calculate derivatives for species creation rates with respect to pressure at constant temperature,...
Definition Kinetics.cpp:451
virtual void getNetRatesOfProgress_ddT(double *drop)
Calculate derivatives for net rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:902
virtual void getRevRatesOfProgress_ddP(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:840
Eigen::SparseMatrix< double > netProductionRates_ddCi()
Calculate derivatives for species net production rates with respect to species concentration at const...
Definition Kinetics.cpp:580
virtual void getRevRatesOfProgress_ddC(double *drop)
Calculate derivatives for reverse rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:854
virtual void getDerivativeSettings(AnyMap &settings) const
Retrieve derivative settings.
Definition Kinetics.h:698
void getNetProductionRates_ddP(double *dwdot)
Calculate derivatives for species net production rates with respect to pressure at constant temperatu...
Definition Kinetics.cpp:559
virtual void getNetRatesOfProgress_ddP(double *drop)
Calculate derivatives for net rates-of-progress with respect to pressure at constant temperature,...
Definition Kinetics.h:913
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddCi()
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
Definition Kinetics.h:891
virtual void getFwdRatesOfProgress_ddT(double *drop)
Calculate derivatives for forward rates-of-progress with respect to temperature at constant pressure,...
Definition Kinetics.h:756
void getDestructionRates_ddC(double *dwdot)
Calculate derivatives for species destruction rates with respect to molar concentration at constant t...
Definition Kinetics.cpp:519
virtual void getNetRatesOfProgress_ddC(double *drop)
Calculate derivatives for net rates-of-progress with respect to molar concentration at constant tempe...
Definition Kinetics.h:927
virtual Eigen::SparseMatrix< double > revRatesOfProgress_ddX()
Calculate derivatives for reverse rates-of-progress with respect to species mole fractions at constan...
Definition Kinetics.h:871
virtual void getFwdRateConstants_ddC(double *dkfwd)
Calculate derivatives for forward rate constants with respect to molar concentration at constant temp...
Definition Kinetics.h:745
virtual void getFwdRatesOfProgress_ddC(double *drop)
Calculate derivatives for forward rates-of-progress with respect to molar concentration at constant t...
Definition Kinetics.h:781
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180