crocoddyl  1.9.0
Contact RObot COntrol by Differential DYnamic programming Library (Crocoddyl)
box-ddp.cpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2019-2021, CNRS-LAAS, University of Edinburgh
5 // Copyright note valid unless otherwise stated in individual files.
6 // All rights reserved.
8 
9 #include <iostream>
10 
11 #include "crocoddyl/core/solvers/box-ddp.hpp"
12 #include "crocoddyl/core/utils/exception.hpp"
13 
14 namespace crocoddyl {
15 
16 SolverBoxDDP::SolverBoxDDP(boost::shared_ptr<ShootingProblem> problem)
17  : SolverDDP(problem), qp_(problem->get_runningModels()[0]->get_nu(), 100, 0.1, 1e-5, 0.) {
18  allocateData();
19 
20  const std::size_t n_alphas = 10;
21  alphas_.resize(n_alphas);
22  for (std::size_t n = 0; n < n_alphas; ++n) {
23  alphas_[n] = 1. / pow(2., static_cast<double>(n));
24  }
25  // Change the default convergence tolerance since the gradient of the Lagrangian is smaller
26  // than an unconstrained OC problem (i.e. gradient = Qu - mu^T * C where mu > 0 and C defines
27  // the inequality matrix that bounds the control); and we don't have access to mu from the
28  // box QP.
29  th_stop_ = 5e-5;
30 }
31 
32 SolverBoxDDP::~SolverBoxDDP() {}
33 
35  START_PROFILER("SolverBoxDDP::resizeData");
37 
38  const std::size_t T = problem_->get_T();
39  const std::vector<boost::shared_ptr<ActionModelAbstract> >& models = problem_->get_runningModels();
40  for (std::size_t t = 0; t < T; ++t) {
41  const boost::shared_ptr<ActionModelAbstract>& model = models[t];
42  const std::size_t nu = model->get_nu();
43  Quu_inv_[t].conservativeResize(nu, nu);
44  du_lb_[t].conservativeResize(nu);
45  du_ub_[t].conservativeResize(nu);
46  }
47  STOP_PROFILER("SolverBoxDDP::resizeData");
48 }
49 
52 
53  const std::size_t T = problem_->get_T();
54  Quu_inv_.resize(T);
55  du_lb_.resize(T);
56  du_ub_.resize(T);
57  const std::vector<boost::shared_ptr<ActionModelAbstract> >& models = problem_->get_runningModels();
58  for (std::size_t t = 0; t < T; ++t) {
59  const boost::shared_ptr<ActionModelAbstract>& model = models[t];
60  const std::size_t nu = model->get_nu();
61  Quu_inv_[t] = Eigen::MatrixXd::Zero(nu, nu);
62  du_lb_[t] = Eigen::VectorXd::Zero(nu);
63  du_ub_[t] = Eigen::VectorXd::Zero(nu);
64  }
65 }
66 
67 void SolverBoxDDP::computeGains(const std::size_t t) {
68  START_PROFILER("SolverBoxDDP::computeGains");
69  const std::size_t nu = problem_->get_runningModels()[t]->get_nu();
70  if (nu > 0) {
71  if (!problem_->get_runningModels()[t]->get_has_control_limits() || !is_feasible_) {
72  // No control limits on this model: Use vanilla DDP
74  return;
75  }
76 
77  du_lb_[t] = problem_->get_runningModels()[t]->get_u_lb() - us_[t];
78  du_ub_[t] = problem_->get_runningModels()[t]->get_u_ub() - us_[t];
79 
80  START_PROFILER("SolverBoxDDP::boxQP");
81  const BoxQPSolution& boxqp_sol = qp_.solve(Quu_[t], Qu_[t], du_lb_[t], du_ub_[t], k_[t]);
82  START_PROFILER("SolverBoxDDP::boxQP");
83 
84  // Compute controls
85  START_PROFILER("SolverBoxDDP::Quu_invproj");
86  Quu_inv_[t].setZero();
87  for (std::size_t i = 0; i < boxqp_sol.free_idx.size(); ++i) {
88  for (std::size_t j = 0; j < boxqp_sol.free_idx.size(); ++j) {
89  Quu_inv_[t](boxqp_sol.free_idx[i], boxqp_sol.free_idx[j]) = boxqp_sol.Hff_inv(i, j);
90  }
91  }
92  STOP_PROFILER("SolverBoxDDP::Quu_invproj");
93  START_PROFILER("SolverBoxDDP::Quu_invproj_Qxu");
94  K_[t].noalias() = Quu_inv_[t] * Qxu_[t].transpose();
95  STOP_PROFILER("SolverBoxDDP::Quu_invproj_Qxu");
96  k_[t] = -boxqp_sol.x;
97 
98  // The box-QP clamped the gradient direction; this is important for accounting
99  // the algorithm advancement (i.e. stopping criteria)
100  START_PROFILER("SolverBoxDDP::Qu_proj");
101  for (std::size_t i = 0; i < boxqp_sol.clamped_idx.size(); ++i) {
102  Qu_[t](boxqp_sol.clamped_idx[i]) = 0.;
103  }
104  STOP_PROFILER("SolverBoxDDP::Qu_proj");
105  }
106  STOP_PROFILER("SolverBoxDDP::computeGains");
107 }
108 
109 void SolverBoxDDP::forwardPass(double steplength) {
110  if (steplength > 1. || steplength < 0.) {
111  throw_pretty("Invalid argument: "
112  << "invalid step length, value is between 0. to 1.");
113  }
114  START_PROFILER("SolverBoxDDP::forwardPass");
115  cost_try_ = 0.;
116  xnext_ = problem_->get_x0();
117  const std::size_t T = problem_->get_T();
118  const std::vector<boost::shared_ptr<ActionModelAbstract> >& models = problem_->get_runningModels();
119  const std::vector<boost::shared_ptr<ActionDataAbstract> >& datas = problem_->get_runningDatas();
120  for (std::size_t t = 0; t < T; ++t) {
121  const boost::shared_ptr<ActionModelAbstract>& m = models[t];
122  const boost::shared_ptr<ActionDataAbstract>& d = datas[t];
123  const std::size_t nu = m->get_nu();
124 
125  xs_try_[t] = xnext_;
126  m->get_state()->diff(xs_[t], xs_try_[t], dx_[t]);
127  if (nu != 0) {
128  us_try_[t].noalias() = us_[t] - k_[t] * steplength - K_[t] * dx_[t];
129  if (m->get_has_control_limits()) { // clamp control
130  us_try_[t] = us_try_[t].cwiseMax(m->get_u_lb()).cwiseMin(m->get_u_ub());
131  }
132  m->calc(d, xs_try_[t], us_try_[t]);
133  } else {
134  m->calc(d, xs_try_[t]);
135  }
136  xnext_ = d->xnext;
137  cost_try_ += d->cost;
138 
139  if (raiseIfNaN(cost_try_)) {
140  STOP_PROFILER("SolverBoxDDP::forwardPass");
141  throw_pretty("forward_error");
142  }
143  if (raiseIfNaN(xnext_.lpNorm<Eigen::Infinity>())) {
144  STOP_PROFILER("SolverBoxDDP::forwardPass");
145  throw_pretty("forward_error");
146  }
147  }
148 
149  const boost::shared_ptr<ActionModelAbstract>& m = problem_->get_terminalModel();
150  const boost::shared_ptr<ActionDataAbstract>& d = problem_->get_terminalData();
151  if ((is_feasible_) || (steplength == 1)) {
152  xs_try_.back() = xnext_;
153  } else {
154  m->get_state()->integrate(xnext_, fs_.back() * (steplength - 1), xs_try_.back());
155  }
156  m->calc(d, xs_try_.back());
157  cost_try_ += d->cost;
158 
159  if (raiseIfNaN(cost_try_)) {
160  STOP_PROFILER("SolverBoxDDP::forwardPass");
161  throw_pretty("forward_error");
162  }
163 }
164 
165 const std::vector<Eigen::MatrixXd>& SolverBoxDDP::get_Quu_inv() const { return Quu_inv_; }
166 
167 } // namespace crocoddyl
virtual void allocateData()
Allocate all the internal data needed for the solver.
Definition: box-ddp.cpp:50
bool is_feasible_
Label that indicates is the iteration is feasible.
virtual void forwardPass(const double steplength)
Run the forward pass or rollout.
Definition: box-ddp.cpp:109
Box QP solution.
Definition: box-qp.hpp:28
std::vector< MatrixXdRowMajor > K_
Feedback gains .
Definition: ddp.hpp:316
virtual void allocateData()
Allocate all the internal data needed for the solver.
Definition: ddp.cpp:376
std::vector< Eigen::VectorXd > xs_try_
State trajectory computed by line-search procedure.
Definition: ddp.hpp:303
virtual void resizeData()
Resizing the solver data.
Definition: ddp.cpp:172
const BoxQPSolution & solve(const Eigen::MatrixXd &H, const Eigen::VectorXd &q, const Eigen::VectorXd &lb, const Eigen::VectorXd &ub, const Eigen::VectorXd &xinit)
Compute the solution of bound-constrained QP based on Newton projection.
Definition: box-qp.cpp:58
std::vector< Eigen::VectorXd > us_
Control trajectory.
std::vector< Eigen::VectorXd > k_
Feed-forward terms .
Definition: ddp.hpp:317
virtual void computeGains(const std::size_t t)
Compute the feedforward and feedback terms using a Cholesky decomposition.
Definition: box-ddp.cpp:67
Eigen::MatrixXd Hff_inv
Inverse of the free space Hessian.
Definition: box-qp.hpp:48
Eigen::VectorXd xnext_
Next state .
Definition: ddp.hpp:319
boost::shared_ptr< ShootingProblem > problem_
optimal control problem
Eigen::VectorXd x
Decision vector.
Definition: box-qp.hpp:49
std::vector< Eigen::VectorXd > xs_
State trajectory.
std::vector< Eigen::VectorXd > dx_
State error during the roll-out/forward-pass (size T)
Definition: ddp.hpp:305
std::vector< Eigen::VectorXd > fs_
Gaps/defects between shooting nodes.
virtual void computeGains(const std::size_t t)
Compute the feedforward and feedback terms using a Cholesky decomposition.
Definition: ddp.cpp:337
double th_stop_
Tolerance for stopping the algorithm.
std::vector< size_t > clamped_idx
Clamped space indexes.
Definition: box-qp.hpp:51
virtual void resizeData()
Resizing the solver data.
Definition: box-ddp.cpp:34
std::vector< Eigen::MatrixXd > Qxu_
Hessian of the Hamiltonian .
Definition: ddp.hpp:312
std::vector< size_t > free_idx
Free space indexes.
Definition: box-qp.hpp:50
std::vector< Eigen::VectorXd > us_try_
Control trajectory computed by line-search procedure.
Definition: ddp.hpp:304
std::vector< Eigen::MatrixXd > Quu_
Hessian of the Hamiltonian .
Definition: ddp.hpp:313
std::vector< Eigen::VectorXd > Qu_
Gradient of the Hamiltonian .
Definition: ddp.hpp:315
double cost_try_
Total cost computed by line-search procedure.
Definition: ddp.hpp:302