alpaqa 1.1.0a1
Nonconvex constrained optimization
Loading...
Searching...
No Matches
panoc.tpp
Go to the documentation of this file.
1#pragma once
2
4
5#include <cassert>
6#include <cmath>
7#include <iomanip>
8#include <iostream>
9#include <stdexcept>
10
14#include <alpaqa/util/print.hpp>
15#include <guanaqo/timed.hpp>
16
17namespace alpaqa {
18
19template <class DirectionProviderT>
21 return "PANOCSolver<" + std::string(direction.get_name()) + ">";
22}
23
24template <class DirectionProviderT>
26 /// [in] Problem description
27 const Problem &problem,
28 /// [in] Solve options
29 const SolveOptions &opts,
30 /// [inout] Decision variable @f$ x @f$
31 rvec x,
32 /// [inout] Lagrange multipliers @f$ y @f$
33 rvec y,
34 /// [in] Constraint weights @f$ \Sigma @f$
35 crvec Σ,
36 /// [out] Slack variable error @f$ g(x) - \Pi_D(g(x) + \Sigma^{-1} y) @f$
37 rvec err_z) -> Stats {
38
39 if (opts.check)
40 problem.check();
41
42 using std::chrono::nanoseconds;
43 auto os = opts.os ? opts.os : this->os;
44 auto start_time = std::chrono::steady_clock::now();
45 Stats s;
46
47 const auto n = problem.get_num_variables();
48 const auto m = problem.get_num_constraints();
49
50 // Represents an iterate in the algorithm, keeping track of some
51 // intermediate values and function evaluations.
52 struct Iterate {
53 vec x; //< Decision variables
54 vec x̂; //< Decision variables after proximal gradient step
55 vec grad_ψ; //< Gradient of cost in x
56 vec grad_ψx̂; //< Gradient of cost in x̂
57 vec p; //< Proximal gradient step in x
58 vec ŷx̂; //< Candidate Lagrange multipliers in x̂
59 real_t ψx = NaN<config_t>; //< Cost in x
60 real_t ψx̂ = NaN<config_t>; //< Cost in x̂
61 real_t γ = NaN<config_t>; //< Step size γ
62 real_t L = NaN<config_t>; //< Lipschitz estimate L
63 real_t pᵀp = NaN<config_t>; //< Norm squared of p
64 real_t grad_ψᵀp = NaN<config_t>; //< Dot product of gradient and p
65 real_t hx̂ = NaN<config_t>; //< Non-smooth function value in x̂
66 bool have_grad_ψx̂ = false;
67
68 // @pre @ref ψx, @ref hx̂ @ref pᵀp, @ref grad_ψᵀp
69 // @return φγ
70 real_t fbe() const { return ψx + hx̂ + pᵀp / (2 * γ) + grad_ψᵀp; }
71
72 Iterate(length_t n, length_t m)
73 : x(n), x̂(n), grad_ψ(n), grad_ψx̂(n), p(n), ŷx̂(m) {}
74 } iterates[2]{{n, m}, {n, m}};
75 Iterate *curr = &iterates[0];
76 Iterate *next = &iterates[1];
77
78 bool need_grad_ψx̂ = Helpers::stop_crit_requires_grad_ψx̂(params.stop_crit);
79 vec work_n(n), work_m(m);
80 vec q(n); // (quasi-)Newton step Hₖ pₖ
81
82 // Helper functions --------------------------------------------------------
83
84 auto qub_violated = [this](const Iterate &i) {
85 real_t margin =
86 (1 + std::abs(i.ψx)) * params.quadratic_upperbound_tolerance_factor;
87 return i.ψx̂ > i.ψx + i.grad_ψᵀp + real_t(0.5) * i.L * i.pᵀp + margin;
88 };
89
90 auto linesearch_violated = [this](const Iterate &curr,
91 const Iterate &next) {
92 if (params.force_linesearch)
93 return false;
94 real_t β = params.linesearch_strictness_factor;
95 real_t σ = β * (1 - curr.γ * curr.L) / (2 * curr.γ);
96 real_t φγ = curr.fbe();
97 real_t margin = (1 + std::abs(φγ)) * params.linesearch_tolerance_factor;
98 return next.fbe() > φγ - σ * curr.pᵀp + margin;
99 };
100
101 // Problem functions -------------------------------------------------------
102
103 auto eval_ψ_grad_ψ = [&problem, &y, &Σ, &work_n, &work_m](Iterate &i) {
104 i.ψx = problem.eval_augmented_lagrangian_and_gradient(
105 i.x, y, Σ, i.grad_ψ, work_n, work_m);
106 };
107 auto eval_prox_grad_step = [&problem](Iterate &i) {
108 i.hx̂ =
109 problem.eval_proximal_gradient_step(i.γ, i.x, i.grad_ψ, i.x̂, i.p);
110 i.pᵀp = i.p.squaredNorm();
111 i.grad_ψᵀp = i.p.dot(i.grad_ψ);
112 };
113 auto eval_ψx̂ = [&problem, &y, &Σ, &work_n, this](Iterate &i) {
114 if (params.eager_gradient_eval)
115 i.ψx̂ = problem.eval_augmented_lagrangian_and_gradient(
116 i.x̂, y, Σ, i.grad_ψx̂, work_n, i.ŷx̂);
117 else
118 i.ψx̂ = problem.eval_augmented_lagrangian(i.x̂, y, Σ, i.ŷx̂);
119 i.have_grad_ψx̂ = params.eager_gradient_eval;
120 };
121 auto eval_grad_ψx̂ = [&problem, &work_n](Iterate &i) {
122 problem.eval_lagrangian_gradient(i.x̂, i.ŷx̂, i.grad_ψx̂, work_n);
123 i.have_grad_ψx̂ = true;
124 };
125
126 // Printing ----------------------------------------------------------------
127
128 std::array<char, 64> print_buf;
129 auto print_real = [this, &print_buf](real_t x) {
130 return float_to_str_vw(print_buf, x, params.print_precision);
131 };
132 auto print_real3 = [&print_buf](real_t x) {
133 return float_to_str_vw(print_buf, x, 3);
134 };
135 auto print_progress_1 = [&print_real, os](unsigned k, real_t φₖ, real_t ψₖ,
136 crvec grad_ψₖ, real_t pₖᵀpₖ,
137 real_t γₖ, real_t εₖ) {
138 if (k == 0)
139 *os << "┌─[PANOC]\n";
140 else
141 *os << "├─ " << std::setw(6) << k << '\n';
142 *os << "│ φγ = " << print_real(φₖ) //
143 << ", ψ = " << print_real(ψₖ) //
144 << ", ‖∇ψ‖ = " << print_real(grad_ψₖ.norm()) //
145 << ", ‖p‖ = " << print_real(std::sqrt(pₖᵀpₖ)) //
146 << ", γ = " << print_real(γₖ) //
147 << ", ε = " << print_real(εₖ) << '\n';
148 };
149 auto print_progress_2 = [&print_real, &print_real3, os](crvec qₖ, real_t τₖ,
150 bool reject) {
151 const char *color = τₖ == 1 ? "\033[0;32m"
152 : τₖ > 0 ? "\033[0;33m"
153 : "\033[0;35m";
154 *os << "│ ‖q‖ = " << print_real(qₖ.norm()) //
155 << ", τ = " << color << print_real3(τₖ) << "\033[0m" //
156 << ", dir update "
157 << (reject ? "\033[0;31mrejected\033[0m"
158 : "\033[0;32maccepted\033[0m") //
159 << std::endl; // Flush for Python buffering
160 };
161 auto print_progress_n = [&](SolverStatus status) {
162 *os << "└─ " << status << " ──"
163 << std::endl; // Flush for Python buffering
164 };
165
166 auto do_progress_cb = [this, &s, &problem, &Σ, &y,
167 &opts](unsigned k, Iterate &it, crvec q, real_t τ,
168 real_t εₖ, SolverStatus status) {
169 if (!progress_cb)
170 return;
172 guanaqo::Timed t{s.time_progress_callback};
173 auto &&grad_ψx̂ =
174 it.have_grad_ψx̂ ? crvec{it.grad_ψx̂} : crvec{null_vec<config_t>};
176 .k = k,
177 .status = status,
178 .x = it.x,
179 .p = it.p,
180 .norm_sq_p = it.pᵀp,
181 .x̂ = it.x̂,
182 .ŷ = it.ŷx̂,
183 .φγ = it.fbe(),
184 .ψ = it.ψx,
185 .grad_ψ = it.grad_ψ,
186 .ψ_hat = it.ψx̂,
187 .grad_ψ_hat = grad_ψx̂,
188 .q = q,
189 .L = it.L,
190 .γ = it.γ,
191 .τ = τ,
192 .ε = εₖ,
193 .Σ = Σ,
194 .y = y,
195 .outer_iter = opts.outer_iter,
196 .problem = &problem,
197 .params = &params,
198 });
199 };
200
201 // Initialization ----------------------------------------------------------
202
203 curr->x = x;
204
205 // Estimate Lipschitz constant ---------------------------------------------
206
207 // Finite difference approximation of ∇²ψ in starting point
208 if (params.Lipschitz.L_0 <= 0) {
210 problem, curr->x, y, Σ, params.Lipschitz.ε, params.Lipschitz.δ,
211 params.L_min, params.L_max,
212 /* in ⟹ out */ curr->ψx, curr->grad_ψ, curr->x̂, next->grad_ψ,
213 work_n, work_m);
214 }
215 // Initial Lipschitz constant provided by the user
216 else {
217 curr->L = params.Lipschitz.L_0;
218 // Calculate ψ(xₖ), ∇ψ(x₀)
219 eval_ψ_grad_ψ(*curr);
220 }
221 if (not std::isfinite(curr->L)) {
223 return s;
224 }
225 curr->γ = params.Lipschitz.Lγ_factor / curr->L;
226
227 // First proximal gradient step --------------------------------------------
228
229 eval_prox_grad_step(*curr);
230 eval_ψx̂(*curr);
231
232 // Quadratic upper bound
233 while (curr->L < params.L_max && qub_violated(*curr)) {
234 curr->γ /= 2;
235 curr->L *= 2;
236 eval_prox_grad_step(*curr);
237 eval_ψx̂(*curr);
239 }
240
241 // Loop data ---------------------------------------------------------------
242
243 unsigned k = 0; // iteration
244 real_t τ = NaN<config_t>; // line search parameter
245 // Keep track of how many successive iterations didn't update the iterate
246 unsigned no_progress = 0;
247
248 // Main PANOC loop
249 // =========================================================================
250
251 ScopedMallocBlocker mb; // Don't allocate in the inner loop
252 while (true) {
253
254 // Check stopping criteria ---------------------------------------------
255
256 // Calculate ∇ψ(x̂ₖ)
257 if (need_grad_ψx̂ && !curr->have_grad_ψx̂)
258 eval_grad_ψx̂(*curr);
259
261 problem, params.stop_crit, curr->p, curr->γ, curr->x, curr->x̂,
262 curr->ŷx̂, curr->grad_ψ, curr->grad_ψx̂, work_n, next->p);
263
264 // Print progress ------------------------------------------------------
265 bool do_print =
266 params.print_interval != 0 && k % params.print_interval == 0;
267 if (do_print)
268 print_progress_1(k, curr->fbe(), curr->ψx, curr->grad_ψ, curr->pᵀp,
269 curr->γ, εₖ);
270
271 // Return solution -----------------------------------------------------
272
273 auto time_elapsed = std::chrono::steady_clock::now() - start_time;
274 auto stop_status = Helpers::check_all_stop_conditions(
275 params, opts, time_elapsed, k, stop_signal, εₖ, no_progress);
276 if (stop_status != SolverStatus::Busy) {
277 do_progress_cb(k, *curr, null_vec<config_t>, -1, εₖ, stop_status);
278 bool do_final_print = params.print_interval != 0;
279 if (!do_print && do_final_print)
280 print_progress_1(k, curr->fbe(), curr->ψx, curr->grad_ψ,
281 curr->pᵀp, curr->γ, εₖ);
282 if (do_print || do_final_print)
283 print_progress_n(stop_status);
284 if (stop_status == SolverStatus::Converged ||
285 stop_status == SolverStatus::Interrupted ||
286 opts.always_overwrite_results) {
287 // Calculate ŷ
288 if (params.eager_gradient_eval)
289 eval_ψx̂(*curr);
290 auto &ŷ = curr->ŷx̂;
291 if (err_z.size() > 0)
292 err_z = (ŷ - y).cwiseQuotient(Σ);
293 x = curr->x̂;
294 y = curr->ŷx̂;
295 }
296 s.iterations = k;
297 s.ε = εₖ;
298 s.elapsed_time = duration_cast<nanoseconds>(time_elapsed);
299 s.status = stop_status;
300 s.final_γ = curr->γ;
301 s.final_ψ = curr->ψx̂;
302 s.final_h = curr->hx̂;
303 s.final_φγ = curr->fbe();
304 return s;
305 }
306
307 // Calculate quasi-Newton step -----------------------------------------
308
309 real_t τ_init = NaN<config_t>;
310 if (k == 0) { // Initialize L-BFGS
312 direction.initialize(problem, y, Σ, curr->γ, curr->x, curr->x̂,
313 curr->p, curr->grad_ψ);
314 τ_init = 0;
315 }
316 if (k > 0 || direction.has_initial_direction()) {
317 τ_init = direction.apply(curr->γ, curr->x, curr->x̂, curr->p,
318 curr->grad_ψ, q)
319 ? 1
320 : 0;
321 // Make sure quasi-Newton step is valid
322 if (τ_init == 1 && not q.allFinite())
323 τ_init = 0;
324 if (τ_init != 1) { // If we computed a quasi-Newton step
326 direction.reset(); // Is there anything else we can do?
327 }
328 }
329
330 // Line search ---------------------------------------------------------
331
332 next->γ = curr->γ;
333 next->L = curr->L;
334 τ = τ_init;
335 real_t τ_prev = -1;
336 bool update_lbfgs_in_linesearch = params.update_direction_in_candidate;
337 bool updated_lbfgs = false;
338 bool dir_rejected = true;
339
340 // xₖ₊₁ = xₖ + pₖ
341 auto take_safe_step = [&] {
342 // Calculate ∇ψ(xₖ₊₁)
343 if (not curr->have_grad_ψx̂)
344 eval_grad_ψx̂(*curr);
345 next->x = curr->x̂; // → safe prox step
346 next->ψx = curr->ψx̂;
347 next->grad_ψ.swap(curr->grad_ψx̂);
348 curr->have_grad_ψx̂ = next->have_grad_ψx̂ = false;
349 };
350
351 // xₖ₊₁ = xₖ + (1-τ) pₖ + τ qₖ
352 auto take_accelerated_step = [&](real_t τ) {
353 if (τ == 1) // → faster quasi-Newton step
354 next->x = curr->x + q;
355 else
356 next->x = curr->x + (1 - τ) * curr->p + τ * q;
357 // Calculate ψ(xₖ₊₁), ∇ψ(xₖ₊₁)
358 eval_ψ_grad_ψ(*next);
359 next->have_grad_ψx̂ = false;
360 };
361
362 while (!stop_signal.stop_requested()) {
363
364 // Recompute step only if τ changed
365 if (τ != τ_prev) {
366 τ != 0 ? take_accelerated_step(τ) : take_safe_step();
367 τ_prev = τ;
368 }
369
370 // If the cost is not finite, or if the quadratic upper bound could
371 // not be satisfied, abandon the direction entirely, don't even
372 // bother backtracking.
373 bool fail = !std::isfinite(next->ψx);
374 fail |= next->L >= params.L_max && !(curr->L >= params.L_max);
375 if (τ > 0 && fail) {
376 // Don't allow a bad accelerated step to destroy the FBS step
377 // size
378 next->L = curr->L;
379 next->γ = curr->γ;
380 // Line search failed
381 τ = 0;
382 direction.reset();
383 // Update the direction in the FB iterate later
384 update_lbfgs_in_linesearch = false;
385 continue;
386 }
387
388 // Calculate x̂ₖ₊₁, ψ(x̂ₖ₊₁)
389 eval_prox_grad_step(*next);
390 eval_ψx̂(*next);
391
392 // Quadratic upper bound step size condition
393 if (next->L < params.L_max && qub_violated(*next)) {
394 next->γ /= 2;
395 next->L *= 2;
396 if (τ > 0)
397 τ = τ_init;
399 // If the step size changes, we need extra care when updating
400 // the direction later
401 update_lbfgs_in_linesearch = false;
402 continue;
403 }
404
405 // Update L-BFGS in candidate (even if we don't accept this point)
406 if (update_lbfgs_in_linesearch && !updated_lbfgs) {
407 s.direction_update_rejected += dir_rejected =
408 not direction.update(curr->γ, next->γ, curr->x, next->x,
409 curr->p, next->p, curr->grad_ψ,
410 next->grad_ψ);
411 update_lbfgs_in_linesearch = false;
412 updated_lbfgs = true;
413 }
414
415 // Line search condition
416 if (τ > 0 && linesearch_violated(*curr, *next)) {
417 τ *= params.linesearch_coefficient_update_factor;
418 if (τ < params.min_linesearch_coefficient)
419 τ = 0;
421 continue;
422 }
423
424 // QUB and line search satisfied (or τ is 0 and L > L_max)
425 break;
426 }
427 // If τ < τ_min the line search failed and we accepted the prox step
428 s.linesearch_failures += (τ == 0 && τ_init > 0);
429 s.τ_1_accepted += τ == 1;
430 s.count_τ += (τ_init > 0);
431 s.sum_τ += τ;
432
433 // Check if we made any progress
434 if (no_progress > 0 || k % params.max_no_progress == 0)
435 no_progress = curr->x == next->x ? no_progress + 1 : 0;
436
437 // Update L-BFGS -------------------------------------------------------
438
439 if (!updated_lbfgs) {
440 if (curr->γ != next->γ) { // Flush L-BFGS if γ changed
441 direction.changed_γ(next->γ, curr->γ);
442 if (params.recompute_last_prox_step_after_stepsize_change) {
443 curr->γ = next->γ;
444 curr->L = next->L;
445 eval_prox_grad_step(*curr);
446 }
447 }
448 s.direction_update_rejected += dir_rejected = not direction.update(
449 curr->γ, next->γ, curr->x, next->x, curr->p, next->p,
450 curr->grad_ψ, next->grad_ψ);
451 }
452
453 // Print ---------------------------------------------------------------
454 do_progress_cb(k, *curr, q, τ, εₖ, SolverStatus::Busy);
455 if (do_print && (k != 0 || direction.has_initial_direction()))
456 print_progress_2(q, τ, dir_rejected);
457
458 // Advance step --------------------------------------------------------
459 std::swap(curr, next);
460 ++k;
461 }
462 throw std::logic_error("[PANOC] loop error");
463}
464
465} // namespace alpaqa
std::string get_name() const
Definition panoc.tpp:20
PANOCStats< config_t > Stats
Definition panoc.hpp:149
std::function< void(const ProgressInfo &)> progress_cb
Definition panoc.hpp:200
Stats operator()(const Problem &problem, const SolveOptions &opts, rvec x, rvec y, crvec Σ, rvec err_z)
Definition panoc.tpp:25
Direction direction
Definition panoc.hpp:204
InnerSolveOptions< config_t > SolveOptions
Definition panoc.hpp:151
guanaqo::AtomicStopSignal stop_signal
Definition panoc.hpp:199
PANOCProgressInfo< config_t > ProgressInfo
Definition panoc.hpp:150
TypeErasedProblem< config_t > Problem
Definition panoc.hpp:146
std::ostream * os
Definition panoc.hpp:205
unsigned stepsize_backtracks
Definition panoc.hpp:99
unsigned direction_update_rejected
Definition panoc.hpp:101
unsigned τ_1_accepted
Definition panoc.hpp:102
SolverStatus
Exit status of a numerical solver such as ALM or PANOC.
@ Interrupted
Solver was interrupted by the user.
@ Converged
Converged and reached given tolerance.
@ NotFinite
Intermediate results were infinite or not-a-number.
std::chrono::nanoseconds time_progress_callback
Definition panoc.hpp:95
std::chrono::nanoseconds elapsed_time
Definition panoc.hpp:94
typename Conf::real_t real_t
Definition config.hpp:86
const rvec< Conf > null_vec
Global empty vector for convenience.
Definition config.hpp:193
constexpr const auto NaN
Definition config.hpp:114
unsigned direction_failures
Definition panoc.hpp:100
unsigned linesearch_backtracks
Definition panoc.hpp:98
typename Conf::length_t length_t
Definition config.hpp:103
typename Conf::rvec rvec
Definition config.hpp:91
typename Conf::crvec crvec
Definition config.hpp:92
unsigned linesearch_failures
Definition panoc.hpp:97
typename Conf::vec vec
Definition config.hpp:88
unsigned iterations
Definition panoc.hpp:96
SolverStatus status
Definition panoc.hpp:92
unsigned count_τ
Definition panoc.hpp:103
static bool stop_crit_requires_grad_ψx̂(PANOCStopCrit crit)
static real_t initial_lipschitz_estimate(const Problem &problem, crvec x, crvec y, crvec Σ, real_t ε, real_t δ, real_t L_min, real_t L_max, real_t &ψ, rvec grad_ψ, rvec work_x, rvec work_grad_ψ, rvec work_n, rvec work_m)
static real_t calc_error_stop_crit(const Problem &problem, PANOCStopCrit crit, crvec pₖ, real_t γ, crvec xₖ, crvec x̂ₖ, crvec ŷₖ, crvec grad_ψₖ, crvec grad_̂ψₖ, rvec work_n1, rvec work_n2)
static SolverStatus check_all_stop_conditions(const ParamsT &params, const InnerSolveOptions< config_t > &opts, DurationT time_elapsed, unsigned iteration, const guanaqo::AtomicStopSignal &stop_signal, real_t εₖ, unsigned no_progress)