Crocoddyl
kkt.hpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2019-2020, LAAS-CNRS, New York University, Max Planck
5 // Gesellschaft,
6 // University of Edinburgh
7 // Copyright note valid unless otherwise stated in individual files.
8 // All rights reserved.
10 
11 #ifndef CROCODDYL_CORE_SOLVERS_KKT_HPP_
12 #define CROCODDYL_CORE_SOLVERS_KKT_HPP_
13 
14 #include <Eigen/Cholesky>
15 #include <Eigen/Dense>
16 
17 #include "crocoddyl/core/solver-base.hpp"
18 
19 namespace crocoddyl {
20 
21 class SolverKKT : public SolverAbstract {
22  public:
23  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
24 
25  explicit SolverKKT(boost::shared_ptr<ShootingProblem> problem);
26  virtual ~SolverKKT();
27 
28  virtual bool solve(
29  const std::vector<Eigen::VectorXd>& init_xs = DEFAULT_VECTOR,
30  const std::vector<Eigen::VectorXd>& init_us = DEFAULT_VECTOR,
31  const std::size_t maxiter = 100, const bool is_feasible = false,
32  const double regInit = 1e-9);
33  virtual void computeDirection(const bool recalc = true);
34  virtual double tryStep(const double steplength = 1);
35  virtual double stoppingCriteria();
36  virtual const Eigen::Vector2d& expectedImprovement();
37 
38  const Eigen::MatrixXd& get_kkt() const;
39  const Eigen::VectorXd& get_kktref() const;
40  const Eigen::VectorXd& get_primaldual() const;
41  const std::vector<Eigen::VectorXd>& get_dxs() const;
42  const std::vector<Eigen::VectorXd>& get_dus() const;
43  const std::vector<Eigen::VectorXd>& get_lambdas() const;
44  std::size_t get_nx() const;
45  std::size_t get_ndx() const;
46  std::size_t get_nu() const;
47 
48  protected:
49  double reg_incfactor_;
50  double reg_decfactor_;
51  double reg_min_;
52  double reg_max_;
53  double cost_try_;
54  std::vector<Eigen::VectorXd> xs_try_;
55  std::vector<Eigen::VectorXd> us_try_;
56 
57  private:
58  double calcDiff();
59  void computePrimalDual();
60  void increaseRegularization();
61  void decreaseRegularization();
62  void allocateData();
63 
64  std::size_t nx_;
65  std::size_t ndx_;
66  std::size_t nu_;
67  std::vector<Eigen::VectorXd> dxs_;
68  std::vector<Eigen::VectorXd> dus_;
69  std::vector<Eigen::VectorXd> lambdas_;
70 
71  // allocate data
72  Eigen::MatrixXd kkt_;
73  Eigen::VectorXd kktref_;
74  Eigen::VectorXd primaldual_;
75  Eigen::VectorXd primal_;
76  Eigen::VectorXd dual_;
77  std::vector<double> alphas_;
78  double th_grad_;
79  bool was_feasible_;
80  Eigen::VectorXd kkt_primal_;
81  Eigen::VectorXd dF;
82 };
83 
84 } // namespace crocoddyl
85 
86 #endif // CROCODDYL_CORE_SOLVERS_KKT_HPP_
Abstract class for optimal control solvers.
Definition: solver-base.hpp:61
virtual bool solve(const std::vector< Eigen::VectorXd > &init_xs=DEFAULT_VECTOR, const std::vector< Eigen::VectorXd > &init_us=DEFAULT_VECTOR, const std::size_t maxiter=100, const bool is_feasible=false, const double regInit=1e-9)
Compute the optimal trajectory as lists of and terms.
Definition: kkt.cpp:36
virtual double tryStep(const double steplength=1)
Try a predefined step length and compute its cost improvement .
Definition: kkt.cpp:119
virtual const Eigen::Vector2d & expectedImprovement()
Return the expected improvement from a given current search direction .
Definition: kkt.cpp:166
virtual double stoppingCriteria()
Return a positive value that quantifies the algorithm termination.
Definition: kkt.cpp:139
virtual void computeDirection(const bool recalc=true)
Compute the search direction for the current guess .
Definition: kkt.cpp:89