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