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