quadprog.cpp 16 KB

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