Cantera  3.1.0a1
CVodesIntegrator.cpp
Go to the documentation of this file.
1 //! @file CVodesIntegrator.cpp
2 
3 // This file is part of Cantera. See License.txt in the top-level directory or
4 // at https://cantera.org/license.txt for license and copyright information.
5 
8 
9 #include <iostream>
10 using namespace std;
11 
12 #include "cantera/numerics/sundials_headers.h"
13 
14 namespace {
15 
16 N_Vector newNVector(size_t N, Cantera::SundialsContext& context)
17 {
18 #if CT_SUNDIALS_VERSION >= 60
19  return N_VNew_Serial(static_cast<sd_size_t>(N), context.get());
20 #else
21  return N_VNew_Serial(static_cast<sd_size_t>(N));
22 #endif
23 }
24 
25 } // end anonymous namespace
26 
27 namespace Cantera
28 {
29 
30 extern "C" {
31  /**
32  * Function called by cvodes to evaluate ydot given y. The CVODE integrator
33  * allows passing in a void* pointer to access external data. This pointer
34  * is cast to a pointer to a instance of class FuncEval. The equations to be
35  * integrated should be specified by deriving a class from FuncEval that
36  * evaluates the desired equations.
37  * @ingroup odeGroup
38  */
39  static int cvodes_rhs(realtype t, N_Vector y, N_Vector ydot, void* f_data)
40  {
41  FuncEval* f = (FuncEval*) f_data;
42  return f->evalNoThrow(t, NV_DATA_S(y), NV_DATA_S(ydot));
43  }
44 
45  //! Function called by CVodes when an error is encountered instead of
46  //! writing to stdout. Here, save the error message provided by CVodes so
47  //! that it can be included in the subsequently raised CanteraError.
48  static void cvodes_err(int error_code, const char* module,
49  const char* function, char* msg, void* eh_data)
50  {
51  CVodesIntegrator* integrator = (CVodesIntegrator*) eh_data;
52  integrator->m_error_message = msg;
53  integrator->m_error_message += "\n";
54  }
55 
56  static int cvodes_prec_setup(realtype t, N_Vector y, N_Vector ydot, booleantype jok,
57  booleantype *jcurPtr, realtype gamma, void *f_data)
58  {
59  FuncEval* f = (FuncEval*) f_data;
60  if (!jok) {
61  *jcurPtr = true; // jacobian data was recomputed
62  return f->preconditioner_setup_nothrow(t, NV_DATA_S(y), gamma);
63  } else {
64  f->updatePreconditioner(gamma); // updates preconditioner with new gamma
65  *jcurPtr = false; // indicates that Jacobian data was not recomputed
66  return 0; // no error because not recomputed
67  }
68  }
69 
70  static int cvodes_prec_solve(realtype t, N_Vector y, N_Vector ydot, N_Vector r,
71  N_Vector z, realtype gamma, realtype delta, int lr,
72  void* f_data)
73  {
74  FuncEval* f = (FuncEval*) f_data;
75  return f->preconditioner_solve_nothrow(NV_DATA_S(r),NV_DATA_S(z));
76  }
77 }
78 
79 CVodesIntegrator::CVodesIntegrator()
80  : m_itol(CV_SS)
81  , m_method(CV_BDF)
82 {
83 }
84 
85 CVodesIntegrator::~CVodesIntegrator()
86 {
87  if (m_cvode_mem) {
88  if (m_np > 0) {
89  CVodeSensFree(m_cvode_mem);
90  }
91  CVodeFree(&m_cvode_mem);
92  }
93 
94  SUNLinSolFree((SUNLinearSolver) m_linsol);
95  SUNMatDestroy((SUNMatrix) m_linsol_matrix);
96 
97  if (m_y) {
98  N_VDestroy_Serial(m_y);
99  }
100  if (m_abstol) {
101  N_VDestroy_Serial(m_abstol);
102  }
103  if (m_dky) {
104  N_VDestroy_Serial(m_dky);
105  }
106  if (m_yS) {
107  #if CT_SUNDIALS_VERSION >= 60
108  N_VDestroyVectorArray(m_yS, static_cast<int>(m_np));
109  #else
110  N_VDestroyVectorArray_Serial(m_yS, static_cast<int>(m_np));
111  #endif
112  }
113 }
114 
115 
116 double& CVodesIntegrator::solution(size_t k)
117 {
118  return NV_Ith_S(m_y, k);
119 }
120 
122 {
123  return NV_DATA_S(m_y);
124 }
125 
126 void CVodesIntegrator::setTolerances(double reltol, size_t n, double* abstol)
127 {
128  m_itol = CV_SV;
129  m_nabs = n;
130  if (n != m_neq) {
131  if (m_abstol) {
132  N_VDestroy_Serial(m_abstol);
133  }
134  m_abstol = newNVector(n, m_sundials_ctx);
135  }
136  for (size_t i=0; i<n; i++) {
137  NV_Ith_S(m_abstol, i) = abstol[i];
138  }
139  m_reltol = reltol;
140 }
141 
142 void CVodesIntegrator::setTolerances(double reltol, double abstol)
143 {
144  m_itol = CV_SS;
145  m_reltol = reltol;
146  m_abstols = abstol;
147 }
148 
149 void CVodesIntegrator::setSensitivityTolerances(double reltol, double abstol)
150 {
151  m_reltolsens = reltol;
152  m_abstolsens = abstol;
153 }
154 
156 {
157  if (t == BDF_Method) {
158  m_method = CV_BDF;
159  } else if (t == Adams_Method) {
160  m_method = CV_ADAMS;
161  } else {
162  throw CanteraError("CVodesIntegrator::setMethod", "unknown method");
163  }
164 }
165 
167 {
168  m_hmax = hmax;
169  if (m_cvode_mem) {
170  CVodeSetMaxStep(m_cvode_mem, hmax);
171  }
172 }
173 
175 {
176  m_hmin = hmin;
177  if (m_cvode_mem) {
178  CVodeSetMinStep(m_cvode_mem, hmin);
179  }
180 }
181 
183 {
184  m_maxsteps = nmax;
185  if (m_cvode_mem) {
186  CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
187  }
188 }
189 
191 {
192  return m_maxsteps;
193 }
194 
196 {
197  m_maxErrTestFails = n;
198  if (m_cvode_mem) {
199  CVodeSetMaxErrTestFails(m_cvode_mem, n);
200  }
201 }
202 
203 void CVodesIntegrator::sensInit(double t0, FuncEval& func)
204 {
205  m_np = func.nparams();
206  m_sens_ok = false;
207 
208  N_Vector y = newNVector(func.neq(), m_sundials_ctx);
209  #if CT_SUNDIALS_VERSION >= 60
210  m_yS = N_VCloneVectorArray(static_cast<int>(m_np), y);
211  #else
212  m_yS = N_VCloneVectorArray_Serial(static_cast<int>(m_np), y);
213  #endif
214  for (size_t n = 0; n < m_np; n++) {
215  N_VConst(0.0, m_yS[n]);
216  }
217  N_VDestroy_Serial(y);
218 
219  int flag = CVodeSensInit(m_cvode_mem, static_cast<int>(m_np),
220  CV_STAGGERED, CVSensRhsFn(0), m_yS);
221  checkError(flag, "sensInit", "CVodeSensInit");
222 
223  vector<double> atol(m_np);
224  for (size_t n = 0; n < m_np; n++) {
225  // This scaling factor is tuned so that reaction and species enthalpy
226  // sensitivities can be computed simultaneously with the same abstol.
227  atol[n] = m_abstolsens / func.m_paramScales[n];
228  }
229  flag = CVodeSensSStolerances(m_cvode_mem, m_reltolsens, atol.data());
230  checkError(flag, "sensInit", "CVodeSensSStolerances");
231 }
232 
234 {
235  m_neq = func.neq();
236  m_t0 = t0;
237  m_time = t0;
238  m_tInteg = t0;
239  m_func = &func;
240  func.clearErrors();
241  // Initialize preconditioner if applied
242  if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
243  m_preconditioner->initialize(m_neq);
244  }
245  if (m_y) {
246  N_VDestroy_Serial(m_y); // free solution vector if already allocated
247  }
248  m_y = newNVector(m_neq, m_sundials_ctx); // allocate solution vector
249  N_VConst(0.0, m_y);
250  if (m_dky) {
251  N_VDestroy_Serial(m_dky); // free derivative vector if already allocated
252  }
253  m_dky = newNVector(m_neq, m_sundials_ctx); // allocate derivative vector
254  N_VConst(0.0, m_dky);
255  // check abs tolerance array size
256  if (m_itol == CV_SV && m_nabs < m_neq) {
257  throw CanteraError("CVodesIntegrator::initialize",
258  "not enough absolute tolerance values specified.");
259  }
260 
261  func.getState(NV_DATA_S(m_y));
262 
263  if (m_cvode_mem) {
264  CVodeFree(&m_cvode_mem);
265  }
266 
267  //! Specify the method and the iteration type. Cantera Defaults:
268  //! CV_BDF - Use BDF methods
269  //! CV_NEWTON - use Newton's method
270  #if CT_SUNDIALS_VERSION < 40
271  m_cvode_mem = CVodeCreate(m_method, CV_NEWTON);
272  #elif CT_SUNDIALS_VERSION < 60
273  m_cvode_mem = CVodeCreate(m_method);
274  #else
275  m_cvode_mem = CVodeCreate(m_method, m_sundials_ctx.get());
276  #endif
277  if (!m_cvode_mem) {
278  throw CanteraError("CVodesIntegrator::initialize",
279  "CVodeCreate failed.");
280  }
281 
282  int flag = CVodeInit(m_cvode_mem, cvodes_rhs, m_t0, m_y);
283  if (flag != CV_SUCCESS) {
284  if (flag == CV_MEM_FAIL) {
285  throw CanteraError("CVodesIntegrator::initialize",
286  "Memory allocation failed.");
287  } else if (flag == CV_ILL_INPUT) {
288  throw CanteraError("CVodesIntegrator::initialize",
289  "Illegal value for CVodeInit input argument.");
290  } else {
291  throw CanteraError("CVodesIntegrator::initialize",
292  "CVodeInit failed.");
293  }
294  }
295  flag = CVodeSetErrHandlerFn(m_cvode_mem, &cvodes_err, this);
296 
297  if (m_itol == CV_SV) {
298  flag = CVodeSVtolerances(m_cvode_mem, m_reltol, m_abstol);
299  checkError(flag, "initialize", "CVodeSVtolerances");
300  } else {
301  flag = CVodeSStolerances(m_cvode_mem, m_reltol, m_abstols);
302  checkError(flag, "initialize", "CVodeSStolerances");
303  }
304 
305  flag = CVodeSetUserData(m_cvode_mem, &func);
306  checkError(flag, "initialize", "CVodeSetUserData");
307 
308  if (func.nparams() > 0) {
309  sensInit(t0, func);
310  flag = CVodeSetSensParams(m_cvode_mem, func.m_sens_params.data(),
311  func.m_paramScales.data(), NULL);
312  checkError(flag, "initialize", "CVodeSetSensParams");
313  }
314  applyOptions();
315 }
316 
317 void CVodesIntegrator::reinitialize(double t0, FuncEval& func)
318 {
319  m_t0 = t0;
320  m_time = t0;
321  m_tInteg = t0;
322  func.getState(NV_DATA_S(m_y));
323  m_func = &func;
324  func.clearErrors();
325  // reinitialize preconditioner if applied
326  if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
327  m_preconditioner->initialize(m_neq);
328  }
329  int result = CVodeReInit(m_cvode_mem, m_t0, m_y);
330  checkError(result, "reinitialize", "CVodeReInit");
331  applyOptions();
332 }
333 
335 {
336  if (m_type == "DENSE") {
337  sd_size_t N = static_cast<sd_size_t>(m_neq);
338  SUNLinSolFree((SUNLinearSolver) m_linsol);
339  SUNMatDestroy((SUNMatrix) m_linsol_matrix);
340  #if CT_SUNDIALS_VERSION >= 60
341  m_linsol_matrix = SUNDenseMatrix(N, N, m_sundials_ctx.get());
342  #else
343  m_linsol_matrix = SUNDenseMatrix(N, N);
344  #endif
345  if (m_linsol_matrix == nullptr) {
346  throw CanteraError("CVodesIntegrator::applyOptions",
347  "Unable to create SUNDenseMatrix of size {0} x {0}", N);
348  }
349  int flag;
350  #if CT_SUNDIALS_VERSION >= 60
351  #if CT_SUNDIALS_USE_LAPACK
352  m_linsol = SUNLinSol_LapackDense(m_y, (SUNMatrix) m_linsol_matrix,
353  m_sundials_ctx.get());
354  #else
355  m_linsol = SUNLinSol_Dense(m_y, (SUNMatrix) m_linsol_matrix,
356  m_sundials_ctx.get());
357  #endif
358  flag = CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
359  (SUNMatrix) m_linsol_matrix);
360  #else
361  #if CT_SUNDIALS_USE_LAPACK
362  m_linsol = SUNLapackDense(m_y, (SUNMatrix) m_linsol_matrix);
363  #else
364  m_linsol = SUNDenseLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
365  #endif
366  flag = CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
367  (SUNMatrix) m_linsol_matrix);
368  #endif
369  if (m_linsol == nullptr) {
370  throw CanteraError("CVodesIntegrator::applyOptions",
371  "Error creating Sundials dense linear solver object");
372  } else if (flag != CV_SUCCESS) {
373  throw CanteraError("CVodesIntegrator::applyOptions",
374  "Error connecting linear solver to CVODES. "
375  "Sundials error code: {}", flag);
376  }
377 
378  // throw preconditioner error for DENSE + NOJAC
379  if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
380  throw CanteraError("CVodesIntegrator::applyOptions",
381  "Preconditioning is not available with the specified problem type.");
382  }
383  } else if (m_type == "DIAG") {
384  CVDiag(m_cvode_mem);
385  // throw preconditioner error for DIAG
386  if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
387  throw CanteraError("CVodesIntegrator::applyOptions",
388  "Preconditioning is not available with the specified problem type.");
389  }
390  } else if (m_type == "GMRES") {
391  #if CT_SUNDIALS_VERSION >= 60
392  m_linsol = SUNLinSol_SPGMR(m_y, PREC_NONE, 0, m_sundials_ctx.get());
393  CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol, nullptr);
394  #elif CT_SUNDIALS_VERSION >= 40
395  m_linsol = SUNLinSol_SPGMR(m_y, PREC_NONE, 0);
396  CVSpilsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol);
397  # else
398  m_linsol = SUNSPGMR(m_y, PREC_NONE, 0);
399  CVSpilsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol);
400  #endif
401  // set preconditioner if used
402  #if CT_SUNDIALS_VERSION >= 40
403  if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
404  SUNLinSol_SPGMRSetPrecType((SUNLinearSolver) m_linsol,
405  static_cast<int>(m_prec_side));
406  CVodeSetPreconditioner(m_cvode_mem, cvodes_prec_setup,
407  cvodes_prec_solve);
408  }
409  #else
410  if (m_prec_side != PreconditionerSide::NO_PRECONDITION) {
411  SUNSPGMRSetPrecType((SUNLinearSolver) m_linsol,
412  static_cast<int>(m_prec_side));
413  CVSpilsSetPreconditioner(m_cvode_mem, cvodes_prec_setup,
414  cvodes_prec_solve);
415  }
416  #endif
417  } else if (m_type == "BAND") {
418  sd_size_t N = static_cast<sd_size_t>(m_neq);
419  sd_size_t nu = m_mupper;
420  sd_size_t nl = m_mlower;
421  SUNLinSolFree((SUNLinearSolver) m_linsol);
422  SUNMatDestroy((SUNMatrix) m_linsol_matrix);
423  #if CT_SUNDIALS_VERSION >= 60
424  m_linsol_matrix = SUNBandMatrix(N, nu, nl, m_sundials_ctx.get());
425  #elif CT_SUNDIALS_VERSION >= 40
426  m_linsol_matrix = SUNBandMatrix(N, nu, nl);
427  #else
428  m_linsol_matrix = SUNBandMatrix(N, nu, nl, nu+nl);
429  #endif
430  if (m_linsol_matrix == nullptr) {
431  throw CanteraError("CVodesIntegrator::applyOptions",
432  "Unable to create SUNBandMatrix of size {} with bandwidths "
433  "{} and {}", N, nu, nl);
434  }
435  #if CT_SUNDIALS_VERSION >= 60
436  #if CT_SUNDIALS_USE_LAPACK
437  m_linsol = SUNLinSol_LapackBand(m_y, (SUNMatrix) m_linsol_matrix,
438  m_sundials_ctx.get());
439  #else
440  m_linsol = SUNLinSol_Band(m_y, (SUNMatrix) m_linsol_matrix,
441  m_sundials_ctx.get());
442  #endif
443  CVodeSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
444  (SUNMatrix) m_linsol_matrix);
445  #else
446  #if CT_SUNDIALS_USE_LAPACK
447  m_linsol = SUNLapackBand(m_y, (SUNMatrix) m_linsol_matrix);
448  #else
449  m_linsol = SUNBandLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
450  #endif
451  CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
452  (SUNMatrix) m_linsol_matrix);
453  #endif
454  } else {
455  throw CanteraError("CVodesIntegrator::applyOptions",
456  "unsupported linear solver flag '{}'", m_type);
457  }
458 
459  if (m_maxord > 0) {
460  int flag = CVodeSetMaxOrd(m_cvode_mem, m_maxord);
461  checkError(flag, "applyOptions", "CVodeSetMaxOrd");
462  }
463  if (m_maxsteps > 0) {
464  CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
465  }
466  if (m_hmax > 0) {
467  CVodeSetMaxStep(m_cvode_mem, m_hmax);
468  }
469  if (m_hmin > 0) {
470  CVodeSetMinStep(m_cvode_mem, m_hmin);
471  }
472  if (m_maxErrTestFails > 0) {
473  CVodeSetMaxErrTestFails(m_cvode_mem, m_maxErrTestFails);
474  }
475 }
476 
478 {
479  if (tout == m_time) {
480  return;
481  } else if (tout < m_time) {
482  throw CanteraError("CVodesIntegrator::integrate",
483  "Cannot integrate backwards in time.\n"
484  "Requested time {} < current time {}",
485  tout, m_time);
486  }
487  int nsteps = 0;
488  while (m_tInteg < tout) {
489  if (nsteps >= m_maxsteps) {
490  throw CanteraError("CVodesIntegrator::integrate",
491  "Maximum number of timesteps ({}) taken without reaching output "
492  "time ({}).\nCurrent integrator time: {}",
493  nsteps, tout, m_tInteg);
494  }
495  int flag = CVode(m_cvode_mem, tout, m_y, &m_tInteg, CV_ONE_STEP);
496  if (flag != CV_SUCCESS) {
497  string f_errs = m_func->getErrors();
498  if (!f_errs.empty()) {
499  f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
500  }
501  throw CanteraError("CVodesIntegrator::integrate",
502  "CVodes error encountered. Error code: {}\n{}\n"
503  "{}"
504  "Components with largest weighted error estimates:\n{}",
505  flag, m_error_message, f_errs, getErrorInfo(10));
506  }
507  nsteps++;
508  }
509  int flag = CVodeGetDky(m_cvode_mem, tout, 0, m_y);
510  checkError(flag, "integrate", "CVodeGetDky");
511  m_time = tout;
512  m_sens_ok = false;
513 }
514 
515 double CVodesIntegrator::step(double tout)
516 {
517  int flag = CVode(m_cvode_mem, tout, m_y, &m_tInteg, CV_ONE_STEP);
518  if (flag != CV_SUCCESS) {
519  string f_errs = m_func->getErrors();
520  if (!f_errs.empty()) {
521  f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
522  }
523  throw CanteraError("CVodesIntegrator::step",
524  "CVodes error encountered. Error code: {}\n{}\n"
525  "{}"
526  "Components with largest weighted error estimates:\n{}",
527  flag, f_errs, m_error_message, getErrorInfo(10));
528 
529  }
530  m_sens_ok = false;
531  m_time = m_tInteg;
532  return m_time;
533 }
534 
535 double* CVodesIntegrator::derivative(double tout, int n)
536 {
537  int flag = CVodeGetDky(m_cvode_mem, tout, n, m_dky);
538  checkError(flag, "derivative", "CVodeGetDky");
539  return NV_DATA_S(m_dky);
540 }
541 
543 {
544  int ord;
545  CVodeGetLastOrder(m_cvode_mem, &ord);
546  return ord;
547 }
548 
550 {
551  long int ne;
552  CVodeGetNumRhsEvals(m_cvode_mem, &ne);
553  return ne;
554 }
555 
557 {
558  // AnyMap to return stats
559  AnyMap stats;
560 
561  // long int linear solver stats provided by CVodes
562  long int steps = 0, rhsEvals = 0, errTestFails = 0, jacEvals = 0, linSetup = 0,
563  linRhsEvals = 0, linIters = 0, linConvFails = 0, precEvals = 0,
564  precSolves = 0, jtSetupEvals = 0, jTimesEvals = 0, nonlinIters = 0,
565  nonlinConvFails = 0, orderReductions = 0;
566  int lastOrder = 0;
567 ;
568  #if CT_SUNDIALS_VERSION >= 40
569  CVodeGetNumSteps(m_cvode_mem, &steps);
570  CVodeGetNumRhsEvals(m_cvode_mem, &rhsEvals);
571  CVodeGetNonlinSolvStats(m_cvode_mem, &nonlinIters, &nonlinConvFails);
572  CVodeGetNumErrTestFails(m_cvode_mem, &errTestFails);
573  CVodeGetLastOrder(m_cvode_mem, &lastOrder);
574  CVodeGetNumStabLimOrderReds(m_cvode_mem, &orderReductions);
575  CVodeGetNumJacEvals(m_cvode_mem, &jacEvals);
576  CVodeGetNumLinRhsEvals(m_cvode_mem, &linRhsEvals);
577  CVodeGetNumLinSolvSetups(m_cvode_mem, &linSetup);
578  CVodeGetNumLinIters(m_cvode_mem, &linIters);
579  CVodeGetNumLinConvFails(m_cvode_mem, &linConvFails);
580  CVodeGetNumPrecEvals(m_cvode_mem, &precEvals);
581  CVodeGetNumPrecSolves(m_cvode_mem, &precSolves);
582  CVodeGetNumJTSetupEvals(m_cvode_mem, &jtSetupEvals);
583  CVodeGetNumJtimesEvals(m_cvode_mem, &jTimesEvals);
584  #else
585  warn_user("CVodesIntegrator::solverStats", "Function not"
586  "supported with sundials versions less than 4.");
587  #endif
588 
589  #if CT_SUNDIALS_VERION >= 60
590  long int stepSolveFails = 0;
591  CVodeGetNumStepSolveFails(m_cvode_mem, &stepSolveFails);
592  stats["step_solve_fails"] = stepSolveFails;
593  #endif
594 
595  stats["steps"] = steps;
596  stats["rhs_evals"] = rhsEvals;
597  stats["nonlinear_iters"] = nonlinIters;
598  stats["nonlinear_conv_fails"] = nonlinConvFails;
599  stats["err_test_fails"] = errTestFails;
600  stats["last_order"] = lastOrder;
601  stats["stab_order_reductions"] = orderReductions;
602 
603  stats["jac_evals"] = jacEvals;
604  stats["lin_solve_setups"] = linSetup;
605  stats["lin_rhs_evals"] = linRhsEvals;
606  stats["lin_iters"] = linIters;
607  stats["lin_conv_fails"] = linConvFails;
608  stats["prec_evals"] = precEvals;
609  stats["prec_solves"] = precSolves;
610  stats["jt_vec_setup_evals"] = jtSetupEvals;
611  stats["jt_vec_prod_evals"] = jTimesEvals;
612  return stats;
613 }
614 
615 double CVodesIntegrator::sensitivity(size_t k, size_t p)
616 {
617  if (m_time == m_t0) {
618  // calls to CVodeGetSens are only allowed after a successful time step.
619  return 0.0;
620  }
621  if (!m_sens_ok && m_np) {
622  int flag = CVodeGetSensDky(m_cvode_mem, m_time, 0, m_yS);
623  checkError(flag, "sensitivity", "CVodeGetSens");
624  m_sens_ok = true;
625  }
626 
627  if (k >= m_neq) {
628  throw CanteraError("CVodesIntegrator::sensitivity",
629  "sensitivity: k out of range ({})", k);
630  }
631  if (p >= m_np) {
632  throw CanteraError("CVodesIntegrator::sensitivity",
633  "sensitivity: p out of range ({})", p);
634  }
635  return NV_Ith_S(m_yS[p],k);
636 }
637 
639 {
640  N_Vector errs = newNVector(m_neq, m_sundials_ctx);
641  N_Vector errw = newNVector(m_neq, m_sundials_ctx);
642  CVodeGetErrWeights(m_cvode_mem, errw);
643  CVodeGetEstLocalErrors(m_cvode_mem, errs);
644 
645  vector<tuple<double, double, size_t>> weightedErrors;
646  for (size_t i=0; i<m_neq; i++) {
647  double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
648  weightedErrors.emplace_back(-abs(err), err, i);
649  }
650  N_VDestroy(errs);
651  N_VDestroy(errw);
652 
653  N = std::min(N, static_cast<int>(m_neq));
654  sort(weightedErrors.begin(), weightedErrors.end());
655  fmt::memory_buffer s;
656  for (int i=0; i<N; i++) {
657  fmt_append(s, "{}: {}\n",
658  get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
659  }
660  return to_string(s);
661 }
662 
663 void CVodesIntegrator::checkError(long flag, const string& ctMethod,
664  const string& cvodesMethod) const
665 {
666  if (flag == CV_SUCCESS) {
667  return;
668  } else if (flag == CV_MEM_NULL) {
669  throw CanteraError("CVodesIntegrator::" + ctMethod,
670  "CVODES integrator is not initialized");
671  } else {
672  const char* flagname = CVodeGetReturnFlagName(flag);
673  throw CanteraError("CVodesIntegrator::" + ctMethod,
674  "{} returned error code {} ({}):\n{}",
675  cvodesMethod, flag, flagname, m_error_message);
676  }
677 }
678 
679 }
A map of string keys to values whose type can vary at runtime.
Definition: AnyMap.h:427
Wrapper class for 'cvodes' integrator from LLNL.
double step(double tout) override
Integrate the system of equations.
void setMaxStepSize(double hmax) override
Set the maximum step size.
double * solution() override
The current value of the solution of the system of equations.
N_Vector m_y
The system state at m_time.
void * m_linsol
Sundials linear solver object.
int nEvals() const override
The number of function evaluations.
double m_time
The current system time, corresponding to m_y.
bool m_sens_ok
Indicates whether the sensitivities stored in m_yS have been updated for at the current integrator ti...
int maxSteps() override
Returns the maximum number of time-steps the integrator can take before reaching the next output time...
SundialsContext m_sundials_ctx
SUNContext object for Sundials>=6.0.
void checkError(long flag, const string &ctMethod, const string &cvodesMethod) const
Check whether a CVODES method indicated an error.
void setSensitivityTolerances(double reltol, double abstol) override
Set the sensitivity error tolerances.
void applyOptions()
Applies user-specified options to the underlying CVODES solver.
int lastOrder() const override
Order used during the last solution step.
void setMinStepSize(double hmin) override
Set the minimum step size.
void integrate(double tout) override
Integrate the system of equations.
void setTolerances(double reltol, size_t n, double *abstol) override
Set error tolerances.
double * derivative(double tout, int n) override
n-th derivative of the output function at time tout.
string m_error_message
Error message information provide by CVodes.
void setMaxErrTestFails(int n) override
Set the maximum permissible number of error test failures.
void setMaxSteps(int nmax) override
Set the maximum number of time-steps the integrator can take before reaching the next output time.
AnyMap solverStats() const override
Get solver stats from integrator.
void * m_linsol_matrix
matrix used by Sundials
double m_tInteg
The latest time reached by the integrator. May be greater than m_time.
void initialize(double t0, FuncEval &func) override
Initialize the integrator for a new problem.
string getErrorInfo(int N)
Returns a string listing the weighted error estimates associated with each solution component.
void setMethod(MethodType t) override
Set the solution method.
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:66
Virtual base class for ODE/DAE right-hand-side function evaluators.
Definition: FuncEval.h:32
vector< double > m_paramScales
Scaling factors for each sensitivity parameter.
Definition: FuncEval.h:184
void clearErrors()
Clear any previously-stored suppressed errors.
Definition: FuncEval.h:174
int evalNoThrow(double t, double *y, double *ydot)
Evaluate the right-hand side using return code to indicate status.
Definition: FuncEval.cpp:7
virtual size_t neq() const =0
Number of equations.
virtual size_t nparams() const
Number of sensitivity parameters.
Definition: FuncEval.h:156
string getErrors() const
Return a string containing the text of any suppressed errors.
Definition: FuncEval.cpp:71
vector< double > m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
Definition: FuncEval.h:176
virtual void getState(double *y)
Fill in the vector y with the current state of the system.
Definition: FuncEval.h:142
shared_ptr< PreconditionerBase > m_preconditioner
Pointer to preconditioner object used in integration which is set by setPreconditioner and initialize...
Definition: Integrator.h:299
PreconditionerSide m_prec_side
Type of preconditioning used in applyOptions.
Definition: Integrator.h:301
A wrapper for managing a SUNContext object, need for Sundials >= 6.0.
void fmt_append(fmt::memory_buffer &b, Args... args)
Versions 6.2.0 and 6.2.1 of fmtlib do not include this define before they include windows....
Definition: fmt.h:29
static int cvodes_rhs(realtype t, N_Vector y, N_Vector ydot, void *f_data)
Function called by cvodes to evaluate ydot given y.
void warn_user(const string &method, const string &msg, const Args &... args)
Print a user warning raised from method as CanteraWarning.
Definition: global.h:267
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:564
MethodType
Specifies the method used to integrate the system of equations.
Definition: Integrator.h:23
@ BDF_Method
Backward Differentiation.
Definition: Integrator.h:24
@ Adams_Method
Adams.
Definition: Integrator.h:25
static void cvodes_err(int error_code, const char *module, const char *function, char *msg, void *eh_data)
Function called by CVodes when an error is encountered instead of writing to stdout.
Contains declarations for string manipulation functions within Cantera.