Crocoddyl
fddp.hpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2019-2022, LAAS-CNRS, University of Edinburgh,
5 // Heriot-Watt University
6 // Copyright note valid unless otherwise stated in individual files.
7 // All rights reserved.
9 
10 #ifndef CROCODDYL_CORE_SOLVERS_FDDP_HPP_
11 #define CROCODDYL_CORE_SOLVERS_FDDP_HPP_
12 
13 #include <Eigen/Cholesky>
14 #include <vector>
15 
16 #include "crocoddyl/core/solvers/ddp.hpp"
17 
18 namespace crocoddyl {
19 
58 class SolverFDDP : public SolverDDP {
59  public:
60  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
61 
67  explicit SolverFDDP(boost::shared_ptr<ShootingProblem> problem);
68  virtual ~SolverFDDP();
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 
90  virtual const Eigen::Vector2d& expectedImprovement();
91 
96  virtual void forwardPass(const double stepLength);
97 
101  double get_th_acceptnegstep() const;
102 
106  void set_th_acceptnegstep(const double th_acceptnegstep);
107 
108  protected:
109  double dg_;
110  double dq_;
111  double dv_;
114 };
115 
116 } // namespace crocoddyl
117 
118 #endif // CROCODDYL_CORE_SOLVERS_FDDP_HPP_
Feasibility-driven Differential Dynamic Programming (FDDP) solver.
Definition: fddp.hpp:58
void set_th_acceptnegstep(const double th_acceptnegstep)
Modify the threshold used for accepting step along ascent direction.
Definition: fddp.cpp:273
double dg_
Internal data for computing the expected improvement.
Definition: fddp.hpp:109
double dv_
Internal data for computing the expected improvement.
Definition: fddp.hpp:111
double th_acceptnegstep_
Definition: fddp.hpp:112
EIGEN_MAKE_ALIGNED_OPERATOR_NEW SolverFDDP(boost::shared_ptr< ShootingProblem > problem)
Initialize the FDDP solver.
Definition: fddp.cpp:19
virtual const Eigen::Vector2d & expectedImprovement()
Return the expected improvement from a given current search direction .
Definition: fddp.cpp:124
double dq_
Internal data for computing the expected improvement.
Definition: fddp.hpp:110
double get_th_acceptnegstep() const
Return the threshold used for accepting step along ascent direction.
Definition: fddp.cpp:271
void updateExpectedImprovement()
Update internal values for computing the expected improvement.
Definition: fddp.cpp:151