26 #ifndef CASADI_SUNDIALS_INTERFACE_HPP
27 #define CASADI_SUNDIALS_INTERFACE_HPP
29 #include <casadi/interfaces/sundials/casadi_sundials_common_export.h>
30 #include "casadi/core/integrator_impl.hpp"
32 #include <nvector/nvector_serial.h>
33 #include <sundials/sundials_dense.h>
34 #include <sundials/sundials_iterative.h>
35 #include <sundials/sundials_types.h>
45 N_Vector v_xz, v_xzdot,
v_q;
54 double *jac_ode_x, *
jac_alg_x, *jac_ode_z, *jac_alg_z;
62 double hinused, hlast,
hcur, tcur;
68 double hinusedB, hlastB,
hcurB, tcurB;
95 double t0,
const std::vector<double>& tout);
107 void init(
const Dict& opts)
override;
110 void set_work(
void* mem,
const double**& arg,
double**& res,
111 casadi_int*& iw,
double*& w)
const override;
114 int init_mem(
void* mem)
const override;
123 int calc_daeF(
SundialsMemory* m,
double t,
const double* x,
const double* z,
124 double* ode,
double* alg)
const;
127 int calc_daeB(
SundialsMemory* m,
double t,
const double* x,
const double* z,
128 const double* adj_ode,
const double* adj_alg,
const double* adj_quad,
129 double* adj_x,
double* adj_z)
const;
132 int calc_quadF(
SundialsMemory* m,
double t,
const double* x,
const double* z,
136 int calc_quadB(
SundialsMemory* m,
double t,
const double* x,
const double* z,
137 const double* adj_ode,
const double* adj_alg,
double* adj_p,
double* adj_u)
const;
140 int calc_jtimesF(
SundialsMemory* m,
double t,
const double* x,
const double* z,
141 const double* fwd_x,
const double* fwd_z,
double* fwd_ode,
double* fwd_alg)
const;
144 int calc_jacF(
SundialsMemory* m,
double t,
const double* x,
const double* z,
145 double* jac_ode_x,
double* jac_alg_x,
double* jac_ode_z,
double* jac_alg_z)
const;
148 Dict get_stats(
void* mem)
const override;
161 const double* adj_x,
const double* adj_z,
const double* adj_q)
const override;
175 casadi_assert_dev(m);
212 enum NewtonScheme {SD_DIRECT, SD_GMRES, SD_BCGSTAB, SD_TFQMR} newton_scheme_;
221 static void printvar(
const std::string&
id,
double v) {
222 uout() <<
id <<
" = " << v << std::endl;
225 static void printvar(
const std::string&
id, N_Vector v) {
226 std::vector<double> tmp(NV_DATA_S(v), NV_DATA_S(v)+NV_LENGTH_S(v));
227 uout() <<
id <<
" = " << tmp << std::endl;
240 std::vector<double> tmp(NV_DATA_S(v), NV_DATA_S(v)+NV_LENGTH_S(v));
Helper class for Serialization.
Internal storage for integrator related data.
Helper class for Serialization.
Linsol linsolF_
Linear solver.
JacFOut
IO conventions for continuous time dynamics.
double get_abstol() const override
Get absolute tolerance.
static SundialsMemory * to_mem(void *mem)
Cast to memory object.
double get_reltol() const override
Get relative tolerance.
std::string linear_solver_
NewtonScheme
Supported iterative solvers in Sundials.
double nonlin_conv_coeff_
Dict linear_solver_options_
bool disable_internal_warnings_
casadi_int max_multistep_order_
casadi_int max_num_steps_
JtimesFIn
IO conventions for continuous time dynamics.
bool second_order_correction_
JtimesFOut
IO conventions for continuous time dynamics.
casadi_int steps_per_checkpoint_
static void printvar(const std::string &id, double v)
static const Options options_
Options.
static void printvar(const std::string &id, N_Vector v)
const Options & get_options() const override
Options.
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
bool is_regular(const std::vector< T > &v)
Checks if array does not contain NaN or Inf.
Options metadata for a class.
Linear solver data (dense) – what is this?
int mem_linsolF
Linear solver memory objects.
int ncheck
number of checkpoints stored so far