1 #ifndef _EIGEN_QUADSOLVE_HPP_ 2 #define _EIGEN_QUADSOLVE_HPP_ 84 #include <Eigen/Dense> 90 template<
typename Scalar>
99 return a1 * std::sqrt(1.0 + t * t);
105 return b1 * std::sqrt(1.0 + t * t);
107 return a1 * std::sqrt(2.0);
112 inline void compute_d(VectorXd &d,
const MatrixXd& J,
const VectorXd& np)
114 d = J.adjoint() * np;
117 inline void update_z(VectorXd& z,
const MatrixXd& J,
const VectorXd& d,
int iq)
119 z = J.rightCols(z.size()-iq) * d.tail(d.size()-iq);
122 inline void update_r(
const MatrixXd& R, VectorXd& r,
const VectorXd& d,
int iq)
124 r.head(iq)= R.topLeftCorner(iq,iq).triangularView<Upper>().solve(d.head(iq));
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);
132 const MatrixXd & CE,
const VectorXd & ce0,
133 const MatrixXd & CI,
const VectorXd & ci0,
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());
141 LLT<MatrixXd,Lower> chol(G.cols());
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();
149 VectorXi A(m + p), A_old(m + p), iai(m + p);
176 J = chol.matrixU().solve(J);
179 print_matrix(
"J", J, n);
192 f_value = 0.5 * g0.dot(x);
194 std::cerr <<
"Unconstrained solution: " << f_value << std::endl;
195 print_vector(
"x", x, n);
200 for (i = 0; i < me; i++)
207 print_matrix(
"R", R, iq);
208 print_vector(
"z", z, n);
209 print_vector(
"r", r, iq);
210 print_vector(
"d", d, n);
216 if (std::abs(z.dot(z)) > std::numeric_limits<double>::epsilon())
217 t2 = (-np.dot(x) - ce0(i)) / z.dot(np);
223 u.head(iq) -= t2 * r.head(iq);
226 f_value += 0.5 * (t2 * t2) * z.dot(np);
238 for (i = 0; i < mi; i++)
243 print_vector(
"x", x, n);
246 for (i = me; i < iq; i++)
256 for (i = 0; i < mi; i++)
259 sum = CI.col(i).dot(x) + ci0(i);
261 psi += std::min(0.0, sum);
264 print_vector(
"s",
s, mi);
268 if (std::abs(psi) <= mi * std::numeric_limits<double>::epsilon() * c1 * c2* 100.0)
276 u_old.head(iq) = u.head(iq);
277 A_old.head(iq) = A.head(iq);
281 for (i = 0; i < mi; i++)
283 if (
s(i) < ss && iai(i) != -1 && iaexcl[i])
303 std::cerr <<
"Trying with constraint " << ip << std::endl;
304 print_vector(
"np", np, n);
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);
327 for (k = me; k < iq; k++)
330 if (r(k) > 0.0 && ((tmp = u(k) / r(k)) < t1) )
337 if (std::abs(z.dot(z)) > std::numeric_limits<double>::epsilon())
338 t2 = -
s(ip) / z.dot(np);
343 t = std::min(t1, t2);
345 std::cerr <<
"Step sizes: " << t <<
" (t1 = " << t1 <<
", t2 = " << t2 <<
") ";
362 u.head(iq) -= t * r.head(iq);
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);
380 f_value += t * z.dot(np) * (0.5 * t + u(iq));
382 u.head(iq) -= t * r.head(iq);
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);
396 std::cerr <<
"Full step has taken " << t << std::endl;
397 print_vector(
"x", x, n);
406 print_matrix(
"R", R, n);
407 print_ivector(
"A", A, iq);
409 for (i = 0; i < m; i++)
411 for (i = 0; i < iq; i++)
423 print_matrix(
"R", R, n);
424 print_ivector(
"A", A, iq);
431 std::cerr <<
"Partial step has taken " << t << std::endl;
432 print_vector(
"x", x, n);
438 print_matrix(
"R", R, n);
439 print_ivector(
"A", A, iq);
442 s(ip) = CI.col(ip).dot(x) + ci0(ip);
445 print_vector(
"s",
s, mi);
451 inline bool add_constraint(MatrixXd& R, MatrixXd& J, VectorXd& d,
int& iq,
double& R_norm)
455 std::cerr <<
"Add constraint " << iq <<
'/';
458 double cc, ss, h, t1, t2, xny;
464 for (j = n - 1; j >= iq + 1; j--)
490 xny = ss / (1.0 + cc);
491 for (k = 0; k < n; k++)
495 J(k,j - 1) = t1 * cc + t2 * ss;
496 J(k,j) = xny * (t1 + J(k,j - 1)) - t2;
504 R.col(iq-1).head(iq) = d.head(iq);
506 std::cerr << iq << std::endl;
509 if (std::abs(d(iq - 1)) <= std::numeric_limits<double>::epsilon() * R_norm)
512 R_norm = std::max<double>(R_norm, std::abs(d(iq - 1)));
517 inline void delete_constraint(MatrixXd& R, MatrixXd& J, VectorXi& A, VectorXd& u,
int p,
int& iq,
int l)
522 std::cerr <<
"Delete constraint " << l <<
' ' << iq;
525 double cc, ss, h, xny, t1, t2;
528 for (i = p; i < iq; i++)
536 for (i = qq; i < iq - 1; i++)
540 R.col(i) = R.col(i+1);
547 for (j = 0; j < iq; j++)
552 std::cerr <<
'/' << iq << std::endl;
558 for (j = qq; j < iq; j++)
577 xny = ss / (1.0 + cc);
578 for (k = j + 1; k < iq; k++)
582 R(j,k) = t1 * cc + t2 * ss;
583 R(j + 1,k) = xny * (t1 + R(j,k)) - t2;
585 for (k = 0; k < n; k++)
589 J(k,j) = t1 * cc + t2 * ss;
590 J(k,j + 1) = xny * (J(k,j) + t1) - t2;
void delete_constraint(MatrixXd &R, MatrixXd &J, VectorXi &A, VectorXd &u, int p, int &iq, int l)
double solve_quadprog(MatrixXd &G, VectorXd &g0, const MatrixXd &CE, const VectorXd &ce0, const MatrixXd &CI, const VectorXd &ci0, VectorXd &x)
void compute_d(VectorXd &d, const MatrixXd &J, const VectorXd &np)
void update_r(const MatrixXd &R, VectorXd &r, const VectorXd &d, int iq)
void update_z(VectorXd &z, const MatrixXd &J, const VectorXd &d, int iq)
bool add_constraint(MatrixXd &R, MatrixXd &J, VectorXd &d, int &iq, double &R_norm)
Scalar distance(Scalar a, Scalar b)