Cantera  2.1.2
IDA_Solver.cpp
Go to the documentation of this file.
1 /**
2  * @file IDA_Solver.cpp
3  */
4 
5 // Copyright 2006 California Institute of Technology
6 
9 
10 #if SUNDIALS_VERSION >= 24
11 #include "sundials/sundials_types.h"
12 #include "sundials/sundials_math.h"
13 #include "ida/ida.h"
14 #include "ida/ida_dense.h"
15 #include "ida/ida_spgmr.h"
16 #include "ida/ida_band.h"
17 #include "nvector/nvector_serial.h"
18 
19 using namespace std;
20 
21 #if SUNDIALS_VERSION < 25
22 typedef int sd_size_t;
23 #else
24 typedef long int sd_size_t;
25 #endif
26 
27 namespace Cantera
28 {
29 
30 /**
31  * A simple class to hold an array of parameter values and a pointer to
32  * an instance of a subclass of ResidEval.
33  */
34 class ResidData
35 {
36 
37 public:
38 
39  ResidData(ResidJacEval* f, IDA_Solver* s, int npar = 0) {
40  m_func = f;
41  m_solver = s;
42  }
43 
44  virtual ~ResidData() {
45  }
46 
47  ResidJacEval* m_func;
48  IDA_Solver* m_solver;
49 };
50 }
51 
52 extern "C" {
53  //! Function called by IDA to evaluate the residual, given y and ydot.
54  /*!
55  * IDA allows passing in a void* pointer to access external data. Instead of requiring the user to provide a
56  * residual function directly to IDA (which would require using
57  * the sundials data types N_Vector, etc.), we define this function as the single function that IDA always calls. The
58  * real evaluation of the residual is done by an instance of a subclass of ResidEval, passed in to this
59  * function as a pointer in the parameters.
60  *
61  * FROM IDA WRITEUP -> What the IDA solver expects as a return flag from its residual routines ------
62  * A IDAResFn res should return a value of 0 if successful, a positive
63  * value if a recoverable error occured (e.g. yy has an illegal value),
64  * or a negative value if a nonrecoverable error occured. In the latter
65  * case, the program halts. If a recoverable error occured, the integrator
66  * will attempt to correct and retry.
67  */
68  static int ida_resid(realtype t, N_Vector y, N_Vector ydot, N_Vector r, void* f_data)
69  {
70  double* ydata = NV_DATA_S(y);
71  double* ydotdata = NV_DATA_S(ydot);
72  double* rdata = NV_DATA_S(r);
74  Cantera::ResidJacEval* f = d->m_func;
75  Cantera::IDA_Solver* s = d->m_solver;
76  double delta_t = s->getCurrentStepFromIDA();
77  // TODO evaluate evalType. Assumed to be Base_ResidEval
78  int retn = 0;
79  int flag = f->evalResidNJ(t, delta_t, ydata, ydotdata, rdata);
80  if (flag < 0) {
81  // This signals to IDA that a nonrecoverable error has occurred.
82  retn = flag;
83  }
84  return retn;
85  }
86 
87  //! Function called by by IDA to evaluate the Jacobian, given y and ydot.
88  /*!
89  * typedef int (*IDADlsDenseJacFn)(sd_size_t N, realtype t, realtype c_j,
90  * N_Vector y, N_Vector yp, N_Vector r,
91  * DlsMat Jac, void *user_data,
92  * N_Vector tmp1, N_Vector tmp2, N_Vector tmp3);
93  *
94  * A IDADlsDenseJacFn should return
95  * 0 if successful,
96  * a positive int if a recoverable error occurred, or
97  * a negative int if a nonrecoverable error occurred.
98  * In the case of a recoverable error return, the integrator will
99  * attempt to recover by reducing the stepsize (which changes cj).
100  */
101  static int ida_jacobian(sd_size_t nrows, realtype t, realtype c_j, N_Vector y, N_Vector ydot, N_Vector r,
102  DlsMat Jac, void* f_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
103  {
104  doublereal* ydata = NV_DATA_S(y);
105  doublereal* ydotdata = NV_DATA_S(ydot);
106  doublereal* rdata = NV_DATA_S(r);
107  Cantera::ResidData* d = (Cantera::ResidData*) f_data;
108  Cantera::ResidJacEval* f = d->m_func;
109  doublereal* const* colPts = Jac->cols;
110  Cantera::IDA_Solver* s = d->m_solver;
111  double delta_t = s->getCurrentStepFromIDA();
112  // printf(" delta_t = %g 1/cj = %g\n", delta_t, 1.0/c_j);
113  f->evalJacobianDP(t, delta_t, c_j, ydata, ydotdata, colPts, rdata);
114  return 0;
115  }
116 }
117 
118 namespace Cantera
119 {
120 
121 IDA_Solver::IDA_Solver(ResidJacEval& f) :
122  DAE_Solver(f),
123  m_ida_mem(0),
124  m_t0(0.0),
125  m_y(0),
126  m_ydot(0),
127  m_id(0),
128  m_constraints(0),
129  m_abstol(0),
130  m_type(0),
131  m_itol(IDA_SS),
132  m_iter(0),
133  m_reltol(1.e-9),
134  m_abstols(1.e-15),
135  m_nabs(0),
136  m_hmax(0.0),
137  m_hmin(0.0),
138  m_h0(0.0),
139  m_maxsteps(20000),
140  m_maxord(0),
141  m_formJac(0),
142  m_tstop(0.0),
143  m_told_old(0.0),
144  m_told(0.0),
145  m_tcurrent(0.0),
146  m_deltat(0.0),
147  m_maxErrTestFails(-1),
148  m_maxNonlinIters(0),
149  m_maxNonlinConvFails(-1),
150  m_setSuppressAlg(0),
151  m_fdata(0),
152  m_mupper(0),
153  m_mlower(0)
154 {
155 }
156 
157 IDA_Solver::~IDA_Solver()
158 {
159  if (m_ida_mem) {
160  IDAFree(&m_ida_mem);
161  }
162  if (m_y) {
163  N_VDestroy_Serial(m_y);
164  }
165  if (m_ydot) {
166  N_VDestroy_Serial(m_ydot);
167  }
168  if (m_abstol) {
169  N_VDestroy_Serial(m_abstol);
170  }
171  if (m_constraints) {
172  N_VDestroy_Serial(m_constraints);
173  }
174  delete m_fdata;
175 }
176 
177 doublereal IDA_Solver::solution(int k) const
178 {
179  return NV_Ith_S(m_y,k);
180 }
181 
182 const doublereal* IDA_Solver::solutionVector() const
183 {
184  return NV_DATA_S(m_y);
185 }
186 
187 doublereal IDA_Solver::derivative(int k) const
188 {
189  return NV_Ith_S(m_ydot,k);
190 }
191 
192 const doublereal* IDA_Solver::derivativeVector() const
193 {
194  return NV_DATA_S(m_ydot);
195 }
196 
197 void IDA_Solver::setTolerances(double reltol, double* abstol)
198 {
199  m_itol = IDA_SV;
200  if (!m_abstol) {
201  m_abstol = N_VNew_Serial(m_neq);
202  }
203  for (int i = 0; i < m_neq; i++) {
204  NV_Ith_S(m_abstol, i) = abstol[i];
205  }
206  m_reltol = reltol;
207  if (m_ida_mem) {
208  int flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
209  if (flag != IDA_SUCCESS) {
210  throw IDA_Err("Memory allocation failed.");
211  }
212  }
213 }
214 
215 void IDA_Solver::setTolerances(doublereal reltol, doublereal abstol)
216 {
217  m_itol = IDA_SS;
218  m_reltol = reltol;
219  m_abstols = abstol;
220  if (m_ida_mem) {
221  int flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
222  if (flag != IDA_SUCCESS) {
223  throw IDA_Err("Memory allocation failed.");
224  }
225  }
226 }
227 
228 void IDA_Solver::setLinearSolverType(int solverType)
229 {
230  m_type = solverType;
231 }
232 
233 void IDA_Solver::setDenseLinearSolver()
234 {
235  setLinearSolverType(0);
236 }
237 
238 void IDA_Solver::setBandedLinearSolver(int m_upper, int m_lower)
239 {
240  m_type = 2;
241  m_upper = m_mupper;
242  m_mlower = m_lower;
243 }
244 
245 void IDA_Solver::setMaxOrder(int n)
246 {
247  m_maxord = n;
248 }
249 
250 void IDA_Solver::setMaxNumSteps(int n)
251 {
252  m_maxsteps = n;
253 }
254 
255 void IDA_Solver::setInitialStepSize(doublereal h0)
256 {
257  m_h0 = h0;
258 }
259 
260 void IDA_Solver::setStopTime(doublereal tstop)
261 {
262  m_tstop = tstop;
263 }
264 
265 doublereal IDA_Solver::getCurrentStepFromIDA()
266 {
267  doublereal hcur;
268  IDAGetCurrentStep(m_ida_mem, &hcur);
269  return hcur;
270 }
271 
272 void IDA_Solver::setJacobianType(int formJac)
273 {
274  m_formJac = formJac;
275  if (m_ida_mem) {
276  if (m_formJac == 1) {
277  int flag = IDADlsSetDenseJacFn(m_ida_mem, ida_jacobian);
278  if (flag != IDA_SUCCESS) {
279  throw IDA_Err("IDADlsSetDenseJacFn failed.");
280  }
281  }
282  }
283 }
284 
285 void IDA_Solver::setMaxErrTestFailures(int maxErrTestFails)
286 {
287  m_maxErrTestFails = maxErrTestFails;
288 }
289 
290 void IDA_Solver::setMaxNonlinIterations(int n)
291 {
292  m_maxNonlinIters = n;
293 }
294 
295 void IDA_Solver::setMaxNonlinConvFailures(int n)
296 {
297  m_maxNonlinConvFails = n;
298 }
299 
300 void IDA_Solver::inclAlgebraicInErrorTest(bool yesno)
301 {
302  if (yesno) {
303  m_setSuppressAlg = 0;
304  } else {
305  m_setSuppressAlg = 1;
306  }
307 }
308 
309 void IDA_Solver::init(doublereal t0)
310 {
311  m_t0 = t0;
312  m_told = t0;
313  m_told_old = t0;
314  m_tcurrent = t0;
315  if (m_y) {
316  N_VDestroy_Serial(m_y);
317  }
318  if (m_ydot) {
319  N_VDestroy_Serial(m_ydot);
320  }
321  if (m_id) {
322  N_VDestroy_Serial(m_id);
323  }
324  if (m_constraints) {
325  N_VDestroy_Serial(m_constraints);
326  }
327 
328  m_y = N_VNew_Serial(m_neq);
329  m_ydot = N_VNew_Serial(m_neq);
330  m_constraints = N_VNew_Serial(m_neq);
331 
332  for (int i=0; i<m_neq; i++) {
333  NV_Ith_S(m_y, i) = 0.0;
334  NV_Ith_S(m_ydot, i) = 0.0;
335  NV_Ith_S(m_constraints, i) = 0.0;
336  }
337 
338  // get the initial conditions
339  m_resid.getInitialConditions(m_t0, NV_DATA_S(m_y), NV_DATA_S(m_ydot));
340 
341  if (m_ida_mem) {
342  IDAFree(&m_ida_mem);
343  }
344 
345  /* Call IDACreate */
346  m_ida_mem = IDACreate();
347 
348  int flag = 0;
349 
350 
351 
352  if (m_itol == IDA_SV) {
353 #if SUNDIALS_VERSION <= 23
354  // vector atol
355  flag = IDAMalloc(m_ida_mem, ida_resid, m_t0, m_y, m_ydot,
356  m_itol, m_reltol, m_abstol);
357  if (flag != IDA_SUCCESS) {
358  if (flag == IDA_MEM_FAIL) {
359  throw IDA_Err("Memory allocation failed.");
360  } else if (flag == IDA_ILL_INPUT) {
361  throw IDA_Err("Illegal value for IDAMalloc input argument.");
362  } else {
363  throw IDA_Err("IDAMalloc failed.");
364  }
365  }
366 
367 #elif SUNDIALS_VERSION >= 24
368  flag = IDAInit(m_ida_mem, ida_resid, m_t0, m_y, m_ydot);
369  if (flag != IDA_SUCCESS) {
370  if (flag == IDA_MEM_FAIL) {
371  throw IDA_Err("Memory allocation failed.");
372  } else if (flag == IDA_ILL_INPUT) {
373  throw IDA_Err("Illegal value for IDAMalloc input argument.");
374  } else {
375  throw IDA_Err("IDAMalloc failed.");
376  }
377  }
378  flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
379  if (flag != IDA_SUCCESS) {
380  throw IDA_Err("Memory allocation failed.");
381  }
382 #endif
383  } else {
384 #if SUNDIALS_VERSION <= 23
385  // scalar atol
386  flag = IDAMalloc(m_ida_mem, ida_resid, m_t0, m_y, m_ydot,
387  m_itol, m_reltol, &m_abstols);
388  if (flag != IDA_SUCCESS) {
389  if (flag == IDA_MEM_FAIL) {
390  throw IDA_Err("Memory allocation failed.");
391  } else if (flag == IDA_ILL_INPUT) {
392  throw IDA_Err("Illegal value for IDAMalloc input argument.");
393  } else {
394  throw IDA_Err("IDAMalloc failed.");
395  }
396  }
397 
398 #elif SUNDIALS_VERSION >= 24
399  flag = IDAInit(m_ida_mem, ida_resid, m_t0, m_y, m_ydot);
400  if (flag != IDA_SUCCESS) {
401  if (flag == IDA_MEM_FAIL) {
402  throw IDA_Err("Memory allocation failed.");
403  } else if (flag == IDA_ILL_INPUT) {
404  throw IDA_Err("Illegal value for IDAMalloc input argument.");
405  } else {
406  throw IDA_Err("IDAMalloc failed.");
407  }
408  }
409  flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
410  if (flag != IDA_SUCCESS) {
411  throw IDA_Err("Memory allocation failed.");
412  }
413 #endif
414  }
415 
416  //-----------------------------------
417  // set the linear solver type
418  //-----------------------------------
419 
420  if (m_type == 1 || m_type == 0) {
421  long int N = m_neq;
422  flag = IDADense(m_ida_mem, N);
423  if (flag) {
424  throw IDA_Err("IDADense failed");
425  }
426  } else if (m_type == 2) {
427  long int N = m_neq;
428  long int nu = m_mupper;
429  long int nl = m_mlower;
430  IDABand(m_ida_mem, N, nu, nl);
431  } else {
432  throw IDA_Err("unsupported linear solver type");
433  }
434 
435  if (m_formJac == 1) {
436  flag = IDADlsSetDenseJacFn(m_ida_mem, ida_jacobian);
437  if (flag != IDA_SUCCESS) {
438  throw IDA_Err("IDADlsSetDenseJacFn failed.");
439  }
440  }
441 
442  // pass a pointer to func in m_data
443  m_fdata = new ResidData(&m_resid, this, m_resid.nparams());
444 #if SUNDIALS_VERSION <= 23
445  flag = IDASetRdata(m_ida_mem, (void*)m_fdata);
446  if (flag != IDA_SUCCESS) {
447  throw IDA_Err("IDASetRdata failed.");
448  }
449 #elif SUNDIALS_VERSION >= 24
450  flag = IDASetUserData(m_ida_mem, (void*)m_fdata);
451  if (flag != IDA_SUCCESS) {
452  throw IDA_Err("IDASetUserData failed.");
453  }
454 #endif
455 
456  // set options
457  if (m_maxord > 0) {
458  flag = IDASetMaxOrd(m_ida_mem, m_maxord);
459  if (flag != IDA_SUCCESS) {
460  throw IDA_Err("IDASetMaxOrd failed.");
461  }
462  }
463  if (m_maxsteps > 0) {
464  flag = IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
465  if (flag != IDA_SUCCESS) {
466  throw IDA_Err("IDASetMaxNumSteps failed.");
467  }
468  }
469  if (m_h0 > 0.0) {
470  flag = IDASetInitStep(m_ida_mem, m_h0);
471  if (flag != IDA_SUCCESS) {
472  throw IDA_Err("IDASetInitStep failed.");
473  }
474  }
475  if (m_tstop > 0.0) {
476  flag = IDASetStopTime(m_ida_mem, m_tstop);
477  if (flag != IDA_SUCCESS) {
478  throw IDA_Err("IDASetStopTime failed.");
479  }
480  }
481  if (m_maxErrTestFails >= 0) {
482  flag = IDASetMaxErrTestFails(m_ida_mem, m_maxErrTestFails);
483  if (flag != IDA_SUCCESS) {
484  throw IDA_Err("IDASetMaxErrTestFails failed.");
485  }
486  }
487  if (m_maxNonlinIters > 0) {
488  flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
489  if (flag != IDA_SUCCESS) {
490  throw IDA_Err("IDASetmaxNonlinIters failed.");
491  }
492  }
493  if (m_maxNonlinConvFails >= 0) {
494  flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
495  if (flag != IDA_SUCCESS) {
496  throw IDA_Err("IDASetMaxConvFails failed.");
497  }
498  }
499  if (m_setSuppressAlg != 0) {
500  flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
501  if (flag != IDA_SUCCESS) {
502  throw IDA_Err("IDASetSuppressAlg failed.");
503  }
504  }
505 }
506 
507 void IDA_Solver::correctInitial_Y_given_Yp(doublereal* y, doublereal* yp, doublereal tout)
508 {
509  int icopt = IDA_Y_INIT;
510  doublereal tout1 = tout;
511  if (tout == 0.0) {
512  double h0 = 1.0E-5;
513  if (m_h0 > 0.0) {
514  h0 = m_h0;
515  }
516  tout1 = m_t0 + h0;
517  }
518 
519  int flag = IDACalcIC(m_ida_mem, icopt, tout1);
520  if (flag != IDA_SUCCESS) {
521  throw IDA_Err("IDACalcIC failed: error = " + int2str(flag));
522  }
523 
524 
525  flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
526  if (flag != IDA_SUCCESS) {
527  throw IDA_Err("IDAGetSolution failed: error = " + int2str(flag));
528  }
529  doublereal* yy = NV_DATA_S(m_y);
530  doublereal* yyp = NV_DATA_S(m_ydot);
531 
532  for (int i = 0; i < m_neq; i++) {
533  y[i] = yy[i];
534  yp[i] = yyp[i];
535  }
536 }
537 
538 void IDA_Solver::correctInitial_YaYp_given_Yd(doublereal* y, doublereal* yp, doublereal tout)
539 {
540  int icopt = IDA_YA_YDP_INIT;
541  doublereal tout1 = tout;
542  if (tout == 0.0) {
543  double h0 = 1.0E-5;
544  if (m_h0 > 0.0) {
545  h0 = m_h0;
546  }
547  tout1 = m_t0 + h0;
548  }
549 
550  int flag = IDACalcIC(m_ida_mem, icopt, tout1);
551  if (flag != IDA_SUCCESS) {
552  throw IDA_Err("IDACalcIC failed: error = " + int2str(flag));
553  }
554 
555 
556  flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
557  if (flag != IDA_SUCCESS) {
558  throw IDA_Err("IDAGetSolution failed: error = " + int2str(flag));
559  }
560  doublereal* yy = NV_DATA_S(m_y);
561  doublereal* yyp = NV_DATA_S(m_ydot);
562 
563  for (int i = 0; i < m_neq; i++) {
564  y[i] = yy[i];
565  yp[i] = yyp[i];
566  }
567 }
568 
569 int IDA_Solver::solve(double tout)
570 {
571  double tretn;
572  int flag;
573  flag = IDASetStopTime(m_ida_mem, tout);
574  if (flag != IDA_SUCCESS) {
575  throw IDA_Err(" IDA error encountered.");
576  }
577  do {
578  if (tout <= m_tcurrent) {
579  throw IDA_Err(" tout <= tcurrent");
580  }
581  m_told_old = m_told;
582  m_told = m_tcurrent;
583  flag = IDASolve(m_ida_mem, tout, &tretn, m_y, m_ydot, IDA_ONE_STEP);
584  if (flag < 0) {
585  throw IDA_Err(" IDA error encountered.");
586  } else if (flag == IDA_TSTOP_RETURN) {
587  // we've reached our goal, and have actually integrated past it
588  } else if (flag == IDA_ROOT_RETURN) {
589  // not sure what to do with this yet
590  } else if (flag == IDA_WARNING) {
591  throw IDA_Err(" IDA Warning encountered.");
592  }
593  m_tcurrent = tretn;
594  m_deltat = m_tcurrent - m_told;
595  } while (tretn < tout);
596 
597  if (flag != IDA_SUCCESS && flag != IDA_TSTOP_RETURN) {
598  throw IDA_Err(" IDA error encountered.");
599  }
600  return flag;
601 }
602 
603 double IDA_Solver::step(double tout)
604 {
605  double t;
606  int flag;
607  if (tout <= m_tcurrent) {
608  throw IDA_Err(" tout <= tcurrent");
609  }
610  m_told_old = m_told;
611  m_told = m_tcurrent;
612  flag = IDASolve(m_ida_mem, tout, &t, m_y, m_ydot, IDA_ONE_STEP);
613  if (flag < 0) {
614  throw IDA_Err(" IDA error encountered.");
615  } else if (flag == IDA_TSTOP_RETURN) {
616  // we've reached our goal, and have actually integrated past it
617  } else if (flag == IDA_ROOT_RETURN) {
618  // not sure what to do with this yet
619  } else if (flag == IDA_WARNING) {
620  throw IDA_Err(" IDA Warning encountered.");
621  }
622  m_tcurrent = t;
623  m_deltat = m_tcurrent - m_told;
624  return t;
625 }
626 
627 doublereal IDA_Solver::getOutputParameter(int flag) const
628 {
629  long int lenrw, leniw;
630  switch (flag) {
631  case REAL_WORKSPACE_SIZE:
632  flag = IDAGetWorkSpace(m_ida_mem, &lenrw, &leniw);
633  return doublereal(lenrw);
634  break;
635  }
636  return 0.0;
637 }
638 
639 }
640 #endif
virtual int evalJacobianDP(const doublereal t, const doublereal delta_t, doublereal cj, const doublereal *const y, const doublereal *const ydot, doublereal *const *jacobianColPts, doublereal *const resid)
Calculate an analytical jacobian and the residual at the current time and values. ...
A simple class to hold an array of parameter values and a pointer to an instance of a subclass of Res...
Definition: IDA_Solver.cpp:34
std::string int2str(const int n, const std::string &fmt)
Convert an int to a string using a format converter.
Definition: stringUtils.cpp:40
virtual int evalResidNJ(const doublereal t, const doublereal delta_t, const doublereal *const y, const doublereal *const ydot, doublereal *const resid, const ResidEval_Type_Enum evalType=Base_ResidEval, const int id_x=-1, const doublereal delta_x=0.0)
Evaluate the residual function.
int solve(DenseMatrix &A, double *b)
Solve Ax = b. Array b is overwritten on exit with x.
static int ida_resid(realtype t, N_Vector y, N_Vector ydot, N_Vector r, void *f_data)
Function called by IDA to evaluate the residual, given y and ydot.
Definition: IDA_Solver.cpp:68
Wrappers for the function evaluators for Nonlinear solvers and Time steppers.
Definition: ResidJacEval.h:51
Contains declarations for string manipulation functions within Cantera.
Header file for class IDA_Solver.
static int ida_jacobian(sd_size_t nrows, realtype t, realtype c_j, N_Vector y, N_Vector ydot, N_Vector r, DlsMat Jac, void *f_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Function called by by IDA to evaluate the Jacobian, given y and ydot.
Definition: IDA_Solver.cpp:101