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