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