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