9 #ifndef __CROCODDYL_CORE_SOLVERS_IPOPT_IPOPT_IFACE_HPP__
10 #define __CROCODDYL_CORE_SOLVERS_IPOPT_IPOPT_IFACE_HPP__
16 #include "crocoddyl/core/optctrl/shooting.hpp"
20 struct IpoptInterfaceData;
64 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
90 virtual bool get_nlp_info(Ipopt::Index& n, Ipopt::Index& m,
91 Ipopt::Index& nnz_jac_g, Ipopt::Index& nnz_h_lag,
92 IndexStyleEnum& index_style);
122 Ipopt::Number* x_u, Ipopt::Index m,
123 Ipopt::Number* g_l, Ipopt::Number* g_u);
159 bool init_z, Ipopt::Number* z_L,
160 Ipopt::Number* z_U, Ipopt::Index m,
161 bool init_lambda, Ipopt::Number* lambda);
180 virtual bool eval_f(Ipopt::Index n,
const Ipopt::Number* x,
bool new_x,
181 Ipopt::Number& obj_value);
200 virtual bool eval_grad_f(Ipopt::Index n,
const Ipopt::Number* x,
bool new_x,
201 Ipopt::Number* grad_f);
220 virtual bool eval_g(Ipopt::Index n,
const Ipopt::Number* x,
bool new_x,
221 Ipopt::Index m, Ipopt::Number* g);
266 virtual bool eval_jac_g(Ipopt::Index n,
const Ipopt::Number* x,
bool new_x,
267 Ipopt::Index m, Ipopt::Index nele_jac,
268 Ipopt::Index* iRow, Ipopt::Index* jCol,
269 Ipopt::Number* values);
325 virtual bool eval_h(Ipopt::Index n,
const Ipopt::Number* x,
bool new_x,
326 Ipopt::Number obj_factor, Ipopt::Index m,
327 const Ipopt::Number* lambda,
bool new_lambda,
328 Ipopt::Index nele_hess, Ipopt::Index* iRow,
329 Ipopt::Index* jCol, Ipopt::Number* values);
377 Ipopt::SolverReturn status, Ipopt::Index n,
const Ipopt::Number* x,
378 const Ipopt::Number* z_L,
const Ipopt::Number* z_U, Ipopt::Index m,
379 const Ipopt::Number* g,
const Ipopt::Number* lambda,
380 Ipopt::Number obj_value,
const Ipopt::IpoptData* ip_data,
381 Ipopt::IpoptCalculatedQuantities* ip_cq);
406 Ipopt::AlgorithmMode mode, Ipopt::Index iter, Ipopt::Number obj_value,
407 Ipopt::Number inf_pr, Ipopt::Number inf_du, Ipopt::Number mu,
408 Ipopt::Number d_norm, Ipopt::Number regularization_size,
409 Ipopt::Number alpha_du, Ipopt::Number alpha_pr, Ipopt::Index ls_trials,
410 const Ipopt::IpoptData* ip_data, Ipopt::IpoptCalculatedQuantities* ip_cq);
417 std::shared_ptr<IpoptInterfaceData>
createData(
const std::size_t nx,
418 const std::size_t ndx,
419 const std::size_t nu);
437 const std::vector<Eigen::VectorXd>&
get_xs()
const;
442 const std::vector<Eigen::VectorXd>&
get_us()
const;
447 const std::shared_ptr<crocoddyl::ShootingProblem>&
get_problem()
const;
449 double get_cost()
const;
454 void set_xs(
const std::vector<Eigen::VectorXd>& xs);
459 void set_us(
const std::vector<Eigen::VectorXd>& us);
462 std::shared_ptr<crocoddyl::ShootingProblem>
464 std::vector<Eigen::VectorXd> xs_;
465 std::vector<Eigen::VectorXd> us_;
466 std::vector<std::size_t> ixu_;
469 std::vector<std::shared_ptr<IpoptInterfaceData>> datas_;
478 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
481 const std::size_t nu)
520 void resize(
const std::size_t nx,
const std::size_t ndx,
521 const std::size_t nu) {
522 x.conservativeResize(nx);
523 xnext.conservativeResize(nx);
524 dx.conservativeResize(ndx);
525 dxnext.conservativeResize(ndx);
526 x_diff.conservativeResize(ndx);
527 u.conservativeResize(nu);
528 Jint_dx.conservativeResize(ndx, ndx);
530 Jdiff_x.conservativeResize(ndx, ndx);
532 Jg_dx.conservativeResize(ndx, ndx);
534 Jg_u.conservativeResize(ndx, ndx);
535 Jg_ic.conservativeResize(ndx, ndx);
537 Ldx.conservativeResize(ndx);
538 Ldxdx.conservativeResize(ndx, ndx);
539 Ldxu.conservativeResize(ndx, nu);
Class for interfacing a crocoddyl::ShootingProblem with IPOPT.
virtual bool get_nlp_info(Ipopt::Index &n, Ipopt::Index &m, Ipopt::Index &nnz_jac_g, Ipopt::Index &nnz_h_lag, IndexStyleEnum &index_style)
Methods to gather information about the NLP.
std::size_t get_nvar() const
Return the total number of optimization variables (states and controls)
virtual bool eval_grad_f(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Number *grad_f)
Method to request the gradient of the objective function.
void set_xs(const std::vector< Eigen::VectorXd > &xs)
Modify the state vector.
virtual void finalize_solution(Ipopt::SolverReturn status, Ipopt::Index n, const Ipopt::Number *x, const Ipopt::Number *z_L, const Ipopt::Number *z_U, Ipopt::Index m, const Ipopt::Number *g, const Ipopt::Number *lambda, Ipopt::Number obj_value, const Ipopt::IpoptData *ip_data, Ipopt::IpoptCalculatedQuantities *ip_cq)
This method is called when the algorithm has finished (successfully or not) so the TNLP can digest th...
virtual bool get_starting_point(Ipopt::Index n, bool init_x, Ipopt::Number *x, bool init_z, Ipopt::Number *z_L, Ipopt::Number *z_U, Ipopt::Index m, bool init_lambda, Ipopt::Number *lambda)
Method to request the starting point before iterating.
virtual bool eval_h(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Number obj_factor, Ipopt::Index m, const Ipopt::Number *lambda, bool new_lambda, Ipopt::Index nele_hess, Ipopt::Index *iRow, Ipopt::Index *jCol, Ipopt::Number *values)
Method to request either the sparsity structure or the values of the Hessian of the Lagrangian.
EIGEN_MAKE_ALIGNED_OPERATOR_NEW IpoptInterface(const std::shared_ptr< crocoddyl::ShootingProblem > &problem)
Initialize the Ipopt interface.
void set_us(const std::vector< Eigen::VectorXd > &us)
Modify the control vector.
virtual bool eval_f(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Number &obj_value)
Method to request the value of the objective function.
virtual bool get_bounds_info(Ipopt::Index n, Ipopt::Number *x_l, Ipopt::Number *x_u, Ipopt::Index m, Ipopt::Number *g_l, Ipopt::Number *g_u)
Method to request bounds on the variables and constraints. T.
const std::vector< Eigen::VectorXd > & get_xs() const
Return the state vector.
std::size_t get_nconst() const
Return the total number of constraints in the NLP.
const std::vector< Eigen::VectorXd > & get_us() const
Return the control vector.
virtual bool eval_g(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Index m, Ipopt::Number *g)
Method to request the constraint values.
bool intermediate_callback(Ipopt::AlgorithmMode mode, Ipopt::Index iter, Ipopt::Number obj_value, Ipopt::Number inf_pr, Ipopt::Number inf_du, Ipopt::Number mu, Ipopt::Number d_norm, Ipopt::Number regularization_size, Ipopt::Number alpha_du, Ipopt::Number alpha_pr, Ipopt::Index ls_trials, const Ipopt::IpoptData *ip_data, Ipopt::IpoptCalculatedQuantities *ip_cq)
Intermediate Callback method for the user.
std::shared_ptr< IpoptInterfaceData > createData(const std::size_t nx, const std::size_t ndx, const std::size_t nu)
Create the data structure to store temporary computations.
const std::shared_ptr< crocoddyl::ShootingProblem > & get_problem() const
Return the crocoddyl::ShootingProblem to be solved.
virtual bool eval_jac_g(Ipopt::Index n, const Ipopt::Number *x, bool new_x, Ipopt::Index m, Ipopt::Index nele_jac, Ipopt::Index *iRow, Ipopt::Index *jCol, Ipopt::Number *values)
Method to request either the sparsity structure or the values of the Jacobian of the constraints.
Eigen::MatrixXd Ldxdx
Hessian of the cost w.r.t dxdx.
Eigen::MatrixXd Jg_dx
Jacobian of the dynamic constraint w.r.t dx.
Eigen::MatrixXd FxJint_dx
Intermediate computation needed for Jg_ic.
Eigen::VectorXd xnext
Integrated state at next node.
Eigen::MatrixXd Jint_dxnext
Jacobian of the sum operation w.r.t dx at next node.
Eigen::MatrixXd Jdiff_x
Jacobian of the diff operation w.r.t the first element.
Eigen::MatrixXd Ldxu
Hessian of the cost w.r.t dxu.
Eigen::VectorXd dxnext
Increment in the tangent space at next node.
Eigen::VectorXd Ldx
Jacobian of the cost w.r.t dx.
Eigen::VectorXd x
Integrated state.
Eigen::VectorXd x_diff
State difference.
Eigen::MatrixXd Jdiff_xnext
Eigen::VectorXd dx
Increment in the tangent space.
Eigen::MatrixXd Jg_ic
Jacobian of the initial condition constraint w.r.t dx.
Eigen::MatrixXd Jg_dxnext
Jacobian of the dynamic constraint w.r.t dxnext.
Eigen::MatrixXd Jg_u
Jacobian of the dynamic constraint w.r.t u.
Eigen::MatrixXd Jint_dx
Jacobian of the sum operation w.r.t dx.
Eigen::VectorXd u
Control.