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