idas_interface.cpp
1 /*
2  * This file is part of CasADi.
3  *
4  * CasADi -- A symbolic framework for dynamic optimization.
5  * Copyright (C) 2010-2023 Joel Andersson, Joris Gillis, Moritz Diehl,
6  * KU Leuven. All rights reserved.
7  * Copyright (C) 2011-2014 Greg Horn
8  *
9  * CasADi is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public
11  * License as published by the Free Software Foundation; either
12  * version 3 of the License, or (at your option) any later version.
13  *
14  * CasADi is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with CasADi; if not, write to the Free Software
21  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
22  *
23  */
24 
25 
26 #include "idas_interface.hpp"
27 #include "casadi/core/casadi_misc.hpp"
28 
29 // Macro for error handling
30 #define THROWING(fcn, ...) \
31 idas_error(CASADI_STR(fcn), fcn(__VA_ARGS__))
32 
33 namespace casadi {
34 
35 extern "C"
36 int CASADI_INTEGRATOR_IDAS_EXPORT
37  casadi_register_integrator_idas(Integrator::Plugin* plugin) {
38  plugin->creator = IdasInterface::creator;
39  plugin->name = "idas";
40  plugin->doc = IdasInterface::meta_doc.c_str();
41  plugin->version = CASADI_VERSION;
42  plugin->options = &IdasInterface::options_;
43  plugin->deserialize = &IdasInterface::deserialize;
44  return 0;
45 }
46 
47 extern "C"
48 void CASADI_INTEGRATOR_IDAS_EXPORT casadi_load_integrator_idas() {
50 }
51 
52 IdasInterface::IdasInterface(const std::string& name, const Function& dae,
53  double t0, const std::vector<double>& tout) : SundialsInterface(name, dae, t0, tout) {
54 }
55 
57  clear_mem();
58 }
59 
62  {{"suppress_algebraic",
63  {OT_BOOL,
64  "Suppress algebraic variables in the error testing"}},
65  {"calc_ic",
66  {OT_BOOL,
67  "Use IDACalcIC to get consistent initial conditions."}},
68  {"constraints",
69  {OT_INTVECTOR,
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."}},
72  {"calc_icB",
73  {OT_BOOL,
74  "Use IDACalcIC to get consistent initial conditions for "
75  "backwards system [default: equal to calc_ic]."}},
76  {"abstolv",
78  "Absolute tolerarance for each component"}},
79  {"max_step_size",
80  {OT_DOUBLE,
81  "Maximim step size"}},
82  {"first_time",
83  {OT_DOUBLE,
84  "First requested time as a fraction of the time interval"}},
85  {"cj_scaling",
86  {OT_BOOL,
87  "IDAS scaling on cj for the user-defined linear solver module"}},
88  {"init_xdot",
90  "Initial values for the state derivatives"}}
91  }
92 };
93 
94 void IdasInterface::init(const Dict& opts) {
95  if (verbose_) casadi_message(name_ + "::init");
96 
97  // Call the base class init
99 
100  // Default options
101  cj_scaling_ = true;
102  calc_ic_ = true;
103  suppress_algebraic_ = false;
104 
105  // Read options
106  for (auto&& op : opts) {
107  if (op.first=="init_xdot") {
108  init_xdot_ = op.second;
109  } else if (op.first=="cj_scaling") {
110  cj_scaling_ = op.second;
111  } else if (op.first=="calc_ic") {
112  calc_ic_ = op.second;
113  } else if (op.first=="suppress_algebraic") {
114  suppress_algebraic_ = op.second;
115  } else if (op.first=="constraints") {
116  y_c_ = op.second;
117  } else if (op.first=="abstolv") {
118  abstolv_ = op.second;
119  }
120  }
121 
122  // Default dependent options
124  first_time_ = tout_.back();
125 
126  // Read dependent options
127  for (auto&& op : opts) {
128  if (op.first=="calc_icB") {
129  calc_icB_ = op.second;
130  } else if (op.first=="first_time") {
131  first_time_ = op.second;
132  }
133  }
134 
135  // Get initial conditions for the state derivatives
136  if (init_xdot_.empty()) {
137  init_xdot_.resize(nx_, 0);
138  } else {
139  casadi_assert(
140  init_xdot_.size()==nx_,
141  "Option \"init_xdot\" has incorrect length. Expecting " + str(nx_) + ", "
142  "but got " + str(init_xdot_.size()) + ". "
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.");
146  }
147 
148  // Constraints
149  casadi_assert(y_c_.size() == nx_+nz_ || y_c_.empty(),
150  "Constraint vector if supplied, must be of length nx+nz, but got "
151  + str(y_c_.size()) + " and nx+nz = " + str(nx_+nz_) + ".");
152 
153  // For Jacobian calculation
154  alloc_w(nx_ + nz_); // casadi_copy_block
155 }
156 
157 int IdasInterface::resF(double t, N_Vector xz, N_Vector xzdot, N_Vector rr, void *user_data) {
158  try {
159  auto m = to_mem(user_data);
160  auto& s = m->self;
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;
163 
164  // Subtract state derivative to get residual
165  casadi_axpy(s.nx_, -1., NV_DATA_S(xzdot), NV_DATA_S(rr));
166  return 0;
167  } catch(std::exception& e) { // non-recoverable error
168  uerr() << "res failed: " << e.what() << std::endl;
169  return -1;
170  }
171 }
172 
173 void IdasInterface::ehfun(int error_code, const char *module, const char *function,
174  char *msg, void *eh_data) {
175  try {
176  //auto m = to_mem(eh_data);
177  //auto& s = m->self;
178  uerr() << msg << std::endl;
179  } catch(std::exception& e) {
180  uerr() << "ehfun failed: " << e.what() << std::endl;
181  }
182 }
183 
184 int IdasInterface::jtimesF(double t, N_Vector xz, N_Vector xzdot, N_Vector rr, N_Vector v,
185  N_Vector Jv, double cj, void *user_data, N_Vector tmp1, N_Vector tmp2) {
186  try {
187  auto m = to_mem(user_data);
188  auto& s = m->self;
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;
192 
193  // Subtract state derivative to get residual
194  casadi_axpy(s.nx_, -cj, NV_DATA_S(v), NV_DATA_S(Jv));
195 
196  return 0;
197  } catch(std::exception& e) { // non-recoverable error
198  uerr() << "jtimesF failed: " << e.what() << std::endl;
199  return -1;
200  }
201 }
202 
203 int IdasInterface::jtimesB(double t, N_Vector xz, N_Vector xzdot, N_Vector rxz,
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) {
206  try {
207  auto m = to_mem(user_data);
208  auto& s = m->self;
209  // The function is linear so the Jacobian-times-vector function is the function itself
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;
213  // Subtract state derivative to get residual
214  casadi_axpy(s.nrx_, cjB, NV_DATA_S(v), NV_DATA_S(Jv));
215 
216  return 0;
217  } catch(std::exception& e) { // non-recoverable error
218  uerr() << "jtimesB failed: " << e.what() << std::endl;
219  return -1;
220  }
221 }
222 
223 int IdasInterface::init_mem(void* mem) const {
224  if (SundialsInterface::init_mem(mem)) return 1;
225  auto m = to_mem(mem);
226 
227  // Create IDAS memory block
228  m->mem = IDACreate();
229  casadi_assert(m->mem!=nullptr, "IDACreate: Creation failed");
230 
231  // Set error handler function
232  THROWING(IDASetErrHandlerFn, m->mem, ehfun, m);
233 
234  // Set user data
235  THROWING(IDASetUserData, m->mem, m);
236 
237  // Allocate n-vectors for ivp
238  m->v_xzdot = N_VNew_Serial(nx_+nz_);
239 
240  // Initialize Idas
241  double t0 = 0;
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");
246 
247  // Include algebraic variables in error testing
248  THROWING(IDASetSuppressAlg, m->mem, suppress_algebraic_);
249 
250  // Maxinum order for the multistep method
251  THROWING(IDASetMaxOrd, m->mem, max_multistep_order_);
252 
253  // Initial step size
254  if (step0_!=0) THROWING(IDASetInitStep, m->mem, step0_);
255 
256  // Set maximum step size
257  if (max_step_size_!=0) THROWING(IDASetMaxStep, m->mem, max_step_size_);
258 
259  // Set constraints
260  if (!y_c_.empty()) {
261  N_Vector domain = N_VNew_Serial(nx_+nz_);
262  std::copy(y_c_.begin(), y_c_.end(), NV_DATA_S(domain));
263 
264  // Pass to IDA
265  int flag = IDASetConstraints(m->mem, domain);
266  casadi_assert_dev(flag==IDA_SUCCESS);
267 
268  // Free the temporary vector
269  N_VDestroy_Serial(domain);
270  }
271 
272  // Maximum order of method
273  if (max_order_) THROWING(IDASetMaxOrd, m->mem, max_order_);
274 
275  // Coeff. in the nonlinear convergence test
276  if (nonlin_conv_coeff_!=0) THROWING(IDASetNonlinConvCoef, m->mem, nonlin_conv_coeff_);
277 
278  // Scaling
279  if (!abstolv_.empty()) {
280  // Vector absolute tolerances
281  N_Vector nv_abstol = N_VNew_Serial(static_cast<long>(abstolv_.size()));
282  std::copy(abstolv_.begin(), abstolv_.end(), NV_DATA_S(nv_abstol));
283  THROWING(IDASVtolerances, m->mem, reltol_, nv_abstol);
284  N_VDestroy_Serial(nv_abstol);
285  } else if (scale_abstol_) {
286  // Scale absolute tolerances with nominal values
287  THROWING(IDASVtolerances, m->mem, reltol_, m->abstolv);
288  } else {
289  // Scalar absolute tolerances
290  THROWING(IDASStolerances, m->mem, reltol_, abstol_);
291  }
292 
293  // Maximum number of steps
294  THROWING(IDASetMaxNumSteps, m->mem, max_num_steps_);
295 
296  // Set algebraic components
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);
300 
301  // Pass this information to IDAS
302  THROWING(IDASetId, m->mem, id);
303 
304  // Delete the allocated memory
305  N_VDestroy_Serial(id);
306 
307  // attach a linear solver
308  if (newton_scheme_==SD_DIRECT) {
309  // Direct scheme
310  IDAMem IDA_mem = IDAMem(m->mem);
311  IDA_mem->ida_lmem = m;
312  IDA_mem->ida_lsetup = lsetupF;
313  IDA_mem->ida_lsolve = lsolveF;
314  IDA_mem->ida_setupNonNull = TRUE;
315  } else {
316  // Iterative scheme
317  switch (newton_scheme_) {
318  case SD_DIRECT: casadi_assert_dev(0);
319  case SD_GMRES: THROWING(IDASpgmr, m->mem, max_krylov_); break;
320  case SD_BCGSTAB: THROWING(IDASpbcg, m->mem, max_krylov_); break;
321  case SD_TFQMR: THROWING(IDASptfqmr, m->mem, max_krylov_); break;
322  }
323  THROWING(IDASpilsSetJacTimesVecFn, m->mem, jtimesF);
324  if (use_precon_) THROWING(IDASpilsSetPreconditioner, m->mem, psetupF, psolveF);
325  }
326 
327  // Quadrature equations
328  if (nq1_ > 0) {
329 
330  // Initialize quadratures in IDAS
331  THROWING(IDAQuadInit, m->mem, rhsQF, m->v_q);
332 
333  // Should the quadrature errors be used for step size control?
334  if (quad_err_con_) {
335  THROWING(IDASetQuadErrCon, m->mem, true);
336 
337  // Quadrature error tolerances
338  // TODO(Joel): vector absolute tolerances
339  THROWING(IDAQuadSStolerances, m->mem, reltol_, abstol_);
340  }
341  }
342 
343  if (verbose_) casadi_message("Attached linear solver");
344 
345  // Adjoint sensitivity problem
346  if (nadj_ > 0) {
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);
350  }
351  if (verbose_) casadi_message("Initialized adjoint sensitivities");
352 
353  // Initialize adjoint sensitivities
354  if (nadj_ > 0) {
355  int interpType = interp_==SD_HERMITE ? IDA_HERMITE : IDA_POLYNOMIAL;
356  THROWING(IDAAdjInit, m->mem, steps_per_checkpoint_, interpType);
357  }
358 
359  m->first_callB = true;
360  return 0;
361 }
362 
363 void IdasInterface::reset(IntegratorMemory* mem, bool first_call) const {
364  if (verbose_) casadi_message(name_ + "::reset");
365  auto m = to_mem(mem);
366 
367  // Reset the base classes
368  SundialsInterface::reset(mem, first_call);
369 
370  // Only reinitialize solver at first call
371  // May want to change this after more testing
372  if (first_call) {
373  // Re-initialize
374  N_VConst(0.0, m->v_xzdot);
375  std::copy(init_xdot_.begin(), init_xdot_.end(), NV_DATA_S(m->v_xzdot));
376 
377  THROWING(IDAReInit, m->mem, m->t, m->v_xz, m->v_xzdot);
378 
379  // Re-initialize quadratures
380  if (nq1_ > 0) THROWING(IDAQuadReInit, m->mem, m->v_q);
381 
382  // Correct initial conditions, if necessary
383  if (calc_ic_) {
384  THROWING(IDACalcIC, m->mem, IDA_YA_YDP_INIT , first_time_);
385  THROWING(IDAGetConsistentIC, m->mem, m->v_xz, m->v_xzdot);
386  }
387 
388  // Re-initialize backward integration
389  if (nadj_ > 0) THROWING(IDAAdjReInit, m->mem);
390  }
391 }
392 
394  auto m = to_mem(mem);
395 
396  // Do not integrate past change in input signals or past the end
397  // The event handling may cause the stop time to become smaller than internal time reached,
398  // in which case the stop time cannot be enforced
399  if (m->t_stop >= m->tcur) {
400  THROWING(IDASetStopTime, m->mem, m->t_stop);
401  }
402 
403  // Integrate, unless already at desired time
404  double ttol = 1e-9; // tolerance
405  if (fabs(m->t - m->t_next) >= ttol) {
406  // Integrate forward ...
407  double tret = m->t;
408  if (nrx_>0) { // ... with taping
409  THROWING(IDASolveF, m->mem, m->t_next, &tret, m->v_xz, m->v_xzdot, IDA_NORMAL, &m->ncheck);
410  } else { // ... without taping
411  THROWING(IDASolve, m->mem, m->t_next, &tret, m->v_xz, m->v_xzdot, IDA_NORMAL);
412  }
413  // Get quadratures
414  if (nq1_ > 0) THROWING(IDAGetQuad, m->mem, &tret, m->v_q);
415  }
416 
417  // Set function outputs
418  casadi_copy(NV_DATA_S(m->v_xz), nx_ + nz_, m->x);
419  casadi_copy(NV_DATA_S(m->v_q), nq_, m->q);
420 
421  // Get stats
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);
425 
426  return 0;
427 }
428 
430  if (verbose_) casadi_message(name_ + "::resetB");
431  auto m = to_mem(mem);
432 
433  // Reset initial guess
434  N_VConst(0.0, m->v_adj_xz);
435 
436  // Reset the base classes
438 
439  // Reset initial guess
440  N_VConst(0.0, m->v_adj_xzdot);
441 }
442 
443 void IdasInterface::z_impulseB(IdasMemory* m, const double* adj_z) const {
444  // Quick return if nothing to propagate
445  if (all_zero(adj_z, nrz_)) return;
446  // We have the following solved nonlinear system of equations:
447  // f_alg(x, z) == 0,
448  // which implicitly defines z as a function of x
449  // Linearize w.r.t. x:
450  // df_alg/dz * dz/dx + df_alg/dx == 0 <=> dz/dx = -inv(df_alg/dz)*df_alg/dx
451  // Want to calculate:
452  // adj_x = (dz/dx)^T * adj_z = -(df_alg/dx)^T * inv((df_alg/dz)^T) * adj_z
453  // = -(df_alg/dx)^T * w,
454  // where
455  // (df_alg/dz)^T * w = adj_z
456  // Augment linear system to get the system we are able to factorize
457  // [(df_ode/dx)^T - cj*I, (df_alg/dx)^T; (df_ode/dz)^T, (df_alg/dz)^T] * [v; w] = [0; adj_z]
458  // (Re)factorize linear system
459  if (psetupF(m->t, m->v_xz, m->v_xzdot, nullptr, m->cj_last, m, nullptr, nullptr, nullptr))
460  casadi_error("Linear system factorization for backwards initial conditions failed");
461  // Right-hand-side for linear system in m->tmp2
462  casadi_clear(m->tmp2, nrx_);
463  casadi_copy(adj_z, nrz_, m->tmp2 + nrx_);
464  // Solve transposed linear system of equations (Note: m->tmp2 not used since rxz null)
465  if (solve_transposed(m, m->t, NV_DATA_S(m->v_xz), nullptr, m->tmp2, m->tmp2)) {
466  casadi_error("Linear system solve for backwards initial conditions failed");
467  }
468  // Calculate: -adj_x = (df_alg/dx)^T * w
469  casadi_clear(m->tmp2, nrx_);
470  if (calc_daeB(m, m->t, NV_DATA_S(m->v_xz), NV_DATA_S(m->v_xz) + nx_,
471  m->tmp2, m->tmp2 + nrx_, nullptr, m->tmp1, m->tmp1 + nrx_)) {
472  casadi_error("Adjoint seed propagation for backwards initial conditions failed");
473  }
474  // Add contribution to backward state
475  casadi_axpy(nrx_, -1., m->tmp1, NV_DATA_S(m->v_adj_xz));
476 }
477 
479  const double* adj_x, const double* adj_z, const double* adj_q) const {
480  auto m = to_mem(mem);
481 
482  // Call method in base class
483  SundialsInterface::impulseB(mem, adj_x, adj_z, adj_q);
484 
485  // Propagate impulse from adj_z to adj_x
486  z_impulseB(m, adj_z);
487 
488  if (m->first_callB) {
489  // Create backward problem
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);
492  THROWING(IDASStolerancesB, m->mem, m->whichB, reltol_, abstol_);
493  THROWING(IDASetUserDataB, m->mem, m->whichB, m);
494  THROWING(IDASetMaxNumStepsB, m->mem, m->whichB, max_num_steps_);
495 
496 
497  // Set algebraic components
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);
503 
504  // attach linear solver
505  if (newton_scheme_==SD_DIRECT) {
506  // Direct scheme
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;
515  } else {
516  // Iterative scheme
517  switch (newton_scheme_) {
518  case SD_DIRECT: casadi_assert_dev(0);
519  case SD_GMRES: THROWING(IDASpgmrB, m->mem, m->whichB, max_krylov_); break;
520  case SD_BCGSTAB: THROWING(IDASpbcgB, m->mem, m->whichB, max_krylov_); break;
521  case SD_TFQMR: THROWING(IDASptfqmrB, m->mem, m->whichB, max_krylov_); break;
522  }
523  THROWING(IDASpilsSetJacTimesVecFnB, m->mem, m->whichB, jtimesB);
524  if (use_precon_) THROWING(IDASpilsSetPreconditionerB, m->mem, m->whichB, psetupB, psolveB);
525  }
526 
527  // Quadratures for the adjoint problem
528  THROWING(IDAQuadInitB, m->mem, m->whichB, rhsQB, m->v_adj_pu);
529  if (quad_err_con_) {
530  THROWING(IDASetQuadErrConB, m->mem, m->whichB, true);
531  THROWING(IDAQuadSStolerancesB, m->mem, m->whichB, reltol_, abstol_);
532  }
533 
534  // Mark initialized
535  m->first_callB = false;
536  } else {
537  // Re-initialize
538  THROWING(IDAReInitB, m->mem, m->whichB, m->t, m->v_adj_xz, m->v_adj_xzdot);
539  if (nrq_ > 0 || nuq_ > 0) {
540  // Workaround (bug in SUNDIALS)
541  // THROWING(IDAQuadReInitB, m->mem, m->whichB[dir], m->rq[dir]);
542  void* memB = IDAGetAdjIDABmem(m->mem, m->whichB);
543  THROWING(IDAQuadReInit, memB, m->v_adj_pu);
544  }
545  }
546 
547  // Correct initial values for the integration if necessary
548  if (calc_icB_ && m->k == nt() - 1) {
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);
551  }
552 }
553 
554 void IdasInterface::retreat(IntegratorMemory* mem, const double* u,
555  double* adj_x, double* adj_p, double* adj_u) const {
556  auto m = to_mem(mem);
557 
558  // Set controls
559  casadi_copy(u, nu_, m->u);
560 
561  // Integrate, unless already at desired time
562  if (m->t_next < m->t) {
563  double tret = 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);
566  if (nrq_ > 0 || nuq_ > 0) {
567  THROWING(IDAGetQuadB, m->mem, m->whichB, &tret, m->v_adj_pu);
568  }
569  // Interpolate to get current state
570  THROWING(IDAGetAdjY, m->mem, m->t_next, m->v_xz, m->v_xzdot);
571  }
572 
573  // Save outputs
574  casadi_copy(NV_DATA_S(m->v_adj_xz), nrx_, adj_x);
575  casadi_copy(NV_DATA_S(m->v_adj_pu), nrq_, adj_p);
576  casadi_copy(NV_DATA_S(m->v_adj_pu) + nrq_, nuq_, adj_u);
577 
578  // Get stats
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);
586 
587  // Add offset corresponding to counters that were set to zero at reinitializations
588  add_offsets(m);
589 }
590 
591 void IdasInterface::idas_error(const char* module, int flag) {
592  // Successfull return or warning
593  if (flag>=IDA_SUCCESS) return;
594  // Construct error message
595  char* flagname = IDAGetReturnFlagName(flag);
596  std::stringstream ss;
597  ss << module << " returned \"" << flagname << "\". Consult IDAS documentation.";
598  free(flagname); // NOLINT
599  casadi_error(ss.str());
600 }
601 
602 int IdasInterface::rhsQF(double t, N_Vector xz, N_Vector xzdot, N_Vector qdot, void *user_data) {
603  try {
604  auto m = to_mem(user_data);
605  auto& s = m->self;
606  if (s.calc_quadF(m, t, NV_DATA_S(xz), NV_DATA_S(xz) + s.nx_, NV_DATA_S(qdot))) return 1;
607 
608  return 0;
609  } catch(std::exception& e) { // non-recoverable error
610  uerr() << "rhsQ failed: " << e.what() << std::endl;
611  return -1;
612  }
613 }
614 
615 int IdasInterface::resB(double t, N_Vector xz, N_Vector xzdot, N_Vector rxz,
616  N_Vector rxzdot, N_Vector rr, void *user_data) {
617  try {
618  auto m = to_mem(user_data);
619  auto& s = m->self;
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;
623 
624  // Subtract state derivative to get residual
625  casadi_axpy(s.nrx_, 1., NV_DATA_S(rxzdot), NV_DATA_S(rr));
626 
627  return 0;
628  } catch(std::exception& e) { // non-recoverable error
629  uerr() << "resB failed: " << e.what() << std::endl;
630  return -1;
631  }
632 }
633 
634 int IdasInterface::rhsQB(double t, N_Vector xz, N_Vector xzdot, N_Vector rxz,
635  N_Vector rxzdot, N_Vector ruqdot, void *user_data) {
636  try {
637  auto m = to_mem(user_data);
638  auto& s = m->self;
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;
642 
643  // Negate (note definition of g)
644  casadi_scal(s.nrq_ + s.nuq_, -1., NV_DATA_S(ruqdot));
645 
646  return 0;
647  } catch(std::exception& e) { // non-recoverable error
648  uerr() << "resQB failed: " << e.what() << std::endl;
649  return -1;
650  }
651 }
652 
653 int IdasInterface::psolveF(double t, N_Vector xz, N_Vector xzdot, N_Vector rr,
654  N_Vector rvec, N_Vector zvec, double cj, double delta, void *user_data, N_Vector tmp) {
655  try {
656  auto m = to_mem(user_data);
657  auto& s = m->self;
658 
659  // Get right-hand sides in m->tmp1, ordered by sensitivity directions
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) {
664  casadi_copy(vx + d * s.nx1_, s.nx1_, v_it);
665  v_it += s.nx1_;
666  casadi_copy(vz + d * s.nz1_, s.nz1_, v_it);
667  v_it += s.nz1_;
668  }
669 
670  // Solve for undifferentiated right-hand-side, save to output
671  if (s.linsolF_.solve(m->jacF, m->tmp1, 1, false, m->mem_linsolF))
672  return 1;
673  vx = NV_DATA_S(zvec); // possibly different from rvec
674  vz = vx + s.nx_;
675  casadi_copy(m->tmp1, s.nx1_, vx);
676  casadi_copy(m->tmp1 + s.nx1_, s.nz1_, vz);
677 
678  // Sensitivity equations
679  if (s.nfwd_ > 0) {
680  // Second order correction
681  if (s.second_order_correction_) {
682  // The outputs will double as seeds for jtimesF
683  casadi_clear(vx + s.nx1_, s.nx_ - s.nx1_);
684  casadi_clear(vz + s.nz1_, s.nz_ - s.nz1_);
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;
687 
688  // Subtract m->tmp2 (reordered) from m->tmp1
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);
692  v_it += s.nx1_;
693  casadi_axpy(s.nz1_, -1., m->tmp2 + s.nx_ + d*s.nz1_, v_it);
694  v_it += s.nz1_;
695  }
696  }
697 
698  // Solve for sensitivity right-hand-sides
699  if (s.linsolF_.solve(m->jacF, m->tmp1 + s.nx1_ + s.nz1_, s.nfwd_,
700  false, m->mem_linsolF)) return 1;
701 
702  // Save to output, reordered
703  v_it = m->tmp1 + s.nx1_ + s.nz1_;
704  for (int d = 1; d <= s.nfwd_; ++d) {
705  casadi_copy(v_it, s.nx1_, vx + d * s.nx1_);
706  v_it += s.nx1_;
707  casadi_copy(v_it, s.nz1_, vz + d * s.nz1_);
708  v_it += s.nz1_;
709  }
710  }
711 
712  return 0;
713  } catch(std::exception& e) { // non-recoverable error
714  uerr() << "psolve failed: " << e.what() << std::endl;
715  return -1;
716  }
717 }
718 
719 int IdasInterface::solve_transposed(IdasMemory* m, double t, const double* xz, const double* rxz,
720  const double* rhs, double* sol) const {
721  // Get right-hand sides in m->tmp1, ordered by sensitivity directions
722  double* v_it = m->tmp1;
723  for (int d = 0; d <= nfwd_; ++d) {
724  for (int a = 0; a < nadj_; ++a) {
725  casadi_copy(rhs + (d * nadj_ + a) * nrx1_, nrx1_, v_it);
726  v_it += nrx1_;
727  casadi_copy(rhs + nrx_ + (d * nadj_ + a) * nrz1_, nrz1_, v_it);
728  v_it += nrz1_;
729  }
730  }
731 
732  // Solve for undifferentiated right-hand-side, save to output
733  if (linsolF_.solve(m->jacF, m->tmp1, nadj_, true, m->mem_linsolF)) return 1;
734  for (int a = 0; a < nadj_; ++a) {
735  casadi_copy(m->tmp1 + a * (nrx1_ + nrz1_), nrx1_, sol + a * nrx1_);
736  casadi_copy(m->tmp1 + a * (nrx1_ + nrz1_) + nrx1_, nrz1_, sol + nrx_ + a * nrz1_);
737  }
738 
739  // Sensitivity equations
740  if (nfwd_ > 0) {
741  // Second order correction
742  if (second_order_correction_ && rxz) {
743  // The outputs will double as seeds for calc_daeB
744  casadi_clear(sol + nrx1_ * nadj_, nrx_ - nrx1_ * nadj_);
745  casadi_clear(sol + nrx_ + nrz1_ * nadj_, nrz_ - nrz1_ * nadj_);
746 
747  // Get second-order-correction, save to m->tmp2
748  if (calc_daeB(m, t, xz, xz + nx_, sol, sol + nrx_, nullptr,
749  m->tmp2, m->tmp2 + nrx_)) return 1;
750 
751  // Subtract m->tmp2 (reordered) from m->tmp1
752  v_it = m->tmp1 + (nrx1_ + nrz1_) * nadj_;
753  for (int d = 1; d <= nfwd_; ++d) {
754  for (int a = 0; a < nadj_; ++a) {
755  casadi_axpy(nrx1_, -1., m->tmp2 + nrx1_ * (d * nadj_ + a), v_it);
756  v_it += nrx1_;
757  casadi_axpy(nrz1_, -1., m->tmp2 + nrx_ + nrz1_ * (d * nadj_ + a), v_it);
758  v_it += nrz1_;
759  }
760  }
761  }
762 
763  // Solve for sensitivity right-hand-sides
764  if (linsolF_.solve(m->jacF, m->tmp1 + nrx1_ * nadj_ + nrz1_ * nadj_,
765  nadj_ * nfwd_, true, m->mem_linsolF)) return 1;
766 
767  // Save to output, reordered
768  v_it = m->tmp1 + (nrx1_ + nrz1_) * nadj_;
769  for (int d = 1; d <= nfwd_; ++d) {
770  for (int a = 0; a < nadj_; ++a) {
771  casadi_copy(v_it, nrx1_, sol + nrx1_ * (d * nadj_ + a));
772  v_it += nrx1_;
773  casadi_copy(v_it, nrz1_, sol + nrx_ + nrz1_ * (d * nadj_ + a));
774  v_it += nrz1_;
775  }
776  }
777  }
778 
779  return 0;
780 }
781 
782 int IdasInterface::psolveB(double t, N_Vector xz, N_Vector xzdot, N_Vector xzB,
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) {
786  try {
787  auto m = to_mem(user_data);
788  auto& s = m->self;
789  return s.solve_transposed(m, t, NV_DATA_S(xz), NV_DATA_S(xzB),
790  NV_DATA_S(rvecB), NV_DATA_S(zvecB));
791 
792  } catch(std::exception& e) { // non-recoverable error
793  uerr() << "psolveB failed: " << e.what() << std::endl;
794  return -1;
795  }
796 }
797 
798 template<typename T1>
799 void casadi_copy_block(const T1* x, const casadi_int* sp_x, T1* y, const casadi_int* sp_y,
800  casadi_int r_begin, casadi_int c_begin, T1* w) {
801  // x and y should be distinct
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;
804  nrow_x = sp_x[0];
805  ncol_x = sp_x[1];
806  colind_x = sp_x+2; row_x = sp_x + 2 + ncol_x+1;
807  ncol_y = sp_y[1];
808  colind_y = sp_y+2; row_y = sp_y + 2 + ncol_y+1;
809  // End of the rows to be copied
810  r_end = r_begin + nrow_x;
811  // w will correspond to a column of x, initialize to zero
812  casadi_clear(w, nrow_x);
813  // Loop over columns in x
814  for (i_x = 0; i_x < ncol_x; ++i_x) {
815  // Corresponding row in y
816  i_y = i_x + c_begin;
817  // Copy x entries to w
818  for (el=colind_x[i_x]; el<colind_x[i_x + 1]; ++el) w[row_x[el]] = x[el];
819  // Copy entries to y, if in interval
820  for (el=colind_y[i_y]; el<colind_y[i_y + 1]; ++el) {
821  j = row_y[el];
822  if (j >= r_begin && j < r_end) y[el] = w[j - r_begin];
823  }
824  // Restore w
825  for (el=colind_x[i_x]; el<colind_x[i_x + 1]; ++el) w[row_x[el]] = 0;
826  }
827 }
828 
829 int IdasInterface::psetupF(double t, N_Vector xz, N_Vector xzdot, N_Vector rr,
830  double cj, void* user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3) {
831  try {
832  auto m = to_mem(user_data);
833  auto& s = m->self;
834 
835  // Sparsity patterns
836  const Function& jacF = s.get_function("jacF");
837  const Sparsity& sp_jac_ode_x = jacF.sparsity_out(JACF_ODE_X);
838  const Sparsity& sp_jac_alg_x = jacF.sparsity_out(JACF_ALG_X);
839  const Sparsity& sp_jac_ode_z = jacF.sparsity_out(JACF_ODE_Z);
840  const Sparsity& sp_jac_alg_z = jacF.sparsity_out(JACF_ALG_Z);
841  const Sparsity& sp_jacF = s.linsolF_.sparsity();
842 
843  // Calculate Jacobian blocks
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;
846 
847  // Copy to jacF structure
848  casadi_int nx_jac = sp_jac_ode_x.size1(); // excludes sensitivity equations
849  casadi_copy_block(m->jac_ode_x, sp_jac_ode_x, m->jacF, sp_jacF, 0, 0, m->w);
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);
853 
854  // Shift diagonal corresponding to jac_ode_x
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;
859  }
860  }
861 
862  // Factorize the linear system
863  if (s.linsolF_.nfact(m->jacF, m->mem_linsolF)) return 1;
864  m->cj_last = cj;
865 
866  return 0;
867  } catch(std::exception& e) { // non-recoverable error
868  uerr() << "psetup failed: " << e.what() << std::endl;
869  return -1;
870  }
871 }
872 
873 int IdasInterface::psetupB(double t, N_Vector xz, N_Vector xzdot, N_Vector rxz, N_Vector rxzdot,
874  N_Vector rresval, double cj, void *user_data, N_Vector tmp1B, N_Vector tmp2B, N_Vector tmp3B) {
875  try {
876  // We use the same linear solver for the forward problem as for the backward problem
877  return psetupF(t, xz, nullptr, nullptr, -cj, user_data, tmp1B, tmp2B, tmp3B);
878 
879  } catch(std::exception& e) { // non-recoverable error
880  uerr() << "psetupB failed: " << e.what() << std::endl;
881  return -1;
882  }
883 }
884 
885 int IdasInterface::lsetupF(IDAMem IDA_mem, N_Vector xz, N_Vector xzdot, N_Vector resp,
886  N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3) {
887  // Current time
888  double t = IDA_mem->ida_tn;
889 
890  // Multiple of df_dydot to be added to the matrix
891  double cj = IDA_mem->ida_cj;
892 
893  // Call the preconditioner setup function (which sets up the linear solver)
894  return psetupF(t, xz, xzdot, nullptr, cj, IDA_mem->ida_lmem,
895  vtemp1, vtemp1, vtemp3);
896 }
897 
898 int IdasInterface::lsetupB(IDAMem IDA_mem, N_Vector xzB, N_Vector xzdotB, N_Vector respB,
899  N_Vector vtemp1B, N_Vector vtemp2B, N_Vector vtemp3B) {
900  try {
901  auto m = to_mem(IDA_mem->ida_lmem);
902  //auto& s = m->self;
903  IDAadjMem IDAADJ_mem;
904  //IDABMem IDAB_mem;
905 
906  // Current time
907  double t = IDA_mem->ida_tn; // TODO(Joel): is this correct?
908  // Multiple of df_dydot to be added to the matrix
909  double cj = IDA_mem->ida_cj;
910 
911  IDA_mem = static_cast<IDAMem>(IDA_mem->ida_user_data);
912 
913  IDAADJ_mem = IDA_mem->ida_adj_mem;
914  //IDAB_mem = IDAADJ_mem->ia_bckpbCrt;
915 
916  // Get FORWARD solution from interpolation.
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,
919  nullptr, nullptr);
920  if (flag != IDA_SUCCESS) casadi_error("Could not interpolate forward states");
921  }
922  // Call the preconditioner setup function (which sets up the linear solver)
923  return psetupB(t, IDAADJ_mem->ia_yyTmp, IDAADJ_mem->ia_ypTmp,
924  xzB, xzdotB, nullptr, cj, static_cast<void*>(m), vtemp1B, vtemp1B, vtemp3B);
925 
926  } catch(std::exception& e) { // non-recoverable error
927  uerr() << "lsetupB failed: " << e.what() << std::endl;
928  return -1;
929  }
930 }
931 
932 int IdasInterface::lsolveF(IDAMem IDA_mem, N_Vector b, N_Vector weight, N_Vector xz,
933  N_Vector xzdot, N_Vector rr) {
934  try {
935  auto m = to_mem(IDA_mem->ida_lmem);
936  auto& s = m->self;
937 
938  // Current time
939  double t = IDA_mem->ida_tn;
940 
941  // Multiple of df_dydot to be added to the matrix
942  double cj = IDA_mem->ida_cj;
943 
944  // Accuracy
945  double delta = 0.0;
946 
947  // Call the preconditioner solve function (which solves the linear system)
948  int flag = psolveF(t, xz, xzdot, rr, b, b, cj,
949  delta, static_cast<void*>(m), nullptr);
950  if (flag) return flag;
951 
952  // Scale the correction to account for change in cj
953  if (s.cj_scaling_) {
954  double cjratio = IDA_mem->ida_cjratio;
955  if (cjratio != 1.0) N_VScale(2.0/(1.0 + cjratio), b, b);
956  }
957 
958  return 0;
959  } catch(std::exception& e) { // non-recoverable error
960  uerr() << "lsolve failed: " << e.what() << std::endl;
961  return -1;
962  }
963 }
964 
965 int IdasInterface::lsolveB(IDAMem IDA_mem, N_Vector b, N_Vector weight, N_Vector xzB,
966  N_Vector xzdotB, N_Vector rrB) {
967  try {
968  auto m = to_mem(IDA_mem->ida_lmem);
969  auto& s = m->self;
970  IDAadjMem IDAADJ_mem;
971  //IDABMem IDAB_mem;
972  int flag;
973 
974  // Current time
975  double t = IDA_mem->ida_tn; // TODO(Joel): is this correct?
976  // Multiple of df_dydot to be added to the matrix
977  double cj = IDA_mem->ida_cj;
978  double cjratio = IDA_mem->ida_cjratio;
979 
980  IDA_mem = (IDAMem) IDA_mem->ida_user_data;
981  IDAADJ_mem = IDA_mem->ida_adj_mem;
982  //IDAB_mem = IDAADJ_mem->ia_bckpbCrt;
983 
984  // Get FORWARD solution from interpolation.
985  if (IDAADJ_mem->ia_noInterp==FALSE) {
986  flag = IDAADJ_mem->ia_getY(IDA_mem, t, IDAADJ_mem->ia_yyTmp, IDAADJ_mem->ia_ypTmp,
987  nullptr, nullptr);
988  if (flag != IDA_SUCCESS) casadi_error("Could not interpolate forward states");
989  }
990 
991  // Accuracy
992  double delta = 0.0;
993 
994  // Call the preconditioner solve function (which solves the linear system)
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;
998 
999  // Scale the correction to account for change in cj
1000  if (s.cj_scaling_) {
1001  if (cjratio != 1.0) N_VScale(2.0/(1.0 + cjratio), b, b);
1002  }
1003  return 0;
1004  } catch(std::exception& e) { // non-recoverable error
1005  uerr() << "lsolveB failed: " << e.what() << std::endl;
1006  return -1;
1007  }
1008 }
1009 
1011  this->mem = nullptr;
1012  this->v_xzdot = nullptr;
1013  this->v_adj_xzdot = nullptr;
1014  this->cj_last = nan;
1015 
1016  // Reset checkpoints counter
1017  this->ncheck = 0;
1018 }
1019 
1021  if (this->mem) IDAFree(&this->mem);
1022  if (this->v_xzdot) N_VDestroy_Serial(this->v_xzdot);
1023  if (this->v_adj_xzdot) N_VDestroy_Serial(this->v_adj_xzdot);
1024  if (this->mem_linsolF >= 0) self.linsolF_.release(this->mem_linsolF);
1025 }
1026 
1028  int version = s.version("IdasInterface", 1, 2);
1029  s.unpack("IdasInterface::cj_scaling", cj_scaling_);
1030  s.unpack("IdasInterface::calc_ic", calc_ic_);
1031  s.unpack("IdasInterface::calc_icB", calc_icB_);
1032  s.unpack("IdasInterface::suppress_algebraic", suppress_algebraic_);
1033  s.unpack("IdasInterface::abstolv", abstolv_);
1034  s.unpack("IdasInterface::first_time", first_time_);
1035  s.unpack("IdasInterface::init_xdot", init_xdot_);
1036 
1037  if (version>=2) {
1038  s.unpack("IdasInterface::max_step_size", max_step_size_);
1039  s.unpack("IdasInterface::y_c", y_c_);
1040  } else {
1041  max_step_size_ = 0;
1042  }
1043 }
1044 
1047  s.version("IdasInterface", 2);
1048  s.pack("IdasInterface::cj_scaling", cj_scaling_);
1049  s.pack("IdasInterface::calc_ic", calc_ic_);
1050  s.pack("IdasInterface::calc_icB", calc_icB_);
1051  s.pack("IdasInterface::suppress_algebraic", suppress_algebraic_);
1052  s.pack("IdasInterface::abstolv", abstolv_);
1053  s.pack("IdasInterface::first_time", first_time_);
1054  s.pack("IdasInterface::init_xdot", init_xdot_);
1055  s.pack("IdasInterface::max_step_size", max_step_size_);
1056  s.pack("IdasInterface::y_c", y_c_);
1057 }
1058 
1059 } // namespace casadi
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.
Function object.
Definition: function.hpp:60
std::vector< std::string > get_function() const
Get a list of all functions.
Definition: function.cpp:1844
const Sparsity & sparsity_out(casadi_int ind) const
Get sparsity of a given output.
Definition: function.cpp:1031
'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.
double t0_
Initial time.
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
Definition: linsol.cpp:73
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.
General sparsity class.
Definition: sparsity.hpp:106
casadi_int size1() const
Get the number of rows.
Definition: sparsity.cpp:124
const casadi_int * row() const
Get a reference to row-vector,.
Definition: sparsity.cpp:164
const casadi_int * colind() const
Get a reference to the colindex of all column element (see class description)
Definition: sparsity.cpp:168
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.
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.
int init_mem(void *mem) const override
Initalize memory block.
void add_offsets(SundialsMemory *m) const
Add stats offsets to stats.
void init(const Dict &opts) override
Initialize.
void resetB(IntegratorMemory *mem) const override
Reset the backward problem and take time to tf.
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.
The casadi namespace.
Definition: archiver.cpp:28
int CASADI_INTEGRATOR_IDAS_EXPORT casadi_register_integrator_idas(Integrator::Plugin *plugin)
std::ostream & uerr()
void CASADI_INTEGRATOR_IDAS_EXPORT casadi_load_integrator_idas()
void casadi_copy(const T1 *x, casadi_int n, T1 *y)
COPY: y <-x.
@ OT_INTVECTOR
@ OT_DOUBLEVECTOR
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.
Definition: calculus.hpp:53
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.
~IdasMemory()
Destructor.
double cj_last
cj used in the last factorization
IdasMemory(const IdasInterface &s)
Constructor.
Options metadata for a class.
Definition: options.hpp:40
int mem_linsolF
Linear solver memory objects.
int ncheck
number of checkpoints stored so far