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