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