conic.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 "conic_impl.hpp"
27 #include "nlpsol_impl.hpp"
28 #include "filesystem_impl.hpp"
29 
30 namespace casadi {
31 
32  bool has_conic(const std::string& name) {
33  return Conic::has_plugin(name);
34  }
35 
36  void load_conic(const std::string& name) {
37  Conic::load_plugin(name);
38  }
39 
40  std::string doc_conic(const std::string& name) {
41  return Conic::getPlugin(name).doc;
42  }
43 
44  Function conic(const std::string& name, const std::string& solver,
45  const SpDict& qp, const Dict& opts) {
46  return Function::create(Conic::instantiate(name, solver, qp), opts);
47  }
48 
49  void conic_debug(const Function& f, const std::string &filename) {
50  std::ofstream file;
52  conic_debug(f, file);
53  }
54 
55  void conic_debug(const Function& f, std::ostream &file) {
56  casadi_assert_dev(!f.is_null());
57  const Conic* n = f.get<Conic>();
58  return n->generateNativeCode(file);
59  }
60 
61  std::vector<std::string> conic_in() {
62  std::vector<std::string> ret(conic_n_in());
63  for (size_t i=0; i<ret.size(); ++i) ret[i]=conic_in(i);
64  return ret;
65  }
66 
67  std::vector<std::string> conic_out() {
68  std::vector<std::string> ret(conic_n_out());
69  for (size_t i=0; i<ret.size(); ++i) ret[i]=conic_out(i);
70  return ret;
71  }
72 
73  std::string conic_in(casadi_int ind) {
74  switch (static_cast<ConicInput>(ind)) {
75  case CONIC_H: return "h";
76  case CONIC_G: return "g";
77  case CONIC_A: return "a";
78  case CONIC_Q: return "q";
79  case CONIC_P: return "p";
80  case CONIC_LBA: return "lba";
81  case CONIC_UBA: return "uba";
82  case CONIC_LBX: return "lbx";
83  case CONIC_UBX: return "ubx";
84  case CONIC_X0: return "x0";
85  case CONIC_LAM_X0: return "lam_x0";
86  case CONIC_LAM_A0: return "lam_a0";
87  case CONIC_NUM_IN: break;
88  }
89  return std::string();
90  }
91 
92  std::string conic_out(casadi_int ind) {
93  switch (static_cast<ConicOutput>(ind)) {
94  case CONIC_X: return "x";
95  case CONIC_COST: return "cost";
96  case CONIC_LAM_A: return "lam_a";
97  case CONIC_LAM_X: return "lam_x";
98  case CONIC_NUM_OUT: break;
99  }
100  return std::string();
101  }
102 
103  casadi_int conic_n_in() {
104  return CONIC_NUM_IN;
105  }
106 
107  casadi_int conic_n_out() {
108  return CONIC_NUM_OUT;
109  }
110 
111  template<typename M>
112  Function qpsol_nlp(const std::string& name, const std::string& solver,
113  const std::map<std::string, M>& qp, const Dict& opts) {
114  // We have: minimize f(x) = 1/2 * x' H x + c'x
115  // subject to lbx <= x <= ubx
116  // lbg <= g(x) = A x + b <= ubg
117  // h(x) >=0 (psd)
118 
119  // Extract 'expand' option
120  bool expand = false;
121  Dict opt = opts;
122  extract_from_dict_inplace(opt, "expand", expand);
123  bool postpone_expand = false;
124  extract_from_dict_inplace(opt, "postpone_expand", postpone_expand);
125  bool error_on_fail = get_from_dict(opts, "error_on_fail", true);
126 
127  if (expand && !postpone_expand && M::type_name()=="MX") {
128  Function f = Function("f", qp, {"x", "p"}, {"f", "g", "h"});
129  std::vector<SX> arg = f.sx_in();
130  std::vector<SX> res = f(arg);
131  SXDict qp_mod;
132  for (casadi_int i=0;i<f.n_in();++i) qp_mod[f.name_in(i)] = arg[i];
133  for (casadi_int i=0;i<f.n_out();++i) qp_mod[f.name_out(i)] = res[i];
134  return qpsol_nlp(name, solver, qp_mod, opt);
135  }
136 
137  M x, p, f, g, h;
138  for (auto&& i : qp) {
139  if (i.first=="x") {
140  x = i.second;
141  } else if (i.first=="p") {
142  p = i.second;
143  } else if (i.first=="f") {
144  f = i.second;
145  } else if (i.first=="g") {
146  g = i.second;
147  } else if (i.first=="h") {
148  h = i.second;
149  } else {
150  casadi_error("No such field: " + i.first);
151  }
152  }
153 
154  if (f.is_empty()) f = 0;
155  if (g.is_empty()) g = M(0, 1);
156 
157  // Dimension checks
158  casadi_assert(g.is_dense() && g.is_vector(),
159  "Expected a dense vector 'g', but got " + g.dim() + ".");
160 
161  casadi_assert(f.is_dense() && f.is_scalar(),
162  "Expected a dense scalar 'f', but got " + f.dim() + ".");
163 
164  casadi_assert(x.is_dense() && x.is_vector(),
165  "Expected a dense vector 'x', but got " + x.dim() + ".");
166 
167  casadi_assert(h.is_square(),
168  "Expected a symmetric matrix 'h', but got " + h.dim() + ".");
169 
170  if (g.is_empty(true)) g = M(0, 1); // workaround
171 
172  // Gradient of the objective: gf == Hx + g
173  M gf = M::gradient(f, x);
174 
175  // Identify the linear term in the objective
176  M c = substitute(gf, x, M::zeros(x.sparsity()));
177 
178  // Identify the quadratic term in the objective
179  M H = M::jacobian(gf, x, {{"symmetric", true}});
180 
181  // Identify constant term in the objective
182  Function r("constant_qp", {x, p}, {substitute(f, x, M::zeros(x.sparsity()))});
183 
184  // Identify the constant term in the constraints
185  M b = substitute(g, x, M::zeros(x.sparsity()));
186 
187  // Identify the linear term in the constraints
188  M A = M::jacobian(g, x);
189 
190  // Identify the constant term in the psd constraints
191  M P = substitute(h, x, M::zeros(x.sparsity()));
192 
193  // Identify the linear term in the psd constraints
194  M Q = M::jacobian(h, x);
195 
196  // Create a function for calculating the required matrices vectors
197  Function prob(name + "_qp", {x, p}, {H, c, A, b, Q, P},
198  {"x", "p"}, {"H", "c", "A", "b", "Q", "P"});
199  if (expand && postpone_expand) prob = prob.expand();
200 
201  // Make sure that the problem is sound
202  casadi_assert(!prob.has_free(), "Cannot create '" + prob.name() + "' "
203  "since " + str(prob.get_free()) + " are free.");
204 
205  // Create the QP solver
206  Function conic_f = conic(name + "_qpsol", solver,
207  {{"h", H.sparsity()}, {"a", A.sparsity()},
208  {"p", P.sparsity()}, {"q", Q.sparsity()}}, opt);
209 
210  // Create an MXFunction with the right signature
211  std::vector<MX> ret_in(NLPSOL_NUM_IN);
212  ret_in[NLPSOL_X0] = MX::sym("x0", x.sparsity());
213  ret_in[NLPSOL_P] = MX::sym("p", p.sparsity());
214  ret_in[NLPSOL_LBX] = MX::sym("lbx", x.sparsity());
215  ret_in[NLPSOL_UBX] = MX::sym("ubx", x.sparsity());
216  ret_in[NLPSOL_LBG] = MX::sym("lbg", g.sparsity());
217  ret_in[NLPSOL_UBG] = MX::sym("ubg", g.sparsity());
218  ret_in[NLPSOL_LAM_X0] = MX::sym("lam_x0", x.sparsity());
219  ret_in[NLPSOL_LAM_G0] = MX::sym("lam_g0", g.sparsity());
220  std::vector<MX> ret_out(NLPSOL_NUM_OUT);
221 
222 
223  std::vector<MX> v(NL_NUM_IN);
224  v[NL_X] = ret_in[NLPSOL_X0];
225  v[NL_P] = ret_in[NLPSOL_P];
226  // Evaluate constant part of objective
227  MX rv = r(v)[0];
228  // Get expressions for the QP matrices and vectors
229  v = prob(v);
230 
231  // Call the QP solver
232  std::vector<MX> w(CONIC_NUM_IN);
233  w[CONIC_H] = v.at(0);
234  w[CONIC_G] = v.at(1);
235  w[CONIC_A] = v.at(2);
236  w[CONIC_Q] = v.at(4);
237  w[CONIC_P] = v.at(5);
238  w[CONIC_LBX] = ret_in[NLPSOL_LBX];
239  w[CONIC_UBX] = ret_in[NLPSOL_UBX];
240  w[CONIC_LBA] = ret_in[NLPSOL_LBG] - v.at(3);
241  w[CONIC_UBA] = ret_in[NLPSOL_UBG] - v.at(3);
242  w[CONIC_X0] = ret_in[NLPSOL_X0];
243  w[CONIC_LAM_X0] = ret_in[NLPSOL_LAM_X0];
244  w[CONIC_LAM_A0] = ret_in[NLPSOL_LAM_G0];
245  w = conic_f(w);
246 
247  // Get expressions for the solution
248  ret_out[NLPSOL_X] = reshape(w[CONIC_X], x.size());
249  ret_out[NLPSOL_F] = rv + w[CONIC_COST];
250  ret_out[NLPSOL_G] = reshape(mtimes(v.at(2), w[CONIC_X]), g.size()) + v.at(3);
251  ret_out[NLPSOL_LAM_X] = reshape(w[CONIC_LAM_X], x.size());
252  ret_out[NLPSOL_LAM_G] = reshape(w[CONIC_LAM_A], g.size());
253  ret_out[NLPSOL_LAM_P] = MX::nan(p.sparsity());
254 
255  Dict fun_opts;
256  fun_opts["default_in"] = nlpsol_default_in();
257  fun_opts["error_on_fail"] = error_on_fail;
258 
259  return Function(name, ret_in, ret_out, nlpsol_in(), nlpsol_out(), fun_opts);
260  }
261 
262  Function qpsol(const std::string& name, const std::string& solver,
263  const SXDict& qp, const Dict& opts) {
264  return qpsol_nlp(name, solver, qp, opts);
265  }
266 
267  Function qpsol(const std::string& name, const std::string& solver,
268  const MXDict& qp, const Dict& opts) {
269  return qpsol_nlp(name, solver, qp, opts);
270  }
271 
272  // Constructor
273  Conic::Conic(const std::string& name, const std::map<std::string, Sparsity> &st)
274  : FunctionInternal(name) {
275 
276  // Set default options
277  error_on_fail_ = true;
278 
279  P_ = Sparsity(0, 0);
280  for (auto i=st.begin(); i!=st.end(); ++i) {
281  if (i->first=="a") {
282  A_ = i->second;
283  } else if (i->first=="h") {
284  H_ = i->second;
285  } else if (i->first=="q") {
286  Q_ = i->second;
287  } else if (i->first=="p") {
288  P_ = i->second;
289  } else {
290  casadi_error("Unrecognized field in QP structure: " + str(i->first));
291  }
292  }
293 
294  // We need either A or H
295  casadi_assert(!A_.is_null() || !H_.is_null(),
296  "Cannot determine dimension");
297 
298  // Generate A or H
299  if (A_.is_null()) {
300  A_ = Sparsity(0, H_.size2());
301  } else if (H_.is_null()) {
302  H_ = Sparsity(A_.size2(), A_.size2());
303  } else {
304  // Consistency check
305  casadi_assert(A_.size2()==H_.size2(),
306  "Got incompatible dimensions.\n"
307  "min x'Hx + G'x s.t. LBA <= Ax <= UBA :\n"
308  "H: " + H_.dim() + " - A: " + A_.dim() + "\n"
309  "We need: H.size2()==A.size2()");
310  }
311 
312  casadi_assert(H_.is_symmetric(),
313  "Got incompatible dimensions. min x'Hx + G'x\n"
314  "H: " + H_.dim() +
315  "We need H square & symmetric");
316 
317  nx_ = A_.size2();
318  na_ = A_.size1();
319 
320  // Check psd constraints (linear part)
321  if (Q_.is_null()) {
322  Q_ = Sparsity(0, 0);
323  np_ = 0;
324  } else {
325  casadi_assert(Q_.size2()==nx_,
326  "Got incompatible dimensions.\n"
327  "Q: " + Q_.dim() +
328  "We need the product Qx to exist.");
329  np_ = static_cast<casadi_int>(sqrt(static_cast<double>(Q_.size1())));
330  casadi_assert(np_*np_==Q_.size1(),
331  "Got incompatible dimensions.\n"
332  "Q: " + Q_.dim() +
333  "We need Q.size1() to have an integer square root.");
334 
335  Sparsity qsum = reshape(sum2(Q_), np_, np_);
336 
337  casadi_assert(qsum.is_symmetric(),
338  "Got incompatible dimensions.");
339  }
340 
341  // Check psd constraints (constant part)
342  if (P_.is_null()) P_ = Sparsity(np_, np_);
343 
344  casadi_assert(P_.is_symmetric(),
345  "Got incompatible dimensions.\n"
346  "P: " + P_.dim() +
347  "We need P square & symmetric.");
348 
349  casadi_assert(P_.size1()==np_,
350  "Got incompatible dimensions.\n"
351  "P: " + P_.dim() +
352  "We need P " + str(np_) + "-by-" + str(np_) + ".");
353 
354  }
355 
357  switch (static_cast<ConicInput>(i)) {
358  case CONIC_X0:
359  case CONIC_G:
360  case CONIC_LBX:
361  case CONIC_UBX:
362  case CONIC_LAM_X0:
363  return get_sparsity_out(CONIC_X);
364  case CONIC_Q:
365  return Q_;
366  case CONIC_P:
367  return P_;
368  case CONIC_LBA:
369  case CONIC_UBA:
370  case CONIC_LAM_A0:
372  case CONIC_A:
373  return A_;
374  case CONIC_H:
375  return H_;
376  case CONIC_NUM_IN: break;
377  }
378  return Sparsity();
379  }
380 
382  switch (static_cast<ConicOutput>(i)) {
383  case CONIC_COST:
384  return Sparsity::scalar();
385  case CONIC_X:
386  case CONIC_LAM_X:
387  return Sparsity::dense(nx_, 1);
388  case CONIC_LAM_A:
389  return Sparsity::dense(na_, 1);
390  case CONIC_NUM_OUT: break;
391  }
392  return Sparsity();
393  }
394 
397  {{"discrete",
398  {OT_BOOLVECTOR,
399  "Indicates which of the variables are discrete, i.e. integer-valued"}},
400  {"equality",
401  {OT_BOOLVECTOR,
402  "Indicate an upfront hint which of the constraints are equalities. "
403  "Some solvers may be able to exploit this knowledge. "
404  "When true, the corresponding lower and upper bounds are assumed equal. "
405  "When false, the corresponding bounds may be equal or different."}},
406  {"print_problem",
407  {OT_BOOL,
408  "Print a numeric description of the problem"}}
409  }
410  };
411 
412  void Conic::init(const Dict& opts) {
413  // Call the init method of the base class
415 
416  print_problem_ = false;
417 
418  // Read options
419  for (auto&& op : opts) {
420  if (op.first=="discrete") {
421  discrete_ = op.second;
422  } else if (op.first=="equality") {
423  equality_ = op.second;
424  } else if (op.first=="print_problem") {
425  print_problem_ = op.second;
426  }
427  }
428 
429  // Check options
430  if (!discrete_.empty()) {
431  casadi_assert(discrete_.size()==nx_, "\"discrete\" option has wrong length");
432  if (std::find(discrete_.begin(), discrete_.end(), true)!=discrete_.end()) {
433  casadi_assert(integer_support(),
434  "Discrete variables require a solver with integer support");
435  }
436  }
437 
438  if (!equality_.empty()) {
439  casadi_assert(equality_.size()==na_, "\"equality\" option has wrong length. "
440  "Expected " + str(na_) + " elements, but got " +
441  str(equality_.size()) + " instead.");
442  }
443 
444  casadi_assert(np_==0 || psd_support(),
445  "Selected solver does not support psd constraints.");
446 
447  set_qp_prob();
448  }
449 
451  int Conic::init_mem(void* mem) const {
452  if (ProtoFunction::init_mem(mem)) return 1;
453 
454  return 0;
455  }
456 
458  void Conic::set_work(void* mem, const double**& arg, double**& res,
459  casadi_int*& iw, double*& w) const {
460 
461  auto m = static_cast<ConicMemory*>(mem);
462 
463  casadi_qp_data<double>& d_qp = m->d_qp;
464  d_qp.h = arg[CONIC_H];
465  d_qp.g = arg[CONIC_G];
466  d_qp.a = arg[CONIC_A];
467  d_qp.lbx = arg[CONIC_LBX];
468  d_qp.ubx = arg[CONIC_UBX];
469  d_qp.uba = arg[CONIC_UBA];
470  d_qp.lba = arg[CONIC_LBA];
471  d_qp.x0 = arg[CONIC_X0];
472  d_qp.lam_x0 = arg[CONIC_LAM_X0];
473  d_qp.lam_a0 = arg[CONIC_LAM_A0];
474  d_qp.x = res[CONIC_X];
475  d_qp.lam_x = res[CONIC_LAM_X];
476  d_qp.lam_a = res[CONIC_LAM_A];
477  d_qp.f = res[CONIC_COST];
478 
479  // Problem has not been solved at this point
480  d_qp.success = false;
482  d_qp.iter_count = -1;
483  }
484 
486  }
487 
488  void Conic::check_inputs(const double* lbx, const double* ubx,
489  const double* lba, const double* uba) const {
490  for (casadi_int i=0; i<nx_; ++i) {
491  double lb = lbx ? lbx[i] : 0., ub = ubx ? ubx[i] : 0.;
492  casadi_assert(lb <= ub && lb!=inf && ub!=-inf,
493  "Ill-posed problem detected: "
494  "LBX[" + str(i) + "] <= UBX[" + str(i) + "] was violated. "
495  "Got LBX[" + str(i) + "]=" + str(lb) + " and UBX[" + str(i) + "] = " + str(ub) + ".");
496  }
497  for (casadi_int i=0; i<na_; ++i) {
498  double lb = lba ? lba[i] : 0., ub = uba ? uba[i] : 0.;
499  casadi_assert(lb <= ub && lb!=inf && ub!=-inf,
500  "Ill-posed problem detected: "
501  "LBA[" + str(i) + "] <= UBA[" + str(i) + "] was violated. "
502  "Got LBA[" + str(i) + "] = " + str(lb) + " and UBA[" + str(i) + "] = " + str(ub) + ".");
503  }
504  }
505 
506  void Conic::generateNativeCode(std::ostream& file) const {
507  casadi_error("generateNativeCode not defined for class " + class_name());
508  }
509 
510  std::map<std::string, Conic::Plugin> Conic::solvers_;
511 
512 #ifdef CASADI_WITH_THREADSAFE_SYMBOLICS
513  std::mutex Conic::mutex_solvers_;
514 #endif // CASADI_WITH_THREADSAFE_SYMBOLICS
515 
516  const std::string Conic::infix_ = "conic";
517 
518  double Conic::get_default_in(casadi_int ind) const {
519  switch (ind) {
520  case CONIC_LBX:
521  case CONIC_LBA:
522  return -std::numeric_limits<double>::infinity();
523  case CONIC_UBX:
524  case CONIC_UBA:
525  return std::numeric_limits<double>::infinity();
526  default:
527  return 0;
528  }
529  }
530 
532  eval(const double** arg, double** res, casadi_int* iw, double* w, void* mem) const {
533  if (print_problem_) {
534  uout() << "H:";
535  DM::print_dense(uout(), H_, arg[CONIC_H], false);
536  uout() << std::endl;
537  uout() << "G:" << std::vector<double>(arg[CONIC_G], arg[CONIC_G]+nx_) << std::endl;
538  uout() << "A:";
539  DM::print_dense(uout(), A_, arg[CONIC_A], false);
540  uout() << std::endl;
541  uout() << "lba:" << std::vector<double>(arg[CONIC_LBA], arg[CONIC_LBA]+na_) << std::endl;
542  uout() << "uba:" << std::vector<double>(arg[CONIC_UBA], arg[CONIC_UBA]+na_) << std::endl;
543  uout() << "lbx:" << std::vector<double>(arg[CONIC_LBX], arg[CONIC_LBX]+nx_) << std::endl;
544  uout() << "ubx:" << std::vector<double>(arg[CONIC_UBX], arg[CONIC_UBX]+nx_) << std::endl;
545  }
546  auto m = static_cast<ConicMemory*>(mem);
547 
548  if (inputs_check_) {
549  check_inputs(arg[CONIC_LBX], arg[CONIC_UBX], arg[CONIC_LBA], arg[CONIC_UBA]);
550  }
551 
552  setup(mem, arg, res, iw, w);
553 
554  int ret = solve(arg, res, iw, w, mem);
555 
556  if (m->d_qp.success) m->d_qp.unified_return_status = SOLVER_RET_SUCCESS;
557 
558  if (error_on_fail_ && !m->d_qp.success)
559  casadi_error("conic process failed. "
560  "Set 'error_on_fail' option to false to ignore this error.");
561  return ret;
562  }
563 
564  std::vector<std::string> conic_options(const std::string& name) {
565  return Conic::plugin_options(name).all();
566  }
567 
568  std::string conic_option_type(const std::string& name, const std::string& op) {
569  return Conic::plugin_options(name).type(op);
570  }
571 
572  std::string conic_option_info(const std::string& name, const std::string& op) {
573  return Conic::plugin_options(name).info(op);
574  }
575 
576  bool Conic::is_a(const std::string& type, bool recursive) const {
577  return type=="Conic" || (recursive && FunctionInternal::is_a(type, recursive));
578  }
579 
581 
582  Sparsity qsum = reshape(sum2(Q_), np_, np_);
583 
584  // Block detection
585  Sparsity aggregate = qsum+P_;
586 
587  std::vector<casadi_int> p;
588  casadi_int nb = aggregate.scc(p, mem.r);
589 
590  std::string pattern_message = "Pattern not recognised";
591 
592  casadi_assert(p==range(p.size()), pattern_message);
593 
594  const casadi_int* row = aggregate.row();
595  const casadi_int* colind = aggregate.colind();
596 
597  // Check fishbone-structure
598  for (casadi_int i=0;i<nb;++i) {
599  casadi_int block_size = mem.r[i+1]-mem.r[i];
600  // number of nonzeros in column mem.r[i+1]-1
601  casadi_int nz = colind[mem.r[i+1]]-colind[mem.r[i+1]-1];
602  // Last column of block should be dense
603  casadi_assert(nz==block_size, pattern_message);
604  for (casadi_int k=0;k<block_size-1;++k) {
605  casadi_assert(colind[mem.r[i]+k+1]-colind[mem.r[i]+k], pattern_message);
606  casadi_assert(*(row++)==k+mem.r[i], pattern_message);
607  casadi_assert(*(row++)==mem.r[i]+block_size-1, pattern_message);
608  }
609 
610  for (casadi_int k=0;k<block_size;++k)
611  casadi_assert(*(row++)==k+mem.r[i], pattern_message);
612  }
613 
633  /*
634 
635  Aggregate pattern:
636 
637  x (x)
638  x(x)
639  xx(x)
640  x (x)
641  x(x)
642  xx(x)
643 
644  We are interested in the parts in parenthesis (target).
645 
646  Find out which rows in Q correspond to those targets
647  */
648 
649  // Lookup vector for target start and end
650  std::vector<casadi_int> target_start(nb), target_stop(nb);
651  for (casadi_int i=0;i<nb;++i) {
652  target_start[i] = (mem.r[i+1]-1)*aggregate.size1()+mem.r[i];
653  target_stop[i] = target_start[i]+mem.r[i+1]-mem.r[i];
654  }
655 
656  // Collect the nonzero indices in Q that that correspond to the target area
657  std::vector<casadi_int> q_nz;
658  // Triplet form for map_Q sparsity
659  std::vector<casadi_int> q_row, q_col;
660 
661  // Loop over Q's columns (decision variables)
662  for (casadi_int j=0; j<Q_.size2(); ++j) {
663  casadi_int block_index = 0;
664  // Loop over Q's rows
665  for (casadi_int k=Q_.colind(j); k<Q_.colind(j+1); ++k) {
666  casadi_int i = Q_.row(k);
667 
668  // Increment block_index if i runs ahead
669  while (i>target_stop[block_index] && block_index<nb-1) block_index++;
670 
671  if (i>=target_start[block_index] && i<target_stop[block_index]) {
672  // Got a nonzero in the target region
673  q_nz.push_back(k);
674  q_row.push_back(mem.r[block_index]+i-target_start[block_index]);
675  q_col.push_back(j);
676  }
677  }
678  }
679 
680  mem.map_Q = IM::triplet(q_row, q_col, q_nz, mem.r[nb], nx_);
681 
682  // Add the [X1;Z1;X2;Z2;...] part
683  mem.map_Q = horzcat(mem.map_Q, -IM::eye(mem.r[nb])).T();
684 
685  // Get maximum nonzero count of any column
686  casadi_int max_nnz = 0;
687  for (casadi_int i=0;i<mem.map_Q.size2();++i) {
688  max_nnz = std::max(max_nnz, mem.map_Q.colind(i+1)-mem.map_Q.colind(i));
689  }
690 
691  // ind/val size needs to cover max nonzero count
692  mem.indval_size = std::max(nx_, max_nnz);
693 
694  // Collect the indices for the P target area
695  mem.map_P.resize(mem.r[nb], -1);
696  for (casadi_int i=0;i<nb;++i) {
697  for (casadi_int k=P_.colind(mem.r[i+1]-1); k<P_.colind(mem.r[i+1]); ++k) {
698  casadi_int r = P_.row(k);
699  mem.map_P[r] = k;
700  }
701  }
702 
703  // ind/val size needs to cover blocksize
704  for (casadi_int i=0;i<nb;++i)
705  mem.indval_size = std::max(mem.indval_size, mem.r[i+1]-mem.r[i]);
706 
707  // Get the transpose and mapping
708  mem.AT = A_.transpose(mem.A_mapping);
709  }
710 
711  Dict Conic::get_stats(void* mem) const {
712  Dict stats = FunctionInternal::get_stats(mem);
713  auto m = static_cast<ConicMemory*>(mem);
714 
715  stats["success"] = m->d_qp.success;
716  stats["unified_return_status"] = string_from_UnifiedReturnStatus(m->d_qp.unified_return_status);
717  stats["iter_count"] = m->d_qp.iter_count;
718  return stats;
719  }
720 
722  s.pack("Conic::SDPToSOCPMem::r", m.r);
723  s.pack("Conic::SDPToSOCPMem::AT", m.AT);
724  s.pack("Conic::SDPToSOCPMem::A_mapping", m.A_mapping);
725  s.pack("Conic::SDPToSOCPMem::map_Q", m.map_Q);
726  s.pack("Conic::SDPToSOCPMem::map_P", m.map_P);
727  s.pack("Conic::SDPToSOCPMem::indval_size", m.indval_size);
728  }
730  s.unpack("Conic::SDPToSOCPMem::r", m.r);
731  s.unpack("Conic::SDPToSOCPMem::AT", m.AT);
732  s.unpack("Conic::SDPToSOCPMem::A_mapping", m.A_mapping);
733  s.unpack("Conic::SDPToSOCPMem::map_Q", m.map_Q);
734  s.unpack("Conic::SDPToSOCPMem::map_P", m.map_P);
735  s.unpack("Conic::SDPToSOCPMem::indval_size", m.indval_size);
736  }
737 
740 
741  s.version("Conic", 3);
742  s.pack("Conic::discrete", discrete_);
743  s.pack("Conic::equality", equality_);
744  s.pack("Conic::print_problem", print_problem_);
745  s.pack("Conic::H", H_);
746  s.pack("Conic::A", A_);
747  s.pack("Conic::Q", Q_);
748  s.pack("Conic::P", P_);
749  s.pack("Conic::nx", nx_);
750  s.pack("Conic::na", na_);
751  s.pack("Conic::np", np_);
752  }
753 
757  }
758 
761  }
762 
764  int version = s.version("Conic", 1, 3);
765  s.unpack("Conic::discrete", discrete_);
766  if (version>=3) {
767  s.unpack("Conic::equality", equality_);
768  }
769  s.unpack("Conic::print_problem", print_problem_);
770  if (version==1) {
771  s.unpack("Conic::error_on_fail", error_on_fail_);
772  }
773  s.unpack("Conic::H", H_);
774  s.unpack("Conic::A", A_);
775  set_qp_prob();
776  s.unpack("Conic::Q", Q_);
777  s.unpack("Conic::P", P_);
778  s.unpack("Conic::nx", nx_);
779  s.unpack("Conic::na", na_);
780  s.unpack("Conic::np", np_);
781  }
782 
783  void Conic::set_qp_prob() {
784  p_qp_.sp_a = A_;
785  p_qp_.sp_h = H_;
786  casadi_qp_setup(&p_qp_);
787  }
788 
791  g.local("d_qp", "struct casadi_qp_data");
792  g.local("p_qp", "struct casadi_qp_prob");
793 
794  g << "d_qp.prob = &p_qp;\n";
795  g << "p_qp.sp_a = " << g.sparsity(A_) << ";\n";
796  g << "p_qp.sp_h = " << g.sparsity(H_) << ";\n";
797  g << "casadi_qp_setup(&p_qp);\n";
798  g << "casadi_qp_init(&d_qp, &iw, &w);\n";
799 
800 
801  g << "d_qp.h = arg[" << CONIC_H << "];\n";
802  g << "d_qp.g = arg[" << CONIC_G << "];\n";
803  g << "d_qp.a = arg[" << CONIC_A << "];\n";
804  g << "d_qp.lbx = arg[" << CONIC_LBX << "];\n";
805  g << "d_qp.ubx = arg[" << CONIC_UBX << "];\n";
806  g << "d_qp.lba = arg[" << CONIC_LBA << "];\n";
807  g << "d_qp.uba = arg[" << CONIC_UBA << "];\n";
808  g << "d_qp.x0 = arg[" << CONIC_X0 << "];\n";
809  g << "d_qp.lam_x0 = arg[" << CONIC_LAM_X0 << "];\n";
810  g << "d_qp.lam_a0 = arg[" << CONIC_LAM_A0 << "];\n";
811 
812  g << "d_qp.f = res[" << CONIC_COST << "];\n";
813  g << "d_qp.x = res[" << CONIC_X << "];\n";
814  g << "d_qp.lam_x = res[" << CONIC_LAM_X << "];\n";
815  g << "d_qp.lam_a = res[" << CONIC_LAM_A << "];\n";
816  }
817 
818 } // namespace casadi
Helper class for C code generation.
void local(const std::string &name, const std::string &type, const std::string &ref="")
Declare a local variable.
std::string sparsity(const Sparsity &sp, bool canonical=true)
void add_auxiliary(Auxiliary f, const std::vector< std::string > &inst={"casadi_real"})
Add a built-in auxiliary function.
Internal class.
Definition: conic_impl.hpp:44
static const Options options_
Options.
Definition: conic_impl.hpp:83
casadi_int nx_
Number of decision variables.
Definition: conic_impl.hpp:169
bool is_a(const std::string &type, bool recursive) const override
Check if the function is of a particular type.
Definition: conic.cpp:576
double get_default_in(casadi_int ind) const override
Get default input value.
Definition: conic.cpp:518
virtual int solve(const double **arg, double **res, casadi_int *iw, double *w, void *mem) const =0
Solve the QP.
int init_mem(void *mem) const override
Initalize memory block.
Definition: conic.cpp:451
casadi_int np_
The shape of psd constraint matrix.
Definition: conic_impl.hpp:175
Sparsity get_sparsity_in(casadi_int i) override
Sparsities of function inputs and outputs.
Definition: conic.cpp:356
casadi_int na_
The number of constraints (counting both equality and inequality) == A.size1()
Definition: conic_impl.hpp:172
static std::map< std::string, Plugin > solvers_
Collection of solvers.
Definition: conic_impl.hpp:123
virtual void generateNativeCode(std::ostream &file) const
Definition: conic.cpp:506
virtual void check_inputs(const double *lbx, const double *ubx, const double *lba, const double *uba) const
Check if the numerical values of the supplied bounds make sense.
Definition: conic.cpp:488
Sparsity H_
Problem structure.
Definition: conic_impl.hpp:166
void init(const Dict &opts) override
Initialize.
Definition: conic.cpp:412
void deserialize(DeserializingStream &s, SDPToSOCPMem &m)
Definition: conic.cpp:729
std::vector< bool > discrete_
Options.
Definition: conic_impl.hpp:161
std::vector< bool > equality_
Definition: conic_impl.hpp:162
Conic(const std::string &name, const std::map< std::string, Sparsity > &st)
Definition: conic.cpp:273
Sparsity get_sparsity_out(casadi_int i) override
Sparsities of function inputs and outputs.
Definition: conic.cpp:381
~Conic() override=0
Definition: conic.cpp:485
void serialize_body(SerializingStream &s) const override
Serialize an object without type information.
Definition: conic.cpp:738
void set_work(void *mem, const double **&arg, double **&res, casadi_int *&iw, double *&w) const override
Set the (persistent) work vectors.
Definition: conic.cpp:458
virtual bool psd_support() const
Can psd constraints be treated.
Definition: conic_impl.hpp:149
Dict get_stats(void *mem) const override
Get all statistics.
Definition: conic.cpp:711
static const std::string infix_
Infix.
Definition: conic_impl.hpp:130
casadi_qp_prob< double > p_qp_
Definition: conic_impl.hpp:47
void serialize_type(SerializingStream &s) const override
Serialize type information.
Definition: conic.cpp:754
void qp_codegen_body(CodeGenerator &g) const
Generate code for the function body.
Definition: conic.cpp:789
void sdp_to_socp_init(SDPToSOCPMem &mem) const
SDP to SOCP conversion initialization.
Definition: conic.cpp:580
virtual bool integer_support() const
Can discrete variables be treated.
Definition: conic_impl.hpp:146
int eval(const double **arg, double **res, casadi_int *iw, double *w, void *mem) const final
Solve the QP.
Definition: conic.cpp:532
void serialize(SerializingStream &s, const SDPToSOCPMem &m) const
Definition: conic.cpp:721
Helper class for Serialization.
void unpack(Sparsity &e)
Reconstruct an object from the input stream.
void version(const std::string &name, int v)
static void open(std::ofstream &, const std::string &path, std::ios_base::openmode mode=std::ios_base::out)
Definition: filesystem.cpp:115
Internal class for Function.
Dict get_stats(void *mem) const override
Get all statistics.
void init(const Dict &opts) override
Initialize.
void serialize_body(SerializingStream &s) const override
Serialize an object without type information.
virtual bool is_a(const std::string &type, bool recursive) const
Check if the function is of a particular type.
bool inputs_check_
Errors are thrown if numerical values of inputs look bad.
static const Options options_
Options.
void serialize_type(SerializingStream &s) const override
Serialize type information.
void setup(void *mem, const double **arg, double **res, casadi_int *iw, double *w) const
Set the (persistent and temporary) work vectors.
static std::string string_from_UnifiedReturnStatus(UnifiedReturnStatus status)
Function object.
Definition: function.hpp:60
FunctionInternal * get() const
Definition: function.cpp:353
Function expand() const
Expand a function to SX.
Definition: function.cpp:308
const std::vector< std::string > & name_in() const
Get input scheme.
Definition: function.cpp:961
static Function create(FunctionInternal *node)
Create from node.
Definition: function.cpp:336
const SX sx_in(casadi_int iind) const
Get symbolic primitives equivalent to the input expressions.
Definition: function.cpp:1552
casadi_int n_out() const
Get the number of function outputs.
Definition: function.cpp:823
casadi_int n_in() const
Get the number of function inputs.
Definition: function.cpp:819
const std::vector< std::string > & name_out() const
Get output scheme.
Definition: function.cpp:965
casadi_int size2() const
Get the second dimension (i.e. number of columns)
static MX sym(const std::string &name, casadi_int nrow=1, casadi_int ncol=1)
Create an nrow-by-ncol symbolic primitive.
const casadi_int * colind() const
Get the sparsity pattern. See the Sparsity class for details.
bool is_null() const
Is a null pointer?
MX - Matrix expression.
Definition: mx.hpp:92
static MX nan(const Sparsity &sp)
create a matrix with all nan
Definition: mx.cpp:576
Matrix< Scalar > T() const
Transpose the matrix.
static Matrix< casadi_int > triplet(const std::vector< casadi_int > &row, const std::vector< casadi_int > &col, const Matrix< casadi_int > &d)
Construct a sparse matrix from triplet form.
void print_dense(std::ostream &stream, bool truncate=true) const
Print dense matrix-stype.
static Matrix< casadi_int > eye(casadi_int n)
create an n-by-n identity matrix
static bool has_plugin(const std::string &pname, bool verbose=false)
Check if a plugin is available or can be loaded.
static Conic * instantiate(const std::string &fname, const std::string &pname, Problem problem)
void serialize_type(SerializingStream &s) const
Serialize type information.
static const Options & plugin_options(const std::string &pname)
Get the plugin options.
static Plugin & getPlugin(const std::string &pname)
Load and get the creator function.
static ProtoFunction * deserialize(DeserializingStream &s)
Deserialize with type disambiguation.
static Plugin load_plugin(const std::string &pname, bool register_plugin=true, bool needs_lock=true)
Load a plugin dynamically.
Base class for FunctionInternal and LinsolInternal.
bool error_on_fail_
Throw an exception on failure?
virtual int init_mem(void *mem) const
Initalize memory block.
Helper class for Serialization.
void version(const std::string &name, int v)
void pack(const Sparsity &e)
Serializes an object to the output stream.
virtual std::string class_name() const =0
Readable name of the internal class.
General sparsity class.
Definition: sparsity.hpp:106
casadi_int size1() const
Get the number of rows.
Definition: sparsity.cpp:124
casadi_int scc(std::vector< casadi_int > &index, std::vector< casadi_int > &offset) const
Find the strongly connected components of the bigraph defined by the sparsity pattern.
Definition: sparsity.cpp:703
std::string dim(bool with_nz=false) const
Get the dimension as a string.
Definition: sparsity.cpp:587
static Sparsity dense(casadi_int nrow, casadi_int ncol=1)
Create a dense rectangular sparsity pattern *.
Definition: sparsity.cpp:1012
Sparsity transpose(std::vector< casadi_int > &mapping, bool invert_mapping=false) const
Transpose the matrix and get the reordering of the non-zero entries.
Definition: sparsity.cpp:390
casadi_int size2() const
Get the number of columns.
Definition: sparsity.cpp:128
const casadi_int * row() const
Get a reference to row-vector,.
Definition: sparsity.cpp:164
static Sparsity scalar(bool dense_scalar=true)
Create a scalar sparsity pattern *.
Definition: sparsity.hpp:153
const casadi_int * colind() const
Get a reference to the colindex of all column element (see class description)
Definition: sparsity.cpp:168
bool is_symmetric() const
Is symmetric?
Definition: sparsity.cpp:317
Function qpsol(const std::string &name, const std::string &solver, const SXDict &qp, const Dict &opts)
Definition: conic.cpp:262
std::string conic_option_info(const std::string &name, const std::string &op)
Get documentation for a particular option.
Definition: conic.cpp:572
Function conic(const std::string &name, const std::string &solver, const SpDict &qp, const Dict &opts)
Definition: conic.cpp:44
std::string conic_option_type(const std::string &name, const std::string &op)
Get type info for a particular option.
Definition: conic.cpp:568
bool has_conic(const std::string &name)
Check if a particular plugin is available.
Definition: conic.cpp:32
std::string doc_conic(const std::string &name)
Get the documentation string for a plugin.
Definition: conic.cpp:40
std::vector< std::string > conic_options(const std::string &name)
Get all options for a plugin.
Definition: conic.cpp:564
std::vector< std::string > conic_out()
Get QP solver output scheme of QP solvers.
Definition: conic.cpp:67
casadi_int conic_n_in()
Get the number of QP solver inputs.
Definition: conic.cpp:103
casadi_int conic_n_out()
Get the number of QP solver outputs.
Definition: conic.cpp:107
void load_conic(const std::string &name)
Explicitly load a plugin dynamically.
Definition: conic.cpp:36
std::vector< std::string > conic_in()
Get input scheme of QP solvers.
Definition: conic.cpp:61
void conic_debug(const Function &f, const std::string &filename)
Definition: conic.cpp:49
std::vector< std::string > nlpsol_in()
Get input scheme of NLP solvers.
Definition: nlpsol.cpp:200
std::vector< std::string > nlpsol_out()
Get NLP solver output scheme of NLP solvers.
Definition: nlpsol.cpp:206
double nlpsol_default_in(casadi_int ind)
Default input for an NLP solver.
Definition: nlpsol.cpp:212
The casadi namespace.
Definition: archiver.cpp:28
@ NLPSOL_P
Value of fixed parameters (np x 1)
Definition: nlpsol.hpp:198
@ NLPSOL_UBX
Decision variables upper bound (nx x 1), default +inf.
Definition: nlpsol.hpp:202
@ NLPSOL_X0
Decision variables, initial guess (nx x 1)
Definition: nlpsol.hpp:196
@ NLPSOL_LAM_G0
Lagrange multipliers for bounds on G, initial guess (ng x 1)
Definition: nlpsol.hpp:210
@ NLPSOL_UBG
Constraints upper bound (ng x 1), default +inf.
Definition: nlpsol.hpp:206
@ NLPSOL_LAM_X0
Lagrange multipliers for bounds on X, initial guess (nx x 1)
Definition: nlpsol.hpp:208
@ NLPSOL_NUM_IN
Definition: nlpsol.hpp:211
@ NLPSOL_LBG
Constraints lower bound (ng x 1), default -inf.
Definition: nlpsol.hpp:204
@ NLPSOL_LBX
Decision variables lower bound (nx x 1), default -inf.
Definition: nlpsol.hpp:200
std::map< std::string, MX > MXDict
Definition: mx.hpp:1009
std::vector< casadi_int > range(casadi_int start, casadi_int stop, casadi_int step, casadi_int len)
Range function.
@ NLPSOL_G
Constraints function at the optimal solution (ng x 1)
Definition: nlpsol.hpp:221
@ NLPSOL_X
Decision variables at the optimal solution (nx x 1)
Definition: nlpsol.hpp:217
@ NLPSOL_NUM_OUT
Definition: nlpsol.hpp:228
@ NLPSOL_LAM_P
Lagrange multipliers for bounds on P at the solution (np x 1)
Definition: nlpsol.hpp:227
@ NLPSOL_F
Cost function value at the optimal solution (1 x 1)
Definition: nlpsol.hpp:219
@ NLPSOL_LAM_G
Lagrange multipliers for bounds on G at the solution (ng x 1)
Definition: nlpsol.hpp:225
@ NLPSOL_LAM_X
Lagrange multipliers for bounds on X at the solution (nx x 1)
Definition: nlpsol.hpp:223
T get_from_dict(const std::map< std::string, T > &d, const std::string &key, const T &default_value)
@ NL_X
Decision variable.
Definition: nlpsol.hpp:170
@ NL_P
Fixed parameter.
Definition: nlpsol.hpp:172
@ NL_NUM_IN
Number of NLP inputs.
Definition: nlpsol.hpp:174
ConicInput
Input arguments of a QP problem.
Definition: conic.hpp:170
@ CONIC_UBA
dense, (nc x 1)
Definition: conic.hpp:181
@ CONIC_X0
dense, (n x 1)
Definition: conic.hpp:187
@ CONIC_A
The matrix A: sparse, (nc x n) - product with x must be dense.
Definition: conic.hpp:177
@ CONIC_G
The vector g: dense, (n x 1)
Definition: conic.hpp:175
@ CONIC_Q
The matrix Q: sparse symmetric, (np^2 x n)
Definition: conic.hpp:193
@ CONIC_LBA
dense, (nc x 1)
Definition: conic.hpp:179
@ CONIC_UBX
dense, (n x 1)
Definition: conic.hpp:185
@ CONIC_H
Definition: conic.hpp:173
@ CONIC_LAM_A0
dense
Definition: conic.hpp:191
@ CONIC_LBX
dense, (n x 1)
Definition: conic.hpp:183
@ CONIC_NUM_IN
Definition: conic.hpp:196
@ CONIC_LAM_X0
dense
Definition: conic.hpp:189
@ CONIC_P
The matrix P: sparse symmetric, (np x np)
Definition: conic.hpp:195
@ OT_BOOLVECTOR
std::map< std::string, SX > SXDict
Definition: sx_fwd.hpp:40
void extract_from_dict_inplace(Dict &d, const std::string &key, T &value)
std::string str(const T &v)
String representation, any type.
GenericType::Dict Dict
C++ equivalent of Python's dict or MATLAB's struct.
const double inf
infinity
Definition: calculus.hpp:50
std::map< std::string, Sparsity > SpDict
Definition: sparsity.hpp:1502
std::ostream & uout()
Function qpsol_nlp(const std::string &name, const std::string &solver, const std::map< std::string, M > &qp, const Dict &opts)
Definition: conic.cpp:112
std::string filename(const std::string &path)
Definition: ghc.cpp:55
@ SOLVER_RET_SUCCESS
@ SOLVER_RET_UNKNOWN
ConicOutput
Output arguments of an QP Solver.
Definition: conic.hpp:199
@ CONIC_X
The primal solution.
Definition: conic.hpp:201
@ CONIC_LAM_A
The dual solution corresponding to linear bounds.
Definition: conic.hpp:205
@ CONIC_COST
The optimal cost.
Definition: conic.hpp:203
@ CONIC_LAM_X
The dual solution corresponding to simple bounds.
Definition: conic.hpp:207
@ CONIC_NUM_OUT
Definition: conic.hpp:208
casadi_qp_data< double > d_qp
Definition: conic_impl.hpp:39
SDP to SOCP conversion memory.
Definition: conic_impl.hpp:178
std::vector< casadi_int > r
Definition: conic_impl.hpp:180
std::vector< casadi_int > A_mapping
Definition: conic_impl.hpp:184
std::vector< casadi_int > map_P
Definition: conic_impl.hpp:190
Options metadata for a class.
Definition: options.hpp:40
std::string type(const std::string &name) const
Definition: options.cpp:289
std::vector< std::string > all() const
Definition: options.cpp:283
std::string info(const std::string &name) const
Definition: options.cpp:295
const T1 * h
Definition: casadi_qp.hpp:64
const T1 * lba
Definition: casadi_qp.hpp:64
casadi_int iter_count
Definition: casadi_qp.hpp:61
const T1 * lbx
Definition: casadi_qp.hpp:64
const T1 * uba
Definition: casadi_qp.hpp:64
const T1 * ubx
Definition: casadi_qp.hpp:64
const T1 * lam_x0
Definition: casadi_qp.hpp:64
const T1 * x0
Definition: casadi_qp.hpp:64
UnifiedReturnStatus unified_return_status
Definition: casadi_qp.hpp:57
const T1 * g
Definition: casadi_qp.hpp:64
const T1 * lam_a0
Definition: casadi_qp.hpp:64
const T1 * a
Definition: casadi_qp.hpp:64
const casadi_int * sp_h
Definition: casadi_qp.hpp:31
const casadi_int * sp_a
Definition: casadi_qp.hpp:31