Cantera  2.1.2
solveProb.cpp
Go to the documentation of this file.
1 /**
2  * @file: solveProb.cpp Implicit solver for nonlinear problems
3  */
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 
13 #include "cantera/base/clockWC.h"
16 
17 using namespace std;
18 namespace Cantera
19 {
20 
21 /***************************************************************************
22  * STATIC ROUTINES DEFINED IN THIS FILE
23  ***************************************************************************/
24 
25 static doublereal calcWeightedNorm(const doublereal [], const doublereal dx[], size_t);
26 
27 solveProb::solveProb(ResidEval* resid) :
28  m_residFunc(resid),
29  m_neq(0),
30  m_atol(0),
31  m_rtol(1.0E-4),
32  m_maxstep(1000),
33  m_ioflag(0)
34 {
36 
37  // Dimension solution vector
38  size_t dim1 = std::max<size_t>(1, m_neq);
39 
40  m_atol.resize(dim1, 1.0E-9);
41  m_netProductionRatesSave.resize(dim1, 0.0);
42  m_numEqn1.resize(dim1, 0.0);
43  m_numEqn2.resize(dim1, 0.0);
44  m_CSolnSave.resize(dim1, 0.0);
45  m_CSolnSP.resize(dim1, 0.0);
46  m_CSolnSPInit.resize(dim1, 0.0);
47  m_CSolnSPOld.resize(dim1, 0.0);
48  m_wtResid.resize(dim1, 0.0);
49  m_wtSpecies.resize(dim1, 0.0);
50  m_resid.resize(dim1, 0.0);
51  m_ipiv.resize(dim1, 0);
52  m_topBounds.resize(dim1, 1.0);
53  m_botBounds.resize(dim1, 0.0);
54 
55  m_Jac.resize(dim1, dim1, 0.0);
56  m_JacCol.resize(dim1, 0);
57  for (size_t k = 0; k < dim1; k++) {
58  m_JacCol[k] = m_Jac.ptrColumn(k);
59  }
60 
61 }
62 
63 solveProb::~solveProb()
64 {
65 }
66 
67 int solveProb::solve(int ifunc, doublereal time_scale,
68  doublereal reltol)
69 {
70  /*
71  * The following calculation is a Newton's method to get the surface fractions
72  * of the surface and bulk species by requiring that the surface species
73  * production rate = 0 and that the bulk fractions are proportional to their
74  * production rates.
75  */
76  doublereal EXTRA_ACCURACY = 0.001;
77  if (ifunc == SOLVEPROB_JACOBIAN) {
78  EXTRA_ACCURACY *= 0.001;
79  }
80  int info = 0;
81  size_t label_t = npos; /* Species IDs for time control */
82  size_t label_d; /* Species IDs for damping control */
83  size_t label_t_old = npos;
84  doublereal label_factor = 1.0;
85  int iter=0; // iteration number on numlinear solver
86  int iter_max=1000; // maximum number of nonlinear iterations
87  int nrhs=1;
88  doublereal deltaT = 1.0E-10; // Delta time step
89  doublereal damp=1.0, tmp;
90  // Weighted L2 norm of the residual. Currently, this is only
91  // used for IO purposes. It doesn't control convergence.
92  // Therefore, it is turned off when DEBUG_SOLVEPROB isn't defined.
93  doublereal resid_norm;
94  doublereal inv_t = 0.0;
95  doublereal t_real = 0.0, update_norm = 1.0E6;
96 
97  bool do_time = false, not_converged = true;
98 
99 #ifdef DEBUG_SOLVEPROB
100 #ifdef DEBUG_SOLVEPROB_TIME
101  doublereal t1;
102 #endif
103 #else
104  if (m_ioflag > 1) {
105  m_ioflag = 1;
106  }
107 #endif
108 
109 #ifdef DEBUG_SOLVEPROB
110 #ifdef DEBUG_SOLVEPROB_TIME
111  Cantera::clockWC wc;
112  if (m_ioflag) {
113  t1 = wc.secondsWC();
114  }
115 #endif
116 #endif
117 
118  /*
119  * Set the initial value of the do_time parameter
120  */
121  if (ifunc == SOLVEPROB_INITIALIZE || ifunc == SOLVEPROB_TRANSIENT) {
122  do_time = true;
123  }
124 
125  /*
126  * upload the initial conditions
127  */
129 
130  /*
131  * Store the initial guess in the soln vector,
132  * CSolnSP, and in an separate vector CSolnSPInit.
133  */
134  std::copy(m_CSolnSP.begin(), m_CSolnSP.end(), m_CSolnSPInit.begin());
135 
136 
137 
138  if (m_ioflag) {
139  print_header(m_ioflag, ifunc, time_scale, reltol,
141  }
142 
143  /*
144  * Quick return when there isn't a surface problem to solve
145  */
146  if (m_neq == 0) {
147  not_converged = false;
148  update_norm = 0.0;
149  }
150 
151  /* ------------------------------------------------------------------
152  * Start of Newton's method
153  * ------------------------------------------------------------------
154  */
155  while (not_converged && iter < iter_max) {
156  iter++;
157  /*
158  * Store previous iteration's solution in the old solution vector
159  */
160  std::copy(m_CSolnSP.begin(), m_CSolnSP.end(), m_CSolnSPOld.begin());
161 
162  /*
163  * Evaluate the largest surface species for each surface phase every
164  * 5 iterations.
165  */
166  // if (iter%5 == 4) {
167  // evalSurfLarge(DATA_PTR(m_CSolnSP));
168  // }
169 
170  /*
171  * Calculate the value of the time step
172  * - heuristics to stop large oscillations in deltaT
173  */
174  if (do_time) {
175  /* don't hurry increase in time step at the same time as damping */
176  if (damp < 1.0) {
177  label_factor = 1.0;
178  }
180  &label_t, &label_t_old, &label_factor, m_ioflag);
181  if (iter < 10) {
182  inv_t = tmp;
183  } else if (tmp > 2.0*inv_t) {
184  inv_t = 2.0*inv_t;
185  } else {
186  inv_t = tmp;
187  }
188 
189  /*
190  * Check end condition
191  */
192 
193  if (ifunc == SOLVEPROB_TRANSIENT) {
194  tmp = t_real + 1.0/inv_t;
195  if (tmp > time_scale) {
196  inv_t = 1.0/(time_scale - t_real);
197  }
198  }
199  } else {
200  /* make steady state calc a step of 1 million seconds to
201  prevent singular jacobians for some pathological cases */
202  inv_t = 1.0e-6;
203  }
204  deltaT = 1.0/inv_t;
205 
206  /*
207  * Call the routine to numerically evaluation the jacobian
208  * and residual for the current iteration.
209  */
211  DATA_PTR(m_CSolnSPOld), do_time, deltaT);
212 
213  /*
214  * Calculate the weights. Make sure the calculation is carried
215  * out on the first iteration.
216  */
217  if (iter%4 == 1) {
220  }
221 
222  /*
223  * Find the weighted norm of the residual
224  */
225  resid_norm = calcWeightedNorm(DATA_PTR(m_wtResid), DATA_PTR(m_resid), m_neq);
226 
227 #ifdef DEBUG_SOLVEPROB
228  if (m_ioflag > 1) {
229  printIterationHeader(m_ioflag, damp, inv_t, t_real, iter, do_time);
230  /*
231  * Print out the residual and jacobian
232  */
233  printResJac(m_ioflag, m_neq, m_Jac, DATA_PTR(m_resid),
234  DATA_PTR(m_wtResid), resid_norm);
235  }
236 #endif
237 
238  /*
239  * Solve Linear system (with LAPACK). The solution is in resid[]
240  */
241 
242  ct_dgetrf(m_neq, m_neq, m_JacCol[0], m_neq, DATA_PTR(m_ipiv), info);
243  if (info==0) {
244  ct_dgetrs(ctlapack::NoTranspose, m_neq, nrhs, m_JacCol[0],
246  info);
247  }
248  /*
249  * Force convergence if residual is small to avoid
250  * "nan" results from the linear solve.
251  */
252  else {
253  if (m_ioflag) {
254  printf("solveSurfSS: Zero pivot, assuming converged: %g (%d)\n",
255  resid_norm, info);
256  }
257  for (size_t jcol = 0; jcol < m_neq; jcol++) {
258  m_resid[jcol] = 0.0;
259  }
260 
261  /* print out some helpful info */
262  if (m_ioflag > 1) {
263  printf("-----\n");
264  printf("solveSurfProb: iter %d t_real %g delta_t %g\n\n",
265  iter,t_real, 1.0/inv_t);
266  printf("solveSurfProb: init guess, current concentration,"
267  "and prod rate:\n");
268 
269  printf("-----\n");
270  }
271  if (do_time) {
272  t_real += time_scale;
273  }
274 #ifdef DEBUG_SOLVEPROB
275  if (m_ioflag) {
276  printf("\nResidual is small, forcing convergence!\n");
277  }
278 #endif
279  }
280 
281  /*
282  * Calculate the Damping factor needed to keep all unknowns
283  * between 0 and 1, and not allow too large a change (factor of 2)
284  * in any unknown.
285  */
286 
287 
288  damp = calc_damping(DATA_PTR(m_CSolnSP), DATA_PTR(m_resid), m_neq, &label_d);
289 
290 
291  /*
292  * Calculate the weighted norm of the update vector
293  * Here, resid is the delta of the solution, in concentration
294  * units.
295  */
296  update_norm = calcWeightedNorm(DATA_PTR(m_wtSpecies),
298  /*
299  * Update the solution vector and real time
300  * Crop the concentrations to zero.
301  */
302  for (size_t irow = 0; irow < m_neq; irow++) {
303  m_CSolnSP[irow] -= damp * m_resid[irow];
304  }
305 
306 
307  if (do_time) {
308  t_real += damp/inv_t;
309  }
310 
311  if (m_ioflag) {
312  printIteration(m_ioflag, damp, label_d, label_t, inv_t, t_real, iter,
313  update_norm, resid_norm,
316  DATA_PTR(m_wtSpecies), m_neq, do_time);
317  }
318 
319  if (ifunc == SOLVEPROB_TRANSIENT) {
320  not_converged = (t_real < time_scale);
321  } else {
322  if (do_time) {
323  if (t_real > time_scale ||
324  (resid_norm < 1.0e-7 &&
325  update_norm*time_scale/t_real < EXTRA_ACCURACY)) {
326  do_time = false;
327 #ifdef DEBUG_SOLVEPROB
328  if (m_ioflag > 1) {
329  printf("\t\tSwitching to steady solve.\n");
330  }
331 #endif
332  }
333  } else {
334  not_converged = ((update_norm > EXTRA_ACCURACY) ||
335  (resid_norm > EXTRA_ACCURACY));
336  }
337  }
338  } /* End of Newton's Method while statement */
339 
340  /*
341  * End Newton's method. If not converged, print error message and
342  * recalculate sdot's at equal site fractions.
343  */
344  if (not_converged) {
345  if (m_ioflag) {
346  printf("#$#$#$# Error in solveProb $#$#$#$ \n");
347  printf("Newton iter on surface species did not converge, "
348  "update_norm = %e \n", update_norm);
349  printf("Continuing anyway\n");
350  }
351  }
352 #ifdef DEBUG_SOLVEPROB
353 #ifdef DEBUG_SOLVEPROB_TIME
354  if (m_ioflag) {
355  printf("\nEnd of solve, time used: %e\n", wc.secondsWC()-t1);
356  }
357 #endif
358 #endif
359 
360  /*
361  * Decide on what to return in the solution vector
362  * - right now, will always return the last solution
363  * no matter how bad
364  */
365  if (m_ioflag) {
367  false, deltaT);
368  resid_norm = calcWeightedNorm(DATA_PTR(m_wtResid),
370  printFinal(m_ioflag, damp, label_d, label_t, inv_t, t_real, iter,
371  update_norm, resid_norm, DATA_PTR(m_netProductionRatesSave),
374  DATA_PTR(m_wtResid), m_neq, do_time);
375  }
376 
377  /*
378  * Return with the appropriate flag
379  */
380  if (update_norm > 1.0) {
381  return -1;
382  }
383  return 0;
384 }
385 
386 void solveProb::reportState(doublereal* const CSolnSP) const
387 {
388  std::copy(m_CSolnSP.begin(), m_CSolnSP.end(), CSolnSP);
389 }
390 
391 void solveProb::fun_eval(doublereal* const resid, const doublereal* const CSoln,
392  const doublereal* const CSolnOld, const bool do_time,
393  const doublereal deltaT)
394 {
395  /*
396  * This routine uses the m_numEqn1 and m_netProductionRatesSave vectors
397  * as temporary internal storage.
398  */
399  if (do_time) {
400  m_residFunc->evalSimpleTD(0.0, CSoln, CSolnOld, deltaT, resid);
401  } else {
402  m_residFunc->evalSS(0.0, CSoln, resid);
403  }
404 }
405 
406 void solveProb::resjac_eval(std::vector<doublereal*> &JacCol,
407  doublereal resid[], doublereal CSoln[],
408  const doublereal CSolnOld[], const bool do_time,
409  const doublereal deltaT)
410 {
411  doublereal dc, cSave, sd;
412  doublereal* col_j;
413  /*
414  * Calculate the residual
415  */
416  fun_eval(resid, CSoln, CSolnOld, do_time, deltaT);
417  /*
418  * Now we will look over the columns perturbing each unknown.
419  */
420 
421  for (size_t kCol = 0; kCol < m_neq; kCol++) {
422  cSave = CSoln[kCol];
423  sd = fabs(cSave) + fabs(CSoln[kCol]) + m_atol[kCol] * 1.0E6;
424  if (sd < 1.0E-200) {
425  sd = 1.0E-4;
426  }
427  dc = std::max(1.0E-11 * sd, fabs(cSave) * 1.0E-6);
428  CSoln[kCol] += dc;
429  // Use the m_numEqn2 vector as temporary internal storage.
430  fun_eval(DATA_PTR(m_numEqn2), CSoln, CSolnOld, do_time, deltaT);
431  col_j = JacCol[kCol];
432  for (size_t i = 0; i < m_neq; i++) {
433  col_j[i] = (m_numEqn2[i] - resid[i])/dc;
434  }
435  CSoln[kCol] = cSave;
436  }
437 
438 }
439 
440 doublereal solveProb::calc_damping(doublereal x[], doublereal dxneg[], size_t dim, size_t* label)
441 {
442  const doublereal APPROACH = 0.50;
443  doublereal damp = 1.0, xnew, xtop, xbot;
444  static doublereal damp_old = 1.0;
445  *label = npos;
446 
447  for (size_t i = 0; i < dim; i++) {
448  doublereal topBounds = m_topBounds[i];
449  doublereal botBounds = m_botBounds[i];
450  /*
451  * Calculate the new suggested new value of x[i]
452  */
453  double delta_x = - dxneg[i];
454  xnew = x[i] - damp * dxneg[i];
455 
456  /*
457  * Calculate the allowed maximum and minimum values of x[i]
458  * - Only going to allow x[i] to converge to the top and bottom bounds by a
459  * single order of magnitude at one time
460  */
461  bool canCrossOrigin = false;
462  if (topBounds > 0.0 && botBounds < 0.0) {
463  canCrossOrigin = true;
464  }
465 
466  xtop = topBounds - 0.1 * fabs(topBounds - x[i]);
467 
468  xbot = botBounds + 0.1 * fabs(x[i] - botBounds);
469 
470  if (xnew > xtop) {
471  damp = - APPROACH * (xtop - x[i]) / dxneg[i];
472  *label = i;
473  } else if (xnew < xbot) {
474  damp = APPROACH * (x[i] - xbot) / dxneg[i];
475  *label = i;
476  }
477  // else if (fabs(xnew) > 2.0*MAX(fabs(x[i]), 1.0E-10)) {
478  // damp = 0.5 * MAX(fabs(x[i]), 1.0E-9)/ fabs(xnew);
479  // *label = i;
480  // }
481  double denom = fabs(x[i]) + 1.0E5 * m_atol[i];
482  if ((fabs(delta_x) / denom) > 0.3) {
483  double newdamp = 0.3 * denom / fabs(delta_x);
484  if (canCrossOrigin) {
485  if (xnew * x[i] < 0.0) {
486  if (fabs(x[i]) < 1.0E8 * m_atol[i]) {
487  newdamp = 2.0 * fabs(x[i]) / fabs(delta_x);
488  }
489  }
490  }
491  damp = std::min(damp, newdamp);
492  }
493 
494  }
495 
496  /*
497  * Only allow the damping parameter to increase by a factor of three each
498  * iteration. Heuristic to avoid oscillations in the value of damp
499  */
500  if (damp > damp_old*3) {
501  damp = damp_old*3;
502  *label = npos;
503  }
504 
505  /*
506  * Save old value of the damping parameter for use
507  * in subsequent calls.
508  */
509  damp_old = damp;
510  return damp;
511 
512 }
513 
514 /*
515  * This function calculates the norm of an update, dx[],
516  * based on the weighted values of x.
517  */
518 static doublereal calcWeightedNorm(const doublereal wtX[], const doublereal dx[], size_t dim)
519 {
520  doublereal norm = 0.0;
521  doublereal tmp;
522  if (dim == 0) {
523  return 0.0;
524  }
525  for (size_t i = 0; i < dim; i++) {
526  tmp = dx[i] / wtX[i];
527  norm += tmp * tmp;
528  }
529  return sqrt(norm/dim);
530 }
531 
532 void solveProb::calcWeights(doublereal wtSpecies[], doublereal wtResid[],
533  const doublereal CSoln[])
534 {
535  /*
536  * First calculate the weighting factor
537  */
538 
539  for (size_t k = 0; k < m_neq; k++) {
540  wtSpecies[k] = m_atol[k] + m_rtol * fabs(CSoln[k]);
541  }
542  /*
543  * Now do the residual Weights. Since we have the Jacobian, we
544  * will use it to generate a number based on the what a significant
545  * change in a solution variable does to each residual.
546  * This is a row sum scale operation.
547  */
548  for (size_t k = 0; k < m_neq; k++) {
549  wtResid[k] = 0.0;
550  for (size_t jcol = 0; jcol < m_neq; jcol++) {
551  wtResid[k] += fabs(m_Jac(k,jcol) * wtSpecies[jcol]);
552  }
553  }
554 }
555 
556 doublereal solveProb::
557 calc_t(doublereal netProdRateSolnSP[], doublereal Csoln[],
558  size_t* label, size_t* label_old, doublereal* label_factor, int ioflag)
559 {
560  doublereal tmp, inv_timeScale=0.0;
561  for (size_t k = 0; k < m_neq; k++) {
562  if (Csoln[k] <= 1.0E-10) {
563  tmp = 1.0E-10;
564  } else {
565  tmp = Csoln[k];
566  }
567  tmp = fabs(netProdRateSolnSP[k]/ tmp);
568 
569 
570  if (netProdRateSolnSP[k]> 0.0) {
571  tmp /= 100.;
572  }
573  if (tmp > inv_timeScale) {
574  inv_timeScale = tmp;
575  *label = k;
576  }
577  }
578 
579  /*
580  * Increase time step exponentially as same species repeatedly
581  * controls time step
582  */
583  if (*label == *label_old) {
584  *label_factor *= 1.5;
585  } else {
586  *label_old = *label;
587  *label_factor = 1.0;
588  }
589  inv_timeScale = inv_timeScale / *label_factor;
590 #ifdef DEBUG_SOLVEPROB
591  if (ioflag > 1) {
592  if (*label_factor > 1.0) {
593  printf("Delta_t increase due to repeated controlling species = %e\n",
594  *label_factor);
595  }
596  int kkin = m_kinSpecIndex[*label];
597 
598  string sn = " "
599  printf("calc_t: spec=%d(%s) sf=%e pr=%e dt=%e\n",
600  *label, sn.c_str(), XMolSolnSP[*label],
601  netProdRateSolnSP[*label], 1.0/inv_timeScale);
602  }
603 #endif
604 
605  return inv_timeScale;
606 
607 }
608 
609 void solveProb::setBounds(const doublereal botBounds[], const doublereal topBounds[])
610 {
611  for (size_t k = 0; k < m_neq; k++) {
612  m_botBounds[k] = botBounds[k];
613  m_topBounds[k] = topBounds[k];
614  }
615 }
616 
617 #ifdef DEBUG_SOLVEPROB
618 void solveProb::printResJac(int ioflag, int neq, const Array2D& Jac,
619  doublereal resid[], doublereal wtRes[],
620  doublereal norm)
621 {
622 
623 }
624 #endif
625 
626 void solveProb::print_header(int ioflag, int ifunc, doublereal time_scale,
627  doublereal reltol,
628  doublereal netProdRate[])
629 {
630  int damping = 1;
631  if (ioflag) {
632  printf("\n================================ SOLVEPROB CALL SETUP "
633  "========================================\n");
634  if (ifunc == SOLVEPROB_INITIALIZE) {
635  printf("\n SOLVEPROB Called with Initialization turned on\n");
636  printf(" Time scale input = %9.3e\n", time_scale);
637  } else if (ifunc == SOLVEPROB_RESIDUAL) {
638  printf("\n SOLVEPROB Called to calculate steady state residual\n");
639  printf(" from a good initial guess\n");
640  } else if (ifunc == SOLVEPROB_JACOBIAN) {
641  printf("\n SOLVEPROB Called to calculate steady state jacobian\n");
642  printf(" from a good initial guess\n");
643  } else if (ifunc == SOLVEPROB_TRANSIENT) {
644  printf("\n SOLVEPROB Called to integrate surface in time\n");
645  printf(" for a total of %9.3e sec\n", time_scale);
646  } else {
647  fprintf(stderr,"Unknown ifunc flag = %d\n", ifunc);
648  exit(EXIT_FAILURE);
649  }
650 
651 
652 
653  if (damping) {
654  printf(" Damping is ON \n");
655  } else {
656  printf(" Damping is OFF \n");
657  }
658 
659  printf(" Reltol = %9.3e, Abstol = %9.3e\n", reltol, m_atol[0]);
660  }
661 
662  /*
663  * Print out the initial guess
664  */
665 #ifdef DEBUG_SOLVEPROB
666  if (ioflag > 1) {
667  printf("\n================================ INITIAL GUESS "
668  "========================================\n");
669  int kindexSP = 0;
670  for (int isp = 0; isp < m_numSurfPhases; isp++) {
671  InterfaceKinetics* m_kin = m_objects[isp];
672  int surfIndex = m_kin->surfacePhaseIndex();
673  int nPhases = m_kin->nPhases();
674  m_kin->getNetProductionRates(netProdRate);
675  updateMFKinSpecies(XMolKinSpecies, isp);
676 
677  printf("\n IntefaceKinetics Object # %d\n\n", isp);
678 
679  printf("\t Number of Phases = %d\n", nPhases);
680  printf("\t Phase:SpecName Prod_Rate MoleFraction kindexSP\n");
681  printf("\t -------------------------------------------------------"
682  "----------\n");
683 
684  int kspindex = 0;
685  bool inSurfacePhase = false;
686  for (int ip = 0; ip < nPhases; ip++) {
687  if (ip == surfIndex) {
688  inSurfacePhase = true;
689  } else {
690  inSurfacePhase = false;
691  }
692  ThermoPhase& THref = m_kin->thermo(ip);
693  int nsp = THref.nSpecies();
694  string pname = THref.id();
695  for (int k = 0; k < nsp; k++) {
696  string sname = THref.speciesName(k);
697  string cname = pname + ":" + sname;
698  if (inSurfacePhase) {
699  printf("\t %-24s %10.3e %10.3e %d\n", cname.c_str(),
700  netProdRate[kspindex], XMolKinSpecies[kspindex],
701  kindexSP);
702  kindexSP++;
703  } else {
704  printf("\t %-24s %10.3e %10.3e\n", cname.c_str(),
705  netProdRate[kspindex], XMolKinSpecies[kspindex]);
706  }
707  kspindex++;
708  }
709  }
710  printf("=========================================================="
711  "=================================\n");
712  }
713  }
714 #endif
715  if (ioflag == 1) {
716  printf("\n\n\t Iter Time Del_t Damp DelX "
717  " Resid Name-Time Name-Damp\n");
718  printf("\t -----------------------------------------------"
719  "------------------------------------\n");
720  }
721 }
722 //================================================================================================
723 void solveProb::printIteration(int ioflag, doublereal damp, size_t label_d,
724  size_t label_t,
725  doublereal inv_t, doublereal t_real, int iter,
726  doublereal update_norm, doublereal resid_norm,
727  doublereal netProdRate[], doublereal CSolnSP[],
728  doublereal resid[],
729  doublereal wtSpecies[], size_t dim, bool do_time)
730 {
731  size_t i, k;
732  string nm;
733  if (ioflag == 1) {
734 
735  printf("\t%6d ", iter);
736  if (do_time) {
737  printf("%9.4e %9.4e ", t_real, 1.0/inv_t);
738  } else
739  for (i = 0; i < 22; i++) {
740  printf(" ");
741  }
742  if (damp < 1.0) {
743  printf("%9.4e ", damp);
744  } else
745  for (i = 0; i < 11; i++) {
746  printf(" ");
747  }
748  printf("%9.4e %9.4e", update_norm, resid_norm);
749  if (do_time) {
750  k = label_t;
751  printf(" %s", int2str(k).c_str());
752  } else {
753  for (i = 0; i < 16; i++) {
754  printf(" ");
755  }
756  }
757  if (label_d != npos) {
758  k = label_d;
759  printf(" %s", int2str(k).c_str());
760  }
761  printf("\n");
762  }
763 #ifdef DEBUG_SOLVEPROB
764  else if (ioflag > 1) {
765 
766  updateMFSolnSP(XMolSolnSP);
767  printf("\n\t Weighted norm of update = %10.4e\n", update_norm);
768 
769  printf("\t Name Prod_Rate XMol Conc "
770  " Conc_Old wtConc");
771  if (damp < 1.0) {
772  printf(" UnDamped_Conc");
773  }
774  printf("\n");
775  printf("\t---------------------------------------------------------"
776  "-----------------------------\n");
777  int kindexSP = 0;
778  for (int isp = 0; isp < m_numSurfPhases; isp++) {
779  int nsp = m_nSpeciesSurfPhase[isp];
780  InterfaceKinetics* m_kin = m_objects[isp];
781  //int surfPhaseIndex = m_kinObjPhaseIDSurfPhase[isp];
783  for (int k = 0; k < nsp; k++, kindexSP++) {
784  int kspIndex = m_kinSpecIndex[kindexSP];
785  nm = m_kin->kineticsSpeciesName(kspIndex);
786  printf("\t%-16s %10.3e %10.3e %10.3e %10.3e %10.3e ",
787  nm.c_str(),
788  m_numEqn1[kspIndex],
789  XMolSolnSP[kindexSP],
790  CSolnSP[kindexSP], CSolnSP[kindexSP]+damp*resid[kindexSP],
791  wtSpecies[kindexSP]);
792  if (damp < 1.0) {
793  printf("%10.4e ", CSolnSP[kindexSP]+(damp-1.0)*resid[kindexSP]);
794  if (label_d == kindexSP) {
795  printf(" Damp ");
796  }
797  }
798  if (label_t == kindexSP) {
799  printf(" Tctrl");
800  }
801  printf("\n");
802  }
803 
804  }
805 
806  printf("\t--------------------------------------------------------"
807  "------------------------------\n");
808  }
809 #endif
810 }
811 
812 void solveProb::printFinal(int ioflag, doublereal damp, size_t label_d, size_t label_t,
813  doublereal inv_t, doublereal t_real, int iter,
814  doublereal update_norm, doublereal resid_norm,
815  doublereal netProdRateKinSpecies[], const doublereal CSolnSP[],
816  const doublereal resid[],
817  const doublereal wtSpecies[], const doublereal wtRes[],
818  size_t dim, bool do_time)
819 {
820  size_t i, k;
821  string nm;
822  if (ioflag == 1) {
823 
824  printf("\tFIN%3d ", iter);
825  if (do_time) {
826  printf("%9.4e %9.4e ", t_real, 1.0/inv_t);
827  } else
828  for (i = 0; i < 22; i++) {
829  printf(" ");
830  }
831  if (damp < 1.0) {
832  printf("%9.4e ", damp);
833  } else
834  for (i = 0; i < 11; i++) {
835  printf(" ");
836  }
837  printf("%9.4e %9.4e", update_norm, resid_norm);
838  if (do_time) {
839  k = label_t;
840  printf(" %s", int2str(k).c_str());
841  } else {
842  for (i = 0; i < 16; i++) {
843  printf(" ");
844  }
845  }
846  if (label_d != npos) {
847  k = label_d;
848 
849  printf(" %s", int2str(k).c_str());
850  }
851  printf(" -- success\n");
852  }
853 #ifdef DEBUG_SOLVEPROB
854  else if (ioflag > 1) {
855 
856 
857  printf("\n================================== FINAL RESULT ========="
858  "==================================================\n");
859 
860  printf("\n Weighted norm of solution update = %10.4e\n", update_norm);
861  printf(" Weighted norm of residual update = %10.4e\n\n", resid_norm);
862 
863  printf(" Name Prod_Rate XMol Conc "
864  " wtConc Resid Resid/wtResid wtResid");
865  if (damp < 1.0) {
866  printf(" UnDamped_Conc");
867  }
868  printf("\n");
869  printf("---------------------------------------------------------------"
870  "---------------------------------------------\n");
871 
872  for (int k = 0; k < m_neq; k++, k++) {
873  printf("%-16s %10.3e %10.3e %10.3e %10.3e %10.3e %10.3e %10.3e",
874  nm.c_str(),
875  m_numEqn1[k],
876  XMolSolnSP[k],
877  CSolnSP[k],
878  wtSpecies[k],
879  resid[k],
880  resid[k]/wtRes[k], wtRes[k]);
881  if (damp < 1.0) {
882  printf("%10.4e ", CSolnSP[k]+(damp-1.0)*resid[k]);
883  if (label_d == k) {
884  printf(" Damp ");
885  }
886  }
887  if (label_t == k) {
888  printf(" Tctrl");
889  }
890  printf("\n");
891  }
892 
893  printf("\n");
894  printf("==============================================================="
895  "============================================\n\n");
896  }
897 #endif
898 }
899 
900 #ifdef DEBUG_SOLVEPROB
901 void solveProb::
902 printIterationHeader(int ioflag, doublereal damp,doublereal inv_t, doublereal t_real,
903  int iter, bool do_time)
904 {
905  if (ioflag > 1) {
906  printf("\n===============================Iteration %5d "
907  "=================================\n", iter);
908  if (do_time) {
909  printf(" Transient step with: Real Time_n-1 = %10.4e sec,", t_real);
910  printf(" Time_n = %10.4e sec\n", t_real + 1.0/inv_t);
911  printf(" Delta t = %10.4e sec", 1.0/inv_t);
912  } else {
913  printf(" Steady Solve ");
914  }
915  if (damp < 1.0) {
916  printf(", Damping value = %10.4e\n", damp);
917  } else {
918  printf("\n");
919  }
920  }
921 }
922 #endif
923 
924 void solveProb::setAtol(const doublereal atol[])
925 {
926  for (size_t k = 0; k < m_neq; k++, k++) {
927  m_atol[k] = atol[k];
928  }
929 }
930 
931 void solveProb::setAtolConst(const doublereal atolconst)
932 {
933  for (size_t k = 0; k < m_neq; k++, k++) {
934  m_atol[k] = atolconst;
935  }
936 }
937 
938 }
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:40
virtual void setBounds(const doublereal botBounds[], const doublereal topBounds[])
Set the bottom and top bounds on the solution vector.
Definition: solveProb.cpp:609
const int SOLVEPROB_INITIALIZE
Solution Methods.
Definition: solveProb.h:52
virtual void getNetProductionRates(doublereal *net)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
void resize(size_t n, size_t m, doublereal v=0.0)
Resize the array, and fill the new entries with 'v'.
Definition: Array.h:128
Declarations for a simple class that implements an Ansi C wall clock timer (see Cantera::clockWC).
virtual doublereal calc_t(doublereal netProdRateSolnSP[], doublereal Csoln[], size_t *label, size_t *label_old, doublereal *label_factor, int ioflag)
Calculate a conservative delta T to use in a pseudo-steady state algorithm.
Definition: solveProb.cpp:557
thermo_t & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism...
Definition: Kinetics.h:288
virtual doublereal calc_damping(doublereal x[], doublereal dxneg[], size_t dim, size_t *label)
This function calculates a damping factor for the Newton iteration update vector, dxneg...
Definition: solveProb.cpp:440
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:173
std::string kineticsSpeciesName(size_t k) const
Return the name of the kth species in the kinetics manager.
Definition: Kinetics.cpp:186
vector_fp m_wtSpecies
Weights for the species concentrations norm calculation.
Definition: solveProb.h:403
The class provides the wall clock timer in seconds.
Definition: clockWC.h:51
int solve(int ifunc, doublereal time_scale, doublereal reltol)
Main routine that actually calculates the pseudo steady state of the surface problem.
Definition: solveProb.cpp:67
virtual void printIteration(int ioflag, doublereal damp, size_t label_d, size_t label_t, doublereal inv_t, doublereal t_real, int iter, doublereal update_norm, doublereal resid_norm, doublereal netProdRate[], doublereal CSolnSP[], doublereal resid[], doublereal wtSpecies[], size_t dim, bool do_time)
Printing routine that gets called after every iteration.
Definition: solveProb.cpp:723
doublereal * ptrColumn(size_t j)
Return a pointer to the top of column j, columns are contiguous in memory.
Definition: Array.h:363
size_t m_neq
Total number of equations to solve in the implicit problem.
Definition: solveProb.h:349
vector_fp m_numEqn1
Temporary vector with length MAX(1, m_neq)
Definition: solveProb.h:367
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:33
ResidEval * m_residFunc
residual function pointer to be solved.
Definition: solveProb.h:343
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:101
vector_fp m_botBounds
Bottom bounds for the solution vector.
Definition: solveProb.h:447
vector_fp m_CSolnSPOld
Saved solution vector at the old time step.
Definition: solveProb.h:391
size_t surfacePhaseIndex()
This returns the integer index of the phase which has ThermoPhase type cSurf.
Definition: Kinetics.h:263
virtual void reportState(doublereal *const CSoln) const
Report the current state of the solution.
Definition: solveProb.cpp:386
virtual int getInitialConditions(const doublereal t0, doublereal *const y, doublereal *const ydot)
Fill in the initial conditions.
Definition: ResidEval.h:122
Virtual base class for DAE residual function evaluators.
Definition: ResidEval.h:33
size_t nPhases() const
The number of phases participating in the reaction mechanism.
Definition: Kinetics.h:228
A kinetics manager for heterogeneous reaction mechanisms.
double secondsWC()
Returns the wall clock time in seconds since the last reset.
Definition: clockWC.cpp:35
doublereal m_rtol
m_rtol is the relative error tolerance.
Definition: solveProb.h:355
vector_fp m_atol
m_atol is the absolute tolerance in real units.
Definition: solveProb.h:352
vector_fp m_resid
Residual for the surface problem.
Definition: solveProb.h:416
std::string id() const
Return the string id for the phase.
Definition: Phase.cpp:119
vector_fp m_CSolnSP
Solution vector.
Definition: solveProb.h:379
virtual int nEquations() const =0
Return the number of equations in the equation system.
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:252
Array2D m_Jac
Jacobian.
Definition: solveProb.h:435
virtual void fun_eval(doublereal *const resid, const doublereal *const CSolnSP, const doublereal *const CSolnSPOld, const bool do_time, const doublereal deltaT)
Main Function evaluation.
Definition: solveProb.cpp:391
virtual 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: solveProb.cpp:406
std::vector< doublereal * > m_JacCol
Vector of pointers to the top of the columns of the jacobians.
Definition: solveProb.h:429
vector_fp m_CSolnSPInit
Saved initial solution vector.
Definition: solveProb.h:385
Header file for implicit nonlinear solver with the option of a pseudotransient (see Numerical Utiliti...
Contains declarations for string manipulation functions within Cantera.
virtual void calcWeights(doublereal wtSpecies[], doublereal wtResid[], const doublereal CSolnSP[])
Calculate the solution and residual weights.
Definition: solveProb.cpp:532
#define DATA_PTR(vec)
Creates a pointer to the start of the raw data for a vector.
Definition: ct_defs.h:36
vector_fp m_netProductionRatesSave
Temporary vector with length MAX(1, m_neq)
Definition: solveProb.h:364
vector_int m_ipiv
pivots
Definition: solveProb.h:422
virtual void print_header(int ioflag, int ifunc, doublereal time_scale, doublereal reltol, doublereal netProdRate[])
Printing routine that gets called at the start of every invocation.
Definition: solveProb.cpp:626
vector_fp m_CSolnSave
Temporary vector with length MAX(1, m_neq)
Definition: solveProb.h:373
vector_fp m_numEqn2
Temporary vector with length MAX(1, m_neq)
Definition: solveProb.h:370
virtual void printFinal(int ioflag, doublereal damp, size_t label_d, size_t label_t, doublereal inv_t, doublereal t_real, int iter, doublereal update_norm, doublereal resid_norm, doublereal netProdRateKinSpecies[], const doublereal CSolnSP[], const doublereal resid[], const doublereal wtSpecies[], const doublereal wtRes[], size_t dim, bool do_time)
Print a summary of the solution.
Definition: solveProb.cpp:812
vector_fp m_wtResid
Weights for the residual norm calculation.
Definition: solveProb.h:397
vector_fp m_topBounds
Top bounds for the solution vector.
Definition: solveProb.h:441
std::string speciesName(size_t k) const
Name of the species with index k.
Definition: Phase.cpp:246