Crocoddyl
ipopt-iface.hpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2022-2023, IRI: CSIC-UPC, Heriot-Watt University
5 // Copyright note valid unless otherwise stated in individual files.
6 // All rights reserved.
8 
9 #ifndef __CROCODDYL_CORE_SOLVERS_IPOPT_IPOPT_IFACE_HPP__
10 #define __CROCODDYL_CORE_SOLVERS_IPOPT_IPOPT_IFACE_HPP__
11 
12 #define HAVE_CSTDDEF
13 #include <IpTNLP.hpp>
14 #undef HAVE_CSTDDEF
15 
16 #include "crocoddyl/core/mathbase.hpp"
17 #include "crocoddyl/core/optctrl/shooting.hpp"
18 
19 namespace crocoddyl {
20 
21 struct IpoptInterfaceData;
22 
63 class IpoptInterface : public Ipopt::TNLP {
64  public:
65  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
66 
72  IpoptInterface(const boost::shared_ptr<crocoddyl::ShootingProblem>& problem);
73 
74  virtual ~IpoptInterface();
75 
91  virtual bool get_nlp_info(Ipopt::Index& n, Ipopt::Index& m,
92  Ipopt::Index& nnz_jac_g, Ipopt::Index& nnz_h_lag,
93  IndexStyleEnum& index_style);
94 
122  virtual bool get_bounds_info(Ipopt::Index n, Ipopt::Number* x_l,
123  Ipopt::Number* x_u, Ipopt::Index m,
124  Ipopt::Number* g_l, Ipopt::Number* g_u);
125 
158  // [TNLP_get_starting_point]
159  virtual bool get_starting_point(Ipopt::Index n, bool init_x, Ipopt::Number* x,
160  bool init_z, Ipopt::Number* z_L,
161  Ipopt::Number* z_U, Ipopt::Index m,
162  bool init_lambda, Ipopt::Number* lambda);
163 
181  virtual bool eval_f(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
182  Ipopt::Number& obj_value);
183 
201  virtual bool eval_grad_f(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
202  Ipopt::Number* grad_f);
203 
221  virtual bool eval_g(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
222  Ipopt::Index m, Ipopt::Number* g);
223 
267  virtual bool eval_jac_g(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
268  Ipopt::Index m, Ipopt::Index nele_jac,
269  Ipopt::Index* iRow, Ipopt::Index* jCol,
270  Ipopt::Number* values);
271 
326  virtual bool eval_h(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
327  Ipopt::Number obj_factor, Ipopt::Index m,
328  const Ipopt::Number* lambda, bool new_lambda,
329  Ipopt::Index nele_hess, Ipopt::Index* iRow,
330  Ipopt::Index* jCol, Ipopt::Number* values);
331 
377  virtual void finalize_solution(
378  Ipopt::SolverReturn status, Ipopt::Index n, const Ipopt::Number* x,
379  const Ipopt::Number* z_L, const Ipopt::Number* z_U, Ipopt::Index m,
380  const Ipopt::Number* g, const Ipopt::Number* lambda,
381  Ipopt::Number obj_value, const Ipopt::IpoptData* ip_data,
382  Ipopt::IpoptCalculatedQuantities* ip_cq);
383 
407  Ipopt::AlgorithmMode mode, Ipopt::Index iter, Ipopt::Number obj_value,
408  Ipopt::Number inf_pr, Ipopt::Number inf_du, Ipopt::Number mu,
409  Ipopt::Number d_norm, Ipopt::Number regularization_size,
410  Ipopt::Number alpha_du, Ipopt::Number alpha_pr, Ipopt::Index ls_trials,
411  const Ipopt::IpoptData* ip_data, Ipopt::IpoptCalculatedQuantities* ip_cq);
412 
418  boost::shared_ptr<IpoptInterfaceData> createData(const std::size_t nx,
419  const std::size_t ndx,
420  const std::size_t nu);
421 
422  void resizeData();
423 
428  std::size_t get_nvar() const;
429 
433  std::size_t get_nconst() const;
434 
438  const std::vector<Eigen::VectorXd>& get_xs() const;
439 
443  const std::vector<Eigen::VectorXd>& get_us() const;
444 
448  const boost::shared_ptr<crocoddyl::ShootingProblem>& get_problem() const;
449 
450  double get_cost() const;
451 
455  void set_xs(const std::vector<Eigen::VectorXd>& xs);
456 
460  void set_us(const std::vector<Eigen::VectorXd>& us);
461 
462  private:
463  boost::shared_ptr<crocoddyl::ShootingProblem>
464  problem_;
465  std::vector<Eigen::VectorXd> xs_;
466  std::vector<Eigen::VectorXd> us_;
467  std::vector<std::size_t> ixu_;
468  std::size_t nvar_;
469  std::size_t nconst_;
470  std::vector<boost::shared_ptr<IpoptInterfaceData>>
471  datas_;
472  double cost_;
473 
475 
476  IpoptInterface& operator=(const IpoptInterface&);
477 };
478 
480  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
481 
482  IpoptInterfaceData(const std::size_t nx, const std::size_t ndx,
483  const std::size_t nu)
484  : x(nx),
485  xnext(nx),
486  dx(ndx),
487  dxnext(ndx),
488  x_diff(ndx),
489  u(nu),
490  Jint_dx(ndx, ndx),
491  Jint_dxnext(ndx, ndx),
492  Jdiff_x(ndx, ndx),
493  Jdiff_xnext(ndx, ndx),
494  Jg_dx(ndx, ndx),
495  Jg_dxnext(ndx, ndx),
496  Jg_u(ndx, ndx),
497  Jg_ic(ndx, ndx),
498  FxJint_dx(ndx, ndx),
499  Ldx(ndx),
500  Ldxdx(ndx, ndx),
501  Ldxu(ndx, nu) {
502  x.setZero();
503  xnext.setZero();
504  dx.setZero();
505  dxnext.setZero();
506  x_diff.setZero();
507  u.setZero();
508  Jint_dx.setZero();
509  Jint_dxnext.setZero();
510  Jdiff_x.setZero();
511  Jdiff_xnext.setZero();
512  Jg_dx.setZero();
513  Jg_dxnext.setZero();
514  Jg_u.setZero();
515  Jg_ic.setZero();
516  FxJint_dx.setZero();
517  Ldx.setZero();
518  Ldxdx.setZero();
519  Ldxu.setZero();
520  }
521 
522  void resize(const std::size_t nx, const std::size_t ndx,
523  const std::size_t nu) {
524  x.conservativeResize(nx);
525  xnext.conservativeResize(nx);
526  dx.conservativeResize(ndx);
527  dxnext.conservativeResize(ndx);
528  x_diff.conservativeResize(ndx);
529  u.conservativeResize(nu);
530  Jint_dx.conservativeResize(ndx, ndx);
531  Jint_dxnext.conservativeResize(ndx, ndx);
532  Jdiff_x.conservativeResize(ndx, ndx);
533  Jdiff_xnext.conservativeResize(ndx, ndx);
534  Jg_dx.conservativeResize(ndx, ndx);
535  Jg_dxnext.conservativeResize(ndx, ndx);
536  Jg_u.conservativeResize(ndx, ndx);
537  Jg_ic.conservativeResize(ndx, ndx);
538  FxJint_dx.conservativeResize(ndx, ndx);
539  Ldx.conservativeResize(ndx);
540  Ldxdx.conservativeResize(ndx, ndx);
541  Ldxu.conservativeResize(ndx, nu);
542  }
543 
544  Eigen::VectorXd x;
545  Eigen::VectorXd xnext;
546  Eigen::VectorXd dx;
547  Eigen::VectorXd dxnext;
548  Eigen::VectorXd x_diff;
549  Eigen::VectorXd u;
550  Eigen::MatrixXd Jint_dx;
551  Eigen::MatrixXd
553  Eigen::MatrixXd
555  Eigen::MatrixXd Jdiff_xnext;
557  Eigen::MatrixXd Jg_dx;
558  Eigen::MatrixXd
560  Eigen::MatrixXd Jg_u;
561  Eigen::MatrixXd
563  Eigen::MatrixXd FxJint_dx;
564  Eigen::VectorXd Ldx;
565  Eigen::MatrixXd Ldxdx;
566  Eigen::MatrixXd Ldxu;
567 };
568 
569 } // namespace crocoddyl
570 
571 #endif
Class for interfacing a crocoddyl::ShootingProblem with IPOPT.
Definition: ipopt-iface.hpp:63
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.
Definition: ipopt-iface.cpp:91
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.
EIGEN_MAKE_ALIGNED_OPERATOR_NEW IpoptInterface(const boost::shared_ptr< crocoddyl::ShootingProblem > &problem)
Initialize the Ipopt interface.
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.
void set_us(const std::vector< Eigen::VectorXd > &us)
Modify the control vector.
const boost::shared_ptr< crocoddyl::ShootingProblem > & get_problem() const
Return the crocoddyl::ShootingProblem to be solved.
boost::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.
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.
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::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.