Crocoddyl
 
Loading...
Searching...
No Matches
fddp.cpp
1
2// 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#ifdef CROCODDYL_WITH_MULTITHREADING
11#include <omp.h>
12#endif // CROCODDYL_WITH_MULTITHREADING
13
14#include "crocoddyl/core/solvers/fddp.hpp"
15#include "crocoddyl/core/utils/exception.hpp"
16
17namespace crocoddyl {
18
19SolverFDDP::SolverFDDP(std::shared_ptr<ShootingProblem> problem)
20 : SolverDDP(problem), dg_(0), dq_(0), dv_(0), th_acceptnegstep_(2) {}
21
22SolverFDDP::~SolverFDDP() {}
23
24bool SolverFDDP::solve(const std::vector<Eigen::VectorXd>& init_xs,
25 const std::vector<Eigen::VectorXd>& init_us,
26 const std::size_t maxiter, const bool is_feasible,
27 const double init_reg) {
28 START_PROFILER("SolverFDDP::solve");
29 if (problem_->is_updated()) {
30 resizeData();
31 }
32 xs_try_[0] =
33 problem_->get_x0(); // it is needed in case that init_xs[0] is infeasible
34 setCandidate(init_xs, init_us, is_feasible);
35
36 if (std::isnan(init_reg)) {
37 preg_ = reg_min_;
38 dreg_ = reg_min_;
39 } else {
40 preg_ = init_reg;
41 dreg_ = init_reg;
42 }
43 was_feasible_ = false;
44
45 bool recalcDiff = true;
46 for (iter_ = 0; iter_ < maxiter; ++iter_) {
47 while (true) {
48 try {
49 computeDirection(recalcDiff);
50 } catch (std::exception& e) {
51 recalcDiff = false;
52 increaseRegularization();
53 if (preg_ == reg_max_) {
54 return false;
55 } else {
56 continue;
57 }
58 }
59 break;
60 }
62
63 // We need to recalculate the derivatives when the step length passes
64 recalcDiff = false;
65 for (std::vector<double>::const_iterator it = alphas_.begin();
66 it != alphas_.end(); ++it) {
67 steplength_ = *it;
68
69 try {
70 dV_ = tryStep(steplength_);
71 } catch (std::exception& e) {
72 continue;
73 }
75 dVexp_ = steplength_ * (d_[0] + 0.5 * steplength_ * d_[1]);
76
77 if (dVexp_ >= 0) { // descend direction
78 if (std::abs(d_[0]) < th_grad_ || dV_ > th_acceptstep_ * dVexp_) {
79 was_feasible_ = is_feasible_;
80 setCandidate(xs_try_, us_try_, (was_feasible_) || (steplength_ == 1));
81 cost_ = cost_try_;
82 recalcDiff = true;
83 break;
84 }
85 } else { // reducing the gaps by allowing a small increment in the cost
86 // value
87 if (!is_feasible_ && dV_ > th_acceptnegstep_ * dVexp_) {
88 was_feasible_ = is_feasible_;
89 setCandidate(xs_try_, us_try_, (was_feasible_) || (steplength_ == 1));
90 cost_ = cost_try_;
91 recalcDiff = true;
92 break;
93 }
94 }
95 }
96
97 if (steplength_ > th_stepdec_) {
98 decreaseRegularization();
99 }
100 if (steplength_ <= th_stepinc_) {
101 increaseRegularization();
102 if (preg_ == reg_max_) {
103 STOP_PROFILER("SolverFDDP::solve");
104 return false;
105 }
106 }
107 stoppingCriteria();
108
109 const std::size_t n_callbacks = callbacks_.size();
110 for (std::size_t c = 0; c < n_callbacks; ++c) {
111 CallbackAbstract& callback = *callbacks_[c];
112 callback(*this);
113 }
114
115 if (was_feasible_ && stop_ < th_stop_) {
116 STOP_PROFILER("SolverFDDP::solve");
117 return true;
118 }
119 }
120 STOP_PROFILER("SolverFDDP::solve");
121 return false;
122}
123
124const Eigen::Vector2d& SolverFDDP::expectedImprovement() {
125 dv_ = 0;
126 const std::size_t T = this->problem_->get_T();
127 if (!is_feasible_) {
128 // NB: The dimension of vectors xs_try_ and xs_ are T+1, whereas the
129 // dimension of dx_ is T. Here, we are re-using the final element of dx_ for
130 // the computation of the difference at the terminal node. Using the access
131 // iterator back() this re-use of the final element is fine. Cf. the
132 // discussion at https://github.com/loco-3d/crocoddyl/issues/1022
133 problem_->get_terminalModel()->get_state()->diff(xs_try_.back(), xs_.back(),
134 dx_.back());
135 fTVxx_p_.noalias() = Vxx_.back() * dx_.back();
136 dv_ -= fs_.back().dot(fTVxx_p_);
137 const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
138 problem_->get_runningModels();
139
140 for (std::size_t t = 0; t < T; ++t) {
141 models[t]->get_state()->diff(xs_try_[t], xs_[t], dx_[t]);
142 fTVxx_p_.noalias() = Vxx_[t] * dx_[t];
143 dv_ -= fs_[t].dot(fTVxx_p_);
144 }
145 }
146 d_[0] = dg_ + dv_;
147 d_[1] = dq_ - 2 * dv_;
148 return d_;
149}
150
152 dg_ = 0;
153 dq_ = 0;
154 const std::size_t T = this->problem_->get_T();
155 if (!is_feasible_) {
156 dg_ -= Vx_.back().dot(fs_.back());
157 fTVxx_p_.noalias() = Vxx_.back() * fs_.back();
158 dq_ += fs_.back().dot(fTVxx_p_);
159 }
160 const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
161 problem_->get_runningModels();
162 for (std::size_t t = 0; t < T; ++t) {
163 const std::size_t nu = models[t]->get_nu();
164 if (nu != 0) {
165 dg_ += Qu_[t].dot(k_[t]);
166 dq_ -= k_[t].dot(Quuk_[t]);
167 }
168 if (!is_feasible_) {
169 dg_ -= Vx_[t].dot(fs_[t]);
170 fTVxx_p_.noalias() = Vxx_[t] * fs_[t];
171 dq_ += fs_[t].dot(fTVxx_p_);
172 }
173 }
174}
175
176void SolverFDDP::forwardPass(const double steplength) {
177 if (steplength > 1. || steplength < 0.) {
178 throw_pretty("Invalid argument: "
179 << "invalid step length, value is between 0. to 1.");
180 }
181 START_PROFILER("SolverFDDP::forwardPass");
182 cost_try_ = 0.;
183 xnext_ = problem_->get_x0();
184 const std::size_t T = problem_->get_T();
185 const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
186 problem_->get_runningModels();
187 const std::vector<std::shared_ptr<ActionDataAbstract> >& datas =
188 problem_->get_runningDatas();
189 if ((is_feasible_) || (steplength == 1)) {
190 for (std::size_t t = 0; t < T; ++t) {
191 const std::shared_ptr<ActionModelAbstract>& m = models[t];
192 const std::shared_ptr<ActionDataAbstract>& d = datas[t];
193 const std::size_t nu = m->get_nu();
194
195 xs_try_[t] = xnext_;
196 m->get_state()->diff(xs_[t], xs_try_[t], dx_[t]);
197 if (nu != 0) {
198 us_try_[t].noalias() = us_[t] - k_[t] * steplength - K_[t] * dx_[t];
199 m->calc(d, xs_try_[t], us_try_[t]);
200 } else {
201 m->calc(d, xs_try_[t]);
202 }
203 xnext_ = d->xnext;
204 cost_try_ += d->cost;
205
206 if (raiseIfNaN(cost_try_)) {
207 STOP_PROFILER("SolverFDDP::forwardPass");
208 throw_pretty("forward_error");
209 }
210 if (raiseIfNaN(xnext_.lpNorm<Eigen::Infinity>())) {
211 STOP_PROFILER("SolverFDDP::forwardPass");
212 throw_pretty("forward_error");
213 }
214 }
215
216 const std::shared_ptr<ActionModelAbstract>& m =
217 problem_->get_terminalModel();
218 const std::shared_ptr<ActionDataAbstract>& d = problem_->get_terminalData();
219 xs_try_.back() = xnext_;
220 m->calc(d, xs_try_.back());
221 cost_try_ += d->cost;
222
223 if (raiseIfNaN(cost_try_)) {
224 STOP_PROFILER("SolverFDDP::forwardPass");
225 throw_pretty("forward_error");
226 }
227 } else {
228 for (std::size_t t = 0; t < T; ++t) {
229 const std::shared_ptr<ActionModelAbstract>& m = models[t];
230 const std::shared_ptr<ActionDataAbstract>& d = datas[t];
231 const std::size_t nu = m->get_nu();
232 m->get_state()->integrate(xnext_, fs_[t] * (steplength - 1), xs_try_[t]);
233 m->get_state()->diff(xs_[t], xs_try_[t], dx_[t]);
234 if (nu != 0) {
235 us_try_[t].noalias() = us_[t] - k_[t] * steplength - K_[t] * dx_[t];
236 m->calc(d, xs_try_[t], us_try_[t]);
237 } else {
238 m->calc(d, xs_try_[t]);
239 }
240 xnext_ = d->xnext;
241 cost_try_ += d->cost;
242
243 if (raiseIfNaN(cost_try_)) {
244 STOP_PROFILER("SolverFDDP::forwardPass");
245 throw_pretty("forward_error");
246 }
247 if (raiseIfNaN(xnext_.lpNorm<Eigen::Infinity>())) {
248 STOP_PROFILER("SolverFDDP::forwardPass");
249 throw_pretty("forward_error");
250 }
251 }
252
253 const std::shared_ptr<ActionModelAbstract>& m =
254 problem_->get_terminalModel();
255 const std::shared_ptr<ActionDataAbstract>& d = problem_->get_terminalData();
256 m->get_state()->integrate(xnext_, fs_.back() * (steplength - 1),
257 xs_try_.back());
258 m->calc(d, xs_try_.back());
259 cost_try_ += d->cost;
260
261 if (raiseIfNaN(cost_try_)) {
262 STOP_PROFILER("SolverFDDP::forwardPass");
263 throw_pretty("forward_error");
264 }
265 }
266 STOP_PROFILER("SolverFDDP::forwardPass");
267}
268
270
271void SolverFDDP::set_th_acceptnegstep(const double th_acceptnegstep) {
272 if (0. > th_acceptnegstep) {
273 throw_pretty(
274 "Invalid argument: " << "th_acceptnegstep value has to be positive.");
275 }
276 th_acceptnegstep_ = th_acceptnegstep;
277}
278
279} // namespace crocoddyl
void set_th_acceptnegstep(const double th_acceptnegstep)
Modify the threshold used for accepting step along ascent direction.
Definition fddp.cpp:271
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
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:269
EIGEN_MAKE_ALIGNED_OPERATOR_NEW SolverFDDP(std::shared_ptr< ShootingProblem > problem)
Initialize the FDDP solver.
Definition fddp.cpp:19
void updateExpectedImprovement()
Update internal values for computing the expected improvement.
Definition fddp.cpp:151