GNU Octave  9.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
__ode15__.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 2016-2024 The Octave Project Developers
4 //
5 // See the file COPYRIGHT.md in the top-level directory of this
6 // distribution or <https://octave.org/copyright/>.
7 //
8 // This file is part of Octave.
9 //
10 // Octave is free software: you can redistribute it and/or modify it
11 // under the terms of the GNU General Public License as published by
12 // the Free Software Foundation, either version 3 of the License, or
13 // (at your option) any later version.
14 //
15 // Octave is distributed in the hope that it will be useful, but
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 // GNU General Public License for more details.
19 //
20 // You should have received a copy of the GNU General Public License
21 // along with Octave; see the file COPYING. If not, see
22 // <https://www.gnu.org/licenses/>.
23 //
24 ////////////////////////////////////////////////////////////////////////
25 
26 #if defined (HAVE_CONFIG_H)
27 # include "config.h"
28 #endif
29 
30 #include "dColVector.h"
31 #include "dMatrix.h"
32 #include "dSparse.h"
33 #include "f77-fcn.h"
34 #include "lo-utils.h"
35 
36 #include "Cell.h"
37 #include "defun-dld.h"
38 #include "error.h"
39 #include "errwarn.h"
40 #include "interpreter-private.h"
41 #include "interpreter.h"
42 #include "oct-map.h"
43 #include "ov.h"
44 #include "ovl.h"
45 #include "pager.h"
46 
47 #if defined (HAVE_SUNDIALS)
48 
49 # if defined (HAVE_NVECTOR_NVECTOR_SERIAL_H)
50 # include <nvector/nvector_serial.h>
51 # endif
52 
53 # if defined (HAVE_IDA_IDA_H)
54 # include <ida/ida.h>
55 # elif defined (HAVE_IDA_H)
56 # include <ida.h>
57 # endif
58 # if defined (HAVE_IDA_IDA_DIRECT_H)
59 # include <ida/ida_direct.h>
60 # elif defined (HAVE_IDA_DIRECT_H)
61 # include <ida_direct.h>
62 # endif
63 
64 # if defined (HAVE_SUNLINSOL_SUNLINSOL_DENSE_H)
65 # include <sunlinsol/sunlinsol_dense.h>
66 # endif
67 
68 # if defined (HAVE_SUNLINSOL_SUNLINSOL_KLU_H)
69 # if defined (HAVE_KLU_H)
70 # include <klu.h>
71 # endif
72 # if defined (HAVE_KLU_KLU_H)
73 # include <klu/klu.h>
74 # endif
75 # if defined (HAVE_SUITESPARSE_KLU_H)
76 # include <suitesparse/klu.h>
77 # endif
78 # if defined (HAVE_UFPARSE_KLU_H)
79 # include <ufsparse/klu.h>
80 # endif
81 # include <sunlinsol/sunlinsol_klu.h>
82 # endif
83 
84 #endif
85 
87 
88 #if defined (HAVE_SUNDIALS)
89 
90 # if ! defined (HAVE_IDASETJACFN) && defined (HAVE_IDADLSSETJACFN)
91 static inline int
92 IDASetJacFn (void *ida_mem, IDADlsJacFn jac)
93 {
94  return IDADlsSetJacFn (ida_mem, jac);
95 }
96 # endif
97 
98 # if ! defined (HAVE_IDASETLINEARSOLVER) && defined (HAVE_IDADLSSETLINEARSOLVER)
99 static inline int
100 IDASetLinearSolver (void *ida_mem, SUNLinearSolver LS, SUNMatrix A)
101 {
102  return IDADlsSetLinearSolver (ida_mem, LS, A);
103 }
104 # endif
105 
106 # if ! defined (HAVE_SUNLINSOL_DENSE) && defined (HAVE_SUNDENSELINEARSOLVER)
107 static inline SUNLinearSolver
108 SUNLinSol_Dense (N_Vector y, SUNMatrix A)
109 {
110  return SUNDenseLinearSolver (y, A);
111 }
112 # endif
113 
114 # if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
115 # if ! defined (HAVE_SUNLINSOL_KLU) && defined (HAVE_SUNKLU)
116 static inline SUNLinearSolver
117 SUNLinSol_KLU (N_Vector y, SUNMatrix A)
118 {
119  return SUNKLU (y, A);
120 }
121 # endif
122 # endif
123 
124 static inline OCTAVE_SUNREALTYPE *
125 nv_data_s (N_Vector& v)
126 {
127 # if defined (HAVE_PRAGMA_GCC_DIAGNOSTIC)
128  // Disable warning from GCC about old-style casts in Sundials
129  // macro expansions. Do this in a function so that this
130  // diagnostic may still be enabled for the rest of the file.
131 # pragma GCC diagnostic push
132 # pragma GCC diagnostic ignored "-Wold-style-cast"
133 # endif
134 
135  return NV_DATA_S (v);
136 
137 # if defined (HAVE_PRAGMA_GCC_DIAGNOSTIC)
138  // Restore prevailing warning state for remainder of the file.
139 # pragma GCC diagnostic pop
140 # endif
141 }
142 
143 class IDA
144 {
145 public:
146 
147  typedef
148  ColumnVector (*DAERHSFuncIDA) (const ColumnVector& x,
149  const ColumnVector& xdot,
150  OCTAVE_SUNREALTYPE t,
151  const octave_value& idaf);
152 
153  typedef
154  Matrix (*DAEJacFuncDense) (const ColumnVector& x,
155  const ColumnVector& xdot, OCTAVE_SUNREALTYPE t,
156  OCTAVE_SUNREALTYPE cj, const octave_value& idaj);
157 
158  typedef
159  SparseMatrix (*DAEJacFuncSparse) (const ColumnVector& x,
160  const ColumnVector& xdot,
161  OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
162  const octave_value& idaj);
163 
164  typedef
165  Matrix (*DAEJacCellDense) (Matrix *dfdy, Matrix *dfdyp,
166  OCTAVE_SUNREALTYPE cj);
167 
168  typedef
169  SparseMatrix (*DAEJacCellSparse) (SparseMatrix *dfdy,
170  SparseMatrix *dfdyp, OCTAVE_SUNREALTYPE cj);
171 
172  //Default
173  IDA ()
174  : m_t0 (0.0), m_y0 (), m_yp0 (), m_havejac (false), m_havejacfcn (false),
175  m_havejacsparse (false), m_mem (nullptr), m_num (), m_ida_fcn (),
176  m_ida_jac (), m_dfdy (nullptr), m_dfdyp (nullptr), m_spdfdy (nullptr),
177  m_spdfdyp (nullptr), m_fcn (nullptr), m_jacfcn (nullptr),
178  m_jacspfcn (nullptr), m_jacdcell (nullptr), m_jacspcell (nullptr),
179  m_sunJacMatrix (nullptr), m_sunLinearSolver (nullptr)
180  { }
181 
182 
183  IDA (OCTAVE_SUNREALTYPE t, ColumnVector y, ColumnVector yp,
184  const octave_value& ida_fcn, DAERHSFuncIDA daefun)
185  : m_t0 (t), m_y0 (y), m_yp0 (yp), m_havejac (false), m_havejacfcn (false),
186  m_havejacsparse (false), m_mem (nullptr), m_num (), m_ida_fcn (ida_fcn),
187  m_ida_jac (), m_dfdy (nullptr), m_dfdyp (nullptr), m_spdfdy (nullptr),
188  m_spdfdyp (nullptr), m_fcn (daefun), m_jacfcn (nullptr),
189  m_jacspfcn (nullptr), m_jacdcell (nullptr), m_jacspcell (nullptr),
190  m_sunJacMatrix (nullptr), m_sunLinearSolver (nullptr)
191  { }
192 
193  OCTAVE_DISABLE_COPY_MOVE (IDA)
194 
195  ~IDA ()
196  {
197  IDAFree (&m_mem);
198  SUNLinSolFree (m_sunLinearSolver);
199  SUNMatDestroy (m_sunJacMatrix);
200 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
201  SUNContext_Free (&m_sunContext);
202 # endif
203  }
204 
205  IDA&
206  set_jacobian (const octave_value& jac, DAEJacFuncDense j)
207  {
208  m_jacfcn = j;
209  m_ida_jac = jac;
210  m_havejac = true;
211  m_havejacfcn = true;
212  m_havejacsparse = false;
213 
214  return *this;
215  }
216 
217  IDA&
218  set_jacobian (const octave_value& jac, DAEJacFuncSparse j)
219  {
220  m_jacspfcn = j;
221  m_ida_jac = jac;
222  m_havejac = true;
223  m_havejacfcn = true;
224  m_havejacsparse = true;
225 
226  return *this;
227  }
228 
229  IDA&
230  set_jacobian (Matrix *dy, Matrix *dyp, DAEJacCellDense j)
231  {
232  m_jacdcell = j;
233  m_dfdy = dy;
234  m_dfdyp = dyp;
235  m_havejac = true;
236  m_havejacfcn = false;
237  m_havejacsparse = false;
238 
239  return *this;
240  }
241 
242  IDA&
243  set_jacobian (SparseMatrix *dy, SparseMatrix *dyp,
244  DAEJacCellSparse j)
245  {
246  m_jacspcell = j;
247  m_spdfdy = dy;
248  m_spdfdyp = dyp;
249  m_havejac = true;
250  m_havejacfcn = false;
251  m_havejacsparse = true;
252 
253  return *this;
254  }
255 
256  void set_userdata ();
257 
258  void initialize ();
259 
260  static ColumnVector NVecToCol (N_Vector& v, octave_f77_int_type n);
261 
262 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
263  N_Vector ColToNVec (const ColumnVector& data, octave_f77_int_type n);
264 # else
265  static N_Vector ColToNVec (const ColumnVector& data, octave_f77_int_type n);
266 # endif
267 
268  void
269  set_up (const ColumnVector& y);
270 
271  void
272  set_tolerance (ColumnVector& abstol, OCTAVE_SUNREALTYPE reltol);
273 
274  void
275  set_tolerance (OCTAVE_SUNREALTYPE abstol, OCTAVE_SUNREALTYPE reltol);
276 
277  static int
278  resfun (OCTAVE_SUNREALTYPE t, N_Vector yy, N_Vector yyp,
279  N_Vector rr, void *user_data);
280 
281  void
282  resfun_impl (OCTAVE_SUNREALTYPE t, N_Vector& yy,
283  N_Vector& yyp, N_Vector& rr);
284  static int
285  jacdense (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj, N_Vector yy,
286  N_Vector yyp, N_Vector, SUNMatrix JJ, void *user_data, N_Vector,
287  N_Vector, N_Vector)
288  {
289  IDA *self = static_cast <IDA *> (user_data);
290  self->jacdense_impl (t, cj, yy, yyp, JJ);
291  return 0;
292  }
293 
294  void
295  jacdense_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
296  N_Vector& yy, N_Vector& yyp, SUNMatrix& JJ);
297 
298 # if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
299  static int
300  jacsparse (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj, N_Vector yy,
301  N_Vector yyp, N_Vector, SUNMatrix Jac, void *user_data, N_Vector,
302  N_Vector, N_Vector)
303  {
304  IDA *self = static_cast <IDA *> (user_data);
305  self->jacsparse_impl (t, cj, yy, yyp, Jac);
306  return 0;
307  }
308 
309  void
310  jacsparse_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
311  N_Vector& yy, N_Vector& yyp, SUNMatrix& Jac);
312 # endif
313 
314  void set_maxstep (OCTAVE_SUNREALTYPE maxstep);
315 
316  void set_initialstep (OCTAVE_SUNREALTYPE initialstep);
317 
318  bool
319  interpolate (octave_idx_type& cont, Matrix& output, ColumnVector& tout,
320  int refine, OCTAVE_SUNREALTYPE tend, bool haveoutputfcn,
321  bool haveoutputsel, const octave_value& output_fcn,
322  ColumnVector& outputsel, bool haveeventfunction,
323  const octave_value& event_fcn, ColumnVector& te,
324  Matrix& ye, ColumnVector& ie, ColumnVector& oldval,
325  ColumnVector& oldisterminal, ColumnVector& olddir,
326  octave_idx_type& temp, ColumnVector& yold,
327  const octave_idx_type num_event_args);
328 
329  bool
330  outputfun (const octave_value& output_fcn, bool haveoutputsel,
331  const ColumnVector& output, OCTAVE_SUNREALTYPE tout,
332  OCTAVE_SUNREALTYPE tend, ColumnVector& outputsel,
333  const std::string& flag);
334 
335 
336  bool
337  event (const octave_value& event_fcn,
338  ColumnVector& te, Matrix& ye, ColumnVector& ie,
339  OCTAVE_SUNREALTYPE tsol, const ColumnVector& y,
340  const std::string& flag, const ColumnVector& yp, ColumnVector& oldval,
341  ColumnVector& oldisterminal, ColumnVector& olddir,
342  octave_idx_type cont, octave_idx_type& temp, OCTAVE_SUNREALTYPE told,
343  ColumnVector& yold,
344  const octave_idx_type num_event_args);
345 
346  void set_maxorder (int maxorder);
347 
349  integrate (const octave_idx_type numt, const ColumnVector& tt,
350  const ColumnVector& y0, const ColumnVector& yp0,
351  const int refine, bool haverefine, bool haveoutputfcn,
352  const octave_value& output_fcn, bool haveoutputsel,
353  ColumnVector& outputsel, bool haveeventfunction,
354  const octave_value& event_fcn,
355  const octave_idx_type num_event_args);
356 
357  void print_stat ();
358 
359 private:
360 
361  OCTAVE_SUNREALTYPE m_t0;
362  ColumnVector m_y0;
363  ColumnVector m_yp0;
364  bool m_havejac;
365  bool m_havejacfcn;
366  bool m_havejacsparse;
367  void *m_mem;
368  octave_f77_int_type m_num;
369  octave_value m_ida_fcn;
370  octave_value m_ida_jac;
371  Matrix *m_dfdy;
372  Matrix *m_dfdyp;
373  SparseMatrix *m_spdfdy;
374  SparseMatrix *m_spdfdyp;
375  DAERHSFuncIDA m_fcn;
376  DAEJacFuncDense m_jacfcn;
377  DAEJacFuncSparse m_jacspfcn;
378  DAEJacCellDense m_jacdcell;
379  DAEJacCellSparse m_jacspcell;
380 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
381  SUNContext m_sunContext;
382 # endif
383  SUNMatrix m_sunJacMatrix;
384  SUNLinearSolver m_sunLinearSolver;
385 };
386 
387 int
388 IDA::resfun (OCTAVE_SUNREALTYPE t, N_Vector yy, N_Vector yyp, N_Vector rr,
389  void *user_data)
390 {
391  IDA *self = static_cast <IDA *> (user_data);
392  self->resfun_impl (t, yy, yyp, rr);
393  return 0;
394 }
395 
396 void
397 IDA::resfun_impl (OCTAVE_SUNREALTYPE t, N_Vector& yy,
398  N_Vector& yyp, N_Vector& rr)
399 {
400  ColumnVector y = IDA::NVecToCol (yy, m_num);
401 
402  ColumnVector yp = IDA::NVecToCol (yyp, m_num);
403 
404  ColumnVector res = (*m_fcn) (y, yp, t, m_ida_fcn);
405 
406  OCTAVE_SUNREALTYPE *puntrr = nv_data_s (rr);
407 
408  for (octave_idx_type i = 0; i < m_num; i++)
409  puntrr[i] = res(i);
410 }
411 
412 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
413 # define OCTAVE_SUNCONTEXT , m_sunContext
414 # else
415 # define OCTAVE_SUNCONTEXT
416 # endif
417 
418 void
419 IDA::set_up (const ColumnVector& y)
420 {
421  N_Vector yy = ColToNVec (y, m_num);
422  octave::unwind_action act ([&yy] () { N_VDestroy_Serial (yy); });
423 
424  if (m_havejacsparse)
425  {
426 # if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
427 # if defined (HAVE_SUNSPARSEMATRIX_REALLOCATE)
428  // Initially allocate memory for 0 entries. We will reallocate when we
429  // get the Jacobian matrix from the user and know the actual number of
430  // entries.
431  m_sunJacMatrix = SUNSparseMatrix (m_num, m_num, 0, CSC_MAT
432  OCTAVE_SUNCONTEXT);
433 # else
434  octave_f77_int_type max_elems;
435  if (math::int_multiply_overflow (m_num, m_num, &max_elems))
436  error ("Unable to allocate memory for sparse Jacobian");
437 
438  m_sunJacMatrix = SUNSparseMatrix (m_num, m_num, max_elems, CSC_MAT
439  OCTAVE_SUNCONTEXT);
440 # endif
441  if (! m_sunJacMatrix)
442  error ("Unable to create sparse Jacobian for Sundials");
443 
444  m_sunLinearSolver = SUNLinSol_KLU (yy, m_sunJacMatrix
445  OCTAVE_SUNCONTEXT);
446  if (! m_sunLinearSolver)
447  error ("Unable to create KLU sparse solver");
448 
449  if (IDASetLinearSolver (m_mem, m_sunLinearSolver, m_sunJacMatrix))
450  error ("Unable to set sparse linear solver");
451 
452  IDASetJacFn (m_mem, IDA::jacsparse);
453 
454 # else
455  error ("SUNDIALS SUNLINSOL KLU was unavailable or disabled when "
456  "Octave was built");
457 
458 # endif
459 
460  }
461  else
462  {
463  m_sunJacMatrix = SUNDenseMatrix (m_num, m_num OCTAVE_SUNCONTEXT);
464  if (! m_sunJacMatrix)
465  error ("Unable to create dense Jacobian for Sundials");
466 
467  m_sunLinearSolver = SUNLinSol_Dense (yy, m_sunJacMatrix
468  OCTAVE_SUNCONTEXT);
469  if (! m_sunLinearSolver)
470  error ("Unable to create dense linear solver");
471 
472  if (IDASetLinearSolver (m_mem, m_sunLinearSolver, m_sunJacMatrix))
473  error ("Unable to set dense linear solver");
474 
475  if (m_havejac && IDASetJacFn (m_mem, IDA::jacdense) != 0)
476  error ("Unable to set dense Jacobian function");
477  }
478 }
479 
480 void
481 IDA::jacdense_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
482  N_Vector& yy, N_Vector& yyp, SUNMatrix& JJ)
483 
484 {
485  octave_f77_int_type Neq = NV_LENGTH_S (yy);
486 
487  ColumnVector y = NVecToCol (yy, Neq);
488 
489  ColumnVector yp = NVecToCol (yyp, Neq);
490 
491  Matrix jac;
492 
493  if (m_havejacfcn)
494  jac = (*m_jacfcn) (y, yp, t, cj, m_ida_jac);
495  else
496  jac = (*m_jacdcell) (m_dfdy, m_dfdyp, cj);
497 
498  octave_f77_int_type num_jac = to_f77_int (jac.numel ());
499  std::copy (jac.fortran_vec (),
500  jac.fortran_vec () + num_jac,
501  SUNDenseMatrix_Data (JJ));
502 }
503 
504 # if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
505 void
506 IDA::jacsparse_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj, N_Vector& yy,
507  N_Vector& yyp, SUNMatrix& Jac)
508 
509 {
510  ColumnVector y = NVecToCol (yy, m_num);
511 
512  ColumnVector yp = NVecToCol (yyp, m_num);
513 
514  SparseMatrix jac;
515 
516  if (m_havejacfcn)
517  jac = (*m_jacspfcn) (y, yp, t, cj, m_ida_jac);
518  else
519  jac = (*m_jacspcell) (m_spdfdy, m_spdfdyp, cj);
520 
521 # if defined (HAVE_SUNSPARSEMATRIX_REALLOCATE)
522  octave_f77_int_type nnz = to_f77_int (jac.nnz ());
523  if (nnz > SUNSparseMatrix_NNZ (Jac))
524  {
525  // Allocate memory for sparse Jacobian defined in user function.
526  // This will always be required at least once since we set the number
527  // of nonzero elements to zero initially.
528  if (SUNSparseMatrix_Reallocate (Jac, nnz))
529  error ("Unable to allocate sufficient memory for IDA sparse matrix");
530  }
531 # endif
532 
533  SUNMatZero_Sparse (Jac);
534  // We have to use "sunindextype *" here but still need to check that
535  // conversion of each element to "octave_f77_int_type" is save.
536  sunindextype *colptrs = SUNSparseMatrix_IndexPointers (Jac);
537  sunindextype *rowvals = SUNSparseMatrix_IndexValues (Jac);
538 
539  for (octave_f77_int_type i = 0; i < m_num + 1; i++)
540  colptrs[i] = to_f77_int (jac.cidx (i));
541 
542  double *d = SUNSparseMatrix_Data (Jac);
543  for (octave_f77_int_type i = 0; i < to_f77_int (jac.nnz ()); i++)
544  {
545  rowvals[i] = to_f77_int (jac.ridx (i));
546  d[i] = jac.data (i);
547  }
548 }
549 # endif
550 
552 IDA::NVecToCol (N_Vector& v, octave_f77_int_type n)
553 {
554  ColumnVector data (n);
555  OCTAVE_SUNREALTYPE *punt = nv_data_s (v);
556 
557  for (octave_f77_int_type i = 0; i < n; i++)
558  data(i) = punt[i];
559 
560  return data;
561 }
562 
563 N_Vector
564 IDA::ColToNVec (const ColumnVector& data, octave_f77_int_type n)
565 {
566  N_Vector v = N_VNew_Serial (n OCTAVE_SUNCONTEXT);
567 
568  OCTAVE_SUNREALTYPE *punt = nv_data_s (v);
569 
570  for (octave_f77_int_type i = 0; i < n; i++)
571  punt[i] = data(i);
572 
573  return v;
574 }
575 
576 void
577 IDA::set_userdata ()
578 {
579  void *userdata = this;
580 
581  if (IDASetUserData (m_mem, userdata) != 0)
582  error ("User data not set");
583 }
584 
585 void
586 IDA::initialize ()
587 {
588  m_num = to_f77_int (m_y0.numel ());
589 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
590 # if ! defined (SUN_COMM_NULL)
591 # define SUN_COMM_NULL nullptr
592 # endif
593  if (SUNContext_Create (SUN_COMM_NULL, &m_sunContext) < 0)
594  error ("__ode15__: unable to create context for SUNDIALS");
595  m_mem = IDACreate (m_sunContext);
596 # else
597  m_mem = IDACreate ();
598 # endif
599 
600  N_Vector yy = ColToNVec (m_y0, m_num);
601  N_Vector yyp = ColToNVec (m_yp0, m_num);
602  octave::unwind_action act ([&yy, &yyp] ()
603  { N_VDestroy_Serial (yy);
604  N_VDestroy_Serial (yyp);
605  });
606 
607  IDA::set_userdata ();
608 
609  if (IDAInit (m_mem, IDA::resfun, m_t0, yy, yyp) != 0)
610  error ("IDA not initialized");
611 }
612 
613 void
614 IDA::set_tolerance (ColumnVector& abstol, OCTAVE_SUNREALTYPE reltol)
615 {
616  N_Vector abs_tol = ColToNVec (abstol, m_num);
617  octave::unwind_action act ([&abs_tol] () { N_VDestroy_Serial (abs_tol); });
618 
619  if (IDASVtolerances (m_mem, reltol, abs_tol) != 0)
620  error ("IDA: Tolerance not set");
621 }
622 
623 void
624 IDA::set_tolerance (OCTAVE_SUNREALTYPE abstol, OCTAVE_SUNREALTYPE reltol)
625 {
626  if (IDASStolerances (m_mem, reltol, abstol) != 0)
627  error ("IDA: Tolerance not set");
628 }
629 
631 IDA::integrate (const octave_idx_type numt, const ColumnVector& tspan,
632  const ColumnVector& y, const ColumnVector& yp,
633  const int refine, bool haverefine, bool haveoutputfcn,
634  const octave_value& output_fcn, bool haveoutputsel,
635  ColumnVector& outputsel, bool haveeventfunction,
636  const octave_value& event_fcn,
637  const octave_idx_type num_event_args)
638 {
639  // Set up output
640  ColumnVector tout, yout (m_num), ypout (m_num), ysel (outputsel.numel ());
641  ColumnVector ie, te, oldval, oldisterminal, olddir;
642  Matrix output, ye;
643  octave_idx_type cont = 0, temp = 0;
644  bool status = false;
645  std::string string = "";
646  ColumnVector yold = y;
647 
648  OCTAVE_SUNREALTYPE tsol = tspan(0);
649  OCTAVE_SUNREALTYPE tend = tspan(numt-1);
650 
651  N_Vector yyp = ColToNVec (yp, m_num);
652  N_Vector yy = ColToNVec (y, m_num);
653  octave::unwind_action act ([&yyp, &yy] ()
654  { N_VDestroy_Serial (yyp);
655  N_VDestroy_Serial (yy);
656  });
657 
658  // Initialize OutputFcn
659  if (haveoutputfcn)
660  status = IDA::outputfun (output_fcn, haveoutputsel, y,
661  tsol, tend, outputsel, "init");
662 
663  // Initialize Events
664  if (haveeventfunction)
665  status = IDA::event (event_fcn, te, ye, ie, tsol, y,
666  "init", yp, oldval, oldisterminal,
667  olddir, cont, temp, tsol, yold, num_event_args);
668 
669  if (numt > 2)
670  {
671  // First output value
672  tout.resize (numt);
673  tout(0) = tsol;
674  output.resize (numt, m_num);
675 
676  for (octave_idx_type i = 0; i < m_num; i++)
677  output.elem (0, i) = y.elem (i);
678 
679  //Main loop
680  for (octave_idx_type j = 1; j < numt && status == 0; j++)
681  {
682  // IDANORMAL already interpolates tspan(j)
683 
684  if (IDASolve (m_mem, tspan (j), &tsol, yy, yyp, IDA_NORMAL) != 0)
685  error ("IDASolve failed");
686 
687  yout = NVecToCol (yy, m_num);
688  ypout = NVecToCol (yyp, m_num);
689  tout(j) = tsol;
690 
691  for (octave_idx_type i = 0; i < m_num; i++)
692  output.elem (j, i) = yout.elem (i);
693 
694  if (haveoutputfcn)
695  status = IDA::outputfun (output_fcn, haveoutputsel, yout, tsol,
696  tend, outputsel, string);
697 
698  if (haveeventfunction)
699  status = IDA::event (event_fcn, te, ye, ie, tout(j), yout,
700  string, ypout, oldval, oldisterminal,
701  olddir, j, temp, tout(j-1), yold,
702  num_event_args);
703 
704  // If integration is stopped, return only the reached steps
705  if (status == 1)
706  {
707  output.resize (j + 1, m_num);
708  tout.resize (j + 1);
709  }
710 
711  }
712  }
713  else // numel (tspan) == 2
714  {
715  // First output value
716  tout.resize (1);
717  tout(0) = tsol;
718  output.resize (1, m_num);
719 
720  for (octave_idx_type i = 0; i < m_num; i++)
721  output.elem (0, i) = y.elem (i);
722 
723  bool posdirection = (tend > tsol);
724 
725  //main loop
726  while (((posdirection == 1 && tsol < tend)
727  || (posdirection == 0 && tsol > tend))
728  && status == 0)
729  {
730  if (IDASolve (m_mem, tend, &tsol, yy, yyp, IDA_ONE_STEP) != 0)
731  error ("IDASolve failed");
732 
733  if (haverefine)
734  status = IDA::interpolate (cont, output, tout, refine, tend,
735  haveoutputfcn, haveoutputsel,
736  output_fcn, outputsel,
737  haveeventfunction, event_fcn, te,
738  ye, ie, oldval, oldisterminal,
739  olddir, temp, yold,
740  num_event_args);
741 
742  ypout = NVecToCol (yyp, m_num);
743  cont += 1;
744  output.resize (cont + 1, m_num); // This may be not efficient
745  tout.resize (cont + 1);
746  tout(cont) = tsol;
747  yout = NVecToCol (yy, m_num);
748 
749  for (octave_idx_type i = 0; i < m_num; i++)
750  output.elem (cont, i) = yout.elem (i);
751 
752  if (haveoutputfcn && ! haverefine && tout(cont) < tend)
753  status = IDA::outputfun (output_fcn, haveoutputsel, yout, tsol,
754  tend, outputsel, string);
755 
756  if (haveeventfunction && ! haverefine && tout(cont) < tend)
757  status = IDA::event (event_fcn, te, ye, ie, tout(cont), yout,
758  string, ypout, oldval, oldisterminal,
759  olddir, cont, temp, tout(cont-1), yold,
760  num_event_args);
761  }
762 
763  if (status == 0)
764  {
765  // Interpolate in tend
766  N_Vector dky = N_VNew_Serial (m_num OCTAVE_SUNCONTEXT);
767  octave::unwind_action act2 ([&dky] () { N_VDestroy_Serial (dky); });
768 
769  if (IDAGetDky (m_mem, tend, 0, dky) != 0)
770  error ("IDA failed to interpolate y");
771 
772  tout(cont) = tend;
773  yout = NVecToCol (dky, m_num);
774 
775  for (octave_idx_type i = 0; i < m_num; i++)
776  output.elem (cont, i) = yout.elem (i);
777 
778  // Plot final value
779  if (haveoutputfcn)
780  {
781  status = IDA::outputfun (output_fcn, haveoutputsel, yout,
782  tend, tend, outputsel, string);
783 
784  // Events during last step
785  if (haveeventfunction)
786  status = IDA::event (event_fcn, te, ye, ie, tend, yout,
787  string, ypout, oldval, oldisterminal,
788  olddir, cont, temp, tout(cont-1),
789  yold, num_event_args);
790  }
791  }
792 
793  // Cleanup plotter
794  status = IDA::outputfun (output_fcn, haveoutputsel, yout, tend, tend,
795  outputsel, "done");
796  }
797 
798  // Index of Events (ie) variable must use 1-based indexing
799  return ovl (tout, output, te, ye, ie + 1.0);
800 }
801 
802 bool
803 IDA::event (const octave_value& event_fcn,
804  ColumnVector& te, Matrix& ye, ColumnVector& ie,
805  OCTAVE_SUNREALTYPE tsol, const ColumnVector& y,
806  const std::string& flag, const ColumnVector& yp,
807  ColumnVector& oldval, ColumnVector& oldisterminal,
808  ColumnVector& olddir, octave_idx_type cont, octave_idx_type& temp,
809  OCTAVE_SUNREALTYPE told, ColumnVector& yold,
810  const octave_idx_type num_event_args)
811 {
812  bool status = false;
813 
814  octave_value_list args;
815  if (num_event_args == 2)
816  args = ovl (tsol, y);
817  else
818  args = ovl (tsol, y, yp);
819 
820  // cont is the number of steps reached by the solver
821  // temp is the number of events registered
822 
823  interpreter& interp = __get_interpreter__ ();
824 
825  if (flag == "init")
826  {
827  octave_value_list output = interp.feval (event_fcn, args, 3);
828  oldval = output(0).vector_value ();
829  oldisterminal = output(1).vector_value ();
830  olddir = output(2).vector_value ();
831  }
832  else if (flag == "")
833  {
834  ColumnVector index (0);
835  octave_value_list output = interp.feval (event_fcn, args, 3);
836  ColumnVector val = output(0).vector_value ();
837  ColumnVector isterminal = output(1).vector_value ();
838  ColumnVector dir = output(2).vector_value ();
839 
840  // Get the index of the changed values
841  for (octave_idx_type i = 0; i < val.numel (); i++)
842  {
843  // Check for sign change and whether a rising / falling edge
844  // either passes through zero or detaches from zero (bug #59063)
845  if ((dir(i) != -1
846  && ((val(i) >= 0 && oldval(i) < 0)
847  || (val(i) > 0 && oldval(i) <= 0))) // increasing
848  || (dir(i) != 1
849  && ((val(i) <= 0 && oldval(i) > 0)
850  || (val(i) < 0 && oldval(i) >= 0)))) // decreasing
851  {
852  index.resize (index.numel () + 1);
853  index (index.numel () - 1) = i;
854  }
855  }
856 
857  if (cont == 1 && index.numel () > 0) // Events in first step
858  {
859  temp = 1; // register only the first event
860  te.resize (1);
861  ye.resize (1, m_num);
862  ie.resize (1);
863 
864  // Linear interpolation
865  ie(0) = index(0);
866  te(0) = tsol - val (index(0)) * (tsol - told)
867  / (val (index(0)) - oldval (index(0)));
868 
869  ColumnVector ytemp
870  = y - ((tsol - te(0)) * (y - yold) / (tsol - told));
871 
872  for (octave_idx_type i = 0; i < m_num; i++)
873  ye.elem (0, i) = ytemp.elem (i);
874 
875  }
876  else if (index.numel () > 0)
877  // Not first step: register all events and test
878  // if stop integration or not
879  {
880  te.resize (temp + index.numel ());
881  ye.resize (temp + index.numel (), m_num);
882  ie.resize (temp + index.numel ());
883 
884  for (octave_idx_type i = 0; i < index.numel (); i++)
885  {
886 
887  if (isterminal (index(i)) == 1)
888  status = 1; // Stop integration
889 
890  // Linear interpolation
891  ie(temp+i) = index(i);
892  te(temp+i) = tsol - val(index(i)) * (tsol - told)
893  / (val(index(i)) - oldval(index(i)));
894 
895  ColumnVector ytemp
896  = y - (tsol - te (temp + i)) * (y - yold) / (tsol - told);
897 
898  for (octave_idx_type j = 0; j < m_num; j++)
899  ye.elem (temp + i, j) = ytemp.elem (j);
900 
901  }
902 
903  temp += index.numel ();
904  }
905 
906  // Update variables
907  yold = y;
908  told = tsol;
909  olddir = dir;
910  oldval = val;
911  oldisterminal = isterminal;
912  }
913 
914  return status;
915 }
916 
917 bool
918 IDA::interpolate (octave_idx_type& cont, Matrix& output, ColumnVector& tout,
919  int refine, OCTAVE_SUNREALTYPE tend, bool haveoutputfcn,
920  bool haveoutputsel, const octave_value& output_fcn,
921  ColumnVector& outputsel, bool haveeventfunction,
922  const octave_value& event_fcn, ColumnVector& te,
923  Matrix& ye, ColumnVector& ie, ColumnVector& oldval,
924  ColumnVector& oldisterminal, ColumnVector& olddir,
925  octave_idx_type& temp, ColumnVector& yold,
926  const octave_idx_type num_event_args)
927 {
928  OCTAVE_SUNREALTYPE h = 0, tcur = 0;
929  bool status = false;
930 
931  N_Vector dky = N_VNew_Serial (m_num OCTAVE_SUNCONTEXT);
932  N_Vector dkyp = N_VNew_Serial (m_num OCTAVE_SUNCONTEXT);
933  octave::unwind_action act ([&dky, &dkyp] ()
934  { N_VDestroy_Serial (dky);
935  N_VDestroy_Serial (dkyp);
936  });
937 
938  ColumnVector yout (m_num);
939  ColumnVector ypout (m_num);
940  std::string string = "";
941 
942  if (IDAGetLastStep (m_mem, &h) != 0)
943  error ("IDA failed to return last step");
944 
945  if (IDAGetCurrentTime (m_mem, &tcur) != 0)
946  error ("IDA failed to return the current time");
947 
948  OCTAVE_SUNREALTYPE tin = tcur - h;
949 
950  OCTAVE_SUNREALTYPE step = h / refine;
951 
952  for (octave_idx_type i = 1;
953  i < refine && tin + step * i < tend && status == 0;
954  i++)
955  {
956  if (IDAGetDky (m_mem, tin + step*i, 0, dky) != 0)
957  error ("IDA failed to interpolate y");
958 
959  if (IDAGetDky (m_mem, tin + step*i, 1, dkyp) != 0)
960  error ("IDA failed to interpolate yp");
961 
962  cont += 1;
963  output.resize (cont + 1, m_num);
964  tout.resize (cont + 1);
965 
966  tout(cont) = tin + step * i;
967  yout = NVecToCol (dky, m_num);
968  ypout = NVecToCol (dkyp, m_num);
969 
970  for (octave_idx_type j = 0; j < m_num; j++)
971  output.elem (cont, j) = yout.elem (j);
972 
973  if (haveoutputfcn)
974  status = IDA::outputfun (output_fcn, haveoutputsel, yout,
975  tout(cont), tend, outputsel, "");
976 
977  if (haveeventfunction)
978  status = IDA::event (event_fcn, te, ye, ie, tout(cont),
979  yout, string, ypout, oldval,
980  oldisterminal, olddir, cont, temp,
981  tout(cont-1), yold, num_event_args);
982  }
983 
984  return status;
985 }
986 
987 bool
988 IDA::outputfun (const octave_value& output_fcn, bool haveoutputsel,
989  const ColumnVector& yout, OCTAVE_SUNREALTYPE tsol,
990  OCTAVE_SUNREALTYPE tend, ColumnVector& outputsel,
991  const std::string& flag)
992 {
993  bool status = false;
994 
995  octave_value_list output;
996  output(2) = flag;
997 
998  ColumnVector ysel (outputsel.numel ());
999  if (haveoutputsel)
1000  {
1001  for (octave_idx_type i = 0; i < outputsel.numel (); i++)
1002  ysel(i) = yout(outputsel(i));
1003 
1004  output(1) = ysel;
1005  }
1006  else
1007  output(1) = yout;
1008 
1009  interpreter& interp = __get_interpreter__ ();
1010 
1011  if (flag == "init")
1012  {
1013  ColumnVector toutput (2);
1014  toutput(0) = tsol;
1015  toutput(1) = tend;
1016  output(0) = toutput;
1017 
1018  interp.feval (output_fcn, output, 0);
1019  }
1020  else if (flag == "")
1021  {
1022  output(0) = tsol;
1023  octave_value_list val = interp.feval (output_fcn, output, 1);
1024  status = val(0).bool_value ();
1025  }
1026  else
1027  {
1028  // Cleanup plotter
1029  output(0) = tend;
1030  interp.feval (output_fcn, output, 0);
1031  }
1032 
1033  return status;
1034 }
1035 
1036 void
1037 IDA::set_maxstep (OCTAVE_SUNREALTYPE maxstep)
1038 {
1039  if (IDASetMaxStep (m_mem, maxstep) != 0)
1040  error ("IDA: Max Step not set");
1041 }
1042 
1043 void
1044 IDA::set_initialstep (OCTAVE_SUNREALTYPE initialstep)
1045 {
1046  if (IDASetInitStep (m_mem, initialstep) != 0)
1047  error ("IDA: Initial Step not set");
1048 }
1049 
1050 void
1051 IDA::set_maxorder (int maxorder)
1052 {
1053  if (IDASetMaxOrd (m_mem, maxorder) != 0)
1054  error ("IDA: Max Order not set");
1055 }
1056 
1057 void
1058 IDA::print_stat ()
1059 {
1060  long int nsteps = 0, netfails = 0, nrevals = 0;
1061 
1062  if (IDAGetNumSteps (m_mem, &nsteps) != 0)
1063  error ("IDA failed to return the number of internal steps");
1064 
1065  if (IDAGetNumErrTestFails (m_mem, &netfails) != 0)
1066  error ("IDA failed to return the number of internal errors");
1067 
1068  if (IDAGetNumResEvals (m_mem, &nrevals) != 0)
1069  error ("IDA failed to return the number of residual evaluations");
1070 
1071  octave_stdout << nsteps << " successful steps\n";
1072  octave_stdout << netfails << " failed attempts\n";
1073  octave_stdout << nrevals << " function evaluations\n";
1074  // octave_stdout << " partial derivatives\n";
1075  // octave_stdout << " LU decompositions\n";
1076  // octave_stdout << " solutions of linear systems\n";
1077 }
1078 
1079 static ColumnVector
1080 ida_user_function (const ColumnVector& x, const ColumnVector& xdot,
1081  double t, const octave_value& ida_fc)
1082 {
1083  octave_value_list tmp;
1084 
1085  try
1086  {
1087  interpreter& interp = __get_interpreter__ ();
1088 
1089  tmp = interp.feval (ida_fc, ovl (t, x, xdot), 1);
1090  }
1091  catch (execution_exception& ee)
1092  {
1093  err_user_supplied_eval (ee, "__ode15__");
1094  }
1095 
1096  return tmp(0).vector_value ();
1097 }
1098 
1099 static Matrix
1100 ida_dense_jac (const ColumnVector& x, const ColumnVector& xdot,
1101  double t, double cj, const octave_value& ida_jc)
1102 {
1103  octave_value_list tmp;
1104 
1105  try
1106  {
1107  interpreter& interp = __get_interpreter__ ();
1108 
1109  tmp = interp.feval (ida_jc, ovl (t, x, xdot), 2);
1110  }
1111  catch (execution_exception& ee)
1112  {
1113  err_user_supplied_eval (ee, "__ode15__");
1114  }
1115 
1116  return tmp(0).matrix_value () + cj * tmp(1).matrix_value ();
1117 }
1118 
1119 static SparseMatrix
1120 ida_sparse_jac (const ColumnVector& x, const ColumnVector& xdot,
1121  double t, double cj, const octave_value& ida_jc)
1122 {
1123  octave_value_list tmp;
1124 
1125  try
1126  {
1127  interpreter& interp = __get_interpreter__ ();
1128 
1129  tmp = interp.feval (ida_jc, ovl (t, x, xdot), 2);
1130  }
1131  catch (execution_exception& ee)
1132  {
1133  err_user_supplied_eval (ee, "__ode15__");
1134  }
1135 
1136  return tmp(0).sparse_matrix_value () + cj * tmp(1).sparse_matrix_value ();
1137 }
1138 
1139 static Matrix
1140 ida_dense_cell_jac (Matrix *dfdy, Matrix *dfdyp, double cj)
1141 {
1142  return (*dfdy) + cj * (*dfdyp);
1143 }
1144 
1145 static SparseMatrix
1146 ida_sparse_cell_jac (SparseMatrix *spdfdy, SparseMatrix *spdfdyp,
1147  double cj)
1148 {
1149  return (*spdfdy) + cj * (*spdfdyp);
1150 }
1151 
1152 static octave_value_list
1153 do_ode15 (const octave_value& ida_fcn,
1154  const ColumnVector& tspan,
1155  const octave_idx_type numt,
1156  const OCTAVE_SUNREALTYPE t0,
1157  const ColumnVector& y0,
1158  const ColumnVector& yp0,
1159  const octave_scalar_map& options,
1160  const octave_idx_type num_event_args)
1161 {
1162  octave_value_list retval;
1163 
1164  // Create object
1165  IDA dae (t0, y0, yp0, ida_fcn, ida_user_function);
1166 
1167  // Set Jacobian
1168  bool havejac = options.getfield ("havejac").bool_value ();
1169 
1170  bool havejacsparse = options.getfield ("havejacsparse").bool_value ();
1171 
1172  bool havejacfcn = options.getfield ("havejacfcn").bool_value ();
1173 
1174  Matrix ida_dfdy, ida_dfdyp;
1175  SparseMatrix ida_spdfdy, ida_spdfdyp;
1176 
1177  if (havejac)
1178  {
1179  if (havejacfcn)
1180  {
1181  octave_value ida_jac = options.getfield ("Jacobian");
1182 
1183  if (havejacsparse)
1184  dae.set_jacobian (ida_jac, ida_sparse_jac);
1185  else
1186  dae.set_jacobian (ida_jac, ida_dense_jac);
1187  }
1188  else
1189  {
1190  Cell jaccell = options.getfield ("Jacobian").cell_value ();
1191 
1192  if (havejacsparse)
1193  {
1194  ida_spdfdy = jaccell(0).sparse_matrix_value ();
1195  ida_spdfdyp = jaccell(1).sparse_matrix_value ();
1196 
1197  dae.set_jacobian (&ida_spdfdy, &ida_spdfdyp,
1198  ida_sparse_cell_jac);
1199  }
1200  else
1201  {
1202  ida_dfdy = jaccell(0).matrix_value ();
1203  ida_dfdyp = jaccell(1).matrix_value ();
1204 
1205  dae.set_jacobian (&ida_dfdy, &ida_dfdyp, ida_dense_cell_jac);
1206  }
1207  }
1208  }
1209 
1210  // Initialize IDA
1211  dae.initialize ();
1212 
1213  // Set tolerances
1214  OCTAVE_SUNREALTYPE rel_tol = options.getfield ("RelTol").double_value ();
1215 
1216  bool haveabstolvec = options.getfield ("haveabstolvec").bool_value ();
1217 
1218  if (haveabstolvec)
1219  {
1220  ColumnVector abs_tol = options.getfield ("AbsTol").vector_value ();
1221 
1222  dae.set_tolerance (abs_tol, rel_tol);
1223  }
1224  else
1225  {
1226  OCTAVE_SUNREALTYPE abs_tol = options.getfield ("AbsTol").double_value ();
1227 
1228  dae.set_tolerance (abs_tol, rel_tol);
1229  }
1230 
1231  //Set max step
1232  OCTAVE_SUNREALTYPE maxstep = options.getfield ("MaxStep").double_value ();
1233 
1234  dae.set_maxstep (maxstep);
1235 
1236  //Set initial step
1237  if (! options.getfield ("InitialStep").isempty ())
1238  {
1239  OCTAVE_SUNREALTYPE initialstep = options.getfield ("InitialStep").double_value ();
1240 
1241  dae.set_initialstep (initialstep);
1242  }
1243 
1244  //Set max order FIXME: it doesn't work
1245  int maxorder = options.getfield ("MaxOrder").int_value ();
1246 
1247  dae.set_maxorder (maxorder);
1248 
1249  //Set Refine
1250  const int refine = options.getfield ("Refine").int_value ();
1251 
1252  bool haverefine = (refine > 1);
1253 
1254  octave_value output_fcn;
1255  ColumnVector outputsel;
1256 
1257  // OutputFcn
1258  bool haveoutputfunction
1259  = options.getfield ("haveoutputfunction").bool_value ();
1260 
1261  if (haveoutputfunction)
1262  output_fcn = options.getfield ("OutputFcn");
1263 
1264  // OutputSel
1265  bool haveoutputsel = options.getfield ("haveoutputselection").bool_value ();
1266 
1267  if (haveoutputsel)
1268  outputsel = options.getfield ("OutputSel").vector_value ();
1269 
1270  octave_value event_fcn;
1271 
1272  // Events
1273  bool haveeventfunction
1274  = options.getfield ("haveeventfunction").bool_value ();
1275 
1276  if (haveeventfunction)
1277  event_fcn = options.getfield ("Events");
1278 
1279  // Set up linear solver
1280  dae.set_up (y0);
1281 
1282  // Integrate
1283  retval = dae.integrate (numt, tspan, y0, yp0, refine,
1284  haverefine, haveoutputfunction,
1285  output_fcn, haveoutputsel, outputsel,
1286  haveeventfunction, event_fcn, num_event_args);
1287 
1288  // Statistics
1289  bool havestats = options.getfield ("havestats").bool_value ();
1290 
1291  if (havestats)
1292  dae.print_stat ();
1293 
1294  return retval;
1295 }
1296 
1297 #endif
1298 
1299 DEFUN_DLD (__ode15__, args, ,
1300  doc: /* -*- texinfo -*-
1301 @deftypefn {} {@var{t}, @var{y} =} __ode15__ (@var{fcn}, @var{tspan}, @var{y0}, @var{yp0}, @var{ode_opt}, @var{num_event_args})
1302 Undocumented internal function.
1303 @end deftypefn */)
1304 {
1305 
1306 #if defined (HAVE_SUNDIALS)
1307 
1308  // Check number of parameters
1309  if (args.length () != 6)
1310  print_usage ();
1311 
1312  // Check ODE function
1313  octave_value ida_fcn = args(0);
1314 
1315  if (! ida_fcn.is_function_handle ())
1316  error ("__ode15__: FCN must be a function handle");
1317 
1318  // Check input tspan
1319  ColumnVector tspan
1320  = args(1).xvector_value ("__ode15__: TRANGE must be a vector of numbers");
1321 
1322  octave_idx_type numt = tspan.numel ();
1323 
1324  OCTAVE_SUNREALTYPE t0 = tspan(0);
1325 
1326  if (numt < 2)
1327  error ("__ode15__: TRANGE must contain at least 2 elements");
1328  else if (tspan.issorted () == UNSORTED || tspan(0) == tspan(numt - 1))
1329  error ("__ode15__: TRANGE must be strictly monotonic");
1330 
1331  // input y0 and yp0
1332  ColumnVector y0
1333  = args(2).xvector_value ("__ode15__: initial state Y0 must be a vector");
1334 
1335  ColumnVector yp0
1336  = args(3).xvector_value ("__ode15__: initial state YP0 must be a vector");
1337 
1338 
1339  if (y0.numel () != yp0.numel ())
1340  error ("__ode15__: initial state Y0 and YP0 must have the same length");
1341  else if (y0.numel () < 1)
1342  error ("__ode15__: initial state YP0 must be a vector or a scalar");
1343 
1344 
1345  if (! args(4).isstruct ())
1346  error ("__ode15__: ODE_OPT argument must be a structure");
1347 
1348  octave_scalar_map options
1349  = args(4).xscalar_map_value ("__ode15__: ODE_OPT argument must be a scalar structure");
1350 
1351  // Provided number of arguments in the ode callback function
1352  octave_idx_type num_event_args
1353  = args(5).xidx_type_value ("__ode15__: NUM_EVENT_ARGS must be an integer");
1354 
1355  if (num_event_args != 2 && num_event_args != 3)
1356  error ("__ode15__: number of input arguments in event callback must be 2 or 3");
1357 
1358  return do_ode15 (ida_fcn, tspan, numt, t0, y0, yp0, options, num_event_args);
1359 
1360 #else
1361 
1362  octave_unused_parameter (args);
1363 
1364  err_disabled_feature ("__ode15__", "sundials_ida, sundials_nvecserial");
1365 
1366 #endif
1367 }
1368 
1369 /*
1370 ## No test needed for internal helper function.
1371 %!assert (1)
1372 */
1373 
1374 OCTAVE_END_NAMESPACE(octave)
T & elem(octave_idx_type n)
Size of the specified dimension.
Definition: Array.h:562
sortmode issorted(sortmode mode=UNSORTED) const
Ordering is auto-detected or can be specified.
Definition: Array-base.cc:2047
octave_idx_type numel() const
Number of elements in the array.
Definition: Array.h:414
Definition: Cell.h:43
void resize(octave_idx_type n, const double &rfv=0)
Definition: dColVector.h:114
Definition: dMatrix.h:42
void resize(octave_idx_type nr, octave_idx_type nc, double rfv=0)
Definition: dMatrix.h:158
octave_value_list feval(const char *name, const octave_value_list &args=octave_value_list(), int nargout=0)
Evaluate an Octave function (built-in or interpreted) and return the list of result values.
octave_value getfield(const std::string &key) const
Definition: oct-map.cc:183
void resize(octave_idx_type n, const octave_value &rfv=octave_value())
Definition: ovl.h:117
octave_idx_type length() const
Definition: ovl.h:113
bool is_function_handle() const
Definition: ov.h:768
bool bool_value(bool warn=false) const
Definition: ov.h:885
Cell cell_value() const
int int_value(bool req_int=false, bool frc_str_conv=false) const
Definition: ov.h:812
bool isempty() const
Definition: ov.h:601
Array< double > vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
double double_value(bool frc_str_conv=false) const
Definition: ov.h:841
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
#define DEFUN_DLD(name, args_name, nargout_name, doc)
Macro to define an at run time dynamically loadable builtin function.
Definition: defun-dld.h:61
void print_usage(void)
Definition: defun-int.h:72
void() error(const char *fmt,...)
Definition: error.cc:988
void err_disabled_feature(const std::string &fcn, const std::string &feature, const std::string &pkg)
Definition: errwarn.cc:53
void err_user_supplied_eval(const char *name)
Definition: errwarn.cc:152
interpreter & __get_interpreter__()
F77_RET_T const F77_INT F77_CMPLX * A
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
F77_RET_T const F77_DBLE * x
bool int_multiply_overflow(int a, int b, int *r)
Definition: lo-utils.cc:499
octave_idx_type n
Definition: mx-inlines.cc:761
@ UNSORTED
Definition: oct-sort.h:97
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
Definition: ovl.h:219
#define octave_stdout
Definition: pager.h:309