Cantera  2.0
CVodesIntegrator.cpp
Go to the documentation of this file.
1 /**
2  * @file CVodesIntegrator.cpp
3  *
4  */
5 
6 // Copyright 2001 California Institute of Technology
7 #include "cantera/base/config.h"
8 
9 #include "CVodesIntegrator.h"
11 
12 #include <iostream>
13 using namespace std;
14 
15 #if SUNDIALS_VERSION == 22
16 
17 #include "sundials_types.h"
18 #include "sundials_math.h"
19 #include "cvodes.h"
20 #include "cvodes_dense.h"
21 #include "cvodes_diag.h"
22 #include "cvodes_spgmr.h"
23 #include "cvodes_band.h"
24 #include "nvector_serial.h"
25 
26 #else
27 
28 #if SUNDIALS_VERSION >= 23
29 #include "sundials/sundials_types.h"
30 #include "sundials/sundials_math.h"
31 #include "sundials/sundials_nvector.h"
32 #include "nvector/nvector_serial.h"
33 #include "cvodes/cvodes.h"
34 #include "cvodes/cvodes_dense.h"
35 #include "cvodes/cvodes_diag.h"
36 #include "cvodes/cvodes_spgmr.h"
37 #include "cvodes/cvodes_band.h"
38 
39 #else
40 #error unsupported Sundials version!
41 #endif
42 
43 #if SUNDIALS_VERSION >= 24
44 #define CV_SS 1
45 #define CV_SV 2
46 #endif
47 
48 #endif
49 
50 inline static N_Vector nv(void* x)
51 {
52  return reinterpret_cast<N_Vector>(x);
53 }
54 
55 namespace Cantera
56 {
57 
58 class FuncData
59 {
60 public:
61  FuncData(FuncEval* f, int npar = 0) {
62  m_pars.resize(npar, 1.0);
63  m_func = f;
64  }
65  virtual ~FuncData() {}
66  vector_fp m_pars;
67  FuncEval* m_func;
68 };
69 }
70 
71 
72 extern "C" {
73 
74  /**
75  * Function called by cvodes to evaluate ydot given y. The cvode
76  * integrator allows passing in a void* pointer to access
77  * external data. This pointer is cast to a pointer to a instance
78  * of class FuncEval. The equations to be integrated should be
79  * specified by deriving a class from FuncEval that evaluates the
80  * desired equations.
81  * @ingroup odeGroup
82  */
83  static int cvodes_rhs(realtype t, N_Vector y, N_Vector ydot,
84  void* f_data)
85  {
86  double* ydata = NV_DATA_S(y); //N_VDATA(y);
87  double* ydotdata = NV_DATA_S(ydot); //N_VDATA(ydot);
88  Cantera::FuncData* d = (Cantera::FuncData*)f_data;
89  Cantera::FuncEval* f = d->m_func;
90  //try {
91  if (d->m_pars.size() == 0) {
92  f->eval(t, ydata, ydotdata, NULL);
93  } else {
94  f->eval(t, ydata, ydotdata, DATA_PTR(d->m_pars));
95  }
96  //}
97  //catch (...) {
98  //Cantera::showErrors();
99  //Cantera::error("Teminating execution");
100  //}
101  return 0;
102  }
103 
104 }
105 
106 namespace Cantera
107 {
108 
109 
110 /**
111  * Constructor. Default settings: dense jacobian, no user-supplied
112  * Jacobian function, Newton iteration.
113  */
114 CVodesIntegrator::CVodesIntegrator() :
115  m_neq(0),
116  m_cvode_mem(0),
117  m_t0(0.0),
118  m_y(0),
119  m_abstol(0),
120  m_type(DENSE+NOJAC),
121  m_itol(CV_SS),
122  m_method(CV_BDF),
123  m_iter(CV_NEWTON),
124  m_maxord(0),
125  m_reltol(1.e-9),
126  m_abstols(1.e-15),
127  m_reltolsens(1.0e-5),
128  m_abstolsens(1.0e-4),
129  m_nabs(0),
130  m_hmax(0.0),
131  m_maxsteps(20000),
132  m_fdata(0),
133  m_np(0),
134  m_mupper(0), m_mlower(0)
135 {
136  //m_ropt.resize(OPT_SIZE,0.0);
137  //m_iopt = new long[OPT_SIZE];
138  //fill(m_iopt, m_iopt+OPT_SIZE,0);
139 }
140 
141 
142 /// Destructor.
143 CVodesIntegrator::~CVodesIntegrator()
144 {
145  if (m_cvode_mem) {
146  if (m_np > 0) {
147  CVodeSensFree(m_cvode_mem);
148  }
149  CVodeFree(&m_cvode_mem);
150  }
151  if (m_y) {
152  N_VDestroy_Serial(nv(m_y));
153  }
154  if (m_abstol) {
155  N_VDestroy_Serial(nv(m_abstol));
156  }
157  delete m_fdata;
158 
159  //delete[] m_iopt;
160 }
161 
162 double& CVodesIntegrator::solution(size_t k)
163 {
164  return NV_Ith_S(nv(m_y),k);
165 }
166 
167 double* CVodesIntegrator::solution()
168 {
169  return NV_DATA_S(nv(m_y));
170 }
171 
172 void CVodesIntegrator::setTolerances(double reltol, size_t n, double* abstol)
173 {
174  m_itol = CV_SV;
175  m_nabs = n;
176  if (n != m_neq) {
177  if (m_abstol) {
178  N_VDestroy_Serial(nv(m_abstol));
179  }
180  m_abstol = reinterpret_cast<void*>(N_VNew_Serial(n));
181  }
182  for (size_t i=0; i<n; i++) {
183  NV_Ith_S(nv(m_abstol), i) = abstol[i];
184  }
185  m_reltol = reltol;
186 }
187 
188 void CVodesIntegrator::setTolerances(double reltol, double abstol)
189 {
190  m_itol = CV_SS;
191  m_reltol = reltol;
192  m_abstols = abstol;
193 }
194 
195 void CVodesIntegrator::setSensitivityTolerances(double reltol, double abstol)
196 {
197  m_reltolsens = reltol;
198  m_abstolsens = abstol;
199 }
200 
201 void CVodesIntegrator::setProblemType(int probtype)
202 {
203  m_type = probtype;
204 }
205 
206 void CVodesIntegrator::setMethod(MethodType t)
207 {
208  if (t == BDF_Method) {
209  m_method = CV_BDF;
210  } else if (t == Adams_Method) {
211  m_method = CV_ADAMS;
212  } else {
213  throw CVodesErr("unknown method");
214  }
215 }
216 
217 void CVodesIntegrator::setMaxStepSize(doublereal hmax)
218 {
219  m_hmax = hmax;
220  if (m_cvode_mem) {
221  CVodeSetMaxStep(m_cvode_mem, hmax);
222  }
223 }
224 
225 void CVodesIntegrator::setMinStepSize(doublereal hmin)
226 {
227  m_hmin = hmin;
228  if (m_cvode_mem) {
229  CVodeSetMinStep(m_cvode_mem, hmin);
230  }
231 }
232 
233 void CVodesIntegrator::setMaxSteps(int nmax)
234 {
235  m_maxsteps = nmax;
236  if (m_cvode_mem) {
237  CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
238  }
239 }
240 
241 void CVodesIntegrator::setIterator(IterType t)
242 {
243  if (t == Newton_Iter) {
244  m_iter = CV_NEWTON;
245  } else if (t == Functional_Iter) {
246  m_iter = CV_FUNCTIONAL;
247  } else {
248  throw CVodesErr("unknown iterator");
249  }
250 }
251 
252 void CVodesIntegrator::sensInit(double t0, FuncEval& func)
253 {
254  m_np = func.nparams();
255  size_t nv = func.neq();
256 
257  doublereal* data;
258  N_Vector y;
259  y = N_VNew_Serial(nv);
260  m_yS = N_VCloneVectorArray_Serial(m_np, y);
261  for (size_t n = 0; n < m_np; n++) {
262  data = NV_DATA_S(m_yS[n]);
263  for (size_t j = 0; j < nv; j++) {
264  data[j] =0.0;
265  }
266  }
267 
268  int flag;
269 
270 #if SUNDIALS_VERSION <= 23
271  flag = CVodeSensMalloc(m_cvode_mem, m_np, CV_STAGGERED, m_yS);
272  if (flag != CV_SUCCESS) {
273  throw CVodesErr("Error in CVodeSensMalloc");
274  }
275  vector_fp atol(m_np, m_abstolsens);
276  double rtol = m_reltolsens;
277  flag = CVodeSetSensTolerances(m_cvode_mem, CV_SS, rtol, DATA_PTR(atol));
278 #elif SUNDIALS_VERSION >= 24
279  flag = CVodeSensInit(m_cvode_mem, m_np, CV_STAGGERED,
280  CVSensRhsFn(0), m_yS);
281 
282  if (flag != CV_SUCCESS) {
283  throw CVodesErr("Error in CVodeSensMalloc");
284  }
285  vector_fp atol(m_np, m_abstolsens);
286  double rtol = m_reltolsens;
287  flag = CVodeSensSStolerances(m_cvode_mem, rtol, DATA_PTR(atol));
288 #endif
289 
290 }
291 
292 void CVodesIntegrator::initialize(double t0, FuncEval& func)
293 {
294  m_neq = func.neq();
295  m_t0 = t0;
296 
297  if (m_y) {
298  N_VDestroy_Serial(nv(m_y)); // free solution vector if already allocated
299  }
300  m_y = reinterpret_cast<void*>(N_VNew_Serial(m_neq)); // allocate solution vector
301  for (size_t i = 0; i < m_neq; i++) {
302  NV_Ith_S(nv(m_y), i) = 0.0;
303  }
304  // check abs tolerance array size
305  if (m_itol == CV_SV && m_nabs < m_neq) {
306  throw CVodesErr("not enough absolute tolerance values specified.");
307  }
308 
309  func.getInitialConditions(m_t0, m_neq, NV_DATA_S(nv(m_y)));
310 
311  if (m_cvode_mem) {
312  CVodeFree(&m_cvode_mem);
313  }
314 
315  /*
316  * Specify the method and the iteration type:
317  * Cantera Defaults:
318  * CV_BDF - Use BDF methods
319  * CV_NEWTON - use newton's method
320  */
321  m_cvode_mem = CVodeCreate(m_method, m_iter);
322  if (!m_cvode_mem) {
323  throw CVodesErr("CVodeCreate failed.");
324  }
325 
326  int flag = 0;
327 #if SUNDIALS_VERSION <= 23
328  if (m_itol == CV_SV) {
329  // vector atol
330  flag = CVodeMalloc(m_cvode_mem, cvodes_rhs, m_t0, nv(m_y), m_itol,
331  m_reltol, nv(m_abstol));
332  } else {
333  // scalar atol
334  flag = CVodeMalloc(m_cvode_mem, cvodes_rhs, m_t0, nv(m_y), m_itol,
335  m_reltol, &m_abstols);
336  }
337  if (flag != CV_SUCCESS) {
338  if (flag == CV_MEM_FAIL) {
339  throw CVodesErr("Memory allocation failed.");
340  } else if (flag == CV_ILL_INPUT) {
341  throw CVodesErr("Illegal value for CVodeMalloc input argument.");
342  } else {
343  throw CVodesErr("CVodeMalloc failed.");
344  }
345  }
346 #elif SUNDIALS_VERSION >= 24
347 
348  flag = CVodeInit(m_cvode_mem, cvodes_rhs, m_t0, nv(m_y));
349  if (flag != CV_SUCCESS) {
350  if (flag == CV_MEM_FAIL) {
351  throw CVodesErr("Memory allocation failed.");
352  } else if (flag == CV_ILL_INPUT) {
353  throw CVodesErr("Illegal value for CVodeInit input argument.");
354  } else {
355  throw CVodesErr("CVodeInit failed.");
356  }
357  }
358 
359  if (m_itol == CV_SV) {
360  flag = CVodeSVtolerances(m_cvode_mem, m_reltol, nv(m_abstol));
361  } else {
362  flag = CVodeSStolerances(m_cvode_mem, m_reltol, m_abstols);
363  }
364  if (flag != CV_SUCCESS) {
365  if (flag == CV_MEM_FAIL) {
366  throw CVodesErr("Memory allocation failed.");
367  } else if (flag == CV_ILL_INPUT) {
368  throw CVodesErr("Illegal value for CVodeInit input argument.");
369  } else {
370  throw CVodesErr("CVodeInit failed.");
371  }
372  }
373 #endif
374 
375  if (m_type == DENSE + NOJAC) {
376  long int N = m_neq;
377  CVDense(m_cvode_mem, N);
378  } else if (m_type == DIAG) {
379  CVDiag(m_cvode_mem);
380  } else if (m_type == GMRES) {
381  CVSpgmr(m_cvode_mem, PREC_NONE, 0);
382  } else if (m_type == BAND + NOJAC) {
383  long int N = m_neq;
384  long int nu = m_mupper;
385  long int nl = m_mlower;
386  CVBand(m_cvode_mem, N, nu, nl);
387  } else {
388  throw CVodesErr("unsupported option");
389  }
390 
391  // pass a pointer to func in m_data
392  delete m_fdata;
393  m_fdata = new FuncData(&func, func.nparams());
394 
395  //m_data = (void*)&func;
396 #if SUNDIALS_VERSION <= 23
397  flag = CVodeSetFdata(m_cvode_mem, (void*)m_fdata);
398  if (flag != CV_SUCCESS) {
399  throw CVodesErr("CVodeSetFdata failed.");
400  }
401 #elif SUNDIALS_VERSION >= 24
402  flag = CVodeSetUserData(m_cvode_mem, (void*)m_fdata);
403  if (flag != CV_SUCCESS) {
404  throw CVodesErr("CVodeSetUserData failed.");
405  }
406 #endif
407  if (func.nparams() > 0) {
408  sensInit(t0, func);
409  flag = CVodeSetSensParams(m_cvode_mem, DATA_PTR(m_fdata->m_pars),
410  NULL, NULL);
411  }
412 
413  // set options
414  if (m_maxord > 0) {
415  flag = CVodeSetMaxOrd(m_cvode_mem, m_maxord);
416  }
417  if (m_maxsteps > 0) {
418  flag = CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
419  }
420  if (m_hmax > 0) {
421  flag = CVodeSetMaxStep(m_cvode_mem, m_hmax);
422  }
423 }
424 
425 
426 void CVodesIntegrator::reinitialize(double t0, FuncEval& func)
427 {
428  m_t0 = t0;
429  //try {
430  func.getInitialConditions(m_t0, m_neq, NV_DATA_S(nv(m_y)));
431  //}
432  //catch (CanteraError) {
433  //showErrors();
434  //error("Teminating execution");
435  //}
436 
437  int result;
438 
439 #if SUNDIALS_VERSION <= 23
440  if (m_itol == CV_SV) {
441  result = CVodeReInit(m_cvode_mem, cvodes_rhs, m_t0, nv(m_y),
442  m_itol, m_reltol,
443  nv(m_abstol));
444  } else {
445  result = CVodeReInit(m_cvode_mem, cvodes_rhs, m_t0, nv(m_y),
446  m_itol, m_reltol,
447  &m_abstols);
448  }
449  if (result != CV_SUCCESS) {
450  throw CVodesErr("CVodeReInit failed. result = "+int2str(result));
451  }
452 #elif SUNDIALS_VERSION >= 24
453  result = CVodeReInit(m_cvode_mem, m_t0, nv(m_y));
454  if (result != CV_SUCCESS) {
455  throw CVodesErr("CVodeReInit failed. result = "+int2str(result));
456  }
457 #endif
458 
459  if (m_type == DENSE + NOJAC) {
460  long int N = m_neq;
461  CVDense(m_cvode_mem, N);
462  } else if (m_type == DIAG) {
463  CVDiag(m_cvode_mem);
464  } else if (m_type == BAND + NOJAC) {
465  long int N = m_neq;
466  long int nu = m_mupper;
467  long int nl = m_mlower;
468  CVBand(m_cvode_mem, N, nu, nl);
469  } else if (m_type == GMRES) {
470  CVSpgmr(m_cvode_mem, PREC_NONE, 0);
471  } else {
472  throw CVodesErr("unsupported option");
473  }
474 
475 
476  // set options
477  if (m_maxord > 0) {
478  CVodeSetMaxOrd(m_cvode_mem, m_maxord);
479  }
480  if (m_maxsteps > 0) {
481  CVodeSetMaxNumSteps(m_cvode_mem, m_maxsteps);
482  }
483  if (m_hmax > 0) {
484  CVodeSetMaxStep(m_cvode_mem, m_hmax);
485  }
486 }
487 
488 void CVodesIntegrator::integrate(double tout)
489 {
490  double t;
491  int flag;
492  flag = CVode(m_cvode_mem, tout, nv(m_y), &t, CV_NORMAL);
493  if (flag != CV_SUCCESS) {
494  throw CVodesErr(" CVodes error encountered. Error code: " + int2str(flag));
495  }
496 #if SUNDIALS_VERSION <= 23
497  if (m_np > 0) {
498  CVodeGetSens(m_cvode_mem, tout, m_yS);
499  }
500 #elif SUNDIALS_VERSION >= 24
501  double tretn;
502  if (m_np > 0) {
503  CVodeGetSens(m_cvode_mem, &tretn, m_yS);
504  if (fabs(tretn - tout) > 1.0E-5) {
505  throw CVodesErr("Time of Sensitivities different than time of tout");
506  }
507  }
508 #endif
509 }
510 
511 double CVodesIntegrator::step(double tout)
512 {
513  double t;
514  int flag;
515  flag = CVode(m_cvode_mem, tout, nv(m_y), &t, CV_ONE_STEP);
516  if (flag != CV_SUCCESS) {
517  throw CVodesErr(" CVodes error encountered. Error code: " + int2str(flag));
518  }
519  return t;
520 }
521 
522 int CVodesIntegrator::nEvals() const
523 {
524  long int ne;
525  CVodeGetNumRhsEvals(m_cvode_mem, &ne);
526  return ne;
527  //return m_iopt[NFE];
528 }
529 
530 double CVodesIntegrator::sensitivity(size_t k, size_t p)
531 {
532  if (k >= m_neq) {
533  throw CVodesErr("sensitivity: k out of range ("+int2str(p)+")");
534  }
535  if (p >= m_np) {
536  throw CVodesErr("sensitivity: p out of range ("+int2str(p)+")");
537  }
538  return NV_Ith_S(m_yS[p],k);
539 }
540 }
541 
542