00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031 #include <iostream>
00032 #include <algorithm>
00033 #include <cmath>
00034 #include <limits>
00035 #include <sstream>
00036 #include <stdexcept>
00037 #include "QuadProg.h"
00038
00039 #ifdef TRACE_SOLVER
00040 #undef TRACE_SOLVER
00041 #endif
00042
00043 namespace QuadProgPP{
00044
00045 void compute_d(Vector<double>& d, const Matrix<double>& J, const Vector<double>& np);
00046 void update_z(Vector<double>& z, const Matrix<double>& J, const Vector<double>& d, int iq);
00047 void update_r(const Matrix<double>& R, Vector<double>& r, const Vector<double>& d, int iq);
00048 bool add_constraint(Matrix<double>& R, Matrix<double>& J, Vector<double>& d, int& iq, double& rnorm);
00049 void delete_constraint(Matrix<double>& R, Matrix<double>& J, Vector<int>& A, Vector<double>& u, int n, int p, int& iq, int l);
00050
00051
00052
00053 void cholesky_decomposition(Matrix<double>& A);
00054 void cholesky_solve(const Matrix<double>& L, Vector<double>& x, const Vector<double>& b);
00055 void forward_elimination(const Matrix<double>& L, Vector<double>& y, const Vector<double>& b);
00056 void backward_elimination(const Matrix<double>& U, Vector<double>& x, const Vector<double>& y);
00057
00058
00059
00060 double scalar_product(const Vector<double>& x, const Vector<double>& y);
00061 double distance(double a, double b);
00062
00063
00064 void print_matrix(const char* name, const Matrix<double>& A, int n = -1, int m = -1);
00065
00066 template<typename T>
00067 void print_vector(const char* name, const Vector<T>& v, int n = -1);
00068
00069
00070
00071 double QuadProg::solve_quadprog(Matrix<double>& G, Vector<double>& g0,
00072 const Matrix<double>& CE, const Vector<double>& ce0,
00073 const Matrix<double>& CI, const Vector<double>& ci0,
00074 Vector<double>& x)
00075 {
00076 std::ostringstream msg;
00077 {
00078
00079
00080 unsigned mx = std::numeric_limits<int>::max();
00081 if(G.ncols() >= mx || G.nrows() >= mx ||
00082 CE.nrows() >= mx || CE.ncols() >= mx ||
00083 CI.nrows() >= mx || CI.ncols() >= mx ||
00084 ci0.size() >= mx || ce0.size() >= mx || g0.size() >= mx){
00085 msg << "The dimensions of one of the input matrices or vectors were "
00086 << "too large." << std::endl
00087 << "The maximum allowable size for inputs to solve_quadprog is:"
00088 << mx << std::endl;
00089 throw std::logic_error(msg.str());
00090 }
00091 }
00092 int n = G.ncols(), p = CE.ncols(), m = CI.ncols();
00093 if ((int)G.nrows() != n)
00094 {
00095 msg << "The matrix G is not a square matrix (" << G.nrows() << " x "
00096 << G.ncols() << ")";
00097 throw std::logic_error(msg.str());
00098 }
00099 if ((int)CE.nrows() != n)
00100 {
00101 msg << "The matrix CE is incompatible (incorrect number of rows "
00102 << CE.nrows() << " , expecting " << n << ")";
00103 throw std::logic_error(msg.str());
00104 }
00105 if ((int)ce0.size() != p)
00106 {
00107 msg << "The vector ce0 is incompatible (incorrect dimension "
00108 << ce0.size() << ", expecting " << p << ")";
00109 throw std::logic_error(msg.str());
00110 }
00111 if ((int)CI.nrows() != n)
00112 {
00113 msg << "The matrix CI is incompatible (incorrect number of rows "
00114 << CI.nrows() << " , expecting " << n << ")";
00115 throw std::logic_error(msg.str());
00116 }
00117 if ((int)ci0.size() != m)
00118 {
00119 msg << "The vector ci0 is incompatible (incorrect dimension "
00120 << ci0.size() << ", expecting " << m << ")";
00121 throw std::logic_error(msg.str());
00122 }
00123 x.resize(n);
00124 register int i, j, k, l;
00125 int ip;
00126 Matrix<double> R(n, n), J(n, n);
00127 Vector<double> s(m + p), z(n), r(m + p), d(n), np(n), u(m + p), x_old(n), u_old(m + p);
00128 double f_value, psi, c1, c2, sum, ss, R_norm;
00129 double inf;
00130 if (std::numeric_limits<double>::has_infinity)
00131 inf = std::numeric_limits<double>::infinity();
00132 else
00133 inf = 1.0E300;
00134 double t, t1, t2;
00135
00136 Vector<int> A(m + p), A_old(m + p), iai(m + p);
00137 int q, iq, iter = 0;
00138 Vector<bool> iaexcl(m + p);
00139
00140
00141
00142 q = 0;
00143 #ifdef TRACE_SOLVER
00144 std::cout << std::endl << "Starting solve_quadprog" << std::endl;
00145 print_matrix("G", G);
00146 print_vector("g0", g0);
00147 print_matrix("CE", CE);
00148 print_vector("ce0", ce0);
00149 print_matrix("CI", CI);
00150 print_vector("ci0", ci0);
00151 #endif
00152
00153
00154
00155
00156
00157
00158 c1 = 0.0;
00159 for (i = 0; i < n; i++)
00160 {
00161 c1 += G[i][i];
00162 }
00163
00164 cholesky_decomposition(G);
00165 #ifdef TRACE_SOLVER
00166 print_matrix("G", G);
00167 #endif
00168
00169 for (i = 0; i < n; i++)
00170 {
00171 d[i] = 0.0;
00172 for (j = 0; j < n; j++)
00173 R[i][j] = 0.0;
00174 }
00175 R_norm = 1.0;
00176
00177
00178 c2 = 0.0;
00179 for (i = 0; i < n; i++)
00180 {
00181 d[i] = 1.0;
00182 forward_elimination(G, z, d);
00183 for (j = 0; j < n; j++)
00184 J[i][j] = z[j];
00185 c2 += z[i];
00186 d[i] = 0.0;
00187 }
00188 #ifdef TRACE_SOLVER
00189 print_matrix("J", J);
00190 #endif
00191
00192
00193
00194
00195
00196
00197
00198
00199 cholesky_solve(G, x, g0);
00200 for (i = 0; i < n; i++)
00201 x[i] = -x[i];
00202
00203 f_value = 0.5 * scalar_product(g0, x);
00204 #ifdef TRACE_SOLVER
00205 std::cout << "Unconstrained solution: " << f_value << std::endl;
00206 print_vector("x", x);
00207 #endif
00208
00209
00210 iq = 0;
00211 for (i = 0; i < p; i++)
00212 {
00213 for (j = 0; j < n; j++)
00214 np[j] = CE[j][i];
00215 compute_d(d, J, np);
00216 update_z(z, J, d, iq);
00217 update_r(R, r, d, iq);
00218 #ifdef TRACE_SOLVER
00219 print_matrix("R", R, n, iq);
00220 print_vector("z", z);
00221 print_vector("r", r, iq);
00222 print_vector("d", d);
00223 #endif
00224
00225
00226
00227 t2 = 0.0;
00228 if (fabs(scalar_product(z, z)) > std::numeric_limits<double>::epsilon())
00229 t2 = (-scalar_product(np, x) - ce0[i]) / scalar_product(z, np);
00230
00231
00232 for (k = 0; k < n; k++)
00233 x[k] += t2 * z[k];
00234
00235
00236 u[iq] = t2;
00237 for (k = 0; k < iq; k++)
00238 u[k] -= t2 * r[k];
00239
00240
00241 f_value += 0.5 * (t2 * t2) * scalar_product(z, np);
00242 A[i] = -i - 1;
00243
00244 if (!add_constraint(R, J, d, iq, R_norm))
00245 {
00246
00247 throw std::runtime_error("Constraints are linearly dependent");
00248 return f_value;
00249 }
00250 }
00251
00252
00253 for (i = 0; i < m; i++)
00254 iai[i] = i;
00255
00256 l1: iter++;
00257 #ifdef TRACE_SOLVER
00258 print_vector("x", x);
00259 #endif
00260
00261 for (i = p; i < iq; i++)
00262 {
00263 ip = A[i];
00264 iai[ip] = -1;
00265 }
00266
00267
00268 ss = 0.0;
00269 psi = 0.0;
00270 ip = 0;
00271 for (i = 0; i < m; i++)
00272 {
00273 iaexcl[i] = true;
00274 sum = 0.0;
00275 for (j = 0; j < n; j++)
00276 sum += CI[j][i] * x[j];
00277 sum += ci0[i];
00278 s[i] = sum;
00279 psi += std::min(0.0, sum);
00280 }
00281 #ifdef TRACE_SOLVER
00282 print_vector("s", s, m);
00283 #endif
00284
00285
00286 if (fabs(psi) <= m * std::numeric_limits<double>::epsilon() * c1 * c2* 100.0)
00287 {
00288
00289 q = iq;
00290
00291 return f_value;
00292 }
00293
00294
00295 for (i = 0; i < iq; i++)
00296 {
00297 u_old[i] = u[i];
00298 A_old[i] = A[i];
00299 }
00300
00301 for (i = 0; i < n; i++)
00302 x_old[i] = x[i];
00303
00304 l2:
00305 for (i = 0; i < m; i++)
00306 {
00307 if (s[i] < ss && iai[i] != -1 && iaexcl[i])
00308 {
00309 ss = s[i];
00310 ip = i;
00311 }
00312 }
00313 if (ss >= 0.0)
00314 {
00315 q = iq;
00316
00317 return f_value;
00318 }
00319
00320
00321 for (i = 0; i < n; i++)
00322 np[i] = CI[i][ip];
00323
00324 u[iq] = 0.0;
00325
00326 A[iq] = ip;
00327
00328 #ifdef TRACE_SOLVER
00329 std::cout << "Trying with constraint " << ip << std::endl;
00330 print_vector("np", np);
00331 #endif
00332
00333 l2a:
00334
00335 compute_d(d, J, np);
00336 update_z(z, J, d, iq);
00337
00338 update_r(R, r, d, iq);
00339 #ifdef TRACE_SOLVER
00340 std::cout << "Step direction z" << std::endl;
00341 print_vector("z", z);
00342 print_vector("r", r, iq + 1);
00343 print_vector("u", u, iq + 1);
00344 print_vector("d", d);
00345 print_vector("A", A, iq + 1);
00346 #endif
00347
00348
00349 l = 0;
00350
00351 t1 = inf;
00352
00353 for (k = p; k < iq; k++)
00354 {
00355 if (r[k] > 0.0)
00356 {
00357 if (u[k] / r[k] < t1)
00358 {
00359 t1 = u[k] / r[k];
00360 l = A[k];
00361 }
00362 }
00363 }
00364
00365 if (fabs(scalar_product(z, z)) > std::numeric_limits<double>::epsilon())
00366 t2 = -s[ip] / scalar_product(z, np);
00367 else
00368 t2 = inf;
00369
00370
00371 t = std::min(t1, t2);
00372 #ifdef TRACE_SOLVER
00373 std::cout << "Step sizes: " << t << " (t1 = " << t1 << ", t2 = " << t2 << ") ";
00374 #endif
00375
00376
00377
00378
00379 if (t >= inf)
00380 {
00381
00382
00383 q = iq;
00384 return inf;
00385 }
00386
00387 if (t2 >= inf)
00388 {
00389
00390 for (k = 0; k < iq; k++)
00391 u[k] -= t * r[k];
00392 u[iq] += t;
00393 iai[l] = l;
00394 delete_constraint(R, J, A, u, n, p, iq, l);
00395 #ifdef TRACE_SOLVER
00396 std::cout << " in dual space: "
00397 << f_value << std::endl;
00398 print_vector("x", x);
00399 print_vector("z", z);
00400 print_vector("A", A, iq + 1);
00401 #endif
00402 goto l2a;
00403 }
00404
00405
00406
00407
00408 for (k = 0; k < n; k++)
00409 x[k] += t * z[k];
00410
00411 f_value += t * scalar_product(z, np) * (0.5 * t + u[iq]);
00412
00413 for (k = 0; k < iq; k++)
00414 u[k] -= t * r[k];
00415 u[iq] += t;
00416 #ifdef TRACE_SOLVER
00417 std::cout << " in both spaces: "
00418 << f_value << std::endl;
00419 print_vector("x", x);
00420 print_vector("u", u, iq + 1);
00421 print_vector("r", r, iq + 1);
00422 print_vector("A", A, iq + 1);
00423 #endif
00424
00425 if (fabs(t - t2) < std::numeric_limits<double>::epsilon())
00426 {
00427 #ifdef TRACE_SOLVER
00428 std::cout << "Full step has taken " << t << std::endl;
00429 print_vector("x", x);
00430 #endif
00431
00432
00433 if (!add_constraint(R, J, d, iq, R_norm))
00434 {
00435 iaexcl[ip] = false;
00436 delete_constraint(R, J, A, u, n, p, iq, ip);
00437 #ifdef TRACE_SOLVER
00438 print_matrix("R", R);
00439 print_vector("A", A, iq);
00440 print_vector("iai", iai);
00441 #endif
00442 for (i = 0; i < m; i++)
00443 iai[i] = i;
00444 for (i = p; i < iq; i++)
00445 {
00446 A[i] = A_old[i];
00447 u[i] = u_old[i];
00448 iai[A[i]] = -1;
00449 }
00450 for (i = 0; i < n; i++)
00451 x[i] = x_old[i];
00452 goto l2;
00453 }
00454 else
00455 iai[ip] = -1;
00456 #ifdef TRACE_SOLVER
00457 print_matrix("R", R);
00458 print_vector("A", A, iq);
00459 print_vector("iai", iai);
00460 #endif
00461 goto l1;
00462 }
00463
00464
00465 #ifdef TRACE_SOLVER
00466 std::cout << "Partial step has taken " << t << std::endl;
00467 print_vector("x", x);
00468 #endif
00469
00470 iai[l] = l;
00471 delete_constraint(R, J, A, u, n, p, iq, l);
00472 #ifdef TRACE_SOLVER
00473 print_matrix("R", R);
00474 print_vector("A", A, iq);
00475 #endif
00476
00477
00478 sum = 0.0;
00479 for (k = 0; k < n; k++)
00480 sum += CI[k][ip] * x[k];
00481 s[ip] = sum + ci0[ip];
00482
00483 #ifdef TRACE_SOLVER
00484 print_vector("s", s, m);
00485 #endif
00486 goto l2a;
00487 }
00488
00489 inline void compute_d(Vector<double>& d, const Matrix<double>& J, const Vector<double>& np)
00490 {
00491 register int i, j, n = d.size();
00492 register double sum;
00493
00494
00495 for (i = 0; i < n; i++)
00496 {
00497 sum = 0.0;
00498 for (j = 0; j < n; j++)
00499 sum += J[j][i] * np[j];
00500 d[i] = sum;
00501 }
00502 }
00503
00504 inline void update_z(Vector<double>& z, const Matrix<double>& J, const Vector<double>& d, int iq)
00505 {
00506 register int i, j, n = z.size();
00507
00508
00509 for (i = 0; i < n; i++)
00510 {
00511 z[i] = 0.0;
00512 for (j = iq; j < n; j++)
00513 z[i] += J[i][j] * d[j];
00514 }
00515 }
00516
00517 inline void update_r(const Matrix<double>& R, Vector<double>& r, const Vector<double>& d, int iq)
00518 {
00519 register int i, j;
00520 register double sum;
00521
00522
00523 for (i = iq - 1; i >= 0; i--)
00524 {
00525 sum = 0.0;
00526 for (j = i + 1; j < iq; j++)
00527 sum += R[i][j] * r[j];
00528 r[i] = (d[i] - sum) / R[i][i];
00529 }
00530 }
00531
00532 bool add_constraint(Matrix<double>& R, Matrix<double>& J, Vector<double>& d, int& iq, double& R_norm)
00533 {
00534 int n = d.size();
00535 #ifdef TRACE_SOLVER
00536 std::cout << "Add constraint " << iq << '/';
00537 #endif
00538 register int i, j, k;
00539 double cc, ss, h, t1, t2, xny;
00540
00541
00542
00543
00544
00545 for (j = n - 1; j >= iq + 1; j--)
00546 {
00547
00548
00549
00550
00551
00552
00553
00554
00555 cc = d[j - 1];
00556 ss = d[j];
00557 h = distance(cc, ss);
00558 if (fabs(h) < std::numeric_limits<double>::epsilon())
00559 continue;
00560 d[j] = 0.0;
00561 ss = ss / h;
00562 cc = cc / h;
00563 if (cc < 0.0)
00564 {
00565 cc = -cc;
00566 ss = -ss;
00567 d[j - 1] = -h;
00568 }
00569 else
00570 d[j - 1] = h;
00571 xny = ss / (1.0 + cc);
00572 for (k = 0; k < n; k++)
00573 {
00574 t1 = J[k][j - 1];
00575 t2 = J[k][j];
00576 J[k][j - 1] = t1 * cc + t2 * ss;
00577 J[k][j] = xny * (t1 + J[k][j - 1]) - t2;
00578 }
00579 }
00580
00581 iq++;
00582
00583
00584
00585 for (i = 0; i < iq; i++)
00586 R[i][iq - 1] = d[i];
00587 #ifdef TRACE_SOLVER
00588 std::cout << iq << std::endl;
00589 print_matrix("R", R, iq, iq);
00590 print_matrix("J", J);
00591 print_vector("d", d, iq);
00592 #endif
00593
00594 if (fabs(d[iq - 1]) <= std::numeric_limits<double>::epsilon() * R_norm)
00595 {
00596
00597 return false;
00598 }
00599 R_norm = std::max<double>(R_norm, fabs(d[iq - 1]));
00600 return true;
00601 }
00602
00603 void delete_constraint(Matrix<double>& R, Matrix<double>& J, Vector<int>& A, Vector<double>& u, int n, int p, int& iq, int l)
00604 {
00605 #ifdef TRACE_SOLVER
00606 std::cout << "Delete constraint " << l << ' ' << iq;
00607 #endif
00608 register int i, j, k, qq = -1;
00609 double cc, ss, h, xny, t1, t2;
00610
00611
00612 for (i = p; i < iq; i++)
00613 if (A[i] == l)
00614 {
00615 qq = i;
00616 break;
00617 }
00618
00619
00620 for (i = qq; i < iq - 1; i++)
00621 {
00622 A[i] = A[i + 1];
00623 u[i] = u[i + 1];
00624 for (j = 0; j < n; j++)
00625 R[j][i] = R[j][i + 1];
00626 }
00627
00628 A[iq - 1] = A[iq];
00629 u[iq - 1] = u[iq];
00630 A[iq] = 0;
00631 u[iq] = 0.0;
00632 for (j = 0; j < iq; j++)
00633 R[j][iq - 1] = 0.0;
00634
00635 iq--;
00636 #ifdef TRACE_SOLVER
00637 std::cout << '/' << iq << std::endl;
00638 #endif
00639
00640 if (iq == 0)
00641 return;
00642
00643 for (j = qq; j < iq; j++)
00644 {
00645 cc = R[j][j];
00646 ss = R[j + 1][j];
00647 h = distance(cc, ss);
00648 if (fabs(h) < std::numeric_limits<double>::epsilon())
00649 continue;
00650 cc = cc / h;
00651 ss = ss / h;
00652 R[j + 1][j] = 0.0;
00653 if (cc < 0.0)
00654 {
00655 R[j][j] = -h;
00656 cc = -cc;
00657 ss = -ss;
00658 }
00659 else
00660 R[j][j] = h;
00661
00662 xny = ss / (1.0 + cc);
00663 for (k = j + 1; k < iq; k++)
00664 {
00665 t1 = R[j][k];
00666 t2 = R[j + 1][k];
00667 R[j][k] = t1 * cc + t2 * ss;
00668 R[j + 1][k] = xny * (t1 + R[j][k]) - t2;
00669 }
00670 for (k = 0; k < n; k++)
00671 {
00672 t1 = J[k][j];
00673 t2 = J[k][j + 1];
00674 J[k][j] = t1 * cc + t2 * ss;
00675 J[k][j + 1] = xny * (J[k][j] + t1) - t2;
00676 }
00677 }
00678 }
00679
00680 inline double distance(double a, double b)
00681 {
00682 register double a1, b1, t;
00683 a1 = fabs(a);
00684 b1 = fabs(b);
00685 if (a1 > b1)
00686 {
00687 t = (b1 / a1);
00688 return a1 * ::std::sqrt(1.0 + t * t);
00689 }
00690 else
00691 if (b1 > a1)
00692 {
00693 t = (a1 / b1);
00694 return b1 * ::std::sqrt(1.0 + t * t);
00695 }
00696 return a1 * ::std::sqrt(2.0);
00697 }
00698
00699
00700 inline double scalar_product(const Vector<double>& x, const Vector<double>& y)
00701 {
00702 register int i, n = x.size();
00703 register double sum;
00704
00705 sum = 0.0;
00706 for (i = 0; i < n; i++)
00707 sum += x[i] * y[i];
00708 return sum;
00709 }
00710
00711 void cholesky_decomposition(Matrix<double>& A)
00712 {
00713 register int i, j, k, n = A.nrows();
00714 register double sum;
00715
00716 for (i = 0; i < n; i++)
00717 {
00718 for (j = i; j < n; j++)
00719 {
00720 sum = A[i][j];
00721 for (k = i - 1; k >= 0; k--)
00722 sum -= A[i][k]*A[j][k];
00723 if (i == j)
00724 {
00725 if (sum <= 0.0)
00726 {
00727 std::ostringstream os;
00728
00729 print_matrix("A", A);
00730 os << "Error in cholesky decomposition, sum: " << sum;
00731 throw std::logic_error(os.str());
00732 exit(-1);
00733 }
00734 A[i][i] = ::std::sqrt(sum);
00735 }
00736 else
00737 A[j][i] = sum / A[i][i];
00738 }
00739 for (k = i + 1; k < n; k++)
00740 A[i][k] = A[k][i];
00741 }
00742 }
00743
00744 void cholesky_solve(const Matrix<double>& L, Vector<double>& x, const Vector<double>& b)
00745 {
00746 int n = L.nrows();
00747 Vector<double> y(n);
00748
00749
00750 forward_elimination(L, y, b);
00751
00752 backward_elimination(L, x, y);
00753 }
00754
00755 inline void forward_elimination(const Matrix<double>& L, Vector<double>& y, const Vector<double>& b)
00756 {
00757 register int i, j, n = L.nrows();
00758
00759 y[0] = b[0] / L[0][0];
00760 for (i = 1; i < n; i++)
00761 {
00762 y[i] = b[i];
00763 for (j = 0; j < i; j++)
00764 y[i] -= L[i][j] * y[j];
00765 y[i] = y[i] / L[i][i];
00766 }
00767 }
00768
00769 inline void backward_elimination(const Matrix<double>& U, Vector<double>& x, const Vector<double>& y)
00770 {
00771 register int i, j, n = U.nrows();
00772
00773 x[n - 1] = y[n - 1] / U[n - 1][n - 1];
00774 for (i = n - 2; i >= 0; i--)
00775 {
00776 x[i] = y[i];
00777 for (j = i + 1; j < n; j++)
00778 x[i] -= U[i][j] * x[j];
00779 x[i] = x[i] / U[i][i];
00780 }
00781 }
00782
00783 void print_matrix(const char* name, const Matrix<double>& A, int n, int m)
00784 {
00785 std::ostringstream s;
00786 std::string t;
00787 if (n == -1)
00788 n = A.nrows();
00789 if (m == -1)
00790 m = A.ncols();
00791
00792 s << name << ": " << std::endl;
00793 for (int i = 0; i < n; i++)
00794 {
00795 s << " ";
00796 for (int j = 0; j < m; j++)
00797 s << A[i][j] << ", ";
00798 s << std::endl;
00799 }
00800 t = s.str();
00801 t = t.substr(0, t.size() - 3);
00802
00803 std::cout << t << std::endl;
00804 }
00805
00806 template<typename T>
00807 void print_vector(const char* name, const Vector<T>& v, int n)
00808 {
00809 std::ostringstream s;
00810 std::string t;
00811 if (n == -1)
00812 n = v.size();
00813
00814 s << name << ": " << std::endl << " ";
00815 for (int i = 0; i < n; i++)
00816 {
00817 s << v[i] << ", ";
00818 }
00819 t = s.str();
00820 t = t.substr(0, t.size() - 2);
00821
00822 std::cout << t << std::endl;
00823 }
00824 }