Crocoddyl
ddp.hpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2019-2022, LAAS-CNRS, University of Edinburgh,
5 // University of Oxford, Heriot-Watt University
6 // Copyright note valid unless otherwise stated in individual files.
7 // All rights reserved.
9 
10 #ifndef CROCODDYL_CORE_SOLVERS_DDP_HPP_
11 #define CROCODDYL_CORE_SOLVERS_DDP_HPP_
12 
13 #include "crocoddyl/core/solver-base.hpp"
14 #include "crocoddyl/core/utils/deprecate.hpp"
15 
16 namespace crocoddyl {
17 
52 class SolverDDP : public SolverAbstract {
53  public:
54  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
55 
56  typedef typename MathBaseTpl<double>::MatrixXsRowMajor MatrixXdRowMajor;
57 
63  explicit SolverDDP(std::shared_ptr<ShootingProblem> problem);
64  virtual ~SolverDDP();
65 
66  virtual bool solve(
67  const std::vector<Eigen::VectorXd>& init_xs = DEFAULT_VECTOR,
68  const std::vector<Eigen::VectorXd>& init_us = DEFAULT_VECTOR,
69  const std::size_t maxiter = 100, const bool is_feasible = false,
70  const double init_reg = NAN);
71  virtual void computeDirection(const bool recalc = true);
72  virtual double tryStep(const double steplength = 1);
73  virtual double stoppingCriteria();
74  virtual const Eigen::Vector2d& expectedImprovement();
75  virtual void resizeData();
76 
86  virtual double calcDiff();
87 
120  virtual void backwardPass();
121 
135  virtual void forwardPass(const double stepLength);
136 
145  virtual void computeActionValueFunction(
146  const std::size_t t, const std::shared_ptr<ActionModelAbstract>& model,
147  const std::shared_ptr<ActionDataAbstract>& data);
148 
158  virtual void computeValueFunction(
159  const std::size_t t, const std::shared_ptr<ActionModelAbstract>& model);
160 
175  virtual void computeGains(const std::size_t t);
176 
181  void increaseRegularization();
182 
187  void decreaseRegularization();
188 
192  virtual void allocateData();
193 
197  double get_reg_incfactor() const;
198 
202  double get_reg_decfactor() const;
203 
207  DEPRECATED("Use get_reg_incfactor() or get_reg_decfactor()",
208  double get_regfactor() const;)
209 
213  double get_reg_min() const;
214  DEPRECATED("Use get_reg_min()", double get_regmin() const);
215 
219  double get_reg_max() const;
220  DEPRECATED("Use get_reg_max()", double get_regmax() const);
221 
225  const std::vector<double>& get_alphas() const;
226 
230  double get_th_stepdec() const;
231 
235  double get_th_stepinc() const;
236 
241  double get_th_grad() const;
242 
246  const std::vector<Eigen::MatrixXd>& get_Vxx() const;
247 
251  const std::vector<Eigen::VectorXd>& get_Vx() const;
252 
257  const std::vector<Eigen::MatrixXd>& get_Qxx() const;
258 
263  const std::vector<Eigen::MatrixXd>& get_Qxu() const;
264 
269  const std::vector<Eigen::MatrixXd>& get_Quu() const;
270 
275  const std::vector<Eigen::VectorXd>& get_Qx() const;
276 
281  const std::vector<Eigen::VectorXd>& get_Qu() const;
282 
286  const std::vector<MatrixXdRowMajor>& get_K() const;
287 
291  const std::vector<Eigen::VectorXd>& get_k() const;
292 
296  void set_reg_incfactor(const double reg_factor);
297 
301  void set_reg_decfactor(const double reg_factor);
302 
306  DEPRECATED("Use set_reg_incfactor() or set_reg_decfactor()",
307  void set_regfactor(const double reg_factor);)
308 
312  void set_reg_min(const double regmin);
313  DEPRECATED("Use set_reg_min()", void set_regmin(const double regmin));
314 
318  void set_reg_max(const double regmax);
319  DEPRECATED("Use set_reg_max()", void set_regmax(const double regmax));
320 
324  void set_alphas(const std::vector<double>& alphas);
325 
329  void set_th_stepdec(const double th_step);
330 
334  void set_th_stepinc(const double th_step);
335 
340  void set_th_grad(const double th_grad);
341 
342  protected:
343  double reg_incfactor_;
345  double reg_decfactor_;
347  double reg_min_;
348  double reg_max_;
349 
350  double cost_try_;
351  std::vector<Eigen::VectorXd>
352  xs_try_;
353  std::vector<Eigen::VectorXd>
354  us_try_;
355  std::vector<Eigen::VectorXd>
356  dx_;
357 
358  // allocate data
359  std::vector<Eigen::MatrixXd>
360  Vxx_;
361  Eigen::MatrixXd
362  Vxx_tmp_;
363  std::vector<Eigen::VectorXd>
364  Vx_;
365  std::vector<Eigen::MatrixXd>
366  Qxx_;
367  std::vector<Eigen::MatrixXd>
368  Qxu_;
369  std::vector<Eigen::MatrixXd>
370  Quu_;
371  std::vector<Eigen::VectorXd>
372  Qx_;
373  std::vector<Eigen::VectorXd>
374  Qu_;
375  std::vector<MatrixXdRowMajor> K_;
376  std::vector<Eigen::VectorXd> k_;
377 
378  Eigen::VectorXd xnext_;
379  MatrixXdRowMajor FxTVxx_p_;
381  std::vector<MatrixXdRowMajor>
382  FuTVxx_p_;
385  Eigen::VectorXd fTVxx_p_;
387  std::vector<Eigen::LLT<Eigen::MatrixXd> > Quu_llt_;
388  std::vector<Eigen::VectorXd>
389  Quuk_;
391  std::vector<double>
392  alphas_;
393  double th_grad_;
395  double
396  th_stepdec_;
397  double
398  th_stepinc_;
399 };
400 
401 } // namespace crocoddyl
402 
403 #endif // CROCODDYL_CORE_SOLVERS_DDP_HPP_