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 # if defined (HAVE_UFPARSE_KLU_H)
79 # include <ufsparse/klu.h>
81 # include <sunlinsol/sunlinsol_klu.h>
88 #if defined (HAVE_SUNDIALS)
90 # if ! defined (HAVE_IDASETJACFN) && defined (HAVE_IDADLSSETJACFN)
92 IDASetJacFn (
void *ida_mem, IDADlsJacFn jac)
94 return IDADlsSetJacFn (ida_mem, jac);
98 # if ! defined (HAVE_IDASETLINEARSOLVER) && defined (HAVE_IDADLSSETLINEARSOLVER)
100 IDASetLinearSolver (
void *ida_mem, SUNLinearSolver LS, SUNMatrix
A)
102 return IDADlsSetLinearSolver (ida_mem, LS,
A);
106 # if ! defined (HAVE_SUNLINSOL_DENSE) && defined (HAVE_SUNDENSELINEARSOLVER)
107 static inline SUNLinearSolver
108 SUNLinSol_Dense (N_Vector y, SUNMatrix
A)
110 return SUNDenseLinearSolver (y,
A);
114 # if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
115 # if ! defined (HAVE_SUNLINSOL_KLU) && defined (HAVE_SUNKLU)
116 static inline SUNLinearSolver
117 SUNLinSol_KLU (N_Vector y, SUNMatrix
A)
119 return SUNKLU (y,
A);
124 static inline OCTAVE_SUNREALTYPE *
125 nv_data_s (N_Vector& v)
127 # if defined (HAVE_PRAGMA_GCC_DIAGNOSTIC)
131 # pragma GCC diagnostic push
132 # pragma GCC diagnostic ignored "-Wold-style-cast"
135 return NV_DATA_S (v);
137 # if defined (HAVE_PRAGMA_GCC_DIAGNOSTIC)
139 # pragma GCC diagnostic pop
150 OCTAVE_SUNREALTYPE t,
161 OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
166 OCTAVE_SUNREALTYPE cj);
174 : m_t0 (0.0), m_y0 (), m_yp0 (), m_havejac (false), m_havejacfcn (false),
175 m_havejacsparse (false), m_mem (nullptr), m_num (), m_ida_fcn (),
176 m_ida_jac (), m_dfdy (nullptr), m_dfdyp (nullptr), m_spdfdy (nullptr),
177 m_spdfdyp (nullptr), m_fcn (nullptr), m_jacfcn (nullptr),
178 m_jacspfcn (nullptr), m_jacdcell (nullptr), m_jacspcell (nullptr),
179 m_sunJacMatrix (nullptr), m_sunLinearSolver (nullptr)
185 : m_t0 (t), m_y0 (y), m_yp0 (yp), m_havejac (false), m_havejacfcn (false),
186 m_havejacsparse (false), m_mem (nullptr), m_num (), m_ida_fcn (ida_fcn),
187 m_ida_jac (), m_dfdy (nullptr), m_dfdyp (nullptr), m_spdfdy (nullptr),
188 m_spdfdyp (nullptr), m_fcn (daefun), m_jacfcn (nullptr),
189 m_jacspfcn (nullptr), m_jacdcell (nullptr), m_jacspcell (nullptr),
190 m_sunJacMatrix (nullptr), m_sunLinearSolver (nullptr)
193 OCTAVE_DISABLE_COPY_MOVE (IDA)
198 SUNLinSolFree (m_sunLinearSolver);
199 SUNMatDestroy (m_sunJacMatrix);
200 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
201 SUNContext_Free (&m_sunContext);
206 set_jacobian (
const octave_value& jac, DAEJacFuncDense j)
212 m_havejacsparse =
false;
218 set_jacobian (
const octave_value& jac, DAEJacFuncSparse j)
224 m_havejacsparse =
true;
230 set_jacobian (
Matrix *dy,
Matrix *dyp, DAEJacCellDense j)
236 m_havejacfcn =
false;
237 m_havejacsparse =
false;
250 m_havejacfcn =
false;
251 m_havejacsparse =
true;
256 void set_userdata ();
260 static ColumnVector NVecToCol (N_Vector& v, octave_f77_int_type
n);
262 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
263 N_Vector ColToNVec (
const ColumnVector& data, octave_f77_int_type
n);
265 static N_Vector ColToNVec (
const ColumnVector& data, octave_f77_int_type
n);
272 set_tolerance (
ColumnVector& abstol, OCTAVE_SUNREALTYPE reltol);
275 set_tolerance (OCTAVE_SUNREALTYPE abstol, OCTAVE_SUNREALTYPE reltol);
278 resfun (OCTAVE_SUNREALTYPE t, N_Vector yy, N_Vector yyp,
279 N_Vector rr,
void *user_data);
282 resfun_impl (OCTAVE_SUNREALTYPE t, N_Vector& yy,
283 N_Vector& yyp, N_Vector& rr);
285 jacdense (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj, N_Vector yy,
286 N_Vector yyp, N_Vector, SUNMatrix JJ,
void *user_data, N_Vector,
289 IDA *
self =
static_cast <IDA *
> (user_data);
290 self->jacdense_impl (t, cj, yy, yyp, JJ);
295 jacdense_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
296 N_Vector& yy, N_Vector& yyp, SUNMatrix& JJ);
298 # if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
300 jacsparse (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj, N_Vector yy,
301 N_Vector yyp, N_Vector, SUNMatrix Jac,
void *user_data, N_Vector,
304 IDA *
self =
static_cast <IDA *
> (user_data);
305 self->jacsparse_impl (t, cj, yy, yyp, Jac);
310 jacsparse_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
311 N_Vector& yy, N_Vector& yyp, SUNMatrix& Jac);
314 void set_maxstep (OCTAVE_SUNREALTYPE maxstep);
316 void set_initialstep (OCTAVE_SUNREALTYPE initialstep);
320 int refine, OCTAVE_SUNREALTYPE tend,
bool haveoutputfcn,
330 outputfun (
const octave_value& output_fcn,
bool haveoutputsel,
333 const std::string& flag);
346 void set_maxorder (
int maxorder);
351 const int refine,
bool haverefine,
bool haveoutputfcn,
361 OCTAVE_SUNREALTYPE m_t0;
366 bool m_havejacsparse;
368 octave_f77_int_type m_num;
376 DAEJacFuncDense m_jacfcn;
377 DAEJacFuncSparse m_jacspfcn;
378 DAEJacCellDense m_jacdcell;
379 DAEJacCellSparse m_jacspcell;
380 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
381 SUNContext m_sunContext;
383 SUNMatrix m_sunJacMatrix;
384 SUNLinearSolver m_sunLinearSolver;
388 IDA::resfun (OCTAVE_SUNREALTYPE t, N_Vector yy, N_Vector yyp, N_Vector rr,
391 IDA *
self =
static_cast <IDA *
> (user_data);
392 self->resfun_impl (t, yy, yyp, rr);
397 IDA::resfun_impl (OCTAVE_SUNREALTYPE t, N_Vector& yy,
398 N_Vector& yyp, N_Vector& rr)
406 OCTAVE_SUNREALTYPE *puntrr = nv_data_s (rr);
412 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
413 # define OCTAVE_SUNCONTEXT , m_sunContext
415 # define OCTAVE_SUNCONTEXT
421 N_Vector yy = ColToNVec (y, m_num);
422 octave::unwind_action act ([&yy] () { N_VDestroy_Serial (yy); });
426 # if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
427 # if defined (HAVE_SUNSPARSEMATRIX_REALLOCATE)
431 m_sunJacMatrix = SUNSparseMatrix (m_num, m_num, 0, CSC_MAT
434 octave_f77_int_type max_elems;
436 error (
"Unable to allocate memory for sparse Jacobian");
438 m_sunJacMatrix = SUNSparseMatrix (m_num, m_num, max_elems, CSC_MAT
441 if (! m_sunJacMatrix)
442 error (
"Unable to create sparse Jacobian for Sundials");
444 m_sunLinearSolver = SUNLinSol_KLU (yy, m_sunJacMatrix
446 if (! m_sunLinearSolver)
447 error (
"Unable to create KLU sparse solver");
449 if (IDASetLinearSolver (m_mem, m_sunLinearSolver, m_sunJacMatrix))
450 error (
"Unable to set sparse linear solver");
452 IDASetJacFn (m_mem, IDA::jacsparse);
455 error (
"SUNDIALS SUNLINSOL KLU was unavailable or disabled when "
463 m_sunJacMatrix = SUNDenseMatrix (m_num, m_num OCTAVE_SUNCONTEXT);
464 if (! m_sunJacMatrix)
465 error (
"Unable to create dense Jacobian for Sundials");
467 m_sunLinearSolver = SUNLinSol_Dense (yy, m_sunJacMatrix
469 if (! m_sunLinearSolver)
470 error (
"Unable to create dense linear solver");
472 if (IDASetLinearSolver (m_mem, m_sunLinearSolver, m_sunJacMatrix))
473 error (
"Unable to set dense linear solver");
475 if (m_havejac && IDASetJacFn (m_mem, IDA::jacdense) != 0)
476 error (
"Unable to set dense Jacobian function");
481 IDA::jacdense_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj,
482 N_Vector& yy, N_Vector& yyp, SUNMatrix& JJ)
485 octave_f77_int_type Neq = NV_LENGTH_S (yy);
494 jac = (*m_jacfcn) (y, yp, t, cj, m_ida_jac);
496 jac = (*m_jacdcell) (m_dfdy, m_dfdyp, cj);
498 octave_f77_int_type num_jac = to_f77_int (jac.numel ());
499 std::copy (jac.fortran_vec (),
500 jac.fortran_vec () + num_jac,
501 SUNDenseMatrix_Data (JJ));
504 # if defined (HAVE_SUNDIALS_SUNLINSOL_KLU)
506 IDA::jacsparse_impl (OCTAVE_SUNREALTYPE t, OCTAVE_SUNREALTYPE cj, N_Vector& yy,
507 N_Vector& yyp, SUNMatrix& Jac)
517 jac = (*m_jacspfcn) (y, yp, t, cj, m_ida_jac);
519 jac = (*m_jacspcell) (m_spdfdy, m_spdfdyp, cj);
521 # if defined (HAVE_SUNSPARSEMATRIX_REALLOCATE)
522 octave_f77_int_type nnz = to_f77_int (jac.nnz ());
523 if (nnz > SUNSparseMatrix_NNZ (Jac))
528 if (SUNSparseMatrix_Reallocate (Jac, nnz))
529 error (
"Unable to allocate sufficient memory for IDA sparse matrix");
533 SUNMatZero_Sparse (Jac);
536 sunindextype *colptrs = SUNSparseMatrix_IndexPointers (Jac);
537 sunindextype *rowvals = SUNSparseMatrix_IndexValues (Jac);
539 for (octave_f77_int_type i = 0; i < m_num + 1; i++)
540 colptrs[i] = to_f77_int (jac.cidx (i));
542 double *
d = SUNSparseMatrix_Data (Jac);
543 for (octave_f77_int_type i = 0; i < to_f77_int (jac.nnz ()); i++)
545 rowvals[i] = to_f77_int (jac.ridx (i));
552 IDA::NVecToCol (N_Vector& v, octave_f77_int_type
n)
555 OCTAVE_SUNREALTYPE *punt = nv_data_s (v);
557 for (octave_f77_int_type i = 0; i <
n; i++)
564 IDA::ColToNVec (
const ColumnVector& data, octave_f77_int_type
n)
566 N_Vector v = N_VNew_Serial (
n OCTAVE_SUNCONTEXT);
568 OCTAVE_SUNREALTYPE *punt = nv_data_s (v);
570 for (octave_f77_int_type i = 0; i <
n; i++)
579 void *userdata =
this;
581 if (IDASetUserData (m_mem, userdata) != 0)
582 error (
"User data not set");
588 m_num = to_f77_int (m_y0.numel ());
589 # if defined (HAVE_SUNDIALS_SUNCONTEXT)
590 # if ! defined (SUN_COMM_NULL)
591 # define SUN_COMM_NULL nullptr
593 if (SUNContext_Create (SUN_COMM_NULL, &m_sunContext) < 0)
594 error (
"__ode15__: unable to create context for SUNDIALS");
595 m_mem = IDACreate (m_sunContext);
597 m_mem = IDACreate ();
600 N_Vector yy = ColToNVec (m_y0, m_num);
601 N_Vector yyp = ColToNVec (m_yp0, m_num);
602 octave::unwind_action act ([&yy, &yyp] ()
603 { N_VDestroy_Serial (yy);
604 N_VDestroy_Serial (yyp);
607 IDA::set_userdata ();
609 if (IDAInit (m_mem, IDA::resfun, m_t0, yy, yyp) != 0)
610 error (
"IDA not initialized");
614 IDA::set_tolerance (
ColumnVector& abstol, OCTAVE_SUNREALTYPE reltol)
616 N_Vector abs_tol = ColToNVec (abstol, m_num);
617 octave::unwind_action act ([&abs_tol] () { N_VDestroy_Serial (abs_tol); });
619 if (IDASVtolerances (m_mem, reltol, abs_tol) != 0)
620 error (
"IDA: Tolerance not set");
624 IDA::set_tolerance (OCTAVE_SUNREALTYPE abstol, OCTAVE_SUNREALTYPE reltol)
626 if (IDASStolerances (m_mem, reltol, abstol) != 0)
627 error (
"IDA: Tolerance not set");
633 const int refine,
bool haverefine,
bool haveoutputfcn,
645 std::string
string =
"";
648 OCTAVE_SUNREALTYPE tsol = tspan(0);
649 OCTAVE_SUNREALTYPE tend = tspan(numt-1);
651 N_Vector yyp = ColToNVec (yp, m_num);
652 N_Vector yy = ColToNVec (y, m_num);
653 octave::unwind_action act ([&yyp, &yy] ()
654 { N_VDestroy_Serial (yyp);
655 N_VDestroy_Serial (yy);
660 status = IDA::outputfun (output_fcn, haveoutputsel, y,
661 tsol, tend, outputsel,
"init");
664 if (haveeventfunction)
665 status = IDA::event (event_fcn, te, ye, ie, tsol, y,
666 "init", yp, oldval, oldisterminal,
667 olddir, cont, temp, tsol, yold, num_event_args);
674 output.
resize (numt, m_num);
684 if (IDASolve (m_mem, tspan (j), &tsol, yy, yyp, IDA_NORMAL) != 0)
685 error (
"IDASolve failed");
687 yout = NVecToCol (yy, m_num);
688 ypout = NVecToCol (yyp, m_num);
692 output.
elem (j, i) = yout.elem (i);
695 status = IDA::outputfun (output_fcn, haveoutputsel, yout, tsol,
696 tend, outputsel,
string);
698 if (haveeventfunction)
699 status = IDA::event (event_fcn, te, ye, ie, tout(j), yout,
700 string, ypout, oldval, oldisterminal,
701 olddir, j, temp, tout(j-1), yold,
707 output.
resize (j + 1, m_num);
723 bool posdirection = (tend > tsol);
726 while (((posdirection == 1 && tsol < tend)
727 || (posdirection == 0 && tsol > tend))
730 if (IDASolve (m_mem, tend, &tsol, yy, yyp, IDA_ONE_STEP) != 0)
731 error (
"IDASolve failed");
734 status = IDA::interpolate (cont, output, tout, refine, tend,
735 haveoutputfcn, haveoutputsel,
736 output_fcn, outputsel,
737 haveeventfunction, event_fcn, te,
738 ye, ie, oldval, oldisterminal,
742 ypout = NVecToCol (yyp, m_num);
744 output.
resize (cont + 1, m_num);
747 yout = NVecToCol (yy, m_num);
750 output.
elem (cont, i) = yout.elem (i);
752 if (haveoutputfcn && ! haverefine && tout(cont) < tend)
753 status = IDA::outputfun (output_fcn, haveoutputsel, yout, tsol,
754 tend, outputsel,
string);
756 if (haveeventfunction && ! haverefine && tout(cont) < tend)
757 status = IDA::event (event_fcn, te, ye, ie, tout(cont), yout,
758 string, ypout, oldval, oldisterminal,
759 olddir, cont, temp, tout(cont-1), yold,
766 N_Vector dky = N_VNew_Serial (m_num OCTAVE_SUNCONTEXT);
767 octave::unwind_action act2 ([&dky] () { N_VDestroy_Serial (dky); });
769 if (IDAGetDky (m_mem, tend, 0, dky) != 0)
770 error (
"IDA failed to interpolate y");
773 yout = NVecToCol (dky, m_num);
776 output.
elem (cont, i) = yout.elem (i);
781 status = IDA::outputfun (output_fcn, haveoutputsel, yout,
782 tend, tend, outputsel,
string);
785 if (haveeventfunction)
786 status = IDA::event (event_fcn, te, ye, ie, tend, yout,
787 string, ypout, oldval, oldisterminal,
788 olddir, cont, temp, tout(cont-1),
789 yold, num_event_args);
794 status = IDA::outputfun (output_fcn, haveoutputsel, yout, tend, tend,
799 return ovl (tout, output, te, ye, ie + 1.0);
815 if (num_event_args == 2)
816 args =
ovl (tsol, y);
818 args =
ovl (tsol, y, yp);
828 oldval = output(0).vector_value ();
829 oldisterminal = output(1).vector_value ();
830 olddir = output(2).vector_value ();
846 && ((val(i) >= 0 && oldval(i) < 0)
847 || (val(i) > 0 && oldval(i) <= 0)))
849 && ((val(i) <= 0 && oldval(i) > 0)
850 || (val(i) < 0 && oldval(i) >= 0))))
852 index.
resize (index.numel () + 1);
853 index (index.numel () - 1) = i;
857 if (cont == 1 && index.numel () > 0)
866 te(0) = tsol - val (index(0)) * (tsol - told)
867 / (val (index(0)) - oldval (index(0)));
870 = y - ((tsol - te(0)) * (y - yold) / (tsol - told));
873 ye.
elem (0, i) = ytemp.elem (i);
876 else if (index.numel () > 0)
880 te.
resize (temp + index.numel ());
881 ye.
resize (temp + index.numel (), m_num);
882 ie.
resize (temp + index.numel ());
887 if (isterminal (index(i)) == 1)
891 ie(temp+i) = index(i);
892 te(temp+i) = tsol - val(index(i)) * (tsol - told)
893 / (val(index(i)) - oldval(index(i)));
896 = y - (tsol - te (temp + i)) * (y - yold) / (tsol - told);
899 ye.
elem (temp + i, j) = ytemp.elem (j);
903 temp += index.numel ();
911 oldisterminal = isterminal;
919 int refine, OCTAVE_SUNREALTYPE tend,
bool haveoutputfcn,
928 OCTAVE_SUNREALTYPE h = 0, tcur = 0;
931 N_Vector dky = N_VNew_Serial (m_num OCTAVE_SUNCONTEXT);
932 N_Vector dkyp = N_VNew_Serial (m_num OCTAVE_SUNCONTEXT);
933 octave::unwind_action act ([&dky, &dkyp] ()
934 { N_VDestroy_Serial (dky);
935 N_VDestroy_Serial (dkyp);
940 std::string
string =
"";
942 if (IDAGetLastStep (m_mem, &h) != 0)
943 error (
"IDA failed to return last step");
945 if (IDAGetCurrentTime (m_mem, &tcur) != 0)
946 error (
"IDA failed to return the current time");
948 OCTAVE_SUNREALTYPE tin = tcur - h;
950 OCTAVE_SUNREALTYPE step = h / refine;
953 i < refine && tin + step * i < tend && status == 0;
956 if (IDAGetDky (m_mem, tin + step*i, 0, dky) != 0)
957 error (
"IDA failed to interpolate y");
959 if (IDAGetDky (m_mem, tin + step*i, 1, dkyp) != 0)
960 error (
"IDA failed to interpolate yp");
963 output.
resize (cont + 1, m_num);
966 tout(cont) = tin + step * i;
967 yout = NVecToCol (dky, m_num);
968 ypout = NVecToCol (dkyp, m_num);
971 output.
elem (cont, j) = yout.elem (j);
974 status = IDA::outputfun (output_fcn, haveoutputsel, yout,
975 tout(cont), tend, outputsel,
"");
977 if (haveeventfunction)
978 status = IDA::event (event_fcn, te, ye, ie, tout(cont),
979 yout,
string, ypout, oldval,
980 oldisterminal, olddir, cont, temp,
981 tout(cont-1), yold, num_event_args);
988 IDA::outputfun (
const octave_value& output_fcn,
bool haveoutputsel,
991 const std::string& flag)
1002 ysel(i) = yout(outputsel(i));
1016 output(0) = toutput;
1018 interp.
feval (output_fcn, output, 0);
1020 else if (flag ==
"")
1024 status = val(0).bool_value ();
1030 interp.
feval (output_fcn, output, 0);
1037 IDA::set_maxstep (OCTAVE_SUNREALTYPE maxstep)
1039 if (IDASetMaxStep (m_mem, maxstep) != 0)
1040 error (
"IDA: Max Step not set");
1044 IDA::set_initialstep (OCTAVE_SUNREALTYPE initialstep)
1046 if (IDASetInitStep (m_mem, initialstep) != 0)
1047 error (
"IDA: Initial Step not set");
1051 IDA::set_maxorder (
int maxorder)
1053 if (IDASetMaxOrd (m_mem, maxorder) != 0)
1054 error (
"IDA: Max Order not set");
1060 long int nsteps = 0, netfails = 0, nrevals = 0;
1062 if (IDAGetNumSteps (m_mem, &nsteps) != 0)
1063 error (
"IDA failed to return the number of internal steps");
1065 if (IDAGetNumErrTestFails (m_mem, &netfails) != 0)
1066 error (
"IDA failed to return the number of internal errors");
1068 if (IDAGetNumResEvals (m_mem, &nrevals) != 0)
1069 error (
"IDA failed to return the number of residual evaluations");
1089 tmp = interp.
feval (ida_fc,
ovl (t,
x, xdot), 1);
1096 return tmp(0).vector_value ();
1109 tmp = interp.
feval (ida_jc,
ovl (t,
x, xdot), 2);
1116 return tmp(0).matrix_value () + cj * tmp(1).matrix_value ();
1129 tmp = interp.
feval (ida_jc,
ovl (t,
x, xdot), 2);
1136 return tmp(0).sparse_matrix_value () + cj * tmp(1).sparse_matrix_value ();
1140 ida_dense_cell_jac (
Matrix *dfdy,
Matrix *dfdyp,
double cj)
1142 return (*dfdy) + cj * (*dfdyp);
1149 return (*spdfdy) + cj * (*spdfdyp);
1156 const OCTAVE_SUNREALTYPE t0,
1165 IDA dae (t0, y0, yp0, ida_fcn, ida_user_function);
1174 Matrix ida_dfdy, ida_dfdyp;
1184 dae.set_jacobian (ida_jac, ida_sparse_jac);
1186 dae.set_jacobian (ida_jac, ida_dense_jac);
1194 ida_spdfdy = jaccell(0).sparse_matrix_value ();
1195 ida_spdfdyp = jaccell(1).sparse_matrix_value ();
1197 dae.set_jacobian (&ida_spdfdy, &ida_spdfdyp,
1198 ida_sparse_cell_jac);
1202 ida_dfdy = jaccell(0).matrix_value ();
1203 ida_dfdyp = jaccell(1).matrix_value ();
1205 dae.set_jacobian (&ida_dfdy, &ida_dfdyp, ida_dense_cell_jac);
1222 dae.set_tolerance (abs_tol, rel_tol);
1228 dae.set_tolerance (abs_tol, rel_tol);
1234 dae.set_maxstep (maxstep);
1241 dae.set_initialstep (initialstep);
1247 dae.set_maxorder (maxorder);
1252 bool haverefine = (refine > 1);
1258 bool haveoutputfunction
1261 if (haveoutputfunction)
1262 output_fcn = options.
getfield (
"OutputFcn");
1273 bool haveeventfunction
1276 if (haveeventfunction)
1277 event_fcn = options.
getfield (
"Events");
1283 retval = dae.integrate (numt, tspan, y0, yp0, refine,
1284 haverefine, haveoutputfunction,
1285 output_fcn, haveoutputsel, outputsel,
1286 haveeventfunction, event_fcn, num_event_args);
1306 #if defined (HAVE_SUNDIALS)
1316 error (
"__ode15__: FCN must be a function handle");
1320 = args(1).xvector_value (
"__ode15__: TRANGE must be a vector of numbers");
1324 OCTAVE_SUNREALTYPE t0 = tspan(0);
1327 error (
"__ode15__: TRANGE must contain at least 2 elements");
1329 error (
"__ode15__: TRANGE must be strictly monotonic");
1333 = args(2).xvector_value (
"__ode15__: initial state Y0 must be a vector");
1336 = args(3).xvector_value (
"__ode15__: initial state YP0 must be a vector");
1340 error (
"__ode15__: initial state Y0 and YP0 must have the same length");
1341 else if (y0.
numel () < 1)
1342 error (
"__ode15__: initial state YP0 must be a vector or a scalar");
1345 if (! args(4).isstruct ())
1346 error (
"__ode15__: ODE_OPT argument must be a structure");
1349 = args(4).xscalar_map_value (
"__ode15__: ODE_OPT argument must be a scalar structure");
1353 = args(5).xidx_type_value (
"__ode15__: NUM_EVENT_ARGS must be an integer");
1355 if (num_event_args != 2 && num_event_args != 3)
1356 error (
"__ode15__: number of input arguments in event callback must be 2 or 3");
1358 return do_ode15 (ida_fcn, tspan, numt, t0, y0, yp0, options, num_event_args);
1362 octave_unused_parameter (args);
1374 OCTAVE_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
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
F77_RET_T const F77_DBLE * x
bool int_multiply_overflow(int a, int b, int *r)
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.