mrgsolve
lsoda_functions.h
Go to the documentation of this file.
1 #include <R.h>
2 typedef odeproblem* dtype;
3 
4 #ifndef LSODA_FUNCTIONS_H
5 #define LSODA_FUNCTIONS_H
6 #define ETA 2.2204460492503131e-16
7 
8 void LSODA::hmax_(const double value) {
9  rworks[2] = value;
10  if(value !=0) iopt = 1;
11 }
12 
13 void LSODA::hmin_(const double value) {
14  rworks[3] = value;
15  if(value !=0) iopt = 1;
16 }
17 
18 void LSODA::maxsteps_(const int value) {
19  Maxsteps = value;
20  iworks[3] = value;
21  if(value !=0) iopt=1;
22 }
23 
24 void LSODA::ixpr_(const int value) {
25  iworks[2] = value;
26  if(value!=0) iopt = 1;
27 }
28 
29 void LSODA::mxhnil_(const int value) {
30  iworks[4] = value;
31  if(value!=0) iopt = 1;
32 }
33 
34 bool LSODA::abs_compare(double a, double b)
35 {
36  return (std::abs(a) < std::abs(b));
37 }
38 
39 /* Purpose : Find largest component of double vector dx */
40 size_t LSODA::idamax1(const vector<double> &dx, const size_t n, const size_t offset = 0)
41 {
42 
43  size_t v = 0, vmax = 0;
44  size_t idmax = 1;
45  for (size_t i = 1; i <= n; ++i)
46  {
47  v = abs(dx[i + offset]);
48  if (v > vmax)
49  {
50  vmax = v;
51  idmax = i;
52  }
53  }
54  return idmax;
55 
56  // Following has failed with seg-fault. Probably issue with STL.
57  // return std::max_element( dx.begin()+1+offset, dx.begin()+1+n, LSODA::abs_compare) - dx.begin() - offset;
58 }
59 
60 /* Purpose : scalar vector multiplication
61  dx = da * dx
62 */
63 void LSODA::dscal1(const double da, vector<double> &dx, const size_t n,
64  const size_t offset = 0)
65 {
66  std::transform(dx.begin() + 1 + offset, dx.end(), dx.begin() + 1 + offset,
67  [&da](double x) -> double { return da * x; });
68 }
69 
70 /* Purpose : Inner product dx . dy */
71 double LSODA::ddot1(const vector<double> &a, const vector<double> &b,
72  const size_t n, const size_t offsetA = 0,
73  const size_t offsetB = 0)
74 {
75  double sum = 0.0;
76  for (size_t i = 1; i <= n; ++i)
77  sum += a[i + offsetA] * b[i + offsetB];
78  return sum;
79 }
80 
81 void LSODA::daxpy1(const double da, const vector<double> &dx,
82  vector<double> &dy, const size_t n, const size_t offsetX = 0,
83  const size_t offsetY = 0)
84 {
85 
86  for (size_t i = 1; i <= n; ++i)
87  dy[i + offsetY] = da * dx[i + offsetX] + dy[i + offsetY];
88 }
89 
90 // See BLAS documentation. The first argument has been changed to vector.
91 void LSODA::dgesl1(const vector<vector<double>> &a, const size_t n,
92  vector<int> &ipvt, vector<double> &b, const size_t job)
93 {
94  size_t k, j;
95  double t;
96 
97  /*
98  Job = 0, solve a * x = b.
99  */
100  if (job == 0)
101  {
102  /*
103  First solve L * y = b.
104  */
105  for (k = 1; k <= n; ++k)
106  {
107  t = ddot1(a[k], b, k - 1);
108  b[k] = (b[k] - t) / a[k][k];
109  }
110  /*
111  Now solve U * x = y.
112  */
113  for (k = n - 1; k >= 1; k--)
114  {
115  b[k] = b[k] + ddot1(a[k], b, n - k, k, k);
116  j = ipvt[k];
117  if (j != k)
118  {
119  t = b[j];
120  b[j] = b[k];
121  b[k] = t;
122  }
123  }
124  return;
125  }
126  /*
127  Job = nonzero, solve Transpose(a) * x = b.
128 
129  First solve Transpose(U) * y = b.
130  */
131  for (k = 1; k <= n - 1; ++k)
132  {
133  j = ipvt[k];
134  t = b[j];
135  if (j != k)
136  {
137  b[j] = b[k];
138  b[k] = t;
139  }
140  daxpy1(t, a[k], b, n - k, k, k);
141  }
142  /*
143  Now solve Transpose(L) * x = y.
144  */
145  for (k = n; k >= 1; k--)
146  {
147  b[k] = b[k] / a[k][k];
148  t = -b[k];
149  daxpy1(t, a[k], b, k - 1);
150  }
151 }
152 
153 // See BLAS documentation. All double* has been changed to std::vector .
154 void LSODA::dgefa1(vector<vector<double>> &a, const size_t n, vector<int> &ipvt,
155  size_t *const info)
156 {
157  size_t j = 0, k = 0, i = 0;
158  double t = 0.0;
159 
160  /* Gaussian elimination with partial pivoting. */
161 
162  *info = 0;
163  for (k = 1; k <= n - 1; ++k)
164  {
165  /*
166  Find j = pivot index. Note that a[k]+k-1 is the address of
167  the 0-th element of the row vector whose 1st element is a[k][k].
168  */
169  j = idamax1(a[k], n - k + 1, k - 1) + k - 1;
170  ipvt[k] = j;
171  /*
172  Zero pivot implies this row already triangularized.
173  */
174  if (a[k][j] == 0.)
175  {
176  *info = k;
177  continue;
178  }
179  /*
180  Interchange if necessary.
181  */
182  if (j != k)
183  {
184  t = a[k][j];
185  a[k][j] = a[k][k];
186  a[k][k] = t;
187  }
188  /*
189  Compute multipliers.
190  */
191  t = -1. / a[k][k];
192  dscal1(t, a[k], n - k, k);
193 
194  /*
195  Column elimination with row indexing.
196  */
197  for (i = k + 1; i <= n; ++i)
198  {
199  t = a[i][j];
200  if (j != k)
201  {
202  a[i][j] = a[i][k];
203  a[i][k] = t;
204  }
205  daxpy1(t, a[k], a[i], n - k, k, k);
206  }
207  } /* end k-loop */
208 
209  ipvt[n] = n;
210  if (a[n][n] == 0.)
211  *info = n;
212 }
213 
214 /* Terminate lsoda due to illegal input. */
215 void LSODA::terminate(int *istate)
216 {
217  if (illin == 5) {
218  // cerr << "[lsoda] repeated occurrence of illegal input. run aborted.. "
219  // "apparent infinite loop."
220  // << endl;
221  REprintf("[lsoda] repeated occurrence of illegal input. run aborted;\n");
222  REprintf(" apparent infinite loop.\n");
223  } else
224  {
225  illin++;
226  *istate = -3;
227  }
228  return;
229 }
230 
231 /* Terminate lsoda due to various error conditions. */
232 void LSODA::terminate2(vector<double> &y, double *t)
233 {
234  for (size_t i = 1; i <= n; ++i)
235  y[i] = yh_[1][i];
236  *t = tn_;
237  illin = 0;
238  return;
239 }
240 
241 /*
242  The following block handles all successful returns from lsoda.
243  If itask != 1, y is loaded from yh_ and t is set accordingly.
244  *Istate is set to 2, the illegal input counter is zeroed, and the
245  optional outputs are loaded into the work arrays before returning.
246 */
247 
248 void LSODA::successreturn(vector<double> &y, double *t, int itask, int ihit,
249  double tcrit, int *istate)
250 {
251  for (size_t i = 1; i <= n; ++i)
252  y[i] = yh_[1][i];
253  *t = tn_;
254  if (itask == 4 || itask == 5)
255  if (ihit)
256  *t = tcrit;
257  *istate = 2;
258  illin = 0;
259 }
260 
261 /*
262 c references..
263 c 1. alan c. hindmarsh, odepack, a systematized collection of ode
264 c solvers, in scientific computing, r. s. stepleman et al. (eds.),
265 c north-holland, amsterdam, 1983, pp. 55-64.
266 c 2. linda r. petzold, automatic selection of methods for solving
267 c stiff and nonstiff systems of ordinary differential equations,
268 c siam j. sci. stat. comput. 4 (1983), pp. 136-148.
269 c-----------------------------------------------------------------------
270 */
271 
272 void LSODA::stoda(const size_t neq, vector<double> &y, LSODA_ODE_SYSTEM_TYPE f,
273  dtype _data)
274 {
275  assert(neq + 1 == y.size());
276 
277  size_t corflag = 0, orderflag = 0;
278  size_t i = 0, i1 = 0, j = 0, m = 0, ncf = 0;
279  double del = 0.0, delp = 0.0, dsm = 0.0, dup = 0.0, exup = 0.0, r = 0.0,
280  rh = 0.0, rhup = 0.0, told = 0.0;
281  double pdh = 0.0, pnorm = 0.0;
282 
283  /*
284  stoda performs one step of the integration of an initial value
285  problem for a system of ordinary differential equations.
286  Note.. stoda is independent of the value of the iteration method
287  indicator miter, when this is != 0, and hence is independent
288  of the type of chord method used, or the Jacobian structure.
289  Communication with stoda is done with the following variables:
290 
291  jstart = an integer used for input only, with the following
292  values and meanings:
293 
294  0 perform the first step,
295  > 0 take a new step continuing from the last,
296  -1 take the next step with a new value of h_,
297  n, meth_, miter, and/or matrix parameters.
298  -2 take the next step with a new value of h_,
299  but with other inputs unchanged.
300 
301  kflag = a completion code with the following meanings:
302 
303  0 the step was successful,
304  -1 the requested error could not be achieved,
305  -2 corrector convergence could not be achieved,
306  -3 fatal error in prja or solsy.
307 
308  miter = corrector iteration method:
309 
310  0 functional iteration,
311  >0 a chord method corresponding to jacobian type jt.
312 
313  */
314  kflag = 0;
315  told = tn_;
316  ncf = 0;
317  ierpj = 0;
318  iersl = 0;
319  jcur = 0;
320  delp = 0.;
321 
322  /*
323  On the first call, the order is set to 1, and other variables are
324  initialized. rmax is the maximum ratio by which h_ can be increased
325  in a single step. It is initially 1.e4 to compensate for the small
326  initial h_, but then is normally equal to 10. If a filure occurs
327  (in corrector convergence or error test), rmax is set at 2 for
328  the next increase.
329  cfode is called to get the needed coefficients for both methods.
330  */
331  if (jstart == 0)
332  {
333  lmax = maxord + 1;
334  nq = 1;
335  l = 2;
336  ialth = 2;
337  rmax = 10000.;
338  rc = 0.;
339  el0 = 1.;
340  crate = 0.7;
341  hold = h_;
342  nslp = 0;
343  ipup = miter;
344  /*
345  Initialize switching parameters. meth_ = 1 is assumed initially.
346  */
347  icount = 20;
348  irflag = 0;
349  pdest = 0.;
350  pdlast = 0.;
351  ratio = 5.;
352  cfode(2);
353  for (i = 1; i <= 5; ++i)
354  cm2[i] = tesco[i][2] * elco[i][i + 1];
355  cfode(1);
356  for (i = 1; i <= 12; ++i)
357  cm1[i] = tesco[i][2] * elco[i][i + 1];
358  resetcoeff();
359  } /* end if ( jstart == 0 ) */
360  /*
361  The following block handles preliminaries needed when jstart = -1.
362  ipup is set to miter to force a matrix update.
363  If an order increase is about to be considered ( ialth = 1 ),
364  ialth is reset to 2 to postpone consideration one more step.
365  If the caller has changed meth_, cfode is called to reset
366  the coefficients of the method.
367  If h_ is to be changed, yh_ must be rescaled.
368  If h_ or meth_ is being changed, ialth is reset to l = nq + 1
369  to prevent further changes in h_ for that many steps.
370  */
371  if (jstart == -1)
372  {
373  ipup = miter;
374  lmax = maxord + 1;
375  if (ialth == 1)
376  ialth = 2;
377  if (meth_ != mused)
378  {
379  cfode(meth_);
380  ialth = l;
381  resetcoeff();
382  }
383  if (h_ != hold)
384  {
385  rh = h_ / hold;
386  h_ = hold;
387  scaleh(&rh, &pdh);
388  }
389  } /* if ( jstart == -1 ) */
390  if (jstart == -2)
391  {
392  if (h_ != hold)
393  {
394  rh = h_ / hold;
395  h_ = hold;
396  scaleh(&rh, &pdh);
397  }
398  } /* if ( jstart == -2 ) */
399 
400  /*
401  Prediction.
402  This section computes the predicted values by effectively
403  multiplying the yh_ array by the pascal triangle matrix.
404  rc is the ratio of new to old values of the coefficient h_ * el[1].
405  When rc differs from 1 by more than ccmax, ipup is set to miter
406  to force pjac to be called, if a jacobian is involved.
407  In any case, prja is called at least every msbp steps.
408  */
409  while (1)
410  {
411  while (1)
412  {
413  if (fabs(rc - 1.) > ccmax)
414  ipup = miter;
415  if (nst >= nslp + msbp)
416  ipup = miter;
417  tn_ += h_;
418  for (size_t j = nq; j >= 1; j--)
419  for (size_t i1 = j; i1 <= nq; ++i1)
420  for (i = 1; i <= n; ++i)
421  yh_[i1][i] += yh_[i1 + 1][i];
422 
423  pnorm = vmnorm(n, yh_[1], ewt);
424  correction(neq, y, f, &corflag, pnorm, &del, &delp, &told, &ncf, &rh, &m,
425  _data);
426  if (corflag == 0)
427  break;
428  if (corflag == 1)
429  {
430  rh = max(rh, hmin / fabs(h_));
431  scaleh(&rh, &pdh);
432  continue;
433  }
434  if (corflag == 2)
435  {
436  kflag = -2;
437  hold = h_;
438  jstart = 1;
439  return;
440  }
441  } /* end inner while ( corrector loop ) */
442 
443  /*
444  The corrector has converged. jcur is set to 0
445  to signal that the Jacobian involved may need updating later.
446  The local error test is done now.
447  */
448  jcur = 0;
449  if (m == 0)
450  dsm = del / tesco[nq][2];
451  if (m > 0)
452  dsm = vmnorm(n, acor, ewt) / tesco[nq][2];
453 
454  if (dsm <= 1.)
455  {
456  /*
457  After a successful step, update the yh_ array.
458  Decrease icount by 1, and if it is -1, consider switching methods.
459  If a method switch is made, reset various parameters,
460  rescale the yh_ array, and exit. If there is no switch,
461  consider changing h_ if ialth = 1. Otherwise decrease ialth by 1.
462  If ialth is then 1 and nq < maxord, then acor is saved for
463  use in a possible order increase on the next step.
464  If a change in h_ is considered, an increase or decrease in order
465  by one is considered also. A change in h_ is made only if it is by
466  a factor of at least 1.1. If not, ialth is set to 3 to prevent
467  testing for that many steps.
468  */
469  kflag = 0;
470  nst++;
471  hu = h_;
472  nqu = nq;
473  mused = meth_;
474  for (size_t j = 1; j <= l; ++j)
475  {
476  r = el[j];
477  for (i = 1; i <= n; ++i)
478  yh_[j][i] += r * acor[i];
479  }
480  icount--;
481  if (icount < 0)
482  {
483  methodswitch(dsm, pnorm, &pdh, &rh);
484  if (meth_ != mused)
485  {
486  rh = max(rh, hmin / fabs(h_));
487  scaleh(&rh, &pdh);
488  rmax = 10.;
489  endstoda();
490  break;
491  }
492  }
493  /*
494  No method switch is being made. Do the usual step/order selection.
495  */
496  ialth--;
497  if (ialth == 0)
498  {
499  rhup = 0.;
500  if (l != lmax)
501  {
502  for (i = 1; i <= n; ++i)
503  savf[i] = acor[i] - yh_[lmax][i];
504  dup = vmnorm(n, savf, ewt) / tesco[nq][3];
505  exup = 1. / (double)(l + 1);
506  rhup = 1. / (1.4 * pow(dup, exup) + 0.0000014);
507  }
508 
509  orderswitch(&rhup, dsm, &pdh, &rh, &orderflag);
510 
511  /*
512  No change in h_ or nq.
513  */
514  if (orderflag == 0)
515  {
516  endstoda();
517  break;
518  }
519  /*
520  h_ is changed, but not nq.
521  */
522  if (orderflag == 1)
523  {
524  rh = max(rh, hmin / fabs(h_));
525  scaleh(&rh, &pdh);
526  rmax = 10.;
527  endstoda();
528  break;
529  }
530  /*
531  both nq and h_ are changed.
532  */
533  if (orderflag == 2)
534  {
535  resetcoeff();
536  rh = max(rh, hmin / fabs(h_));
537  scaleh(&rh, &pdh);
538  rmax = 10.;
539  endstoda();
540  break;
541  }
542  } /* end if ( ialth == 0 ) */
543  if (ialth > 1 || l == lmax)
544  {
545  endstoda();
546  break;
547  }
548 
549  for (size_t i = 1; i <= n; ++i)
550  yh_[lmax][i] = acor[i];
551 
552  endstoda();
553  break;
554  }
555  /* end if ( dsm <= 1. ) */
556  /*
557  The error test failed. kflag keeps track of multiple failures.
558  Restore tn_ and the yh_ array to their previous values, and prepare
559  to try the step again. Compute the optimum step size for this or
560  one lower. After 2 or more failures, h_ is forced to decrease
561  by a factor of 0.2 or less.
562  */
563  else
564  {
565  kflag--;
566  tn_ = told;
567  for (j = nq; j >= 1; j--)
568  {
569  for (i1 = j; i1 <= nq; ++i1)
570  for (i = 1; i <= n; ++i)
571  yh_[i1][i] -= yh_[i1 + 1][i];
572  }
573  rmax = 2.;
574  if (fabs(h_) <= hmin * 1.00001)
575  {
576  kflag = -1;
577  hold = h_;
578  jstart = 1;
579  break;
580  }
581  if (kflag > -3)
582  {
583  rhup = 0.;
584  orderswitch(&rhup, dsm, &pdh, &rh, &orderflag);
585  if (orderflag == 1 || orderflag == 0)
586  {
587  if (orderflag == 0)
588  rh = min(rh, 0.2);
589  rh = max(rh, hmin / fabs(h_));
590  scaleh(&rh, &pdh);
591  }
592  if (orderflag == 2)
593  {
594  resetcoeff();
595  rh = max(rh, hmin / fabs(h_));
596  scaleh(&rh, &pdh);
597  }
598  continue;
599  }
600  /* if ( kflag > -3 ) */
601  /*
602  Control reaches this section if 3 or more failures have occurred.
603  If 10 failures have occurred, exit with kflag = -1.
604  It is assumed that the derivatives that have accumulated in the
605  yh_ array have errors of the wrong order. Hence the first
606  derivative is recomputed, and the order is set to 1. Then
607  h_ is reduced by a factor of 10, and the step is retried,
608  until it succeeds or h_ reaches hmin.
609  */
610  else
611  {
612  if (kflag == -10)
613  {
614  kflag = -1;
615  hold = h_;
616  jstart = 1;
617  break;
618  }
619  else
620  {
621  rh = 0.1;
622  rh = max(hmin / fabs(h_), rh);
623  h_ *= rh;
624  for (i = 1; i <= n; ++i)
625  y[i] = yh_[1][i];
626  (*f)(tn_, &y[1], &savf[1], _data);
627  nfe++;
628  for (i = 1; i <= n; ++i)
629  yh_[2][i] = h_ * savf[i];
630  ipup = miter;
631  ialth = 5;
632  if (nq == 1)
633  continue;
634  nq = 1;
635  l = 2;
636  resetcoeff();
637  continue;
638  }
639  } /* end else -- kflag <= -3 */
640  } /* end error failure handling */
641  } /* end outer while */
642 
643 } /* end stoda */
644 
645 void LSODA::ewset(const vector<double> &ycur)
646 {
647  switch (itol_)
648  {
649  case 1:
650  for (size_t i = 1; i <= n; ++i)
651  ewt[i] = rtol_[1] * fabs(ycur[i]) + atol_[1];
652  break;
653  case 2:
654  for (size_t i = 1; i <= n; ++i)
655  ewt[i] = rtol_[1] * fabs(ycur[i]) + atol_[i];
656  break;
657  case 3:
658  for (size_t i = 1; i <= n; ++i)
659  ewt[i] = rtol_[i] * fabs(ycur[i]) + atol_[1];
660  break;
661  case 4:
662  for (size_t i = 1; i <= n; ++i)
663  ewt[i] = rtol_[i] * fabs(ycur[i]) + atol_[i];
664  break;
665  }
666 
667 } /* end ewset */
668 
669 /*
670  Intdy computes interpolated values of the k-th derivative of the
671  dependent variable vector y, and stores it in dky. This routine
672  is called within the package with k = 0 and *t = tout, but may
673  also be called by the user for any k up to the current order.
674  ( See detailed instructions in the usage documentation. )
675 
676  The computed values in dky are gotten by interpolation using the
677  Nordsieck history array yh_. This array corresponds uniquely to a
678  vector-valued polynomial of degree nqcur or less, and dky is set
679  to the k-th derivative of this polynomial at t.
680  The formula for dky is
681 
682  q
683  dky[i] = sum c[k][j] * ( t - tn_ )^(j-k) * h_^(-j) * yh_[j+1][i]
684  j=k
685 
686  where c[k][j] = j*(j-1)*...*(j-k+1), q = nqcur, tn_ = tcur, h_ = hcur.
687  The quantities nq = nqcur, l = nq+1, n = neq, tn_, and h_ are declared
688  static globally. The above sum is done in reverse order.
689  *iflag is returned negative if either k or t is out of bounds.
690 */
691 void LSODA::intdy(double t, int k, vector<double> &dky, int *iflag)
692 {
693  int ic, jp1 = 0;
694  double c, r, s, tp;
695 
696  *iflag = 0;
697  if (k < 0 || k > (int)nq)
698  {
699  //fprintf(stderr, "[intdy] k = %d illegal\n", k);
700  REprintf("[intdy] k = %d illegal.\n", k);
701  *iflag = -1;
702  return;
703  }
704  tp = tn_ - hu - 100. * ETA * (tn_ + hu);
705  if ((t - tp) * (t - tn_) > 0.)
706  {
707  // fprintf(stderr,
708  // "intdy -- t = %g illegal. t not in interval tcur - hu to tcur\n",
709  // t);
710  REprintf("[intdy] t = %g illegal. t not in interval tcur - hu to tcur.\n\n",t);
711  *iflag = -2;
712  return;
713  }
714  s = (t - tn_) / h_;
715  ic = 1;
716  for (size_t jj = l - k; jj <= nq; ++jj)
717  ic *= jj;
718  c = (double)ic;
719  for (size_t i = 1; i <= n; ++i)
720  dky[i] = c * yh_[l][i];
721 
722  for (int j = nq - 1; j >= k; j--)
723  {
724  jp1 = j + 1;
725  ic = 1;
726  for (int jj = jp1 - k; jj <= j; ++jj)
727  ic *= jj;
728  c = (double)ic;
729 
730  for (size_t i = 1; i <= n; ++i)
731  dky[i] = c * yh_[jp1][i] + s * dky[i];
732  }
733  if (k == 0)
734  return;
735  r = pow(h_, (double)(-k));
736 
737  for (size_t i = 1; i <= n; ++i)
738  dky[i] *= r;
739 
740 } /* end intdy */
741 
742 void LSODA::cfode(int meth_)
743 {
744  int i, nq, nqm1, nqp1;
745  double agamq, fnq, fnqm1, pc[13], pint, ragq, rqfac, rq1fac, tsign, xpin;
746  /*
747  cfode is called by the integrator routine to set coefficients
748  needed there. The coefficients for the current method, as
749  given by the value of meth_, are set for all orders and saved.
750  The maximum order assumed here is 12 if meth_ = 1 and 5 if meth_ = 2.
751  ( A smaller value of the maximum order is also allowed. )
752  cfode is called once at the beginning of the problem, and
753  is not called again unless and until meth_ is changed.
754 
755  The elco array contains the basic method coefficients.
756  The coefficients el[i], 1 < i < nq+1, for the method of
757  order nq are stored in elco[nq][i]. They are given by a generating
758  polynomial, i.e.,
759 
760  l(x) = el[1] + el[2]*x + ... + el[nq+1]*x^nq.
761 
762  For the implicit Adams method, l(x) is given by
763 
764  dl/dx = (x+1)*(x+2)*...*(x+nq-1)/factorial(nq-1), l(-1) = 0.
765 
766  For the bdf methods, l(x) is given by
767 
768  l(x) = (x+1)*(x+2)*...*(x+nq)/k,
769 
770  where k = factorial(nq)*(1+1/2+...+1/nq).
771 
772  The tesco array contains test constants used for the
773  local error test and the selection of step size and/or order.
774  At order nq, tesco[nq][k] is used for the selection of step
775  size at order nq-1 if k = 1, at order nq if k = 2, and at order
776  nq+1 if k = 3.
777  */
778  if (meth_ == 1)
779  {
780  elco[1][1] = 1.;
781  elco[1][2] = 1.;
782  tesco[1][1] = 0.;
783  tesco[1][2] = 2.;
784  tesco[2][1] = 1.;
785  tesco[12][3] = 0.;
786  pc[1] = 1.;
787  rqfac = 1.;
788  for (nq = 2; nq <= 12; nq++)
789  {
790  /*
791  The pc array will contain the coefficients of the polynomial
792 
793  p(x) = (x+1)*(x+2)*...*(x+nq-1).
794 
795  Initially, p(x) = 1.
796  */
797  rq1fac = rqfac;
798  rqfac = rqfac / (double)nq;
799  nqm1 = nq - 1;
800  fnqm1 = (double)nqm1;
801  nqp1 = nq + 1;
802  /*
803  Form coefficients of p(x)*(x+nq-1).
804  */
805  pc[nq] = 0.;
806  for (i = nq; i >= 2; i--)
807  pc[i] = pc[i - 1] + fnqm1 * pc[i];
808  pc[1] = fnqm1 * pc[1];
809  /*
810  Compute integral, -1 to 0, of p(x) and x*p(x).
811  */
812  pint = pc[1];
813  xpin = pc[1] / 2.;
814  tsign = 1.;
815  for (i = 2; i <= nq; ++i)
816  {
817  tsign = -tsign;
818  pint += tsign * pc[i] / (double)i;
819  xpin += tsign * pc[i] / (double)(i + 1);
820  }
821  /*
822  Store coefficients in elco and tesco.
823  */
824  elco[nq][1] = pint * rq1fac;
825  elco[nq][2] = 1.;
826  for (i = 2; i <= nq; ++i)
827  elco[nq][i + 1] = rq1fac * pc[i] / (double)i;
828  agamq = rqfac * xpin;
829  ragq = 1. / agamq;
830  tesco[nq][2] = ragq;
831  if (nq < 12)
832  tesco[nqp1][1] = ragq * rqfac / (double)nqp1;
833  tesco[nqm1][3] = ragq;
834  } /* end for */
835  return;
836  } /* end if ( meth_ == 1 ) */
837 
838  /* meth_ = 2. */
839  pc[1] = 1.;
840  rq1fac = 1.;
841 
842  /*
843  The pc array will contain the coefficients of the polynomial
844  p(x) = (x+1)*(x+2)*...*(x+nq).
845  Initially, p(x) = 1.
846  */
847  for (nq = 1; nq <= 5; ++nq)
848  {
849  fnq = (double)nq;
850  nqp1 = nq + 1;
851  /*
852  Form coefficients of p(x)*(x+nq).
853  */
854  pc[nqp1] = 0.;
855  for (i = nq + 1; i >= 2; i--)
856  pc[i] = pc[i - 1] + fnq * pc[i];
857  pc[1] *= fnq;
858  /*
859  Store coefficients in elco and tesco.
860  */
861  for (i = 1; i <= nqp1; ++i)
862  elco[nq][i] = pc[i] / pc[2];
863  elco[nq][2] = 1.;
864  tesco[nq][1] = rq1fac;
865  tesco[nq][2] = ((double)nqp1) / elco[nq][1];
866  tesco[nq][3] = ((double)(nq + 2)) / elco[nq][1];
867  rq1fac /= fnq;
868  }
869  return;
870 
871 } /* end cfode */
872 
873 void LSODA::scaleh(double *rh, double *pdh)
874 {
875  double r;
876  /*
877  If h_ is being changed, the h_ ratio rh is checked against rmax, hmin,
878  and hmxi, and the yh_ array is rescaled. ialth is set to l = nq + 1
879  to prevent a change of h_ for that many steps, unless forced by a
880  convergence or error test failure.
881  */
882  *rh = min(*rh, rmax);
883  *rh = *rh / max(1., fabs(h_) * hmxi * *rh);
884  /*
885  If meth_ = 1, also restrict the new step size by the stability region.
886  If this reduces h_, set irflag to 1 so that if there are roundoff
887  problems later, we can assume that is the cause of the trouble.
888  */
889  if (meth_ == 1)
890  {
891  irflag = 0;
892  *pdh = max(fabs(h_) * pdlast, 0.000001);
893  if ((*rh * *pdh * 1.00001) >= sm1[nq])
894  {
895  *rh = sm1[nq] / *pdh;
896  irflag = 1;
897  }
898  }
899  r = 1.;
900  for (size_t j = 2; j <= l; ++j)
901  {
902  r *= *rh;
903  for (size_t i = 1; i <= n; ++i)
904  yh_[j][i] *= r;
905  }
906  h_ *= *rh;
907  rc *= *rh;
908  ialth = l;
909 
910 } /* end scaleh */
911 
912 void LSODA::prja(const size_t neq, vector<double> &y, LSODA_ODE_SYSTEM_TYPE f,
913  dtype _data)
914 {
915  size_t i = 0, ier = 0, j = 0;
916  double fac = 0.0, hl0 = 0.0, r = 0.0, r0 = 0.0, yj = 0.0;
917  /*
918  prja is called by stoda to compute and process the matrix
919  P = I - h_ * el[1] * J, where J is an approximation to the Jacobian.
920  Here J is computed by finite differencing.
921  J, scaled by -h_ * el[1], is stored in wm_. Then the norm of J ( the
922  matrix norm consistent with the weighted max-norm on vectors given
923  by vmnorm ) is computed, and J is overwritten by P. P is then
924  subjected to LU decomposition in preparation for later solution
925  of linear systems with p as coefficient matrix. This is done
926  by dgefa1 if miter = 2, and by dgbfa if miter = 5.
927  */
928  nje++;
929  ierpj = 0;
930  jcur = 1;
931  hl0 = h_ * el0;
932  /*
933  If miter = 2, make n calls to f to approximate J.
934  */
935  if (miter != 2)
936  {
937  //fprintf(stderr, "[prja] miter != 2\n");
938  REprintf("[prja] miter != 2.\n");
939  return;
940  }
941  if (miter == 2)
942  {
943  fac = vmnorm(n, savf, ewt);
944  r0 = 1000. * fabs(h_) * ETA * ((double)n) * fac;
945  if (r0 == 0.)
946  r0 = 1.;
947  for (j = 1; j <= n; ++j)
948  {
949  yj = y[j];
950  r = max(sqrteta * fabs(yj), r0 / ewt[j]);
951  y[j] += r;
952  fac = -hl0 / r;
953  (*f)(tn_, &y[1], &acor[1], _data);
954  for (i = 1; i <= n; ++i)
955  wm_[i][j] = (acor[i] - savf[i]) * fac;
956  y[j] = yj;
957  }
958  nfe += n;
959  /*
960  Compute norm of Jacobian.
961  */
962  pdnorm = fnorm(n, wm_, ewt) / fabs(hl0);
963  /*
964  Add identity matrix.
965  */
966  for (i = 1; i <= n; ++i)
967  wm_[i][i] += 1.;
968  /*
969  Do LU decomposition on P.
970  */
971  dgefa1(wm_, n, ipvt, &ier);
972  if (ier != 0)
973  ierpj = 1;
974  return;
975  }
976 } /* end prja */
977 
978 /*
979  This function routine computes the weighted max-norm
980  of the vector of length n contained in the array v, with weights
981  contained in the array w of length n.
982 
983  vmnorm = max( i = 1, ..., n ) fabs( v[i] ) * w[i].
984 */
985 double LSODA::vmnorm(const size_t n, const vector<double> &v,
986  const vector<double> &w)
987 {
988  double vm = 0.;
989  for (size_t i = 1; i <= n; ++i)
990  vm = max(vm, fabs(v[i]) * w[i]);
991  return vm;
992 }
993 
994 double LSODA::fnorm(int n, const vector<vector<double>> &a,
995  const vector<double> &w)
996 
997 /*
998  This subroutine computes the norm of a full n by n matrix,
999  stored in the array a, that is consistent with the weighted max-norm
1000  on vectors, with weights stored in the array w.
1001 
1002  fnorm = max(i=1,...,n) ( w[i] * sum(j=1,...,n) fabs( a[i][j] ) / w[j] )
1003 */
1004 
1005 {
1006  double an = 0, sum = 0;
1007 
1008  for (size_t i = 1; i <= (size_t)n; ++i)
1009  {
1010  sum = 0.;
1011  for (size_t j = 1; j <= (size_t)n; ++j)
1012  sum += fabs(a[i][j]) / w[j];
1013  an = max(an, sum * w[i]);
1014  }
1015  return an;
1016 }
1017 
1018 /*
1019  *corflag = 0 : corrector converged,
1020 1 : step size to be reduced, redo prediction,
1021 2 : corrector cannot converge, failure flag.
1022 */
1023 void LSODA::correction(const size_t neq, vector<double> &y,
1024  LSODA_ODE_SYSTEM_TYPE f, size_t *corflag, double pnorm,
1025  double *del, double *delp, double *told,size_t *ncf,
1026  double *rh, size_t *m, dtype _data)
1027 {
1028  double rm = 0.0, rate = 0.0, dcon = 0.0;
1029 
1030  /*
1031  Up to maxcor corrector iterations are taken. A convergence test is
1032  made on the r.m.s. norm of each correction, weighted by the error
1033  weight vector ewt. The sum of the corrections is accumulated in the
1034  vector acor[i]. The yh_ array is not altered in the corrector loop.
1035  */
1036 
1037  *m = 0;
1038  *corflag = 0;
1039  *del = 0.;
1040 
1041  for (size_t i = 1; i <= n; ++i)
1042  y[i] = yh_[1][i];
1043 
1044  (*f)(tn_, &y[1], &savf[1], _data);
1045 
1046  nfe++;
1047  /*
1048  If indicated, the matrix P = I - h_ * el[1] * J is reevaluated and
1049  preprocessed before starting the corrector iteration. ipup is set
1050  to 0 as an indicator that this has been done.
1051  */
1052  while (1)
1053  {
1054  if (*m == 0)
1055  {
1056  if (ipup > 0)
1057  {
1058  prja(neq, y, f, _data);
1059  ipup = 0;
1060  rc = 1.;
1061  nslp = nst;
1062  crate = 0.7;
1063  if (ierpj != 0)
1064  {
1065  corfailure(told, rh, ncf, corflag);
1066  return;
1067  }
1068  }
1069  for (size_t i = 1; i <= n; ++i)
1070  acor[i] = 0.;
1071  } /* end if ( *m == 0 ) */
1072  if (miter == 0)
1073  {
1074  /*
1075  In case of functional iteration, update y directly from
1076  the result of the last function evaluation.
1077  */
1078  for (size_t i = 1; i <= n; ++i)
1079  {
1080  savf[i] = h_ * savf[i] - yh_[2][i];
1081  y[i] = savf[i] - acor[i];
1082  }
1083  *del = vmnorm(n, y, ewt);
1084  for (size_t i = 1; i <= n; ++i)
1085  {
1086  y[i] = yh_[1][i] + el[1] * savf[i];
1087  acor[i] = savf[i];
1088  }
1089  }
1090  /* end functional iteration */
1091  /*
1092  In the case of the chord method, compute the corrector error,
1093  and solve the linear system with that as right-hand side and
1094  P as coefficient matrix.
1095  */
1096  else
1097  {
1098  for (size_t i = 1; i <= n; ++i)
1099  y[i] = h_ * savf[i] - (yh_[2][i] + acor[i]);
1100 
1101  solsy(y);
1102  *del = vmnorm(n, y, ewt);
1103 
1104  for (size_t i = 1; i <= n; ++i)
1105  {
1106  acor[i] += y[i];
1107  y[i] = yh_[1][i] + el[1] * acor[i];
1108  }
1109  } /* end chord method */
1110  /*
1111  Test for convergence. If *m > 0, an estimate of the convergence
1112  rate constant is stored in crate, and this is used in the test.
1113 
1114  We first check for a change of iterates that is the size of
1115  roundoff error. If this occurs, the iteration has converged, and a
1116  new rate estimate is not formed.
1117  In all other cases, force at least two iterations to estimate a
1118  local Lipschitz constant estimate for Adams method.
1119  On convergence, form pdest = local maximum Lipschitz constant
1120  estimate. pdlast is the most recent nonzero estimate.
1121  */
1122  if (*del <= 100. * pnorm * ETA)
1123  break;
1124  if (*m != 0 || meth_ != 1)
1125  {
1126  if (*m != 0)
1127  {
1128  rm = 1024.0;
1129  if (*del <= (1024. * *delp))
1130  rm = *del / *delp;
1131  rate = max(rate, rm);
1132  crate = max(0.2 * crate, rm);
1133  }
1134  dcon = *del * min(1., 1.5 * crate) / (tesco[nq][2] * conit);
1135  if (dcon <= 1.)
1136  {
1137  pdest = max(pdest, rate / fabs(h_ * el[1]));
1138  if (pdest != 0.)
1139  pdlast = pdest;
1140  break;
1141  }
1142  }
1143  /*
1144  The corrector iteration failed to converge.
1145  If miter != 0 and the Jacobian is out of date, prja is called for
1146  the next try. Otherwise the yh_ array is retracted to its values
1147  before prediction, and h_ is reduced, if possible. If h_ cannot be
1148  reduced or mxncf failures have occured, exit with corflag = 2.
1149  */
1150  (*m)++;
1151  if (*m == maxcor || (*m >= 2 && *del > 2. * *delp))
1152  {
1153  if (miter == 0 || jcur == 1)
1154  {
1155  corfailure(told, rh, ncf, corflag);
1156  return;
1157  }
1158  ipup = miter;
1159  /*
1160  Restart corrector if Jacobian is recomputed.
1161  */
1162  *m = 0;
1163  rate = 0.;
1164  *del = 0.;
1165  for (size_t i = 1; i <= n; ++i)
1166  y[i] = yh_[1][i];
1167 
1168  (*f)(tn_, &y[1], &savf[1], _data);
1169 
1170  nfe++;
1171  }
1172  /*
1173  Iterate corrector.
1174  */
1175  else
1176  {
1177  *delp = *del;
1178  (*f)(tn_, &y[1], &savf[1], _data);
1179  nfe++;
1180  }
1181  } /* end while */
1182 } /* end correction */
1183 
1184 void LSODA::corfailure(double *told, double *rh, size_t *ncf, size_t *corflag)
1185 {
1186  (*ncf)++;
1187  rmax = 2.;
1188  tn_ = *told;
1189  for (size_t j = nq; j >= 1; j--)
1190  for (size_t i1 = j; i1 <= nq; ++i1)
1191  for (size_t i = 1; i <= n; ++i)
1192  yh_[i1][i] -= yh_[i1 + 1][i];
1193 
1194  if (fabs(h_) <= hmin * 1.00001 || *ncf == mxncf)
1195  {
1196  *corflag = 2;
1197  return;
1198  }
1199  *corflag = 1;
1200  *rh = 0.25;
1201  ipup = miter;
1202 }
1203 
1204 /*
1205  This routine manages the solution of the linear system arising from
1206  a chord iteration. It is called if miter != 0.
1207  If miter is 2, it calls dgesl1 to accomplish this.
1208  If miter is 5, it calls dgbsl.
1209 
1210  y = the right-hand side vector on input, and the solution vector
1211  on output.
1212 */
1213 void LSODA::solsy(vector<double> &y)
1214 {
1215  iersl = 0;
1216  if (miter != 2)
1217  {
1218  //printf("solsy -- miter != 2\n");
1219  REprintf("solsy -- miter != 2.\n");
1220  return;
1221  }
1222  if (miter == 2)
1223  dgesl1(wm_, n, ipvt, y, 0);
1224  return;
1225 }
1226 
1227 void LSODA::methodswitch(double dsm, double pnorm, double *pdh, double *rh)
1228 {
1229  int lm1, lm1p1, lm2, lm2p1, nqm1, nqm2;
1230  double rh1, rh2, rh1it, exm2, dm2, exm1, dm1, alpha, exsm;
1231 
1232  /*
1233  We are current using an Adams method. Consider switching to bdf.
1234  If the current order is greater than 5, assume the problem is
1235  not stiff, and skip this section.
1236  If the Lipschitz constant and error estimate are not polluted
1237  by roundoff, perform the usual test.
1238  Otherwise, switch to the bdf methods if the last step was
1239  restricted to insure stability ( irflag = 1 ), and stay with Adams
1240  method if not. When switching to bdf with polluted error estimates,
1241  in the absence of other information, double the step size.
1242 
1243  When the estimates are ok, we make the usual test by computing
1244  the step size we could have (ideally) used on this step,
1245  with the current (Adams) method, and also that for the bdf.
1246  If nq > mxords, we consider changing to order mxords on switching.
1247  Compare the two step sizes to decide whether to switch.
1248  The step size advantage must be at least ratio = 5 to switch.
1249  */
1250  if (meth_ == 1)
1251  {
1252  if (nq > 5)
1253  return;
1254  if (dsm <= (100. * pnorm * ETA) || pdest == 0.)
1255  {
1256  if (irflag == 0)
1257  return;
1258  rh2 = 2.;
1259  nqm2 = min(nq, mxords);
1260  }
1261  else
1262  {
1263  exsm = 1. / (double)l;
1264  rh1 = 1. / (1.2 * pow(dsm, exsm) + 0.0000012);
1265  rh1it = 2. * rh1;
1266  *pdh = pdlast * fabs(h_);
1267  if ((*pdh * rh1) > 0.00001)
1268  rh1it = sm1[nq] / *pdh;
1269  rh1 = min(rh1, rh1it);
1270  if (nq > mxords)
1271  {
1272  nqm2 = mxords;
1273  lm2 = mxords + 1;
1274  exm2 = 1. / (double)lm2;
1275  lm2p1 = lm2 + 1;
1276  dm2 = vmnorm(n, yh_[lm2p1], ewt) / cm2[mxords];
1277  rh2 = 1. / (1.2 * pow(dm2, exm2) + 0.0000012);
1278  }
1279  else
1280  {
1281  dm2 = dsm * (cm1[nq] / cm2[nq]);
1282  rh2 = 1. / (1.2 * pow(dm2, exsm) + 0.0000012);
1283  nqm2 = nq;
1284  }
1285  if (rh2 < ratio * rh1)
1286  return;
1287  }
1288  /*
1289  The method switch test passed. Reset relevant quantities for bdf.
1290  */
1291  *rh = rh2;
1292  icount = 20;
1293  meth_ = 2;
1294  miter = jtyp;
1295  pdlast = 0.;
1296  nq = nqm2;
1297  l = nq + 1;
1298  return;
1299  } /* end if ( meth_ == 1 ) */
1300 
1301  /*
1302  We are currently using a bdf method, considering switching to Adams.
1303  Compute the step size we could have (ideally) used on this step,
1304  with the current (bdf) method, and also that for the Adams.
1305  If nq > mxordn, we consider changing to order mxordn on switching.
1306  Compare the two step sizes to decide whether to switch.
1307  The step size advantage must be at least 5/ratio = 1 to switch.
1308  If the step size for Adams would be so small as to cause
1309  roundoff pollution, we stay with bdf.
1310  */
1311  exsm = 1. / (double)l;
1312  if (mxordn < nq)
1313  {
1314  nqm1 = mxordn;
1315  lm1 = mxordn + 1;
1316  exm1 = 1. / (double)lm1;
1317  lm1p1 = lm1 + 1;
1318  dm1 = vmnorm(n, yh_[lm1p1], ewt) / cm1[mxordn];
1319  rh1 = 1. / (1.2 * pow(dm1, exm1) + 0.0000012);
1320  }
1321  else
1322  {
1323  dm1 = dsm * (cm2[nq] / cm1[nq]);
1324  rh1 = 1. / (1.2 * pow(dm1, exsm) + 0.0000012);
1325  nqm1 = nq;
1326  exm1 = exsm;
1327  }
1328  rh1it = 2. * rh1;
1329  *pdh = pdnorm * fabs(h_);
1330  if ((*pdh * rh1) > 0.00001)
1331  rh1it = sm1[nqm1] / *pdh;
1332  rh1 = min(rh1, rh1it);
1333  rh2 = 1. / (1.2 * pow(dsm, exsm) + 0.0000012);
1334  if ((rh1 * ratio) < (5. * rh2))
1335  return;
1336  alpha = max(0.001, rh1);
1337  dm1 *= pow(alpha, exm1);
1338  if (dm1 <= 1000. * ETA * pnorm)
1339  return;
1340  /*
1341  The switch test passed. Reset relevant quantities for Adams.
1342  */
1343  *rh = rh1;
1344  icount = 20;
1345  meth_ = 1;
1346  miter = 0;
1347  pdlast = 0.;
1348  nq = nqm1;
1349  l = nq + 1;
1350 } /* end methodswitch */
1351 
1352 /*
1353  This routine returns from stoda to lsoda. Hence freevectors() is
1354  not executed.
1355 */
1357 {
1358  double r = 1. / tesco[nqu][2];
1359  for (size_t i = 1; i <= n; ++i)
1360  acor[i] *= r;
1361  hold = h_;
1362  jstart = 1;
1363 }
1364 
1365 /*
1366  Regardless of the success or failure of the step, factors
1367  rhdn, rhsm, and rhup are computed, by which h_ could be multiplied
1368  at order nq - 1, order nq, or order nq + 1, respectively.
1369  In the case of a failure, rhup = 0. to avoid an order increase.
1370  The largest of these is determined and the new order chosen
1371  accordingly. If the order is to be increased, we compute one
1372  additional scaled derivative.
1373 
1374  orderflag = 0 : no change in h_ or nq,
1375  1 : change in h_ but not nq,
1376  2 : change in both h_ and nq.
1377 */
1378 void LSODA::orderswitch(double *rhup, double dsm, double *pdh, double *rh,
1379  size_t *orderflag)
1380 {
1381  size_t newq = 0;
1382  double exsm, rhdn, rhsm, ddn, exdn, r;
1383 
1384  *orderflag = 0;
1385 
1386  exsm = 1. / (double)l;
1387  rhsm = 1. / (1.2 * pow(dsm, exsm) + 0.0000012);
1388 
1389  rhdn = 0.;
1390  if (nq != 1)
1391  {
1392  ddn = vmnorm(n, yh_[l], ewt) / tesco[nq][1];
1393  exdn = 1. / (double)nq;
1394  rhdn = 1. / (1.3 * pow(ddn, exdn) + 0.0000013);
1395  }
1396  /*
1397  If meth_ = 1, limit rh accordinfg to the stability region also.
1398  */
1399  if (meth_ == 1)
1400  {
1401  *pdh = max(fabs(h_) * pdlast, 0.000001);
1402  if (l < lmax)
1403  *rhup = min(*rhup, sm1[l] / *pdh);
1404  rhsm = min(rhsm, sm1[nq] / *pdh);
1405  if (nq > 1)
1406  rhdn = min(rhdn, sm1[nq - 1] / *pdh);
1407  pdest = 0.;
1408  }
1409  if (rhsm >= *rhup)
1410  {
1411  if (rhsm >= rhdn)
1412  {
1413  newq = nq;
1414  *rh = rhsm;
1415  }
1416  else
1417  {
1418  newq = nq - 1;
1419  *rh = rhdn;
1420  if (kflag < 0 && *rh > 1.)
1421  *rh = 1.;
1422  }
1423  }
1424  else
1425  {
1426  if (*rhup <= rhdn)
1427  {
1428  newq = nq - 1;
1429  *rh = rhdn;
1430  if (kflag < 0 && *rh > 1.)
1431  *rh = 1.;
1432  }
1433  else
1434  {
1435  *rh = *rhup;
1436  if (*rh >= 1.1)
1437  {
1438  r = el[l] / (double)l;
1439  nq = l;
1440  l = nq + 1;
1441  for (size_t i = 1; i <= n; ++i)
1442  yh_[l][i] = acor[i] * r;
1443 
1444  *orderflag = 2;
1445  return;
1446  }
1447  else
1448  {
1449  ialth = 3;
1450  return;
1451  }
1452  }
1453  }
1454  /*
1455  If meth_ = 1 and h_ is restricted by stability, bypass 10 percent test.
1456  */
1457  if (1 == meth_)
1458  {
1459  if ((*rh * *pdh * 1.00001) < sm1[newq])
1460  if (kflag == 0 && *rh < 1.1)
1461  {
1462  ialth = 3;
1463  return;
1464  }
1465  }
1466  else
1467  {
1468  if (kflag == 0 && *rh < 1.1)
1469  {
1470  ialth = 3;
1471  return;
1472  }
1473  }
1474  if (kflag <= -2)
1475  *rh = min(*rh, 0.2);
1476  /*
1477  If there is a change of order, reset nq, l, and the coefficients.
1478  In any case h_ is reset according to rh and the yh_ array is rescaled.
1479  Then exit or redo the step.
1480  */
1481  if (newq == nq)
1482  {
1483  *orderflag = 1;
1484  return;
1485  }
1486  nq = newq;
1487  l = nq + 1;
1488  *orderflag = 2;
1489 
1490 } /* end orderswitch */
1491 
1493 /*
1494  The el vector and related constants are reset
1495  whenever the order nq is changed, or at the start of the problem.
1496 */
1497 {
1498  array<double, 14> ep1;
1499 
1500  ep1 = elco[nq];
1501  for (size_t i = 1; i <= l; ++i)
1502  el[i] = ep1[i];
1503  rc = rc * el[1] / el0;
1504  el0 = el[1];
1505  conit = 0.5 / (double)(nq + 2);
1506 }
1507 
1508 #endif
LSODA::daxpy1
void daxpy1(const double da, const vector< double > &dx, vector< double > &dy, const size_t n, const size_t offsetX, const size_t offsetY)
Definition: lsoda_functions.h:81
LSODA::ipvt
vector< int > ipvt
Definition: LSODA.h:159
LSODA::nslp
size_t nslp
Definition: LSODA.h:177
LSODA::orderswitch
void orderswitch(double *rhup, double dsm, double *pdh, double *rh, size_t *orderflag)
Definition: lsoda_functions.h:1378
LSODA::nje
size_t nje
Definition: LSODA.h:167
LSODA::jstart
int jstart
Definition: LSODA.h:136
LSODA::l
size_t l
Definition: LSODA.h:162
LSODA::endstoda
void endstoda(void)
Definition: lsoda_functions.h:1356
ETA
#define ETA
Definition: lsoda_functions.h:6
LSODA::ewset
void ewset(const vector< double > &ycur)
Definition: lsoda_functions.h:645
LSODA::itol_
int itol_
Definition: LSODA.h:182
LSODA::irflag
int irflag
Definition: LSODA.h:179
LSODA::h_
double h_
Definition: LSODA.h:171
LSODA::sqrteta
double sqrteta
Definition: LSODA.h:139
LSODA::ixpr_
void ixpr_(const int value)
Definition: lsoda_functions.h:24
LSODA::el0
double el0
Definition: LSODA.h:171
LSODA::mxords
size_t mxords
Definition: LSODA.h:164
LSODA::nst
size_t nst
Definition: LSODA.h:167
LSODA::mused
size_t mused
Definition: LSODA.h:164
LSODA::pdest
double pdest
Definition: LSODA.h:178
LSODA::hmin
double hmin
Definition: LSODA.h:172
LSODA::intdy
void intdy(double t, int k, vector< double > &dky, int *iflag)
Definition: lsoda_functions.h:691
LSODA::resetcoeff
void resetcoeff(void)
Definition: lsoda_functions.h:1492
LSODA::msbp
size_t msbp
Definition: LSODA.h:162
LSODA::ratio
double ratio
Definition: LSODA.h:178
LSODA::correction
void correction(const size_t neq, vector< double > &y, LSODA_ODE_SYSTEM_TYPE f, size_t *corflag, double pnorm, double *del, double *delp, double *told, size_t *ncf, double *rh, size_t *m, dtype _data)
Definition: lsoda_functions.h:1023
LSODA::rc
double rc
Definition: LSODA.h:172
LSODA::ialth
size_t ialth
Definition: LSODA.h:176
LSODA::mxordn
size_t mxordn
Definition: LSODA.h:164
LSODA::sm1
array< double, 13 > sm1
Definition: LSODA.h:145
LSODA::mxhnil_
void mxhnil_(const int value)
Definition: lsoda_functions.h:29
LSODA::kflag
int kflag
Definition: LSODA.h:136
LSODA::dgesl1
void dgesl1(const vector< vector< double >> &a, const size_t n, vector< int > &ipvt, vector< double > &b, const size_t job)
Definition: lsoda_functions.h:91
LSODA_ODE_SYSTEM_TYPE
void(* LSODA_ODE_SYSTEM_TYPE)(double t, double *y, double *dydt, dtype _data)
Definition: LSODA.h:35
LSODA::hmin_
void hmin_(const double value)
Definition: lsoda_functions.h:13
LSODA::Maxsteps
int Maxsteps
Definition: LSODA.h:55
LSODA::stoda
void stoda(const size_t neq, vector< double > &y, LSODA_ODE_SYSTEM_TYPE f, dtype _data)
Definition: lsoda_functions.h:272
LSODA::elco
array< array< double, 14 >, 13 > elco
Definition: LSODA.h:151
LSODA::vmnorm
double vmnorm(const size_t n, const vector< double > &v, const vector< double > &w)
Definition: lsoda_functions.h:985
LSODA::hmax_
void hmax_(const double value)
Definition: lsoda_functions.h:8
LSODA::pdlast
double pdlast
Definition: LSODA.h:178
LSODA::tesco
array< array< double, 4 >, 13 > tesco
Definition: LSODA.h:152
LSODA::acor
vector< double > acor
Definition: LSODA.h:156
LSODA::maxord
size_t maxord
Definition: LSODA.h:162
LSODA::crate
double crate
Definition: LSODA.h:174
LSODA::idamax1
size_t idamax1(const vector< double > &dx, const size_t n, const size_t offset)
Definition: lsoda_functions.h:40
dtype
odeproblem * dtype
Definition: lsoda_functions.h:2
LSODA::terminate2
void terminate2(vector< double > &y, double *t)
Definition: lsoda_functions.h:232
LSODA::maxsteps_
void maxsteps_(const int value)
Definition: lsoda_functions.h:18
LSODA::iworks
array< int, 7 > iworks
Definition: LSODA.h:64
LSODA::ierpj
size_t ierpj
Definition: LSODA.h:161
LSODA::savf
vector< double > savf
Definition: LSODA.h:155
LSODA::ewt
vector< double > ewt
Definition: LSODA.h:154
LSODA::cfode
void cfode(int meth_)
Definition: lsoda_functions.h:742
LSODA::n
size_t n
Definition: LSODA.h:167
LSODA::mxncf
size_t mxncf
Definition: LSODA.h:162
LSODA::nqu
size_t nqu
Definition: LSODA.h:167
LSODA::successreturn
void successreturn(vector< double > &y, double *t, int itask, int ihit, double tcrit, int *istate)
Definition: lsoda_functions.h:248
LSODA::nfe
size_t nfe
Definition: LSODA.h:167
LSODA::methodswitch
void methodswitch(double dsm, double pnorm, double *pdh, double *rh)
Definition: lsoda_functions.h:1227
LSODA::abs_compare
static bool abs_compare(double a, double b)
Definition: lsoda_functions.h:34
LSODA::el
array< double, 14 > el
Definition: LSODA.h:147
LSODA::lmax
size_t lmax
Definition: LSODA.h:176
LSODA::dscal1
void dscal1(const double da, vector< double > &dx, const size_t n, const size_t offset)
Definition: lsoda_functions.h:63
LSODA::scaleh
void scaleh(double *rh, double *pdh)
Definition: lsoda_functions.h:873
LSODA::jtyp
size_t jtyp
Definition: LSODA.h:164
LSODA::wm_
vector< vector< double > > wm_
Definition: LSODA.h:158
LSODA::meth_
size_t meth_
Definition: LSODA.h:165
LSODA::fnorm
double fnorm(int n, const vector< vector< double >> &a, const vector< double > &w)
Definition: lsoda_functions.h:994
LSODA::yh_
vector< vector< double > > yh_
Definition: LSODA.h:157
LSODA::cm2
array< double, 6 > cm2
Definition: LSODA.h:149
LSODA::tn_
double tn_
Definition: LSODA.h:172
LSODA::miter
size_t miter
Definition: LSODA.h:162
LSODA::corfailure
void corfailure(double *told, double *rh, size_t *ncf, size_t *corflag)
Definition: lsoda_functions.h:1184
LSODA::ipup
size_t ipup
Definition: LSODA.h:176
LSODA::rmax
double rmax
Definition: LSODA.h:174
LSODA::pdnorm
double pdnorm
Definition: LSODA.h:173
LSODA::iersl
size_t iersl
Definition: LSODA.h:161
LSODA::ccmax
double ccmax
Definition: LSODA.h:171
LSODA::conit
double conit
Definition: LSODA.h:174
LSODA::dgefa1
void dgefa1(vector< vector< double >> &a, const size_t n, vector< int > &ipvt, size_t *const info)
Definition: lsoda_functions.h:154
LSODA::rtol_
std::vector< double > rtol_
Definition: LSODA.h:183
LSODA::terminate
void terminate(int *istate)
Definition: lsoda_functions.h:215
LSODA::iopt
int iopt
Definition: LSODA.h:50
LSODA::prja
void prja(const size_t neq, vector< double > &y, LSODA_ODE_SYSTEM_TYPE f, dtype _data)
Definition: lsoda_functions.h:912
LSODA::cm1
array< double, 13 > cm1
Definition: LSODA.h:148
LSODA::atol_
std::vector< double > atol_
Definition: LSODA.h:184
LSODA::hu
double hu
Definition: LSODA.h:172
LSODA::jcur
size_t jcur
Definition: LSODA.h:161
LSODA::maxcor
size_t maxcor
Definition: LSODA.h:162
LSODA::hmxi
double hmxi
Definition: LSODA.h:172
LSODA::rworks
array< double, 4 > rworks
Definition: LSODA.h:65
LSODA::solsy
void solsy(vector< double > &y)
Definition: lsoda_functions.h:1213
LSODA::hold
double hold
Definition: LSODA.h:174
LSODA::illin
size_t illin
Definition: LSODA.h:161
odeproblem
Definition: odeproblem.h:94
LSODA::icount
int icount
Definition: LSODA.h:179
LSODA::itask
int itask
Definition: LSODA.h:52
LSODA::nq
size_t nq
Definition: LSODA.h:167
LSODA::ddot1
double ddot1(const vector< double > &a, const vector< double > &b, const size_t n, const size_t offsetA, const size_t offsetB)
Definition: lsoda_functions.h:71