Crocoddyl
box-qp.cpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2019-2022, University of Edinburgh, Heriot-Watt University
5 // Copyright note valid unless otherwise stated in individual files.
6 // All rights reserved.
8 
9 #include "crocoddyl/core/solvers/box-qp.hpp"
10 
11 #include <iostream>
12 
13 #include "crocoddyl/core/utils/exception.hpp"
14 
15 namespace crocoddyl {
16 
17 BoxQP::BoxQP(const std::size_t nx, const std::size_t maxiter,
18  const double th_acceptstep, const double th_grad, const double reg)
19  : nx_(nx),
20  maxiter_(maxiter),
21  th_acceptstep_(th_acceptstep),
22  th_grad_(th_grad),
23  reg_(reg),
24  fold_(0.),
25  fnew_(0.),
26  x_(nx),
27  xnew_(nx),
28  g_(nx),
29  dx_(nx),
30  xo_(nx),
31  dxo_(nx),
32  qo_(nx),
33  Ho_(nx, nx) {
34  // Check if values have a proper range
35  if (0. >= th_acceptstep && th_acceptstep >= 0.5) {
36  std::cerr << "Warning: th_acceptstep value should between 0 and 0.5"
37  << std::endl;
38  }
39  if (0. > th_grad) {
40  std::cerr << "Warning: th_grad value has to be positive." << std::endl;
41  }
42  if (0. > reg) {
43  std::cerr << "Warning: reg value has to be positive." << std::endl;
44  }
45 
46  // Initialized the values of vectors
47  x_.setZero();
48  xnew_.setZero();
49  g_.setZero();
50  dx_.setZero();
51  xo_.setZero();
52  dxo_.setZero();
53  qo_.setZero();
54  Ho_.setZero();
55 
56  // Reserve the space and compute alphas
57  solution_.x = Eigen::VectorXd::Zero(nx);
58  solution_.clamped_idx.reserve(nx_);
59  solution_.free_idx.reserve(nx_);
60  const std::size_t n_alphas_ = 10;
61  alphas_.resize(n_alphas_);
62  for (std::size_t n = 0; n < n_alphas_; ++n) {
63  alphas_[n] = 1. / pow(2., static_cast<double>(n));
64  }
65 }
66 
68 
69 const BoxQPSolution& BoxQP::solve(const Eigen::MatrixXd& H,
70  const Eigen::VectorXd& q,
71  const Eigen::VectorXd& lb,
72  const Eigen::VectorXd& ub,
73  const Eigen::VectorXd& xinit) {
74  if (static_cast<std::size_t>(H.rows()) != nx_ ||
75  static_cast<std::size_t>(H.cols()) != nx_) {
76  throw_pretty("Invalid argument: "
77  << "H has wrong dimension (it should be " +
78  std::to_string(nx_) + "," + std::to_string(nx_) + ")");
79  }
80  if (static_cast<std::size_t>(q.size()) != nx_) {
81  throw_pretty("Invalid argument: "
82  << "q has wrong dimension (it should be " +
83  std::to_string(nx_) + ")");
84  }
85  if (static_cast<std::size_t>(lb.size()) != nx_) {
86  throw_pretty("Invalid argument: "
87  << "lb has wrong dimension (it should be " +
88  std::to_string(nx_) + ")");
89  }
90  if (static_cast<std::size_t>(ub.size()) != nx_) {
91  throw_pretty("Invalid argument: "
92  << "ub has wrong dimension (it should be " +
93  std::to_string(nx_) + ")");
94  }
95  if (static_cast<std::size_t>(xinit.size()) != nx_) {
96  throw_pretty("Invalid argument: "
97  << "xinit has wrong dimension (it should be " +
98  std::to_string(nx_) + ")");
99  }
100 
101  // We need to enforce feasible warm-starting of the algorithm
102  for (std::size_t i = 0; i < nx_; ++i) {
103  x_(i) = std::max(std::min(xinit(i), ub(i)), lb(i));
104  }
105 
106  // Start the numerical iterations
107  for (std::size_t k = 0; k < maxiter_; ++k) {
108  solution_.clamped_idx.clear();
109  solution_.free_idx.clear();
110  // Compute the Cauchy point and active set
111  g_ = q;
112  g_.noalias() += H * x_;
113  for (std::size_t j = 0; j < nx_; ++j) {
114  const double gj = g_(j);
115  const double xj = x_(j);
116  const double lbj = lb(j);
117  const double ubj = ub(j);
118  if ((xj == lbj && gj > 0.) || (xj == ubj && gj < 0.)) {
119  solution_.clamped_idx.push_back(j);
120  } else {
121  solution_.free_idx.push_back(j);
122  }
123  }
124 
125  // Compute the search direction as Newton step along the free space
126  nf_ = solution_.free_idx.size();
127  nc_ = solution_.clamped_idx.size();
128  Eigen::VectorBlock<Eigen::VectorXd> xf = xo_.head(nf_);
129  Eigen::VectorBlock<Eigen::VectorXd> xc = xo_.tail(nc_);
130  Eigen::VectorBlock<Eigen::VectorXd> dxf = dxo_.head(nf_);
131  Eigen::VectorBlock<Eigen::VectorXd> qf = qo_.head(nf_);
132  Eigen::Block<Eigen::MatrixXd> Hff = Ho_.topLeftCorner(nf_, nf_);
133  Eigen::Block<Eigen::MatrixXd> Hfc = Ho_.topRightCorner(nf_, nc_);
134  for (std::size_t i = 0; i < nf_; ++i) {
135  const std::size_t fi = solution_.free_idx[i];
136  qf(i) = q(fi);
137  xf(i) = x_(fi);
138  for (std::size_t j = 0; j < nf_; ++j) {
139  Hff(i, j) = H(fi, solution_.free_idx[j]);
140  }
141  for (std::size_t j = 0; j < nc_; ++j) {
142  const std::size_t cj = solution_.clamped_idx[j];
143  xc(j) = x_(cj);
144  Hfc(i, j) = H(fi, cj);
145  }
146  }
147  if (reg_ != 0.) {
148  Hff.diagonal().array() += reg_;
149  }
150  Hff_inv_llt_.compute(Hff);
151  const Eigen::ComputationInfo& info = Hff_inv_llt_.info();
152  if (info != Eigen::Success) {
153  throw_pretty("backward_error");
154  }
155  solution_.Hff_inv.setIdentity(nf_, nf_);
156  Hff_inv_llt_.solveInPlace(solution_.Hff_inv);
157  qf.noalias() += Hff * xf;
158  if (nc_ != 0) {
159  qf.noalias() += Hfc * xc;
160  }
161  dxf = -qf;
162  Hff_inv_llt_.solveInPlace(dxf);
163  dx_.setZero();
164  for (std::size_t i = 0; i < nf_; ++i) {
165  dx_(solution_.free_idx[i]) = dxf(i);
166  g_(solution_.free_idx[i]) = -qf(i);
167  }
168 
169  // Try different step lengths
170  fold_ = 0.5 * x_.dot(H * x_) + q.dot(x_);
171  for (std::vector<double>::const_iterator it = alphas_.begin();
172  it != alphas_.end(); ++it) {
173  double steplength = *it;
174  for (std::size_t i = 0; i < nx_; ++i) {
175  xnew_(i) =
176  std::max(std::min(x_(i) + steplength * dx_(i), ub(i)), lb(i));
177  }
178  fnew_ = 0.5 * xnew_.dot(H * xnew_) + q.dot(xnew_);
179  if (fold_ - fnew_ > th_acceptstep_ * g_.dot(x_ - xnew_)) {
180  x_ = xnew_;
181  break;
182  }
183  }
184 
185  // Check convergence
186  if (qf.lpNorm<Eigen::Infinity>() <= th_grad_) {
187  solution_.x = x_;
188  return solution_;
189  }
190  }
191  solution_.x = x_;
192  return solution_;
193 }
194 
195 const BoxQPSolution& BoxQP::get_solution() const { return solution_; }
196 
197 std::size_t BoxQP::get_nx() const { return nx_; }
198 
199 std::size_t BoxQP::get_maxiter() const { return maxiter_; }
200 
201 double BoxQP::get_th_acceptstep() const { return th_acceptstep_; }
202 
203 double BoxQP::get_th_grad() const { return th_grad_; }
204 
205 double BoxQP::get_reg() const { return reg_; }
206 
207 const std::vector<double>& BoxQP::get_alphas() const { return alphas_; }
208 
209 void BoxQP::set_nx(const std::size_t nx) {
210  nx_ = nx;
211  x_.conservativeResize(nx);
212  xnew_.conservativeResize(nx);
213  g_.conservativeResize(nx);
214  dx_.conservativeResize(nx);
215  xo_.conservativeResize(nx);
216  dxo_.conservativeResize(nx);
217  qo_.conservativeResize(nx);
218  Ho_.conservativeResize(nx, nx);
219 }
220 
221 void BoxQP::set_maxiter(const std::size_t maxiter) { maxiter_ = maxiter; }
222 
223 void BoxQP::set_th_acceptstep(const double th_acceptstep) {
224  if (0. >= th_acceptstep && th_acceptstep >= 0.5) {
225  throw_pretty("Invalid argument: "
226  << "th_acceptstep value should between 0 and 0.5");
227  }
228  th_acceptstep_ = th_acceptstep;
229 }
230 
231 void BoxQP::set_th_grad(const double th_grad) {
232  if (0. > th_grad) {
233  throw_pretty("Invalid argument: "
234  << "th_grad value has to be positive.");
235  }
236  th_grad_ = th_grad;
237 }
238 
239 void BoxQP::set_reg(const double reg) {
240  if (0. > reg) {
241  throw_pretty("Invalid argument: "
242  << "reg value has to be positive.");
243  }
244  reg_ = reg;
245 }
246 
247 void BoxQP::set_alphas(const std::vector<double>& alphas) {
248  double prev_alpha = alphas[0];
249  if (prev_alpha != 1.) {
250  std::cerr << "Warning: alpha[0] should be 1" << std::endl;
251  }
252  for (std::size_t i = 1; i < alphas.size(); ++i) {
253  double alpha = alphas[i];
254  if (0. >= alpha) {
255  throw_pretty("Invalid argument: "
256  << "alpha values has to be positive.");
257  }
258  if (alpha >= prev_alpha) {
259  throw_pretty("Invalid argument: "
260  << "alpha values are monotonously decreasing.");
261  }
262  prev_alpha = alpha;
263  }
264  alphas_ = alphas;
265 }
266 
267 } // namespace crocoddyl
void set_reg(const double reg)
Modify the regularization value.
Definition: box-qp.cpp:239
void set_th_grad(const double th_grad)
Modify the gradient tolerance threshold.
Definition: box-qp.cpp:231
const BoxQPSolution & get_solution() const
Return the stored solution.
Definition: box-qp.cpp:195
double get_th_grad() const
Return the gradient tolerance threshold.
Definition: box-qp.cpp:203
const BoxQPSolution & solve(const Eigen::MatrixXd &H, const Eigen::VectorXd &q, const Eigen::VectorXd &lb, const Eigen::VectorXd &ub, const Eigen::VectorXd &xinit)
Compute the solution of bound-constrained QP based on Newton projection.
Definition: box-qp.cpp:69
std::size_t get_maxiter() const
Return the maximum allowed number of iterations.
Definition: box-qp.cpp:199
void set_maxiter(const std::size_t maxiter)
Modify the maximum allowed number of iterations.
Definition: box-qp.cpp:221
~BoxQP()
Destroy the Projected-Newton QP solver.
Definition: box-qp.cpp:67
void set_nx(const std::size_t nx)
Modify the decision vector dimension.
Definition: box-qp.cpp:209
EIGEN_MAKE_ALIGNED_OPERATOR_NEW BoxQP(const std::size_t nx, const std::size_t maxiter=100, const double th_acceptstep=0.1, const double th_grad=1e-9, const double reg=1e-9)
Initialize the Projected-Newton QP for bound constraints.
Definition: box-qp.cpp:17
const std::vector< double > & get_alphas() const
Return the stack of step lengths using by the line-search procedure.
Definition: box-qp.cpp:207
void set_alphas(const std::vector< double > &alphas)
Modify the stack of step lengths using by the line-search procedure.
Definition: box-qp.cpp:247
void set_th_acceptstep(const double th_acceptstep)
Modify the acceptance step threshold.
Definition: box-qp.cpp:223
double get_th_acceptstep() const
Return the acceptance step threshold.
Definition: box-qp.cpp:201
double get_reg() const
Return the regularization value.
Definition: box-qp.cpp:205
std::size_t get_nx() const
Return the decision vector dimension.
Definition: box-qp.cpp:197
Box QP solution.
Definition: box-qp.hpp:29
std::vector< size_t > free_idx
Free space indexes.
Definition: box-qp.hpp:52
Eigen::MatrixXd Hff_inv
Inverse of the free space Hessian.
Definition: box-qp.hpp:50
Eigen::VectorXd x
Decision vector.
Definition: box-qp.hpp:51
std::vector< size_t > clamped_idx
Clamped space indexes.
Definition: box-qp.hpp:53