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