quadprog.cpp 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622
  1. #include "quadprog.h"
  2. #include "../matlab_format.h"
  3. #include <vector>
  4. #include <iostream>
  5. #include <cstdio>
  6. /*
  7. FILE eiquadprog.hh
  8. NOTE: this is a modified of uQuadProg++ package, working with Eigen data structures.
  9. uQuadProg++ is itself a port made by Angelo Furfaro of QuadProg++ originally developed by
  10. Luca Di Gaspero, working with ublas data structures.
  11. The quadprog_solve() function implements the algorithm of Goldfarb and Idnani
  12. for the solution of a (convex) Quadratic Programming problem
  13. by means of a dual method.
  14. The problem is in the form:
  15. min 0.5 * x G x + g0 x
  16. s.t.
  17. CE^T x + ce0 = 0
  18. CI^T x + ci0 >= 0
  19. The matrix and vectors dimensions are as follows:
  20. G: n * n
  21. g0: n
  22. CE: n * p
  23. ce0: p
  24. CI: n * m
  25. ci0: m
  26. x: n
  27. The function will return the cost of the solution written in the x vector or
  28. std::numeric_limits::infinity() if the problem is infeasible. In the latter case
  29. the value of the x vector is not correct.
  30. References: D. Goldfarb, A. Idnani. A numerically stable dual method for solving
  31. strictly convex quadratic programs. Mathematical Programming 27 (1983) pp. 1-33.
  32. Notes:
  33. 1. pay attention in setting up the vectors ce0 and ci0.
  34. If the constraints of your problem are specified in the form
  35. A^T x = b and C^T x >= d, then you should set ce0 = -b and ci0 = -d.
  36. 2. The matrix G is modified within the function since it is used to compute
  37. the G = L^T L cholesky factorization for further computations inside the function.
  38. If you need the original matrix G you should make a copy of it and pass the copy
  39. to the function.
  40. The author will be grateful if the researchers using this software will
  41. acknowledge the contribution of this modified function and of Di Gaspero's
  42. original version in their research papers.
  43. LICENSE
  44. Copyright (2010) Gael Guennebaud
  45. Copyright (2008) Angelo Furfaro
  46. Copyright (2006) Luca Di Gaspero
  47. This file is a porting of QuadProg++ routine, originally developed
  48. by Luca Di Gaspero, exploiting uBlas data structures for vectors and
  49. matrices instead of native C++ array.
  50. uquadprog is free software; you can redistribute it and/or modify
  51. it under the terms of the GNU General Public License as published by
  52. the Free Software Foundation; either version 2 of the License, or
  53. (at your option) any later version.
  54. uquadprog is distributed in the hope that it will be useful,
  55. but WITHOUT ANY WARRANTY; without even the implied warranty of
  56. MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  57. GNU General Public License for more details.
  58. You should have received a copy of the GNU General Public License
  59. along with uquadprog; if not, write to the Free Software
  60. Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
  61. */
  62. #include <Eigen/Dense>
  63. IGL_INLINE bool igl::copyleft::quadprog(
  64. const Eigen::MatrixXd & G,
  65. const Eigen::VectorXd & g0,
  66. const Eigen::MatrixXd & CE,
  67. const Eigen::VectorXd & ce0,
  68. const Eigen::MatrixXd & CI,
  69. const Eigen::VectorXd & ci0,
  70. Eigen::VectorXd& x)
  71. {
  72. using namespace Eigen;
  73. typedef double Scalar;
  74. #ifdef TRACE_SOLVER
  75. const auto print_ivector= [](const char* name, const Eigen::MatrixXi & A, int /*n*/)
  76. {
  77. std::cout<<igl::matlab_format(A,name)<<std::endl;
  78. };
  79. const auto print_matrix = [](const char* name, const Eigen::MatrixXd & A, int /*n*/)
  80. {
  81. std::cout<<igl::matlab_format(A,name)<<std::endl;
  82. };
  83. const auto print_vector = [](const char* name, const Eigen::VectorXd & v, int /*n*/)
  84. {
  85. std::cout<<igl::matlab_format(v,name)<<std::endl;
  86. };
  87. #endif
  88. const auto distance = [](Scalar a, Scalar b)->Scalar
  89. {
  90. Scalar a1, b1, t;
  91. a1 = std::abs(a);
  92. b1 = std::abs(b);
  93. if (a1 > b1)
  94. {
  95. t = (b1 / a1);
  96. return a1 * std::sqrt(1.0 + t * t);
  97. }
  98. else
  99. if (b1 > a1)
  100. {
  101. t = (a1 / b1);
  102. return b1 * std::sqrt(1.0 + t * t);
  103. }
  104. return a1 * std::sqrt(2.0);
  105. };
  106. const auto compute_d = [](VectorXd &d, const MatrixXd& J, const VectorXd& np)
  107. {
  108. d = J.adjoint() * np;
  109. };
  110. const auto update_z =
  111. [](VectorXd& z, const MatrixXd& J, const VectorXd& d, int iq)
  112. {
  113. z = J.rightCols(z.size()-iq) * d.tail(d.size()-iq);
  114. };
  115. const auto update_r =
  116. [](const MatrixXd& R, VectorXd& r, const VectorXd& d, int iq)
  117. {
  118. r.head(iq) =
  119. R.topLeftCorner(iq,iq).triangularView<Upper>().solve(d.head(iq));
  120. };
  121. const auto add_constraint = [&distance](
  122. MatrixXd& R,
  123. MatrixXd& J,
  124. VectorXd& d,
  125. int& iq,
  126. double& R_norm)->bool
  127. {
  128. int n=J.rows();
  129. #ifdef TRACE_SOLVER
  130. std::cerr << "Add constraint " << iq << '/';
  131. #endif
  132. int j, k;
  133. double cc, ss, h, t1, t2, xny;
  134. /* we have to find the Givens rotation which will reduce the element
  135. d(j) to zero.
  136. if it is already zero we don't have to do anything, except of
  137. decreasing j */
  138. for (j = n - 1; j >= iq + 1; j--)
  139. {
  140. /* The Givens rotation is done with the matrix (cc cs, cs -cc).
  141. If cc is one, then element (j) of d is zero compared with element
  142. (j - 1). Hence we don't have to do anything.
  143. If cc is zero, then we just have to switch column (j) and column (j - 1)
  144. of J. Since we only switch columns in J, we have to be careful how we
  145. update d depending on the sign of gs.
  146. Otherwise we have to apply the Givens rotation to these columns.
  147. The i - 1 element of d has to be updated to h. */
  148. cc = d(j - 1);
  149. ss = d(j);
  150. h = distance(cc, ss);
  151. if (h == 0.0)
  152. continue;
  153. d(j) = 0.0;
  154. ss = ss / h;
  155. cc = cc / h;
  156. if (cc < 0.0)
  157. {
  158. cc = -cc;
  159. ss = -ss;
  160. d(j - 1) = -h;
  161. }
  162. else
  163. d(j - 1) = h;
  164. xny = ss / (1.0 + cc);
  165. for (k = 0; k < n; k++)
  166. {
  167. t1 = J(k,j - 1);
  168. t2 = J(k,j);
  169. J(k,j - 1) = t1 * cc + t2 * ss;
  170. J(k,j) = xny * (t1 + J(k,j - 1)) - t2;
  171. }
  172. }
  173. /* update the number of constraints added*/
  174. iq++;
  175. /* To update R we have to put the iq components of the d vector
  176. into column iq - 1 of R
  177. */
  178. R.col(iq-1).head(iq) = d.head(iq);
  179. #ifdef TRACE_SOLVER
  180. std::cerr << iq << std::endl;
  181. #endif
  182. if (std::abs(d(iq - 1)) <= std::numeric_limits<double>::epsilon() * R_norm)
  183. {
  184. // problem degenerate
  185. return false;
  186. }
  187. R_norm = std::max<double>(R_norm, std::abs(d(iq - 1)));
  188. return true;
  189. };
  190. const auto delete_constraint = [&distance](
  191. MatrixXd& R,
  192. MatrixXd& J,
  193. VectorXi& A,
  194. VectorXd& u,
  195. int p,
  196. int& iq,
  197. int l)
  198. {
  199. int n = R.rows();
  200. #ifdef TRACE_SOLVER
  201. std::cerr << "Delete constraint " << l << ' ' << iq;
  202. #endif
  203. int i, j, k, qq = -1;
  204. double cc, ss, h, xny, t1, t2;
  205. /* Find the index qq for active constraint l to be removed */
  206. for (i = p; i < iq; i++)
  207. if (A(i) == l)
  208. {
  209. qq = i;
  210. break;
  211. }
  212. /* remove the constraint from the active set and the duals */
  213. for (i = qq; i < iq - 1; i++)
  214. {
  215. A(i) = A(i + 1);
  216. u(i) = u(i + 1);
  217. R.col(i) = R.col(i+1);
  218. }
  219. A(iq - 1) = A(iq);
  220. u(iq - 1) = u(iq);
  221. A(iq) = 0;
  222. u(iq) = 0.0;
  223. for (j = 0; j < iq; j++)
  224. R(j,iq - 1) = 0.0;
  225. /* constraint has been fully removed */
  226. iq--;
  227. #ifdef TRACE_SOLVER
  228. std::cerr << '/' << iq << std::endl;
  229. #endif
  230. if (iq == 0)
  231. return;
  232. for (j = qq; j < iq; j++)
  233. {
  234. cc = R(j,j);
  235. ss = R(j + 1,j);
  236. h = distance(cc, ss);
  237. if (h == 0.0)
  238. continue;
  239. cc = cc / h;
  240. ss = ss / h;
  241. R(j + 1,j) = 0.0;
  242. if (cc < 0.0)
  243. {
  244. R(j,j) = -h;
  245. cc = -cc;
  246. ss = -ss;
  247. }
  248. else
  249. R(j,j) = h;
  250. xny = ss / (1.0 + cc);
  251. for (k = j + 1; k < iq; k++)
  252. {
  253. t1 = R(j,k);
  254. t2 = R(j + 1,k);
  255. R(j,k) = t1 * cc + t2 * ss;
  256. R(j + 1,k) = xny * (t1 + R(j,k)) - t2;
  257. }
  258. for (k = 0; k < n; k++)
  259. {
  260. t1 = J(k,j);
  261. t2 = J(k,j + 1);
  262. J(k,j) = t1 * cc + t2 * ss;
  263. J(k,j + 1) = xny * (J(k,j) + t1) - t2;
  264. }
  265. }
  266. };
  267. int i, k, l; /* indices */
  268. int ip, me, mi;
  269. int n=g0.size(); int p=ce0.size(); int m=ci0.size();
  270. MatrixXd R(G.rows(),G.cols()), J(G.rows(),G.cols());
  271. LLT<MatrixXd,Lower> chol(G.cols());
  272. VectorXd s(m+p), z(n), r(m + p), d(n), np(n), u(m + p);
  273. VectorXd x_old(n), u_old(m + p);
  274. #ifdef TRACE_SOLVER
  275. double f_value;
  276. #endif
  277. double psi, c1, c2, sum, ss, R_norm;
  278. const double inf = std::numeric_limits<double>::infinity();
  279. double t, t1, t2; /* t is the step length, which is the minimum of the partial step length t1
  280. * and the full step length t2 */
  281. VectorXi A(m + p), A_old(m + p), iai(m + p);
  282. int iq;
  283. std::vector<bool> iaexcl(m + p);
  284. me = p; /* number of equality constraints */
  285. mi = m; /* number of inequality constraints */
  286. /*
  287. * Preprocessing phase
  288. */
  289. /* compute the trace of the original matrix G */
  290. c1 = G.trace();
  291. /* decompose the matrix G in the form LL^T */
  292. chol.compute(G);
  293. /* initialize the matrix R */
  294. d.setZero();
  295. R.setZero();
  296. R_norm = 1.0; /* this variable will hold the norm of the matrix R */
  297. /* compute the inverse of the factorized matrix G^-1, this is the initial value for H */
  298. // J = L^-T
  299. J.setIdentity();
  300. J = chol.matrixU().solve(J);
  301. c2 = J.trace();
  302. #ifdef TRACE_SOLVER
  303. print_matrix("J", J, n);
  304. #endif
  305. /* c1 * c2 is an estimate for cond(G) */
  306. /*
  307. * Find the unconstrained minimizer of the quadratic form 0.5 * x G x + g0 x
  308. * this is a feasible point in the dual space
  309. * x = G^-1 * g0
  310. */
  311. x = chol.solve(g0);
  312. x = -x;
  313. /* and compute the current solution value */
  314. #ifdef TRACE_SOLVER
  315. f_value = 0.5 * g0.dot(x);
  316. std::cerr << "Unconstrained solution: " << f_value << std::endl;
  317. print_vector("x", x, n);
  318. #endif
  319. /* Add equality constraints to the working set A */
  320. iq = 0;
  321. for (i = 0; i < me; i++)
  322. {
  323. np = CE.col(i);
  324. compute_d(d, J, np);
  325. update_z(z, J, d, iq);
  326. update_r(R, r, d, iq);
  327. #ifdef TRACE_SOLVER
  328. print_matrix("R", R, iq);
  329. print_vector("z", z, n);
  330. print_vector("r", r, iq);
  331. print_vector("d", d, n);
  332. #endif
  333. /* compute full step length t2: i.e., the minimum step in primal space s.t. the contraint
  334. becomes feasible */
  335. t2 = 0.0;
  336. if (std::abs(z.dot(z)) > std::numeric_limits<double>::epsilon()) // i.e. z != 0
  337. t2 = (-np.dot(x) - ce0(i)) / z.dot(np);
  338. x += t2 * z;
  339. /* set u = u+ */
  340. u(iq) = t2;
  341. u.head(iq) -= t2 * r.head(iq);
  342. /* compute the new solution value */
  343. #ifdef TRACE_SOLVER
  344. f_value += 0.5 * (t2 * t2) * z.dot(np);
  345. #endif
  346. A(i) = -i - 1;
  347. if (!add_constraint(R, J, d, iq, R_norm))
  348. {
  349. // FIXME: it should raise an error
  350. // Equality constraints are linearly dependent
  351. return false;
  352. }
  353. }
  354. /* set iai = K \ A */
  355. for (i = 0; i < mi; i++)
  356. iai(i) = i;
  357. l1:
  358. #ifdef TRACE_SOLVER
  359. print_vector("x", x, n);
  360. #endif
  361. /* step 1: choose a violated constraint */
  362. for (i = me; i < iq; i++)
  363. {
  364. ip = A(i);
  365. iai(ip) = -1;
  366. }
  367. /* compute s(x) = ci^T * x + ci0 for all elements of K \ A */
  368. ss = 0.0;
  369. psi = 0.0; /* this value will contain the sum of all infeasibilities */
  370. ip = 0; /* ip will be the index of the chosen violated constraint */
  371. for (i = 0; i < mi; i++)
  372. {
  373. iaexcl[i] = true;
  374. sum = CI.col(i).dot(x) + ci0(i);
  375. s(i) = sum;
  376. psi += std::min(0.0, sum);
  377. }
  378. #ifdef TRACE_SOLVER
  379. print_vector("s", s, mi);
  380. #endif
  381. if (std::abs(psi) <= mi * std::numeric_limits<double>::epsilon() * c1 * c2* 100.0)
  382. {
  383. /* numerically there are not infeasibilities anymore */
  384. return true;
  385. }
  386. /* save old values for u, x and A */
  387. u_old.head(iq) = u.head(iq);
  388. A_old.head(iq) = A.head(iq);
  389. x_old = x;
  390. l2: /* Step 2: check for feasibility and determine a new S-pair */
  391. for (i = 0; i < mi; i++)
  392. {
  393. if (s(i) < ss && iai(i) != -1 && iaexcl[i])
  394. {
  395. ss = s(i);
  396. ip = i;
  397. }
  398. }
  399. if (ss >= 0.0)
  400. {
  401. return true;
  402. }
  403. /* set np = n(ip) */
  404. np = CI.col(ip);
  405. /* set u = (u 0)^T */
  406. u(iq) = 0.0;
  407. /* add ip to the active set A */
  408. A(iq) = ip;
  409. #ifdef TRACE_SOLVER
  410. std::cerr << "Trying with constraint " << ip << std::endl;
  411. print_vector("np", np, n);
  412. #endif
  413. l2a:/* Step 2a: determine step direction */
  414. /* compute z = H np: the step direction in the primal space (through J, see the paper) */
  415. compute_d(d, J, np);
  416. update_z(z, J, d, iq);
  417. /* compute N* np (if q > 0): the negative of the step direction in the dual space */
  418. update_r(R, r, d, iq);
  419. #ifdef TRACE_SOLVER
  420. std::cerr << "Step direction z" << std::endl;
  421. print_vector("z", z, n);
  422. print_vector("r", r, iq + 1);
  423. print_vector("u", u, iq + 1);
  424. print_vector("d", d, n);
  425. print_ivector("A", A, iq + 1);
  426. #endif
  427. /* Step 2b: compute step length */
  428. l = 0;
  429. /* Compute t1: partial step length (maximum step in dual space without violating dual feasibility */
  430. t1 = inf; /* +inf */
  431. /* find the index l s.t. it reaches the minimum of u+(x) / r */
  432. for (k = me; k < iq; k++)
  433. {
  434. double tmp;
  435. if (r(k) > 0.0 && ((tmp = u(k) / r(k)) < t1) )
  436. {
  437. t1 = tmp;
  438. l = A(k);
  439. }
  440. }
  441. /* Compute t2: full step length (minimum step in primal space such that the constraint ip becomes feasible */
  442. if (std::abs(z.dot(z)) > std::numeric_limits<double>::epsilon()) // i.e. z != 0
  443. t2 = -s(ip) / z.dot(np);
  444. else
  445. t2 = inf; /* +inf */
  446. /* the step is chosen as the minimum of t1 and t2 */
  447. t = std::min(t1, t2);
  448. #ifdef TRACE_SOLVER
  449. std::cerr << "Step sizes: " << t << " (t1 = " << t1 << ", t2 = " << t2 << ") ";
  450. #endif
  451. /* Step 2c: determine new S-pair and take step: */
  452. /* case (i): no step in primal or dual space */
  453. if (t >= inf)
  454. {
  455. /* QPP is infeasible */
  456. // FIXME: unbounded to raise
  457. return false;
  458. }
  459. /* case (ii): step in dual space */
  460. if (t2 >= inf)
  461. {
  462. /* set u = u + t * [-r 1) and drop constraint l from the active set A */
  463. u.head(iq) -= t * r.head(iq);
  464. u(iq) += t;
  465. iai(l) = l;
  466. delete_constraint(R, J, A, u, p, iq, l);
  467. #ifdef TRACE_SOLVER
  468. std::cerr << " in dual space: "
  469. << f_value << std::endl;
  470. print_vector("x", x, n);
  471. print_vector("z", z, n);
  472. print_ivector("A", A, iq + 1);
  473. #endif
  474. goto l2a;
  475. }
  476. /* case (iii): step in primal and dual space */
  477. x += t * z;
  478. /* update the solution value */
  479. #ifdef TRACE_SOLVER
  480. f_value += t * z.dot(np) * (0.5 * t + u(iq));
  481. #endif
  482. u.head(iq) -= t * r.head(iq);
  483. u(iq) += t;
  484. #ifdef TRACE_SOLVER
  485. std::cerr << " in both spaces: "
  486. << f_value << std::endl;
  487. print_vector("x", x, n);
  488. print_vector("u", u, iq + 1);
  489. print_vector("r", r, iq + 1);
  490. print_ivector("A", A, iq + 1);
  491. #endif
  492. if (t == t2)
  493. {
  494. #ifdef TRACE_SOLVER
  495. std::cerr << "Full step has taken " << t << std::endl;
  496. print_vector("x", x, n);
  497. #endif
  498. /* full step has taken */
  499. /* add constraint ip to the active set*/
  500. if (!add_constraint(R, J, d, iq, R_norm))
  501. {
  502. iaexcl[ip] = false;
  503. delete_constraint(R, J, A, u, p, iq, ip);
  504. #ifdef TRACE_SOLVER
  505. print_matrix("R", R, n);
  506. print_ivector("A", A, iq);
  507. #endif
  508. for (i = 0; i < m; i++)
  509. iai(i) = i;
  510. for (i = 0; i < iq; i++)
  511. {
  512. A(i) = A_old(i);
  513. iai(A(i)) = -1;
  514. u(i) = u_old(i);
  515. }
  516. x = x_old;
  517. goto l2; /* go to step 2 */
  518. }
  519. else
  520. iai(ip) = -1;
  521. #ifdef TRACE_SOLVER
  522. print_matrix("R", R, n);
  523. print_ivector("A", A, iq);
  524. #endif
  525. goto l1;
  526. }
  527. /* a patial step has taken */
  528. #ifdef TRACE_SOLVER
  529. std::cerr << "Partial step has taken " << t << std::endl;
  530. print_vector("x", x, n);
  531. #endif
  532. /* drop constraint l */
  533. iai(l) = l;
  534. delete_constraint(R, J, A, u, p, iq, l);
  535. #ifdef TRACE_SOLVER
  536. print_matrix("R", R, n);
  537. print_ivector("A", A, iq);
  538. #endif
  539. s(ip) = CI.col(ip).dot(x) + ci0(ip);
  540. #ifdef TRACE_SOLVER
  541. print_vector("s", s, mi);
  542. #endif
  543. goto l2a;
  544. }