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