quadprog.cpp 15 KB

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