GNU Octave  4.4.1
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
DASPK.cc
Go to the documentation of this file.
1 /*
2 
3 Copyright (C) 1996-2018 John W. Eaton
4 
5 This file is part of Octave.
6 
7 Octave is free software: you can redistribute it and/or modify it
8 under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11 
12 Octave is distributed in the hope that it will be useful, but
13 WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16 
17 You should have received a copy of the GNU General Public License
18 along with Octave; see the file COPYING. If not, see
19 <https://www.gnu.org/licenses/>.
20 
21 */
22 
23 #if defined (HAVE_CONFIG_H)
24 # include "config.h"
25 #endif
26 
27 #include <sstream>
28 
29 #include "DASPK.h"
30 #include "dMatrix.h"
31 #include "f77-fcn.h"
32 #include "lo-error.h"
33 #include "quit.h"
34 
35 typedef F77_INT (*daspk_fcn_ptr) (const double&, const double*, const double*,
36  const double&, double*, F77_INT&, double*,
37  F77_INT*);
38 
39 typedef F77_INT (*daspk_jac_ptr) (const double&, const double*, const double*,
40  double*, const double&, double*, F77_INT*);
41 
42 typedef F77_INT (*daspk_psol_ptr) (const F77_INT&, const double&,
43  const double*, const double*,
44  const double*, const double&,
45  const double*, double*, F77_INT*,
46  double*, const double&, F77_INT&,
47  double*, F77_INT*);
48 
49 extern "C"
50 {
51  F77_RET_T
52  F77_FUNC (ddaspk, DDASPK) (daspk_fcn_ptr, const F77_INT&, F77_DBLE&,
53  F77_DBLE*, F77_DBLE*, F77_DBLE&, const F77_INT*,
54  const F77_DBLE*, const F77_DBLE*, F77_INT&,
55  F77_DBLE*, const F77_INT&, F77_INT*,
56  const F77_INT&, const F77_DBLE*, const F77_INT*,
58 }
59 
62 static F77_INT nn;
63 
64 static F77_INT
65 ddaspk_f (const double& time, const double *state, const double *deriv,
66  const double&, double *delta, F77_INT& ires, double *, F77_INT *)
67 {
68  BEGIN_INTERRUPT_WITH_EXCEPTIONS;
69 
70  ColumnVector tmp_deriv (nn);
71  ColumnVector tmp_state (nn);
72  ColumnVector tmp_delta (nn);
73 
74  for (F77_INT i = 0; i < nn; i++)
75  {
76  tmp_deriv.elem (i) = deriv[i];
77  tmp_state.elem (i) = state[i];
78  }
79 
80  octave_idx_type tmp_ires = ires;
81 
82  tmp_delta = user_fun (tmp_state, tmp_deriv, time, tmp_ires);
83 
84  ires = octave::to_f77_int (tmp_ires);
85 
86  if (ires >= 0)
87  {
88  if (tmp_delta.isempty ())
89  ires = -2;
90  else
91  {
92  for (F77_INT i = 0; i < nn; i++)
93  delta[i] = tmp_delta.elem (i);
94  }
95  }
96 
97  END_INTERRUPT_WITH_EXCEPTIONS;
98 
99  return 0;
100 }
101 
102 //NEQ, T, Y, YPRIME, SAVR, WK, CJ, WGHT,
103 //C WP, IWP, B, EPLIN, IER, RPAR, IPAR)
104 
105 static F77_INT
106 ddaspk_psol (const F77_INT&, const double&, const double *,
107  const double *, const double *, const double&,
108  const double *, double *, F77_INT *, double *,
109  const double&, F77_INT&, double *, F77_INT*)
110 {
111  BEGIN_INTERRUPT_WITH_EXCEPTIONS;
112 
113  (*current_liboctave_error_handler) ("daspk: PSOL is not implemented");
114 
115  END_INTERRUPT_WITH_EXCEPTIONS;
116 
117  return 0;
118 }
119 
120 static F77_INT
121 ddaspk_j (const double& time, const double *state, const double *deriv,
122  double *pd, const double& cj, double *, F77_INT *)
123 {
124  BEGIN_INTERRUPT_WITH_EXCEPTIONS;
125 
126  // FIXME: would be nice to avoid copying the data.
127 
128  ColumnVector tmp_state (nn);
129  ColumnVector tmp_deriv (nn);
130 
131  for (F77_INT i = 0; i < nn; i++)
132  {
133  tmp_deriv.elem (i) = deriv[i];
134  tmp_state.elem (i) = state[i];
135  }
136 
137  Matrix tmp_pd = user_jac (tmp_state, tmp_deriv, time, cj);
138 
139  for (F77_INT j = 0; j < nn; j++)
140  for (F77_INT i = 0; i < nn; i++)
141  pd[nn * j + i] = tmp_pd.elem (i, j);
142 
143  END_INTERRUPT_WITH_EXCEPTIONS;
144 
145  return 0;
146 }
147 
149 DASPK::do_integrate (double tout)
150 {
151  // FIXME: should handle all this option stuff just once for each new problem.
152 
154 
155  if (! initialized || restart || DAEFunc::reset || DASPK_options::reset)
156  {
157  integration_error = false;
158 
159  initialized = true;
160 
161  info.resize (dim_vector (20, 1));
162 
163  for (F77_INT i = 0; i < 20; i++)
164  info(i) = 0;
165 
166  F77_INT n = octave::to_f77_int (size ());
167 
168  nn = n;
169 
170  info(0) = 0;
171 
172  if (stop_time_set)
173  {
174  rwork(0) = stop_time;
175  info(3) = 1;
176  }
177  else
178  info(3) = 0;
179 
180  // DAEFunc
181 
184 
185  if (user_fun)
186  {
187  octave_idx_type ires = 0;
188 
189  ColumnVector res = (*user_fun) (x, xdot, t, ires);
190 
191  if (res.numel () != x.numel ())
192  {
193  // FIXME: Should this be a warning?
194  (*current_liboctave_error_handler)
195  ("daspk: inconsistent sizes for state and residual vectors");
196 
197  integration_error = true;
198  return retval;
199  }
200  }
201  else
202  {
203  // FIXME: Should this be a warning?
204  (*current_liboctave_error_handler)
205  ("daspk: no user supplied RHS subroutine!");
206 
207  integration_error = true;
208  return retval;
209  }
210 
211  info(4) = (user_jac ? 1 : 0);
212 
213  DAEFunc::reset = false;
214 
215  octave_idx_type eiq = enforce_inequality_constraints ();
216  octave_idx_type ccic = compute_consistent_initial_condition ();
217  octave_idx_type eavfet = exclude_algebraic_variables_from_error_test ();
218 
219  liw = 40 + n;
220  if (eiq == 1 || eiq == 3)
221  liw += n;
222  if (ccic == 1 || eavfet == 1)
223  liw += n;
224 
225  lrw = 50 + 9*n + n*n;
226  if (eavfet == 1)
227  lrw += n;
228 
229  iwork.resize (dim_vector (liw, 1));
230  rwork.resize (dim_vector (lrw, 1));
231 
232  // DASPK_options
233 
234  abs_tol = absolute_tolerance ();
235  rel_tol = relative_tolerance ();
236 
237  F77_INT abs_tol_len = octave::to_f77_int (abs_tol.numel ());
238  F77_INT rel_tol_len = octave::to_f77_int (rel_tol.numel ());
239 
240  if (abs_tol_len == 1 && rel_tol_len == 1)
241  {
242  info(1) = 0;
243  }
244  else if (abs_tol_len == n && rel_tol_len == n)
245  {
246  info(1) = 1;
247  }
248  else
249  {
250  // FIXME: Should this be a warning?
251  (*current_liboctave_error_handler)
252  ("daspk: inconsistent sizes for tolerance arrays");
253 
254  integration_error = true;
255  return retval;
256  }
257 
258  double hmax = maximum_step_size ();
259  if (hmax >= 0.0)
260  {
261  rwork(1) = hmax;
262  info(6) = 1;
263  }
264  else
265  info(6) = 0;
266 
267  double h0 = initial_step_size ();
268  if (h0 >= 0.0)
269  {
270  rwork(2) = h0;
271  info(7) = 1;
272  }
273  else
274  info(7) = 0;
275 
276  octave_idx_type maxord = maximum_order ();
277  if (maxord >= 0)
278  {
279  if (maxord > 0 && maxord < 6)
280  {
281  info(8) = 1;
282  iwork(2) = octave::to_f77_int (maxord);
283  }
284  else
285  {
286  // FIXME: Should this be a warning?
287  (*current_liboctave_error_handler)
288  ("daspk: invalid value for maximum order");
289  integration_error = true;
290  return retval;
291  }
292  }
293 
294  switch (eiq)
295  {
296  case 1:
297  case 3:
298  {
299  Array<octave_idx_type> ict = inequality_constraint_types ();
300 
301  F77_INT ict_nel = octave::to_f77_int (ict.numel ());
302 
303  if (ict_nel == n)
304  {
305  for (F77_INT i = 0; i < n; i++)
306  {
307  F77_INT val = octave::to_f77_int (ict(i));
308  if (val < -2 || val > 2)
309  {
310  // FIXME: Should this be a warning?
311  (*current_liboctave_error_handler)
312  ("daspk: invalid value for inequality constraint type");
313  integration_error = true;
314  return retval;
315  }
316  iwork(40+i) = val;
317  }
318  }
319  else
320  {
321  // FIXME: Should this be a warning?
322  (*current_liboctave_error_handler)
323  ("daspk: inequality constraint types size mismatch");
324  integration_error = true;
325  return retval;
326  }
327  }
328  OCTAVE_FALLTHROUGH;
329 
330  case 0:
331  case 2:
332  info(9) = octave::to_f77_int (eiq);
333  break;
334 
335  default:
336  // FIXME: Should this be a warning?
337  (*current_liboctave_error_handler)
338  ("daspk: invalid value for enforce inequality constraints option");
339  integration_error = true;
340  return retval;
341  }
342 
343  if (ccic)
344  {
345  if (ccic == 1)
346  {
347  // FIXME: this code is duplicated below.
348 
349  Array<octave_idx_type> av = algebraic_variables ();
350 
351  F77_INT av_nel = octave::to_f77_int (av.numel ());
352 
353  if (av_nel == n)
354  {
355  F77_INT lid;
356  if (eiq == 0 || eiq == 2)
357  lid = 40;
358  else if (eiq == 1 || eiq == 3)
359  lid = 40 + n;
360  else
362  ("daspk: invalid value for eiq: %d", eiq);
363 
364  for (F77_INT i = 0; i < n; i++)
365  iwork(lid+i) = (av(i) ? -1 : 1);
366  }
367  else
368  {
369  // FIXME: Should this be a warning?
370  (*current_liboctave_error_handler)
371  ("daspk: algebraic variables size mismatch");
372  integration_error = true;
373  return retval;
374  }
375  }
376  else if (ccic != 2)
377  {
378  // FIXME: Should this be a warning?
379  (*current_liboctave_error_handler)
380  ("daspk: invalid value for compute consistent initial condition option");
381  integration_error = true;
382  return retval;
383  }
384 
385  info(10) = octave::to_f77_int (ccic);
386  }
387 
388  if (eavfet)
389  {
390  info(15) = 1;
391 
392  // FIXME: this code is duplicated above.
393 
394  Array<octave_idx_type> av = algebraic_variables ();
395 
396  F77_INT av_nel = octave::to_f77_int (av.numel ());
397 
398  if (av_nel == n)
399  {
400  F77_INT lid;
401  if (eiq == 0 || eiq == 2)
402  lid = 40;
403  else if (eiq == 1 || eiq == 3)
404  lid = 40 + n;
405  else
407  ("daspk: invalid value for eiq: %d", eiq);
408 
409  for (F77_INT i = 0; i < n; i++)
410  iwork(lid+i) = (av(i) ? -1 : 1);
411  }
412  }
413 
414  if (use_initial_condition_heuristics ())
415  {
416  Array<double> ich = initial_condition_heuristics ();
417 
418  if (ich.numel () == 6)
419  {
420  iwork(31) = octave::to_f77_int (octave::math::nint_big (ich(0)));
421  iwork(32) = octave::to_f77_int (octave::math::nint_big (ich(1)));
422  iwork(33) = octave::to_f77_int (octave::math::nint_big (ich(2)));
423  iwork(34) = octave::to_f77_int (octave::math::nint_big (ich(3)));
424 
425  rwork(13) = ich(4);
426  rwork(14) = ich(5);
427  }
428  else
429  {
430  // FIXME: Should this be a warning?
431  (*current_liboctave_error_handler)
432  ("daspk: invalid initial condition heuristics option");
433  integration_error = true;
434  return retval;
435  }
436 
437  info(16) = 1;
438  }
439 
440  octave_idx_type pici = print_initial_condition_info ();
441  switch (pici)
442  {
443  case 0:
444  case 1:
445  case 2:
446  info(17) = octave::to_f77_int (pici);
447  break;
448 
449  default:
450  // FIXME: Should this be a warning?
451  (*current_liboctave_error_handler)
452  ("daspk: invalid value for print initial condition info option");
453  integration_error = true;
454  return retval;
455  break;
456  }
457 
458  DASPK_options::reset = false;
459 
460  restart = false;
461  }
462 
463  double *px = x.fortran_vec ();
464  double *pxdot = xdot.fortran_vec ();
465 
466  F77_INT *pinfo = info.fortran_vec ();
467 
468  double *prel_tol = rel_tol.fortran_vec ();
469  double *pabs_tol = abs_tol.fortran_vec ();
470 
471  double *prwork = rwork.fortran_vec ();
472  F77_INT *piwork = iwork.fortran_vec ();
473 
474  double *dummy = nullptr;
475  F77_INT *idummy = nullptr;
476 
477  F77_INT tmp_istate = octave::to_f77_int (istate);
478 
479  F77_XFCN (ddaspk, DDASPK, (ddaspk_f, nn, t, px, pxdot, tout, pinfo,
480  prel_tol, pabs_tol, tmp_istate, prwork, lrw,
481  piwork, liw, dummy, idummy, ddaspk_j,
482  ddaspk_psol));
483 
484  istate = tmp_istate;
485 
486  switch (istate)
487  {
488  case 1: // A step was successfully taken in intermediate-output
489  // mode. The code has not yet reached TOUT.
490  case 2: // The integration to TSTOP was successfully completed
491  // (T=TSTOP) by stepping exactly to TSTOP.
492  case 3: // The integration to TOUT was successfully completed
493  // (T=TOUT) by stepping past TOUT. Y(*) is obtained by
494  // interpolation. YPRIME(*) is obtained by interpolation.
495  case 4: // The initial condition calculation, with
496  // INFO(11) > 0, was successful, and INFO(14) = 1.
497  // No integration steps were taken, and the solution
498  // is not considered to have been started.
499  retval = x;
500  t = tout;
501  break;
502 
503  case -1: // A large amount of work has been expended. (~500 steps).
504  case -2: // The error tolerances are too stringent.
505  case -3: // The local error test cannot be satisfied because you
506  // specified a zero component in ATOL and the
507  // corresponding computed solution component is zero.
508  // Thus, a pure relative error test is impossible for
509  // this component.
510  case -6: // DDASPK had repeated error test failures on the last
511  // attempted step.
512  case -7: // The corrector could not converge.
513  case -8: // The matrix of partial derivatives is singular.
514  case -9: // The corrector could not converge. There were repeated
515  // error test failures in this step.
516  case -10: // The corrector could not converge because IRES was
517  // equal to minus one.
518  case -11: // IRES equal to -2 was encountered and control is being
519  // returned to the calling program.
520  case -12: // DDASPK failed to compute the initial YPRIME.
521  case -13: // Unrecoverable error encountered inside user's
522  // PSOL routine, and control is being returned to
523  // the calling program.
524  case -14: // The Krylov linear system solver could not
525  // achieve convergence.
526  case -33: // The code has encountered trouble from which it cannot
527  // recover. A message is printed explaining the trouble
528  // and control is returned to the calling program. For
529  // example, this occurs when invalid input is detected.
530  integration_error = true;
531  break;
532 
533  default:
534  integration_error = true;
535  (*current_liboctave_error_handler)
536  ("unrecognized value of istate (= %d) returned from ddaspk", istate);
537  break;
538  }
539 
540  return retval;
541 }
542 
543 Matrix
545 {
546  Matrix dummy;
547  return integrate (tout, dummy);
548 }
549 
550 Matrix
551 DASPK::integrate (const ColumnVector& tout, Matrix& xdot_out)
552 {
553  Matrix retval;
554 
555  octave_idx_type n_out = tout.numel ();
556  F77_INT n = octave::to_f77_int (size ());
557 
558  if (n_out > 0 && n > 0)
559  {
560  retval.resize (n_out, n);
561  xdot_out.resize (n_out, n);
562 
563  for (F77_INT i = 0; i < n; i++)
564  {
565  retval.elem (0, i) = x.elem (i);
566  xdot_out.elem (0, i) = xdot.elem (i);
567  }
568 
569  for (octave_idx_type j = 1; j < n_out; j++)
570  {
571  ColumnVector x_next = do_integrate (tout.elem (j));
572 
573  if (integration_error)
574  return retval;
575 
576  for (F77_INT i = 0; i < n; i++)
577  {
578  retval.elem (j, i) = x_next.elem (i);
579  xdot_out.elem (j, i) = xdot.elem (i);
580  }
581  }
582  }
583 
584  return retval;
585 }
586 
587 Matrix
588 DASPK::do_integrate (const ColumnVector& tout, const ColumnVector& tcrit)
589 {
590  Matrix dummy;
591  return integrate (tout, dummy, tcrit);
592 }
593 
594 Matrix
595 DASPK::integrate (const ColumnVector& tout, Matrix& xdot_out,
596  const ColumnVector& tcrit)
597 {
598  Matrix retval;
599 
600  octave_idx_type n_out = tout.numel ();
601  F77_INT n = octave::to_f77_int (size ());
602 
603  if (n_out > 0 && n > 0)
604  {
605  retval.resize (n_out, n);
606  xdot_out.resize (n_out, n);
607 
608  for (F77_INT i = 0; i < n; i++)
609  {
610  retval.elem (0, i) = x.elem (i);
611  xdot_out.elem (0, i) = xdot.elem (i);
612  }
613 
614  octave_idx_type n_crit = tcrit.numel ();
615 
616  if (n_crit > 0)
617  {
618  octave_idx_type i_crit = 0;
619  octave_idx_type i_out = 1;
620  double next_crit = tcrit.elem (0);
621  double next_out;
622  while (i_out < n_out)
623  {
624  bool do_restart = false;
625 
626  next_out = tout.elem (i_out);
627  if (i_crit < n_crit)
628  next_crit = tcrit.elem (i_crit);
629 
630  bool save_output;
631  double t_out;
632 
633  if (next_crit == next_out)
634  {
635  set_stop_time (next_crit);
636  t_out = next_out;
637  save_output = true;
638  i_out++;
639  i_crit++;
640  do_restart = true;
641  }
642  else if (next_crit < next_out)
643  {
644  if (i_crit < n_crit)
645  {
646  set_stop_time (next_crit);
647  t_out = next_crit;
648  save_output = false;
649  i_crit++;
650  do_restart = true;
651  }
652  else
653  {
654  clear_stop_time ();
655  t_out = next_out;
656  save_output = true;
657  i_out++;
658  }
659  }
660  else
661  {
662  set_stop_time (next_crit);
663  t_out = next_out;
664  save_output = true;
665  i_out++;
666  }
667 
668  ColumnVector x_next = do_integrate (t_out);
669 
670  if (integration_error)
671  return retval;
672 
673  if (save_output)
674  {
675  for (F77_INT i = 0; i < n; i++)
676  {
677  retval.elem (i_out-1, i) = x_next.elem (i);
678  xdot_out.elem (i_out-1, i) = xdot.elem (i);
679  }
680  }
681 
682  if (do_restart)
683  force_restart ();
684  }
685  }
686  else
687  {
688  retval = integrate (tout, xdot_out);
689 
690  if (integration_error)
691  return retval;
692  }
693  }
694 
695  return retval;
696 }
697 
700 {
702 
703  std::ostringstream buf;
704  buf << t;
705  std::string t_curr = buf.str ();
706 
707  switch (istate)
708  {
709  case 1:
710  retval = "a step was successfully taken in intermediate-output mode.";
711  break;
712 
713  case 2:
714  retval = "integration completed by stepping exactly to TOUT";
715  break;
716 
717  case 3:
718  retval = "integration to tout completed by stepping past TOUT";
719  break;
720 
721  case 4:
722  retval = "initial condition calculation completed successfully";
723  break;
724 
725  case -1:
726  retval = "a large amount of work has been expended (t =" + t_curr + ')';
727  break;
728 
729  case -2:
730  retval = "the error tolerances are too stringent";
731  break;
732 
733  case -3:
734  retval = "error weight became zero during problem. (t = " + t_curr +
735  "; solution component i vanished, and atol or atol(i) == 0)";
736  break;
737 
738  case -6:
739  retval = "repeated error test failures on the last attempted step (t = "
740  + t_curr + ')';
741  break;
742 
743  case -7:
744  retval = "the corrector could not converge (t = " + t_curr + ')';
745  break;
746 
747  case -8:
748  retval = "the matrix of partial derivatives is singular (t = " + t_curr +
749  ')';
750  break;
751 
752  case -9:
753  retval = "the corrector could not converge (t = " + t_curr +
754  "; repeated test failures)";
755  break;
756 
757  case -10:
758  retval = "corrector could not converge because IRES was -1 (t = "
759  + t_curr + ')';
760  break;
761 
762  case -11:
763  retval = "return requested in user-supplied function (t = " + t_curr +
764  ')';
765  break;
766 
767  case -12:
768  retval = "failed to compute consistent initial conditions";
769  break;
770 
771  case -13:
772  retval = "unrecoverable error encountered inside user's PSOL function (t = "
773  + t_curr + ')';
774  break;
775 
776  case -14:
777  retval = "the Krylov linear system solver failed to converge (t = "
778  + t_curr + ')';
779  break;
780 
781  case -33:
782  retval = "unrecoverable error (see printed message)";
783  break;
784 
785  default:
786  retval = "unknown error state";
787  break;
788  }
789 
790  return retval;
791 }
static DAEFunc::DAERHSFunc user_fun
Definition: DASPK.cc:60
F77_INT(* daspk_jac_ptr)(const double &, const double *, const double *, double *, const double &, double *, F77_INT *)
Definition: DASPK.cc:39
ColumnVector xdot
Definition: base-dae.h:77
ColumnVector do_integrate(double t)
Definition: DASPK.cc:149
Matrix integrate(const ColumnVector &tout, Matrix &xdot_out)
Definition: DASPK.cc:551
void resize(octave_idx_type nr, octave_idx_type nc, double rfv=0)
Definition: dMatrix.h:148
octave_f77_int_type lrw
Definition: DASPK.h:75
bool isempty(void) const
Definition: Array.h:565
identity matrix If supplied two scalar respectively For allows like xample val
Definition: data.cc:4986
std::string error_message(void) const
Definition: DASPK.cc:699
OCTAVE_NORETURN liboctave_error_handler current_liboctave_error_handler
Definition: lo-error.c:38
static DAEFunc::DAEJacFunc user_jac
Definition: DASPK.cc:61
const T * fortran_vec(void) const
Definition: Array.h:584
bool integration_error
Definition: base-de.h:115
static F77_INT nn
Definition: DASPK.cc:62
octave_idx_type istate
Definition: base-de.h:117
T & elem(octave_idx_type n)
Definition: Array.h:488
F77_INT(* daspk_psol_ptr)(const F77_INT &, const double &, const double *, const double *, const double *, const double &, const double *, double *, F77_INT *, double *, const double &, F77_INT &, double *, F77_INT *)
Definition: DASPK.cc:42
#define F77_XFCN(f, F, args)
Definition: f77-fcn.h:41
subroutine DDASPK(RES, NEQ, T, Y, YPRIME, TOUT, INFO, RTOL, ATOL, IDID, RWORK, LRW, IWORK, LIW, RPAR, IPAR, JAC, PSOL)
Definition: ddaspk.f:7
octave_f77_int_type liw
Definition: DASPK.h:74
F77_RET_T(F77_CONST_CHAR_ARG_DECL, F77_CONST_CHAR_ARG_DECL, F77_CONST_CHAR_ARG_DECL, const F77_INT &, const F77_INT &, const F77_INT &, F77_INT &, F77_INT &, F77_DBLE *, const F77_INT &, F77_DBLE *, const F77_INT &, F77_DBLE *, F77_DBLE *, F77_DBLE *, const F77_INT &, F77_DBLE *, const F77_INT &, F77_DBLE *, const F77_INT &, F77_DBLE *, F77_INT *, F77_INT &F77_CHAR_ARG_LEN_DECL F77_CHAR_ARG_LEN_DECL F77_CHAR_ARG_LEN_DECL)
static F77_INT ddaspk_f(const double &time, const double *state, const double *deriv, const double &, double *delta, F77_INT &ires, double *, F77_INT *)
Definition: DASPK.cc:65
octave_value resize(const dim_vector &dv, bool fill=false) const
Definition: ov.h:511
Array< double > rwork
Definition: DASPK.h:80
bool reset
Definition: DAEFunc.h:101
octave_idx_type size(void) const
Definition: base-de.h:76
Matrix(* DAEJacFunc)(const ColumnVector &x, const ColumnVector &xdot, double t, double cj)
Definition: DAEFunc.h:44
virtual void force_restart(void)
Definition: base-de.h:95
void set_stop_time(double tt)
Definition: base-de.h:82
double t
Definition: base-de.h:107
DAEJacFunc jacobian_function(void) const
Definition: DAEFunc.h:82
F77_INT(* daspk_fcn_ptr)(const double &, const double *, const double *, const double &, double *, F77_INT &, double *, F77_INT *)
Definition: DASPK.cc:35
void resize(const dim_vector &dv, const T &rfv)
Resizing (with fill).
Definition: Array.cc:1010
octave_value retval
Definition: data.cc:6246
Definition: dMatrix.h:36
void clear_stop_time(void)
Definition: base-de.h:89
static uint32_t state[624]
Definition: randmtzig.cc:183
static F77_INT ddaspk_j(const double &time, const double *state, const double *deriv, double *pd, const double &cj, double *, F77_INT *)
Definition: DASPK.cc:121
Array< octave_f77_int_type > iwork
Definition: DASPK.h:78
double stop_time
Definition: base-de.h:109
Array< octave_f77_int_type > info
Definition: DASPK.h:77
ColumnVector x
Definition: base-de.h:105
F77_RET_T F77_FUNC(ddaspk, DDASPK)(daspk_fcn_ptr
octave_f77_int_type F77_INT
Definition: f77-fcn.h:305
Array< double > abs_tol
Definition: DASPK.h:82
octave_idx_type nint_big(double x)
Definition: lo-mappers.cc:182
static F77_INT ddaspk_psol(const F77_INT &, const double &, const double *, const double *, const double *, const double &, const double *, double *, F77_INT *, double *, const double &, F77_INT &, double *, F77_INT *)
Definition: DASPK.cc:106
DAERHSFunc function(void) const
Definition: DAEFunc.h:73
bool initialized
Definition: DASPK.h:72
ColumnVector(* DAERHSFunc)(const ColumnVector &x, const ColumnVector &xdot, double t, octave_idx_type &ires)
Definition: DAEFunc.h:36
for i
Definition: data.cc:5264
octave_idx_type numel(void) const
Number of elements in the array.
Definition: Array.h:366
bool restart
Definition: base-de.h:113
Vector representing the dimensions (size) of an Array.
Definition: dim-vector.h:87
If this string is the system will ring the terminal sometimes it is useful to be able to print the original representation of the string
Definition: utils.cc:888
bool stop_time_set
Definition: base-de.h:111
Array< double > rel_tol
Definition: DASPK.h:83