alpaqa cmake-targets
Nonconvex constrained optimization
Loading...
Searching...
No Matches
panoc-ocp.tpp
Go to the documentation of this file.
10#include <alpaqa/util/timed.hpp>
11#include <concepts>
12#include <iomanip>
13#include <iostream>
14#include <numeric>
15#include <stdexcept>
16#include <type_traits>
17
18namespace alpaqa {
19
20template <Config Conf>
22 return detail::extract_u(*problem, xu);
23}
24
25template <Config Conf>
27 return detail::extract_x(*problem, xu);
28}
29
30template <Config Conf>
32 return detail::extract_u(*problem, x̂u);
33}
34
35template <Config Conf>
37 return detail::extract_x(*problem, x̂u);
38}
39
40template <Config Conf>
41std::string PANOCOCPSolver<Conf>::get_name() const {
42 return "PANOCOCPSolver<" + std::string(config_t::get_name()) + '>';
43}
44
45template <Config Conf>
47 /// [in] Problem description
48 const Problem &problem,
49 /// [in] Solve options
50 const SolveOptions &opts,
51 /// [inout] Decision variable @f$ u @f$
52 rvec u,
53 /// [inout] Lagrange multipliers @f$ y @f$
54 rvec y,
55 /// [in] Penalty factors @f$ \mu @f$
56 crvec μ,
57 /// [out] Slack variable error @f$ c(x) - \Pi_D(c(x) + \mu^{-1} y) @f$
58 rvec err_z) -> Stats {
59
60 if (opts.check)
61 problem.check();
62
63 using std::chrono::nanoseconds;
64 auto os = opts.os ? opts.os : this->os;
65 auto start_time = std::chrono::steady_clock::now();
66 Stats s;
67
68 const auto N = problem.get_N();
69 const auto nu = problem.get_nu();
70 const auto nx = problem.get_nx();
71 const auto nc = problem.get_nc();
72 const auto nc_N = problem.get_nc_N();
73 const auto n = nu * N;
74
75 bool enable_lbfgs = params.gn_interval != 1;
76
77 // Allocate storage --------------------------------------------------------
78
79 // TODO: the L-BFGS objects and vectors allocate on each iteration of ALM,
80 // and there are more vectors than strictly necessary.
81
83 auto &vars = eval.vars;
86 LQRFactor lqr{{.N = N, .nx = nx, .nu = nu}};
87 LBFGS<config_t> lbfgs{params.lbfgs_params, enable_lbfgs ? n : 0};
88 mat jacs = vars.create_AB();
89 vec qr = vars.create_qr();
90
91 vec q(n); // Newton step, including states
95
96 // Workspace storage
97 vec work_2x(nx * 2);
98
99 // ALM
100 assert(μ.size() == nc * N + nc_N);
101 assert(y.size() == nc * N + nc_N);
102
103 // Functions for accessing the LQR matrices and index sets
104 auto ABk = [&](index_t i) -> crmat { return vars.ABk(jacs, i); };
105 auto Qk = [&](rvec storage) {
106 return [&, storage](index_t k) {
107 return [&, k](rmat out) {
109 return eval.Qk(storage, y, μ, D, D_N, k, out);
110 };
111 };
112 };
113 auto Rk = [&](rvec storage) {
114 return [&, storage](index_t k) {
115 return [&, k](crindexvec mask, rmat out) {
117 return eval.Rk(storage, k, mask, out);
118 };
119 };
120 };
121 auto Sk = [&](rvec storage) {
122 return [&, storage](index_t k) {
123 return [&, k](crindexvec mask, rmat out) {
125 return eval.Sk(storage, k, mask, out);
126 };
127 };
128 };
129 auto Rk_prod = [&](rvec storage) {
130 return [&, storage](index_t k) {
131 return [&, k](crindexvec mask_J, crindexvec mask_K, crvec v,
132 rvec out) {
134 return eval.Rk_prod(storage, k, mask_J, mask_K, v, out);
135 };
136 };
137 };
138 auto Sk_prod = [&](rvec storage) {
139 return [&, storage](index_t k) {
140 return [&, k](crindexvec mask_K, crvec v, rvec out) {
142 return eval.Sk_prod(storage, k, mask_K, v, out);
143 };
144 };
145 };
146 auto mut_qrk = [&](index_t k) -> rvec { return vars.qrk(qr, k); };
147 auto mut_q_N = [&]() -> rvec { return vars.qk(qr, N); };
148 auto qk = [&](index_t k) -> crvec { return vars.qk(qr, k); };
149 auto rk = [&](index_t k) -> crvec { return vars.rk(qr, k); };
150 auto uk_eq = [&](index_t k) -> crvec { return q.segment(k * nu, nu); };
151 auto Jk = [&](index_t k) -> crindexvec { return J.indices(k); };
152 auto Kk = [&](index_t k) -> crindexvec { return J.compl_indices(k); };
153
154 // Iterates ----------------------------------------------------------------
155
156 // Represents an iterate in the algorithm, keeping track of some
157 // intermediate values and function evaluations.
158 struct Iterate {
159 vec xu; //< Inputs u interleaved with states x
160 vec ; //< Inputs u interleaved with states x after prox grad
161 vec grad_ψ; //< Gradient of cost in u
162 vec p; //< Proximal gradient step in u
163 vec u; //< Inputs u (used for L-BFGS only)
164 real_t ψu = NaN<config_t>; //< Cost in u
165 real_t ψû = NaN<config_t>; //< Cost in û
166 real_t γ = NaN<config_t>; //< Step size γ
167 real_t L = NaN<config_t>; //< Lipschitz estimate L
168 real_t pᵀp = NaN<config_t>; //< Norm squared of p
169 real_t grad_ψᵀp = NaN<config_t>; //< Dot product of gradient and p
170
171 // @pre @ref ψu, @ref pᵀp, @pre grad_ψᵀp
172 // @return φγ
173 real_t fbe() const { return ψu + pᵀp / (2 * γ) + grad_ψᵀp; }
174
176 : xu{vars.create()}, {vars.create()}, grad_ψ{vars.N * vars.nu()},
177 p{vars.N * vars.nu()}, u{enable_lbfgs ? vars.N * vars.nu() : 0} {}
178 } iterates[2]{
179 {vars, enable_lbfgs},
180 {vars, enable_lbfgs},
181 };
182 Iterate *curr = &iterates[0];
183 Iterate *next = &iterates[1];
184
185 // Helper functions --------------------------------------------------------
186
187 auto eval_proj_grad_step_box = [&U](real_t γ, crvec x, crvec grad_ψ, rvec x̂,
188 rvec p) {
189 using binary_real_f = real_t (*)(real_t, real_t);
190 p = (-γ * grad_ψ)
191 .binaryExpr(U.lowerbound - x, binary_real_f(std::fmax))
192 .binaryExpr(U.upperbound - x, binary_real_f(std::fmin));
193 x̂ = x + p;
194 };
195
196 auto eval_prox_impl = [&](real_t γ, crvec xu, crvec grad_ψ, rvec x̂u,
197 rvec p) {
199 real_t pᵀp = 0;
201 for (index_t t = 0; t < N; ++t) {
202 auto &&grad_ψ_t = grad_ψ.segment(t * nu, nu);
203 auto &&p_t = p.segment(t * nu, nu);
204 eval_proj_grad_step_box(γ, vars.uk(xu, t), grad_ψ_t,
205 /* in ⟹ out */ vars.uk(x̂u, t), p_t);
206 // Calculate ∇ψ(x)ᵀp and ‖p‖²
207 pᵀp += p_t.squaredNorm();
208 grad_ψᵀp += grad_ψ_t.dot(p_t);
209 }
210 return std::make_tuple(pᵀp, grad_ψᵀp);
211 };
212
213 auto calc_error_stop_crit = [this, &eval_prox_impl](
216 rvec work_p) {
217 switch (params.stop_crit) {
219 return vec_util::norm_inf(pₖ);
220 }
222 return std::sqrt(pₖᵀpₖ);
223 }
227 }
229 auto [pTp, gTp] =
231 return std::sqrt(pTp);
232 }
234 return vec_util::norm_inf(pₖ) / γ;
235 }
237 return std::sqrt(pₖᵀpₖ) / γ;
238 }
243 default:
244 throw std::invalid_argument("Unsupported stopping criterion");
245 }
246 };
247
248 auto check_all_stop_conditions =
249 [this, &opts](
250 /// [in] Time elapsed since the start of the algorithm
251 auto time_elapsed,
252 /// [in] The current iteration number
253 unsigned iteration,
254 /// [in] Tolerance of the current iterate
256 /// [in] The number of successive iterations no progress was made
257 unsigned no_progress) {
258 auto max_time = params.max_time;
259 if (opts.max_time)
260 max_time = std::min(max_time, *opts.max_time);
261 auto tolerance = opts.tolerance > 0 ? opts.tolerance : real_t(1e-8);
262 bool out_of_time = time_elapsed > max_time;
263 bool out_of_iter = iteration == params.max_iter;
264 bool interrupted = stop_signal.stop_requested();
265 bool not_finite = not std::isfinite(εₖ);
266 bool conv = εₖ <= tolerance;
267 bool max_no_progress = no_progress > params.max_no_progress;
272 : max_no_progress ? SolverStatus::NoProgress
275 };
276
277 auto assign_interleave_xu = [&vars](crvec u, rvec xu) {
278 detail::assign_interleave_xu(vars, u, xu);
279 };
280 auto assign_extract_u = [&vars](crvec xu, rvec u) {
281 detail::assign_extract_u(vars, xu, u);
282 };
283
284 auto write_solution = [&](Iterate &it) {
285 // Update multipliers and constraint error
286 if (nc > 0 || nc_N > 0) {
287 for (index_t t = 0; t < N; ++t) {
288 auto ct = vars.ck(it.xû, t);
289 auto yt = y.segment(nc * t, nc);
290 auto μt = μ.segment(nc * t, nc);
291 auto ζ = ct + μt.asDiagonal().inverse() * yt;
292 auto et = err_z.segment(nc * t, nc);
293 et = projecting_difference(ζ, D);
294 et -= μt.asDiagonal().inverse() * yt;
295 yt += μt.asDiagonal() * et;
296 }
297 auto ct = vars.ck(it.xû, N);
298 auto yt = y.segment(nc * N, nc_N);
299 auto μt = μ.segment(nc * N, nc_N);
300 auto ζ = ct + μt.asDiagonal().inverse() * yt;
301 auto et = err_z.segment(nc * N, nc_N);
302 et = projecting_difference(ζ, D_N);
303 et -= μt.asDiagonal().inverse() * yt;
304 yt += μt.asDiagonal() * et;
305 }
306 assign_extract_u(it.xû, u);
307 };
308
309 // @pre @ref Iterate::γ, @ref Iterate::xu, @ref Iterate::grad_ψ
310 // @post @ref Iterate::xû, @ref Iterate::p, @ref Iterate::pᵀp,
311 // @ref Iterate::grad_ψᵀp
312 auto eval_prox = [&](Iterate &i) {
313 std::tie(i.pᵀp, i.grad_ψᵀp) =
314 eval_prox_impl(i.γ, i.xu, i.grad_ψ, i.xû, i.p);
315 };
316
317 // @pre @ref Iterate::xu
318 // @post @ref Iterate::ψu
319 auto eval_forward = [&](Iterate &i) {
321 i.ψu = eval.forward(i.xu, D, D_N, μ, y);
322 };
323 // @pre @ref Iterate::xû
324 // @post @ref Iterate::ψû
325 auto eval_forward_hat = [&](Iterate &i) {
327 i.ψû = eval.forward(i.xû, D, D_N, μ, y);
328 };
329
330 // @pre @ref Iterate::xu
331 // @post @ref Iterate::grad_ψ, q, q_N
332 auto eval_backward = [&](Iterate &i) {
334 eval.backward(i.xu, i.grad_ψ, mut_qrk, mut_q_N, D, D_N, μ, y);
335 };
336
337 auto qub_violated = [this](const Iterate &i) {
338 real_t margin =
339 (1 + std::abs(i.ψu)) * params.quadratic_upperbound_tolerance_factor;
340 return i.ψû > i.ψu + i.grad_ψᵀp + real_t(0.5) * i.L * i.pᵀp + margin;
341 };
342
343 auto linesearch_violated = [this](const Iterate &curr,
344 const Iterate &next) {
345 real_t β = params.linesearch_strictness_factor;
346 real_t σ = β * (1 - curr.γ * curr.L) / (2 * curr.γ);
347 real_t φγ = curr.fbe();
348 real_t margin = (1 + std::abs(φγ)) * params.linesearch_tolerance_factor;
349 return next.fbe() > φγ - σ * curr.pᵀp + margin;
350 };
351
352 auto initial_lipschitz_estimate =
353 [&](
354 /// Iterate, updates xu, ψ, grad_ψ, have_jacobians, L
355 Iterate *it,
356 /// [in] Finite difference step size relative to x
357 real_t ε,
358 /// [in] Minimum absolute finite difference step size
359 real_t δ,
360 /// [in] Minimum allowed Lipschitz estimate.
361 real_t L_min,
362 /// [in] Maximum allowed Lipschitz estimate.
363 real_t L_max,
364 /// Workspace with the same dimensions as xu, with x_init
366 /// Workspace with the same dimensions as grad_ψ
368 // Calculate ψ(x₀), ∇ψ(x₀)
369 eval_forward(*it);
370 eval_backward(*it);
371 // Select a small step h for finite differences
372 auto h = it->grad_ψ.unaryExpr([&](real_t g) {
373 return g > 0 ? std::max(g * ε, δ) : std::min(g * ε, -δ);
374 });
375 real_t norm_h = h.norm();
376 // work_xu = xu - h
377 for (index_t t = 0; t < N; ++t)
378 vars.uk(work_xu, t) =
379 vars.uk(it->xu, t) - h.segment(t * nu, nu);
380
381 { // Calculate ψ(x₀ - h)
383 eval.forward_simulate(work_xu); // needed for backwards sweep
384 }
385 { // Calculate ∇ψ(x₀ + h)
387 eval.backward(work_xu, work_grad_ψ, mut_qrk, mut_q_N, D, D_N, μ,
388 y);
389 }
390 // Estimate Lipschitz constant using finite differences
391 it->L = (work_grad_ψ - it->grad_ψ).norm() / norm_h;
392 it->L = std::clamp(it->L, L_min, L_max);
393 };
394
395 // Printing ----------------------------------------------------------------
396
397 std::array<char, 64> print_buf;
398 auto print_real = [&](real_t x) {
399 return float_to_str_vw(print_buf, x, params.print_precision);
400 };
401 auto print_real3 = [&](real_t x) {
402 return float_to_str_vw(print_buf, x, 3);
403 };
404 auto print_progress_1 = [&](unsigned k, real_t φₖ, real_t ψₖ, crvec grad_ψₖ,
406 if (k == 0)
407 *os << "┌─[PANOCOCP]\n";
408 else
409 *os << "├─ " << std::setw(6) << k << '\n';
410 *os << "│ φγ = " << print_real(φₖ) //
411 << ", ψ = " << print_real(ψₖ) //
412 << ", ‖∇ψ‖ = " << print_real(grad_ψₖ.norm()) //
413 << ", ‖p‖ = " << print_real(std::sqrt(pₖᵀpₖ)) //
414 << ", γ = " << print_real(γₖ) //
415 << ", ε = " << print_real(εₖ) << '\n';
416 };
417 auto print_progress_2 = [&](crvec qₖ, real_t τₖ, bool did_gn, length_t nJ,
418 real_t min_rcond, bool reject) {
419 const char *color = τₖ == 1 ? "\033[0;32m"
420 : τₖ > 0 ? "\033[0;33m"
421 : "\033[0;35m";
422 *os << "│ ‖q‖ = " << print_real(qₖ.norm()) //
423 << ", #J = " << std::setw(7 + params.print_precision) << nJ //
424 << ", cond = " << print_real3(real_t(1) / min_rcond) //
425 << ", τ = " << color << print_real3(τₖ) << "\033[0m" //
426 << ", " << (did_gn ? "GN" : "L-BFGS") //
427 << ", dir update "
428 << (reject ? "\033[0;31mrejected\033[0m"
429 : "\033[0;32maccepted\033[0m") //
430 << std::endl; // Flush for Python buffering
431 };
432 auto print_progress_n = [&](SolverStatus status) {
433 *os << "└─ " << status << " ──"
434 << std::endl; // Flush for Python buffering
435 };
436
437 auto do_progress_cb = [this, &s, &problem, &lqr,
438 &opts](unsigned k, Iterate &curr, crvec q, real_t τ,
439 real_t εₖ, bool did_gn, index_t nJ,
440 SolverStatus status) {
441 if (!progress_cb)
442 return;
445 progress_cb({
446 .k = k,
447 .status = status,
448 .xu = curr.xu,
449 .p = curr.p,
450 .norm_sq_p = curr.pᵀp,
451 .x̂u = curr.xû,
452 .φγ = curr.fbe(),
453 .ψ = curr.ψu,
454 .grad_ψ = curr.grad_ψ,
455 .ψ_hat = curr.ψû,
456 .q = q,
457 .gn = did_gn,
458 .nJ = nJ,
459 .lqr_min_rcond = lqr.min_rcond,
460 .L = curr.L,
461 .γ = curr.γ,
462 .τ = status == SolverStatus::Busy ? τ : NaN<config_t>,
463 .ε = εₖ,
464 .outer_iter = opts.outer_iter,
465 .problem = &problem,
466 .params = &params,
467 });
468 };
469
470 // Initialize inputs and initial state (do not simulate states yet) --------
471
472 assign_interleave_xu(u, curr->xu); // initial guess
473 problem.get_x_init(curr->xu.topRows(nx)); // initial state
474 curr->xû.topRows(nx) = curr->xu.topRows(nx); // initial state
475 next->xu.topRows(nx) = curr->xu.topRows(nx); // initial state
476 next->xû.topRows(nx) = curr->xu.topRows(nx); // initial state
477 if (enable_lbfgs)
478 curr->u = u;
479
480 problem.get_U(U); // input box constraints
481 problem.get_D(D); // general constraints
482 problem.get_D_N(D_N); // general terminal constraints
483
484 bool do_gn_step = params.gn_interval > 0 and !params.disable_acceleration;
485 bool did_gn = false;
486
487 // Make sure that we don't allocate any memory in the inner loop
489
490 // Estimate Lipschitz constant ---------------------------------------------
491
492 // Finite difference approximation of ∇²ψ in starting point
493 if (params.Lipschitz.L_0 <= 0) {
494 initial_lipschitz_estimate(curr, params.Lipschitz.ε, params.Lipschitz.δ,
495 params.L_min, params.L_max, next->xu,
496 next->grad_ψ);
497 }
498 // Initial Lipschitz constant provided by the user
499 else {
500 curr->L = params.Lipschitz.L_0;
501 // Calculate ψ(x₀), ∇ψ(x₀)
504 }
505 if (not std::isfinite(curr->L)) {
507 return s;
508 }
509 curr->γ = params.Lipschitz.Lγ_factor / curr->L;
510
511 // First proximal gradient step --------------------------------------------
512
513 eval_prox(*curr);
515
516 // Quadratic upper bound
517 while (curr->L < params.L_max && qub_violated(*curr)) {
518 curr->γ /= 2;
519 curr->L *= 2;
520 eval_prox(*curr);
523 }
524
525 unsigned k = 0;
527 length_t nJ = -1;
528
529 // Keep track of how many successive iterations didn't update the iterate
530 unsigned no_progress = 0;
531
532 // Main PANOC loop
533 // =========================================================================
534 while (true) {
535
536 // Check stop condition ------------------------------------------------
537
538 real_t εₖ = calc_error_stop_crit(curr->γ, curr->xu, curr->grad_ψ,
539 curr->p, curr->pᵀp, next->xû, next->p);
540
541 // Print progress ------------------------------------------------------
542 bool do_print =
543 params.print_interval != 0 && k % params.print_interval == 0;
544 if (do_print)
545 print_progress_1(k, curr->fbe(), curr->ψu, curr->grad_ψ, curr->pᵀp,
546 curr->γ, εₖ);
547
548 // Return solution -----------------------------------------------------
549
550 auto time_elapsed = std::chrono::steady_clock::now() - start_time;
551 auto stop_status =
552 check_all_stop_conditions(time_elapsed, k, εₖ, no_progress);
554 do_progress_cb(k, *curr, null_vec<config_t>, -1, εₖ, false, 0,
556 bool do_final_print = params.print_interval != 0;
557 if (!do_print && do_final_print)
558 print_progress_1(k, curr->fbe(), curr->ψu, curr->grad_ψ,
559 curr->pᵀp, curr->γ, εₖ);
564 opts.always_overwrite_results) {
566 }
567 s.iterations = k;
568 s.ε = εₖ;
572 s.final_γ = curr->γ;
573 s.final_ψ = curr->ψû;
574 s.final_h = 0; // only box constraints
575 s.final_φγ = curr->fbe();
576 return s;
577 }
578
579 // Calculate Gauss-Newton step -----------------------------------------
580
581 real_t τ_init = 1;
583 if (params.disable_acceleration) {
584 τ_init = 0;
585 } else if (do_gn_step) {
586 auto is_constr_inactive = [&](index_t t, index_t i) {
587 real_t ui = vars.uk(curr->xu, t)(i);
588 // Gradient descent step.
589 real_t gs = ui - curr->γ * curr->grad_ψ(t * nu + i);
590 // Check whether the box constraints are active for this index.
591 bool active_lb = gs <= U.lowerbound(i);
592 bool active_ub = gs >= U.upperbound(i);
593 if (active_ub) {
594 q(nu * t + i) = U.upperbound(i) - ui;
595 return false;
596 } else if (active_lb) {
597 q(nu * t + i) = U.lowerbound(i) - ui;
598 return false;
599 } else { // Store inactive indices
600 return true;
601 }
602 };
603 { // Find active indices J
605 J.update(is_constr_inactive);
606 nJ = J.sizes().sum();
607 }
608 { // evaluate the Jacobians
610 for (index_t t = 0; t < N; ++t)
611 problem.eval_jac_f(t, vars.xk(curr->xu, t),
612 vars.uk(curr->xu, t), vars.ABk(jacs, t));
613 }
614 { // LQR factor
616 lqr.factor_masked(ABk, Qk(curr->xu), Rk(curr->xu), Sk(curr->xu),
617 Rk_prod(curr->xu), Sk_prod(curr->xu), qk, rk,
618 uk_eq, Jk, Kk, params.lqr_factor_cholesky);
619 }
620 { // LQR solve
622 lqr.solve_masked(ABk, Jk, q, work_2x);
623 }
624 } else {
625 if (!enable_lbfgs)
626 throw std::logic_error("enable_lbfgs");
627
628 // Find inactive indices J
629 auto is_constr_inactive = [&](index_t t, index_t i) {
630 real_t ui = vars.uk(curr->xu, t)(i);
631 real_t grad_i = curr->grad_ψ(t * nu + i);
632 // Gradient descent step.
633 real_t gs = ui - curr->γ * grad_i;
634 // Check whether the box constraints are active for this index.
635 bool active_lb = gs <= U.lowerbound(i);
636 bool active_ub = gs >= U.upperbound(i);
637 if (active_ub || active_lb) {
638 q(t * nu + i) = curr->p(t * nu + i);
639 return false;
640 } else { // Store inactive indices
641 q(t * nu + i) = -grad_i;
642 return true;
643 }
644 };
645
646 auto J_idx = J.indices();
647 nJ = 0;
648 {
650 for (index_t t = 0; t < N; ++t)
651 for (index_t i = 0; i < nu; ++i)
652 if (is_constr_inactive(t, i))
653 J_idx(nJ++) = t * nu + i;
654 }
655 auto J_lbfgs = J_idx.topRows(nJ);
656
657 // If all indices are inactive, we can use standard L-BFGS,
658 // if there are active indices, we need the specialized version
659 // that only applies L-BFGS to the inactive indices
660 bool success = [&] {
662 return lbfgs.apply_masked(q, curr->γ, J_lbfgs);
663 }();
664 // If L-BFGS application failed, qₖ(J) still contains
665 // -∇ψ(x)(J) - HqK(J) or -∇ψ(x)(J), which is not a valid step.
666 if (not success)
667 τ_init = 0;
668 }
669
670 // Make sure quasi-Newton step is valid
671 if (not q.allFinite()) {
672 τ_init = 0;
673 // Is there anything we can do?
674 if (not did_gn)
675 lbfgs.reset();
676 }
677 s.lbfgs_failures += (τ_init == 0 && k > 0);
678
679 bool do_next_gn = params.gn_interval > 0 &&
680 ((k + 1) % params.gn_interval) == 0 &&
681 !params.disable_acceleration;
682 do_gn_step = do_next_gn || (do_gn_step && params.gn_sticky);
683
684 // Line search ---------------------------------------------------------
685
686 next->γ = curr->γ;
687 next->L = curr->L;
688 τ = τ_init;
689 real_t τ_prev = -1;
690 bool dir_rejected = true;
691
692 // xₖ₊₁ = xₖ + pₖ
693 auto take_safe_step = [&] {
694 next->xu = curr->xû;
695 next->ψu = curr->ψû;
696 // Calculate ∇ψ(xₖ₊₁)
698 };
699
700 // xₖ₊₁ = xₖ + (1-τ) pₖ + τ qₖ
701 auto take_accelerated_step = [&](real_t τ) {
702 if (τ == 1) {
703 for (index_t t = 0; t < N; ++t)
704 vars.uk(next->xu, t) =
705 vars.uk(curr->xu, t) + q.segment(t * nu, nu);
706 } else {
708 for (index_t t = 0; t < N; ++t)
709 vars.uk(next->xu, t) =
710 vars.uk(curr->xu, t) +
711 (1 - τ) * curr->p.segment(t * nu, nu) +
712 τ * q.segment(t * nu, nu);
713 }
714 // Calculate ψ(xₖ₊₁), ∇ψ(xₖ₊₁)
715 eval_forward(*next); // Not necessary for DDP
717 };
718
719 // Backtracking line search loop
720 while (!stop_signal.stop_requested()) {
721
722 // Recompute step only if τ changed
723 if (τ != τ_prev) {
724 τ != 0 ? take_accelerated_step(τ) : take_safe_step();
725 τ_prev = τ;
726 }
727
728 // If the cost is not finite, or if the quadratic upper bound could
729 // not be satisfied, abandon the direction entirely, don't even
730 // bother backtracking.
731 bool fail = next->L >= params.L_max || !std::isfinite(next->ψu);
732 if (τ > 0 && fail) {
733 // Don't allow a bad accelerated step to destroy the FBS step
734 // size
735 next->L = curr->L;
736 next->γ = curr->γ;
737 // Line search failed
738 τ = 0;
739 if (enable_lbfgs)
740 lbfgs.reset();
741 continue;
742 }
743
744 // Calculate x̂ₖ₊₁, ψ(x̂ₖ₊₁)
745 eval_prox(*next);
747
748 // Quadratic upper bound step size condition
749 if (next->L < params.L_max && qub_violated(*next)) {
750 next->γ /= 2;
751 next->L *= 2;
752 if (τ > 0)
753 τ = τ_init;
755 continue;
756 }
757
758 // Line search condition
759 if (τ > 0 && linesearch_violated(*curr, *next)) {
760 τ /= 2;
761 if (τ < params.min_linesearch_coefficient)
762 τ = 0;
764 continue;
765 }
766
767 // QUB and line search satisfied (or τ is 0 and L > L_max)
768 break;
769 }
770
771 // If τ < τ_min the line search failed and we accepted the prox step
772 s.linesearch_failures += (τ == 0 && τ_init > 0);
773 s.τ_1_accepted += τ == 1;
774 s.count_τ += (τ_init > 0);
775 s.sum_τ += τ;
776
777 // Check if we made any progress
778 if (no_progress > 0 || k % params.max_no_progress == 0)
779 no_progress = curr->xu == next->xu ? no_progress + 1 : 0;
780
781 // Update L-BFGS -------------------------------------------------------
782
783 if (enable_lbfgs) {
784 const bool force = true;
785 assign_extract_u(next->xu, next->u);
786 bool reset_because_gn = did_gn && params.reset_lbfgs_on_gn_step;
787 if (reset_because_gn || curr->γ != next->γ) {
788 lbfgs.reset();
789 }
790 if (!reset_because_gn) { // TODO: this may be too restrictive
792 s.lbfgs_rejected += dir_rejected = not lbfgs.update(
793 curr->u, next->u, curr->grad_ψ, next->grad_ψ,
795 }
796 }
797
798 // Print ---------------------------------------------------------------
800 if (do_print && (k != 0 || did_gn))
801 print_progress_2(q, τ, did_gn, nJ, lqr.min_rcond, dir_rejected);
802
803 // Advance step --------------------------------------------------------
804 std::swap(curr, next);
805 ++k;
806 }
807 throw std::logic_error("[PANOC] loop error");
808}
809
810} // namespace alpaqa
Limited memory Broyden–Fletcher–Goldfarb–Shanno (L-BFGS) algorithm.
Definition lbfgs.hpp:117
std::string get_name() const
Definition panoc-ocp.tpp:41
Stats operator()(const Problem &problem, const SolveOptions &opts, rvec u, rvec y, crvec μ, rvec err_z)
Definition panoc-ocp.tpp:46
void assign_extract_u(const OCPVariables< Conf > &dim, crvec< Conf > storage, rvec< Conf > u)
Definition ocp-vars.hpp:460
vec< Conf > extract_u(const TypeErasedControlProblem< Conf > &problem, crvec< Conf > xu)
Definition ocp-vars.hpp:474
vec< Conf > extract_x(const TypeErasedControlProblem< Conf > &problem, crvec< Conf > xu)
Definition ocp-vars.hpp:482
void assign_interleave_xu(const OCPVariables< Conf > &dim, crvec< Conf > u, rvec< Conf > storage)
Definition ocp-vars.hpp:445
auto norm_inf(const Eigen::MatrixBase< Derived > &v)
Get the maximum or infinity-norm of the given vector.
Definition config.hpp:163
typename Conf::mat mat
Definition config.hpp:71
std::chrono::nanoseconds time_jacobians
@ LBFGSBpp
The stopping criterion used by LBFGS++, see https://lbfgspp.statr.me/doc/classLBFGSpp_1_1LBFGSBParam....
@ ProjGradUnitNorm
∞-norm of the projected gradient with unit step size:
@ ProjGradNorm
∞-norm of the projected gradient with step size γ:
@ Ipopt
The stopping criterion used by Ipopt, see https://link.springer.com/article/10.1007/s10107-004-0559-y...
@ FPRNorm2
2-norm of fixed point residual:
@ ProjGradNorm2
2-norm of the projected gradient with step size γ:
@ ApproxKKT
Find an ε-approximate KKT point in the ∞-norm:
@ FPRNorm
∞-norm of fixed point residual:
@ ApproxKKT2
Find an ε-approximate KKT point in the 2-norm:
@ ProjGradUnitNorm2
2-norm of the projected gradient with unit step size:
typename Conf::crmat crmat
Definition config.hpp:75
std::chrono::nanoseconds time_hessians
typename Conf::rmat rmat
Definition config.hpp:74
std::chrono::nanoseconds time_lbfgs_apply
std::chrono::nanoseconds time_lbfgs_indices
SolverStatus
Exit status of a numerical solver such as ALM or PANOC.
@ Interrupted
Solver was interrupted by the user.
@ MaxTime
Maximum allowed execution time exceeded.
@ NoProgress
No progress was made in the last iteration.
@ MaxIter
Maximum number of iterations exceeded.
@ Busy
In progress.
@ Converged
Converged and reached given tolerance.
@ NotFinite
Intermediate results were infinite or not-a-number.
std::chrono::nanoseconds time_progress_callback
std::chrono::nanoseconds time_indices
std::chrono::nanoseconds elapsed_time
std::chrono::nanoseconds time_backward
typename Conf::real_t real_t
Definition config.hpp:65
std::chrono::nanoseconds time_lqr_factor
unsigned linesearch_backtracks
typename Conf::index_t index_t
Definition config.hpp:77
typename Conf::length_t length_t
Definition config.hpp:76
constexpr const auto inf
Definition config.hpp:85
typename Conf::rvec rvec
Definition config.hpp:69
std::string_view float_to_str_vw(auto &buf, double value, int precision=std::numeric_limits< double >::max_digits10)
Definition print.tpp:39
typename Conf::crvec crvec
Definition config.hpp:70
std::chrono::nanoseconds time_prox
std::chrono::nanoseconds time_forward
typename Conf::vec vec
Definition config.hpp:66
std::chrono::nanoseconds time_lbfgs_update
std::chrono::nanoseconds time_lqr_solve
typename Conf::crindexvec crindexvec
Definition config.hpp:80
vec create() const
Definition ocp-vars.hpp:65
length_t nu() const
Definition ocp-vars.hpp:57
static Box NaN(length_t n)
Definition box.hpp:19