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