26 #ifndef CASADI_KINSOL_INTERFACE_HPP
27 #define CASADI_KINSOL_INTERFACE_HPP
29 #include <casadi/interfaces/sundials/casadi_rootfinder_kinsol_export.h>
30 #include "casadi/core/rootfinder_impl.hpp"
31 #include <nvector/nvector_serial.h>
32 #include <sundials/sundials_dense.h>
33 #include <sundials/sundials_types.h>
34 #include <kinsol/kinsol.h>
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>
41 #include <kinsol/kinsol_spils_impl.h>
55 class KinsolInterface;
58 struct KinsolMemory :
public RootfinderMemory {
60 const KinsolInterface&
self;
63 KinsolMemory(
const KinsolInterface& s);
83 class KinsolInterface :
public Rootfinder {
86 explicit KinsolInterface(
const std::string& name,
const Function& f);
89 ~KinsolInterface()
override;
92 static Rootfinder* creator(
const std::string& name,
const Function& f) {
93 return new KinsolInterface(name, f);
98 static const Options options_;
99 const Options& get_options()
const override {
return options_;}
103 void init(
const Dict& opts)
override;
106 int solve(
void* mem)
const override;
109 const char* plugin_name()
const override {
return "kinsol";}
112 std::string class_name()
const override {
return "KinsolInterface";}
115 N_Vector u_scale_, f_scale_;
118 casadi_int strategy_;
121 bool disable_internal_warnings_;
124 casadi_int max_iter_;
127 casadi_int print_level_;
133 enum LinsolType { DENSE, BANDED, ITERATIVE, USER_DEFINED};
134 LinsolType linear_solver_type_;
137 casadi_int upper_bandwidth_, lower_bandwidth_;
143 enum IterativeSolver { GMRES, BCGSTAB, TFQMR};
144 IterativeSolver iterative_solver_;
147 bool use_preconditioner_;
159 void kinsol_error(
const std::string& module,
int flag,
bool fatal=
true)
const;
162 static const std::string meta_doc;
165 void* alloc_mem()
const override {
return new KinsolMemory(*
this);}
168 int init_mem(
void* mem)
const override;
171 void free_mem(
void *mem)
const override {
delete static_cast<KinsolMemory*
>(mem);}
174 void set_work(
void* mem,
const double**& arg,
double**& res,
175 casadi_int*& iw,
double*& w)
const override;
178 static KinsolMemory* to_mem(
void *mem) {
179 KinsolMemory* m =
static_cast<KinsolMemory*
>(mem);
180 casadi_assert_dev(m);
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;
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);
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);
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.