Cantera  2.1.2
BEulerInt.cpp
Go to the documentation of this file.
1 /**
2  * @file BEulerInt.cpp
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 
14 using namespace std;
15 
16 namespace Cantera
17 {
18 
19 BEulerErr::BEulerErr(const std::string& msg) :
20  CanteraError("BEulerInt", msg)
21 {
22 }
23 
25  m_iter(Newton_Iter),
26  m_method(BEulerVarStep),
27  m_jacFormMethod(BEULER_JAC_NUM),
28  m_rowScaling(true),
29  m_colScaling(false),
30  m_matrixConditioning(false),
31  m_itol(0),
32  m_reltol(1.e-4),
33  m_abstols(1.e-10),
34  m_abstol(0),
35  m_ewt(0),
36  m_hmax(0.0),
37  m_maxord(0),
38  m_time_step_num(0),
39  m_time_step_attempts(0),
40  m_max_time_step_attempts(11000000),
41  m_numInitialConstantDeltaTSteps(0),
42  m_failure_counter(0),
43  m_min_newt_its(0),
44  m_printSolnStepInterval(1),
45  m_printSolnNumberToTout(1),
46  m_printSolnFirstSteps(0),
47  m_dumpJacobians(false),
48  m_neq(0),
49  m_y_n(0),
50  m_y_nm1(0),
51  m_y_pred_n(0),
52  m_ydot_n(0),
53  m_ydot_nm1(0),
54  m_t0(0.0),
55  m_time_final(0.0),
56  time_n(0.0),
57  time_nm1(0.0),
58  time_nm2(0.0),
59  delta_t_n(0.0),
60  delta_t_nm1(0.0),
61  delta_t_nm2(0.0),
62  delta_t_np1(1.0E-8),
63  delta_t_max(1.0E300),
64  m_resid(0),
65  m_residWts(0),
66  m_wksp(0),
67  m_func(0),
68  m_rowScales(0),
69  m_colScales(0),
70  tdjac_ptr(0),
71  m_print_flag(3),
72  m_nfe(0),
73  m_nJacEval(0),
74  m_numTotalNewtIts(0),
75  m_numTotalLinearSolves(0),
76  m_numTotalConvFails(0),
77  m_numTotalTruncFails(0),
78  num_failures(0)
79 {
80 
81 }
82 
83 BEulerInt::~BEulerInt()
84 {
85  delete tdjac_ptr;
86 }
87 
88 void BEulerInt::setTolerances(double reltol, size_t n, double* abstol)
89 {
90  m_itol = 1;
91  m_abstol.resize(m_neq);
92  if (static_cast<int>(n) != m_neq) {
93  printf("ERROR n is wrong\n");
94  exit(-1);
95  }
96  for (int i = 0; i < m_neq; i++) {
97  m_abstol[i] = abstol[i];
98  }
99  m_reltol = reltol;
100 }
101 
102 void BEulerInt::setTolerances(double reltol, double abstol)
103 {
104  m_itol = 0;
105  m_reltol = reltol;
106  m_abstols = abstol;
107 }
108 
109 void BEulerInt::setProblemType(int jacFormMethod)
110 {
111  m_jacFormMethod = jacFormMethod;
112 }
113 
114 void BEulerInt::setMethodBEMT(BEulerMethodType t)
115 {
116  m_method = t;
117 }
118 
119 void BEulerInt::setMaxStep(doublereal hmax)
120 {
121  m_hmax = hmax;
122 }
123 
124 void BEulerInt::setMaxNumTimeSteps(int maxNumTimeSteps)
125 {
126  m_max_time_step_attempts = maxNumTimeSteps;
127 }
128 
129 void BEulerInt::setNumInitialConstantDeltaTSteps(int num)
130 {
132 }
133 
134 void BEulerInt::setPrintSolnOptions(int printSolnStepInterval,
135  int printSolnNumberToTout,
136  int printSolnFirstSteps,
137  bool dumpJacobians)
138 {
139  m_printSolnStepInterval = printSolnStepInterval;
140  m_printSolnNumberToTout = printSolnNumberToTout;
141  m_printSolnFirstSteps = printSolnFirstSteps;
142  m_dumpJacobians = dumpJacobians;
143 }
144 
146 {
147  m_iter = t;
148 }
149 
150 void BEulerInt::setNonLinOptions(int min_newt_its, bool matrixConditioning,
151  bool colScaling, bool rowScaling)
152 {
153  m_min_newt_its = min_newt_its;
154  m_matrixConditioning = matrixConditioning;
155  m_colScaling = colScaling;
156  m_rowScaling = rowScaling;
157  if (m_colScaling && m_colScales.empty()) {
158  m_colScales.assign(m_neq, 1.0);
159  }
160  if (m_rowScaling && m_rowScales.empty()) {
161  m_rowScales.assign(m_neq, 1.0);
162  }
163 }
164 
165 void BEulerInt::setInitialTimeStep(double deltaT)
166 {
167  delta_t_np1 = deltaT;
168 }
169 
170 void BEulerInt::setPrintFlag(int print_flag)
171 {
172  m_print_flag = print_flag;
173 }
174 
176 {
177  m_neq = func.nEquations();
178  m_t0 = t0;
179  internalMalloc();
180 
181  /*
182  * Get the initial conditions.
183  */
184  func.getInitialConditions(m_t0, &m_y_n[0], &m_ydot_n[0]);
185 
186  // Store a pointer to the residual routine in the object
187  m_func = &func;
188 
189  /*
190  * Initialize the various time counters in the object
191  */
192  time_n = t0;
193  time_nm1 = time_n;
194  time_nm2 = time_nm1;
195  delta_t_n = 0.0;
196  delta_t_nm1 = 0.0;
197 }
198 
199 void BEulerInt::reinitializeRJE(double t0, ResidJacEval& func)
200 {
201  m_neq = func.nEquations();
202  m_t0 = t0;
203  internalMalloc();
204  /*
205  * At the initial time, get the initial conditions and time and store
206  * them into internal storage in the object, my[].
207  */
208  m_t0 = t0;
209  func.getInitialConditions(m_t0, &m_y_n[0], &m_ydot_n[0]);
210  /**
211  * Set up the internal weights that are used for testing convergence
212  */
213  setSolnWeights();
214 
215  // Store a pointer to the function
216  m_func = &func;
217 
218 }
219 
220 double BEulerInt::getPrintTime(double time_current)
221 {
222  double tnext;
223  if (m_printSolnNumberToTout > 0) {
224  double dt = (m_time_final - m_t0) / m_printSolnNumberToTout;
225  for (int i = 0; i <= m_printSolnNumberToTout; i++) {
226  tnext = m_t0 + dt * i;
227  if (tnext >= time_current) {
228  return tnext;
229  }
230  }
231  }
232  return 1.0E300;
233 }
234 
235 int BEulerInt::nEvals() const
236 {
237  return m_nfe;
238 }
239 
241 {
242  m_ewt.assign(m_neq, 0.0);
243  m_y_n.assign(m_neq, 0.0);
244  m_y_nm1.assign(m_neq, 0.0);
245  m_y_pred_n.assign(m_neq, 0.0);
246  m_ydot_n.assign(m_neq, 0.0);
247  m_ydot_nm1.assign(m_neq, 0.0);
248  m_resid.assign(m_neq, 0.0);
249  m_residWts.assign(m_neq, 0.0);
250  m_wksp.assign(m_neq, 0.0);
251  if (m_rowScaling) {
252  m_rowScales.assign(m_neq, 1.0);
253  }
254  if (m_colScaling) {
255  m_colScales.assign(m_neq, 1.0);
256  }
258 }
259 
261 {
262  int i;
263  if (m_itol == 1) {
264  /*
265  * Adjust the atol vector if we are using vector
266  * atol conditions.
267  */
268  // m_func->adjustAtol(m_abstol);
269 
270  for (i = 0; i < m_neq; i++) {
271  m_ewt[i] = m_abstol[i] + m_reltol * 0.5 *
272  (fabs(m_y_n[i]) + fabs(m_y_pred_n[i]));
273  }
274  } else {
275  for (i = 0; i < m_neq; i++) {
276  m_ewt[i] = m_abstols + m_reltol * 0.5 *
277  (fabs(m_y_n[i]) + fabs(m_y_pred_n[i]));
278  }
279  }
280 }
281 
283 {
284  m_func->calcSolnScales(time_n, &m_y_n[0], &m_y_nm1[0], &m_colScales[0]);
285 }
286 
288 {
289  /*
290  * We compute residual weights here, which we define as the L_0 norm
291  * of the Jacobian Matrix, weighted by the solution weights.
292  * This is the proper way to guage the magnitude of residuals. However,
293  * it does need the evaluation of the jacobian, and the implementation
294  * below is slow, but doesn't take up much memory.
295  *
296  * Here a small weighting indicates that the change in solution is
297  * very sensitive to that equation.
298  */
299  int i, j;
300  double* data = &(*(jac.begin()));
301  double value;
302  for (i = 0; i < m_neq; i++) {
303  m_residWts[i] = fabs(data[i] * m_ewt[0]);
304  for (j = 1; j < m_neq; j++) {
305  value = fabs(data[j*m_neq + i] * m_ewt[j]);
306  m_residWts[i] = std::max(m_residWts[i], value);
307  }
308  }
309 }
310 
311 double BEulerInt::filterNewStep(double timeCurrent, double* y_current, double* ydot_current)
312 {
313  return 0.0;
314 }
315 
316 static void print_line(const char* str, int n)
317 {
318  for (int i = 0; i < n; i++) {
319  printf("%s", str);
320  }
321  printf("\n");
322 }
323 
324 /*
325  * Print out for relevant time step information
326  */
327 static void print_time_step1(int order, int n_time_step, double time,
328  double delta_t_n, double delta_t_nm1,
329  bool step_failed, int num_failures)
330 {
331  const char* string = 0;
332  if (order == 0) {
333  string = "Backward Euler";
334  } else if (order == 1) {
335  string = "Forward/Backward Euler";
336  } else if (order == 2) {
337  string = "Adams-Bashforth/TR";
338  }
339  printf("\n");
340  print_line("=", 80);
341  printf("\nStart of Time Step: %5d Time_n = %9.5g Time_nm1 = %9.5g\n",
342  n_time_step, time, time - delta_t_n);
343  printf("\tIntegration method = %s\n", string);
344  if (step_failed) {
345  printf("\tPreviously attempted step was a failure\n");
346  }
347  if (delta_t_n > delta_t_nm1) {
348  string = "(Increased from previous iteration)";
349  } else if (delta_t_n < delta_t_nm1) {
350  string = "(Decreased from previous iteration)";
351  } else {
352  string = "(same as previous iteration)";
353  }
354  printf("\tdelta_t_n = %8.5e %s", delta_t_n, string);
355  if (num_failures > 0) {
356  printf("\t(Bad_History Failure Counter = %d)", num_failures);
357  }
358  printf("\n\tdelta_t_nm1 = %8.5e\n", delta_t_nm1);
359 }
360 
361 /*
362  * Print out for relevant time step information
363  */
364 static void print_time_step2(int time_step_num, int order,
365  double time, double time_error_factor,
366  double delta_t_n, double delta_t_np1)
367 {
368  printf("\tTime Step Number %5d was a success: time = %10g\n", time_step_num,
369  time);
370  printf("\t\tEstimated Error\n");
371  printf("\t\t-------------------- = %8.5e\n", time_error_factor);
372  printf("\t\tTolerated Error\n\n");
373  printf("\t- Recommended next delta_t (not counting history) = %g\n",
374  delta_t_np1);
375  printf("\n");
376  print_line("=", 80);
377  printf("\n");
378 }
379 
380 /*
381  * Print Out descriptive information on why the current step failed
382  */
383 static void print_time_fail(bool convFailure, int time_step_num,
384  double time, double delta_t_n,
385  double delta_t_np1, double time_error_factor)
386 {
387  printf("\n");
388  print_line("=", 80);
389  if (convFailure) {
390  printf("\tTime Step Number %5d experienced a convergence "
391  "failure\n", time_step_num);
392  printf("\tin the non-linear or linear solver\n");
393  printf("\t\tValue of time at failed step = %g\n", time);
394  printf("\t\tdelta_t of the failed step = %g\n",
395  delta_t_n);
396  printf("\t\tSuggested value of delta_t to try next = %g\n",
397  delta_t_np1);
398  } else {
399  printf("\tTime Step Number %5d experienced a truncation error "
400  "failure!\n", time_step_num);
401  printf("\t\tValue of time at failed step = %g\n", time);
402  printf("\t\tdelta_t of the failed step = %g\n",
403  delta_t_n);
404  printf("\t\tSuggested value of delta_t to try next = %g\n",
405  delta_t_np1);
406  printf("\t\tCalculated truncation error factor = %g\n",
407  time_error_factor);
408  }
409  printf("\n");
410  print_line("=", 80);
411 }
412 
413 /*
414  * Print out the final results and counters
415  */
416 static void print_final(double time, int step_failed,
417  int time_step_num, int num_newt_its,
418  int total_linear_solves, int numConvFails,
419  int numTruncFails, int nfe, int nJacEval)
420 {
421  printf("\n");
422  print_line("=", 80);
423  printf("TIME INTEGRATION ROUTINE HAS FINISHED: ");
424  if (step_failed) {
425  printf(" IT WAS A FAILURE\n");
426  } else {
427  printf(" IT WAS A SUCCESS\n");
428  }
429  printf("\tEnding time = %g\n", time);
430  printf("\tNumber of time steps = %d\n", time_step_num);
431  printf("\tNumber of newt its = %d\n", num_newt_its);
432  printf("\tNumber of linear solves = %d\n", total_linear_solves);
433  printf("\tNumber of convergence failures= %d\n", numConvFails);
434  printf("\tNumber of TimeTruncErr fails = %d\n", numTruncFails);
435  printf("\tNumber of Function evals = %d\n", nfe);
436  printf("\tNumber of Jacobian evals/solvs= %d\n", nJacEval);
437  printf("\n");
438  print_line("=", 80);
439 }
440 
441 /*
442  * Header info for one line comment about a time step
443  */
444 static void print_lvl1_Header(int nTimes)
445 {
446  printf("\n");
447  if (nTimes) {
448  print_line("-", 80);
449  }
450  printf("time Time Time Time ");
451  if (nTimes == 0) {
452  printf(" START");
453  } else {
454  printf(" (continued)");
455  }
456  printf("\n");
457 
458  printf("step (sec) step Newt Aztc bktr trunc ");
459  printf("\n");
460 
461  printf(" No. Rslt size Its Its stps error |");
462  printf(" comment");
463  printf("\n");
464  print_line("-", 80);
465 }
466 
467 /*
468  * One line entry about time step
469  * rslt -> 4 letter code
470  */
471 static void print_lvl1_summary(
472  int time_step_num, double time, const char* rslt, double delta_t_n,
473  int newt_its, int aztec_its, int bktr_stps, double time_error_factor,
474  const char* comment)
475 {
476  printf("%6d %11.6g %4s %10.4g %4d %4d %4d %11.4g",
477  time_step_num, time, rslt, delta_t_n, newt_its, aztec_its,
478  bktr_stps, time_error_factor);
479  if (comment) {
480  printf(" | %s", comment);
481  }
482  printf("\n");
483 }
484 
485 /*
486  * This routine subtracts 2 numbers. If the difference is less
487  * than 1.0E-14 times the magnitude of the smallest number,
488  * then diff returns an exact zero.
489  * It also returns an exact zero if the difference is less than
490  * 1.0E-300.
491  *
492  * returns: a - b
493  *
494  * This routine is used in numerical differencing schemes in order
495  * to avoid roundoff errors resulting in creating Jacobian terms.
496  * Note: This is a slow routine. However, jacobian errors may cause
497  * loss of convergence. Therefore, in practice this routine
498  * has proved cost-effective.
499  */
500 double subtractRD(double a, double b)
501 {
502  double diff = a - b;
503  double d = std::min(fabs(a), fabs(b));
504  d *= 1.0E-14;
505  double ad = fabs(diff);
506  if (ad < 1.0E-300) {
507  diff = 0.0;
508  }
509  if (ad < d) {
510  diff = 0.0;
511  }
512  return diff;
513 }
514 
515 void BEulerInt::beuler_jac(GeneralMatrix& J, double* const f,
516  double time_curr, double CJ,
517  double* const y,
518  double* const ydot,
519  int num_newt_its)
520 {
521  int i, j;
522  double* col_j;
523  double ysave, ydotsave, dy;
524  /**
525  * Clear the factor flag
526  */
527  J.clearFactorFlag();
528 
529 
530  if (m_jacFormMethod & BEULER_JAC_ANAL) {
531  /********************************************************************
532  * Call the function to get a jacobian.
533  */
534  m_func->evalJacobian(time_curr, delta_t_n, CJ, y, ydot, J, f);
535 #ifdef DEBUG_HKM
536  //double dddd = J(89, 89);
537  //checkFinite(dddd);
538 #endif
539  m_nJacEval++;
540  m_nfe++;
541  } else {
542  /*******************************************************************
543  * Generic algorithm to calculate a numerical Jacobian
544  */
545  /*
546  * Calculate the current value of the rhs given the
547  * current conditions.
548  */
549 
550  m_func->evalResidNJ(time_curr, delta_t_n, y, ydot, f, JacBase_ResidEval);
551  m_nfe++;
552  m_nJacEval++;
553 
554 
555  /*
556  * Malloc a vector and call the function object to return a set of
557  * deltaY's that are appropriate for calculating the numerical
558  * derivative.
559  */
560  vector_fp dyVector(m_neq);
561  m_func->calcDeltaSolnVariables(time_curr, y, &m_y_nm1[0], &dyVector[0],
562  &m_ewt[0]);
563 #ifdef DEBUG_HKM
564  bool print_NumJac = false;
565  if (print_NumJac) {
566  FILE* idy = fopen("NumJac.csv", "w");
567  fprintf(idy, "Unk m_ewt y "
568  "dyVector ResN\n");
569  for (int iii = 0; iii < m_neq; iii++) {
570  fprintf(idy, " %4d %16.8e %16.8e %16.8e %16.8e \n",
571  iii, m_ewt[iii], y[iii], dyVector[iii], f[iii]);
572  }
573  fclose(idy);
574  }
575 #endif
576  /*
577  * Loop over the variables, formulating a numerical derivative
578  * of the dense matrix.
579  * For the delta in the variable, we will use a variety of approaches
580  * The original approach was to use the error tolerance amount.
581  * This may not be the best approach, as it could be overly large in
582  * some instances and overly small in others.
583  * We will first protect from being overly small, by using the usual
584  * sqrt of machine precision approach, i.e., 1.0E-7,
585  * to bound the lower limit of the delta.
586  */
587  for (j = 0; j < m_neq; j++) {
588 
589 
590  /*
591  * Get a pointer into the column of the matrix
592  */
593 
594 
595  col_j = (double*) J.ptrColumn(j);
596  ysave = y[j];
597  dy = dyVector[j];
598  //dy = fmaxx(1.0E-6 * m_ewt[j], fabs(ysave)*1.0E-7);
599 
600  y[j] = ysave + dy;
601  dy = y[j] - ysave;
602  ydotsave = ydot[j];
603  ydot[j] += dy * CJ;
604  /*
605  * Call the functon
606  */
607 
608 
609  m_func->evalResidNJ(time_curr, delta_t_n, y, ydot, &m_wksp[0],
610  JacDelta_ResidEval, j, dy);
611  m_nfe++;
612  double diff;
613  for (i = 0; i < m_neq; i++) {
614  diff = subtractRD(m_wksp[i], f[i]);
615  col_j[i] = diff / dy;
616  //col_j[i] = (m_wksp[i] - f[i])/dy;
617  }
618 
619  y[j] = ysave;
620  ydot[j] = ydotsave;
621 
622  }
623  }
624 
625 
626 }
627 
628 void BEulerInt::calc_y_pred(int order)
629 {
630  int i;
631  double c1, c2;
632  switch (order) {
633  case 0:
634  case 1:
635  c1 = delta_t_n;
636  for (i = 0; i < m_neq; i++) {
637  m_y_pred_n[i] = m_y_n[i] + c1 * m_ydot_n[i];
638  }
639  break;
640  case 2:
641  c1 = delta_t_n * (2.0 + delta_t_n / delta_t_nm1) / 2.0;
642  c2 = (delta_t_n * delta_t_n) / (delta_t_nm1 * 2.0);
643  for (i = 0; i < m_neq; i++) {
644  m_y_pred_n[i] = m_y_n[i] + c1 * m_ydot_n[i] - c2 * m_ydot_nm1[i];
645  }
646  break;
647  }
648 
649  /*
650  * Filter the predictions.
651  */
652  m_func->filterSolnPrediction(time_n, &m_y_pred_n[0]);
653 
654 }
655 
656 void BEulerInt::
657 calc_ydot(int order, double* y_curr, double* ydot_curr)
658 {
659  int i;
660  double c1;
661  switch (order) {
662  case 0:
663  case 1: /* First order forward Euler/backward Euler */
664  c1 = 1.0 / delta_t_n;
665  for (i = 0; i < m_neq; i++) {
666  ydot_curr[i] = c1 * (y_curr[i] - m_y_nm1[i]);
667  }
668  return;
669  case 2: /* Second order Adams-Bashforth / Trapezoidal Rule */
670  c1 = 2.0 / delta_t_n;
671  for (i = 0; i < m_neq; i++) {
672  ydot_curr[i] = c1 * (y_curr[i] - m_y_nm1[i]) - m_ydot_nm1[i];
673  }
674  return;
675  }
676 }
677 
679 {
680  int i;
681  double rel_norm, error;
682 #ifdef DEBUG_HKM
683 #define NUM_ENTRIES 5
684  if (m_print_flag > 2) {
685  int imax[NUM_ENTRIES], j, jnum;
686  double dmax;
687  bool used;
688  printf("\t\ttime step truncation error contributors:\n");
689  printf("\t\t I entry actual predicted "
690  " weight ydot\n");
691  printf("\t\t");
692  print_line("-", 70);
693  for (j = 0; j < NUM_ENTRIES; j++) {
694  imax[j] = -1;
695  }
696  for (jnum = 0; jnum < NUM_ENTRIES; jnum++) {
697  dmax = -1.0;
698  for (i = 0; i < m_neq; i++) {
699  used = false;
700  for (j = 0; j < jnum; j++) {
701  if (imax[j] == i) {
702  used = true;
703  }
704  }
705  if (!used) {
706  error = (m_y_n[i] - m_y_pred_n[i]) / m_ewt[i];
707  rel_norm = sqrt(error * error);
708  if (rel_norm > dmax) {
709  imax[jnum] = i;
710  dmax = rel_norm;
711  }
712  }
713  }
714  if (imax[jnum] >= 0) {
715  i = imax[jnum];
716  printf("\t\t%4d %12.4e %12.4e %12.4e %12.4e %12.4e\n",
717  i, dmax, m_y_n[i], m_y_pred_n[i], m_ewt[i], m_ydot_n[i]);
718  }
719  }
720  printf("\t\t");
721  print_line("-", 70);
722  }
723 #endif
724  rel_norm = 0.0;
725  for (i = 0; i < m_neq; i++) {
726  error = (m_y_n[i] - m_y_pred_n[i]) / m_ewt[i];
727  rel_norm += (error * error);
728  }
729  return sqrt(rel_norm / m_neq);
730 }
731 
732 double BEulerInt::time_step_control(int order, double time_error_factor)
733 {
734  double factor = 0.0, power = 0.0, delta_t;
735  const char* yo = "time_step_control";
736 
737  /*
738  * Special case time_error_factor so that zeroes don't cause a problem.
739  */
740  time_error_factor = std::max(1.0E-50, time_error_factor);
741 
742  /*
743  * Calculate the factor for the change in magnitude of time step.
744  */
745  switch (order) {
746  case 1:
747  factor = 1.0/(2.0 *(time_error_factor));
748  power = 0.5;
749  break;
750  case 2:
751  factor = 1.0/(3.0 * (1.0 + delta_t_nm1 / delta_t_n)
752  * (time_error_factor));
753  power = 0.3333333333333333;
754  }
755  factor = pow(factor, power);
756  if (factor < 0.5) {
757  if (m_print_flag > 1) {
758  printf("\t%s: WARNING - Current time step will be chucked\n", yo);
759  printf("\t\tdue to a time step truncation error failure.\n");
760  }
761  delta_t = - 0.5 * delta_t_n;
762  } else {
763  factor = std::min(factor, 1.5);
764  delta_t = factor * delta_t_n;
765  }
766  return delta_t;
767 }
768 
769 double BEulerInt::integrateRJE(double tout, double time_init)
770 {
771  double time_current;
772  bool weAreNotFinished = true;
773  m_time_final = tout;
774  int flag = SUCCESS;
775  /**
776  * Initialize the time step number to zero. step will increment so that
777  * the first time step is number 1
778  */
779  m_time_step_num = 0;
780 
781 
782  /*
783  * Do the integration a step at a time
784  */
785  int istep = 0;
786  int printStep = 0;
787  bool doPrintSoln = false;
788  time_current = time_init;
789  time_n = time_init;
790  time_nm1 = time_init;
791  time_nm2 = time_init;
792  m_func->evalTimeTrackingEqns(time_current, 0.0, &m_y_n[0], &m_ydot_n[0]);
793  double print_time = getPrintTime(time_current);
794  if (print_time == time_current) {
795  m_func->writeSolution(4, time_current, delta_t_n,
796  istep, &m_y_n[0], &m_ydot_n[0]);
797  }
798  /*
799  * We print out column headers here for the case of
800  */
801  if (m_print_flag == 1) {
802  print_lvl1_Header(0);
803  }
804  /*
805  * Call a different user routine at the end of each step,
806  * that will probably print to a file.
807  */
808  m_func->user_out2(0, time_current, 0.0, &m_y_n[0], &m_ydot_n[0]);
809 
810  do {
811 
812  print_time = getPrintTime(time_current);
813  if (print_time >= tout) {
814  print_time = tout;
815  }
816 
817  /************************************************************
818  * Step the solution
819  */
820  time_current = step(tout);
821  istep++;
822  printStep++;
823  /***********************************************************/
824  if (time_current < 0.0) {
825  if (time_current == -1234.) {
826  time_current = 0.0;
827  } else {
828  time_current = -time_current;
829  }
830  flag = FAILURE;
831  }
832 
833  if (flag != FAILURE) {
834  bool retn =
835  m_func->evalStoppingCritera(time_current, delta_t_n,
836  &m_y_n[0], &m_ydot_n[0]);
837  if (retn) {
838  weAreNotFinished = false;
839  doPrintSoln = true;
840  }
841  }
842 
843  /*
844  * determine conditional printing of soln
845  */
846  if (time_current >= print_time) {
847  doPrintSoln = true;
848  }
849  if (m_printSolnStepInterval == printStep) {
850  doPrintSoln = true;
851  }
852  if (m_printSolnFirstSteps > istep) {
853  doPrintSoln = true;
854  }
855 
856  /*
857  * Evaluate time integrated quantities that are calculated at the
858  * end of every successful time step.
859  */
860  if (flag != FAILURE) {
861  m_func->evalTimeTrackingEqns(time_current, delta_t_n,
862  &m_y_n[0], &m_ydot_n[0]);
863  }
864 
865  /*
866  * Call the printout routine.
867  */
868  if (doPrintSoln) {
869  m_func->writeSolution(1, time_current, delta_t_n,
870  istep, &m_y_n[0], &m_ydot_n[0]);
871  printStep = 0;
872  doPrintSoln = false;
873  if (m_print_flag == 1) {
874  print_lvl1_Header(1);
875  }
876  }
877  /*
878  * Call a different user routine at the end of each step,
879  * that will probably print to a file.
880  */
881  if (flag == FAILURE) {
882  m_func->user_out2(-1, time_current, delta_t_n, &m_y_n[0], &m_ydot_n[0]);
883  } else {
884  m_func->user_out2(1, time_current, delta_t_n, &m_y_n[0], &m_ydot_n[0]);
885  }
886 
887  } while (time_current < tout &&
888  m_time_step_attempts < m_max_time_step_attempts &&
889  flag == SUCCESS && weAreNotFinished);
890 
891  /*
892  * Check current time against the max solution time.
893  */
894  if (time_current >= tout) {
895  printf("Simulation completed time integration in %d time steps\n",
897  printf("Final Time: %e\n\n", time_current);
898  } else if (m_time_step_attempts >= m_max_time_step_attempts) {
899  printf("Simulation ran into time step attempt limit in"
900  "%d time steps\n",
902  printf("Final Time: %e\n\n", time_current);
903  } else if (flag == FAILURE) {
904  printf("ERROR: time stepper failed at time = %g\n", time_current);
905  }
906 
907  /*
908  * Print out the final results and counters.
909  */
910  print_final(time_n, flag, m_time_step_num, m_numTotalNewtIts,
913 
914  /*
915  * Call a different user routine at the end of each step,
916  * that will probably print to a file.
917  */
918  m_func->user_out2(2, time_current, delta_t_n, &m_y_n[0], &m_ydot_n[0]);
919 
920 
921  if (flag != SUCCESS) {
922  throw BEulerErr(" BEuler error encountered.");
923  }
924  return time_current;
925 }
926 
927 double BEulerInt::step(double t_max)
928 {
929  double CJ;
930  bool step_failed = false;
931  bool giveUp = false;
932  bool convFailure = false;
933  const char* rslt;
934  double time_error_factor = 0.0;
935  double normFilter = 0.0;
936  int numTSFailures = 0;
937  int bktr_stps = 0;
938  int nonlinearloglevel = m_print_flag;
939  int num_newt_its = 0;
940  int aztec_its = 0;
941  string comment;
942  /*
943  * Increment the time counter - May have to be taken back,
944  * if time step is found to be faulty.
945  */
946  m_time_step_num++;
947 
948  /**
949  * Loop here until we achieve a successful step or we set the giveUp
950  * flag indicating that repeated errors have occurred.
951  */
952  do {
953  m_time_step_attempts++;
954  comment.clear();
955 
956  /*
957  * Possibly adjust the delta_t_n value for this time step from the
958  * recommended delta_t_np1 value determined in the previous step
959  * due to maximum time step constraints or other occurences,
960  * known to happen at a given time.
961  */
962  if ((time_n + delta_t_np1) >= t_max) {
963  delta_t_np1 =t_max - time_n;
964  }
965 
966  if (delta_t_np1 >= delta_t_max) {
967  delta_t_np1 = delta_t_max;
968  }
969 
970  /*
971  * Increment the delta_t counters and the time for the current
972  * time step.
973  */
974 
975  delta_t_nm2 = delta_t_nm1;
976  delta_t_nm1 = delta_t_n;
977  delta_t_n = delta_t_np1;
978  time_n += delta_t_n;
979 
980  /*
981  * Determine the integration order of the current step.
982  *
983  * Special case for start-up of time integration procedure
984  * First time step = Do a predictor step as we
985  * have recently added an initial
986  * ydot input option. And, setting ydot=0
987  * is equivalent to not doing a
988  * predictor step.
989  * Second step = If 2nd order method, do a first order
990  * step for this time-step, only.
991  *
992  * If 2nd order method with a constant time step, the
993  * first and second steps are 1/10 the specified step, and
994  * the third step is 8/10 the specified step. This reduces
995  * the error asociated with using lower order
996  * integration on the first two steps. (RCS 11-6-97)
997  *
998  * If the previous time step failed for one reason or another,
999  * do a linear step. It's more robust.
1000  */
1001  if (m_time_step_num == 1) {
1002  m_order = 1; /* Backward Euler */
1003  } else if (m_time_step_num == 2) {
1004  m_order = 1; /* Forward/Backward Euler */
1005  } else if (step_failed) {
1006  m_order = 1; /* Forward/Backward Euler */
1007  } else if (m_time_step_num > 2) {
1008  m_order = 1; /* Specified
1009  Predictor/Corrector
1010  - not implemented */
1011  }
1012 
1013  /*
1014  * Print out an initial statement about the step.
1015  */
1016  if (m_print_flag > 1) {
1017  print_time_step1(m_order, m_time_step_num, time_n, delta_t_n,
1018  delta_t_nm1, step_failed, m_failure_counter);
1019  }
1020 
1021  /*
1022  * Calculate the predicted solution, m_y_pred_n, for the current
1023  * time step.
1024  */
1026 
1027  /*
1028  * HKM - Commented this out. I may need it for particles later.
1029  * If Solution bounds checking is turned on, we need to crop the
1030  * predicted solution to make sure bounds are enforced
1031  *
1032  *
1033  * cropNorm = 0.0;
1034  * if (Cur_Realm->Realm_Nonlinear.Constraint_Backtracking_Flag ==
1035  * Constraint_Backtrack_Enable) {
1036  * cropNorm = cropPredictor(mesh, x_pred_n, abs_time_error,
1037  * m_reltol);
1038  */
1039 
1040  /*
1041  * Save the old solution, before overwriting with the new solution
1042  * - use
1043  */
1044  m_y_nm1 = m_y_n;
1045 
1046  /*
1047  * Use the predicted value as the initial guess for the corrector
1048  * loop, for
1049  * every step other than the first step.
1050  */
1051  if (m_order > 0) {
1052  m_y_n = m_y_pred_n;
1053  }
1054 
1055  /*
1056  * Save the old time derivative, if necessary, before it is
1057  * overwritten.
1058  * This overwrites ydot_nm1, losing information from the previous time
1059  * step.
1060  */
1061  m_ydot_nm1 = m_ydot_n;
1062 
1063  /*
1064  * Calculate the new time derivative, ydot_n, that is consistent
1065  * with the
1066  * initial guess for the corrected solution vector.
1067  *
1068  */
1069  calc_ydot(m_order, &m_y_n[0], &m_ydot_n[0]);
1070 
1071  /*
1072  * Calculate CJ, the coefficient for the jacobian corresponding to the
1073  * derivative of the residual wrt to the acceleration vector.
1074  */
1075  if (m_order < 2) {
1076  CJ = 1.0 / delta_t_n;
1077  } else {
1078  CJ = 2.0 / delta_t_n;
1079  }
1080 
1081  /*
1082  * Calculate a new Solution Error Weighting vector
1083  */
1084  setSolnWeights();
1085 
1086  /*
1087  * Solve the system of equations at the current time step.
1088  * Note - x_corr_n and x_dot_n are considered to be updated,
1089  * on return from this solution.
1090  */
1091  int ierror = solve_nonlinear_problem(&m_y_n[0], &m_ydot_n[0],
1092  CJ, time_n, *tdjac_ptr, num_newt_its,
1093  aztec_its, bktr_stps,
1094  nonlinearloglevel);
1095  /*
1096  * Set the appropriate flags if a convergence failure is detected.
1097  */
1098  if (ierror < 0) { /* Step failed */
1099  convFailure = true;
1100  step_failed = true;
1101  rslt = "fail";
1103  m_failure_counter +=3;
1104  if (m_print_flag > 1) {
1105  printf("\tStep is Rejected, nonlinear problem didn't converge,"
1106  "ierror = %d\n", ierror);
1107  }
1108  } else { /* Step succeeded */
1109  convFailure = false;
1110  step_failed = false;
1111  rslt = "done";
1112 
1113  /*
1114  * Apply a filter to a new successful step
1115  */
1116  normFilter = filterNewStep(time_n, &m_y_n[0], &m_ydot_n[0]);
1117  if (normFilter > 1.0) {
1118  convFailure = true;
1119  step_failed = true;
1120  rslt = "filt";
1121  if (m_print_flag > 1) {
1122  printf("\tStep is Rejected, too large filter adjustment = %g\n",
1123  normFilter);
1124  }
1125  } else if (normFilter > 0.0) {
1126  if (normFilter > 0.3) {
1127  if (m_print_flag > 1) {
1128  printf("\tStep was filtered, norm = %g, next "
1129  "time step adjusted\n", normFilter);
1130  }
1131  } else {
1132  if (m_print_flag > 1) {
1133  printf("\tStep was filtered, norm = %g\n", normFilter);
1134  }
1135  }
1136  }
1137  }
1138 
1139  /*
1140  * Calculate the time step truncation error for the current step.
1141  */
1142  if (!step_failed) {
1143  time_error_factor = time_error_norm();
1144  } else {
1145  time_error_factor = 1000.;
1146  }
1147 
1148  /*
1149  * Dynamic time step control- delta_t_n, delta_t_nm1 are set here.
1150  */
1151  if (step_failed) {
1152  /*
1153  * For convergence failures, decrease the step-size by a factor of
1154  * 4 and try again.
1155  */
1156  delta_t_np1 = 0.25 * delta_t_n;
1157  } else if (m_method == BEulerVarStep) {
1158 
1159  /*
1160  * If we are doing a predictor/corrector method, and we are
1161  * past a certain number of time steps given by the input file
1162  * then either correct the DeltaT for the next time step or
1163  *
1164  */
1165  if ((m_order > 0) &&
1167  delta_t_np1 = time_step_control(m_order, time_error_factor);
1168  if (normFilter > 0.1) {
1169  if (delta_t_np1 > delta_t_n) {
1170  delta_t_np1 = delta_t_n;
1171  }
1172  }
1173 
1174  /*
1175  * Check for Current time step failing due to violation of
1176  * time step
1177  * truncation bounds.
1178  */
1179  if (delta_t_np1 < 0.0) {
1181  step_failed = true;
1182  delta_t_np1 = -delta_t_np1;
1183  m_failure_counter += 2;
1184  comment += "TIME TRUNC FAILURE";
1185  rslt = "TRNC";
1186  }
1187 
1188  /*
1189  * Prevent churning of the time step by not increasing the
1190  * time step,
1191  * if the recent "History" of the time step behavior is still bad
1192  */
1193  else if (m_failure_counter > 0) {
1194  delta_t_np1 = std::min(delta_t_np1, delta_t_n);
1195  }
1196  } else {
1197  delta_t_np1 = delta_t_n;
1198  }
1199 
1200  /* Decrease time step if a lot of Newton Iterations are
1201  * taken.
1202  * The idea being if more or less Newton iteration are taken
1203  * than the target number of iterations, then adjust the time
1204  * step downwards so that the target number of iterations or lower
1205  * is achieved. This
1206  * should prevent step failure by too many Newton iterations because
1207  * the time step becomes too large. CCO
1208  * hkm -> put in num_new_its min of 3 because the time step
1209  * was being altered even when num_newt_its == 1
1210  */
1211  int max_Newton_steps = 10000;
1212  int target_num_iter = 5;
1213  if (num_newt_its > 3000 && !step_failed) {
1214  if (max_Newton_steps != target_num_iter) {
1215  double iter_diff = num_newt_its - target_num_iter;
1216  double iter_adjust_zone = max_Newton_steps - target_num_iter;
1217  double target_time_step = delta_t_n
1218  *(1.0 - iter_diff*fabs(iter_diff)/
1219  ((2.0*iter_adjust_zone*iter_adjust_zone)));
1220  target_time_step = std::max(0.5*delta_t_n, target_time_step);
1221  if (target_time_step < delta_t_np1) {
1222  printf("\tNext time step will be decreased from %g to %g"
1223  " because of new its restraint\n",
1224  delta_t_np1, target_time_step);
1225  delta_t_np1 = target_time_step;
1226  }
1227  }
1228  }
1229 
1230 
1231  }
1232 
1233  /*
1234  * The final loop in the time stepping algorithm depends on whether the
1235  * current step was a success or not.
1236  */
1237  if (step_failed) {
1238  /*
1239  * Increment the counter indicating the number of consecutive
1240  * failures
1241  */
1242  numTSFailures++;
1243  /*
1244  * Print out a statement about the failure of the time step.
1245  */
1246  if (m_print_flag > 1) {
1247  print_time_fail(convFailure, m_time_step_num, time_n, delta_t_n,
1248  delta_t_np1, time_error_factor);
1249  } else if (m_print_flag == 1) {
1250  print_lvl1_summary(m_time_step_num, time_n, rslt, delta_t_n,
1251  num_newt_its, aztec_its, bktr_stps,
1252  time_error_factor,
1253  comment.c_str());
1254  }
1255 
1256  /*
1257  * Change time step counters back to the previous step before
1258  * the failed
1259  * time step occurred.
1260  */
1261  time_n -= delta_t_n;
1262  delta_t_n = delta_t_nm1;
1263  delta_t_nm1 = delta_t_nm2;
1264 
1265  /*
1266  * Replace old solution vector and time derivative solution vector.
1267  */
1268  m_y_n = m_y_nm1;
1269  m_ydot_n = m_ydot_nm1;
1270  /*
1271  * Decide whether to bail on the whole loop
1272  */
1273  if (numTSFailures > 35) {
1274  giveUp = true;
1275  }
1276  }
1277 
1278  /*
1279  * Do processing for a successful step.
1280  */
1281  else {
1282 
1283  /*
1284  * Decrement the number of consequative failure counter.
1285  */
1286  m_failure_counter = std::max(0, m_failure_counter-1);
1287 
1288  /*
1289  * Print out final results of a successfull time step.
1290  */
1291  if (m_print_flag > 1) {
1292  print_time_step2(m_time_step_num, m_order, time_n, time_error_factor,
1293  delta_t_n, delta_t_np1);
1294  } else if (m_print_flag == 1) {
1295  print_lvl1_summary(m_time_step_num, time_n, " ", delta_t_n,
1296  num_newt_its, aztec_its, bktr_stps, time_error_factor,
1297  comment.c_str());
1298  }
1299 
1300  /*
1301  * Output information at the end of every successful time step, if
1302  * requested.
1303  *
1304  * fill in
1305  */
1306 
1307 
1308  }
1309  } while (step_failed && !giveUp);
1310 
1311  /*
1312  * Send back the overall result of the time step.
1313  */
1314  if (step_failed) {
1315  if (time_n == 0.0) {
1316  return -1234.0;
1317  }
1318  return -time_n;
1319  }
1320  return time_n;
1321 }
1322 
1323 //-----------------------------------------------------------
1324 // Constants
1325 //-----------------------------------------------------------
1326 
1327 const double DampFactor = 4;
1328 const int NDAMP = 10;
1329 
1330 //-----------------------------------------------------------
1331 // MultiNewton methods
1332 //-----------------------------------------------------------
1333 
1334 double BEulerInt::soln_error_norm(const double* const delta_y,
1335  bool printLargest)
1336 {
1337  int i;
1338  double sum_norm = 0.0, error;
1339  for (i = 0; i < m_neq; i++) {
1340  error = delta_y[i] / m_ewt[i];
1341  sum_norm += (error * error);
1342  }
1343  sum_norm = sqrt(sum_norm / m_neq);
1344  if (printLargest) {
1345  const int num_entries = 8;
1346  double dmax1, normContrib;
1347  int j;
1348  vector_int imax(num_entries, -1);
1349  printf("\t\tPrintout of Largest Contributors to norm "
1350  "of value (%g)\n", sum_norm);
1351  printf("\t\t I ysoln deltaY weightY "
1352  "Error_Norm**2\n");
1353  printf("\t\t ");
1354  print_line("-", 80);
1355  for (int jnum = 0; jnum < num_entries; jnum++) {
1356  dmax1 = -1.0;
1357  for (i = 0; i < m_neq; i++) {
1358  bool used = false;
1359  for (j = 0; j < jnum; j++) {
1360  if (imax[j] == i) {
1361  used = true;
1362  }
1363  }
1364  if (!used) {
1365  error = delta_y[i] / m_ewt[i];
1366  normContrib = sqrt(error * error);
1367  if (normContrib > dmax1) {
1368  imax[jnum] = i;
1369  dmax1 = normContrib;
1370  }
1371  }
1372  }
1373  i = imax[jnum];
1374  if (i >= 0) {
1375  printf("\t\t %4d %12.4e %12.4e %12.4e %12.4e\n",
1376  i, m_y_n[i], delta_y[i], m_ewt[i], dmax1);
1377  }
1378  }
1379  printf("\t\t ");
1380  print_line("-", 80);
1381  }
1382  return sum_norm;
1383 }
1384 #ifdef DEBUG_HKM_JAC
1385 SquareMatrix jacBack();
1386 #endif
1387 
1388 void BEulerInt::doNewtonSolve(double time_curr, double* y_curr,
1389  double* ydot_curr, double* delta_y,
1390  GeneralMatrix& jac, int loglevel)
1391 {
1392  int irow, jcol;
1393 
1394  m_func->evalResidNJ(time_curr, delta_t_n, y_curr,
1395  ydot_curr, delta_y, Base_ResidEval);
1396  m_nfe++;
1397  int sz = m_func->nEquations();
1398  for (int n = 0; n < sz; n++) {
1399  delta_y[n] = -delta_y[n];
1400  }
1401 
1402 
1403  /*
1404  * Column scaling -> We scale the columns of the Jacobian
1405  * by the nominal important change in the solution vector
1406  */
1407  if (m_colScaling) {
1408  if (!jac.factored()) {
1409  /*
1410  * Go get new scales
1411  */
1412  setColumnScales();
1413 
1414  /*
1415  * Scale the new Jacobian
1416  */
1417  double* jptr = &(*(jac.begin()));
1418  for (jcol = 0; jcol < m_neq; jcol++) {
1419  for (irow = 0; irow < m_neq; irow++) {
1420  *jptr *= m_colScales[jcol];
1421  jptr++;
1422  }
1423  }
1424  }
1425  }
1426 
1427  if (m_matrixConditioning) {
1428  if (jac.factored()) {
1429  m_func->matrixConditioning(0, sz, delta_y);
1430  } else {
1431  double* jptr = &(*(jac.begin()));
1432  m_func->matrixConditioning(jptr, sz, delta_y);
1433  }
1434  }
1435 
1436  /*
1437  * row sum scaling -> Note, this is an unequivocal success
1438  * at keeping the small numbers well balanced and
1439  * nonnegative.
1440  */
1441  if (m_rowScaling) {
1442  if (! jac.factored()) {
1443  /*
1444  * Ok, this is ugly. jac.begin() returns an vector<double> iterator
1445  * to the first data location.
1446  * Then &(*()) reverts it to a double *.
1447  */
1448  double* jptr = &(*(jac.begin()));
1449  for (irow = 0; irow < m_neq; irow++) {
1450  m_rowScales[irow] = 0.0;
1451  }
1452  for (jcol = 0; jcol < m_neq; jcol++) {
1453  for (irow = 0; irow < m_neq; irow++) {
1454  m_rowScales[irow] += fabs(*jptr);
1455  jptr++;
1456  }
1457  }
1458 
1459  jptr = &(*(jac.begin()));
1460  for (jcol = 0; jcol < m_neq; jcol++) {
1461  for (irow = 0; irow < m_neq; irow++) {
1462  *jptr /= m_rowScales[irow];
1463  jptr++;
1464  }
1465  }
1466  }
1467  for (irow = 0; irow < m_neq; irow++) {
1468  delta_y[irow] /= m_rowScales[irow];
1469  }
1470  }
1471 
1472 #ifdef DEBUG_HKM_JAC
1473  bool printJacContributions = false;
1474  if (m_time_step_num > 304) {
1475  printJacContributions = false;
1476  }
1477  int focusRow = 10;
1478  int numRows = 2;
1479  double RRow[2];
1480  bool freshJac = true;
1481  RRow[0] = delta_y[focusRow];
1482  RRow[1] = delta_y[focusRow+1];
1483  double Pcutoff = 1.0E-70;
1484  if (!jac.m_factored) {
1485  jacBack = jac;
1486  } else {
1487  freshJac = false;
1488  }
1489 #endif
1490  /*
1491  * Solve the system -> This also involves inverting the
1492  * matrix
1493  */
1494  (void) jac.solve(delta_y);
1495 
1496 
1497  /*
1498  * reverse the column scaling if there was any.
1499  */
1500  if (m_colScaling) {
1501  for (irow = 0; irow < m_neq; irow++) {
1502  delta_y[irow] *= m_colScales[irow];
1503  }
1504  }
1505 
1506 #ifdef DEBUG_HKM_JAC
1507  if (printJacContributions) {
1508  for (int iNum = 0; iNum < numRows; iNum++) {
1509  if (iNum > 0) {
1510  focusRow++;
1511  }
1512  double dsum = 0.0;
1513  vector_fp& Jdata = jacBack.data();
1514  double dRow = Jdata[m_neq * focusRow + focusRow];
1515  printf("\n Details on delta_Y for row %d \n", focusRow);
1516  printf(" Value before = %15.5e, delta = %15.5e,"
1517  "value after = %15.5e\n", y_curr[focusRow],
1518  delta_y[focusRow],
1519  y_curr[focusRow] + delta_y[focusRow]);
1520  if (!freshJac) {
1521  printf(" Old Jacobian\n");
1522  }
1523  printf(" col delta_y aij "
1524  "contrib \n");
1525  printf("--------------------------------------------------"
1526  "---------------------------------------------\n");
1527  printf(" Res(%d) %15.5e %15.5e %15.5e (Res = %g)\n",
1528  focusRow, delta_y[focusRow],
1529  dRow, RRow[iNum] / dRow, RRow[iNum]);
1530  dsum += RRow[iNum] / dRow;
1531  for (int ii = 0; ii < m_neq; ii++) {
1532  if (ii != focusRow) {
1533  double aij = Jdata[m_neq * ii + focusRow];
1534  double contrib = aij * delta_y[ii] * (-1.0) / dRow;
1535  dsum += contrib;
1536  if (fabs(contrib) > Pcutoff) {
1537  printf("%6d %15.5e %15.5e %15.5e\n", ii,
1538  delta_y[ii] , aij, contrib);
1539  }
1540  }
1541  }
1542  printf("--------------------------------------------------"
1543  "---------------------------------------------\n");
1544  printf(" %15.5e %15.5e\n",
1545  delta_y[focusRow], dsum);
1546  }
1547  }
1548 
1549 #endif
1550 
1552 }
1553 
1554 double BEulerInt::boundStep(const double* const y,
1555  const double* const step0, int loglevel)
1556 {
1557  int i, i_lower = -1, ifbd = 0, i_fbd = 0;
1558  double fbound = 1.0, f_lowbounds = 1.0, f_delta_bounds = 1.0;
1559  double ff, y_new, ff_alt;
1560  for (i = 0; i < m_neq; i++) {
1561  y_new = y[i] + step0[i];
1562  if ((y_new < (-0.01 * m_ewt[i])) && y[i] >= 0.0) {
1563  ff = 0.9 * (y[i] / (y[i] - y_new));
1564  if (ff < f_lowbounds) {
1565  f_lowbounds = ff;
1566  i_lower = i;
1567  }
1568  }
1569  /*
1570  * Now do a delta bounds
1571  * Increase variables by a factor of 2 only
1572  * decrease variables by a factor of 5 only
1573  */
1574  ff = 1.0;
1575  if ((fabs(y_new) > 2.0 * fabs(y[i])) &&
1576  (fabs(y_new-y[i]) > m_ewt[i])) {
1577  ff = fabs(y[i]/(y_new - y[i]));
1578  ff_alt = fabs(m_ewt[i] / (y_new - y[i]));
1579  ff = std::max(ff, ff_alt);
1580  ifbd = 1;
1581  }
1582  if ((fabs(5.0 * y_new) < fabs(y[i])) &&
1583  (fabs(y_new - y[i]) > m_ewt[i])) {
1584  ff = y[i]/(y_new-y[i]) * (1.0 - 5.0)/5.0;
1585  ff_alt = fabs(m_ewt[i] / (y_new - y[i]));
1586  ff = std::max(ff, ff_alt);
1587  ifbd = 0;
1588  }
1589  if (ff < f_delta_bounds) {
1590  f_delta_bounds = ff;
1591  i_fbd = ifbd;
1592  }
1593  f_delta_bounds = std::min(f_delta_bounds, ff);
1594  }
1595  fbound = std::min(f_lowbounds, f_delta_bounds);
1596  /*
1597  * Report on any corrections
1598  */
1599  if (loglevel > 1) {
1600  if (fbound != 1.0) {
1601  if (f_lowbounds < f_delta_bounds) {
1602  printf("\t\tboundStep: Variable %d causing lower bounds "
1603  "damping of %g\n",
1604  i_lower, f_lowbounds);
1605  } else {
1606  if (ifbd) {
1607  printf("\t\tboundStep: Decrease of Variable %d causing "
1608  "delta damping of %g\n",
1609  i_fbd, f_delta_bounds);
1610  } else {
1611  printf("\t\tboundStep: Increase of variable %d causing"
1612  "delta damping of %g\n",
1613  i_fbd, f_delta_bounds);
1614  }
1615  }
1616  }
1617  }
1618  return fbound;
1619 }
1620 
1621 int BEulerInt::dampStep(double time_curr, const double* y0,
1622  const double* ydot0, const double* step0,
1623  double* y1, double* ydot1, double* step1,
1624  double& s1, GeneralMatrix& jac,
1625  int& loglevel, bool writetitle,
1626  int& num_backtracks)
1627 {
1628  // Compute the weighted norm of the undamped step size step0
1629  double s0 = soln_error_norm(step0);
1630 
1631  // Compute the multiplier to keep all components in bounds
1632  // A value of one indicates that there is no limitation
1633  // on the current step size in the nonlinear method due to
1634  // bounds constraints (either negative values of delta
1635  // bounds constraints.
1636  double fbound = boundStep(y0, step0, loglevel);
1637 
1638  // if fbound is very small, then y0 is already close to the
1639  // boundary and step0 points out of the allowed domain. In
1640  // this case, the Newton algorithm fails, so return an error
1641  // condition.
1642  if (fbound < 1.e-10) {
1643  if (loglevel > 1) {
1644  printf("\t\t\tdampStep: At limits.\n");
1645  }
1646  return -3;
1647  }
1648 
1649  //--------------------------------------------
1650  // Attempt damped step
1651  //--------------------------------------------
1652 
1653  // damping coefficient starts at 1.0
1654  double damp = 1.0;
1655  int j, m;
1656  double ff;
1657  num_backtracks = 0;
1658  for (m = 0; m < NDAMP; m++) {
1659 
1660  ff = fbound*damp;
1661 
1662  // step the solution by the damped step size
1663  /*
1664  * Whenever we update the solution, we must also always
1665  * update the time derivative.
1666  */
1667  for (j = 0; j < m_neq; j++) {
1668  y1[j] = y0[j] + ff*step0[j];
1669  // HKM setting intermediate y's to zero was a tossup.
1670  // slightly different, equivalent results
1671  //#ifdef DEBUG_HKM
1672  // y1[j] = MAX(0.0, y1[j]);
1673  //#endif
1674  }
1675  calc_ydot(m_order, y1, ydot1);
1676 
1677  // compute the next undamped step, step1[], that would result
1678  // if y1[] were accepted.
1679 
1680  doNewtonSolve(time_curr, y1, ydot1, step1, jac, loglevel);
1681 
1682 #ifdef DEBUG_HKM
1683  for (j = 0; j < m_neq; j++) {
1684  checkFinite(step1[j]);
1685  checkFinite(y1[j]);
1686  }
1687 #endif
1688  // compute the weighted norm of step1
1689  s1 = soln_error_norm(step1);
1690 
1691  // write log information
1692  if (loglevel > 3) {
1693  print_solnDelta_norm_contrib((const double*) step0,
1694  "DeltaSolnTrial",
1695  (const double*) step1,
1696  "DeltaSolnTrialTest",
1697  "dampNewt: Important Entries for "
1698  "Weighted Soln Updates:",
1699  y0, y1, ff, 5);
1700  }
1701  if (loglevel > 1) {
1702  printf("\t\t\tdampNewt: s0 = %g, s1 = %g, fbound = %g,"
1703  "damp = %g\n", s0, s1, fbound, damp);
1704  }
1705 #ifdef DEBUG_HKM
1706  if (loglevel > 2) {
1707  if (s1 > 1.00000001 * s0 && s1 > 1.0E-5) {
1708  printf("WARNING: Possible Jacobian Problem "
1709  "-> turning on more debugging for this step!!!\n");
1710  print_solnDelta_norm_contrib((const double*) step0,
1711  "DeltaSolnTrial",
1712  (const double*) step1,
1713  "DeltaSolnTrialTest",
1714  "dampNewt: Important Entries for "
1715  "Weighted Soln Updates:",
1716  y0, y1, ff, 5);
1717  loglevel = 4;
1718  }
1719  }
1720 #endif
1721 
1722  // if the norm of s1 is less than the norm of s0, then
1723  // accept this damping coefficient. Also accept it if this
1724  // step would result in a converged solution. Otherwise,
1725  // decrease the damping coefficient and try again.
1726 
1727  if (s1 < 1.0E-5 || s1 < s0) {
1728  if (loglevel > 2) {
1729  if (s1 > s0) {
1730  if (s1 > 1.0) {
1731  printf("\t\t\tdampStep: current trial step and damping"
1732  " coefficient accepted because test step < 1\n");
1733  printf("\t\t\t s1 = %g, s0 = %g\n", s1, s0);
1734  }
1735  }
1736  }
1737  break;
1738  } else {
1739  if (loglevel > 1) {
1740  printf("\t\t\tdampStep: current step rejected: (s1 = %g > "
1741  "s0 = %g)", s1, s0);
1742  if (m < (NDAMP-1)) {
1743  printf(" Decreasing damping factor and retrying");
1744  } else {
1745  printf(" Giving up!!!");
1746  }
1747  printf("\n");
1748  }
1749  }
1750  num_backtracks++;
1751  damp /= DampFactor;
1752  }
1753 
1754  // If a damping coefficient was found, return 1 if the
1755  // solution after stepping by the damped step would represent
1756  // a converged solution, and return 0 otherwise. If no damping
1757  // coefficient could be found, return -2.
1758  if (m < NDAMP) {
1759  if (s1 > 1.0) {
1760  return 0;
1761  } else {
1762  return 1;
1763  }
1764  } else {
1765  if (s1 < 0.5 && (s0 < 0.5)) {
1766  return 1;
1767  }
1768  if (s1 < 1.0) {
1769  return 0;
1770  }
1771  return -2;
1772  }
1773 }
1774 
1775 int BEulerInt::solve_nonlinear_problem(double* const y_comm,
1776  double* const ydot_comm, double CJ,
1777  double time_curr,
1778  GeneralMatrix& jac,
1779  int& num_newt_its,
1780  int& num_linear_solves,
1781  int& num_backtracks,
1782  int loglevel)
1783 {
1784  int m = 0;
1785  bool forceNewJac = false;
1786  double s1=1.e30;
1787 
1788  vector_fp y_curr(y_comm, y_comm + m_neq);
1789  vector_fp ydot_curr(ydot_comm, ydot_comm + m_neq);
1790  vector_fp stp(m_neq, 0.0);
1791  vector_fp stp1(m_neq, 0.0);
1792  vector_fp y_new(m_neq, 0.0);
1793  vector_fp ydot_new(m_neq, 0.0);
1794 
1795  bool frst = true;
1796  num_newt_its = 0;
1797  num_linear_solves = - m_numTotalLinearSolves;
1798  num_backtracks = 0;
1799  int i_backtracks;
1800 
1801  while (1 > 0) {
1802 
1803  /*
1804  * Increment Newton Solve counter
1805  */
1807  num_newt_its++;
1808 
1809 
1810  if (loglevel > 1) {
1811  printf("\t\tSolve_Nonlinear_Problem: iteration %d:\n",
1812  num_newt_its);
1813  }
1814 
1815  // Check whether the Jacobian should be re-evaluated.
1816 
1817  forceNewJac = true;
1818 
1819  if (forceNewJac) {
1820  if (loglevel > 1) {
1821  printf("\t\t\tGetting a new Jacobian and solving system\n");
1822  }
1823  beuler_jac(jac, &m_resid[0], time_curr, CJ, &y_curr[0], &ydot_curr[0],
1824  num_newt_its);
1825  } else {
1826  if (loglevel > 1) {
1827  printf("\t\t\tSolving system with old jacobian\n");
1828  }
1829  }
1830 
1831  // compute the undamped Newton step
1832  doNewtonSolve(time_curr, &y_curr[0], &ydot_curr[0], &stp[0], jac, loglevel);
1833 
1834  // damp the Newton step
1835  m = dampStep(time_curr, &y_curr[0], &ydot_curr[0], &stp[0], &y_new[0], &ydot_new[0],
1836  &stp1[0], s1, jac, loglevel, frst, i_backtracks);
1837  frst = false;
1838  num_backtracks += i_backtracks;
1839 
1840  /*
1841  * Impose the minimum number of newton iterations critera
1842  */
1843  if (num_newt_its < m_min_newt_its) {
1844  if (m == 1) {
1845  m = 0;
1846  }
1847  }
1848  /*
1849  * Impose max newton iteration
1850  */
1851  if (num_newt_its > 20) {
1852  m = -1;
1853  if (loglevel > 1) {
1854  printf("\t\t\tDampnewton unsuccessful (max newts exceeded) sfinal = %g\n", s1);
1855  }
1856  }
1857 
1858  if (loglevel > 1) {
1859  if (m == 1) {
1860  printf("\t\t\tDampNewton iteration successful, nonlin "
1861  "converged sfinal = %g\n", s1);
1862  } else if (m == 0) {
1863  printf("\t\t\tDampNewton iteration successful, get new"
1864  "direction, sfinal = %g\n", s1);
1865  } else {
1866  printf("\t\t\tDampnewton unsuccessful sfinal = %g\n", s1);
1867  }
1868  }
1869 
1870  // If we are converged, then let's use the best solution possible
1871  // for an end result. We did a resolve in dampStep(). Let's update
1872  // the solution to reflect that.
1873  // HKM 5/16 -> Took this out, since if the last step was a
1874  // damped step, then adding stp1[j] is undamped, and
1875  // may lead to oscillations. It kind of defeats the
1876  // purpose of dampStep() anyway.
1877  // if (m == 1) {
1878  // for (int j = 0; j < m_neq; j++) {
1879  // y_new[j] += stp1[j];
1880  // HKM setting intermediate y's to zero was a tossup.
1881  // slightly different, equivalent results
1882  // #ifdef DEBUG_HKM
1883  // y_new[j] = MAX(0.0, y_new[j]);
1884  // #endif
1885  // }
1886  // }
1887 
1888  bool m_filterIntermediate = false;
1889  if (m_filterIntermediate) {
1890  if (m == 0) {
1891  (void) filterNewStep(time_n, &y_new[0], &ydot_new[0]);
1892  }
1893  }
1894  // Exchange new for curr solutions
1895  if (m == 0 || m == 1) {
1896  y_curr = y_new;
1897  calc_ydot(m_order, &y_curr[0], &ydot_curr[0]);
1898  }
1899 
1900  // convergence
1901  if (m == 1) {
1902  goto done;
1903  }
1904 
1905  // If dampStep fails, first try a new Jacobian if an old
1906  // one was being used. If it was a new Jacobian, then
1907  // return -1 to signify failure.
1908  else if (m < 0) {
1909  goto done;
1910  }
1911  }
1912 
1913 done:
1914  // Copy into the return vectors
1915  copy(y_curr.begin(), y_curr.end(), y_comm);
1916  copy(ydot_curr.begin(), ydot_curr.end(), ydot_comm);
1917  // Increment counters
1918  num_linear_solves += m_numTotalLinearSolves;
1919 
1920  double time_elapsed = 0.0;
1921  if (loglevel > 1) {
1922  if (m == 1) {
1923  printf("\t\tNonlinear problem solved successfully in "
1924  "%d its, time elapsed = %g sec\n",
1925  num_newt_its, time_elapsed);
1926  }
1927  }
1928  return m;
1929 }
1930 
1931 void BEulerInt::
1932 print_solnDelta_norm_contrib(const double* const solnDelta0,
1933  const char* const s0,
1934  const double* const solnDelta1,
1935  const char* const s1,
1936  const char* const title,
1937  const double* const y0,
1938  const double* const y1,
1939  double damp,
1940  int num_entries)
1941 {
1942  int i, j, jnum;
1943  bool used;
1944  double dmax0, dmax1, error, rel_norm;
1945  printf("\t\t%s currentDamp = %g\n", title, damp);
1946  printf("\t\t I ysoln %10s ysolnTrial "
1947  "%10s weight relSoln0 relSoln1\n", s0, s1);
1948  vector_int imax(num_entries, -1);
1949  printf("\t\t ");
1950  print_line("-", 90);
1951  for (jnum = 0; jnum < num_entries; jnum++) {
1952  dmax1 = -1.0;
1953  for (i = 0; i < m_neq; i++) {
1954  used = false;
1955  for (j = 0; j < jnum; j++) {
1956  if (imax[j] == i) {
1957  used = true;
1958  }
1959  }
1960  if (!used) {
1961  error = solnDelta0[i] / m_ewt[i];
1962  rel_norm = sqrt(error * error);
1963  error = solnDelta1[i] / m_ewt[i];
1964  rel_norm += sqrt(error * error);
1965  if (rel_norm > dmax1) {
1966  imax[jnum] = i;
1967  dmax1 = rel_norm;
1968  }
1969  }
1970  }
1971  if (imax[jnum] >= 0) {
1972  i = imax[jnum];
1973  error = solnDelta0[i] / m_ewt[i];
1974  dmax0 = sqrt(error * error);
1975  error = solnDelta1[i] / m_ewt[i];
1976  dmax1 = sqrt(error * error);
1977  printf("\t\t %4d %12.4e %12.4e %12.4e %12.4e "
1978  "%12.4e %12.4e %12.4e\n",
1979  i, y0[i], solnDelta0[i], y1[i],
1980  solnDelta1[i], m_ewt[i], dmax0, dmax1);
1981  }
1982  }
1983  printf("\t\t ");
1984  print_line("-", 90);
1985 }
1986 
1987 } // End of namespace Cantera
virtual void setPrintSolnOptions(int printSolnStepInterval, int printSolnNumberToTout, int printSolnFirstSteps=0, bool dumpJacobians=false)
This routine controls when the solution is printed.
Definition: BEulerInt.cpp:134
double m_time_final
Final time.
Definition: BEulerInt.h:521
virtual doublereal * ptrColumn(size_t j)=0
Return a pointer to the top of column j, columns are assumed to be contiguous in memory.
virtual void user_out2(const int ifunc, const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot)
This function may be used to create output at various points in the execution of an application...
int m_numTotalTruncFails
Total Number of time truncation error failures.
Definition: BEulerInt.h:577
double time_error_norm()
Definition: BEulerInt.cpp:678
virtual void setProblemType(int probtype)
Set the problem type.
Definition: BEulerInt.cpp:109
int m_printSolnNumberToTout
Number of evenly spaced printouts of the solution If zero, there is no printout from this option defa...
Definition: BEulerInt.h:494
double getPrintTime(double time_current)
Get the next time to print out.
Definition: BEulerInt.cpp:220
virtual void setIterator(IterType t)
Set the linear iterator.
Definition: BEulerInt.cpp:145
int m_printSolnStepInterval
Step Interval at which to print out the solution default = 1; If set to zero, there is no printout...
Definition: BEulerInt.h:487
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
Definition: checkFinite.cpp:33
virtual int nEquations() const
Return the number of equations in the equation system.
int m_numTotalNewtIts
Number of total newton iterations.
Definition: BEulerInt.h:568
int m_neq
Number of equations in the ode integrator.
Definition: BEulerInt.h:507
virtual double soln_error_norm(const double *const, bool printLargest=false)
Calculate the solution error norm.
Definition: BEulerInt.cpp:1334
void doNewtonSolve(double, double *, double *, double *, GeneralMatrix &, int)
Compute the undamped Newton step.
Definition: BEulerInt.cpp:1388
int m_failure_counter
Failure Counter -> keeps track of the number of consecutive failures.
Definition: BEulerInt.h:475
BEulerMethodType m_method
MethodType is used to specify how the time step is to be chosen.
Definition: BEulerInt.h:407
bool m_rowScaling
m_rowScaling is a boolean.
Definition: BEulerInt.h:417
virtual int evalResidNJ(const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot, doublereal *const resid, const ResidEval_Type_Enum evalType=Base_ResidEval, const int id_x=-1, const doublereal delta_x=0.0)
Evaluate the residual function.
virtual int evalJacobian(const doublereal t, const doublereal delta_t, doublereal cj, const doublereal *const y, const doublereal *const ydot, GeneralMatrix &J, doublereal *const resid)
Calculate an analytical jacobian and the residual at the current time and values. ...
virtual vector_fp::iterator begin()=0
Return an iterator pointing to the first element.
void beuler_jac(GeneralMatrix &J, double *const f, double, double, double *const, double *const, int)
Definition: BEulerInt.cpp:515
double m_reltol
Relative time truncation error tolerances.
Definition: BEulerInt.h:437
double m_t0
Initial time at the start of the integration.
Definition: BEulerInt.h:518
const double DampFactor
Dampfactor is the factor by which the damping factor is reduced by when a reduction in step length is...
Definition: BEulerInt.cpp:1327
virtual void writeSolution(int ievent, const double time, const double deltaT, const int time_step_num, const double *y, const double *ydot)
Write out to a file or to standard output the current solution.
Definition: ResidEval.h:138
Generic matrix.
Definition: GeneralMatrix.h:23
Base residual calculation for the Jacobian calculation.
Definition: ResidJacEval.h:28
virtual int getInitialConditions(const doublereal t0, doublereal *const y, doublereal *const ydot)
Fill in the initial conditions.
bool m_matrixConditioning
m_matrixConditioning is a boolean.
Definition: BEulerInt.h:429
int m_min_newt_its
Minimum Number of Newton Iterations per nonlinear step. default = 0.
Definition: BEulerInt.h:478
int m_print_flag
Determines the level of printing for each time step.
Definition: BEulerInt.h:553
double subtractRD(double a, double b)
This routine subtracts two numbers for one another.
Definition: BEulerInt.cpp:500
vector_fp m_abstol
Vector of absolute time truncation error tolerance when not uniform for all variables.
Definition: BEulerInt.h:448
IterType m_iter
IterType is used to specify how the nonlinear equations are to be relaxed at each time step...
Definition: BEulerInt.h:400
Wrappers for the function evaluators for Nonlinear solvers and Time steppers.
Definition: ResidJacEval.h:51
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:167
double delta_t_max
Maximum permissible time step.
Definition: BEulerInt.h:532
virtual int calcDeltaSolnVariables(const doublereal t, const doublereal *const y, const doublereal *const ydot, doublereal *const delta_y, const doublereal *const solnWeights=0)
Return a vector of delta y's for calculation of the numerical Jacobian.
int solve_nonlinear_problem(double *const y_comm, double *const ydot_comm, double CJ, double time_curr, GeneralMatrix &jac, int &num_newt_its, int &num_linear_solves, int &num_backtracks, int loglevel)
Solve a nonlinear system.
Definition: BEulerInt.cpp:1775
int m_numTotalLinearSolves
Total number of linear iterations.
Definition: BEulerInt.h:571
GeneralMatrix * tdjac_ptr
Pointer to the jacobian representing the time dependent problem.
Definition: BEulerInt.h:542
virtual void setSolnWeights()
Set the solution weights.
Definition: BEulerInt.cpp:260
virtual int matrixConditioning(doublereal *const matrix, const int nrows, doublereal *const rhs)
Multiply the matrix by another matrix that leads to better conditioning.
int m_printSolnFirstSteps
Number of initial steps that the solution is printed out. default = 0.
Definition: BEulerInt.h:497
virtual void setColumnScales()
Set the column scaling vector at the current time.
Definition: BEulerInt.cpp:282
int m_numInitialConstantDeltaTSteps
Number of initial time steps to take where the time truncation error tolerances are not checked...
Definition: BEulerInt.h:472
void calc_y_pred(int)
Definition: BEulerInt.cpp:628
void error(const std::string &msg)
Write an error message and quit.
Definition: global.cpp:71
int m_order
Current integration order.
Definition: BEulerInt.h:460
int m_nfe
Number of function evaluations.
Definition: BEulerInt.h:559
A class for full (non-sparse) matrices with Fortran-compatible data storage.
Definition: SquareMatrix.h:26
virtual void initializeRJE(double t0, ResidJacEval &func)
Find the initial conditions for y and ydot.
Definition: BEulerInt.cpp:175
virtual void setTolerances(double reltol, size_t n, double *abstol)
Set or reset the number of equations.
Definition: BEulerInt.cpp:88
virtual bool factored() const =0
true if the current factorization is up to date with the matrix
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:68
const int NDAMP
Number of damping steps that are carried out before the solution is deemed a failure.
Definition: BEulerInt.cpp:1328
double m_abstols
Absolute time truncation error tolerances, when uniform for all variables.
Definition: BEulerInt.h:443
virtual void clearFactorFlag()=0
clear the factored flag
double time_step_control(int m_order, double time_error_factor)
Definition: BEulerInt.cpp:732
virtual int evalTimeTrackingEqns(const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot)
Evaluate the time tracking equations, if any.
int dampStep(double, const double *, const double *, const double *, double *, double *, double *, double &, GeneralMatrix &, int &, bool, int &)
Definition: BEulerInt.cpp:1621
int m_max_time_step_attempts
Max time steps allowed.
Definition: BEulerInt.h:467
virtual void calcSolnScales(const doublereal t, const doublereal *const y, const doublereal *const y_old, doublereal *const yScales)
Returns a vector of column scale factors that can be used to column scale Jacobians.
double filterNewStep(double, double *, double *)
Filter a new step.
Definition: BEulerInt.cpp:311
int m_jacFormMethod
m_jacFormMethod determines how a matrix is formed.
Definition: BEulerInt.h:411
int m_numTotalConvFails
Total number of convergence failures.
Definition: BEulerInt.h:574
virtual int nEvals() const
Return the total number of function evaluations.
Definition: BEulerInt.cpp:235
double boundStep(const double *const y, const double *const step0, int loglevel)
Bound the Newton step while relaxing the solution.
Definition: BEulerInt.cpp:1554
bool m_colScaling
m_colScaling is a boolean.
Definition: BEulerInt.h:422
IterType
Specifies the method used for iteration.
Definition: Integrator.h:41
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
void internalMalloc()
Internal routine that sets up the fixed length storage based on the size of the problem to solve...
Definition: BEulerInt.cpp:240
bool m_dumpJacobians
Dump Jacobians to disk. default false.
Definition: BEulerInt.h:500
virtual int solve(doublereal *b)=0
Solves the Ax = b system returning x in the b spot.
virtual doublereal step(double tout)
Integrate the system of equations.
Definition: BEulerInt.cpp:927
Delta residual calculation for the Jacbobian calculation.
Definition: ResidJacEval.h:30
Newton Iteration.
Definition: Integrator.h:43
int m_nJacEval
Number of Jacobian Evaluations and factorization steps (they are the same)
Definition: BEulerInt.h:565
virtual doublereal filterSolnPrediction(const doublereal t, doublereal *const y)
Filter the solution predictions.
void setNonLinOptions(int min_newt_its=0, bool matrixConditioning=false, bool colScaling=false, bool rowScaling=true)
Set the options for the nonlinear method.
Definition: BEulerInt.cpp:150
void calc_ydot(int, double *, double *)
Definition: BEulerInt.cpp:657
double m_hmax
Maximum step size.
Definition: BEulerInt.h:454
int m_itol
If m_itol =1 then each component has an individual value of atol.
Definition: BEulerInt.h:434
virtual bool evalStoppingCritera(const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot)
Evaluate any stopping criteria other than a final time limit.
Base residual calculation for the time-stepping function.
Definition: ResidJacEval.h:26
vector_fp m_ewt
Error Weights. This is a surprisingly important quantity.
Definition: BEulerInt.h:451
int m_time_step_num
Time step number.
Definition: BEulerInt.h:463
void computeResidWts(GeneralMatrix &jac)
Compute Residual Weights.
Definition: BEulerInt.cpp:287