26 #include "idas_interface.hpp"
27 #include "casadi/core/casadi_misc.hpp"
30 #define THROWING(fcn, ...) \
31 idas_error(CASADI_STR(fcn), fcn(__VA_ARGS__))
36 int CASADI_INTEGRATOR_IDAS_EXPORT
39 plugin->name =
"idas";
41 plugin->version = CASADI_VERSION;
53 double t0,
const std::vector<double>& tout) :
SundialsInterface(name, dae, t0, tout) {
62 {{
"suppress_algebraic",
64 "Suppress algebraic variables in the error testing"}},
67 "Use IDACalcIC to get consistent initial conditions."}},
70 "Constrain the solution y=[x,z]. 0 (default): no constraint on yi, "
71 "1: yi >= 0.0, -1: yi <= 0.0, 2: yi > 0.0, -2: yi < 0.0."}},
74 "Use IDACalcIC to get consistent initial conditions for "
75 "backwards system [default: equal to calc_ic]."}},
78 "Absolute tolerarance for each component"}},
81 "Maximim step size"}},
84 "First requested time as a fraction of the time interval"}},
87 "IDAS scaling on cj for the user-defined linear solver module"}},
90 "Initial values for the state derivatives"}}
106 for (
auto&& op : opts) {
107 if (op.first==
"init_xdot") {
109 }
else if (op.first==
"cj_scaling") {
111 }
else if (op.first==
"calc_ic") {
113 }
else if (op.first==
"suppress_algebraic") {
115 }
else if (op.first==
"constraints") {
117 }
else if (op.first==
"abstolv") {
127 for (
auto&& op : opts) {
128 if (op.first==
"calc_icB") {
130 }
else if (op.first==
"first_time") {
141 "Option \"init_xdot\" has incorrect length. Expecting " +
str(
nx_) +
", "
143 "Note that this message may actually be generated by the augmented integrator. "
144 "In that case, make use of the 'augmented_options' options "
145 "to correct 'init_xdot' for the augmented integrator.");
150 "Constraint vector if supplied, must be of length nx+nz, but got "
159 auto m =
to_mem(user_data);
161 if (s.calc_daeF(m, t, NV_DATA_S(xz), NV_DATA_S(xz) + s.nx_,
162 NV_DATA_S(rr), NV_DATA_S(rr) + s.nx_))
return 1;
165 casadi_axpy(s.nx_, -1., NV_DATA_S(xzdot), NV_DATA_S(rr));
167 }
catch(std::exception& e) {
168 uerr() <<
"res failed: " << e.what() << std::endl;
174 char *msg,
void *eh_data) {
178 uerr() << msg << std::endl;
179 }
catch(std::exception& e) {
180 uerr() <<
"ehfun failed: " << e.what() << std::endl;
185 N_Vector Jv,
double cj,
void *user_data, N_Vector tmp1, N_Vector tmp2) {
187 auto m =
to_mem(user_data);
189 if (s.calc_jtimesF(m, t, NV_DATA_S(xz), NV_DATA_S(xz) + s.nx_,
190 NV_DATA_S(v), NV_DATA_S(v) + s.nx_,
191 NV_DATA_S(Jv), NV_DATA_S(Jv) + s.nx_))
return 1;
194 casadi_axpy(s.nx_, -cj, NV_DATA_S(v), NV_DATA_S(Jv));
197 }
catch(std::exception& e) {
198 uerr() <<
"jtimesF failed: " << e.what() << std::endl;
204 N_Vector rxzdot, N_Vector resvalB, N_Vector v, N_Vector Jv,
205 double cjB,
void *user_data, N_Vector tmp1B, N_Vector tmp2B) {
207 auto m =
to_mem(user_data);
210 if (s.calc_daeB(m, t, NV_DATA_S(xz), NV_DATA_S(xz) + s.nx_,
211 NV_DATA_S(v), NV_DATA_S(v) + s.nrx_,
nullptr,
212 NV_DATA_S(Jv), NV_DATA_S(Jv) + s.nrx_))
return 1;
214 casadi_axpy(s.nrx_, cjB, NV_DATA_S(v), NV_DATA_S(Jv));
217 }
catch(std::exception& e) {
218 uerr() <<
"jtimesB failed: " << e.what() << std::endl;
228 m->mem = IDACreate();
229 casadi_assert(m->mem!=
nullptr,
"IDACreate: Creation failed");
232 THROWING(IDASetErrHandlerFn, m->mem,
ehfun, m);
235 THROWING(IDASetUserData, m->mem, m);
238 m->v_xzdot = N_VNew_Serial(
nx_+
nz_);
242 N_VConst(0.0, m->v_xz);
243 N_VConst(0.0, m->v_xzdot);
244 IDAInit(m->mem,
resF, t0, m->v_xz, m->v_xzdot);
245 if (
verbose_) casadi_message(
"IDA initialized");
254 if (
step0_!=0) THROWING(IDASetInitStep, m->mem,
step0_);
261 N_Vector domain = N_VNew_Serial(
nx_+
nz_);
262 std::copy(
y_c_.begin(),
y_c_.end(), NV_DATA_S(domain));
265 int flag = IDASetConstraints(m->mem, domain);
266 casadi_assert_dev(flag==IDA_SUCCESS);
269 N_VDestroy_Serial(domain);
281 N_Vector nv_abstol = N_VNew_Serial(
static_cast<long>(
abstolv_.size()));
283 THROWING(IDASVtolerances, m->mem,
reltol_, nv_abstol);
284 N_VDestroy_Serial(nv_abstol);
287 THROWING(IDASVtolerances, m->mem,
reltol_, m->abstolv);
297 N_Vector
id = N_VNew_Serial(
nx_+
nz_);
298 std::fill_n(NV_DATA_S(
id),
nx_, 1);
299 std::fill_n(NV_DATA_S(
id)+
nx_,
nz_, 0);
302 THROWING(IDASetId, m->mem,
id);
305 N_VDestroy_Serial(
id);
310 IDAMem IDA_mem = IDAMem(m->mem);
311 IDA_mem->ida_lmem = m;
314 IDA_mem->ida_setupNonNull = TRUE;
323 THROWING(IDASpilsSetJacTimesVecFn, m->mem,
jtimesF);
331 THROWING(IDAQuadInit, m->mem,
rhsQF, m->v_q);
335 THROWING(IDASetQuadErrCon, m->mem,
true);
343 if (
verbose_) casadi_message(
"Attached linear solver");
347 m->v_adj_xzdot = N_VNew_Serial(
nrx_+
nrz_);
348 N_VConst(0.0, m->v_adj_xz);
349 N_VConst(0.0, m->v_adj_xzdot);
351 if (
verbose_) casadi_message(
"Initialized adjoint sensitivities");
359 m->first_callB =
true;
374 N_VConst(0.0, m->v_xzdot);
377 THROWING(IDAReInit, m->mem, m->t, m->v_xz, m->v_xzdot);
380 if (
nq1_ > 0) THROWING(IDAQuadReInit, m->mem, m->v_q);
384 THROWING(IDACalcIC, m->mem, IDA_YA_YDP_INIT ,
first_time_);
385 THROWING(IDAGetConsistentIC, m->mem, m->v_xz, m->v_xzdot);
389 if (
nadj_ > 0) THROWING(IDAAdjReInit, m->mem);
399 if (m->t_stop >= m->tcur) {
400 THROWING(IDASetStopTime, m->mem, m->t_stop);
405 if (fabs(m->t - m->t_next) >= ttol) {
409 THROWING(IDASolveF, m->mem, m->t_next, &tret, m->v_xz, m->v_xzdot, IDA_NORMAL, &m->ncheck);
411 THROWING(IDASolve, m->mem, m->t_next, &tret, m->v_xz, m->v_xzdot, IDA_NORMAL);
414 if (
nq1_ > 0) THROWING(IDAGetQuad, m->mem, &tret, m->v_q);
422 THROWING(IDAGetIntegratorStats, m->mem, &m->nsteps, &m->nfevals, &m->nlinsetups,
423 &m->netfails, &m->qlast, &m->qcur, &m->hinused, &m->hlast, &m->hcur, &m->tcur);
424 THROWING(IDAGetNonlinSolvStats, m->mem, &m->nniters, &m->nncfails);
434 N_VConst(0.0, m->v_adj_xz);
440 N_VConst(0.0, m->v_adj_xzdot);
460 casadi_error(
"Linear system factorization for backwards initial conditions failed");
466 casadi_error(
"Linear system solve for backwards initial conditions failed");
472 casadi_error(
"Adjoint seed propagation for backwards initial conditions failed");
479 const double* adj_x,
const double* adj_z,
const double* adj_q)
const {
488 if (m->first_callB) {
490 THROWING(IDACreateB, m->mem, &m->whichB);
491 THROWING(IDAInitB, m->mem, m->whichB,
resB, m->t, m->v_adj_xz, m->v_adj_xzdot);
493 THROWING(IDASetUserDataB, m->mem, m->whichB, m);
498 N_Vector
id = N_VNew_Serial(
nrx_+
nrz_);
499 std::fill_n(NV_DATA_S(
id),
nrx_, 1);
500 std::fill_n(NV_DATA_S(
id)+
nrx_,
nrz_, 0);
501 THROWING(IDASetIdB, m->mem, m->whichB,
id);
502 N_VDestroy_Serial(
id);
507 IDAMem IDA_mem = IDAMem(m->mem);
508 IDAadjMem IDAADJ_mem = IDA_mem->ida_adj_mem;
509 IDABMem IDAB_mem = IDAADJ_mem->IDAB_mem;
510 IDAB_mem->ida_lmem = m;
511 IDAB_mem->IDA_mem->ida_lmem = m;
512 IDAB_mem->IDA_mem->ida_lsetup =
lsetupB;
513 IDAB_mem->IDA_mem->ida_lsolve =
lsolveB;
514 IDAB_mem->IDA_mem->ida_setupNonNull = TRUE;
523 THROWING(IDASpilsSetJacTimesVecFnB, m->mem, m->whichB,
jtimesB);
528 THROWING(IDAQuadInitB, m->mem, m->whichB,
rhsQB, m->v_adj_pu);
530 THROWING(IDASetQuadErrConB, m->mem, m->whichB,
true);
531 THROWING(IDAQuadSStolerancesB, m->mem, m->whichB,
reltol_,
abstol_);
535 m->first_callB =
false;
538 THROWING(IDAReInitB, m->mem, m->whichB, m->t, m->v_adj_xz, m->v_adj_xzdot);
542 void* memB = IDAGetAdjIDABmem(m->mem, m->whichB);
543 THROWING(IDAQuadReInit, memB, m->v_adj_pu);
549 THROWING(IDACalcICB, m->mem, m->whichB,
t0_, m->v_xz, m->v_xzdot);
550 THROWING(IDAGetConsistentICB, m->mem, m->whichB, m->v_adj_xz, m->v_adj_xzdot);
555 double* adj_x,
double* adj_p,
double* adj_u)
const {
562 if (m->t_next < m->t) {
564 THROWING(IDASolveB, m->mem, m->t_next, IDA_NORMAL);
565 THROWING(IDAGetB, m->mem, m->whichB, &tret, m->v_adj_xz, m->v_adj_xzdot);
567 THROWING(IDAGetQuadB, m->mem, m->whichB, &tret, m->v_adj_pu);
570 THROWING(IDAGetAdjY, m->mem, m->t_next, m->v_xz, m->v_xzdot);
579 IDAMem IDA_mem = IDAMem(m->mem);
580 IDAadjMem IDAADJ_mem = IDA_mem->ida_adj_mem;
581 IDABMem IDAB_mem = IDAADJ_mem->IDAB_mem;
582 THROWING(IDAGetIntegratorStats, IDAB_mem->IDA_mem, &m->nstepsB, &m->nfevalsB,
583 &m->nlinsetupsB, &m->netfailsB, &m->qlastB, &m->qcurB, &m->hinusedB,
584 &m->hlastB, &m->hcurB, &m->tcurB);
585 THROWING(IDAGetNonlinSolvStats, IDAB_mem->IDA_mem, &m->nnitersB, &m->nncfailsB);
593 if (flag>=IDA_SUCCESS)
return;
595 char* flagname = IDAGetReturnFlagName(flag);
596 std::stringstream ss;
597 ss << module <<
" returned \"" << flagname <<
"\". Consult IDAS documentation.";
599 casadi_error(ss.str());
604 auto m =
to_mem(user_data);
606 if (s.calc_quadF(m, t, NV_DATA_S(xz), NV_DATA_S(xz) + s.nx_, NV_DATA_S(qdot)))
return 1;
609 }
catch(std::exception& e) {
610 uerr() <<
"rhsQ failed: " << e.what() << std::endl;
616 N_Vector rxzdot, N_Vector rr,
void *user_data) {
618 auto m =
to_mem(user_data);
620 if (s.calc_daeB(m, t, NV_DATA_S(xz), NV_DATA_S(xz) + s.nx_,
621 NV_DATA_S(rxz), NV_DATA_S(rxz) + s.nrx_, m->adj_q,
622 NV_DATA_S(rr), NV_DATA_S(rr) + s.nrx_))
return 1;
625 casadi_axpy(s.nrx_, 1., NV_DATA_S(rxzdot), NV_DATA_S(rr));
628 }
catch(std::exception& e) {
629 uerr() <<
"resB failed: " << e.what() << std::endl;
635 N_Vector rxzdot, N_Vector ruqdot,
void *user_data) {
637 auto m =
to_mem(user_data);
639 if (s.calc_quadB(m, t, NV_DATA_S(xz), NV_DATA_S(xz) + s.nx_,
640 NV_DATA_S(rxz), NV_DATA_S(rxz) + s.nrx_,
641 NV_DATA_S(ruqdot), NV_DATA_S(ruqdot) + s.nrq_))
return 1;
644 casadi_scal(s.nrq_ + s.nuq_, -1., NV_DATA_S(ruqdot));
647 }
catch(std::exception& e) {
648 uerr() <<
"resQB failed: " << e.what() << std::endl;
654 N_Vector rvec, N_Vector zvec,
double cj,
double delta,
void *user_data, N_Vector tmp) {
656 auto m =
to_mem(user_data);
660 double* vx = NV_DATA_S(rvec);
661 double* vz = vx + s.nx_;
662 double* v_it = m->tmp1;
663 for (
int d = 0; d <= s.nfwd_; ++d) {
671 if (s.linsolF_.solve(m->jacF, m->tmp1, 1,
false, m->mem_linsolF))
673 vx = NV_DATA_S(zvec);
681 if (s.second_order_correction_) {
685 if (s.calc_jtimesF(m, t, NV_DATA_S(xz), NV_DATA_S(xz) + s.nx_,
686 vx, vz, m->tmp2, m->tmp2 + s.nx_))
return 1;
689 v_it = m->tmp1 + s.nx1_ + s.nz1_;
690 for (
int d = 1; d <= s.nfwd_; ++d) {
691 casadi_axpy(s.nx1_, -1., m->tmp2 + d*s.nx1_, v_it);
693 casadi_axpy(s.nz1_, -1., m->tmp2 + s.nx_ + d*s.nz1_, v_it);
699 if (s.linsolF_.solve(m->jacF, m->tmp1 + s.nx1_ + s.nz1_, s.nfwd_,
700 false, m->mem_linsolF))
return 1;
703 v_it = m->tmp1 + s.nx1_ + s.nz1_;
704 for (
int d = 1; d <= s.nfwd_; ++d) {
713 }
catch(std::exception& e) {
714 uerr() <<
"psolve failed: " << e.what() << std::endl;
720 const double* rhs,
double* sol)
const {
722 double* v_it = m->
tmp1;
723 for (
int d = 0; d <=
nfwd_; ++d) {
724 for (
int a = 0; a <
nadj_; ++a) {
734 for (
int a = 0; a <
nadj_; ++a) {
753 for (
int d = 1; d <=
nfwd_; ++d) {
754 for (
int a = 0; a <
nadj_; ++a) {
769 for (
int d = 1; d <=
nfwd_; ++d) {
770 for (
int a = 0; a <
nadj_; ++a) {
783 N_Vector xzdotB, N_Vector resvalB, N_Vector rvecB,
784 N_Vector zvecB,
double cjB,
double deltaB,
785 void *user_data, N_Vector tmpB) {
787 auto m =
to_mem(user_data);
789 return s.solve_transposed(m, t, NV_DATA_S(xz), NV_DATA_S(xzB),
790 NV_DATA_S(rvecB), NV_DATA_S(zvecB));
792 }
catch(std::exception& e) {
793 uerr() <<
"psolveB failed: " << e.what() << std::endl;
798 template<
typename T1>
800 casadi_int r_begin, casadi_int c_begin, T1* w) {
802 casadi_int nrow_x, ncol_x, ncol_y, i_x, i_y, j, el, r_end;
803 const casadi_int *colind_x, *row_x, *colind_y, *row_y;
806 colind_x = sp_x+2; row_x = sp_x + 2 + ncol_x+1;
808 colind_y = sp_y+2; row_y = sp_y + 2 + ncol_y+1;
810 r_end = r_begin + nrow_x;
814 for (i_x = 0; i_x < ncol_x; ++i_x) {
818 for (el=colind_x[i_x]; el<colind_x[i_x + 1]; ++el) w[row_x[el]] = x[el];
820 for (el=colind_y[i_y]; el<colind_y[i_y + 1]; ++el) {
822 if (j >= r_begin && j < r_end) y[el] = w[j - r_begin];
825 for (el=colind_x[i_x]; el<colind_x[i_x + 1]; ++el) w[row_x[el]] = 0;
830 double cj,
void* user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3) {
832 auto m =
to_mem(user_data);
841 const Sparsity& sp_jacF = s.linsolF_.sparsity();
844 if (s.calc_jacF(m, t, NV_DATA_S(xz), NV_DATA_S(xz) + s.nx_,
845 m->jac_ode_x, m->jac_alg_x, m->jac_ode_z, m->jac_alg_z))
return 1;
848 casadi_int nx_jac = sp_jac_ode_x.
size1();
850 casadi_copy_block(m->jac_alg_x, sp_jac_alg_x, m->jacF, sp_jacF, nx_jac, 0, m->w);
851 casadi_copy_block(m->jac_ode_z, sp_jac_ode_z, m->jacF, sp_jacF, 0, nx_jac, m->w);
852 casadi_copy_block(m->jac_alg_z, sp_jac_alg_z, m->jacF, sp_jacF, nx_jac, nx_jac, m->w);
855 const casadi_int *colind = sp_jacF.
colind(), *row = sp_jacF.
row();
856 for (casadi_int c = 0; c < nx_jac; ++c) {
857 for (casadi_int k = colind[c]; k < colind[c + 1]; ++k) {
858 if (row[k] == c) m->jacF[k] -= cj;
863 if (s.linsolF_.nfact(m->jacF, m->mem_linsolF))
return 1;
867 }
catch(std::exception& e) {
868 uerr() <<
"psetup failed: " << e.what() << std::endl;
874 N_Vector rresval,
double cj,
void *user_data, N_Vector tmp1B, N_Vector tmp2B, N_Vector tmp3B) {
877 return psetupF(t, xz,
nullptr,
nullptr, -cj, user_data, tmp1B, tmp2B, tmp3B);
879 }
catch(std::exception& e) {
880 uerr() <<
"psetupB failed: " << e.what() << std::endl;
886 N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3) {
888 double t = IDA_mem->ida_tn;
891 double cj = IDA_mem->ida_cj;
894 return psetupF(t, xz, xzdot,
nullptr, cj, IDA_mem->ida_lmem,
895 vtemp1, vtemp1, vtemp3);
899 N_Vector vtemp1B, N_Vector vtemp2B, N_Vector vtemp3B) {
901 auto m =
to_mem(IDA_mem->ida_lmem);
903 IDAadjMem IDAADJ_mem;
907 double t = IDA_mem->ida_tn;
909 double cj = IDA_mem->ida_cj;
911 IDA_mem =
static_cast<IDAMem
>(IDA_mem->ida_user_data);
913 IDAADJ_mem = IDA_mem->ida_adj_mem;
917 if (IDAADJ_mem->ia_noInterp==FALSE) {
918 int flag = IDAADJ_mem->ia_getY(IDA_mem, t, IDAADJ_mem->ia_yyTmp, IDAADJ_mem->ia_ypTmp,
920 if (flag != IDA_SUCCESS) casadi_error(
"Could not interpolate forward states");
923 return psetupB(t, IDAADJ_mem->ia_yyTmp, IDAADJ_mem->ia_ypTmp,
924 xzB, xzdotB,
nullptr, cj,
static_cast<void*
>(m), vtemp1B, vtemp1B, vtemp3B);
926 }
catch(std::exception& e) {
927 uerr() <<
"lsetupB failed: " << e.what() << std::endl;
933 N_Vector xzdot, N_Vector rr) {
935 auto m =
to_mem(IDA_mem->ida_lmem);
939 double t = IDA_mem->ida_tn;
942 double cj = IDA_mem->ida_cj;
948 int flag =
psolveF(t, xz, xzdot, rr, b, b, cj,
949 delta,
static_cast<void*
>(m),
nullptr);
950 if (flag)
return flag;
954 double cjratio = IDA_mem->ida_cjratio;
955 if (cjratio != 1.0) N_VScale(2.0/(1.0 + cjratio), b, b);
959 }
catch(std::exception& e) {
960 uerr() <<
"lsolve failed: " << e.what() << std::endl;
966 N_Vector xzdotB, N_Vector rrB) {
968 auto m =
to_mem(IDA_mem->ida_lmem);
970 IDAadjMem IDAADJ_mem;
975 double t = IDA_mem->ida_tn;
977 double cj = IDA_mem->ida_cj;
978 double cjratio = IDA_mem->ida_cjratio;
980 IDA_mem = (IDAMem) IDA_mem->ida_user_data;
981 IDAADJ_mem = IDA_mem->ida_adj_mem;
985 if (IDAADJ_mem->ia_noInterp==FALSE) {
986 flag = IDAADJ_mem->ia_getY(IDA_mem, t, IDAADJ_mem->ia_yyTmp, IDAADJ_mem->ia_ypTmp,
988 if (flag != IDA_SUCCESS) casadi_error(
"Could not interpolate forward states");
995 flag =
psolveB(t, IDAADJ_mem->ia_yyTmp, IDAADJ_mem->ia_ypTmp, xzB, xzdotB,
996 rrB, b, b, cj, delta,
static_cast<void*
>(m),
nullptr);
997 if (flag)
return flag;
1000 if (s.cj_scaling_) {
1001 if (cjratio != 1.0) N_VScale(2.0/(1.0 + cjratio), b, b);
1004 }
catch(std::exception& e) {
1005 uerr() <<
"lsolveB failed: " << e.what() << std::endl;
1011 this->
mem =
nullptr;
1021 if (this->
mem) IDAFree(&this->
mem);
1028 int version = s.
version(
"IdasInterface", 1, 2);
1047 s.
version(
"IdasInterface", 2);
1056 s.
pack(
"IdasInterface::y_c",
y_c_);
Helper class for Serialization.
void unpack(Sparsity &e)
Reconstruct an object from the input stream.
void version(const std::string &name, int v)
void alloc_w(size_t sz_w, bool persistent=false)
Ensure required length of w field.
std::vector< std::string > get_function() const
Get a list of all functions.
const Sparsity & sparsity_out(casadi_int ind) const
Get sparsity of a given output.
'idas' plugin for Integrator
static int rhsQF(double t, N_Vector xz, N_Vector xzdot, N_Vector qdot, void *user_data)
static int jtimesB(double t, N_Vector xz, N_Vector xzdot, N_Vector xzB, N_Vector xzdotB, N_Vector resvalB, N_Vector vB, N_Vector JvB, double cjB, void *user_data, N_Vector tmp1B, N_Vector tmp2B)
int init_mem(void *mem) const override
Initalize memory block.
static int resF(double t, N_Vector xz, N_Vector xzdot, N_Vector rr, void *user_data)
void retreat(IntegratorMemory *mem, const double *u, double *adj_x, double *adj_p, double *adj_u) const override
Retreat solution in time.
void resetB(IntegratorMemory *mem) const override
Reset the backward problem and take time to tf.
void serialize_body(SerializingStream &s) const override
Serialize an object without type information.
int solve_transposed(IdasMemory *m, double t, const double *xz, const double *rxz, const double *rhs, double *sol) const
Solve transposed linear system.
void impulseB(IntegratorMemory *mem, const double *adj_x, const double *adj_z, const double *adj_q) const override
Introduce an impulse into the backwards integration at the current time.
static void idas_error(const char *module, int flag)
std::vector< double > init_xdot_
static int psetupB(double t, N_Vector xz, N_Vector xzdot, N_Vector rxz, N_Vector rxzdot, N_Vector resvalB, double cjB, void *user_dataB, N_Vector tmp1B, N_Vector tmp2B, N_Vector tmp3B)
static int lsolveB(IDAMem IDA_mem, N_Vector b, N_Vector weight, N_Vector ycur, N_Vector xzdotcur, N_Vector rescur)
static int lsetupF(IDAMem IDA_mem, N_Vector xz, N_Vector xzdot, N_Vector resp, N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
static const std::string meta_doc
A documentation string.
static int lsolveF(IDAMem IDA_mem, N_Vector b, N_Vector weight, N_Vector ycur, N_Vector xzdotcur, N_Vector rescur)
static int resB(double t, N_Vector xz, N_Vector xzdot, N_Vector rxz, N_Vector rxzdot, N_Vector rr, void *user_data)
IdasInterface(const std::string &name, const Function &dae, double t0, const std::vector< double > &tout)
Constructor.
static int rhsQB(double t, N_Vector xz, N_Vector xzdot, N_Vector rxz, N_Vector rxzdot, N_Vector ruqdot, void *user_data)
std::vector< double > abstolv_
static const Options options_
Options.
void reset(IntegratorMemory *mem, bool first_call) const override
Reset the forward solver at the start or after an event.
static int psolveF(double t, N_Vector xz, N_Vector xzdot, N_Vector rr, N_Vector rvec, N_Vector zvec, double cj, double delta, void *user_data, N_Vector tmp)
static Integrator * creator(const std::string &name, const Function &dae, double t0, const std::vector< double > &tout)
Create a new integrator.
void init(const Dict &opts) override
Initialize.
static int psetupF(double t, N_Vector xz, N_Vector xzdot, N_Vector rr, double cj, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
static void ehfun(int error_code, const char *module, const char *function, char *msg, void *eh_data)
~IdasInterface() override
Destructor.
static int lsetupB(IDAMem IDA_mem, N_Vector xz, N_Vector xzdot, N_Vector resp, N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
static IdasMemory * to_mem(void *mem)
Cast to memory object.
static int psolveB(double t, N_Vector xz, N_Vector xzdot, N_Vector rxz, N_Vector rxzdot, N_Vector resvalB, N_Vector rvecB, N_Vector zvecB, double cjB, double deltaB, void *user_dataB, N_Vector tmpB)
void z_impulseB(IdasMemory *m, const double *adj_z) const
Propagate impulse from adj_z to adj_x.
std::vector< casadi_int > y_c_
int advance_noevent(IntegratorMemory *mem) const override
Advance solution in time.
static ProtoFunction * deserialize(DeserializingStream &s)
Deserialize into MX.
static int jtimesF(double t, N_Vector xz, N_Vector xzdot, N_Vector rr, N_Vector v, N_Vector Jv, double cj, void *user_data, N_Vector tmp1, N_Vector tmp2)
casadi_int nfwd_
Number of sensitivities.
casadi_int nrx_
Number of states for the backward integration.
casadi_int nt() const
Number of output times.
static bool all_zero(const double *v, casadi_int n)
Helper function: Vector has only zeros?
std::vector< double > tout_
Output time grid.
casadi_int nu_
Number of controls.
casadi_int nx_
Number of states for the forward integration.
DM solve(const DM &A, const DM &B, bool tr=false) const
static void registerPlugin(const Plugin &plugin, bool needs_lock=true)
Register an integrator in the factory.
bool verbose_
Verbose printout.
void clear_mem()
Clear all memory (called from destructor)
Helper class for Serialization.
void version(const std::string &name, int v)
void pack(const Sparsity &e)
Serializes an object to the output stream.
casadi_int size1() const
Get the number of rows.
const casadi_int * row() const
Get a reference to row-vector,.
const casadi_int * colind() const
Get a reference to the colindex of all column element (see class description)
Linsol linsolF_
Linear solver.
enum casadi::SundialsInterface::InterpType interp_
void impulseB(IntegratorMemory *mem, const double *adj_x, const double *adj_z, const double *adj_q) const override
Introduce an impulse into the backwards integration at the current time.
double nonlin_conv_coeff_
void serialize_body(SerializingStream &s) const override
Serialize an object without type information.
enum casadi::SundialsInterface::NewtonScheme newton_scheme_
void reset(IntegratorMemory *mem, bool first_call) const override
Reset the forward solver at the start or after an event.
casadi_int max_multistep_order_
int init_mem(void *mem) const override
Initalize memory block.
void add_offsets(SundialsMemory *m) const
Add stats offsets to stats.
casadi_int max_num_steps_
bool second_order_correction_
void init(const Dict &opts) override
Initialize.
void resetB(IntegratorMemory *mem) const override
Reset the backward problem and take time to tf.
casadi_int steps_per_checkpoint_
int calc_daeB(SundialsMemory *m, double t, const double *x, const double *z, const double *adj_ode, const double *adj_alg, const double *adj_quad, double *adj_x, double *adj_z) const
static const Options options_
Options.
int CASADI_INTEGRATOR_IDAS_EXPORT casadi_register_integrator_idas(Integrator::Plugin *plugin)
void CASADI_INTEGRATOR_IDAS_EXPORT casadi_load_integrator_idas()
void casadi_copy(const T1 *x, casadi_int n, T1 *y)
COPY: y <-x.
std::string str(const T &v)
String representation, any type.
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
void casadi_copy_block(const T1 *x, const casadi_int *sp_x, T1 *y, const casadi_int *sp_y, casadi_int r_begin, casadi_int c_begin, T1 *w)
const double nan
Not a number.
void casadi_scal(casadi_int n, T1 alpha, T1 *x)
SCAL: x <- alpha*x.
void casadi_axpy(casadi_int n, T1 alpha, const T1 *x, T1 *y)
AXPY: y <- a*x + y.
void casadi_clear(T1 *x, casadi_int n)
CLEAR: x <- 0.
void * mem
Idas memory block.
double cj_last
cj used in the last factorization
IdasMemory(const IdasInterface &s)
Constructor.
Options metadata for a class.
int mem_linsolF
Linear solver memory objects.
int ncheck
number of checkpoints stored so far