9 #include "crocoddyl/core/solvers/box-qp.hpp"
14 const double th_acceptstep,
const double th_grad,
const double reg)
17 th_acceptstep_(th_acceptstep),
31 if (0. >= th_acceptstep && th_acceptstep >= 0.5) {
32 std::cerr <<
"Warning: th_acceptstep value should between 0 and 0.5"
36 std::cerr <<
"Warning: th_grad value has to be positive." << std::endl;
39 std::cerr <<
"Warning: reg value has to be positive." << std::endl;
53 solution_.
x = Eigen::VectorXd::Zero(nx);
56 const std::size_t n_alphas_ = 10;
57 alphas_.resize(n_alphas_);
58 for (std::size_t n = 0; n < n_alphas_; ++n) {
59 alphas_[n] = 1. / pow(2.,
static_cast<double>(n));
66 const Eigen::VectorXd& q,
67 const Eigen::VectorXd& lb,
68 const Eigen::VectorXd& ub,
69 const Eigen::VectorXd& xinit) {
70 if (
static_cast<std::size_t
>(H.rows()) != nx_ ||
71 static_cast<std::size_t
>(H.cols()) != nx_) {
72 throw_pretty(
"Invalid argument: "
73 <<
"H has wrong dimension (it should be " +
74 std::to_string(nx_) +
"," + std::to_string(nx_) +
")");
76 if (
static_cast<std::size_t
>(q.size()) != nx_) {
78 "Invalid argument: " <<
"q has wrong dimension (it should be " +
79 std::to_string(nx_) +
")");
81 if (
static_cast<std::size_t
>(lb.size()) != nx_) {
83 "Invalid argument: " <<
"lb has wrong dimension (it should be " +
84 std::to_string(nx_) +
")");
86 if (
static_cast<std::size_t
>(ub.size()) != nx_) {
88 "Invalid argument: " <<
"ub has wrong dimension (it should be " +
89 std::to_string(nx_) +
")");
91 if (
static_cast<std::size_t
>(xinit.size()) != nx_) {
93 "Invalid argument: " <<
"xinit has wrong dimension (it should be " +
94 std::to_string(nx_) +
")");
98 for (std::size_t i = 0; i < nx_; ++i) {
99 x_(i) = std::max(std::min(xinit(i), ub(i)), lb(i));
103 for (std::size_t k = 0; k < maxiter_; ++k) {
108 g_.noalias() += H * x_;
109 for (std::size_t j = 0; j < nx_; ++j) {
110 const double gj = g_(j);
111 const double xj = x_(j);
112 const double lbj = lb(j);
113 const double ubj = ub(j);
114 if ((xj == lbj && gj > 0.) || (xj == ubj && gj < 0.)) {
124 Eigen::VectorBlock<Eigen::VectorXd> xf = xo_.head(nf_);
125 Eigen::VectorBlock<Eigen::VectorXd> xc = xo_.tail(nc_);
126 Eigen::VectorBlock<Eigen::VectorXd> dxf = dxo_.head(nf_);
127 Eigen::VectorBlock<Eigen::VectorXd> qf = qo_.head(nf_);
128 Eigen::Block<Eigen::MatrixXd> Hff = Ho_.topLeftCorner(nf_, nf_);
129 Eigen::Block<Eigen::MatrixXd> Hfc = Ho_.topRightCorner(nf_, nc_);
130 for (std::size_t i = 0; i < nf_; ++i) {
131 const std::size_t fi = solution_.
free_idx[i];
134 for (std::size_t j = 0; j < nf_; ++j) {
135 Hff(i, j) = H(fi, solution_.
free_idx[j]);
137 for (std::size_t j = 0; j < nc_; ++j) {
140 Hfc(i, j) = H(fi, cj);
144 Hff.diagonal().array() += reg_;
146 Hff_inv_llt_.compute(Hff);
147 const Eigen::ComputationInfo& info = Hff_inv_llt_.info();
148 if (info != Eigen::Success) {
149 throw_pretty(
"backward_error");
151 solution_.
Hff_inv.setIdentity(nf_, nf_);
152 Hff_inv_llt_.solveInPlace(solution_.
Hff_inv);
153 qf.noalias() += Hff * xf;
155 qf.noalias() += Hfc * xc;
158 Hff_inv_llt_.solveInPlace(dxf);
160 for (std::size_t i = 0; i < nf_; ++i) {
161 dx_(solution_.
free_idx[i]) = dxf(i);
166 fold_ = 0.5 * x_.dot(H * x_) + q.dot(x_);
167 for (std::vector<double>::const_iterator it = alphas_.begin();
168 it != alphas_.end(); ++it) {
169 double steplength = *it;
170 for (std::size_t i = 0; i < nx_; ++i) {
172 std::max(std::min(x_(i) + steplength * dx_(i), ub(i)), lb(i));
174 fnew_ = 0.5 * xnew_.dot(H * xnew_) + q.dot(xnew_);
175 if (fold_ - fnew_ > th_acceptstep_ * g_.dot(x_ - xnew_)) {
182 if (qf.lpNorm<Eigen::Infinity>() <= th_grad_) {
207 x_.conservativeResize(nx);
208 xnew_.conservativeResize(nx);
209 g_.conservativeResize(nx);
210 dx_.conservativeResize(nx);
211 xo_.conservativeResize(nx);
212 dxo_.conservativeResize(nx);
213 qo_.conservativeResize(nx);
214 Ho_.conservativeResize(nx, nx);
220 if (0. >= th_acceptstep && th_acceptstep >= 0.5) {
222 "Invalid argument: " <<
"th_acceptstep value should between 0 and 0.5");
224 th_acceptstep_ = th_acceptstep;
229 throw_pretty(
"Invalid argument: " <<
"th_grad value has to be positive.");
236 throw_pretty(
"Invalid argument: " <<
"reg value has to be positive.");
242 double prev_alpha = alphas[0];
243 if (prev_alpha != 1.) {
244 std::cerr <<
"Warning: alpha[0] should be 1" << std::endl;
246 for (std::size_t i = 1; i < alphas.size(); ++i) {
247 double alpha = alphas[i];
249 throw_pretty(
"Invalid argument: " <<
"alpha values has to be positive.");
251 if (alpha >= prev_alpha) {
253 "Invalid argument: " <<
"alpha values are monotonously decreasing.");
void set_reg(const double reg)
Modify the regularization value.
void set_th_grad(const double th_grad)
Modify the gradient tolerance threshold.
const BoxQPSolution & get_solution() const
Return the stored solution.
double get_th_grad() const
Return the gradient tolerance threshold.
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.
std::size_t get_maxiter() const
Return the maximum allowed number of iterations.
void set_maxiter(const std::size_t maxiter)
Modify the maximum allowed number of iterations.
~BoxQP()
Destroy the Projected-Newton QP solver.
void set_nx(const std::size_t nx)
Modify the decision vector dimension.
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.
const std::vector< double > & get_alphas() const
Return the stack of step lengths using by the line-search procedure.
void set_alphas(const std::vector< double > &alphas)
Modify the stack of step lengths using by the line-search procedure.
void set_th_acceptstep(const double th_acceptstep)
Modify the acceptance step threshold.
double get_th_acceptstep() const
Return the acceptance step threshold.
double get_reg() const
Return the regularization value.
std::size_t get_nx() const
Return the decision vector dimension.
std::vector< size_t > free_idx
Free space indexes.
Eigen::MatrixXd Hff_inv
Inverse of the free space Hessian.
Eigen::VectorXd x
Decision vector.
std::vector< size_t > clamped_idx
Clamped space indexes.