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 SundialsMemory : public IntegratorMemory {
44  // N-vectors for the forward integration
45  N_Vector xz, xzdot, q;
46 
47  // N-vectors for the backward integration
48  N_Vector rxz, rxzdot, ruq;
49 
50  // Initialize or reinitialize?
51  bool first_callB;
52 
53  // Parameters
54  double *p, *rp;
55 
56  // Controls
57  double *u;
58 
60  double *jac_ode_x, *jac_alg_x, *jac_ode_z, *jac_alg_z;
61 
62  // Jacobian
63  double *jacF;
64 
66  long nsteps, nfevals, nlinsetups, netfails;
67  int qlast, qcur;
68  double hinused, hlast, hcur, tcur;
69  long nniters, nncfails;
70 
72  long nstepsB, nfevalsB, nlinsetupsB, netfailsB;
73  int qlastB, qcurB;
74  double hinusedB, hlastB, hcurB, tcurB;
75  long nnitersB, nncfailsB;
76 
78  long nstepsB_off, nfevalsB_off, nlinsetupsB_off, netfailsB_off;
79  long nnitersB_off, nncfailsB_off;
80 
81  // Temporaries for [x;z] or [rx;rz]
82  double *v1, *v2;
83 
85  int ncheck;
86 
87  // Absolute tolerance
88  N_Vector abstolv;
89 
91  int mem_linsolF;
92 
94  SundialsMemory();
95 
97  ~SundialsMemory();
98  };
99 
100  class SundialsInterface : public Integrator {
101  public:
103  SundialsInterface(const std::string& name, const Function& dae,
104  double t0, const std::vector<double>& tout);
105 
107  ~SundialsInterface() override=0;
108 
110 
111  static const Options options_;
112  const Options& get_options() const override { return options_;}
114 
116  void init(const Dict& opts) override;
117 
119  void set_work(void* mem, const double**& arg, double**& res,
120  casadi_int*& iw, double*& w) const override;
121 
123  int init_mem(void* mem) const override;
124 
126  double get_reltol() const override { return reltol_;}
127 
129  double get_abstol() const override { return abstol_;}
130 
131  // DAE right-hand-side, forward problem
132  int calc_daeF(SundialsMemory* m, double t, const double* x, const double* z,
133  double* ode, double* alg) const;
134 
135  // DAE right-hand-side, forward problem
136  int calc_daeB(SundialsMemory* m, double t, const double* x, const double* z,
137  const double* rx, const double* rz, const double* rp, double* adj_x, double* adj_z) const;
138 
139  // Quadrature right-hand-side, forward problem
140  int calc_quadF(SundialsMemory* m, double t, const double* x, const double* z,
141  double* quad) const;
142 
143  // Quadrature right-hand-side, backward problem
144  int calc_quadB(SundialsMemory* m, double t, const double* x, const double* z,
145  const double* rx, const double* rz, double* adj_p, double* adj_u) const;
146 
147  // Jacobian of DAE-times-vector function, forward problem
148  int calc_jtimesF(SundialsMemory* m, double t, const double* x, const double* z,
149  const double* fwd_x, const double* fwd_z, double* fwd_ode, double* fwd_alg) const;
150 
151  // Jacobian of DAE right-hand-side function, forward problem
152  int calc_jacF(SundialsMemory* m, double t, const double* x, const double* z,
153  double* jac_ode_x, double* jac_alg_x, double* jac_ode_z, double* jac_alg_z) const;
154 
156  Dict get_stats(void* mem) const override;
157 
159  void print_stats(IntegratorMemory* mem) const override;
160 
162  void reset(IntegratorMemory* mem, const double* u, const double* x,
163  const double* z, const double* p) const override;
164 
166  void resetB(IntegratorMemory* mem) const override;
167 
169  void impulseB(IntegratorMemory* mem,
170  const double* rx, const double* rz, const double* rp) const override;
171 
173  void reset_stats(SundialsMemory* m) const;
174 
176  void save_offsets(SundialsMemory* m) const;
177 
179  void add_offsets(SundialsMemory* m) const;
180 
182  static SundialsMemory* to_mem(void *mem) {
183  SundialsMemory* m = static_cast<SundialsMemory*>(mem);
184  casadi_assert_dev(m);
185  return m;
186  }
187 
189 
190  enum JtimesFIn { JTIMESF_T, JTIMESF_X, JTIMESF_Z, JTIMESF_P, JTIMESF_U, JTIMESF_FWD_X,
191  JTIMESF_FWD_Z, JTIMESF_NUM_IN};
192  enum JtimesFOut { JTIMESF_FWD_ODE, JTIMESF_FWD_ALG, JTIMESF_NUM_OUT};
193  enum JacFOut {JACF_ODE_X, JACF_ALG_X, JACF_ODE_Z, JACF_ALG_Z, JACF_NUM_OUT};
195 
198  double abstol_, reltol_;
199  casadi_int max_num_steps_;
200  bool stop_at_end_;
201  bool quad_err_con_;
202  casadi_int steps_per_checkpoint_;
203  bool disable_internal_warnings_;
204  casadi_int max_multistep_order_;
205  std::string linear_solver_;
206  Dict linear_solver_options_;
207  casadi_int max_krylov_;
208  bool use_precon_;
209  bool second_order_correction_;
210  double step0_;
211  double max_step_size_;
212  double nonlin_conv_coeff_;
213  casadi_int max_order_;
214  bool scale_abstol_;
216 
218  Linsol linsolF_;
219 
221  enum NewtonScheme {SD_DIRECT, SD_GMRES, SD_BCGSTAB, SD_TFQMR} newton_scheme_;
222 
223  // Supported interpolations in Sundials
224  enum InterpType {SD_POLYNOMIAL, SD_HERMITE} interp_;
225 
227  struct LinSolDataDense {};
228 
229  // Print a variable
230  static void printvar(const std::string& id, double v) {
231  uout() << id << " = " << v << std::endl;
232  }
233  // Print an N_Vector
234  static void printvar(const std::string& id, N_Vector v) {
235  std::vector<double> tmp(NV_DATA_S(v), NV_DATA_S(v)+NV_LENGTH_S(v));
236  uout() << id << " = " << tmp << std::endl;
237  }
238 
240  void serialize_body(SerializingStream &s) const override;
241 
242  protected:
244  explicit SundialsInterface(DeserializingStream& s);
245  };
246 
247  // Check if N_Vector is regular
248  inline bool is_regular(N_Vector v) {
249  std::vector<double> tmp(NV_DATA_S(v), NV_DATA_S(v)+NV_LENGTH_S(v));
250  return is_regular(tmp);
251  }
252 
253 } // namespace casadi
254 
256 #endif // CASADI_SUNDIALS_INTERFACE_HPP
The casadi namespace.
CASADI_EXPORT std::ostream & uout()
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.