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/optctrl/shooting.hpp"
17 
18 namespace crocoddyl {
19 
20 struct IpoptInterfaceData;
21 
62 class IpoptInterface : public Ipopt::TNLP {
63  public:
64  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
65 
71  IpoptInterface(const std::shared_ptr<crocoddyl::ShootingProblem>& problem);
72 
73  virtual ~IpoptInterface();
74 
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);
93 
121  virtual bool get_bounds_info(Ipopt::Index n, Ipopt::Number* x_l,
122  Ipopt::Number* x_u, Ipopt::Index m,
123  Ipopt::Number* g_l, Ipopt::Number* g_u);
124 
157  // [TNLP_get_starting_point]
158  virtual bool get_starting_point(Ipopt::Index n, bool init_x, Ipopt::Number* x,
159  bool init_z, Ipopt::Number* z_L,
160  Ipopt::Number* z_U, Ipopt::Index m,
161  bool init_lambda, Ipopt::Number* lambda);
162 
180  virtual bool eval_f(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
181  Ipopt::Number& obj_value);
182 
200  virtual bool eval_grad_f(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
201  Ipopt::Number* grad_f);
202 
220  virtual bool eval_g(Ipopt::Index n, const Ipopt::Number* x, bool new_x,
221  Ipopt::Index m, Ipopt::Number* g);
222 
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);
270 
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);
330 
376  virtual void finalize_solution(
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);
382 
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);
411 
417  std::shared_ptr<IpoptInterfaceData> createData(const std::size_t nx,
418  const std::size_t ndx,
419  const std::size_t nu);
420 
421  void resizeData();
422 
427  std::size_t get_nvar() const;
428 
432  std::size_t get_nconst() const;
433 
437  const std::vector<Eigen::VectorXd>& get_xs() const;
438 
442  const std::vector<Eigen::VectorXd>& get_us() const;
443 
447  const std::shared_ptr<crocoddyl::ShootingProblem>& get_problem() const;
448 
449  double get_cost() const;
450 
454  void set_xs(const std::vector<Eigen::VectorXd>& xs);
455 
459  void set_us(const std::vector<Eigen::VectorXd>& us);
460 
461  private:
462  std::shared_ptr<crocoddyl::ShootingProblem>
463  problem_;
464  std::vector<Eigen::VectorXd> xs_;
465  std::vector<Eigen::VectorXd> us_;
466  std::vector<std::size_t> ixu_;
467  std::size_t nvar_;
468  std::size_t nconst_;
469  std::vector<std::shared_ptr<IpoptInterfaceData>> datas_;
470  double cost_;
471 
473 
474  IpoptInterface& operator=(const IpoptInterface&);
475 };
476 
478  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
479 
480  IpoptInterfaceData(const std::size_t nx, const std::size_t ndx,
481  const std::size_t nu)
482  : x(nx),
483  xnext(nx),
484  dx(ndx),
485  dxnext(ndx),
486  x_diff(ndx),
487  u(nu),
488  Jint_dx(ndx, ndx),
489  Jint_dxnext(ndx, ndx),
490  Jdiff_x(ndx, ndx),
491  Jdiff_xnext(ndx, ndx),
492  Jg_dx(ndx, ndx),
493  Jg_dxnext(ndx, ndx),
494  Jg_u(ndx, ndx),
495  Jg_ic(ndx, ndx),
496  FxJint_dx(ndx, ndx),
497  Ldx(ndx),
498  Ldxdx(ndx, ndx),
499  Ldxu(ndx, nu) {
500  x.setZero();
501  xnext.setZero();
502  dx.setZero();
503  dxnext.setZero();
504  x_diff.setZero();
505  u.setZero();
506  Jint_dx.setZero();
507  Jint_dxnext.setZero();
508  Jdiff_x.setZero();
509  Jdiff_xnext.setZero();
510  Jg_dx.setZero();
511  Jg_dxnext.setZero();
512  Jg_u.setZero();
513  Jg_ic.setZero();
514  FxJint_dx.setZero();
515  Ldx.setZero();
516  Ldxdx.setZero();
517  Ldxu.setZero();
518  }
519 
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);
529  Jint_dxnext.conservativeResize(ndx, ndx);
530  Jdiff_x.conservativeResize(ndx, ndx);
531  Jdiff_xnext.conservativeResize(ndx, ndx);
532  Jg_dx.conservativeResize(ndx, ndx);
533  Jg_dxnext.conservativeResize(ndx, ndx);
534  Jg_u.conservativeResize(ndx, ndx);
535  Jg_ic.conservativeResize(ndx, ndx);
536  FxJint_dx.conservativeResize(ndx, ndx);
537  Ldx.conservativeResize(ndx);
538  Ldxdx.conservativeResize(ndx, ndx);
539  Ldxu.conservativeResize(ndx, nu);
540  }
541 
542  Eigen::VectorXd x;
543  Eigen::VectorXd xnext;
544  Eigen::VectorXd dx;
545  Eigen::VectorXd dxnext;
546  Eigen::VectorXd x_diff;
547  Eigen::VectorXd u;
548  Eigen::MatrixXd Jint_dx;
549  Eigen::MatrixXd
551  Eigen::MatrixXd
553  Eigen::MatrixXd Jdiff_xnext;
555  Eigen::MatrixXd Jg_dx;
556  Eigen::MatrixXd
558  Eigen::MatrixXd Jg_u;
559  Eigen::MatrixXd
561  Eigen::MatrixXd FxJint_dx;
562  Eigen::VectorXd Ldx;
563  Eigen::MatrixXd Ldxdx;
564  Eigen::MatrixXd Ldxu;
565 };
566 
567 } // namespace crocoddyl
568 
569 #endif
Class for interfacing a crocoddyl::ShootingProblem with IPOPT.
Definition: ipopt-iface.hpp:62
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:86
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::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.