quadprog.cpp 15 KB

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