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


eus_qp
Author(s):
autogenerated on Fri May 14 2021 02:51:43