26#if defined (HAVE_CONFIG_H)
47#if defined (HAVE_SUNDIALS)
49# if defined (HAVE_NVECTOR_NVECTOR_SERIAL_H)
50# include <nvector/nvector_serial.h>
53# if defined (HAVE_IDA_IDA_H)
55# elif defined (HAVE_IDA_H)
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>
64# if defined (HAVE_SUNLINSOL_SUNLINSOL_DENSE_H)
65# include <sunlinsol/sunlinsol_dense.h>
68# if defined (HAVE_SUNLINSOL_SUNLINSOL_KLU_H)
69# if defined (HAVE_KLU_H)
72# if defined (HAVE_KLU_KLU_H)
75# if defined (HAVE_SUITESPARSE_KLU_H)
76# include <suitesparse/klu.h>
78# include <sunlinsol/sunlinsol_klu.h>
85#if defined (HAVE_SUNDIALS)
87# if ! defined (HAVE_IDASETJACFN) && defined (HAVE_IDADLSSETJACFN)
89IDASetJacFn (
void *ida_mem, IDADlsJacFn jac)
91 return IDADlsSetJacFn (ida_mem, jac);
95# if ! defined (HAVE_IDASETLINEARSOLVER) && defined (HAVE_IDADLSSETLINEARSOLVER)
97IDASetLinearSolver (
void *ida_mem, SUNLinearSolver LS, SUNMatrix
A)
99 return IDADlsSetLinearSolver (ida_mem, LS,
A);
103# if ! defined (HAVE_SUNLINSOL_DENSE) && defined (HAVE_SUNDENSELINEARSOLVER)
104static inline SUNLinearSolver
105SUNLinSol_Dense (N_Vector y, SUNMatrix
A)
107 return SUNDenseLinearSolver (y,
A);
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)
116 return SUNKLU (y,
A);
121static inline OCTAVE_SUNREALTYPE *
122nv_data_s (N_Vector& v)
124# if defined (HAVE_PRAGMA_GCC_DIAGNOSTIC)
128# pragma GCC diagnostic push
129# pragma GCC diagnostic ignored "-Wold-style-cast"
132 return NV_DATA_S (v);
134# if defined (HAVE_PRAGMA_GCC_DIAGNOSTIC)
136# pragma GCC diagnostic pop
147 OCTAVE_SUNREALTYPE t,
158 OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
163 OCTAVE_SUNREALTYPE cj);
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)
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)
190 OCTAVE_DISABLE_COPY_MOVE (IDA)
195 SUNLinSolFree (m_sunLinearSolver);
196 SUNMatDestroy (m_sunJacMatrix);
197# if defined (HAVE_SUNDIALS_SUNCONTEXT)
198 SUNContext_Free (&m_sunContext);
203 set_jacobian (
const octave_value& jac, DAEJacFuncDense j)
209 m_havejacsparse =
false;
215 set_jacobian (
const octave_value& jac, DAEJacFuncSparse j)
221 m_havejacsparse =
true;
227 set_jacobian (
Matrix *dy,
Matrix *dyp, DAEJacCellDense j)
233 m_havejacfcn =
false;
234 m_havejacsparse =
false;
247 m_havejacfcn =
false;
248 m_havejacsparse =
true;
253 void set_userdata ();
257 static ColumnVector NVecToCol (N_Vector& v, octave_f77_int_type n);
259# if defined (HAVE_SUNDIALS_SUNCONTEXT)
260 N_Vector ColToNVec (
const ColumnVector& data, octave_f77_int_type n);
262 static N_Vector ColToNVec (
const ColumnVector& data, octave_f77_int_type n);
269 set_tolerance (
ColumnVector& abstol, OCTAVE_SUNREALTYPE reltol);
272 set_tolerance (OCTAVE_SUNREALTYPE abstol, OCTAVE_SUNREALTYPE reltol);
275 resfun (OCTAVE_SUNREALTYPE t, N_Vector yy, N_Vector yyp,
276 N_Vector rr,
void *user_data);
279 resfun_impl (OCTAVE_SUNREALTYPE t, N_Vector& yy,
280 N_Vector& yyp, N_Vector& rr);
282 jacdense (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj, N_Vector yy,
283 N_Vector yyp, N_Vector, SUNMatrix JJ,
void *user_data, N_Vector,
286 IDA *self =
static_cast<IDA *
> (user_data);
287 self->jacdense_impl (t, cj, yy, yyp, JJ);
292 jacdense_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
293 N_Vector& yy, N_Vector& yyp, SUNMatrix& JJ);
295# if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
297 jacsparse (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj, N_Vector yy,
298 N_Vector yyp, N_Vector, SUNMatrix Jac,
void *user_data, N_Vector,
301 IDA *self =
static_cast<IDA *
> (user_data);
302 self->jacsparse_impl (t, cj, yy, yyp, Jac);
307 jacsparse_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
308 N_Vector& yy, N_Vector& yyp, SUNMatrix& Jac);
311 void set_maxstep (OCTAVE_SUNREALTYPE maxstep);
313 void set_initialstep (OCTAVE_SUNREALTYPE initialstep);
317 int refine, OCTAVE_SUNREALTYPE tend,
bool haveoutputfcn,
327 outputfun (
const octave_value& output_fcn,
bool haveoutputsel,
330 const std::string& flag);
343 void set_maxorder (
int maxorder);
348 const int refine,
bool haverefine,
bool haveoutputfcn,
358 OCTAVE_SUNREALTYPE m_t0;
363 bool m_havejacsparse;
365 octave_f77_int_type m_num;
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;
380 SUNMatrix m_sunJacMatrix;
381 SUNLinearSolver m_sunLinearSolver;
385IDA::resfun (OCTAVE_SUNREALTYPE t, N_Vector yy, N_Vector yyp, N_Vector rr,
388 IDA *self =
static_cast<IDA *
> (user_data);
389 self->resfun_impl (t, yy, yyp, rr);
394IDA::resfun_impl (OCTAVE_SUNREALTYPE t, N_Vector& yy,
395 N_Vector& yyp, N_Vector& rr)
403 OCTAVE_SUNREALTYPE *puntrr = nv_data_s (rr);
409# if defined (HAVE_SUNDIALS_SUNCONTEXT)
410# define OCTAVE_SUNCONTEXT , m_sunContext
412# define OCTAVE_SUNCONTEXT
418 N_Vector yy = ColToNVec (y, m_num);
419 octave::unwind_action act ([&yy] () { N_VDestroy_Serial (yy); });
423# if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
424# if defined (HAVE_SUNSPARSEMATRIX_REALLOCATE)
428 m_sunJacMatrix = SUNSparseMatrix (m_num, m_num, 0, CSC_MAT
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");
435 m_sunJacMatrix = SUNSparseMatrix (m_num, m_num, max_elems, CSC_MAT
438 if (! m_sunJacMatrix)
439 error (
"Unable to create sparse Jacobian for Sundials");
441 m_sunLinearSolver = SUNLinSol_KLU (yy, m_sunJacMatrix
443 if (! m_sunLinearSolver)
444 error (
"Unable to create KLU sparse solver");
446 if (IDASetLinearSolver (m_mem, m_sunLinearSolver, m_sunJacMatrix))
447 error (
"Unable to set sparse linear solver");
449 IDASetJacFn (m_mem, IDA::jacsparse);
452 error (
"SUNDIALS SUNLINSOL KLU was unavailable or disabled when "
460 m_sunJacMatrix = SUNDenseMatrix (m_num, m_num OCTAVE_SUNCONTEXT);
461 if (! m_sunJacMatrix)
462 error (
"Unable to create dense Jacobian for Sundials");
464 m_sunLinearSolver = SUNLinSol_Dense (yy, m_sunJacMatrix
466 if (! m_sunLinearSolver)
467 error (
"Unable to create dense linear solver");
469 if (IDASetLinearSolver (m_mem, m_sunLinearSolver, m_sunJacMatrix))
470 error (
"Unable to set dense linear solver");
472 if (m_havejac && IDASetJacFn (m_mem, IDA::jacdense) != 0)
473 error (
"Unable to set dense Jacobian function");
478IDA::jacdense_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
479 N_Vector& yy, N_Vector& yyp, SUNMatrix& JJ)
482 octave_f77_int_type Neq = NV_LENGTH_S (yy);
491 jac = (*m_jacfcn) (y, yp, t, cj, m_ida_jac);
493 jac = (*m_jacdcell) (m_dfdy, m_dfdyp, cj);
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));
501# if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
503IDA::jacsparse_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj, N_Vector& yy,
504 N_Vector& yyp, SUNMatrix& Jac)
514 jac = (*m_jacspfcn) (y, yp, t, cj, m_ida_jac);
516 jac = (*m_jacspcell) (m_spdfdy, m_spdfdyp, cj);
518# if defined (HAVE_SUNSPARSEMATRIX_REALLOCATE)
519 octave_f77_int_type nnz = to_f77_int (jac.nnz ());
520 if (nnz > SUNSparseMatrix_NNZ (Jac))
525 if (SUNSparseMatrix_Reallocate (Jac, nnz))
526 error (
"Unable to allocate sufficient memory for IDA sparse matrix");
530 SUNMatZero_Sparse (Jac);
533 sunindextype *colptrs = SUNSparseMatrix_IndexPointers (Jac);
534 sunindextype *rowvals = SUNSparseMatrix_IndexValues (Jac);
536 for (octave_f77_int_type i = 0; i < m_num + 1; i++)
537 colptrs[i] = to_f77_int (jac.cidx (i));
539 double *
d = SUNSparseMatrix_Data (Jac);
540 for (octave_f77_int_type i = 0; i < to_f77_int (jac.nnz ()); i++)
542 rowvals[i] = to_f77_int (jac.ridx (i));
549IDA::NVecToCol (N_Vector& v, octave_f77_int_type n)
552 OCTAVE_SUNREALTYPE *punt = nv_data_s (v);
554 for (octave_f77_int_type i = 0; i < n; i++)
561IDA::ColToNVec (
const ColumnVector& data, octave_f77_int_type n)
563 N_Vector v = N_VNew_Serial (n OCTAVE_SUNCONTEXT);
565 OCTAVE_SUNREALTYPE *punt = nv_data_s (v);
567 for (octave_f77_int_type i = 0; i < n; i++)
576 void *userdata =
this;
578 if (IDASetUserData (m_mem, userdata) != 0)
579 error (
"User data not set");
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
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);
594 m_mem = IDACreate ();
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);
604 IDA::set_userdata ();
606 if (IDAInit (m_mem, IDA::resfun, m_t0, yy, yyp) != 0)
607 error (
"IDA not initialized");
611IDA::set_tolerance (
ColumnVector& abstol, OCTAVE_SUNREALTYPE reltol)
613 N_Vector abs_tol = ColToNVec (abstol, m_num);
614 octave::unwind_action act ([&abs_tol] () { N_VDestroy_Serial (abs_tol); });
616 if (IDASVtolerances (m_mem, reltol, abs_tol) != 0)
617 error (
"IDA: Tolerance not set");
621IDA::set_tolerance (OCTAVE_SUNREALTYPE abstol, OCTAVE_SUNREALTYPE reltol)
623 if (IDASStolerances (m_mem, reltol, abstol) != 0)
624 error (
"IDA: Tolerance not set");
630 const int refine,
bool haverefine,
bool haveoutputfcn,
642 std::string
string =
"";
645 OCTAVE_SUNREALTYPE tsol = tspan(0);
646 OCTAVE_SUNREALTYPE tend = tspan(numt-1);
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);
657 status = IDA::outputfun (output_fcn, haveoutputsel, y,
658 tsol, tend, outputsel,
"init");
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);
671 output.
resize (numt, m_num);
681 if (IDASolve (m_mem, tspan (j), &tsol, yy, yyp, IDA_NORMAL) != 0)
682 error (
"IDASolve failed");
684 yout = NVecToCol (yy, m_num);
685 ypout = NVecToCol (yyp, m_num);
689 output.
elem (j, i) = yout.elem (i);
692 status = IDA::outputfun (output_fcn, haveoutputsel, yout, tsol,
693 tend, outputsel,
string);
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,
704 output.
resize (j + 1, m_num);
720 bool posdirection = (tend > tsol);
723 while (((posdirection == 1 && tsol < tend)
724 || (posdirection == 0 && tsol > tend))
727 if (IDASolve (m_mem, tend, &tsol, yy, yyp, IDA_ONE_STEP) != 0)
728 error (
"IDASolve failed");
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,
739 ypout = NVecToCol (yyp, m_num);
741 output.
resize (cont + 1, m_num);
744 yout = NVecToCol (yy, m_num);
747 output.
elem (cont, i) = yout.elem (i);
749 if (haveoutputfcn && ! haverefine && tout(cont) < tend)
750 status = IDA::outputfun (output_fcn, haveoutputsel, yout, tsol,
751 tend, outputsel,
string);
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,
763 N_Vector dky = N_VNew_Serial (m_num OCTAVE_SUNCONTEXT);
764 octave::unwind_action act2 ([&dky] () { N_VDestroy_Serial (dky); });
766 if (IDAGetDky (m_mem, tend, 0, dky) != 0)
767 error (
"IDA failed to interpolate y");
770 yout = NVecToCol (dky, m_num);
773 output.
elem (cont, i) = yout.elem (i);
778 status = IDA::outputfun (output_fcn, haveoutputsel, yout,
779 tend, tend, outputsel,
string);
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);
791 status = IDA::outputfun (output_fcn, haveoutputsel, yout, tend, tend,
796 return ovl (tout, output, te, ye, ie + 1.0);
812 if (num_event_args == 2)
813 args =
ovl (tsol, y);
815 args =
ovl (tsol, y, yp);
825 oldval = output(0).vector_value ();
826 oldisterminal = output(1).vector_value ();
827 olddir = output(2).vector_value ();
843 && ((val(i) >= 0 && oldval(i) < 0)
844 || (val(i) > 0 && oldval(i) <= 0)))
846 && ((val(i) <= 0 && oldval(i) > 0)
847 || (val(i) < 0 && oldval(i) >= 0))))
849 index.
resize (index.numel () + 1);
850 index (index.numel () - 1) = i;
854 if (cont == 1 && index.numel () > 0)
863 te(0) = tsol - val (index(0)) * (tsol - told)
864 / (val (index(0)) - oldval (index(0)));
867 = y - ((tsol - te(0)) * (y - yold) / (tsol - told));
870 ye.
elem (0, i) = ytemp.elem (i);
873 else if (index.numel () > 0)
877 te.
resize (temp + index.numel ());
878 ye.
resize (temp + index.numel (), m_num);
879 ie.
resize (temp + index.numel ());
884 if (isterminal (index(i)) == 1)
888 ie(temp+i) = index(i);
889 te(temp+i) = tsol - val(index(i)) * (tsol - told)
890 / (val(index(i)) - oldval(index(i)));
893 = y - (tsol - te (temp + i)) * (y - yold) / (tsol - told);
896 ye.
elem (temp + i, j) = ytemp.elem (j);
900 temp += index.numel ();
908 oldisterminal = isterminal;
916 int refine, OCTAVE_SUNREALTYPE tend,
bool haveoutputfcn,
925 OCTAVE_SUNREALTYPE h = 0, tcur = 0;
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);
937 std::string
string =
"";
939 if (IDAGetLastStep (m_mem, &h) != 0)
940 error (
"IDA failed to return last step");
942 if (IDAGetCurrentTime (m_mem, &tcur) != 0)
943 error (
"IDA failed to return the current time");
945 OCTAVE_SUNREALTYPE tin = tcur - h;
947 OCTAVE_SUNREALTYPE step = h / refine;
950 i < refine && tin + step * i < tend && status == 0;
953 if (IDAGetDky (m_mem, tin + step*i, 0, dky) != 0)
954 error (
"IDA failed to interpolate y");
956 if (IDAGetDky (m_mem, tin + step*i, 1, dkyp) != 0)
957 error (
"IDA failed to interpolate yp");
960 output.
resize (cont + 1, m_num);
963 tout(cont) = tin + step * i;
964 yout = NVecToCol (dky, m_num);
965 ypout = NVecToCol (dkyp, m_num);
968 output.
elem (cont, j) = yout.elem (j);
971 status = IDA::outputfun (output_fcn, haveoutputsel, yout,
972 tout(cont), tend, outputsel,
"");
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);
985IDA::outputfun (
const octave_value& output_fcn,
bool haveoutputsel,
988 const std::string& flag)
999 ysel(i) = yout(outputsel(i));
1013 output(0) = toutput;
1015 interp.
feval (output_fcn, output, 0);
1017 else if (flag ==
"")
1021 status = val(0).bool_value ();
1027 interp.
feval (output_fcn, output, 0);
1034IDA::set_maxstep (OCTAVE_SUNREALTYPE maxstep)
1036 if (IDASetMaxStep (m_mem, maxstep) != 0)
1037 error (
"IDA: Max Step not set");
1041IDA::set_initialstep (OCTAVE_SUNREALTYPE initialstep)
1043 if (IDASetInitStep (m_mem, initialstep) != 0)
1044 error (
"IDA: Initial Step not set");
1048IDA::set_maxorder (
int maxorder)
1050 if (IDASetMaxOrd (m_mem, maxorder) != 0)
1051 error (
"IDA: Max Order not set");
1057 long int nsteps = 0, netfails = 0, nrevals = 0;
1059 if (IDAGetNumSteps (m_mem, &nsteps) != 0)
1060 error (
"IDA failed to return the number of internal steps");
1062 if (IDAGetNumErrTestFails (m_mem, &netfails) != 0)
1063 error (
"IDA failed to return the number of internal errors");
1065 if (IDAGetNumResEvals (m_mem, &nrevals) != 0)
1066 error (
"IDA failed to return the number of residual evaluations");
1086 tmp = interp.
feval (ida_fc,
ovl (t,
x, xdot), 1);
1093 return tmp(0).vector_value ();
1106 tmp = interp.
feval (ida_jc,
ovl (t,
x, xdot), 2);
1113 return tmp(0).matrix_value () + cj * tmp(1).matrix_value ();
1126 tmp = interp.
feval (ida_jc,
ovl (t,
x, xdot), 2);
1133 return tmp(0).sparse_matrix_value () + cj * tmp(1).sparse_matrix_value ();
1137ida_dense_cell_jac (
Matrix *dfdy,
Matrix *dfdyp,
double cj)
1139 return (*dfdy) + cj * (*dfdyp);
1146 return (*spdfdy) + cj * (*spdfdyp);
1153 const OCTAVE_SUNREALTYPE t0,
1162 IDA dae (t0, y0, yp0, ida_fcn, ida_user_function);
1171 Matrix ida_dfdy, ida_dfdyp;
1181 dae.set_jacobian (ida_jac, ida_sparse_jac);
1183 dae.set_jacobian (ida_jac, ida_dense_jac);
1191 ida_spdfdy = jaccell(0).sparse_matrix_value ();
1192 ida_spdfdyp = jaccell(1).sparse_matrix_value ();
1194 dae.set_jacobian (&ida_spdfdy, &ida_spdfdyp,
1195 ida_sparse_cell_jac);
1199 ida_dfdy = jaccell(0).matrix_value ();
1200 ida_dfdyp = jaccell(1).matrix_value ();
1202 dae.set_jacobian (&ida_dfdy, &ida_dfdyp, ida_dense_cell_jac);
1219 dae.set_tolerance (abs_tol, rel_tol);
1225 dae.set_tolerance (abs_tol, rel_tol);
1231 dae.set_maxstep (maxstep);
1238 dae.set_initialstep (initialstep);
1244 dae.set_maxorder (maxorder);
1249 bool haverefine = (refine > 1);
1255 bool haveoutputfunction
1258 if (haveoutputfunction)
1259 output_fcn = options.
getfield (
"OutputFcn");
1270 bool haveeventfunction
1273 if (haveeventfunction)
1274 event_fcn = options.
getfield (
"Events");
1280 retval = dae.integrate (numt, tspan, y0, yp0, refine,
1281 haverefine, haveoutputfunction,
1282 output_fcn, haveoutputsel, outputsel,
1283 haveeventfunction, event_fcn, num_event_args);
1303#if defined (HAVE_SUNDIALS)
1313 error (
"__ode15__: FCN must be a function handle");
1317 = args(1).xvector_value (
"__ode15__: TRANGE must be a vector of numbers");
1321 OCTAVE_SUNREALTYPE t0 = tspan(0);
1324 error (
"__ode15__: TRANGE must contain at least 2 elements");
1326 error (
"__ode15__: TRANGE must be strictly monotonic");
1330 = args(2).xvector_value (
"__ode15__: initial state Y0 must be a vector");
1333 = args(3).xvector_value (
"__ode15__: initial state YP0 must be a vector");
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");
1342 if (! args(4).isstruct ())
1343 error (
"__ode15__: ODE_OPT argument must be a structure");
1346 = args(4).xscalar_map_value (
"__ode15__: ODE_OPT argument must be a scalar structure");
1350 = args(5).strict_idx_type_value (
"__ode15__: NUM_EVENT_ARGS must be an integer");
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");
1355 return do_ode15 (ida_fcn, tspan, numt, t0, y0, yp0, options, num_event_args);
1359 octave_unused_parameter (args);
1371OCTAVE_END_NAMESPACE(octave)
T & elem(octave_idx_type n)
Size of the specified dimension.
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.
void resize(octave_idx_type n, const double &rfv=0)
void resize(octave_idx_type nr, octave_idx_type nc, double rfv=0)
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
void resize(octave_idx_type n, const octave_value &rfv=octave_value())
octave_idx_type length() const
bool is_function_handle() const
bool bool_value(bool warn=false) const
int int_value(bool req_int=false, bool frc_str_conv=false) const
Array< double > vector_value(bool frc_str_conv=false, bool frc_vec_conv=false) const
double double_value(bool frc_str_conv=false) const
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.
void error(const char *fmt,...)
void err_disabled_feature(const std::string &fcn, const std::string &feature, const std::string &pkg)
void err_user_supplied_eval(const char *name)
interpreter & __get_interpreter__()
F77_RET_T const F77_INT F77_CMPLX * A
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
F77_RET_T const F77_DBLE * x