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 "crocoddyl/core/solvers/ddp.hpp"
14 
15 namespace crocoddyl {
16 
55 class SolverFDDP : public SolverDDP {
56  public:
57  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
58 
64  explicit SolverFDDP(std::shared_ptr<ShootingProblem> problem);
65  virtual ~SolverFDDP();
66 
67  virtual bool solve(
68  const std::vector<Eigen::VectorXd>& init_xs = DEFAULT_VECTOR,
69  const std::vector<Eigen::VectorXd>& init_us = DEFAULT_VECTOR,
70  const std::size_t maxiter = 100, const bool is_feasible = false,
71  const double init_reg = NAN);
72 
87  virtual const Eigen::Vector2d& expectedImprovement();
88 
93  virtual void forwardPass(const double stepLength);
94 
98  double get_th_acceptnegstep() const;
99 
103  void set_th_acceptnegstep(const double th_acceptnegstep);
104 
105  protected:
106  double dg_;
107  double dq_;
108  double dv_;
111 };
112 
113 } // namespace crocoddyl
114 
115 #endif // CROCODDYL_CORE_SOLVERS_FDDP_HPP_
Feasibility-driven Differential Dynamic Programming (FDDP) solver.
Definition: fddp.hpp:55
void set_th_acceptnegstep(const double th_acceptnegstep)
Modify the threshold used for accepting step along ascent direction.
Definition: fddp.cpp:270
double dg_
Internal data for computing the expected improvement.
Definition: fddp.hpp:106
double dv_
Internal data for computing the expected improvement.
Definition: fddp.hpp:108
double th_acceptnegstep_
Definition: fddp.hpp:109
virtual const Eigen::Vector2d & expectedImprovement()
Return the expected improvement from a given current search direction .
Definition: fddp.cpp:123
double dq_
Internal data for computing the expected improvement.
Definition: fddp.hpp:107
double get_th_acceptnegstep() const
Return the threshold used for accepting step along ascent direction.
Definition: fddp.cpp:268
EIGEN_MAKE_ALIGNED_OPERATOR_NEW SolverFDDP(std::shared_ptr< ShootingProblem > problem)
Initialize the FDDP solver.
Definition: fddp.cpp:18
void updateExpectedImprovement()
Update internal values for computing the expected improvement.
Definition: fddp.cpp:150