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