Cantera  2.4.0
solveSP.h
Go to the documentation of this file.
1 /**
2  * @file solveSP.h Header file for implicit surface problem solver (see \ref
3  * chemkinetics and class \link Cantera::solveSP solveSP\endlink).
4  */
5 
6 // This file is part of Cantera. See License.txt in the top-level directory or
7 // at http://www.cantera.org/license.txt for license and copyright information.
8 
9 #ifndef SOLVESP_H
10 #define SOLVESP_H
11 
14 
15 //! @defgroup solvesp_methods Surface Problem Solver Methods
16 //! @{
17 
18 //! This assumes that the initial guess supplied to the routine is far from
19 //! the correct one. Substantial work plus transient time-stepping is to be
20 //! expected to find a solution.
21 const int SFLUX_INITIALIZE = 1;
22 
23 //! Need to solve the surface problem in order to calculate the surface fluxes
24 //! of gas-phase species. (Can expect a moderate change in the solution
25 //! vector; try to solve the system by direct methods with no damping first,
26 //! then try time-stepping if the first method fails). A "time_scale" supplied
27 //! here is used in the algorithm to determine when to shut off time-stepping.
28 const int SFLUX_RESIDUAL = 2;
29 
30 //! Calculation of the surface problem is due to the need for a numerical
31 //! Jacobian for the gas-problem. The solution is expected to be very close to
32 //! the initial guess, and accuracy is needed because solution variables have
33 //! been perturbed from nominal values to create Jacobian entries.
34 const int SFLUX_JACOBIAN = 3;
35 
36 //! The transient calculation is performed here for an amount of time
37 //! specified by "time_scale". It is not guaranteed to be time-accurate -
38 //! just stable and fairly fast. The solution after del_t time is returned,
39 //! whether it's converged to a steady state or not. This is a poor man's time
40 //! stepping algorithm.
41 const int SFLUX_TRANSIENT = 4;
42 // @}
43 
44 //! @defgroup solvesp_bulkFunc Surface Problem Bulk Phase Mode
45 //! Functionality expected from the bulk phase. This changes the equations
46 //! that will be used to solve for the bulk mole fractions.
47 //! @{
48 
49 //! Deposition of a bulk phase is to be expected. Bulk mole fractions are
50 //! determined from ratios of growth rates of bulk species.
51 const int BULK_DEPOSITION = 1;
52 
53 //! Etching of a bulk phase is to be expected. Bulk mole fractions are assumed
54 //! constant, and given by the initial conditions. This is also used whenever
55 //! the condensed phase is part of the larger solution.
56 const int BULK_ETCH = 2;
57 // @}
58 
59 namespace Cantera
60 {
61 
62 //! Method to solve a pseudo steady state surface problem
63 /*!
64  * The following class handles solving the surface problem. The calculation
65  * uses Newton's method to obtain the surface fractions of the surface and
66  * bulk species by requiring that the surface species production rate = 0 and
67  * that the either the bulk fractions are proportional to their production
68  * rates or they are constants.
69  *
70  * Currently, the bulk mole fractions are treated as constants. Implementation
71  * of their being added to the unknown solution vector is delayed.
72  *
73  * Lets introduce the unknown vector for the "surface problem". The surface
74  * problem is defined as the evaluation of the surface site fractions for
75  * multiple surface phases. The unknown vector will consist of the vector of
76  * surface concentrations for each species in each surface vector. Species
77  * are grouped first by their surface phases
78  *
79  * - C_i_j = Concentration of the ith species in the jth surface phase
80  * - Nj = number of surface species in the jth surface phase
81  *
82  * The unknown solution vector is defined as follows:
83  *
84  * C_i_j | kindexSP
85  * --------- | ----------
86  * C_0_0 | 0
87  * C_1_0 | 1
88  * C_2_0 | 2
89  * . . . | ...
90  * C_N0-1_0 | N0-1
91  * C_0_1 | N0
92  * C_1_1 | N0+1
93  * C_2_1 | N0+2
94  * . . . | ...
95  * C_N1-1_1 | NO+N1-1
96  *
97  * Note there are a couple of different types of species indices floating
98  * around in the formulation of this object.
99  *
100  * kindexSP: This is the species index in the contiguous vector of unknowns
101  * for the surface problem.
102  *
103  * Note, in the future, BULK_DEPOSITION systems will be added, and the
104  * solveSP unknown vector will get more complicated. It will include the mole
105  * fraction and growth rates of specified bulk phases
106  *
107  * Indices which relate to individual kinetics objects use the suffix KSI
108  * (kinetics species index).
109  *
110  * ## Solution Method
111  *
112  * This routine is typically used within a residual calculation in a large code.
113  * It's typically invoked millions of times for large calculations, and it must
114  * work every time. Therefore, requirements demand that it be robust but also
115  * efficient.
116  *
117  * The solution methodology is largely determined by the `ifunc` parameter,
118  * that is input to the solution object. This parameter may have one of the
119  * values defined in @ref solvesp_methods.
120  *
121  * ### Pseudo time stepping algorithm:
122  * The time step is determined from sdot[], so that the time step
123  * doesn't ever change the value of a variable by more than 100%.
124  *
125  * This algorithm does use a damped Newton's method to relax the equations.
126  * Damping is based on a "delta damping" technique. The solution unknowns
127  * are not allowed to vary too much between iterations.
128  *
129  * `EXTRA_ACCURACY`: A constant that is the ratio of the required update norm
130  * in this Newton iteration compared to that in the nonlinear solver. A value
131  * of 0.1 is used so surface species are safely overconverged.
132  *
133  * Functions called:
134  * - `ct_dgetrf` -- First half of LAPACK direct solve of a full Matrix
135  * - `ct_dgetrs` -- Second half of LAPACK direct solve of a full matrix.
136  * Returns solution vector in the right-hand-side vector, resid.
137  */
138 class solveSP
139 {
140 public:
141  //! Constructor for the object
142  /*!
143  * @param surfChemPtr Pointer to the ImplicitSurfChem object that
144  * defines the surface problem to be solved.
145  * @param bulkFunc Integer representing how the bulk phases should be
146  * handled. See @ref solvesp_bulkFunc. Currently,
147  * only the default value of BULK_ETCH is supported.
148  */
149  solveSP(ImplicitSurfChem* surfChemPtr, int bulkFunc = BULK_ETCH);
150 
151  //! Destructor. Deletes the integrator.
152  ~solveSP() {}
153 
154 private:
155  //! Unimplemented private copy constructor
156  solveSP(const solveSP& right);
157 
158  //! Unimplemented private assignment operator
159  solveSP& operator=(const solveSP& right);
160 
161 public:
162  //! Main routine that actually calculates the pseudo steady state
163  //! of the surface problem
164  /*!
165  * The actual converged solution is returned as part of the internal state
166  * of the InterfaceKinetics objects.
167  *
168  * Uses Newton's method to get the surface fractions of the surface and
169  * bulk species by requiring that the surface species production rate = 0
170  * and that the bulk fractions are proportional to their production rates.
171  *
172  * @param ifunc Determines the type of solution algorithm to be used. See
173  * @ref solvesp_methods for possible values.
174  * @param time_scale Time over which to integrate the surface equations,
175  * where applicable
176  * @param TKelvin Temperature (kelvin)
177  * @param PGas Pressure (pascals)
178  * @param reltol Relative tolerance to use
179  * @param abstol absolute tolerance.
180  * @return 1 if the surface problem is successfully solved.
181  * -1 if the surface problem wasn't solved successfully.
182  * Note the actual converged solution is returned as part of the
183  * internal state of the InterfaceKinetics objects.
184  */
185  int solveSurfProb(int ifunc, doublereal time_scale, doublereal TKelvin,
186  doublereal PGas, doublereal reltol, doublereal abstol);
187 
188 private:
189  //! Printing routine that optionally gets called at the start of every
190  //! invocation
191  void print_header(int ioflag, int ifunc, doublereal time_scale,
192  int damping, doublereal reltol, doublereal abstol);
193 
194  //! Printing routine that gets called after every iteration
195  void printIteration(int ioflag, doublereal damp, int label_d, int label_t,
196  doublereal inv_t, doublereal t_real, size_t iter,
197  doublereal update_norm, doublereal resid_norm,
198  bool do_time, bool final=false);
199 
200  //! Calculate a conservative delta T to use in a pseudo-steady state
201  //! algorithm
202  /*!
203  * This routine calculates a pretty conservative 1/del_t based on
204  * MAX_i(sdot_i/(X_i*SDen0)). This probably guarantees diagonal dominance.
205  *
206  * Small surface fractions are allowed to intervene in the del_t
207  * determination, no matter how small. This may be changed.
208  * Now minimum changed to 1.0e-12,
209  *
210  * Maximum time step set to time_scale.
211  *
212  * @param netProdRateSolnSP Output variable. Net production rate of all
213  * of the species in the solution vector.
214  * @param XMolSolnSP output variable. Mole fraction of all of the species
215  * in the solution vector
216  * @param label Output variable. Pointer to the value of the species
217  * index (kindexSP) that is controlling the time step
218  * @param label_old Output variable. Pointer to the value of the species
219  * index (kindexSP) that controlled the time step at the previous
220  * iteration
221  * @param label_factor Output variable. Pointer to the current factor
222  * that is used to indicate the same species is controlling the time
223  * step.
224  * @param ioflag Level of the output requested.
225  * @returns the 1. / delta T to be used on the next step
226  */
227  doublereal calc_t(doublereal netProdRateSolnSP[], doublereal XMolSolnSP[],
228  int* label, int* label_old,
229  doublereal* label_factor, int ioflag);
230 
231  //! Calculate the solution and residual weights
232  /*!
233  * @param wtSpecies Weights to use for the soln unknowns. These are in
234  * concentration units
235  * @param wtResid Weights to sue for the residual unknowns.
236  * @param Jac Jacobian. Row sum scaling is used for the Jacobian
237  * @param CSolnSP Solution vector for the surface problem
238  * @param abstol Absolute error tolerance
239  * @param reltol Relative error tolerance
240  */
241  void calcWeights(doublereal wtSpecies[], doublereal wtResid[],
242  const Array2D& Jac, const doublereal CSolnSP[],
243  const doublereal abstol, const doublereal reltol);
244 
245  /**
246  * Update the surface states of the surface phases.
247  */
248  void updateState(const doublereal* cSurfSpec);
249 
250  //! Update mole fraction vector consisting of unknowns in surface problem
251  /*!
252  * @param XMolSolnSP Vector of mole fractions for the unknowns in the
253  * surface problem.
254  */
255  void updateMFSolnSP(doublereal* XMolSolnSP);
256 
257  //! Update the mole fraction vector for a specific kinetic species vector
258  //! corresponding to one InterfaceKinetics object
259  /*!
260  * @param XMolKinSp Mole fraction vector corresponding to a particular
261  * kinetic species for a single InterfaceKinetics Object
262  * This is a vector over all the species in all of the
263  * phases in the InterfaceKinetics object
264  * @param isp ID of the InterfaceKinetics Object.
265  */
266  void updateMFKinSpecies(doublereal* XMolKinSp, int isp);
267 
268  //! Update the vector that keeps track of the largest species in each
269  //! surface phase.
270  /*!
271  * @param CSolnSP Vector of the current values of the surface concentrations
272  * in all of the surface species.
273  */
274  void evalSurfLarge(const doublereal* CSolnSP);
275 
276  //! Main Function evaluation
277  /*!
278  * @param resid output Vector of residuals, length = m_neq
279  * @param CSolnSP Vector of species concentrations, unknowns in the
280  * problem, length = m_neq
281  * @param CSolnOldSP Old Vector of species concentrations, unknowns in the
282  * problem, length = m_neq
283  * @param do_time Calculate a time dependent residual
284  * @param deltaT Delta time for time dependent problem.
285  */
286  void fun_eval(doublereal* resid, const doublereal* CSolnSP,
287  const doublereal* CSolnOldSP, const bool do_time, const doublereal deltaT);
288 
289  //! Main routine that calculates the current residual and Jacobian
290  /*!
291  * @param jac Jacobian to be evaluated.
292  * @param resid output Vector of residuals, length = m_neq
293  * @param CSolnSP Vector of species concentrations, unknowns in the
294  * problem, length = m_neq. These are tweaked in order
295  * to derive the columns of the Jacobian.
296  * @param CSolnSPOld Old Vector of species concentrations, unknowns in the
297  * problem, length = m_neq
298  * @param do_time Calculate a time dependent residual
299  * @param deltaT Delta time for time dependent problem.
300  */
301  void resjac_eval(DenseMatrix& jac, doublereal* resid,
302  doublereal* CSolnSP,
303  const doublereal* CSolnSPOld, const bool do_time,
304  const doublereal deltaT);
305 
306  //! Pointer to the manager of the implicit surface chemistry problem
307  /*!
308  * This object actually calls the current object. Thus, we are providing a
309  * loop-back functionality here.
310  */
312 
313  //! Vector of interface kinetics objects
314  /*!
315  * Each of these is associated with one and only one surface phase.
316  */
317  std::vector<InterfaceKinetics*> &m_objects;
318 
319  //! Total number of equations to solve in the implicit problem.
320  /*!
321  * Note, this can be zero, and frequently is
322  */
323  size_t m_neq;
324 
325  //! This variable determines how the bulk phases are to be handled
326  /*!
327  * Possible values are given in @ref solvesp_bulkFunc.
328  */
330 
331  //! Number of surface phases in the surface problem
332  /*!
333  * This number is equal to the number of InterfaceKinetics objects
334  * in the problem. (until further noted)
335  */
337 
338  //! Total number of surface species in all surface phases.
339  /*!
340  * This is also the number of equations to solve for m_mode=0 system
341  * It's equal to the sum of the number of species in each of the
342  * m_numSurfPhases.
343  */
345 
346  //! Mapping between the surface phases and the InterfaceKinetics objects
347  /*!
348  * Currently this is defined to be a 1-1 mapping (and probably assumed
349  * in some places)
350  * m_surfKinObjID[i] = i
351  */
352  std::vector<size_t> m_indexKinObjSurfPhase;
353 
354  //! Vector of length number of surface phases containing
355  //! the number of surface species in each phase
356  /*!
357  * Length is equal to the number of surface phases, m_numSurfPhases
358  */
359  std::vector<size_t> m_nSpeciesSurfPhase;
360 
361  //! Vector of surface phase pointers
362  /*!
363  * This is created during the constructor
364  * Length is equal to the number of surface phases, m_numSurfPhases
365  */
366  std::vector<SurfPhase*> m_ptrsSurfPhase;
367 
368  //! Index of the start of the unknowns for each solution phase
369  /*!
370  * i_eqn = m_eqnIndexStartPhase[isp]
371  *
372  * isp is the phase id in the list of phases solved by the
373  * surface problem.
374  *
375  * i_eqn is the equation number of the first unknown in the
376  * solution vector corresponding to isp'th phase.
377  */
378  std::vector<size_t> m_eqnIndexStartSolnPhase;
379 
380  //! Phase ID in the InterfaceKinetics object of the surface phase
381  /*!
382  * For each surface phase, this lists the PhaseId of the
383  * surface phase in the corresponding InterfaceKinetics object
384  *
385  * Length is equal to m_numSurfPhases
386  */
387  std::vector<size_t> m_kinObjPhaseIDSurfPhase;
388 
389  //! Total number of volumetric condensed phases included in the steady state
390  //! problem handled by this routine.
391  /*!
392  * This is equal to or less than the total number of volumetric phases in
393  * all of the InterfaceKinetics objects. We usually do not include bulk
394  * phases. Bulk phases are only included in the calculation when their
395  * domain isn't included in the underlying continuum model conservation
396  * equation system.
397  *
398  * This is equal to 0, for the time being
399  */
401 
402  //! Vector of number of species in the m_numBulkPhases phases.
403  /*!
404  * Length is number of bulk phases
405  */
406  std::vector<size_t> m_numBulkSpecies;
407 
408  //! Total number of species in all bulk phases.
409  /*!
410  * This is also the number of bulk equations to solve when bulk equation
411  * solving is turned on.
412  */
414 
415  //! Vector of bulk phase pointers, length is equal to m_numBulkPhases.
416  std::vector<ThermoPhase*> m_bulkPhasePtrs;
417 
418  //! Index between the equation index and the position in the kinetic
419  //! species array for the appropriate kinetics operator
420  /*!
421  * Length = m_neq.
422  *
423  * ksp = m_kinSpecIndex[ieq]
424  * ksp is the kinetic species index for the ieq'th equation.
425  */
426  std::vector<size_t> m_kinSpecIndex;
427 
428  //! Index between the equation index and the index of the
429  //! InterfaceKinetics object
430  /*!
431  * Length m_neq
432  */
433  std::vector<size_t> m_kinObjIndex;
434 
435  //! Vector containing the indices of the largest species
436  //! in each surface phase
437  /*!
438  * `k = m_spSurfLarge[i]` where `k` is the local species index, i.e., it
439  * varies from 0 to (num species in phase - 1) and `i` is the surface
440  * phase index in the problem. Length is equal to #m_numSurfPhases.
441  */
442  std::vector<size_t> m_spSurfLarge;
443 
444  //! The absolute tolerance in real units. units are (kmol/m2)
445  doublereal m_atol;
446 
447  //! The relative error tolerance.
448  doublereal m_rtol;
449 
450  //! maximum value of the time step. units = seconds
451  doublereal m_maxstep;
452 
453  //! Maximum number of species in any single kinetics operator
454  //! -> also maxed wrt the total # of solution species
456 
457  //! Temporary vector with length equal to max m_maxTotSpecies
459 
460  //! Temporary vector with length equal to max m_maxTotSpecies
462 
463  //! Temporary vector with length equal to max m_maxTotSpecies
465 
466  //! Temporary vector with length equal to max m_maxTotSpecies
468 
469  //! Solution vector. length MAX(1, m_neq)
471 
472  //! Saved initial solution vector. length MAX(1, m_neq)
474 
475  //! Saved solution vector at the old time step. length MAX(1, m_neq)
477 
478  //! Weights for the residual norm calculation. length MAX(1, m_neq)
480 
481  //! Weights for the species concentrations norm calculation
482  /*!
483  * length MAX(1, m_neq)
484  */
486 
487  //! Residual for the surface problem
488  /*!
489  * The residual vector of length "dim" that, that has the value of "sdot"
490  * for surface species. The residuals for the bulk species are a function
491  * of the sdots for all species in the bulk phase. The last residual of
492  * each phase enforces {Sum(fractions) = 1}. After linear solve (dgetrf_ &
493  * dgetrs_), resid holds the update vector.
494  *
495  * length MAX(1, m_neq)
496  */
498 
499  //! Vector of mole fractions. length m_maxTotSpecies
501 
502  //! Jacobian. m_neq by m_neq computed Jacobian matrix for the local
503  //! Newton's method.
505 
506 public:
507  int m_ioflag;
508 };
509 }
510 #endif
void printIteration(int ioflag, doublereal damp, int label_d, int label_t, doublereal inv_t, doublereal t_real, size_t iter, doublereal update_norm, doublereal resid_norm, bool do_time, bool final=false)
Printing routine that gets called after every iteration.
Definition: solveSP.cpp:686
int m_bulkFunc
This variable determines how the bulk phases are to be handled.
Definition: solveSP.h:329
std::vector< InterfaceKinetics * > & m_objects
Vector of interface kinetics objects.
Definition: solveSP.h:317
size_t m_neq
Total number of equations to solve in the implicit problem.
Definition: solveSP.h:323
Method to solve a pseudo steady state surface problem.
Definition: solveSP.h:138
solveSP(ImplicitSurfChem *surfChemPtr, int bulkFunc=BULK_ETCH)
Constructor for the object.
Definition: solveSP.cpp:23
vector_fp m_CSolnSPInit
Saved initial solution vector. length MAX(1, m_neq)
Definition: solveSP.h:473
size_t m_numBulkPhasesSS
Total number of volumetric condensed phases included in the steady state problem handled by this rout...
Definition: solveSP.h:400
void updateMFSolnSP(doublereal *XMolSolnSP)
Update mole fraction vector consisting of unknowns in surface problem.
Definition: solveSP.cpp:305
vector_fp m_CSolnSP
Solution vector. length MAX(1, m_neq)
Definition: solveSP.h:470
const int BULK_ETCH
Etching of a bulk phase is to be expected.
Definition: solveSP.h:56
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:31
void evalSurfLarge(const doublereal *CSolnSP)
Update the vector that keeps track of the largest species in each surface phase.
Definition: solveSP.cpp:322
void resjac_eval(DenseMatrix &jac, doublereal *resid, doublereal *CSolnSP, const doublereal *CSolnSPOld, const bool do_time, const doublereal deltaT)
Main routine that calculates the current residual and Jacobian.
Definition: solveSP.cpp:445
size_t m_numTotSurfSpecies
Total number of surface species in all surface phases.
Definition: solveSP.h:344
void print_header(int ioflag, int ifunc, doublereal time_scale, int damping, doublereal reltol, doublereal abstol)
Printing routine that optionally gets called at the start of every invocation.
Definition: solveSP.cpp:637
int solveSurfProb(int ifunc, doublereal time_scale, doublereal TKelvin, doublereal PGas, doublereal reltol, doublereal abstol)
Main routine that actually calculates the pseudo steady state of the surface problem.
Definition: solveSP.cpp:113
const int SFLUX_INITIALIZE
This assumes that the initial guess supplied to the routine is far from the correct one...
Definition: solveSP.h:21
std::vector< size_t > m_eqnIndexStartSolnPhase
Index of the start of the unknowns for each solution phase.
Definition: solveSP.h:378
vector_fp m_wtResid
Weights for the residual norm calculation. length MAX(1, m_neq)
Definition: solveSP.h:479
const int SFLUX_RESIDUAL
Need to solve the surface problem in order to calculate the surface fluxes of gas-phase species...
Definition: solveSP.h:28
vector_fp m_XMolKinSpecies
Vector of mole fractions. length m_maxTotSpecies.
Definition: solveSP.h:500
doublereal calc_t(doublereal netProdRateSolnSP[], doublereal XMolSolnSP[], int *label, int *label_old, doublereal *label_factor, int ioflag)
Calculate a conservative delta T to use in a pseudo-steady state algorithm.
Definition: solveSP.cpp:591
vector_fp m_numEqn1
Temporary vector with length equal to max m_maxTotSpecies.
Definition: solveSP.h:461
const int BULK_DEPOSITION
Deposition of a bulk phase is to be expected.
Definition: solveSP.h:51
const int SFLUX_JACOBIAN
Calculation of the surface problem is due to the need for a numerical Jacobian for the gas-problem...
Definition: solveSP.h:34
std::vector< size_t > m_spSurfLarge
Vector containing the indices of the largest species in each surface phase.
Definition: solveSP.h:442
vector_fp m_CSolnSave
Temporary vector with length equal to max m_maxTotSpecies.
Definition: solveSP.h:467
std::vector< size_t > m_nSpeciesSurfPhase
Vector of length number of surface phases containing the number of surface species in each phase...
Definition: solveSP.h:359
void updateMFKinSpecies(doublereal *XMolKinSp, int isp)
Update the mole fraction vector for a specific kinetic species vector corresponding to one InterfaceK...
Definition: solveSP.cpp:313
solveSP & operator=(const solveSP &right)
Unimplemented private assignment operator.
vector_fp m_numEqn2
Temporary vector with length equal to max m_maxTotSpecies.
Definition: solveSP.h:464
const int SFLUX_TRANSIENT
The transient calculation is performed here for an amount of time specified by "time_scale".
Definition: solveSP.h:41
std::vector< SurfPhase * > m_ptrsSurfPhase
Vector of surface phase pointers.
Definition: solveSP.h:366
doublereal m_maxstep
maximum value of the time step. units = seconds
Definition: solveSP.h:451
std::vector< size_t > m_kinObjPhaseIDSurfPhase
Phase ID in the InterfaceKinetics object of the surface phase.
Definition: solveSP.h:387
vector_fp m_netProductionRatesSave
Temporary vector with length equal to max m_maxTotSpecies.
Definition: solveSP.h:458
std::vector< size_t > m_kinObjIndex
Index between the equation index and the index of the InterfaceKinetics object.
Definition: solveSP.h:433
Advances the surface coverages of the associated set of SurfacePhase objects in time.
size_t m_numSurfPhases
Number of surface phases in the surface problem.
Definition: solveSP.h:336
std::vector< ThermoPhase * > m_bulkPhasePtrs
Vector of bulk phase pointers, length is equal to m_numBulkPhases.
Definition: solveSP.h:416
ImplicitSurfChem * m_SurfChemPtr
Pointer to the manager of the implicit surface chemistry problem.
Definition: solveSP.h:311
~solveSP()
Destructor. Deletes the integrator.
Definition: solveSP.h:152
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:157
doublereal m_atol
The absolute tolerance in real units. units are (kmol/m2)
Definition: solveSP.h:445
std::vector< size_t > m_numBulkSpecies
Vector of number of species in the m_numBulkPhases phases.
Definition: solveSP.h:406
Headers for the DenseMatrix object, which deals with dense rectangular matrices and description of th...
doublereal m_rtol
The relative error tolerance.
Definition: solveSP.h:448
vector_fp m_CSolnSPOld
Saved solution vector at the old time step. length MAX(1, m_neq)
Definition: solveSP.h:476
vector_fp m_resid
Residual for the surface problem.
Definition: solveSP.h:497
std::vector< size_t > m_kinSpecIndex
Index between the equation index and the position in the kinetic species array for the appropriate ki...
Definition: solveSP.h:426
size_t m_maxTotSpecies
Maximum number of species in any single kinetics operator -> also maxed wrt the total # of solution s...
Definition: solveSP.h:455
DenseMatrix m_Jac
Jacobian.
Definition: solveSP.h:504
std::vector< size_t > m_indexKinObjSurfPhase
Mapping between the surface phases and the InterfaceKinetics objects.
Definition: solveSP.h:352
vector_fp m_wtSpecies
Weights for the species concentrations norm calculation.
Definition: solveSP.h:485
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:8
void fun_eval(doublereal *resid, const doublereal *CSolnSP, const doublereal *CSolnOldSP, const bool do_time, const doublereal deltaT)
Main Function evaluation.
Definition: solveSP.cpp:338
size_t m_numTotBulkSpeciesSS
Total number of species in all bulk phases.
Definition: solveSP.h:413
void updateState(const doublereal *cSurfSpec)
Update the surface states of the surface phases.
Definition: solveSP.cpp:296
A class for full (non-sparse) matrices with Fortran-compatible data storage, which adds matrix operat...
Definition: DenseMatrix.h:50
void calcWeights(doublereal wtSpecies[], doublereal wtResid[], const Array2D &Jac, const doublereal CSolnSP[], const doublereal abstol, const doublereal reltol)
Calculate the solution and residual weights.
Definition: solveSP.cpp:558