30 #ifndef EXOTICA_CORE_BOX_QP_H_
31 #define EXOTICA_CORE_BOX_QP_H_
34 #include <Eigen/Dense>
47 inline BoxQPSolution BoxQP(
const Eigen::MatrixXd& H,
const Eigen::VectorXd& q,
const Eigen::VectorXd& b_low,
const Eigen::VectorXd& b_high,
const Eigen::VectorXd& x_init,
const double th_acceptstep,
const int max_iterations,
const double th_gradient_tolerance,
const double lambda,
bool use_polynomial_linesearch =
true,
bool use_cholesky_factorization =
true)
49 if (lambda < 0.)
ThrowPretty(
"lambda needs to be positive.");
54 const std::size_t nx = x_init.size();
56 Eigen::VectorXd delta_xf(nx),
x = x_init;
57 std::vector<size_t> clamped_idx, free_idx;
58 Eigen::VectorXd grad(nx);
59 Eigen::MatrixXd Hff(nx, nx), Hfc(nx, nx);
61 std::vector<double> alphas_;
62 const std::size_t& n_alphas_ = 10;
63 alphas_.resize(n_alphas_);
64 const Eigen::VectorXd alphas_linear = Eigen::VectorXd::LinSpaced(n_alphas_, 1.0, 0.1);
65 for (std::size_t
n = 0;
n < n_alphas_; ++
n)
67 if (use_polynomial_linesearch)
69 alphas_[
n] = 1. / pow(2.,
static_cast<double>(
n));
73 alphas_[
n] = alphas_linear(
n);
79 for (
int i = 0; i <
x.size(); ++i)
81 x(i) = std::max(std::min(x_init(i), b_high(i)), b_low(i));
87 std::size_t num_free, num_clamped;
88 Eigen::LLT<Eigen::MatrixXd> Hff_inv_llt_;
89 Eigen::VectorXd qf_, xf_, xc_, dxf_, dx_(nx);
91 for (
int k = 0; k < max_iterations; ++k)
98 grad.noalias() += H *
x;
101 for (std::size_t i = 0; i < nx; ++i)
103 if ((
x(i) == b_low(i) && grad(i) > 0.) || (
x(i) == b_high(i) && grad(i) < 0.))
108 num_free = solution.
free_idx.size();
115 if (grad.lpNorm<Eigen::Infinity>() <= th_gradient_tolerance || num_free == 0)
118 if (k == 0 && num_free != 0)
120 Hff.resize(num_free, num_free);
121 for (std::size_t i = 0; i < num_free; ++i)
123 const std::size_t& fi = solution.
free_idx[i];
124 for (std::size_t j = 0; j < num_free; ++j)
126 Hff(i, j) = H(fi, solution.
free_idx[j]);
132 Hff.diagonal().array() += lambda;
136 if (use_cholesky_factorization)
138 Hff_inv_llt_.compute(Hff);
139 const Eigen::ComputationInfo& info = Hff_inv_llt_.info();
140 if (info != Eigen::Success)
142 ThrowPretty(
"Error during Cholesky decomposition of Hff");
144 solution.
Hff_inv.setIdentity(num_free, num_free);
145 Hff_inv_llt_.solveInPlace(solution.
Hff_inv);
149 solution.
Hff_inv = Hff.inverse();
155 solution.
Hff_inv.resize(num_free, num_free);
164 qf_.resize(num_free);
165 xf_.resize(num_free);
166 xc_.resize(num_clamped);
167 dxf_.resize(num_free);
168 Hff.resize(num_free, num_free);
169 Hfc.resize(num_free, num_clamped);
170 for (std::size_t i = 0; i < num_free; ++i)
172 const std::size_t& fi = solution.
free_idx[i];
175 for (std::size_t j = 0; j < num_free; ++j)
177 Hff(i, j) = H(fi, solution.
free_idx[j]);
179 for (std::size_t j = 0; j < num_clamped; ++j)
183 Hfc(i, j) = H(fi, cj);
188 Hff.diagonal().array() += lambda;
192 if (use_cholesky_factorization)
194 Hff_inv_llt_.compute(Hff);
195 const Eigen::ComputationInfo& info = Hff_inv_llt_.info();
196 if (info != Eigen::Success)
198 ThrowPretty(
"Error during Cholesky decomposition of Hff (iter=" << k <<
"):\n"
201 << H <<
"\nnum_free: " << num_free <<
" num_clamped: " << num_clamped <<
" lambda: " << lambda);
203 solution.
Hff_inv.setIdentity(num_free, num_free);
204 Hff_inv_llt_.solveInPlace(solution.
Hff_inv);
208 solution.
Hff_inv = Hff.inverse();
212 if (num_clamped != 0)
214 dxf_.noalias() -= Hfc * xc_;
216 if (use_cholesky_factorization)
218 Hff_inv_llt_.solveInPlace(dxf_);
222 dxf_ = solution.
Hff_inv * dxf_;
226 for (std::size_t i = 0; i < num_free; ++i)
228 dx_(solution.
free_idx[i]) = dxf_(i);
232 Eigen::VectorXd xnew_(nx);
233 fold_ = 0.5 *
x.dot(H *
x) + q.dot(
x);
234 for (std::vector<double>::const_iterator it = alphas_.begin(); it != alphas_.end(); ++it)
236 double steplength = *it;
237 for (std::size_t i = 0; i < nx; ++i)
239 xnew_(i) = std::max(std::min(
x(i) + steplength * dx_(i), b_high(i)), b_low(i));
241 fnew_ = 0.5 * xnew_.dot(H * xnew_) + q.dot(xnew_);
242 if (fold_ - fnew_ > th_acceptstep * grad.dot(
x - xnew_))
249 if (it == alphas_.end() - 1)
262 const Eigen::VectorXd& b_low,
const Eigen::VectorXd& b_high,
263 const Eigen::VectorXd& x_init)
265 constexpr
double epsilon = 1e-5;
266 constexpr
double gamma = 0.1;
267 constexpr
int max_iterations = 100;
268 constexpr
double lambda = 1e-5;
269 return BoxQP(H, q, b_low, b_high, x_init, gamma, max_iterations,
epsilon, lambda,
true,
true);
273 #endif // EXOTICA_CORE_BOX_QP_H_