Crocoddyl
 
Loading...
Searching...
No Matches
box-qp.hpp
1
2// 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#ifndef CROCODDYL_CORE_SOLVERS_BOX_QP_HPP_
10#define CROCODDYL_CORE_SOLVERS_BOX_QP_HPP_
11
12#include "crocoddyl/core/fwd.hpp"
13#include "crocoddyl/core/utils/exception.hpp"
14
15namespace crocoddyl {
16
27 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
28
33
42 BoxQPSolution(const Eigen::MatrixXd& Hff_inv, const Eigen::VectorXd& x,
43 const std::vector<size_t>& free_idx,
44 const std::vector<size_t>& clamped_idx)
46
47 Eigen::MatrixXd Hff_inv;
48 Eigen::VectorXd x;
49 std::vector<size_t> free_idx;
50 std::vector<size_t> clamped_idx;
51};
52
79class BoxQP {
80 public:
81 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
82
93 BoxQP(const std::size_t nx, const std::size_t maxiter = 100,
94 const double th_acceptstep = 0.1, const double th_grad = 1e-9,
95 const double reg = 1e-9);
99 ~BoxQP();
100
112 const BoxQPSolution& solve(const Eigen::MatrixXd& H, const Eigen::VectorXd& q,
113 const Eigen::VectorXd& lb,
114 const Eigen::VectorXd& ub,
115 const Eigen::VectorXd& xinit);
116
120 const BoxQPSolution& get_solution() const;
121
125 std::size_t get_nx() const;
126
130 std::size_t get_maxiter() const;
131
135 double get_th_acceptstep() const;
136
140 double get_th_grad() const;
141
145 double get_reg() const;
146
150 const std::vector<double>& get_alphas() const;
151
155 void set_nx(const std::size_t nx);
156
160 void set_maxiter(const std::size_t maxiter);
161
165 void set_th_acceptstep(const double th_acceptstep);
166
170 void set_th_grad(const double th_grad);
171
175 void set_reg(const double reg);
176
180 void set_alphas(const std::vector<double>& alphas);
181
182 private:
183 std::size_t nx_;
184 BoxQPSolution solution_;
185 std::size_t maxiter_;
186 double th_acceptstep_;
187 double
188 th_grad_;
189 double reg_;
190
191 double fold_;
192 double fnew_;
193 std::size_t nf_;
194 std::size_t nc_;
195 std::vector<double>
196 alphas_;
197 Eigen::VectorXd x_;
198 Eigen::VectorXd xnew_;
199 Eigen::VectorXd g_;
200 Eigen::VectorXd dx_;
201
202 Eigen::VectorXd xo_;
203 Eigen::VectorXd
204 dxo_;
205 Eigen::VectorXd
206 qo_;
207 Eigen::MatrixXd Ho_;
208
209 Eigen::LLT<Eigen::MatrixXd> Hff_inv_llt_;
210};
211
212} // namespace crocoddyl
213
214#endif // CROCODDYL_CORE_SOLVERS_BOX_QP_HPP_
This class implements a Box QP solver based on a Projected Newton method.
Definition box-qp.hpp:79
void set_reg(const double reg)
Modify the regularization value.
Definition box-qp.cpp:234
void set_th_grad(const double th_grad)
Modify the gradient tolerance threshold.
Definition box-qp.cpp:227
const BoxQPSolution & get_solution() const
Return the stored solution.
Definition box-qp.cpp:191
double get_th_grad() const
Return the gradient tolerance threshold.
Definition box-qp.cpp:199
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:65
std::size_t get_maxiter() const
Return the maximum allowed number of iterations.
Definition box-qp.cpp:195
void set_maxiter(const std::size_t maxiter)
Modify the maximum allowed number of iterations.
Definition box-qp.cpp:217
~BoxQP()
Destroy the Projected-Newton QP solver.
Definition box-qp.cpp:63
void set_nx(const std::size_t nx)
Modify the decision vector dimension.
Definition box-qp.cpp:205
const std::vector< double > & get_alphas() const
Return the stack of step lengths using by the line-search procedure.
Definition box-qp.cpp:203
void set_alphas(const std::vector< double > &alphas)
Modify the stack of step lengths using by the line-search procedure.
Definition box-qp.cpp:241
void set_th_acceptstep(const double th_acceptstep)
Modify the acceptance step threshold.
Definition box-qp.cpp:219
double get_th_acceptstep() const
Return the acceptance step threshold.
Definition box-qp.cpp:197
double get_reg() const
Return the regularization value.
Definition box-qp.cpp:201
std::size_t get_nx() const
Return the decision vector dimension.
Definition box-qp.cpp:193
Box QP solution.
Definition box-qp.hpp:26
std::vector< size_t > free_idx
Free space indexes.
Definition box-qp.hpp:49
EIGEN_MAKE_ALIGNED_OPERATOR_NEW BoxQPSolution()
Initialize the QP solution structure.
Definition box-qp.hpp:32
Eigen::MatrixXd Hff_inv
Inverse of the free space Hessian.
Definition box-qp.hpp:47
Eigen::VectorXd x
Decision vector.
Definition box-qp.hpp:48
BoxQPSolution(const Eigen::MatrixXd &Hff_inv, const Eigen::VectorXd &x, const std::vector< size_t > &free_idx, const std::vector< size_t > &clamped_idx)
Initialize the QP solution structure.
Definition box-qp.hpp:42
std::vector< size_t > clamped_idx
Clamped space indexes.
Definition box-qp.hpp:50