kinsol_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_KINSOL_INTERFACE_HPP
27 #define CASADI_KINSOL_INTERFACE_HPP
28 
29 #include <casadi/interfaces/sundials/casadi_rootfinder_kinsol_export.h>
30 #include "casadi/core/rootfinder_impl.hpp"
31 #include <nvector/nvector_serial.h> /* serial N_Vector types, fcts., and macros */
32 #include <sundials/sundials_dense.h> /* definitions DlsMat DENSE_ELEM */
33 #include <sundials/sundials_types.h> /* definition of type double */
34 #include <kinsol/kinsol.h> /* prototypes for CVode fcts. and consts. */
35 #include <kinsol/kinsol_dense.h>
36 #include <kinsol/kinsol_band.h>
37 #include <kinsol/kinsol_spgmr.h>
38 #include <kinsol/kinsol_spbcgs.h>
39 #include <kinsol/kinsol_sptfqmr.h>
40 #include <kinsol/kinsol_impl.h> /* Needed for the provided linear solver */
41 #include <kinsol/kinsol_spils_impl.h> /* Needed for the provided linear solver */
42 #include <ctime>
43 
53 namespace casadi {
55  class KinsolInterface;
56 
57  // Memory
58  struct KinsolMemory : public RootfinderMemory {
60  const KinsolInterface& self;
61 
63  KinsolMemory(const KinsolInterface& s);
64 
66  ~KinsolMemory();
67 
69  void* mem;
70 
72  N_Vector u;
73 
74  // Current Jacobian
75  double* jac;
76  };
77 
83  class KinsolInterface : public Rootfinder {
84  public:
86  explicit KinsolInterface(const std::string& name, const Function& f);
87 
89  ~KinsolInterface() override;
90 
92  static Rootfinder* creator(const std::string& name, const Function& f) {
93  return new KinsolInterface(name, f);
94  }
95 
97 
98  static const Options options_;
99  const Options& get_options() const override { return options_;}
101 
103  void init(const Dict& opts) override;
104 
106  int solve(void* mem) const override;
107 
108  // Get name of the plugin
109  const char* plugin_name() const override { return "kinsol";}
110 
111  // Get name of the class
112  std::string class_name() const override { return "KinsolInterface";}
113 
114  // Scaling
115  N_Vector u_scale_, f_scale_;
116 
118  casadi_int strategy_;
119 
120  // Should KINSOL internal warning messages be ignored
121  bool disable_internal_warnings_;
122 
123  // Maximum number of iterations
124  casadi_int max_iter_;
125 
126  // Print information about iterations
127  casadi_int print_level_;
128 
129  // Use exact Jacobian?
130  bool exact_jac_;
131 
132  // Type of linear solver
133  enum LinsolType { DENSE, BANDED, ITERATIVE, USER_DEFINED};
134  LinsolType linear_solver_type_;
135 
136  // Bandwidth (for banded solvers)
137  casadi_int upper_bandwidth_, lower_bandwidth_;
138 
139  // Krylov subspace size (for iterative solvers)
140  casadi_int maxl_;
141 
142  // Iterative solver
143  enum IterativeSolver { GMRES, BCGSTAB, TFQMR};
144  IterativeSolver iterative_solver_;
145 
146  // Should a preconditioner be used
147  bool use_preconditioner_;
148 
149  // Absolute tolerance
150  double abstol_;
151 
152  // Jacobian times vector function
153  Function jtimes_;
154 
155  // Get jtimes_
156  void get_jtimes();
157 
158  // Raise an error specific to KinSol
159  void kinsol_error(const std::string& module, int flag, bool fatal=true) const;
160 
162  static const std::string meta_doc;
163 
165  void* alloc_mem() const override { return new KinsolMemory(*this);}
166 
168  int init_mem(void* mem) const override;
169 
171  void free_mem(void *mem) const override { delete static_cast<KinsolMemory*>(mem);}
172 
174  void set_work(void* mem, const double**& arg, double**& res,
175  casadi_int*& iw, double*& w) const override;
176 
178  static KinsolMemory* to_mem(void *mem) {
179  KinsolMemory* m = static_cast<KinsolMemory*>(mem);
180  casadi_assert_dev(m);
181  return m;
182  }
183 
185  void func(KinsolMemory& m, N_Vector u, N_Vector fval) const;
186  void djac(KinsolMemory& m, long N, N_Vector u, N_Vector fu,
187  DlsMat J, N_Vector tmp1, N_Vector tmp2) const;
188  void bjac(KinsolMemory& m, long N, long mupper, long mlower, N_Vector u, N_Vector fu, DlsMat J,
189  N_Vector tmp1, N_Vector tmp2) const;
190  void jtimes(KinsolMemory& m, N_Vector v, N_Vector Jv, N_Vector u, int* new_u) const;
191  void psetup(KinsolMemory& m, N_Vector u, N_Vector uscale, N_Vector fval, N_Vector fscale,
192  N_Vector tmp1, N_Vector tmp2) const;
193  void psolve(KinsolMemory& m, N_Vector u, N_Vector uscale,
194  N_Vector fval, N_Vector fscale, N_Vector v, N_Vector tmp) const;
195 
197  static int func_wrapper(N_Vector u, N_Vector fval, void *user_data);
198  static int djac_wrapper(long N, N_Vector u, N_Vector fu, DlsMat J, void *user_data,
199  N_Vector tmp1, N_Vector tmp2);
200  static int bjac_wrapper(long N, long mupper, long mlower, N_Vector u, N_Vector fu, DlsMat J,
201  void *user_data, N_Vector tmp1, N_Vector tmp2);
202  static int jtimes_wrapper(N_Vector v, N_Vector Jv, N_Vector u, int* new_u, void *user_data);
203  static int psetup_wrapper(N_Vector u, N_Vector uscale, N_Vector fval, N_Vector fscale,
204  void* user_data, N_Vector tmp1, N_Vector tmp2);
205  static int psolve_wrapper(N_Vector u, N_Vector uscale, N_Vector fval, N_Vector fscale,
206  N_Vector v, void* user_data, N_Vector tmp);
207 
209  static int lsetup(KINMem kin_mem);
210  static int lsolve(KINMem kin_mem, N_Vector x, N_Vector b, double *sJpnorm, double *sFdotJp);
211  static void ehfun(int error_code, const char *module, const char *function,
212  char *msg, void *eh_data);
213  static void ihfun(const char *module, const char *function, char *msg, void *ih_data);
214  };
215 
216 } // namespace casadi
217 
219 #endif // CASADI_KINSOL_INTERFACE_HPP
The casadi namespace.
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.