Cantera  2.2.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 HAS_SUNDIALS
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  f->evalJacobianDP(t, delta_t, c_j, ydata, ydotdata, colPts, rdata);
113  return 0;
114  }
115 }
116 
117 namespace Cantera
118 {
119 
120 IDA_Solver::IDA_Solver(ResidJacEval& f) :
121  DAE_Solver(f),
122  m_ida_mem(0),
123  m_t0(0.0),
124  m_y(0),
125  m_ydot(0),
126  m_id(0),
127  m_constraints(0),
128  m_abstol(0),
129  m_type(0),
130  m_itol(IDA_SS),
131  m_iter(0),
132  m_reltol(1.e-9),
133  m_abstols(1.e-15),
134  m_nabs(0),
135  m_hmax(0.0),
136  m_hmin(0.0),
137  m_h0(0.0),
138  m_maxsteps(20000),
139  m_maxord(0),
140  m_formJac(0),
141  m_tstop(0.0),
142  m_told_old(0.0),
143  m_told(0.0),
144  m_tcurrent(0.0),
145  m_deltat(0.0),
146  m_maxErrTestFails(-1),
147  m_maxNonlinIters(0),
148  m_maxNonlinConvFails(-1),
149  m_setSuppressAlg(0),
150  m_fdata(0),
151  m_mupper(0),
152  m_mlower(0)
153 {
154 }
155 
156 IDA_Solver::~IDA_Solver()
157 {
158  if (m_ida_mem) {
159  IDAFree(&m_ida_mem);
160  }
161  if (m_y) {
162  N_VDestroy_Serial(m_y);
163  }
164  if (m_ydot) {
165  N_VDestroy_Serial(m_ydot);
166  }
167  if (m_abstol) {
168  N_VDestroy_Serial(m_abstol);
169  }
170  if (m_constraints) {
171  N_VDestroy_Serial(m_constraints);
172  }
173  delete m_fdata;
174 }
175 
176 doublereal IDA_Solver::solution(int k) const
177 {
178  return NV_Ith_S(m_y,k);
179 }
180 
181 const doublereal* IDA_Solver::solutionVector() const
182 {
183  return NV_DATA_S(m_y);
184 }
185 
186 doublereal IDA_Solver::derivative(int k) const
187 {
188  return NV_Ith_S(m_ydot,k);
189 }
190 
191 const doublereal* IDA_Solver::derivativeVector() const
192 {
193  return NV_DATA_S(m_ydot);
194 }
195 
196 void IDA_Solver::setTolerances(double reltol, double* abstol)
197 {
198  m_itol = IDA_SV;
199  if (!m_abstol) {
200  m_abstol = N_VNew_Serial(m_neq);
201  }
202  for (int i = 0; i < m_neq; i++) {
203  NV_Ith_S(m_abstol, i) = abstol[i];
204  }
205  m_reltol = reltol;
206  if (m_ida_mem) {
207  int flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
208  if (flag != IDA_SUCCESS) {
209  throw IDA_Err("Memory allocation failed.");
210  }
211  }
212 }
213 
214 void IDA_Solver::setTolerances(doublereal reltol, doublereal abstol)
215 {
216  m_itol = IDA_SS;
217  m_reltol = reltol;
218  m_abstols = abstol;
219  if (m_ida_mem) {
220  int flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
221  if (flag != IDA_SUCCESS) {
222  throw IDA_Err("Memory allocation failed.");
223  }
224  }
225 }
226 
227 void IDA_Solver::setLinearSolverType(int solverType)
228 {
229  m_type = solverType;
230 }
231 
232 void IDA_Solver::setDenseLinearSolver()
233 {
234  setLinearSolverType(0);
235 }
236 
237 void IDA_Solver::setBandedLinearSolver(int m_upper, int m_lower)
238 {
239  m_type = 2;
240  m_upper = m_mupper;
241  m_mlower = m_lower;
242 }
243 
244 void IDA_Solver::setMaxOrder(int n)
245 {
246  m_maxord = n;
247 }
248 
249 void IDA_Solver::setMaxNumSteps(int n)
250 {
251  m_maxsteps = n;
252 }
253 
254 void IDA_Solver::setInitialStepSize(doublereal h0)
255 {
256  m_h0 = h0;
257 }
258 
259 void IDA_Solver::setStopTime(doublereal tstop)
260 {
261  m_tstop = tstop;
262 }
263 
264 doublereal IDA_Solver::getCurrentStepFromIDA()
265 {
266  doublereal hcur;
267  IDAGetCurrentStep(m_ida_mem, &hcur);
268  return hcur;
269 }
270 
271 void IDA_Solver::setJacobianType(int formJac)
272 {
273  m_formJac = formJac;
274  if (m_ida_mem) {
275  if (m_formJac == 1) {
276  int flag = IDADlsSetDenseJacFn(m_ida_mem, ida_jacobian);
277  if (flag != IDA_SUCCESS) {
278  throw IDA_Err("IDADlsSetDenseJacFn failed.");
279  }
280  }
281  }
282 }
283 
284 void IDA_Solver::setMaxErrTestFailures(int maxErrTestFails)
285 {
286  m_maxErrTestFails = maxErrTestFails;
287 }
288 
289 void IDA_Solver::setMaxNonlinIterations(int n)
290 {
291  m_maxNonlinIters = n;
292 }
293 
294 void IDA_Solver::setMaxNonlinConvFailures(int n)
295 {
296  m_maxNonlinConvFails = n;
297 }
298 
299 void IDA_Solver::inclAlgebraicInErrorTest(bool yesno)
300 {
301  if (yesno) {
302  m_setSuppressAlg = 0;
303  } else {
304  m_setSuppressAlg = 1;
305  }
306 }
307 
308 void IDA_Solver::init(doublereal t0)
309 {
310  m_t0 = t0;
311  m_told = t0;
312  m_told_old = t0;
313  m_tcurrent = t0;
314  if (m_y) {
315  N_VDestroy_Serial(m_y);
316  }
317  if (m_ydot) {
318  N_VDestroy_Serial(m_ydot);
319  }
320  if (m_id) {
321  N_VDestroy_Serial(m_id);
322  }
323  if (m_constraints) {
324  N_VDestroy_Serial(m_constraints);
325  }
326 
327  m_y = N_VNew_Serial(m_neq);
328  m_ydot = N_VNew_Serial(m_neq);
329  m_constraints = N_VNew_Serial(m_neq);
330 
331  for (int i=0; i<m_neq; i++) {
332  NV_Ith_S(m_y, i) = 0.0;
333  NV_Ith_S(m_ydot, i) = 0.0;
334  NV_Ith_S(m_constraints, i) = 0.0;
335  }
336 
337  // get the initial conditions
338  m_resid.getInitialConditions(m_t0, NV_DATA_S(m_y), NV_DATA_S(m_ydot));
339 
340  if (m_ida_mem) {
341  IDAFree(&m_ida_mem);
342  }
343 
344  /* Call IDACreate */
345  m_ida_mem = IDACreate();
346 
347  int flag = 0;
348 
349 
350 
351  if (m_itol == IDA_SV) {
352  flag = IDAInit(m_ida_mem, ida_resid, m_t0, m_y, m_ydot);
353  if (flag != IDA_SUCCESS) {
354  if (flag == IDA_MEM_FAIL) {
355  throw IDA_Err("Memory allocation failed.");
356  } else if (flag == IDA_ILL_INPUT) {
357  throw IDA_Err("Illegal value for IDAMalloc input argument.");
358  } else {
359  throw IDA_Err("IDAMalloc failed.");
360  }
361  }
362  flag = IDASVtolerances(m_ida_mem, m_reltol, m_abstol);
363  if (flag != IDA_SUCCESS) {
364  throw IDA_Err("Memory allocation failed.");
365  }
366  } else {
367  flag = IDAInit(m_ida_mem, ida_resid, m_t0, m_y, m_ydot);
368  if (flag != IDA_SUCCESS) {
369  if (flag == IDA_MEM_FAIL) {
370  throw IDA_Err("Memory allocation failed.");
371  } else if (flag == IDA_ILL_INPUT) {
372  throw IDA_Err("Illegal value for IDAMalloc input argument.");
373  } else {
374  throw IDA_Err("IDAMalloc failed.");
375  }
376  }
377  flag = IDASStolerances(m_ida_mem, m_reltol, m_abstols);
378  if (flag != IDA_SUCCESS) {
379  throw IDA_Err("Memory allocation failed.");
380  }
381  }
382 
383  //-----------------------------------
384  // set the linear solver type
385  //-----------------------------------
386 
387  if (m_type == 1 || m_type == 0) {
388  long int N = m_neq;
389  flag = IDADense(m_ida_mem, N);
390  if (flag) {
391  throw IDA_Err("IDADense failed");
392  }
393  } else if (m_type == 2) {
394  long int N = m_neq;
395  long int nu = m_mupper;
396  long int nl = m_mlower;
397  IDABand(m_ida_mem, N, nu, nl);
398  } else {
399  throw IDA_Err("unsupported linear solver type");
400  }
401 
402  if (m_formJac == 1) {
403  flag = IDADlsSetDenseJacFn(m_ida_mem, ida_jacobian);
404  if (flag != IDA_SUCCESS) {
405  throw IDA_Err("IDADlsSetDenseJacFn failed.");
406  }
407  }
408 
409  // pass a pointer to func in m_data
410  m_fdata = new ResidData(&m_resid, this, m_resid.nparams());
411  flag = IDASetUserData(m_ida_mem, (void*)m_fdata);
412  if (flag != IDA_SUCCESS) {
413  throw IDA_Err("IDASetUserData failed.");
414  }
415 
416  // set options
417  if (m_maxord > 0) {
418  flag = IDASetMaxOrd(m_ida_mem, m_maxord);
419  if (flag != IDA_SUCCESS) {
420  throw IDA_Err("IDASetMaxOrd failed.");
421  }
422  }
423  if (m_maxsteps > 0) {
424  flag = IDASetMaxNumSteps(m_ida_mem, m_maxsteps);
425  if (flag != IDA_SUCCESS) {
426  throw IDA_Err("IDASetMaxNumSteps failed.");
427  }
428  }
429  if (m_h0 > 0.0) {
430  flag = IDASetInitStep(m_ida_mem, m_h0);
431  if (flag != IDA_SUCCESS) {
432  throw IDA_Err("IDASetInitStep failed.");
433  }
434  }
435  if (m_tstop > 0.0) {
436  flag = IDASetStopTime(m_ida_mem, m_tstop);
437  if (flag != IDA_SUCCESS) {
438  throw IDA_Err("IDASetStopTime failed.");
439  }
440  }
441  if (m_maxErrTestFails >= 0) {
442  flag = IDASetMaxErrTestFails(m_ida_mem, m_maxErrTestFails);
443  if (flag != IDA_SUCCESS) {
444  throw IDA_Err("IDASetMaxErrTestFails failed.");
445  }
446  }
447  if (m_maxNonlinIters > 0) {
448  flag = IDASetMaxNonlinIters(m_ida_mem, m_maxNonlinIters);
449  if (flag != IDA_SUCCESS) {
450  throw IDA_Err("IDASetmaxNonlinIters failed.");
451  }
452  }
453  if (m_maxNonlinConvFails >= 0) {
454  flag = IDASetMaxConvFails(m_ida_mem, m_maxNonlinConvFails);
455  if (flag != IDA_SUCCESS) {
456  throw IDA_Err("IDASetMaxConvFails failed.");
457  }
458  }
459  if (m_setSuppressAlg != 0) {
460  flag = IDASetSuppressAlg(m_ida_mem, m_setSuppressAlg);
461  if (flag != IDA_SUCCESS) {
462  throw IDA_Err("IDASetSuppressAlg failed.");
463  }
464  }
465 }
466 
467 void IDA_Solver::correctInitial_Y_given_Yp(doublereal* y, doublereal* yp, doublereal tout)
468 {
469  int icopt = IDA_Y_INIT;
470  doublereal tout1 = tout;
471  if (tout == 0.0) {
472  double h0 = 1.0E-5;
473  if (m_h0 > 0.0) {
474  h0 = m_h0;
475  }
476  tout1 = m_t0 + h0;
477  }
478 
479  int flag = IDACalcIC(m_ida_mem, icopt, tout1);
480  if (flag != IDA_SUCCESS) {
481  throw IDA_Err("IDACalcIC failed: error = " + int2str(flag));
482  }
483 
484 
485  flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
486  if (flag != IDA_SUCCESS) {
487  throw IDA_Err("IDAGetSolution failed: error = " + int2str(flag));
488  }
489  doublereal* yy = NV_DATA_S(m_y);
490  doublereal* yyp = NV_DATA_S(m_ydot);
491 
492  for (int i = 0; i < m_neq; i++) {
493  y[i] = yy[i];
494  yp[i] = yyp[i];
495  }
496 }
497 
498 void IDA_Solver::correctInitial_YaYp_given_Yd(doublereal* y, doublereal* yp, doublereal tout)
499 {
500  int icopt = IDA_YA_YDP_INIT;
501  doublereal tout1 = tout;
502  if (tout == 0.0) {
503  double h0 = 1.0E-5;
504  if (m_h0 > 0.0) {
505  h0 = m_h0;
506  }
507  tout1 = m_t0 + h0;
508  }
509 
510  int flag = IDACalcIC(m_ida_mem, icopt, tout1);
511  if (flag != IDA_SUCCESS) {
512  throw IDA_Err("IDACalcIC failed: error = " + int2str(flag));
513  }
514 
515 
516  flag = IDAGetConsistentIC(m_ida_mem, m_y, m_ydot);
517  if (flag != IDA_SUCCESS) {
518  throw IDA_Err("IDAGetSolution failed: error = " + int2str(flag));
519  }
520  doublereal* yy = NV_DATA_S(m_y);
521  doublereal* yyp = NV_DATA_S(m_ydot);
522 
523  for (int i = 0; i < m_neq; i++) {
524  y[i] = yy[i];
525  yp[i] = yyp[i];
526  }
527 }
528 
529 int IDA_Solver::solve(double tout)
530 {
531  double tretn;
532  int flag;
533  flag = IDASetStopTime(m_ida_mem, tout);
534  if (flag != IDA_SUCCESS) {
535  throw IDA_Err(" IDA error encountered.");
536  }
537  do {
538  if (tout <= m_tcurrent) {
539  throw IDA_Err(" tout <= tcurrent");
540  }
541  m_told_old = m_told;
542  m_told = m_tcurrent;
543  flag = IDASolve(m_ida_mem, tout, &tretn, m_y, m_ydot, IDA_ONE_STEP);
544  if (flag < 0) {
545  throw IDA_Err(" IDA error encountered.");
546  } else if (flag == IDA_TSTOP_RETURN) {
547  // we've reached our goal, and have actually integrated past it
548  } else if (flag == IDA_ROOT_RETURN) {
549  // not sure what to do with this yet
550  } else if (flag == IDA_WARNING) {
551  throw IDA_Err(" IDA Warning encountered.");
552  }
553  m_tcurrent = tretn;
554  m_deltat = m_tcurrent - m_told;
555  } while (tretn < tout);
556 
557  if (flag != IDA_SUCCESS && flag != IDA_TSTOP_RETURN) {
558  throw IDA_Err(" IDA error encountered.");
559  }
560  return flag;
561 }
562 
563 double IDA_Solver::step(double tout)
564 {
565  double t;
566  int flag;
567  if (tout <= m_tcurrent) {
568  throw IDA_Err(" tout <= tcurrent");
569  }
570  m_told_old = m_told;
571  m_told = m_tcurrent;
572  flag = IDASolve(m_ida_mem, tout, &t, m_y, m_ydot, IDA_ONE_STEP);
573  if (flag < 0) {
574  throw IDA_Err(" IDA error encountered.");
575  } else if (flag == IDA_TSTOP_RETURN) {
576  // we've reached our goal, and have actually integrated past it
577  } else if (flag == IDA_ROOT_RETURN) {
578  // not sure what to do with this yet
579  } else if (flag == IDA_WARNING) {
580  throw IDA_Err(" IDA Warning encountered.");
581  }
582  m_tcurrent = t;
583  m_deltat = m_tcurrent - m_told;
584  return t;
585 }
586 
587 doublereal IDA_Solver::getOutputParameter(int flag) const
588 {
589  long int lenrw, leniw;
590  switch (flag) {
591  case REAL_WORKSPACE_SIZE:
592  flag = IDAGetWorkSpace(m_ida_mem, &lenrw, &leniw);
593  return doublereal(lenrw);
594  break;
595  }
596  return 0.0;
597 }
598 
599 }
600 #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:39
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, size_t nrhs, size_t ldb)
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