Cantera  3.1.0a1
Loading...
Searching...
No Matches
InterfaceKinetics.h
Go to the documentation of this file.
1/**
2 * @file InterfaceKinetics.h
3 */
4
5// This file is part of Cantera. See License.txt in the top-level directory or
6// at https://cantera.org/license.txt for license and copyright information.
7
8#ifndef CT_IFACEKINETICS_H
9#define CT_IFACEKINETICS_H
10
11#include "Kinetics.h"
12#include "MultiRate.h"
13
14namespace Cantera
15{
16
17// forward declarations
18class SurfPhase;
19class ImplicitSurfChem;
20
21//! A kinetics manager for heterogeneous reaction mechanisms. The reactions are
22//! assumed to occur at a 2D interface between two 3D phases.
23/*!
24 * There are some important additions to the behavior of the kinetics class due
25 * to the presence of multiple phases and a heterogeneous interface. If a
26 * reactant phase doesn't exists, that is, has a mole number of zero, a
27 * heterogeneous reaction can not proceed from reactants to products. Note it
28 * could perhaps proceed from products to reactants if all of the product phases
29 * exist.
30 *
31 * In order to make the determination of whether a phase exists or not actually
32 * involves the specification of additional information to the kinetics object.,
33 * which heretofore has only had access to intrinsic field information about the
34 * phases (for example, temperature, pressure, and mole fraction).
35 *
36 * The extrinsic specification of whether a phase exists or not must be
37 * specified on top of the intrinsic calculation of the reaction rate. This
38 * class carries a set of booleans indicating whether a phase in the
39 * heterogeneous mechanism exists or not.
40 *
41 * Additionally, the class carries a set of booleans around indicating whether a
42 * product phase is stable or not. If a phase is not thermodynamically stable,
43 * it may be the case that a particular reaction in a heterogeneous mechanism
44 * will create a product species in the unstable phase. However, other reactions
45 * in the mechanism will destruct that species. This may cause oscillations in
46 * the formation of the unstable phase from time step to time step within a ODE
47 * solver, in practice. In order to avoid this situation, a set of booleans is
48 * tracked which sets the stability of a phase. If a phase is deemed to be
49 * unstable, then species in that phase will not be allowed to be birthed by the
50 * kinetics operator. Nonexistent phases are deemed to be unstable by default,
51 * but this can be changed.
52 *
53 * @ingroup kineticsmgr
54 */
56{
57public:
58 //! Constructor
59 InterfaceKinetics() = default;
60
61 ~InterfaceKinetics() override;
62
63 void resizeReactions() override;
64
65 string kineticsType() const override {
66 return "surface";
67 }
68
69 //! Set the electric potential in the nth phase
70 /*!
71 * @param n phase Index in this kinetics object.
72 * @param V Electric potential (volts)
73 */
74 void setElectricPotential(int n, double V);
75
76 //! @name Reaction Rates Of Progress
77 //! @{
78
79 //! Equilibrium constant for all reactions including the voltage term
80 /*!
81 * Kc = exp(deltaG/RT)
82 *
83 * where deltaG is the electrochemical potential difference between
84 * products minus reactants.
85 */
86 void getEquilibriumConstants(double* kc) override;
87
88 void getDeltaGibbs(double* deltaG) override;
89
90 void getDeltaElectrochemPotentials(double* deltaM) override;
91 void getDeltaEnthalpy(double* deltaH) override;
92 void getDeltaEntropy(double* deltaS) override;
93
94 void getDeltaSSGibbs(double* deltaG) override;
95 void getDeltaSSEnthalpy(double* deltaH) override;
96 void getDeltaSSEntropy(double* deltaS) override;
97
98 //! @}
99 //! @name Reaction Mechanism Informational Query Routines
100 //! @{
101
102 void getActivityConcentrations(double* const conc) override;
103
104 bool isReversible(size_t i) override {
105 if (std::find(m_revindex.begin(), m_revindex.end(), i)
106 < m_revindex.end()) {
107 return true;
108 } else {
109 return false;
110 }
111 }
112
113 void getFwdRateConstants(double* kfwd) override;
114 void getRevRateConstants(double* krev, bool doIrreversible=false) override;
115
116 //! @}
117 //! @name Reaction Mechanism Construction
118 //! @{
119
120 //! Add a thermo phase to the kinetics manager object.
121 /*!
122 * This must be done before the function init() is called or
123 * before any reactions are input. The lowest dimensional phase, where reactions
124 * occur, must be added first.
125 *
126 * This function calls Kinetics::addThermo). It also sets the following
127 * fields:
128 *
129 * m_phaseExists[]
130 *
131 * @param thermo Reference to the ThermoPhase to be added.
132 */
133 void addThermo(shared_ptr<ThermoPhase> thermo) override;
134
135 void init() override;
136 void resizeSpecies() override;
137 bool addReaction(shared_ptr<Reaction> r, bool resize=true) override;
138 void modifyReaction(size_t i, shared_ptr<Reaction> rNew) override;
139 void setMultiplier(size_t i, double f) override;
140 //! @}
141
142 //! Internal routine that updates the Rates of Progress of the reactions
143 /*!
144 * This is actually the guts of the functionality of the object
145 */
146 void updateROP() override;
147
148 //! Update properties that depend on temperature
149 /*!
150 * Current objects that this function updates:
151 * m_kdata->m_rfn
152 * updateKc();
153 */
154 void _update_rates_T();
155
156 //! Update properties that depend on the electric potential
157 void _update_rates_phi();
158
159 //! Update properties that depend on the species mole fractions and/or
160 //! concentration,
161 /*!
162 * This method fills out the array of generalized concentrations by calling
163 * method getActivityConcentrations for each phase, which classes
164 * representing phases should overload to return the appropriate quantities.
165 */
166 void _update_rates_C();
167
168 //! Advance the surface coverages in time
169 /*!
170 * This method carries out a time-accurate advancement of the
171 * surface coverages for a specified amount of time.
172 *
173 * @f[
174 * \dot {\theta}_k = \dot s_k (\sigma_k / s_0)
175 * @f]
176 *
177 * @param tstep Time value to advance the surface coverages
178 * @param rtol The relative tolerance for the integrator
179 * @param atol The absolute tolerance for the integrator
180 * @param maxStepSize The maximum step-size the integrator is allowed to take.
181 * If zero, this option is disabled.
182 * @param maxSteps The maximum number of time-steps the integrator can take.
183 * If not supplied, uses the default value in CVodeIntegrator (20000).
184 * @param maxErrTestFails the maximum permissible number of error test failures
185 * If not supplied, uses the default value in CVODES (7).
186 */
187 void advanceCoverages(double tstep, double rtol=1.e-7,
188 double atol=1.e-14, double maxStepSize=0,
189 size_t maxSteps=20000, size_t maxErrTestFails=7);
190
191 //! Solve for the pseudo steady-state of the surface problem
192 /*!
193 * This is the same thing as the advanceCoverages() function,
194 * but at infinite times.
195 *
196 * Note, a direct solve is carried out under the hood here,
197 * to reduce the computational time.
198 *
199 * @param ifuncOverride One of the values defined in @ref solvesp_methods.
200 * The default is -1, which means that the program will decide.
201 * @param timeScaleOverride When a pseudo transient is selected this value
202 * can be used to override the default time scale for
203 * integration which is one. When SFLUX_TRANSIENT is used, this
204 * is equal to the time over which the equations are integrated.
205 * When SFLUX_INITIALIZE is used, this is equal to the time used
206 * in the initial transient algorithm, before the equation
207 * system is solved directly.
208 */
209 void solvePseudoSteadyStateProblem(int ifuncOverride = -1,
210 double timeScaleOverride = 1.0);
211
212 void setIOFlag(int ioFlag);
213
214 //! Update the standard state chemical potentials and species equilibrium
215 //! constant entries
216 /*!
217 * Virtual because it is overridden when dealing with experimental open
218 * circuit voltage overrides
219 */
220 virtual void updateMu0();
221
222 //! Update the equilibrium constants and stored electrochemical potentials
223 //! in molar units for all reversible reactions and for all species.
224 /*!
225 * Irreversible reactions have their equilibrium constant set
226 * to zero. For reactions involving charged species the equilibrium
227 * constant is adjusted according to the electrostatic potential.
228 */
229 void updateKc();
230
231 //! Set the existence of a phase in the reaction object
232 /*!
233 * Tell the kinetics object whether a phase in the object exists. This is
234 * actually an extrinsic specification that must be carried out on top of
235 * the intrinsic calculation of the reaction rate. The routine will also
236 * flip the IsStable boolean within the kinetics object as well.
237 *
238 * @param iphase Index of the phase. This is the order within the
239 * internal thermo vector object
240 * @param exists Boolean indicating whether the phase exists or not
241 */
242 void setPhaseExistence(const size_t iphase, const int exists);
243
244 //! Set the stability of a phase in the reaction object
245 /*!
246 * Tell the kinetics object whether a phase in the object is stable.
247 * Species in an unstable phase will not be allowed to have a positive
248 * rate of formation from this kinetics object. This is actually an
249 * extrinsic specification that must be carried out on top of the
250 * intrinsic calculation of the reaction rate.
251 *
252 * While conceptually not needed since kinetics is consistent with thermo
253 * when taken as a whole, in practice it has found to be very useful to
254 * turn off the creation of phases which shouldn't be forming. Typically
255 * this can reduce the oscillations in phase formation and destruction
256 * which are observed.
257 *
258 * @param iphase Index of the phase. This is the order within the
259 * internal thermo vector object
260 * @param isStable Flag indicating whether the phase is stable or not
261 */
262 void setPhaseStability(const size_t iphase, const int isStable);
263
264 //! Gets the phase existence int for the ith phase
265 /*!
266 * @param iphase Phase Id
267 * @return The int specifying whether the kinetics object thinks the phase
268 * exists or not. If it exists, then species in that phase can be
269 * a reactant in reactions.
270 */
271 int phaseExistence(const size_t iphase) const;
272
273 //! Gets the phase stability int for the ith phase
274 /*!
275 * @param iphase Phase Id
276 * @return The int specifying whether the kinetics object thinks the phase
277 * is stable with nonzero mole numbers. If it stable, then the
278 * kinetics object will allow for rates of production of of
279 * species in that phase that are positive.
280 */
281 int phaseStability(const size_t iphase) const;
282
283 //! Gets the interface current for the ith phase
284 /*!
285 * @param iphase Phase Id
286 * @return The double specifying the interface current. The interface current
287 * is useful when charge transfer reactions occur at an interface. It
288 * is defined here as the net positive charge entering the phase
289 * specified by the Phase Id. (Units: A/m^2 for a surface reaction,
290 * A/m for an edge reaction).
291 */
292 double interfaceCurrent(const size_t iphase);
293
294 void setDerivativeSettings(const AnyMap& settings) override;
295 void getDerivativeSettings(AnyMap& settings) const override;
296 Eigen::SparseMatrix<double> fwdRatesOfProgress_ddCi() override;
297 Eigen::SparseMatrix<double> revRatesOfProgress_ddCi() override;
298 Eigen::SparseMatrix<double> netRatesOfProgress_ddCi() override;
299
300protected:
301 //! @name Internal service methods
302 //!
303 //! @note These methods are for internal use, and seek to avoid code duplication
304 //! while evaluating terms used for rate constants, rates of progress, and
305 //! their derivatives.
306 //! @{
307
308
309 //! Multiply rate with inverse equilibrium constant
310 void applyEquilibriumConstants(double* rop);
311
312 //! Process mole fraction derivative
313 //! @param stoich stoichiometry manager
314 //! @param in rate expression used for the derivative calculation
315 //! @return a sparse matrix of derivative contributions for each reaction of
316 //! dimensions nTotalReactions by nTotalSpecies
317 Eigen::SparseMatrix<double> calculateCompositionDerivatives(StoichManagerN& stoich,
318 const vector<double>& in);
319
320 //! Helper function ensuring that all rate derivatives can be calculated
321 //! @param name method name used for error output
322 //! @throw CanteraError if coverage dependence or electrochemical reactions are
323 //! included
324 void assertDerivativesValid(const string& name);
325
326 //! @}
327
328 //! Temporary work vector of length m_kk
329 vector<double> m_grt;
330
331 //! List of reactions numbers which are reversible reactions
332 /*!
333 * This is a vector of reaction numbers. Each reaction in the list is
334 * reversible. Length = number of reversible reactions
335 */
336 vector<size_t> m_revindex;
337
338 bool m_redo_rates = false;
339
340 //! Vector of rate handlers for interface reactions
341 vector<unique_ptr<MultiRateBase>> m_interfaceRates;
342 map<string, size_t> m_interfaceTypes; //!< Rate handler mapping
343
344 //! Vector of irreversible reaction numbers
345 /*!
346 * vector containing the reaction numbers of irreversible reactions.
347 */
348 vector<size_t> m_irrev;
349
350 //! Array of concentrations for each species in the kinetics mechanism
351 /*!
352 * An array of generalized concentrations @f$ C_k @f$ that are defined
353 * such that @f$ a_k = C_k / C^0_k, @f$ where @f$ C^0_k @f$ is a standard
354 * concentration/ These generalized concentrations are used by this
355 * kinetics manager class to compute the forward and reverse rates of
356 * elementary reactions. The "units" for the concentrations of each phase
357 * depend upon the implementation of kinetics within that phase. The order
358 * of the species within the vector is based on the order of listed
359 * ThermoPhase objects in the class, and the order of the species within
360 * each ThermoPhase class.
361 */
362 vector<double> m_conc;
363
364 //! Array of activity concentrations for each species in the kinetics object
365 /*!
366 * An array of activity concentrations @f$ Ca_k @f$ that are defined
367 * such that @f$ a_k = Ca_k / C^0_k, @f$ where @f$ C^0_k @f$ is a standard
368 * concentration. These activity concentrations are used by this
369 * kinetics manager class to compute the forward and reverse rates of
370 * elementary reactions. The "units" for the concentrations of each phase
371 * depend upon the implementation of kinetics within that phase. The order
372 * of the species within the vector is based on the order of listed
373 * ThermoPhase objects in the class, and the order of the species within
374 * each ThermoPhase class.
375 */
376 vector<double> m_actConc;
377
378 //! Vector of standard state chemical potentials for all species
379 /*!
380 * This vector contains a temporary vector of standard state chemical
381 * potentials for all of the species in the kinetics object
382 *
383 * Length = m_kk. Units = J/kmol.
384 */
385 vector<double> m_mu0;
386
387 //! Vector of chemical potentials for all species
388 /*!
389 * This vector contains a vector of chemical potentials for all of the
390 * species in the kinetics object
391 *
392 * Length = m_kk. Units = J/kmol.
393 */
394 vector<double> m_mu;
395
396 //! Vector of standard state electrochemical potentials modified by a
397 //! standard concentration term.
398 /*!
399 * This vector contains a temporary vector of standard state electrochemical
400 * potentials + RTln(Cs) for all of the species in the kinetics object
401 *
402 * In order to get the units correct for the concentration equilibrium
403 * constant, each species needs to have an RT ln(Cs) added to its
404 * contribution to the equilibrium constant Cs is the standard concentration
405 * for the species. Frequently, for solid species, Cs is equal to 1.
406 * However, for gases Cs is P/RT. Length = m_kk. Units = J/kmol.
407 */
408 vector<double> m_mu0_Kc;
409
410 //! Vector of phase electric potentials
411 /*!
412 * Temporary vector containing the potential of each phase in the kinetics
413 * object. length = number of phases. Units = Volts.
414 */
415 vector<double> m_phi;
416
417 //! Pointer to the single surface phase
418 SurfPhase* m_surf = nullptr;
419
420 //! Pointer to the Implicit surface chemistry object
421 /*!
422 * Note this object is owned by this InterfaceKinetics object. It may only
423 * be used to solve this single InterfaceKinetics object's surface problem
424 * uncoupled from other surface phases.
425 */
427
428 bool m_ROP_ok = false;
429
430 //! Current temperature of the data
431 double m_temp = 0.0;
432
433 //! Int flag to indicate that some phases in the kinetics mechanism are
434 //! non-existent.
435 /*!
436 * We change the ROP vectors to make sure that non-existent phases are
437 * treated correctly in the kinetics operator. The value of this is equal
438 * to the number of phases which don't exist.
439 */
441
442 //! Vector of booleans indicating whether phases exist or not
443 /*!
444 * Vector of booleans indicating whether a phase exists or not. We use this
445 * to set the ROP's so that unphysical things don't happen. For example, a
446 * reaction can't go in the forwards direction if a phase in which a
447 * reactant is present doesn't exist. Because InterfaceKinetics deals with
448 * intrinsic quantities only normally, nowhere else is this extrinsic
449 * concept introduced except here.
450 *
451 * length = number of phases in the object. By default all phases exist.
452 */
453 vector<bool> m_phaseExists;
454
455 //! Vector of int indicating whether phases are stable or not
456 /*!
457 * Vector of booleans indicating whether a phase is stable or not under
458 * the current conditions. We use this to set the ROP's so that
459 * unphysical things don't happen.
460 *
461 * length = number of phases in the object. By default all phases are stable.
462 */
463 vector<int> m_phaseIsStable;
464
465 //! Vector of vector of booleans indicating whether a phase participates in
466 //! a reaction as a reactant
467 /*!
468 * m_rxnPhaseIsReactant[j][p] indicates whether a species in phase p
469 * participates in reaction j as a reactant.
470 */
471 vector<vector<bool>> m_rxnPhaseIsReactant;
472
473 //! Vector of vector of booleans indicating whether a phase participates in a
474 //! reaction as a product
475 /*!
476 * m_rxnPhaseIsReactant[j][p] indicates whether a species in phase p
477 * participates in reaction j as a product.
478 */
479 vector<vector<bool>> m_rxnPhaseIsProduct;
480
481 int m_ioFlag = 0;
482
483 //! Number of dimensions of reacting phase (2 for InterfaceKinetics, 1 for
484 //! EdgeKinetics)
485 size_t m_nDim = 2;
486
487 //! Buffers for partial rop results with length nReactions()
488 vector<double> m_rbuf0;
489 vector<double> m_rbuf1;
490
491 //! A flag used to neglect rate coefficient coverage dependence in derivative
492 //! formation
494 //! A flag used to neglect electrochemical contributions in derivative formation
496 //! Relative tolerance used in developing numerical portions of specific derivatives
497 double m_jac_rtol_delta = 1e-8;
498 //! A flag stating if the object uses electrochemistry
500 //! A flag stating if the object has coverage dependent rates
502};
503
504}
505
506#endif
Base class for kinetics managers and also contains the kineticsmgr module documentation (see Kinetics...
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:427
Advances the surface coverages of the associated set of SurfacePhase objects in time.
A kinetics manager for heterogeneous reaction mechanisms.
bool addReaction(shared_ptr< Reaction > r, bool resize=true) override
Add a single reaction to the mechanism.
int phaseExistence(const size_t iphase) const
Gets the phase existence int for the ith phase.
double interfaceCurrent(const size_t iphase)
Gets the interface current for the ith phase.
SurfPhase * m_surf
Pointer to the single surface phase.
vector< vector< bool > > m_rxnPhaseIsReactant
Vector of vector of booleans indicating whether a phase participates in a reaction as a reactant.
string kineticsType() const override
Identifies the Kinetics manager type.
vector< int > m_phaseIsStable
Vector of int indicating whether phases are stable or not.
Eigen::SparseMatrix< double > revRatesOfProgress_ddCi() override
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
size_t m_nDim
Number of dimensions of reacting phase (2 for InterfaceKinetics, 1 for EdgeKinetics)
void getFwdRateConstants(double *kfwd) override
Return the forward rate constants.
vector< bool > m_phaseExists
Vector of booleans indicating whether phases exist or not.
void getDeltaSSGibbs(double *deltaG) override
Return the vector of values for the reaction standard state Gibbs free energy change.
double m_temp
Current temperature of the data.
void setPhaseStability(const size_t iphase, const int isStable)
Set the stability of a phase in the reaction object.
void getDeltaGibbs(double *deltaG) override
Return the vector of values for the reaction Gibbs free energy change.
bool m_jac_skip_coverage_dependence
A flag used to neglect rate coefficient coverage dependence in derivative formation.
map< string, size_t > m_interfaceTypes
Rate handler mapping.
bool isReversible(size_t i) override
True if reaction i has been declared to be reversible.
void _update_rates_phi()
Update properties that depend on the electric potential.
vector< double > m_phi
Vector of phase electric potentials.
void assertDerivativesValid(const string &name)
Helper function ensuring that all rate derivatives can be calculated.
double m_jac_rtol_delta
Relative tolerance used in developing numerical portions of specific derivatives.
void setMultiplier(size_t i, double f) override
Set the multiplier for reaction i to f.
void getDeltaSSEntropy(double *deltaS) override
Return the vector of values for the change in the standard state entropies for each reaction.
void advanceCoverages(double tstep, double rtol=1.e-7, double atol=1.e-14, double maxStepSize=0, size_t maxSteps=20000, size_t maxErrTestFails=7)
Advance the surface coverages in time.
InterfaceKinetics()=default
Constructor.
void getDeltaSSEnthalpy(double *deltaH) override
Return the vector of values for the change in the standard state enthalpies of reaction.
void resizeSpecies() override
Resize arrays with sizes that depend on the total number of species.
vector< double > m_grt
Temporary work vector of length m_kk.
void getActivityConcentrations(double *const conc) override
Get the vector of activity concentrations used in the kinetics object.
void getRevRateConstants(double *krev, bool doIrreversible=false) override
Return the reverse rate constants.
void getEquilibriumConstants(double *kc) override
Equilibrium constant for all reactions including the voltage term.
vector< size_t > m_revindex
List of reactions numbers which are reversible reactions.
bool m_has_electrochemistry
A flag stating if the object uses electrochemistry.
void updateKc()
Update the equilibrium constants and stored electrochemical potentials in molar units for all reversi...
void setPhaseExistence(const size_t iphase, const int exists)
Set the existence of a phase in the reaction object.
void updateROP() override
Internal routine that updates the Rates of Progress of the reactions.
void applyEquilibriumConstants(double *rop)
Multiply rate with inverse equilibrium constant.
void getDerivativeSettings(AnyMap &settings) const override
Retrieve derivative settings.
void getDeltaEnthalpy(double *deltaH) override
Return the vector of values for the reactions change in enthalpy.
Eigen::SparseMatrix< double > fwdRatesOfProgress_ddCi() override
Calculate derivatives for forward rates-of-progress with respect to species concentration at constant...
void modifyReaction(size_t i, shared_ptr< Reaction > rNew) override
Modify the rate expression associated with a reaction.
bool m_jac_skip_electrochemistry
A flag used to neglect electrochemical contributions in derivative formation.
void _update_rates_C()
Update properties that depend on the species mole fractions and/or concentration,.
vector< double > m_rbuf0
Buffers for partial rop results with length nReactions()
int phaseStability(const size_t iphase) const
Gets the phase stability int for the ith phase.
vector< unique_ptr< MultiRateBase > > m_interfaceRates
Vector of rate handlers for interface reactions.
vector< double > m_mu0_Kc
Vector of standard state electrochemical potentials modified by a standard concentration term.
ImplicitSurfChem * m_integrator
Pointer to the Implicit surface chemistry object.
void _update_rates_T()
Update properties that depend on temperature.
vector< size_t > m_irrev
Vector of irreversible reaction numbers.
void solvePseudoSteadyStateProblem(int ifuncOverride=-1, double timeScaleOverride=1.0)
Solve for the pseudo steady-state of the surface problem.
bool m_has_coverage_dependence
A flag stating if the object has coverage dependent rates.
virtual void updateMu0()
Update the standard state chemical potentials and species equilibrium constant entries.
void addThermo(shared_ptr< ThermoPhase > thermo) override
Add a thermo phase to the kinetics manager object.
vector< double > m_mu
Vector of chemical potentials for all species.
vector< double > m_actConc
Array of activity concentrations for each species in the kinetics object.
vector< double > m_mu0
Vector of standard state chemical potentials for all species.
void getDeltaElectrochemPotentials(double *deltaM) override
Return the vector of values for the reaction electrochemical free energy change.
void setElectricPotential(int n, double V)
Set the electric potential in the nth phase.
void init() override
Prepare the class for the addition of reactions, after all phases have been added.
void resizeReactions() override
Finalize Kinetics object and associated objects.
int m_phaseExistsCheck
Int flag to indicate that some phases in the kinetics mechanism are non-existent.
Eigen::SparseMatrix< double > calculateCompositionDerivatives(StoichManagerN &stoich, const vector< double > &in)
Process mole fraction derivative.
Eigen::SparseMatrix< double > netRatesOfProgress_ddCi() override
Calculate derivatives for net rates-of-progress with respect to species concentration at constant tem...
void getDeltaEntropy(double *deltaS) override
Return the vector of values for the reactions change in entropy.
vector< vector< bool > > m_rxnPhaseIsProduct
Vector of vector of booleans indicating whether a phase participates in a reaction as a product.
void setDerivativeSettings(const AnyMap &settings) override
Set/modify derivative settings.
vector< double > m_conc
Array of concentrations for each species in the kinetics mechanism.
Public interface for kinetics managers.
Definition Kinetics.h:125
ThermoPhase & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition Kinetics.h:242
This class handles operations involving the stoichiometric coefficients on one side of a reaction (re...
A simple thermodynamic model for a surface phase, assuming an ideal solution model.
Definition SurfPhase.h:98
Namespace for the Cantera kernel.
Definition AnyMap.cpp:564