sundials_interface.hpp
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 #ifndef CASADI_SUNDIALS_INTERFACE_HPP
27 #define CASADI_SUNDIALS_INTERFACE_HPP
28 
29 #include <casadi/interfaces/sundials/casadi_sundials_common_export.h>
30 #include "casadi/core/integrator_impl.hpp"
31 
32 #include <nvector/nvector_serial.h>
33 #include <sundials/sundials_dense.h>
34 #include <sundials/sundials_iterative.h>
35 #include <sundials/sundials_types.h>
36 
37 #include <ctime>
38 
40 namespace casadi {
41 
42  // IdasMemory
43  struct CASADI_SUNDIALS_COMMON_EXPORT SundialsMemory : public IntegratorMemory {
44  // N-vectors for the forward integration
45  N_Vector v_xz, v_xzdot, v_q;
46 
47  // N-vectors for the backward integration
48  N_Vector v_adj_xz, v_adj_xzdot, v_adj_pu;
49 
50  // Initialize or reinitialize?
52 
54  double *jac_ode_x, *jac_alg_x, *jac_ode_z, *jac_alg_z;
55 
56  // Jacobian
57  double *jacF;
58 
60  long nsteps, nfevals, nlinsetups, netfails;
61  int qlast, qcur;
62  double hinused, hlast, hcur, tcur;
63  long nniters, nncfails;
64 
66  long nstepsB, nfevalsB, nlinsetupsB, netfailsB;
67  int qlastB, qcurB;
68  double hinusedB, hlastB, hcurB, tcurB;
69  long nnitersB, nncfailsB;
70 
72  long nstepsB_off, nfevalsB_off, nlinsetupsB_off, netfailsB_off;
73  long nnitersB_off, nncfailsB_off;
74 
76  int ncheck;
77 
78  // Absolute tolerance
79  N_Vector abstolv;
80 
83 
86 
88  ~SundialsMemory();
89  };
90 
91  class CASADI_SUNDIALS_COMMON_EXPORT SundialsInterface : public Integrator {
92  public:
94  SundialsInterface(const std::string& name, const Function& dae,
95  double t0, const std::vector<double>& tout);
96 
98  ~SundialsInterface() override=0;
99 
101 
102  static const Options options_;
103  const Options& get_options() const override { return options_;}
105 
107  void init(const Dict& opts) override;
108 
110  void set_work(void* mem, const double**& arg, double**& res,
111  casadi_int*& iw, double*& w) const override;
112 
114  int init_mem(void* mem) const override;
115 
117  double get_reltol() const override { return reltol_;}
118 
120  double get_abstol() const override { return abstol_;}
121 
122  // DAE right-hand-side, forward problem
123  int calc_daeF(SundialsMemory* m, double t, const double* x, const double* z,
124  double* ode, double* alg) const;
125 
126  // DAE right-hand-side, forward problem
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;
130 
131  // Quadrature right-hand-side, forward problem
132  int calc_quadF(SundialsMemory* m, double t, const double* x, const double* z,
133  double* quad) const;
134 
135  // Quadrature right-hand-side, backward problem
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;
138 
139  // Jacobian of DAE-times-vector function, forward problem
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;
142 
143  // Jacobian of DAE right-hand-side function, forward problem
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;
146 
148  Dict get_stats(void* mem) const override;
149 
151  void print_stats(IntegratorMemory* mem) const override;
152 
154  void reset(IntegratorMemory* mem, bool first_call) const override;
155 
157  void resetB(IntegratorMemory* mem) const override;
158 
160  void impulseB(IntegratorMemory* mem,
161  const double* adj_x, const double* adj_z, const double* adj_q) const override;
162 
164  void reset_stats(SundialsMemory* m) const;
165 
167  void save_offsets(SundialsMemory* m) const;
168 
170  void add_offsets(SundialsMemory* m) const;
171 
173  static SundialsMemory* to_mem(void *mem) {
174  SundialsMemory* m = static_cast<SundialsMemory*>(mem);
175  casadi_assert_dev(m);
176  return m;
177  }
178 
180 
181  enum JtimesFIn { JTIMESF_T, JTIMESF_X, JTIMESF_Z, JTIMESF_P, JTIMESF_U, JTIMESF_FWD_X,
182  JTIMESF_FWD_Z, JTIMESF_NUM_IN};
183  enum JtimesFOut { JTIMESF_FWD_ODE, JTIMESF_FWD_ALG, JTIMESF_NUM_OUT};
184  enum JacFOut {JACF_ODE_X, JACF_ALG_X, JACF_ODE_Z, JACF_ALG_Z, JACF_NUM_OUT};
186 
189  double abstol_, reltol_;
190  casadi_int max_num_steps_;
196  std::string linear_solver_;
198  casadi_int max_krylov_;
201  double step0_;
204  casadi_int max_order_;
207 
210 
212  enum NewtonScheme {SD_DIRECT, SD_GMRES, SD_BCGSTAB, SD_TFQMR} newton_scheme_;
213 
214  // Supported interpolations in Sundials
215  enum InterpType {SD_POLYNOMIAL, SD_HERMITE} interp_;
216 
218  struct LinSolDataDense {};
219 
220  // Print a variable
221  static void printvar(const std::string& id, double v) {
222  uout() << id << " = " << v << std::endl;
223  }
224  // Print an N_Vector
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;
228  }
229 
231  void serialize_body(SerializingStream &s) const override;
232 
233  protected:
236  };
237 
238  // Check if N_Vector is regular
239  inline bool is_regular(N_Vector v) {
240  std::vector<double> tmp(NV_DATA_S(v), NV_DATA_S(v)+NV_LENGTH_S(v));
241  return is_regular(tmp);
242  }
243 
244 } // namespace casadi
245 
247 #endif // CASADI_SUNDIALS_INTERFACE_HPP
Helper class for Serialization.
Function object.
Definition: function.hpp:60
Internal storage for integrator related data.
Linear solver.
Definition: linsol.hpp:55
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.
NewtonScheme
Supported iterative solvers in Sundials.
JtimesFIn
IO conventions for continuous time dynamics.
JtimesFOut
IO conventions for continuous time dynamics.
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.
The casadi namespace.
Definition: archiver.cpp:28
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.
std::ostream & uout()
Options metadata for a class.
Definition: options.hpp:40
Linear solver data (dense) – what is this?
int mem_linsolF
Linear solver memory objects.
int ncheck
number of checkpoints stored so far