Cantera  2.5.1
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 "sundials/sundials_types.h"
13 #include "sundials/sundials_math.h"
14 #include "sundials/sundials_nvector.h"
15 #include "nvector/nvector_serial.h"
16 #include "cvodes/cvodes.h"
17 #if CT_SUNDIALS_VERSION >= 30
18  #if CT_SUNDIALS_USE_LAPACK
19  #include "sunlinsol/sunlinsol_lapackdense.h"
20  #include "sunlinsol/sunlinsol_lapackband.h"
21  #else
22  #include "sunlinsol/sunlinsol_dense.h"
23  #include "sunlinsol/sunlinsol_band.h"
24  #endif
25  #include "sunlinsol/sunlinsol_spgmr.h"
26  #include "cvodes/cvodes_direct.h"
27  #include "cvodes/cvodes_diag.h"
28  #include "cvodes/cvodes_spils.h"
29 #else
30  #if CT_SUNDIALS_USE_LAPACK
31  #include "cvodes/cvodes_lapack.h"
32  #else
33  #include "cvodes/cvodes_dense.h"
34  #include "cvodes/cvodes_band.h"
35  #endif
36  #include "cvodes/cvodes_diag.h"
37  #include "cvodes/cvodes_spgmr.h"
38 #endif
39 
40 #define CV_SS 1
41 #define CV_SV 2
42 
43 #if CT_SUNDIALS_VERSION < 25
44 typedef int sd_size_t;
45 #else
46 typedef long int sd_size_t;
47 #endif
48 
49 namespace Cantera
50 {
51 
52 extern "C" {
53  /**
54  * Function called by cvodes to evaluate ydot given y. The CVODE integrator
55  * allows passing in a void* pointer to access external data. This pointer
56  * is cast to a pointer to a instance of class FuncEval. The equations to be
57  * integrated should be specified by deriving a class from FuncEval that
58  * evaluates the desired equations.
59  * @ingroup odeGroup
60  */
61  static int cvodes_rhs(realtype t, N_Vector y, N_Vector ydot, void* f_data)
62  {
63  FuncEval* f = (FuncEval*) f_data;
64  return f->eval_nothrow(t, NV_DATA_S(y), NV_DATA_S(ydot));
65  }
66 
67  //! Function called by CVodes when an error is encountered instead of
68  //! writing to stdout. Here, save the error message provided by CVodes so
69  //! that it can be included in the subsequently raised CanteraError.
70  static void cvodes_err(int error_code, const char* module,
71  const char* function, char* msg, void* eh_data)
72  {
73  CVodesIntegrator* integrator = (CVodesIntegrator*) eh_data;
74  integrator->m_error_message = msg;
75  integrator->m_error_message += "\n";
76  }
77 }
78 
79 CVodesIntegrator::CVodesIntegrator() :
80  m_neq(0),
81  m_cvode_mem(0),
82  m_linsol(0),
83  m_linsol_matrix(0),
84  m_func(0),
85  m_t0(0.0),
86  m_y(0),
87  m_abstol(0),
88  m_dky(0),
89  m_type(DENSE+NOJAC),
90  m_itol(CV_SS),
91  m_method(CV_BDF),
92  m_maxord(0),
93  m_reltol(1.e-9),
94  m_abstols(1.e-15),
95  m_reltolsens(1.0e-5),
96  m_abstolsens(1.0e-4),
97  m_nabs(0),
98  m_hmax(0.0),
99  m_hmin(0.0),
100  m_maxsteps(20000),
101  m_maxErrTestFails(0),
102  m_yS(nullptr),
103  m_np(0),
104  m_mupper(0), m_mlower(0),
105  m_sens_ok(false)
106 {
107 }
108 
109 CVodesIntegrator::~CVodesIntegrator()
110 {
111  if (m_cvode_mem) {
112  if (m_np > 0) {
113  CVodeSensFree(m_cvode_mem);
114  }
115  CVodeFree(&m_cvode_mem);
116  }
117 
118  #if CT_SUNDIALS_VERSION >= 30
119  SUNLinSolFree((SUNLinearSolver) m_linsol);
120  SUNMatDestroy((SUNMatrix) m_linsol_matrix);
121  #endif
122 
123  if (m_y) {
124  N_VDestroy_Serial(m_y);
125  }
126  if (m_abstol) {
127  N_VDestroy_Serial(m_abstol);
128  }
129  if (m_dky) {
130  N_VDestroy_Serial(m_dky);
131  }
132  if (m_yS) {
133  N_VDestroyVectorArray_Serial(m_yS, static_cast<sd_size_t>(m_np));
134  }
135 }
136 
137 double& CVodesIntegrator::solution(size_t k)
138 {
139  return NV_Ith_S(m_y, k);
140 }
141 
143 {
144  return NV_DATA_S(m_y);
145 }
146 
147 void CVodesIntegrator::setTolerances(double reltol, size_t n, double* abstol)
148 {
149  m_itol = CV_SV;
150  m_nabs = n;
151  if (n != m_neq) {
152  if (m_abstol) {
153  N_VDestroy_Serial(m_abstol);
154  }
155  m_abstol = N_VNew_Serial(static_cast<sd_size_t>(n));
156  }
157  for (size_t i=0; i<n; i++) {
158  NV_Ith_S(m_abstol, i) = abstol[i];
159  }
160  m_reltol = reltol;
161 }
162 
163 void CVodesIntegrator::setTolerances(double reltol, double abstol)
164 {
165  m_itol = CV_SS;
166  m_reltol = reltol;
167  m_abstols = abstol;
168 }
169 
170 void CVodesIntegrator::setSensitivityTolerances(double reltol, double abstol)
171 {
172  m_reltolsens = reltol;
173  m_abstolsens = abstol;
174 }
175 
177 {
178  m_type = probtype;
179 }
180 
182 {
183  if (t == BDF_Method) {
184  m_method = CV_BDF;
185  } else if (t == Adams_Method) {
186  m_method = CV_ADAMS;
187  } else {
188  throw CanteraError("CVodesIntegrator::setMethod", "unknown method");
189  }
190 }
191 
192 void CVodesIntegrator::setMaxStepSize(doublereal hmax)
193 {
194  m_hmax = hmax;
195  if (m_cvode_mem) {
196  CVodeSetMaxStep(m_cvode_mem, hmax);
197  }
198 }
199 
200 void CVodesIntegrator::setMinStepSize(doublereal hmin)
201 {
202  m_hmin = hmin;
203  if (m_cvode_mem) {
204  CVodeSetMinStep(m_cvode_mem, hmin);
205  }
206 }
207 
209 {
210  m_maxsteps = nmax;
211  if (m_cvode_mem) {
212  CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
213  }
214 }
215 
217 {
218  return m_maxsteps;
219 }
220 
222 {
223  m_maxErrTestFails = n;
224  if (m_cvode_mem) {
225  CVodeSetMaxErrTestFails(m_cvode_mem, n);
226  }
227 }
228 
229 void CVodesIntegrator::sensInit(double t0, FuncEval& func)
230 {
231  m_np = func.nparams();
232  m_sens_ok = false;
233 
234  N_Vector y = N_VNew_Serial(static_cast<sd_size_t>(func.neq()));
235  m_yS = N_VCloneVectorArray_Serial(static_cast<sd_size_t>(m_np), y);
236  for (size_t n = 0; n < m_np; n++) {
237  N_VConst(0.0, m_yS[n]);
238  }
239  N_VDestroy_Serial(y);
240 
241  int flag = CVodeSensInit(m_cvode_mem, static_cast<sd_size_t>(m_np),
242  CV_STAGGERED, CVSensRhsFn(0), m_yS);
243 
244  if (flag != CV_SUCCESS) {
245  throw CanteraError("CVodesIntegrator::sensInit", "Error in CVodeSensInit");
246  }
247  vector_fp atol(m_np);
248  for (size_t n = 0; n < m_np; n++) {
249  // This scaling factor is tuned so that reaction and species enthalpy
250  // sensitivities can be computed simultaneously with the same abstol.
251  atol[n] = m_abstolsens / func.m_paramScales[n];
252  }
253  flag = CVodeSensSStolerances(m_cvode_mem, m_reltolsens, atol.data());
254 }
255 
257 {
258  m_neq = func.neq();
259  m_t0 = t0;
260  m_time = t0;
261  m_func = &func;
262  func.clearErrors();
263 
264  if (m_y) {
265  N_VDestroy_Serial(m_y); // free solution vector if already allocated
266  }
267  m_y = N_VNew_Serial(static_cast<sd_size_t>(m_neq)); // allocate solution vector
268  N_VConst(0.0, m_y);
269  if (m_dky) {
270  N_VDestroy_Serial(m_dky); // free derivative vector if already allocated
271  }
272  m_dky = N_VNew_Serial(static_cast<sd_size_t>(m_neq)); // allocate derivative vector
273  N_VConst(0.0, m_dky);
274  // check abs tolerance array size
275  if (m_itol == CV_SV && m_nabs < m_neq) {
276  throw CanteraError("CVodesIntegrator::initialize",
277  "not enough absolute tolerance values specified.");
278  }
279 
280  func.getState(NV_DATA_S(m_y));
281 
282  if (m_cvode_mem) {
283  CVodeFree(&m_cvode_mem);
284  }
285 
286  //! Specify the method and the iteration type. Cantera Defaults:
287  //! CV_BDF - Use BDF methods
288  //! CV_NEWTON - use Newton's method
289  #if CT_SUNDIALS_VERSION < 40
290  m_cvode_mem = CVodeCreate(m_method, CV_NEWTON);
291  #else
292  m_cvode_mem = CVodeCreate(m_method);
293  #endif
294  if (!m_cvode_mem) {
295  throw CanteraError("CVodesIntegrator::initialize",
296  "CVodeCreate failed.");
297  }
298 
299  int flag = CVodeInit(m_cvode_mem, cvodes_rhs, m_t0, m_y);
300  if (flag != CV_SUCCESS) {
301  if (flag == CV_MEM_FAIL) {
302  throw CanteraError("CVodesIntegrator::initialize",
303  "Memory allocation failed.");
304  } else if (flag == CV_ILL_INPUT) {
305  throw CanteraError("CVodesIntegrator::initialize",
306  "Illegal value for CVodeInit input argument.");
307  } else {
308  throw CanteraError("CVodesIntegrator::initialize",
309  "CVodeInit failed.");
310  }
311  }
312  CVodeSetErrHandlerFn(m_cvode_mem, &cvodes_err, this);
313 
314  if (m_itol == CV_SV) {
315  flag = CVodeSVtolerances(m_cvode_mem, m_reltol, m_abstol);
316  } else {
317  flag = CVodeSStolerances(m_cvode_mem, m_reltol, m_abstols);
318  }
319  if (flag != CV_SUCCESS) {
320  if (flag == CV_MEM_FAIL) {
321  throw CanteraError("CVodesIntegrator::initialize",
322  "Memory allocation failed.");
323  } else if (flag == CV_ILL_INPUT) {
324  throw CanteraError("CVodesIntegrator::initialize",
325  "Illegal value for CVodeInit input argument.");
326  } else {
327  throw CanteraError("CVodesIntegrator::initialize",
328  "CVodeInit failed.");
329  }
330  }
331 
332  flag = CVodeSetUserData(m_cvode_mem, &func);
333  if (flag != CV_SUCCESS) {
334  throw CanteraError("CVodesIntegrator::initialize",
335  "CVodeSetUserData failed.");
336  }
337  if (func.nparams() > 0) {
338  sensInit(t0, func);
339  flag = CVodeSetSensParams(m_cvode_mem, func.m_sens_params.data(),
340  func.m_paramScales.data(), NULL);
341  if (flag != CV_SUCCESS) {
342  throw CanteraError("CVodesIntegrator::initialize",
343  "CVodeSetSensParams failed.");
344  }
345  }
346  applyOptions();
347 }
348 
349 void CVodesIntegrator::reinitialize(double t0, FuncEval& func)
350 {
351  m_t0 = t0;
352  m_time = t0;
353  func.getState(NV_DATA_S(m_y));
354  m_func = &func;
355  func.clearErrors();
356 
357  int result = CVodeReInit(m_cvode_mem, m_t0, m_y);
358  if (result != CV_SUCCESS) {
359  throw CanteraError("CVodesIntegrator::reinitialize",
360  "CVodeReInit failed. result = {}", result);
361  }
362  applyOptions();
363 }
364 
366 {
367  if (m_type == DENSE + NOJAC) {
368  sd_size_t N = static_cast<sd_size_t>(m_neq);
369  #if CT_SUNDIALS_VERSION >= 30
370  SUNLinSolFree((SUNLinearSolver) m_linsol);
371  SUNMatDestroy((SUNMatrix) m_linsol_matrix);
372  m_linsol_matrix = SUNDenseMatrix(N, N);
373  if (m_linsol_matrix == nullptr) {
374  throw CanteraError("CVodesIntegrator::applyOptions",
375  "Unable to create SUNDenseMatrix of size {0} x {0}", N);
376  }
377  #if CT_SUNDIALS_USE_LAPACK
378  m_linsol = SUNLapackDense(m_y, (SUNMatrix) m_linsol_matrix);
379  #else
380  m_linsol = SUNDenseLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
381  #endif
382  CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
383  (SUNMatrix) m_linsol_matrix);
384  #else
385  #if CT_SUNDIALS_USE_LAPACK
386  CVLapackDense(m_cvode_mem, N);
387  #else
388  CVDense(m_cvode_mem, N);
389  #endif
390  #endif
391  } else if (m_type == DIAG) {
392  CVDiag(m_cvode_mem);
393  } else if (m_type == GMRES) {
394  #if CT_SUNDIALS_VERSION >= 30
395  m_linsol = SUNSPGMR(m_y, PREC_NONE, 0);
396  CVSpilsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol);
397  #else
398  CVSpgmr(m_cvode_mem, PREC_NONE, 0);
399  #endif
400  } else if (m_type == BAND + NOJAC) {
401  sd_size_t N = static_cast<sd_size_t>(m_neq);
402  long int nu = m_mupper;
403  long int nl = m_mlower;
404  #if CT_SUNDIALS_VERSION >= 30
405  SUNLinSolFree((SUNLinearSolver) m_linsol);
406  SUNMatDestroy((SUNMatrix) m_linsol_matrix);
407  #if CT_SUNDIALS_VERSION < 40
408  m_linsol_matrix = SUNBandMatrix(N, nu, nl, nu+nl);
409  #else
410  m_linsol_matrix = SUNBandMatrix(N, nu, nl);
411  #endif
412  if (m_linsol_matrix == nullptr) {
413  throw CanteraError("CVodesIntegrator::applyOptions",
414  "Unable to create SUNBandMatrix of size {} with bandwidths "
415  "{} and {}", N, nu, nl);
416  }
417  #if CT_SUNDIALS_USE_LAPACK
418  m_linsol = SUNLapackBand(m_y, (SUNMatrix) m_linsol_matrix);
419  #else
420  m_linsol = SUNBandLinearSolver(m_y, (SUNMatrix) m_linsol_matrix);
421  #endif
422  CVDlsSetLinearSolver(m_cvode_mem, (SUNLinearSolver) m_linsol,
423  (SUNMatrix) m_linsol_matrix);
424  #else
425  #if CT_SUNDIALS_USE_LAPACK
426  CVLapackBand(m_cvode_mem, N, nu, nl);
427  #else
428  CVBand(m_cvode_mem, N, nu, nl);
429  #endif
430  #endif
431  } else {
432  throw CanteraError("CVodesIntegrator::applyOptions",
433  "unsupported option");
434  }
435 
436  if (m_maxord > 0) {
437  CVodeSetMaxOrd(m_cvode_mem, m_maxord);
438  }
439  if (m_maxsteps > 0) {
440  CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
441  }
442  if (m_hmax > 0) {
443  CVodeSetMaxStep(m_cvode_mem, m_hmax);
444  }
445  if (m_hmin > 0) {
446  CVodeSetMinStep(m_cvode_mem, m_hmin);
447  }
448  if (m_maxErrTestFails > 0) {
449  CVodeSetMaxErrTestFails(m_cvode_mem, m_maxErrTestFails);
450  }
451 }
452 
454 {
455  if (tout == m_time) {
456  return;
457  }
458  int flag = CVode(m_cvode_mem, tout, m_y, &m_time, CV_NORMAL);
459  if (flag != CV_SUCCESS) {
460  string f_errs = m_func->getErrors();
461  if (!f_errs.empty()) {
462  f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
463  }
464  throw CanteraError("CVodesIntegrator::integrate",
465  "CVodes error encountered. Error code: {}\n{}\n"
466  "{}"
467  "Components with largest weighted error estimates:\n{}",
468  flag, m_error_message, f_errs, getErrorInfo(10));
469  }
470  m_sens_ok = false;
471 }
472 
473 double CVodesIntegrator::step(double tout)
474 {
475  int flag = CVode(m_cvode_mem, tout, m_y, &m_time, CV_ONE_STEP);
476  if (flag != CV_SUCCESS) {
477  string f_errs = m_func->getErrors();
478  if (!f_errs.empty()) {
479  f_errs = "Exceptions caught during RHS evaluation:\n" + f_errs;
480  }
481  throw CanteraError("CVodesIntegrator::step",
482  "CVodes error encountered. Error code: {}\n{}\n"
483  "{}"
484  "Components with largest weighted error estimates:\n{}",
485  flag, f_errs, m_error_message, getErrorInfo(10));
486 
487  }
488  m_sens_ok = false;
489  return m_time;
490 }
491 
492 double* CVodesIntegrator::derivative(double tout, int n)
493 {
494  int flag = CVodeGetDky(m_cvode_mem, tout, n, m_dky);
495  if (flag != CV_SUCCESS) {
496  string f_errs = m_func->getErrors();
497  if (!f_errs.empty()) {
498  f_errs = "Exceptions caught evaluating derivative:\n" + f_errs;
499  }
500  throw CanteraError("CVodesIntegrator::derivative",
501  "CVodes error encountered. Error code: {}\n{}\n"
502  "{}",
503  flag, m_error_message, f_errs);
504  }
505  return NV_DATA_S(m_dky);
506 }
507 
509 {
510  int ord;
511  CVodeGetLastOrder(m_cvode_mem, &ord);
512  return ord;
513 }
514 
516 {
517  long int ne;
518  CVodeGetNumRhsEvals(m_cvode_mem, &ne);
519  return ne;
520 }
521 
522 double CVodesIntegrator::sensitivity(size_t k, size_t p)
523 {
524  if (m_time == m_t0) {
525  // calls to CVodeGetSens are only allowed after a successful time step.
526  return 0.0;
527  }
528  if (!m_sens_ok && m_np) {
529  int flag = CVodeGetSens(m_cvode_mem, &m_time, m_yS);
530  if (flag != CV_SUCCESS) {
531  throw CanteraError("CVodesIntegrator::sensitivity",
532  "CVodeGetSens failed. Error code: {}", flag);
533  }
534  m_sens_ok = true;
535  }
536 
537  if (k >= m_neq) {
538  throw CanteraError("CVodesIntegrator::sensitivity",
539  "sensitivity: k out of range ({})", k);
540  }
541  if (p >= m_np) {
542  throw CanteraError("CVodesIntegrator::sensitivity",
543  "sensitivity: p out of range ({})", p);
544  }
545  return NV_Ith_S(m_yS[p],k);
546 }
547 
549 {
550  N_Vector errs = N_VNew_Serial(static_cast<sd_size_t>(m_neq));
551  N_Vector errw = N_VNew_Serial(static_cast<sd_size_t>(m_neq));
552  CVodeGetErrWeights(m_cvode_mem, errw);
553  CVodeGetEstLocalErrors(m_cvode_mem, errs);
554 
555  vector<tuple<double, double, size_t> > weightedErrors;
556  for (size_t i=0; i<m_neq; i++) {
557  double err = NV_Ith_S(errs, i) * NV_Ith_S(errw, i);
558  weightedErrors.emplace_back(-abs(err), err, i);
559  }
560  N_VDestroy(errs);
561  N_VDestroy(errw);
562 
563  N = std::min(N, static_cast<int>(m_neq));
564  sort(weightedErrors.begin(), weightedErrors.end());
565  fmt::memory_buffer s;
566  for (int i=0; i<N; i++) {
567  format_to(s, "{}: {}\n",
568  get<2>(weightedErrors[i]), get<1>(weightedErrors[i]));
569  }
570  return to_string(s);
571 }
572 
573 }
Wrapper class for 'cvodes' integrator from LLNL.
virtual int nEvals() const
The number of function evaluations.
virtual int lastOrder() const
Order used during the last solution step.
virtual void setMinStepSize(double hmin)
Set the minimum step size.
void * m_linsol
Sundials linear solver object.
std::string m_error_message
Error message information provide by CVodes.
double m_time
The current integrator time.
bool m_sens_ok
Indicates whether the sensitivities stored in m_yS have been updated for at the current integrator ti...
virtual void initialize(double t0, FuncEval &func)
Initialize the integrator for a new problem.
virtual void integrate(double tout)
Integrate the system of equations.
virtual void setMaxSteps(int nmax)
Set the maximum number of time-steps the integrator can take before reaching the next output time.
virtual doublereal step(double tout)
Integrate the system of equations.
virtual void setSensitivityTolerances(double reltol, double abstol)
Set the sensitivity error tolerances.
virtual double * solution()
The current value of the solution of the system of equations.
virtual void setTolerances(double reltol, size_t n, double *abstol)
Set error tolerances.
virtual int maxSteps()
Returns the maximum number of time-steps the integrator can take before reaching the next output time...
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
void applyOptions()
Applies user-specified options to the underlying CVODES solver.
virtual double * derivative(double tout, int n)
n-th derivative of the output function at time tout.
virtual void setProblemType(int probtype)
Set the problem type.
virtual void setMethod(MethodType t)
Set the solution method.
void * m_linsol_matrix
matrix used by Sundials
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
virtual std::string getErrorInfo(int N)
Returns a string listing the weighted error estimates associated with each solution component.
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
Virtual base class for ODE right-hand-side function evaluators.
Definition: FuncEval.h:27
virtual size_t nparams()
Number of sensitivity parameters.
Definition: FuncEval.h:61
vector_fp m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
Definition: FuncEval.h:81
void clearErrors()
Clear any previously-stored suppressed errors.
Definition: FuncEval.h:79
virtual size_t neq()=0
Number of equations.
std::string getErrors() const
Return a string containing the text of any suppressed errors.
Definition: FuncEval.cpp:45
int eval_nothrow(double t, double *y, double *ydot)
Evaluate the right-hand side using return code to indicate status.
Definition: FuncEval.cpp:12
vector_fp m_paramScales
Scaling factors for each sensitivity parameter.
Definition: FuncEval.h:89
virtual void getState(double *y)
Fill in the vector y with the current state of the system.
Definition: FuncEval.h:53
static int cvodes_rhs(realtype t, N_Vector y, N_Vector ydot, void *f_data)
Function called by cvodes to evaluate ydot given y.
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:180
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:264
MethodType
Specifies the method used to integrate the system of equations.
Definition: Integrator.h:32
@ BDF_Method
Backward Differentiation.
Definition: Integrator.h:33
@ Adams_Method
Adams.
Definition: Integrator.h:34
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.