crocoddyl  1.9.0
Contact RObot COntrol by Differential DYnamic programming Library (Crocoddyl)
fddp.hpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2019-2021, LAAS-CNRS, University of Edinburgh
5 // Copyright note valid unless otherwise stated in individual files.
6 // All rights reserved.
8 
9 #ifndef CROCODDYL_CORE_SOLVERS_FDDP_HPP_
10 #define CROCODDYL_CORE_SOLVERS_FDDP_HPP_
11 
12 #include <Eigen/Cholesky>
13 #include <vector>
14 
15 #include "crocoddyl/core/solvers/ddp.hpp"
16 
17 namespace crocoddyl {
18 
50 class SolverFDDP : public SolverDDP {
51  public:
52  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
53 
59  explicit SolverFDDP(boost::shared_ptr<ShootingProblem> problem);
60  virtual ~SolverFDDP();
61 
62  virtual bool solve(const std::vector<Eigen::VectorXd>& init_xs = DEFAULT_VECTOR,
63  const std::vector<Eigen::VectorXd>& init_us = DEFAULT_VECTOR, const std::size_t maxiter = 100,
64  const bool is_feasible = false, const double regInit = 1e-9);
65 
78  virtual const Eigen::Vector2d& expectedImprovement();
79 
84  virtual void forwardPass(const double stepLength);
85 
89  double get_th_acceptnegstep() const;
90 
94  void set_th_acceptnegstep(const double th_acceptnegstep);
95 
96  protected:
97  double dg_;
98  double dq_;
99  double dv_;
100 
101  private:
102  double th_acceptnegstep_;
103 };
104 
105 } // namespace crocoddyl
106 
107 #endif // CROCODDYL_CORE_SOLVERS_FDDP_HPP_
virtual void forwardPass(const double stepLength)
Run the forward pass or rollout.
Definition: fddp.cpp:166
void updateExpectedImprovement()
Update internal values for computing the expected improvement.
Definition: fddp.cpp:142
virtual bool solve(const std::vector< Eigen::VectorXd > &init_xs=DEFAULT_VECTOR, const std::vector< Eigen::VectorXd > &init_us=DEFAULT_VECTOR, const std::size_t maxiter=100, const bool is_feasible=false, const double regInit=1e-9)
Compute the optimal trajectory as lists of and terms.
Definition: fddp.cpp:23
double dq_
Internal data for computing the expected improvement.
Definition: fddp.hpp:98
void set_th_acceptnegstep(const double th_acceptnegstep)
Modify the threshold used for accepting step along ascent direction.
Definition: fddp.cpp:256
EIGEN_MAKE_ALIGNED_OPERATOR_NEW SolverFDDP(boost::shared_ptr< ShootingProblem > problem)
Initialize the FDDP solver.
Definition: fddp.cpp:18
Feasibility-driven Differential Dynamic Programming (FDDP) solver.
Definition: fddp.hpp:50
double dg_
Internal data for computing the expected improvement.
Definition: fddp.hpp:97
double get_th_acceptnegstep() const
Return the threshold used for accepting step along ascent direction.
Definition: fddp.cpp:254
double dv_
Internal data for computing the expected improvement.
Definition: fddp.hpp:99
Differential Dynamic Programming (DDP) solver.
Definition: ddp.hpp:49
virtual const Eigen::Vector2d & expectedImprovement()
Definition: fddp.cpp:118