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