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