Crocoddyl
 
Loading...
Searching...
No Matches
ddp.hpp
1
2// 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 <Eigen/Cholesky>
14#include <vector>
15
16#include "crocoddyl/core/mathbase.hpp"
17#include "crocoddyl/core/solver-base.hpp"
18#include "crocoddyl/core/utils/deprecate.hpp"
19
20namespace crocoddyl {
21
56class SolverDDP : public SolverAbstract {
57 public:
58 EIGEN_MAKE_ALIGNED_OPERATOR_NEW
59
60 typedef typename MathBaseTpl<double>::MatrixXsRowMajor MatrixXdRowMajor;
61
67 explicit SolverDDP(std::shared_ptr<ShootingProblem> problem);
68 virtual ~SolverDDP();
69
70 virtual bool solve(
71 const std::vector<Eigen::VectorXd>& init_xs = DEFAULT_VECTOR,
72 const std::vector<Eigen::VectorXd>& init_us = DEFAULT_VECTOR,
73 const std::size_t maxiter = 100, const bool is_feasible = false,
74 const double init_reg = NAN);
75 virtual void computeDirection(const bool recalc = true);
76 virtual double tryStep(const double steplength = 1);
77 virtual double stoppingCriteria();
78 virtual const Eigen::Vector2d& expectedImprovement();
79 virtual void resizeData();
80
90 virtual double calcDiff();
91
124 virtual void backwardPass();
125
139 virtual void forwardPass(const double stepLength);
140
149 virtual void computeActionValueFunction(
150 const std::size_t t, const std::shared_ptr<ActionModelAbstract>& model,
151 const std::shared_ptr<ActionDataAbstract>& data);
152
162 virtual void computeValueFunction(
163 const std::size_t t, const std::shared_ptr<ActionModelAbstract>& model);
164
179 virtual void computeGains(const std::size_t t);
180
185 void increaseRegularization();
186
191 void decreaseRegularization();
192
196 virtual void allocateData();
197
201 double get_reg_incfactor() const;
202
206 double get_reg_decfactor() const;
207
211 DEPRECATED("Use get_reg_incfactor() or get_reg_decfactor()",
212 double get_regfactor() const;)
213
217 double get_reg_min() const;
218 DEPRECATED("Use get_reg_min()", double get_regmin() const);
219
223 double get_reg_max() const;
224 DEPRECATED("Use get_reg_max()", double get_regmax() const);
225
229 const std::vector<double>& get_alphas() const;
230
234 double get_th_stepdec() const;
235
239 double get_th_stepinc() const;
240
245 double get_th_grad() const;
246
250 const std::vector<Eigen::MatrixXd>& get_Vxx() const;
251
255 const std::vector<Eigen::VectorXd>& get_Vx() const;
256
261 const std::vector<Eigen::MatrixXd>& get_Qxx() const;
262
267 const std::vector<Eigen::MatrixXd>& get_Qxu() const;
268
273 const std::vector<Eigen::MatrixXd>& get_Quu() const;
274
279 const std::vector<Eigen::VectorXd>& get_Qx() const;
280
285 const std::vector<Eigen::VectorXd>& get_Qu() const;
286
290 const std::vector<MatrixXdRowMajor>& get_K() const;
291
295 const std::vector<Eigen::VectorXd>& get_k() const;
296
300 void set_reg_incfactor(const double reg_factor);
301
305 void set_reg_decfactor(const double reg_factor);
306
310 DEPRECATED("Use set_reg_incfactor() or set_reg_decfactor()",
311 void set_regfactor(const double reg_factor);)
312
316 void set_reg_min(const double regmin);
317 DEPRECATED("Use set_reg_min()", void set_regmin(const double regmin));
318
322 void set_reg_max(const double regmax);
323 DEPRECATED("Use set_reg_max()", void set_regmax(const double regmax));
324
328 void set_alphas(const std::vector<double>& alphas);
329
333 void set_th_stepdec(const double th_step);
334
338 void set_th_stepinc(const double th_step);
339
344 void set_th_grad(const double th_grad);
345
346 protected:
347 double reg_incfactor_;
349 double reg_decfactor_;
351 double reg_min_;
352 double reg_max_;
353
354 double cost_try_;
355 std::vector<Eigen::VectorXd>
356 xs_try_;
357 std::vector<Eigen::VectorXd>
358 us_try_;
359 std::vector<Eigen::VectorXd>
360 dx_;
361
362 // allocate data
363 std::vector<Eigen::MatrixXd>
364 Vxx_;
365 Eigen::MatrixXd
366 Vxx_tmp_;
367 std::vector<Eigen::VectorXd>
368 Vx_;
369 std::vector<Eigen::MatrixXd>
370 Qxx_;
371 std::vector<Eigen::MatrixXd>
372 Qxu_;
373 std::vector<Eigen::MatrixXd>
374 Quu_;
375 std::vector<Eigen::VectorXd>
376 Qx_;
377 std::vector<Eigen::VectorXd>
378 Qu_;
379 std::vector<MatrixXdRowMajor> K_;
380 std::vector<Eigen::VectorXd> k_;
381
382 Eigen::VectorXd xnext_;
383 MatrixXdRowMajor FxTVxx_p_;
385 std::vector<MatrixXdRowMajor>
386 FuTVxx_p_;
389 Eigen::VectorXd fTVxx_p_;
391 std::vector<Eigen::LLT<Eigen::MatrixXd> > Quu_llt_;
392 std::vector<Eigen::VectorXd>
393 Quuk_;
395 std::vector<double>
396 alphas_;
397 double th_grad_;
399 double
400 th_stepdec_;
401 double
402 th_stepinc_;
403};
404
405} // namespace crocoddyl
406
407#endif // CROCODDYL_CORE_SOLVERS_DDP_HPP_