Crocoddyl
 
Loading...
Searching...
No Matches
kkt.hpp
1
2// 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 "crocoddyl/core/solver-base.hpp"
15
16namespace crocoddyl {
17
18class SolverKKT : public SolverAbstract {
19 public:
20 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
21
22 explicit SolverKKT(std::shared_ptr<ShootingProblem> problem);
23 virtual ~SolverKKT();
24
25 virtual bool solve(
26 const std::vector<Eigen::VectorXd>& init_xs = DEFAULT_VECTOR,
27 const std::vector<Eigen::VectorXd>& init_us = DEFAULT_VECTOR,
28 const std::size_t maxiter = 100, const bool is_feasible = false,
29 const double regInit = 1e-9);
30 virtual void computeDirection(const bool recalc = true);
31 virtual double tryStep(const double steplength = 1);
32 virtual double stoppingCriteria();
33 virtual const Eigen::Vector2d& expectedImprovement();
34
35 const Eigen::MatrixXd& get_kkt() const;
36 const Eigen::VectorXd& get_kktref() const;
37 const Eigen::VectorXd& get_primaldual() const;
38 const std::vector<Eigen::VectorXd>& get_dxs() const;
39 const std::vector<Eigen::VectorXd>& get_dus() const;
40 const std::vector<Eigen::VectorXd>& get_lambdas() const;
41 std::size_t get_nx() const;
42 std::size_t get_ndx() const;
43 std::size_t get_nu() const;
44
45 protected:
46 double reg_incfactor_;
47 double reg_decfactor_;
48 double reg_min_;
49 double reg_max_;
50 double cost_try_;
51 std::vector<Eigen::VectorXd> xs_try_;
52 std::vector<Eigen::VectorXd> us_try_;
53
54 private:
55 double calcDiff();
56 void computePrimalDual();
57 void increaseRegularization();
58 void decreaseRegularization();
59 void allocateData();
60
61 std::size_t nx_;
62 std::size_t ndx_;
63 std::size_t nu_;
64 std::vector<Eigen::VectorXd> dxs_;
65 std::vector<Eigen::VectorXd> dus_;
66 std::vector<Eigen::VectorXd> lambdas_;
67
68 // allocate data
69 Eigen::MatrixXd kkt_;
70 Eigen::VectorXd kktref_;
71 Eigen::VectorXd primaldual_;
72 Eigen::VectorXd primal_;
73 Eigen::VectorXd dual_;
74 std::vector<double> alphas_;
75 double th_grad_;
76 bool was_feasible_;
77 Eigen::VectorXd kkt_primal_;
78 Eigen::VectorXd dF;
79};
80
81} // namespace crocoddyl
82
83#endif // CROCODDYL_CORE_SOLVERS_KKT_HPP_
Abstract class for optimal control solvers.
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:165
virtual double stoppingCriteria()
Return a positive value that quantifies the algorithm termination.
Definition kkt.cpp:138
virtual void computeDirection(const bool recalc=true)
Compute the search direction for the current guess .
Definition kkt.cpp:89