Crocoddyl
solver-base.cpp
1 // BSD 3-Clause License
3 //
4 // Copyright (C) 2019-2024, LAAS-CNRS, University of Edinburgh,
5 // Heriot-Watt University, University of Oxford
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/solver-base.hpp"
15 #include "crocoddyl/core/utils/exception.hpp"
16 
17 namespace crocoddyl {
18 
19 SolverAbstract::SolverAbstract(std::shared_ptr<ShootingProblem> problem)
20  : problem_(problem),
21  is_feasible_(false),
22  was_feasible_(false),
23  cost_(0.),
24  merit_(0.),
25  stop_(0.),
26  dV_(0.),
27  dPhi_(0.),
28  dVexp_(0.),
29  dPhiexp_(0.),
30  dfeas_(0.),
31  feas_(0.),
32  ffeas_(0.),
33  gfeas_(0.),
34  hfeas_(0.),
35  ffeas_try_(0.),
36  gfeas_try_(0.),
37  hfeas_try_(0.),
38  preg_(0.),
39  dreg_(0.),
40  steplength_(1.),
41  th_acceptstep_(0.1),
42  th_stop_(1e-9),
43  th_gaptol_(1e-16),
44  feasnorm_(LInf),
45  iter_(0),
46  tmp_feas_(0.) {
47  // Allocate common data
48  const std::size_t ndx = problem_->get_ndx();
49  const std::size_t T = problem_->get_T();
50  const std::size_t ng_T = problem_->get_terminalModel()->get_ng_T();
51  xs_.resize(T + 1);
52  us_.resize(T);
53  fs_.resize(T + 1);
54  g_adj_.resize(T + 1);
55  const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
56  problem_->get_runningModels();
57  for (std::size_t t = 0; t < T; ++t) {
58  const std::shared_ptr<ActionModelAbstract>& model = models[t];
59  const std::size_t nu = model->get_nu();
60  const std::size_t ng = model->get_ng();
61  xs_[t] = model->get_state()->zero();
62  us_[t] = Eigen::VectorXd::Zero(nu);
63  fs_[t] = Eigen::VectorXd::Zero(ndx);
64  g_adj_[t] = Eigen::VectorXd::Zero(ng);
65  }
66  xs_.back() = problem_->get_terminalModel()->get_state()->zero();
67  fs_.back() = Eigen::VectorXd::Zero(ndx);
68  g_adj_.back() = Eigen::VectorXd::Zero(ng_T);
69 }
70 
71 SolverAbstract::~SolverAbstract() {}
72 
74  START_PROFILER("SolverAbstract::resizeData");
75  const std::size_t T = problem_->get_T();
76  const std::size_t ng_T = problem_->get_terminalModel()->get_ng_T();
77  const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
78  problem_->get_runningModels();
79  for (std::size_t t = 0; t < T; ++t) {
80  const std::shared_ptr<ActionModelAbstract>& model = models[t];
81  const std::size_t nu = model->get_nu();
82  const std::size_t ng = model->get_ng();
83  us_[t].conservativeResize(nu);
84  g_adj_[t].conservativeResize(ng);
85  }
86 
87  g_adj_.back().conservativeResize(ng_T);
88 
89  STOP_PROFILER("SolverAbstract::resizeData");
90 }
91 
93  tmp_feas_ = 0.;
94  const std::size_t T = problem_->get_T();
95  const Eigen::VectorXd& x0 = problem_->get_x0();
96  const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
97  problem_->get_runningModels();
98  const std::vector<std::shared_ptr<ActionDataAbstract> >& datas =
99  problem_->get_runningDatas();
100 
101  models[0]->get_state()->diff(xs_[0], x0, fs_[0]);
102 #ifdef CROCODDYL_WITH_MULTITHREADING
103 #pragma omp parallel for num_threads(problem_->get_nthreads())
104 #endif
105  for (std::size_t t = 0; t < T; ++t) {
106  const std::shared_ptr<ActionModelAbstract>& m = models[t];
107  const std::shared_ptr<ActionDataAbstract>& d = datas[t];
108  m->get_state()->diff(xs_[t + 1], d->xnext, fs_[t + 1]);
109  }
110  switch (feasnorm_) {
111  case LInf:
112  tmp_feas_ = std::max(tmp_feas_, fs_[0].lpNorm<Eigen::Infinity>());
113  for (std::size_t t = 0; t < T; ++t) {
114  tmp_feas_ = std::max(tmp_feas_, fs_[t + 1].lpNorm<Eigen::Infinity>());
115  }
116  break;
117  case L1:
118  tmp_feas_ = fs_[0].lpNorm<1>();
119  for (std::size_t t = 0; t < T; ++t) {
120  tmp_feas_ += fs_[t + 1].lpNorm<1>();
121  }
122  break;
123  }
124  return tmp_feas_;
125 }
126 
128  tmp_feas_ = 0.;
129  const std::size_t T = problem_->get_T();
130  const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
131  problem_->get_runningModels();
132  const std::vector<std::shared_ptr<ActionDataAbstract> >& datas =
133  problem_->get_runningDatas();
134 
135  switch (feasnorm_) {
136  case LInf:
137  for (std::size_t t = 0; t < T; ++t) {
138  if (models[t]->get_ng() > 0) {
139  g_adj_[t] = datas[t]
140  ->g.cwiseMax(models[t]->get_g_lb())
141  .cwiseMin(models[t]->get_g_ub());
142  tmp_feas_ = std::max(
143  tmp_feas_, (datas[t]->g - g_adj_[t]).lpNorm<Eigen::Infinity>());
144  }
145  }
146  if (problem_->get_terminalModel()->get_ng_T() > 0) {
147  g_adj_.back() =
148  problem_->get_terminalData()
149  ->g.cwiseMax(problem_->get_terminalModel()->get_g_lb())
150  .cwiseMin(problem_->get_terminalModel()->get_g_ub());
151  tmp_feas_ += (problem_->get_terminalData()->g - g_adj_.back())
152  .lpNorm<Eigen::Infinity>();
153  }
154  break;
155  case L1:
156  for (std::size_t t = 0; t < T; ++t) {
157  if (models[t]->get_ng() > 0) {
158  g_adj_[t] = datas[t]
159  ->g.cwiseMax(models[t]->get_g_lb())
160  .cwiseMin(models[t]->get_g_ub());
161  tmp_feas_ =
162  std::max(tmp_feas_, (datas[t]->g - g_adj_[t]).lpNorm<1>());
163  }
164  }
165  if (problem_->get_terminalModel()->get_ng_T() > 0) {
166  g_adj_.back() =
167  problem_->get_terminalData()
168  ->g.cwiseMax(problem_->get_terminalModel()->get_g_lb())
169  .cwiseMin(problem_->get_terminalModel()->get_g_ub());
170  tmp_feas_ +=
171  (problem_->get_terminalData()->g - g_adj_.back()).lpNorm<1>();
172  }
173  break;
174  }
175  return tmp_feas_;
176 }
177 
179  tmp_feas_ = 0.;
180  const std::size_t T = problem_->get_T();
181  const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
182  problem_->get_runningModels();
183  const std::vector<std::shared_ptr<ActionDataAbstract> >& datas =
184  problem_->get_runningDatas();
185  switch (feasnorm_) {
186  case LInf:
187  for (std::size_t t = 0; t < T; ++t) {
188  if (models[t]->get_nh() > 0) {
189  tmp_feas_ =
190  std::max(tmp_feas_, datas[t]->h.lpNorm<Eigen::Infinity>());
191  }
192  }
193  if (problem_->get_terminalModel()->get_nh_T() > 0) {
194  tmp_feas_ =
195  std::max(tmp_feas_,
196  problem_->get_terminalData()->h.lpNorm<Eigen::Infinity>());
197  }
198  break;
199  case L1:
200  for (std::size_t t = 0; t < T; ++t) {
201  if (models[t]->get_nh() > 0) {
202  tmp_feas_ += datas[t]->h.lpNorm<1>();
203  }
204  }
205  if (problem_->get_terminalModel()->get_nh_T() > 0) {
206  tmp_feas_ += problem_->get_terminalData()->h.lpNorm<1>();
207  }
208  break;
209  }
210  return tmp_feas_;
211 }
212 
213 void SolverAbstract::setCandidate(const std::vector<Eigen::VectorXd>& xs_warm,
214  const std::vector<Eigen::VectorXd>& us_warm,
215  bool is_feasible) {
216  const std::size_t T = problem_->get_T();
217 
218  const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
219  problem_->get_runningModels();
220  if (xs_warm.size() == 0) {
221  for (std::size_t t = 0; t < T; ++t) {
222  const std::shared_ptr<ActionModelAbstract>& model = models[t];
223  xs_[t] = model->get_state()->zero();
224  }
225  xs_.back() = problem_->get_terminalModel()->get_state()->zero();
226  } else {
227  if (xs_warm.size() != T + 1) {
228  throw_pretty("Warm start state vector has wrong dimension, got "
229  << xs_warm.size() << " expecting " << (T + 1));
230  }
231  for (std::size_t t = 0; t < T; ++t) {
232  const std::size_t nx = models[t]->get_state()->get_nx();
233  if (static_cast<std::size_t>(xs_warm[t].size()) != nx) {
234  throw_pretty("Invalid argument: "
235  << "xs_init[" + std::to_string(t) +
236  "] has wrong dimension ("
237  << xs_warm[t].size()
238  << " provided - it should be equal to " +
239  std::to_string(nx) + "). ActionModel: "
240  << *models[t]);
241  }
242  }
243  const std::size_t nx = problem_->get_terminalModel()->get_state()->get_nx();
244  if (static_cast<std::size_t>(xs_warm[T].size()) != nx) {
245  throw_pretty("Invalid argument: "
246  << "xs_init[" + std::to_string(T) +
247  "] (terminal state) has wrong dimension ("
248  << xs_warm[T].size()
249  << " provided - it should be equal to " +
250  std::to_string(nx) + "). ActionModel: "
251  << *problem_->get_terminalModel());
252  }
253  std::copy(xs_warm.begin(), xs_warm.end(), xs_.begin());
254  }
255 
256  if (us_warm.size() == 0) {
257  for (std::size_t t = 0; t < T; ++t) {
258  const std::shared_ptr<ActionModelAbstract>& model = models[t];
259  const std::size_t nu = model->get_nu();
260  us_[t] = Eigen::VectorXd::Zero(nu);
261  }
262  } else {
263  if (us_warm.size() != T) {
264  throw_pretty("Warm start control has wrong dimension, got "
265  << us_warm.size() << " expecting " << T);
266  }
267  for (std::size_t t = 0; t < T; ++t) {
268  const std::shared_ptr<ActionModelAbstract>& model = models[t];
269  const std::size_t nu = model->get_nu();
270  if (static_cast<std::size_t>(us_warm[t].size()) != nu) {
271  throw_pretty("Invalid argument: "
272  << "us_init[" + std::to_string(t) +
273  "] has wrong dimension ("
274  << us_warm[t].size()
275  << " provided - it should be equal to " +
276  std::to_string(nu) + "). ActionModel: "
277  << *model);
278  }
279  }
280  std::copy(us_warm.begin(), us_warm.end(), us_.begin());
281  }
282  is_feasible_ = is_feasible;
283 }
284 
286  const std::vector<std::shared_ptr<CallbackAbstract> >& callbacks) {
287  callbacks_ = callbacks;
288 }
289 
290 const std::vector<std::shared_ptr<CallbackAbstract> >&
292  return callbacks_;
293 }
294 
295 const std::shared_ptr<ShootingProblem>& SolverAbstract::get_problem() const {
296  return problem_;
297 }
298 
299 const std::vector<Eigen::VectorXd>& SolverAbstract::get_xs() const {
300  return xs_;
301 }
302 
303 const std::vector<Eigen::VectorXd>& SolverAbstract::get_us() const {
304  return us_;
305 }
306 
307 const std::vector<Eigen::VectorXd>& SolverAbstract::get_fs() const {
308  return fs_;
309 }
310 
312 
313 double SolverAbstract::get_cost() const { return cost_; }
314 
315 double SolverAbstract::get_merit() const { return merit_; }
316 
317 double SolverAbstract::get_stop() const { return stop_; }
318 
319 const Eigen::Vector2d& SolverAbstract::get_d() const { return d_; }
320 
321 double SolverAbstract::get_dV() const { return dV_; }
322 
323 double SolverAbstract::get_dPhi() const { return dPhi_; }
324 
325 double SolverAbstract::get_dVexp() const { return dVexp_; }
326 
327 double SolverAbstract::get_dPhiexp() const { return dPhiexp_; }
328 
329 double SolverAbstract::get_dfeas() const { return dfeas_; }
330 
331 double SolverAbstract::get_feas() const { return feas_; }
332 
333 double SolverAbstract::get_ffeas() const { return ffeas_; }
334 
335 double SolverAbstract::get_gfeas() const { return gfeas_; }
336 
337 double SolverAbstract::get_hfeas() const { return hfeas_; }
338 
339 double SolverAbstract::get_ffeas_try() const { return ffeas_try_; }
340 
341 double SolverAbstract::get_gfeas_try() const { return gfeas_try_; }
342 
343 double SolverAbstract::get_hfeas_try() const { return hfeas_try_; }
344 
345 double SolverAbstract::get_preg() const { return preg_; }
346 
347 double SolverAbstract::get_dreg() const { return dreg_; }
348 
349 DEPRECATED(
350  "Use get_preg for gettting the primal-dual regularization",
351  double SolverAbstract::get_xreg() const { return preg_; })
352 
353 DEPRECATED(
354  "Use get_preg for gettting the primal-dual regularization",
355  double SolverAbstract::get_ureg() const { return preg_; })
356 
357 double SolverAbstract::get_steplength() const { return steplength_; }
358 
360 
361 double SolverAbstract::get_th_stop() const { return th_stop_; }
362 
363 double SolverAbstract::get_th_gaptol() const { return th_gaptol_; }
364 
365 FeasibilityNorm SolverAbstract::get_feasnorm() const { return feasnorm_; }
366 
367 std::size_t SolverAbstract::get_iter() const { return iter_; }
368 
369 void SolverAbstract::set_xs(const std::vector<Eigen::VectorXd>& xs) {
370  const std::size_t T = problem_->get_T();
371  if (xs.size() != T + 1) {
372  throw_pretty("Invalid argument: " << "xs list has to be of length " +
373  std::to_string(T + 1));
374  }
375 
376  const std::size_t nx = problem_->get_nx();
377  for (std::size_t t = 0; t < T; ++t) {
378  if (static_cast<std::size_t>(xs[t].size()) != nx) {
379  throw_pretty("Invalid argument: "
380  << "xs[" + std::to_string(t) + "] has wrong dimension ("
381  << xs[t].size()
382  << " provided - it should be " + std::to_string(nx) + ")")
383  }
384  }
385  if (static_cast<std::size_t>(xs[T].size()) != nx) {
386  throw_pretty("Invalid argument: "
387  << "xs[" + std::to_string(T) +
388  "] (terminal state) has wrong dimension ("
389  << xs[T].size()
390  << " provided - it should be " + std::to_string(nx) + ")")
391  }
392  xs_ = xs;
393 }
394 
395 void SolverAbstract::set_us(const std::vector<Eigen::VectorXd>& us) {
396  const std::size_t T = problem_->get_T();
397  if (us.size() != T) {
398  throw_pretty("Invalid argument: " << "us list has to be of length " +
399  std::to_string(T));
400  }
401 
402  const std::vector<std::shared_ptr<ActionModelAbstract> >& models =
403  problem_->get_runningModels();
404  for (std::size_t t = 0; t < T; ++t) {
405  const std::shared_ptr<ActionModelAbstract>& model = models[t];
406  const std::size_t nu = model->get_nu();
407  if (static_cast<std::size_t>(us[t].size()) != nu) {
408  throw_pretty("Invalid argument: "
409  << "us[" + std::to_string(t) + "] has wrong dimension ("
410  << us[t].size()
411  << " provided - it should be " + std::to_string(nu) + ")")
412  }
413  }
414  us_ = us;
415 }
416 
417 void SolverAbstract::set_preg(const double preg) {
418  if (preg < 0.) {
419  throw_pretty("Invalid argument: " << "preg value has to be positive.");
420  }
421  preg_ = preg;
422 }
423 
424 void SolverAbstract::set_dreg(const double dreg) {
425  if (dreg < 0.) {
426  throw_pretty("Invalid argument: " << "dreg value has to be positive.");
427  }
428  dreg_ = dreg;
429 }
430 
431 DEPRECATED(
432  "Use set_preg for gettting the primal-variable regularization",
433  void SolverAbstract::set_xreg(const double xreg) {
434  if (xreg < 0.) {
435  throw_pretty("Invalid argument: " << "xreg value has to be positive.");
436  }
437  xreg_ = xreg;
438  preg_ = xreg;
439  })
440 
441 DEPRECATED(
442  "Use set_preg for gettting the primal-variable regularization",
443  void SolverAbstract::set_ureg(const double ureg) {
444  if (ureg < 0.) {
445  throw_pretty("Invalid argument: " << "ureg value has to be positive.");
446  }
447  ureg_ = ureg;
448  preg_ = ureg;
449  })
450 
451 void SolverAbstract::set_th_acceptstep(const double th_acceptstep) {
452  if (0. >= th_acceptstep || th_acceptstep > 1) {
453  throw_pretty(
454  "Invalid argument: " << "th_acceptstep value should between 0 and 1.");
455  }
456  th_acceptstep_ = th_acceptstep;
457 }
458 
459 void SolverAbstract::set_th_stop(const double th_stop) {
460  if (th_stop <= 0.) {
461  throw_pretty("Invalid argument: " << "th_stop value has to higher than 0.");
462  }
463  th_stop_ = th_stop;
464 }
465 
466 void SolverAbstract::set_th_gaptol(const double th_gaptol) {
467  if (0. > th_gaptol) {
468  throw_pretty("Invalid argument: " << "th_gaptol value has to be positive.");
469  }
470  th_gaptol_ = th_gaptol;
471 }
472 
473 void SolverAbstract::set_feasnorm(const FeasibilityNorm feasnorm) {
474  feasnorm_ = feasnorm;
475 }
476 
477 bool raiseIfNaN(const double value) {
478  if (std::isnan(value) || std::isinf(value) || value >= 1e30) {
479  return true;
480  } else {
481  return false;
482  }
483 }
484 
485 } // namespace crocoddyl
double get_cost() const
Return the cost for the current guess.
std::vector< Eigen::VectorXd > g_adj_
Adjusted inequality bound.
double get_dPhi() const
Return the reduction in the merit function .
double get_th_gaptol() const
Return the threshold for accepting a gap as non-zero.
double dVexp_
Expected reduction in the cost function.
std::vector< Eigen::VectorXd > xs_
State trajectory.
std::size_t get_iter() const
Return the number of iterations performed by the solver.
double get_hfeas() const
Return the equality feasibility for the current guess.
void set_th_stop(const double th_stop)
Modify the tolerance for stopping the algorithm.
double stop_
Value computed by stoppingCriteria()
void set_xs(const std::vector< Eigen::VectorXd > &xs)
Modify the state trajectory .
double get_dVexp() const
Return the expected reduction in the cost function .
double dreg_
Current dual-variable regularization value.
double feas_
Total feasibility for the current guess.
bool is_feasible_
Label that indicates is the iteration is feasible.
EIGEN_MAKE_ALIGNED_OPERATOR_NEW SolverAbstract(std::shared_ptr< ShootingProblem > problem)
Initialize the solver.
Definition: solver-base.cpp:19
std::shared_ptr< ShootingProblem > problem_
optimal control problem
std::vector< Eigen::VectorXd > us_
Control trajectory.
double get_dPhiexp() const
Return the expected reduction in the merit function .
double th_acceptstep_
Threshold used for accepting step.
double get_steplength() const
Return the step length .
void set_th_gaptol(const double th_gaptol)
Modify the threshold for accepting a gap as non-zero.
double get_merit() const
Return the merit for the current guess.
double dPhi_
Reduction in the merit function computed by tryStep()
const std::vector< std::shared_ptr< CallbackAbstract > > & getCallbacks() const
Return the list of callback functions using for diagnostic.
double computeInequalityFeasibility()
Compute the feasibility of the inequality constraints for the current guess.
double get_preg() const
Return the primal-variable regularization.
double get_hfeas_try() const
Return the equality feasibility for the current step length.
double th_stop_
Tolerance for stopping the algorithm.
double computeDynamicFeasibility()
Compute the dynamic feasibility for the current guess .
Definition: solver-base.cpp:92
const Eigen::Vector2d & get_d() const
Return the linear and quadratic terms of the expected improvement.
double dPhiexp_
Expected reduction in the merit function.
enum FeasibilityNorm feasnorm_
void setCandidate(const std::vector< Eigen::VectorXd > &xs_warm=DEFAULT_VECTOR, const std::vector< Eigen::VectorXd > &us_warm=DEFAULT_VECTOR, const bool is_feasible=false)
Set the solver candidate trajectories .
double get_th_stop() const
Return the tolerance for stopping the algorithm.
double dfeas_
Reduction in the feasibility.
void set_feasnorm(const FeasibilityNorm feas_norm)
Modify the current norm used for computed the dynamic and constraint feasibility.
double get_ffeas() const
Return the dynamic feasibility for the current guess.
double get_gfeas() const
Return the inequality feasibility for the current guess.
double cost_
Cost for the current guess.
std::vector< std::shared_ptr< CallbackAbstract > > callbacks_
Callback functions.
void set_us(const std::vector< Eigen::VectorXd > &us)
Modify the control trajectory .
void setCallbacks(const std::vector< std::shared_ptr< CallbackAbstract > > &callbacks)
Set a list of callback functions using for the solver diagnostic.
double th_gaptol_
Threshold limit to check non-zero gaps.
std::size_t iter_
Number of iteration performed by the solver.
double get_feas() const
Return the total feasibility for the current guess.
double dV_
Reduction in the cost function computed by tryStep()
double get_stop() const
Return the stopping-criteria value computed by stoppingCriteria()
double get_ffeas_try() const
Return the dynamic feasibility for the current step length.
void set_dreg(const double dreg)
Modify the dual-variable regularization value.
Eigen::Vector2d d_
LQ approximation of the expected improvement.
double get_dV() const
Return the reduction in the cost function .
const std::vector< Eigen::VectorXd > & get_fs() const
Return the dynamic infeasibility .
const std::vector< Eigen::VectorXd > & get_xs() const
Return the state trajectory .
double get_dfeas() const
Return the reduction in the feasibility.
void set_preg(const double preg)
Modify the primal-variable regularization value.
double ffeas_
Feasibility of the dynamic constraints for the current guess.
double get_gfeas_try() const
Return the inequality feasibility for the current step length.
bool get_is_feasible() const
Return the feasibility status of the trajectory.
double preg_
Current primal-variable regularization value.
const std::vector< Eigen::VectorXd > & get_us() const
Return the control trajectory .
double merit_
Merit for the current guess.
virtual void resizeData()
Resizing the solver data.
Definition: solver-base.cpp:73
std::vector< Eigen::VectorXd > fs_
Gaps/defects between shooting nodes.
void set_th_acceptstep(const double th_acceptstep)
Modify the threshold used for accepting step.
const std::shared_ptr< ShootingProblem > & get_problem() const
Return the shooting problem.
double tmp_feas_
Temporal variables used for computed the feasibility.
double get_th_acceptstep() const
Return the threshold used for accepting a step.
FeasibilityNorm get_feasnorm() const
Return the type of norm used to evaluate the dynamic and constraints feasibility.
double get_dreg() const
Return the dual-variable regularization.
double computeEqualityFeasibility()
Compute the feasibility of the equality constraints for the current guess.