knitro_interface.cpp
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 #include "knitro_interface.hpp"
27 #include "casadi/core/casadi_misc.hpp"
28 #include <ctime>
29 #include <cstdio>
30 #include <cstdlib>
31 #include <algorithm>
32 #ifdef CASADI_WITH_THREAD
33 #ifdef CASADI_WITH_THREAD_MINGW
34 #include <mingw.thread.h>
35 #else // CASADI_WITH_THREAD_MINGW
36 #include <thread>
37 #endif // CASADI_WITH_THREAD_MINGW
38 #endif //CASADI_WITH_THREAD
39 
40 namespace casadi {
41 
42  extern "C"
43  int CASADI_NLPSOL_KNITRO_EXPORT
44  casadi_register_nlpsol_knitro(Nlpsol::Plugin* plugin) {
45  plugin->creator = KnitroInterface::creator;
46  plugin->name = "knitro";
47  plugin->doc = KnitroInterface::meta_doc.c_str();
48  plugin->version = CASADI_VERSION;
49  plugin->options = &KnitroInterface::options_;
50  plugin->deserialize = &KnitroInterface::deserialize;
51  return 0;
52  }
53 
54  extern "C"
55  void CASADI_NLPSOL_KNITRO_EXPORT casadi_load_nlpsol_knitro() {
57  }
58 
59  KnitroInterface::KnitroInterface(const std::string& name, const Function& nlp)
60  : Nlpsol(name, nlp) {
61  }
62 
63 
65  clear_mem();
66  }
67 
69  = {{&Nlpsol::options_},
70  {{"knitro",
71  {OT_DICT,
72  "Options to be passed to KNITRO"}},
73  {"options_file",
74  {OT_STRING,
75  "Read options from file (solver specific)"}},
76  {"detect_linear_constraints",
77  {OT_BOOL,
78  "Detect type of constraints"}},
79  {"contype",
80  {OT_INTVECTOR,
81  "Type of constraint"}},
82  {"complem_variables",
84  "List of complementary constraints on simple bounds. "
85  "Pair (i, j) encodes complementarity between the bounds on variable i and variable j."}}
86  }
87  };
88 
89  void KnitroInterface::init(const Dict& opts) {
90  // Call the init method of the base class
91  Nlpsol::init(opts);
92  bool detect_linear_constraints = true;
93  std::vector< std::vector<casadi_int> > complem_variables;
94  options_file_ = "";
95 
96  // Read user options
97  for (auto&& op : opts) {
98  if (op.first=="knitro") {
99  opts_ = op.second;
100  } else if (op.first=="options_file") {
101  options_file_ = op.second.to_string();
102  } else if (op.first=="contype") {
103  contype_ = op.second;
104  } else if (op.first=="detect_linear_constraints") {
105  detect_linear_constraints = op.second;
106  } else if (op.first=="complem_variables") {
107  complem_variables = op.second;
108  }
109  }
110 
111  // Type of constraints, general by default
112  if (contype_.empty()) {
113  contype_.resize(ng_, KN_CONTYPE_GENERAL);
114  if (detect_linear_constraints) {
115  std::vector<bool> nl_g = oracle_.which_depends("x", {"g"}, 2, true);
116  for (casadi_int i=0;i<ng_;++i)
117  contype_[i] = nl_g[i] ? KN_CONTYPE_GENERAL : KN_CONTYPE_LINEAR;
118  }
119  }
120 
121  casadi_assert_dev(contype_.size()==ng_);
122 
123 
124  comp_type_.resize(complem_variables.size(), KN_CCTYPE_VARVAR);
125  comp_i1_.reserve(complem_variables.size());
126  comp_i2_.reserve(complem_variables.size());
127  for (auto && e : complem_variables) {
128  casadi_assert(e.size()==2, "Complementary constraints must come in pairs.");
129  casadi_assert(e[0]>=0, "Invalid variable index.");
130  casadi_assert(e[1]>=0, "Invalid variable index.");
131  casadi_assert(e[0]<nx_, "Invalid variable index.");
132  casadi_assert(e[1]<nx_, "Invalid variable index.");
133  comp_i1_.push_back(e[0]);
134  comp_i2_.push_back(e[1]);
135  }
136 
137  // Setup NLP functions
138  create_function("nlp_fg", {"x", "p"}, {"f", "g"});
139  Function gf_jg_fcn = create_function("nlp_gf_jg", {"x", "p"}, {"grad:f:x", "jac:g:x"});
140  jacg_sp_ = gf_jg_fcn.sparsity_out(1);
141  Function hess_l_fcn = create_function("nlp_hess_l", {"x", "p", "lam:f", "lam:g"},
142  {"triu:hess:gamma:x:x"},
143  {{"gamma", {"f", "g"}}});
144  hesslag_sp_ = hess_l_fcn.sparsity_out(0);
145 
146  unsigned int hc = 0;
147  #ifdef CASADI_WITH_THREAD
148  //may return 0 when not able to detect. If it's the case, then return 8.
149  hc = std::thread::hardware_concurrency();
150  #endif
151  int processor_count = hc ? hc : 8;
152  //Obtain maximum number of threads needed
153  int ms_numthreads = 1;
154  int findiff_numthreads = 1;
155  int numthreads = 1;
156  int mip_numthreads = 1;
157  for (auto&& op : opts_) {
158  if (op.first=="ms_numthreads") {
159  ms_numthreads = op.second;
160  }
161  if (op.first=="findiff_numthreads") {
162  findiff_numthreads = op.second;
163  }
164  if (op.first=="numthreads") {
165  numthreads = op.second;
166  }
167  if (op.first=="mip_numthreads") {
168  mip_numthreads = op.second;
169  }
170  }
171  max_num_threads_ = std::max({processor_count, ms_numthreads,
172  findiff_numthreads, numthreads, mip_numthreads});
173 
174  // Allocate persistent memory
175  alloc_w(nx_, true); // wlbx_
176  alloc_w(nx_, true); // wubx_
177  alloc_w(ng_, true); // wlbg_
178  alloc_w(ng_, true); // wubg_
179  }
180 
181  int KnitroInterface::init_mem(void* mem) const {
182  return Nlpsol::init_mem(mem);
183  //auto m = static_cast<KnitroMemory*>(mem);
184 
185  // Commented out since I have not found out how to change the bounds
186  // Allocate KNITRO memory block
187  // m.kc = KN_new();
188  }
189 
190  void KnitroInterface::set_work(void* mem, const double**& arg, double**& res,
191  casadi_int*& iw, double*& w) const {
192  auto m = static_cast<KnitroMemory*>(mem);
193 
194  // Set work in base classes
195  Nlpsol::set_work(mem, arg, res, iw, w);
196 
197  // Copy inputs to temporary arrays
198  m->wlbx = w; w += nx_;
199  m->wubx = w; w += nx_;
200  m->wlbg = w; w += ng_;
201  m->wubg = w; w += ng_;
202  }
203 
204  int casadi_KN_puts(const char * const str, void * const userParams) {
205  std::string s(str);
206  uout() << s << std::flush;
207  return s.size();
208  }
209 
210  int KnitroInterface::solve(void* mem) const {
211  auto m = static_cast<KnitroMemory*>(mem);
212  auto d_nlp = &m->d_nlp;
213  casadi_int status;
214 
215  // Allocate KNITRO memory block (move back to init!)
216  casadi_assert_dev(m->kc==nullptr);
217  status = KN_new(&m->kc);
218  casadi_assert_dev(m->kc!=nullptr);
219 
220  status = KN_set_puts_callback(m->kc, casadi_KN_puts, nullptr);
221  casadi_assert(status == 0, "KN_set_puts_callback failed");
222 
223  // Jacobian sparsity
224  std::vector<int> Jcol, Jrow;
225  if (!jacg_sp_.is_null()) {
226  assign_vector(jacg_sp_.get_col(), Jcol);
227  assign_vector(jacg_sp_.get_row(), Jrow);
228  }
229 
230  if (!options_file_.empty()) {
231  status = KN_load_param_file(m->kc, options_file_.c_str());
232  casadi_assert(status==0, "KN_load_param_file failed");
233  }
234 
235  // Hessian sparsity
236  casadi_int nnzH = hesslag_sp_.is_null() ? 0 : hesslag_sp_.nnz();
237  std::vector<int> Hcol, Hrow;
238  if (nnzH>0) {
241  status = KN_set_int_param_by_name(m->kc, "hessopt", KN_HESSOPT_EXACT);
242  casadi_assert(status==0, "KN_set_int_param failed");
243  } else {
244  status = KN_set_int_param_by_name(m->kc, "hessopt", KN_HESSOPT_LBFGS);
245  casadi_assert(status==0, "KN_set_int_param failed");
246  }
247 
248  // Pass user set options
249  for (auto&& op : opts_) {
250  int param_id;
251  casadi_assert(KN_get_param_id(m->kc, op.first.c_str(), &param_id)==0,
252  "Unknown parameter '" + op.first + "'.");
253 
254  int param_type;
255  casadi_assert(!KN_get_param_type(m->kc, param_id, &param_type),
256  "Error when setting option '" + op.first + "'.");
257 
258  switch (param_type) {
259  case KN_PARAMTYPE_INTEGER:
260  casadi_assert(!KN_set_int_param(m->kc, param_id, op.second),
261  "Error when setting option '" + op.first + "'.");
262  continue;
263  case KN_PARAMTYPE_FLOAT:
264  casadi_assert(!KN_set_double_param(m->kc, param_id, op.second),
265  "Error when setting option '" + op.first + "'.");
266  continue;
267  case KN_PARAMTYPE_STRING:
268  {
269  std::string str = op.second.to_string();
270  casadi_assert(!KN_set_char_param(m->kc, param_id, str.c_str()),
271  "Error when setting option '" + op.first + "'.");
272  }
273  continue;
274  default:
275  casadi_error("Error when setting option '" + op.first + "'.");
276  }
277  }
278 
279  // "Correct" upper and lower bounds
280  casadi_copy(d_nlp->lbz, nx_, m->wlbx);
281  casadi_copy(d_nlp->ubz, nx_, m->wubx);
282  casadi_copy(d_nlp->lbz+nx_, ng_, m->wlbg);
283  casadi_copy(d_nlp->ubz+nx_, ng_, m->wubg);
284  for (casadi_int i=0; i<nx_; ++i) if (isinf(m->wlbx[i])) m->wlbx[i] = -KN_INFINITY;
285  for (casadi_int i=0; i<nx_; ++i) if (isinf(m->wubx[i])) m->wubx[i] = KN_INFINITY;
286  for (casadi_int i=0; i<ng_; ++i) if (isinf(m->wlbg[i])) m->wlbg[i] = -KN_INFINITY;
287  for (casadi_int i=0; i<ng_; ++i) if (isinf(m->wubg[i])) m->wubg[i] = KN_INFINITY;
288 
289  std::vector<int> xindex(nx_);
290  std::vector<int> gindex(ng_);
291  iota(begin(xindex), end(xindex), 0);
292  iota(begin(gindex), end(gindex), 0);
293 
294  status = KN_add_vars(m->kc, nx_, get_ptr(xindex));
295  casadi_assert(status==0, "KN_add_vars failed");
296  status = KN_set_obj_goal(m->kc, KN_OBJGOAL_MINIMIZE);
297  casadi_assert(status==0, "KN_set_obj_goal failed");
298  status = KN_set_var_lobnds_all(m->kc, m->wlbx);
299  casadi_assert(status==0, "KN_set_var_lobnds failed");
300  status = KN_set_var_upbnds_all(m->kc, m->wubx);
301  casadi_assert(status==0, "KN_set_var_upbnds failed");
302  status = KN_add_cons(m->kc, ng_, get_ptr(gindex));
303  casadi_assert(status==0, "KN_add_cons failed");
304  status = KN_set_con_lobnds_all(m->kc, m->wlbg);
305  casadi_assert(status==0, "KN_set_con_lobnds failed");
306  status = KN_set_con_upbnds_all(m->kc, m->wubg);
307  casadi_assert(status==0, "KN_set_con_upbnds failed");
308  status = KN_set_var_primal_init_values_all(m->kc, d_nlp->z);
309  casadi_assert(status==0, "KN_set_var_primal_init_values failed");
310  if (mi_) {
311  // Types of variables
312  std::vector<int> vtype;
313  vtype.reserve(nx_);
314  for (auto&& e : discrete_) {
315  vtype.push_back(e ? KN_VARTYPE_INTEGER : KN_VARTYPE_CONTINUOUS);
316  }
317  status = KN_set_var_types_all(m->kc, get_ptr(vtype));
318  casadi_assert(status==0, "KN_set_var_types failed");
319  }
320 
321  // Complementarity constraints
322  status = KN_set_compcons(m->kc, comp_i1_.size(),
324  casadi_assert(status==0, "KN_set_compcons failed");
325 
326  // Register callback functions
327  status = KN_add_eval_callback(m->kc, true, ng_, get_ptr(gindex), &callback, &m->cb);
328  casadi_assert(status==0, "KN_add_eval_callback failed");
329 
330  status = KN_set_cb_grad(m->kc, m->cb, nx_, get_ptr(xindex),
331  Jcol.size(), get_ptr(Jrow), get_ptr(Jcol), &callback);
332  casadi_assert(status==0, "KN_set_cb_grad failed");
333 
334  if (nnzH>0) {
335  status = KN_set_cb_hess(m->kc, m->cb, nnzH, get_ptr(Hrow), get_ptr(Hcol), &callback);
336  casadi_assert(status==0, "KN_set_cb_hess failed");
337  }
338 
339  status = KN_set_cb_user_params(m->kc, m->cb, static_cast<void*>(m));
340  casadi_assert(status==0, "KN_set_cb_user_params failed");
341 
342  // NumThreads to 1 to prevent segmentation fault
343  status = KN_set_int_param_by_name(m->kc, "numthreads", 1);
344  casadi_assert(status==0, "KN_set_cb_user_params failed");
345 
346  // Lagrange multipliers
347  std::vector<double> lambda(nx_+ng_);
348 
349  // objective solution
350  double objSol;
351 
352  // Solve NLP
353  status = KN_solve(m->kc);
354  int statusKnitro = static_cast<int>(status);
355 
356  m->return_status = return_codes(status);
357  m->success = status==KN_RC_OPTIMAL_OR_SATISFACTORY ||
358  status==KN_RC_NEAR_OPT;
359  if (status==KN_RC_ITER_LIMIT_FEAS ||
360  status==KN_RC_TIME_LIMIT_FEAS ||
361  status==KN_RC_FEVAL_LIMIT_FEAS ||
362  status==KN_RC_ITER_LIMIT_INFEAS ||
363  status==KN_RC_TIME_LIMIT_INFEAS ||
364  status==KN_RC_FEVAL_LIMIT_INFEAS)
365  m->unified_return_status = SOLVER_RET_LIMITED;
366 
367  // Output optimal cost
368  casadi_int error;
369  error = KN_get_solution(m->kc, &statusKnitro, &objSol, d_nlp->z, get_ptr(lambda));
370  casadi_assert(error == 0, "KN_get_solution failed");
371  // Output dual solution
372  casadi_copy(get_ptr(lambda), ng_, d_nlp->lam + nx_);
373  casadi_copy(get_ptr(lambda)+ng_, nx_, d_nlp->lam);
374 
375  d_nlp->objective = objSol;
376 
377  // Calculate constraints
378  if (ng_>0) {
379  m->arg[0] = d_nlp->z;
380  m->arg[1] = d_nlp->p;
381  m->res[0] = nullptr;
382  m->res[1] = d_nlp->z + nx_;
383  calc_function(m, "nlp_fg");
384  }
385 
386  // Free memory (move to destructor!)
387  status = KN_free(&m->kc);
388  casadi_assert(status == 0, "KN_free failed");
389  m->kc = nullptr;
390  return 0;
391  }
392 
393  int KnitroInterface::callback(KN_context_ptr kc,
394  CB_context_ptr cb,
395  KN_eval_request_ptr const evalRequest,
396  KN_eval_result_ptr const evalResult,
397  void * const userParams) {
398  try {
399  int thread_id = evalRequest->threadID;
400  // Get a pointer to the calling object
401  auto m = static_cast<KnitroMemory*>(userParams);
402  // Get thread local memory
403  auto ml = m->thread_local_mem.at(thread_id);
404  auto d_nlp = &m->d_nlp;
405  const double *x;
406  // Direct to the correct function
407  switch (evalRequest->type) {
408  case KN_RC_EVALFC:
409  double *obj;
410  double *c;
411  x = evalRequest->x;
412  obj = evalResult->obj;
413  c = evalResult->c;
414  ml->arg[0] = x;
415  ml->arg[1] = d_nlp->p;
416  ml->res[0] = obj;
417  ml->res[1] = c;
418  if (m->self.calc_function(m, "nlp_fg", nullptr, thread_id)) return KN_RC_EVAL_ERR;
419  break;
420  case KN_RC_EVALGA:
421  double *objGrad;
422  double *jac;
423  x = evalRequest->x;
424  objGrad = evalResult->objGrad;
425  jac = evalResult->jac;
426  ml->arg[0] = x;
427  ml->arg[1] = d_nlp->p;
428  ml->res[0] = objGrad;
429  ml->res[1] = jac;
430  if (m->self.calc_function(m, "nlp_gf_jg", nullptr, thread_id)) return KN_RC_EVAL_ERR;
431  break;
432  case KN_RC_EVALH_NO_F:
433  case KN_RC_EVALH:
434  const double *lambda;
435  double sigma;
436  double *hess;
437  x = evalRequest->x;
438  hess = evalResult->hess;
439  lambda = evalRequest->lambda;
440  sigma = *(evalRequest->sigma);
441  ml->arg[0] = x;
442  ml->arg[1] = d_nlp->p;
443  ml->arg[2] = &sigma;
444  ml->arg[3] = lambda;
445  ml->res[0] = hess;
446  if (m->self.calc_function(m, "nlp_hess_l", nullptr, thread_id)) {
447  casadi_error("calc_hess_l failed");
448  }
449  break;
450  default:
451  casadi_error("KnitroInterface::callback: unknown method");
452  }
453 
454  return 0;
455  } catch(KeyboardInterruptException& ex) {
456  return KN_RC_USER_TERMINATION;
457  } catch(std::exception& ex) {
458  uerr() << "KnitroInterface::callback caught exception: "
459  << ex.what() << std::endl;
460  return -1;
461  }
462 
463  }
464 
465  const char* KnitroInterface::return_codes(int flag) {
466  switch (flag) {
467  case KN_RC_OPTIMAL_OR_SATISFACTORY: return "KN_RC_OPTIMAL_OR_SATISFACTORY";
468  case KN_RC_NEAR_OPT: return "KN_RC_NEAR_OPT";
469  case KN_RC_FEAS_XTOL: return "KN_RC_FEAS_XTOL";
470  case KN_RC_FEAS_NO_IMPROVE: return "KN_RC_FEAS_NO_IMPROVE";
471  case KN_RC_FEAS_FTOL: return "KN_RC_FEAS_FTOL";
472  case KN_RC_INFEASIBLE: return "KN_RC_INFEASIBLE";
473  case KN_RC_INFEAS_XTOL: return "KN_RC_INFEAS_XTOL";
474  case KN_RC_INFEAS_NO_IMPROVE: return "KN_RC_INFEAS_NO_IMPROVE";
475  case KN_RC_INFEAS_MULTISTART: return "KN_RC_INFEAS_MULTISTART";
476  case KN_RC_INFEAS_CON_BOUNDS: return "KN_RC_INFEAS_CON_BOUNDS";
477  case KN_RC_INFEAS_VAR_BOUNDS: return "KN_RC_INFEAS_VAR_BOUNDS";
478  case KN_RC_UNBOUNDED: return "KN_RC_UNBOUNDED";
479  case KN_RC_ITER_LIMIT_FEAS: return "KN_RC_ITER_LIMIT_FEAS";
480  case KN_RC_TIME_LIMIT_FEAS: return "KN_RC_TIME_LIMIT_FEAS";
481  case KN_RC_FEVAL_LIMIT_FEAS: return "KN_RC_FEVAL_LIMIT_FEAS";
482  case KN_RC_MIP_EXH_FEAS: return "KN_RC_MIP_EXH_FEAS";
483  case KN_RC_MIP_TERM_FEAS: return "KN_RC_MIP_TERM_FEAS";
484  case KN_RC_MIP_SOLVE_LIMIT_FEAS: return "KN_RC_MIP_SOLVE_LIMIT_FEAS";
485  case KN_RC_MIP_NODE_LIMIT_FEAS: return "KN_RC_MIP_NODE_LIMIT_FEAS";
486  case KN_RC_ITER_LIMIT_INFEAS: return "KN_RC_ITER_LIMIT_INFEAS";
487  case KN_RC_TIME_LIMIT_INFEAS: return "KN_RC_TIME_LIMIT_INFEAS";
488  case KN_RC_FEVAL_LIMIT_INFEAS: return "KN_RC_FEVAL_LIMIT_INFEAS";
489  case KN_RC_MIP_EXH_INFEAS: return "KN_RC_MIP_EXH_INFEAS";
490  case KN_RC_MIP_SOLVE_LIMIT_INFEAS: return "KN_RC_MIP_SOLVE_LIMIT_INFEAS";
491  case KN_RC_MIP_NODE_LIMIT_INFEAS: return "KN_RC_MIP_NODE_LIMIT_INFEAS";
492  case KN_RC_CALLBACK_ERR: return "KN_RC_CALLBACK_ERR";
493  case KN_RC_LP_SOLVER_ERR: return "KN_RC_LP_SOLVER_ERR";
494  case KN_RC_EVAL_ERR: return "KN_RC_EVAL_ERR";
495  case KN_RC_OUT_OF_MEMORY: return "KN_RC_OUT_OF_MEMORY";
496  case KN_RC_USER_TERMINATION: return "KN_RC_USER_TERMINATION";
497  case KN_RC_OPEN_FILE_ERR: return "KN_RC_OPEN_FILE_ERR";
498  case KN_RC_BAD_N_OR_F: return "KN_RC_BAD_N_OR_F";
499  case KN_RC_BAD_CONSTRAINT: return "KN_RC_BAD_CONSTRAINT";
500  case KN_RC_BAD_JACOBIAN: return "KN_RC_BAD_JACOBIAN";
501  case KN_RC_BAD_HESSIAN: return "KN_RC_BAD_HESSIAN";
502  case KN_RC_BAD_CON_INDEX: return "KN_RC_BAD_CON_INDEX";
503  case KN_RC_BAD_JAC_INDEX: return "KN_RC_BAD_JAC_INDEX";
504  case KN_RC_BAD_HESS_INDEX: return "KN_RC_BAD_HESS_INDEX";
505  case KN_RC_BAD_CON_BOUNDS: return "KN_RC_BAD_CON_BOUNDS";
506  case KN_RC_BAD_VAR_BOUNDS: return "KN_RC_BAD_VAR_BOUNDS";
507  case KN_RC_ILLEGAL_CALL: return "KN_RC_ILLEGAL_CALL";
508  case KN_RC_BAD_KCPTR: return "KN_RC_BAD_KCPTR";
509  case KN_RC_NULL_POINTER: return "KN_RC_NULL_POINTER";
510  case KN_RC_BAD_INIT_VALUE: return "KN_RC_BAD_INIT_VALUE";
511  case KN_RC_BAD_PARAMINPUT: return "KN_RC_BAD_PARAMINPUT";
512  case KN_RC_LINEAR_SOLVER_ERR: return "KN_RC_LINEAR_SOLVER_ERR";
513  case KN_RC_DERIV_CHECK_FAILED: return "KN_RC_DERIV_CHECK_FAILED";
514  case KN_RC_DERIV_CHECK_TERMINATE: return "KN_RC_DERIV_CHECK_TERMINATE";
515  case KN_RC_INTERNAL_ERROR: return "KN_RC_INTERNAL_ERROR";
516  }
517  return nullptr;
518  }
519 
520  Dict KnitroInterface::get_stats(void* mem) const {
521  Dict stats = Nlpsol::get_stats(mem);
522  auto m = static_cast<KnitroMemory*>(mem);
523  stats["return_status"] = m->return_status;
524 
525  return stats;
526  }
527 
529  int version = s.version("KnitroInterface", 1, 2);
530  s.unpack("KnitroInterface::contype", contype_);
531  s.unpack("KnitroInterface::comp_type", comp_type_);
532  s.unpack("KnitroInterface::comp_i1", comp_i1_);
533  s.unpack("KnitroInterface::comp_i2", comp_i2_);
534  s.unpack("KnitroInterface::opts", opts_);
535  s.unpack("KnitroInterface::jacg_sp", jacg_sp_);
536  s.unpack("KnitroInterface::hesslag_sp", hesslag_sp_);
537  if (version>=2) {
538  s.unpack("KnitroInterface::options_file", options_file_);
539  } else {
540  options_file_ = "";
541  }
542  }
543 
546  s.version("KnitroInterface", 2);
547  s.pack("KnitroInterface::contype", contype_);
548  s.pack("KnitroInterface::comp_type", comp_type_);
549  s.pack("KnitroInterface::comp_i1", comp_i1_);
550  s.pack("KnitroInterface::comp_i2", comp_i2_);
551  s.pack("KnitroInterface::opts", opts_);
552  s.pack("KnitroInterface::jacg_sp", jacg_sp_);
553  s.pack("KnitroInterface::hesslag_sp", hesslag_sp_);
554  s.pack("KnitroInterface::options_file", options_file_);
555  }
556 
557  KnitroMemory::KnitroMemory(const KnitroInterface& self) : self(self) {
558  this->kc = nullptr;
559  }
560 
562  // Currently no persistent memory since KNITRO requires knowledge of nature of bounds
563  // if (this->kc) {
564  // KTR_free(&this->kc);
565  //}
566  }
567 
568 } // namespace casadi
const char * what() const override
Display error.
Definition: exception.hpp:90
Helper class for Serialization.
void unpack(Sparsity &e)
Reconstruct an object from the input stream.
void version(const std::string &name, int v)
void alloc_w(size_t sz_w, bool persistent=false)
Ensure required length of w field.
Function object.
Definition: function.hpp:60
std::vector< bool > which_depends(const std::string &s_in, const std::vector< std::string > &s_out, casadi_int order=1, bool tr=false) const
Which variables enter with some order.
Definition: function.cpp:1826
const Sparsity & sparsity_out(casadi_int ind) const
Get sparsity of a given output.
Definition: function.cpp:1031
bool is_null() const
Is a null pointer?
'knitro' plugin for Nlpsol
std::vector< int > comp_i1_
static int callback(KN_context_ptr kc, CB_context_ptr cb, KN_eval_request_ptr const evalRequest, KN_eval_result_ptr const evalResult, void *const userParams)
static const std::string meta_doc
A documentation string.
Dict get_stats(void *mem) const override
Get all statistics.
static Nlpsol * creator(const std::string &name, const Function &nlp)
Create a new NLP Solver.
void set_work(void *mem, const double **&arg, double **&res, casadi_int *&iw, double *&w) const override
Set the (persistent) work vectors.
void serialize_body(SerializingStream &s) const override
Serialize an object without type information.
std::vector< int > comp_type_
static const Options options_
Options.
static const char * return_codes(int flag)
void init(const Dict &opts) override
Initialize.
int init_mem(void *mem) const override
Initalize memory block.
std::vector< int > comp_i2_
static ProtoFunction * deserialize(DeserializingStream &s)
Deserialize into MX.
std::vector< int > contype_
int solve(void *mem) const override
KnitroInterface(const std::string &name, const Function &nlp)
NLP solver storage class.
Definition: nlpsol_impl.hpp:59
Dict get_stats(void *mem) const override
Get all statistics.
Definition: nlpsol.cpp:1162
static const Options options_
Options.
void init(const Dict &opts) override
Initialize.
Definition: nlpsol.cpp:420
casadi_int ng_
Number of constraints.
Definition: nlpsol_impl.hpp:69
int init_mem(void *mem) const override
Initalize memory block.
Definition: nlpsol.cpp:603
void serialize_body(SerializingStream &s) const override
Serialize an object without type information.
Definition: nlpsol.cpp:1306
std::vector< bool > discrete_
Options.
casadi_int nx_
Number of variables.
Definition: nlpsol_impl.hpp:66
void set_work(void *mem, const double **&arg, double **&res, casadi_int *&iw, double *&w) const override
Set the (persistent) work vectors.
Definition: nlpsol.cpp:795
Function oracle_
Oracle: Used to generate other functions.
Function create_function(const Function &oracle, const std::string &fname, const std::vector< std::string > &s_in, const std::vector< std::string > &s_out, const Function::AuxOut &aux=Function::AuxOut(), const Dict &opts=Dict())
int calc_function(OracleMemory *m, const std::string &fcn, const double *const *arg=nullptr, int thread_id=0) const
static void registerPlugin(const Plugin &plugin, bool needs_lock=true)
Register an integrator in the factory.
void clear_mem()
Clear all memory (called from destructor)
Helper class for Serialization.
void version(const std::string &name, int v)
void pack(const Sparsity &e)
Serializes an object to the output stream.
std::vector< casadi_int > get_col() const
Get the column for each non-zero entry.
Definition: sparsity.cpp:368
casadi_int nnz() const
Get the number of (structural) non-zeros.
Definition: sparsity.cpp:148
std::vector< casadi_int > get_row() const
Get the row for each non-zero entry.
Definition: sparsity.cpp:372
The casadi namespace.
Definition: archiver.cpp:28
int casadi_KN_puts(const char *const str, void *const userParams)
void assign_vector(const std::vector< S > &s, std::vector< D > &d)
std::ostream & uerr()
void casadi_copy(const T1 *x, casadi_int n, T1 *y)
COPY: y <-x.
void CASADI_NLPSOL_KNITRO_EXPORT casadi_load_nlpsol_knitro()
@ OT_INTVECTOR
@ OT_INTVECTORVECTOR
std::string str(const T &v)
String representation, any type.
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
T * get_ptr(std::vector< T > &v)
Get a pointer to the data contained in the vector.
int CASADI_NLPSOL_KNITRO_EXPORT casadi_register_nlpsol_knitro(Nlpsol::Plugin *plugin)
std::ostream & uout()
@ SOLVER_RET_LIMITED
KnitroMemory(const KnitroInterface &self)
Constructor.
casadi_nlpsol_data< double > d_nlp
Definition: nlpsol_impl.hpp:42
Options metadata for a class.
Definition: options.hpp:40
std::vector< LocalOracleMemory * > thread_local_mem