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