Crocoddyl
 
Loading...
Searching...
No Matches
intro.hpp
1
2// BSD 3-Clause License
3//
4// Copyright (C) 2021-2023, Heriot-Watt University, University of Edinburgh
5// Copyright note valid unless otherwise stated in individual files.
6// All rights reserved.
8
9#ifndef CROCODDYL_CORE_SOLVERS_INTRO_HPP_
10#define CROCODDYL_CORE_SOLVERS_INTRO_HPP_
11
12#include "crocoddyl/core/solvers/fddp.hpp"
13
14namespace crocoddyl {
15
16enum EqualitySolverType { LuNull = 0, QrNull, Schur };
17
18class SolverIntro : public SolverFDDP {
19 public:
20 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
21
29 explicit SolverIntro(std::shared_ptr<ShootingProblem> problem);
30 virtual ~SolverIntro();
31
32 virtual bool solve(
33 const std::vector<Eigen::VectorXd>& init_xs = DEFAULT_VECTOR,
34 const std::vector<Eigen::VectorXd>& init_us = DEFAULT_VECTOR,
35 const std::size_t maxiter = 100, const bool is_feasible = false,
36 const double init_reg = NAN);
37 virtual double tryStep(const double step_length = 1);
38 virtual double stoppingCriteria();
39 virtual void resizeData();
40 virtual double calcDiff();
41 virtual void computeValueFunction(
42 const std::size_t t, const std::shared_ptr<ActionModelAbstract>& model);
43 virtual void computeGains(const std::size_t t);
44
48 EqualitySolverType get_equality_solver() const;
49
53 double get_th_feas() const;
54
58 double get_rho() const;
59
64 double get_upsilon() const;
65
69 const std::vector<std::size_t>& get_Hu_rank() const;
70
75 const std::vector<Eigen::MatrixXd>& get_YZ() const;
76
80 const std::vector<Eigen::MatrixXd>& get_Qzz() const;
81
85 const std::vector<Eigen::MatrixXd>& get_Qxz() const;
86
90 const std::vector<Eigen::MatrixXd>& get_Quz() const;
91
95 const std::vector<Eigen::VectorXd>& get_Qz() const;
96
101 const std::vector<Eigen::MatrixXd>& get_Hy() const;
102
107 const std::vector<Eigen::VectorXd>& get_kz() const;
108
112 const std::vector<Eigen::MatrixXd>& get_Kz() const;
113
117 const std::vector<Eigen::VectorXd>& get_ks() const;
118
122 const std::vector<Eigen::MatrixXd>& get_Ks() const;
123
130 bool get_zero_upsilon() const;
131
139 void set_equality_solver(const EqualitySolverType type);
140
144 void set_th_feas(const double th_feas);
145
149 void set_rho(const double rho);
150
157 void set_zero_upsilon(const bool zero_upsilon);
158
159 protected:
160 enum EqualitySolverType
162 double th_feas_;
163 double rho_;
165 double
170
171 std::vector<std::size_t>
173 std::vector<Eigen::MatrixXd> KQuu_tmp_;
174 std::vector<Eigen::MatrixXd>
178 std::vector<Eigen::MatrixXd>
181 std::vector<Eigen::VectorXd>
183 std::vector<Eigen::MatrixXd>
185 std::vector<Eigen::MatrixXd>
187 std::vector<Eigen::MatrixXd>
189 std::vector<Eigen::VectorXd>
191 std::vector<Eigen::MatrixXd>
193 std::vector<Eigen::VectorXd>
195 std::vector<Eigen::MatrixXd>
197 std::vector<Eigen::MatrixXd> QuuinvHuT_;
198 std::vector<Eigen::LLT<Eigen::MatrixXd> > Qzz_llt_;
199 std::vector<Eigen::FullPivLU<Eigen::MatrixXd> >
202 std::vector<Eigen::ColPivHouseholderQR<Eigen::MatrixXd> >
205 std::vector<Eigen::PartialPivLU<Eigen::MatrixXd> >
208};
209
210} // namespace crocoddyl
211
212#endif // CROCODDYL_CORE_SOLVERS_INTRO_HPP_
Feasibility-driven Differential Dynamic Programming (FDDP) solver.
Definition fddp.hpp:55
std::vector< Eigen::MatrixXd > Kz_
Feedback gain in the nullspace of .
Definition intro.hpp:192
std::vector< Eigen::FullPivLU< Eigen::MatrixXd > > Hu_lu_
Definition intro.hpp:200
std::vector< Eigen::MatrixXd > Hy_
Definition intro.hpp:179
double th_feas_
Threshold for switching to feasibility.
Definition intro.hpp:162
std::vector< std::size_t > Hu_rank_
Rank of the control Jacobian of the equality constraints.
Definition intro.hpp:172
std::vector< Eigen::MatrixXd > Qxz_
Hessian of the reduced Hamiltonian .
Definition intro.hpp:186
std::vector< Eigen::VectorXd > Qz_
Jacobian of the reduced Hamiltonian .
Definition intro.hpp:182
double get_th_feas() const
Return the threshold for switching to feasibility.
Definition intro.cpp:405
const std::vector< Eigen::VectorXd > & get_kz() const
Return feedforward term related to the nullspace of .
Definition intro.cpp:435
EqualitySolverType get_equality_solver() const
Return the type of solver used for handling the equality constraints.
Definition intro.cpp:401
const std::vector< Eigen::MatrixXd > & get_YZ() const
Return the rank of control-equality constraints \mathbf{H_u}\f.
Definition intro.cpp:417
const std::vector< Eigen::MatrixXd > & get_Kz() const
Return feedback gain related to the nullspace of .
Definition intro.cpp:437
void set_equality_solver(const EqualitySolverType type)
Modify the type of solver used for handling the equality constraints.
Definition intro.cpp:443
std::vector< Eigen::MatrixXd > Quz_
Hessian of the reduced Hamiltonian .
Definition intro.hpp:188
std::vector< Eigen::MatrixXd > Ks_
Feedback gain related to the equality constraints.
Definition intro.hpp:196
std::vector< Eigen::LLT< Eigen::MatrixXd > > Qzz_llt_
Cholesky LLT solver.
Definition intro.hpp:198
std::vector< Eigen::ColPivHouseholderQR< Eigen::MatrixXd > > Hu_qr_
Definition intro.hpp:203
std::vector< Eigen::MatrixXd > Qzz_
Hessian of the reduced Hamiltonian .
Definition intro.hpp:184
void set_zero_upsilon(const bool zero_upsilon)
Modify the zero-upsilon label.
Definition intro.cpp:456
const std::vector< Eigen::MatrixXd > & get_Hy() const
Return span-projected Jacobian of the equality-constraint with respect to the control.
Definition intro.cpp:433
double get_rho() const
Return the rho parameter used in the merit function.
Definition intro.cpp:407
const std::vector< Eigen::MatrixXd > & get_Qxz() const
Return Hessian of the reduced Hamiltonian .
Definition intro.cpp:423
const std::vector< Eigen::MatrixXd > & get_Ks() const
Return feedback gain related to the equality constraints.
Definition intro.cpp:441
void set_th_feas(const double th_feas)
Modify the threshold for switching to feasibility.
Definition intro.cpp:447
bool get_zero_upsilon() const
Return the zero-upsilon label.
Definition intro.cpp:411
double get_upsilon() const
Return the estimated penalty parameter that balances relative contribution of the cost function and e...
Definition intro.cpp:409
std::vector< Eigen::MatrixXd > YZ_
Definition intro.hpp:175
const std::vector< Eigen::MatrixXd > & get_Quz() const
Return Hessian of the reduced Hamiltonian .
Definition intro.cpp:427
enum EqualitySolverType eq_solver_
Strategy used for handling the equality constraints.
Definition intro.hpp:160
const std::vector< Eigen::MatrixXd > & get_Qzz() const
Return Hessian of the reduced Hamiltonian .
Definition intro.cpp:419
std::vector< Eigen::VectorXd > kz_
Feedforward term in the nullspace of .
Definition intro.hpp:190
const std::vector< Eigen::VectorXd > & get_ks() const
Return feedforward term related to the equality constraints.
Definition intro.cpp:439
const std::vector< Eigen::VectorXd > & get_Qz() const
Return Jacobian of the reduced Hamiltonian .
Definition intro.cpp:431
std::vector< Eigen::PartialPivLU< Eigen::MatrixXd > > Hy_lu_
Definition intro.hpp:206
std::vector< Eigen::VectorXd > ks_
Feedforward term related to the equality constraints.
Definition intro.hpp:194
void set_rho(const double rho)
Modify the rho parameter used in the merit function.
Definition intro.cpp:449