ilqg_solver.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2019, University of Edinburgh
3 // All rights reserved.
4 //
5 // Redistribution and use in source and binary forms, with or without
6 // modification, are permitted provided that the following conditions are met:
7 //
8 // * Redistributions of source code must retain the above copyright notice,
9 // this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright
11 // notice, this list of conditions and the following disclaimer in the
12 // documentation and/or other materials provided with the distribution.
13 // * Neither the name of nor the names of its contributors may be used to
14 // endorse or promote products derived from this software without specific
15 // prior written permission.
16 //
17 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20 // ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21 // LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22 // CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23 // SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24 // INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25 // CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26 // ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27 // POSSIBILITY OF SUCH DAMAGE.
28 //
29 
31 
33 
34 namespace exotica
35 {
36 void ILQGSolver::SpecifyProblem(PlanningProblemPtr pointer)
37 {
38  if (pointer->type() != "exotica::DynamicTimeIndexedShootingProblem")
39  {
40  ThrowNamed("This ILQGSolver can't solve problem of type '" << pointer->type() << "'!");
41  }
42 
43  MotionSolver::SpecifyProblem(pointer);
44  prob_ = std::static_pointer_cast<DynamicTimeIndexedShootingProblem>(pointer);
45  dynamics_solver_ = prob_->GetScene()->GetDynamicsSolver();
46  if (debug_) HIGHLIGHT_NAMED("ILQGSolver", "initialized");
47 }
48 
49 void ILQGSolver::BackwardPass()
50 {
51  constexpr double min_clamp_ = -1e10;
52  constexpr double max_clamp_ = 1e10;
53  const int T = prob_->get_T();
54  const double dt = dynamics_solver_->get_dt();
55  const int NU = prob_->GetScene()->get_num_controls();
56 
57  // Noise terms
58  Eigen::VectorXd big_C_times_little_c = Eigen::VectorXd::Zero(NU, 1);
59  Eigen::MatrixXd big_C_times_big_C = Eigen::MatrixXd::Zero(NU, NU);
60  Eigen::MatrixXd little_c_times_little_c = Eigen::MatrixXd::Zero(1, 1);
61 
62  // Value function and derivatives at the final timestep
63  double s0 = prob_->GetStateCost(T - 1);
64  Eigen::MatrixXd s = prob_->GetStateCostJacobian(T - 1);
65  Eigen::MatrixXd S = prob_->GetStateCostHessian(T - 1);
66 
67  for (int t = T - 2; t > 0; --t)
68  {
69  // eq. 3
70  Eigen::VectorXd x = prob_->get_X(t), u = prob_->get_U(t);
71  dynamics_solver_->ComputeDerivatives(x, u);
72  Eigen::MatrixXd A = dynamics_solver_->get_Fx();
73  Eigen::MatrixXd B = dynamics_solver_->get_Fu();
74 
75  double q0 = dt * (prob_->GetStateCost(t) + prob_->GetControlCost(t));
76  // Aliases from the paper used. These are used with different names in e.g. DDPSolver.
77  Eigen::MatrixXd q = dt * prob_->GetStateCostJacobian(t);
78  Eigen::MatrixXd Q = dt * prob_->GetStateCostHessian(t);
79  Eigen::MatrixXd r = dt * prob_->GetControlCostJacobian(t);
80  Eigen::MatrixXd R = dt * prob_->GetControlCostHessian(t);
81  Eigen::MatrixXd P = dt * prob_->GetStateControlCostHessian();
82 
83  Eigen::MatrixXd g = r + B.transpose() * s;
84  Eigen::MatrixXd G = P + B.transpose() * S * A;
85  Eigen::MatrixXd H = R + B.transpose() * S * B;
86 
87  if (parameters_.IncludeNoiseTerms)
88  {
89  Eigen::MatrixXd F = prob_->get_F(t);
90  for (int i = 0; i < NU; ++i)
91  {
92  Eigen::MatrixXd C = std::sqrt(dt) * prob_->GetControlNoiseJacobian(i);
93  Eigen::VectorXd c = std::sqrt(dt) * F.col(i);
94 
95  big_C_times_little_c = big_C_times_little_c + C.transpose() * S * c;
96  big_C_times_big_C = big_C_times_big_C + C.transpose() * S * C;
97  little_c_times_little_c = little_c_times_little_c + c.transpose() * S * c;
98  }
99 
100  g = g + big_C_times_little_c;
101  H = H + big_C_times_big_C;
102  }
103 
104  // optimal U
105  Eigen::EigenSolver<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>> eig_solver(H);
106  auto d = eig_solver.eigenvalues();
107  auto V = eig_solver.eigenvectors();
108  Eigen::MatrixXcd D = Eigen::MatrixXcd::Zero(d.size(), d.size());
109 
110  for (int i = 0; i < d.size(); ++i)
111  {
112  if (d[i].real() < 0)
113  d[i] = 0;
114  d[i] = 1. / (d[i] + parameters_.RegularizationRate);
115  D(i, i) = d[i];
116  }
117 
118  Eigen::MatrixXd H1 = (V * D * V.transpose()).real();
119 
120  l_gains_[t] = -H1 * g;
121  L_gains_[t] = -H1 * G;
122 
123  // Recursive terms update
124  S = Q + A.transpose() * S * A + L_gains_[t].transpose() * H * L_gains_[t] + L_gains_[t].transpose() * G + G.transpose() * L_gains_[t];
125  s = q + A.transpose() * s + L_gains_[t].transpose() * H * l_gains_[t] +
126  L_gains_[t].transpose() * g + G.transpose() * l_gains_[t];
127  s0 = q0 + s0 + (l_gains_[t].transpose() * H * l_gains_[t] / 2.0 + l_gains_[t].transpose() * g)(0);
128 
129  if (parameters_.IncludeNoiseTerms)
130  {
131  s0 = s0 + 0.5 * little_c_times_little_c(0);
132  }
133 
134  // fix for large values
135  S = S.unaryExpr([min_clamp_, max_clamp_](double x) -> double {
136  return std::min(std::max(x, min_clamp_), max_clamp_);
137  });
138  s = s.unaryExpr([min_clamp_, max_clamp_](double x) -> double {
139  return std::min(std::max(x, min_clamp_), max_clamp_);
140  });
141  s0 = std::min(std::max(s0, min_clamp_), max_clamp_);
142  }
143 }
144 
145 double ILQGSolver::ForwardPass(const double alpha, Eigen::MatrixXdRefConst ref_x, Eigen::MatrixXdRefConst ref_u)
146 {
147  double cost = 0;
148  const int T = prob_->get_T();
149  const Eigen::MatrixXd control_limits = dynamics_solver_->get_control_limits();
150  const double dt = dynamics_solver_->get_dt();
151 
152  // NOTE: Todorov uses the linearized system in forward simulation.
153  // We here use the full system dynamics. YMMV
154  for (int t = 0; t < T - 1; ++t)
155  {
156  Eigen::VectorXd u = ref_u.col(t);
157  // eq. 12
158  Eigen::VectorXd delta_uk = l_gains_[t] +
159  L_gains_[t] * dynamics_solver_->StateDelta(prob_->get_X(t), ref_x.col(t));
160 
161  u.noalias() += alpha * delta_uk;
162  // clamp controls
163  u = u.cwiseMax(control_limits.col(0)).cwiseMin(control_limits.col(1));
164 
165  prob_->Update(u, t);
166  cost += dt * (prob_->GetControlCost(t) + prob_->GetStateCost(t));
167  }
168 
169  // add terminal cost
170  cost += prob_->GetStateCost(T - 1);
171  return cost;
172 }
173 
174 void ILQGSolver::Solve(Eigen::MatrixXd& solution)
175 {
176  if (!prob_) ThrowNamed("Solver has not been initialized!");
177  Timer planning_timer, backward_pass_timer, line_search_timer;
178  // TODO: This is an interesting approach but might give us incorrect results.
179  prob_->DisableStochasticUpdates();
180 
181  const int T = prob_->get_T();
182  const int NU = prob_->GetScene()->get_num_controls();
183  const int NX = prob_->GetScene()->get_num_state();
184  const double dt = dynamics_solver_->get_dt();
185  prob_->ResetCostEvolution(GetNumberOfMaxIterations() + 1);
186  prob_->PreUpdate();
187 
188  double initial_cost = 0;
189  for (int t = 0; t < T - 1; ++t)
190  initial_cost += dt * (prob_->GetControlCost(t) + prob_->GetStateCost(t));
191 
192  // add terminal cost
193  initial_cost += prob_->GetStateCost(T - 1);
194  prob_->SetCostEvolution(0, initial_cost);
195 
196  // initialize Gain matrices
197  l_gains_.assign(T, Eigen::MatrixXd::Zero(NU, 1));
198  L_gains_.assign(T, Eigen::MatrixXd::Zero(NU, NX));
199 
200  // all of the below are not pointers, since we want to copy over
201  // solutions across iterations
202  Eigen::MatrixXd new_U, global_best_U = prob_->get_U();
203  solution.resize(T - 1, NU);
204 
205  if (debug_) HIGHLIGHT_NAMED("ILQGSolver", "Running ILQG solver for max " << parameters_.MaxIterations << " iterations");
206 
207  double last_cost = initial_cost, global_best_cost = initial_cost;
208  int last_best_iteration = 0;
209 
210  for (int iteration = 1; iteration <= GetNumberOfMaxIterations(); ++iteration)
211  {
212  // Check whether user interrupted (Ctrl+C)
213  if (Server::IsRos() && !ros::ok())
214  {
215  if (debug_) HIGHLIGHT("Solving cancelled by user");
216  prob_->termination_criterion = TerminationCriterion::UserDefined;
217  break;
218  }
219 
220  // Backwards pass computes the gains
221  backward_pass_timer.Reset();
222  BackwardPass();
223  if (debug_) HIGHLIGHT_NAMED("ILQGSolver", "Backward pass complete in " << backward_pass_timer.GetDuration());
224  // if (debug_) HIGHLIGHT_NAMED("ILQGSolver", "Backward pass complete in " << backward_pass_timer.GetDuration());
225 
226  line_search_timer.Reset();
227  // forward pass to compute new control trajectory
228  // TODO: Configure line-search space from xml
229  // TODO (Wolf): What you are doing is a forward line-search when we may try a
230  // backtracking line search for a performance improvement later - this just as an aside.
231  const Eigen::VectorXd alpha_space = Eigen::VectorXd::LinSpaced(10, 0.1, 1.0);
232  double current_cost = 0, best_alpha = 0;
233 
234  Eigen::MatrixXd ref_x = prob_->get_X(),
235  ref_u = prob_->get_U();
236  // perform a linear search to find the best rate
237  for (int ai = 0; ai < alpha_space.rows(); ++ai)
238  {
239  double alpha = alpha_space(ai);
240  double cost = ForwardPass(alpha, ref_x, ref_u);
241 
242  if (ai == 0 || (cost < current_cost && !std::isnan(cost)))
243  {
244  current_cost = cost;
245  new_U = prob_->get_U();
246  best_alpha = alpha;
247  }
248  // else if (cost > current_cost)
249  // break;
250  }
251 
252  // finite checks
253  if (!new_U.allFinite() || !std::isfinite(current_cost))
254  {
255  prob_->termination_criterion = TerminationCriterion::Divergence;
256  WARNING_NAMED("ILQGSolver", "Diverged!");
257  return;
258  }
259 
260  if (debug_)
261  {
262  HIGHLIGHT_NAMED("ILQGSolver", "Forward pass complete in " << line_search_timer.GetDuration() << " with cost: " << current_cost << " and alpha " << best_alpha);
263  HIGHLIGHT_NAMED("ILQGSolver", "Final state: " << prob_->get_X(T - 1).transpose());
264  }
265 
266  // copy solutions for next iteration
267  if (global_best_cost > current_cost)
268  {
269  global_best_cost = current_cost;
270  last_best_iteration = iteration;
271  global_best_U = new_U;
272  best_ref_x_ = ref_x;
273  best_ref_u_ = ref_u;
274  }
275 
276  if (iteration - last_best_iteration > parameters_.FunctionTolerancePatience)
277  {
278  if (debug_) HIGHLIGHT_NAMED("ILQGSolver", "Early stopping criterion reached. Time: " << planning_timer.GetDuration());
279  prob_->termination_criterion = TerminationCriterion::FunctionTolerance;
280  break;
281  }
282 
283  if (last_cost - current_cost < parameters_.FunctionTolerance && last_cost - current_cost > 0)
284  {
285  if (debug_) HIGHLIGHT_NAMED("ILQGSolver", "Function tolerance reached. Time: " << planning_timer.GetDuration());
286  prob_->termination_criterion = TerminationCriterion::Divergence;
287  break;
288  }
289 
290  if (debug_ && iteration == parameters_.MaxIterations)
291  {
292  HIGHLIGHT_NAMED("ILQGSolver", "Max iterations reached. Time: " << planning_timer.GetDuration());
293  prob_->termination_criterion = TerminationCriterion::IterationLimit;
294  }
295 
296  last_cost = current_cost;
297  for (int t = 0; t < T - 1; ++t)
298  prob_->Update(new_U.col(t), t);
299  prob_->SetCostEvolution(iteration, current_cost);
300  }
301 
302  // store the best solution found over all iterations
303  for (int t = 0; t < T - 1; ++t)
304  {
305  solution.row(t) = global_best_U.col(t).transpose();
306  prob_->Update(global_best_U.col(t), t);
307  }
308 
309  planning_time_ = planning_timer.GetDuration();
310 
311  // TODO: See note at disable.
312  prob_->EnableStochasticUpdates();
313 }
314 
315 Eigen::VectorXd ILQGSolver::GetFeedbackControl(Eigen::VectorXdRefConst x, int t) const
316 {
317  const Eigen::MatrixXd control_limits = dynamics_solver_->get_control_limits();
318 
319  Eigen::VectorXd delta_uk = l_gains_[t] +
320  L_gains_[t] * dynamics_solver_->StateDelta(x, best_ref_x_.col(t));
321 
322  Eigen::VectorXd u = best_ref_u_.col(t) + delta_uk;
323  return u.cwiseMax(control_limits.col(0)).cwiseMin(control_limits.col(1));
324 }
325 
326 } // namespace exotica
d
#define HIGHLIGHT(x)
#define WARNING_NAMED(name, x)
XmlRpcServer s
#define REGISTER_MOTIONSOLVER_TYPE(TYPE, DERIV)
const Eigen::Ref< const Eigen::MatrixXd > & MatrixXdRefConst
#define ThrowNamed(m)
ROSCPP_DECL bool ok()
const Eigen::Ref< const Eigen::VectorXd > & VectorXdRefConst
#define HIGHLIGHT_NAMED(name, x)
double GetDuration() const
std::shared_ptr< PlanningProblem > PlanningProblemPtr
ScenePtr GetScene() const
const AutoDiffScalar< DerType > & real(const AutoDiffScalar< DerType > &x)
double x
r


exotica_ilqg_solver
Author(s): Traiko Dinev
autogenerated on Sat Apr 10 2021 02:36:28