PANOC-ALM  quadratic-penalty
Nonconvex constrained optimization
structured-panoc-lbfgs.hpp
Go to the documentation of this file.
1 #pragma once
2 
6 
7 #include <cassert>
8 #include <cmath>
9 #include <iomanip>
10 #include <iostream>
11 #include <stdexcept>
12 
13 namespace pa {
14 
15 using std::chrono::duration_cast;
16 using std::chrono::microseconds;
17 
20  const Problem &problem,
22  crvec Σ,
24  real_t ε,
26  bool always_overwrite_results,
28  rvec x,
30  rvec y,
32  rvec err_z) {
33 
34  auto start_time = std::chrono::steady_clock::now();
35  Stats s;
36 
37  const auto n = problem.n;
38  const auto m = problem.m;
39 
40  // Allocate vectors, init L-BFGS -------------------------------------------
41 
42  // TODO: the L-BFGS objects and vectors allocate on each iteration of ALM,
43  // and there are more vectors than strictly necessary.
44 
45  vec xₖ = x, // Value of x at the beginning of the iteration
46  x̂ₖ(n), // Value of x after a projected gradient step
47  xₖ₊₁(n), // xₖ for next iteration
48  x̂ₖ₊₁(n), // x̂ₖ for next iteration
49  ŷx̂ₖ(m), // ŷ(x̂ₖ) = Σ (g(x̂ₖ) - ẑₖ)
50  ŷx̂ₖ₊₁(m), // ŷ(x̂ₖ) for next iteration
51  pₖ(n), // Projected gradient step pₖ = x̂ₖ - xₖ
52  pₖ₊₁(n), // Projected gradient step pₖ₊₁ = x̂ₖ₊₁ - xₖ₊₁
53  qₖ(n), // Newton step Hₖ pₖ
54  grad_ψₖ(n), // ∇ψ(xₖ)
55  grad_̂ψₖ(n), // ∇ψ(x̂ₖ)
56  grad_ψₖ₊₁(n); // ∇ψ(xₖ₊₁)
57 
58  vec work_n(n), work_m(m);
59 
60  vec work_n2(n);
61  vec HqK(n);
62 
63  using indexvec = std::vector<vec::Index>;
64  indexvec J;
65  J.reserve(n);
66  lbfgs.resize(n);
67 
68  // Keep track of how many successive iterations didn't update the iterate
69  unsigned no_progress = 0;
70 
71  // Helper functions --------------------------------------------------------
72 
73  // Wrappers for helper functions that automatically pass along any arguments
74  // that are constant within PANOC (for readability in the main algorithm)
75  auto calc_ψ_ŷ = [&problem, &y, &Σ](crvec x, rvec ŷ) {
76  return detail::calc_ψ_ŷ(problem, x, y, Σ, ŷ);
77  };
78  auto calc_ψ_grad_ψ = [&problem, &y, &Σ, &work_n, &work_m](crvec x,
79  rvec grad_ψ) {
80  return detail::calc_ψ_grad_ψ(problem, x, y, Σ, grad_ψ, work_n, work_m);
81  };
82  auto calc_grad_ψ_from_ŷ = [&problem, &work_n](crvec x, crvec ŷ,
83  rvec grad_ψ) {
84  detail::calc_grad_ψ_from_ŷ(problem, x, ŷ, grad_ψ, work_n);
85  };
86  auto calc_x̂ = [&problem](real_t γ, crvec x, crvec grad_ψ, rvec x̂, rvec p) {
87  detail::calc_x̂(problem, γ, x, grad_ψ, x̂, p);
88  };
89  auto calc_err_z = [&problem, &y, &Σ](crvec x̂, rvec err_z) {
91  };
92  auto descent_lemma = [this, &problem, &y,
93  &Σ](crvec xₖ, real_t ψₖ, crvec grad_ψₖ, rvec x̂ₖ,
94  rvec pₖ, rvec ŷx̂ₖ, real_t &ψx̂ₖ, real_t &pₖᵀpₖ,
95  real_t &grad_ψₖᵀpₖ, real_t &Lₖ, real_t &γₖ) {
96  return detail::descent_lemma(
98  xₖ, ψₖ, grad_ψₖ, y, Σ, x̂ₖ, pₖ, ŷx̂ₖ, ψx̂ₖ, pₖᵀpₖ, grad_ψₖᵀpₖ, Lₖ, γₖ);
99  };
100  auto print_progress = [&](unsigned k, real_t ψₖ, crvec grad_ψₖ,
101  real_t pₖᵀpₖ, real_t γₖ, real_t εₖ) {
102  std::cout << "[PANOC] " << std::setw(6) << k
103  << ": ψ = " << std::setw(13) << ψₖ
104  << ", ‖∇ψ‖ = " << std::setw(13) << grad_ψₖ.norm()
105  << ", ‖p‖ = " << std::setw(13) << std::sqrt(pₖᵀpₖ)
106  << ", γ = " << std::setw(13) << γₖ
107  << ", εₖ = " << std::setw(13) << εₖ << "\r\n";
108  };
109 
110  // Estimate Lipschitz constant ---------------------------------------------
111 
112  real_t ψₖ, Lₖ;
113  // Finite difference approximation of ∇²ψ in starting point
114  if (params.Lipschitz.L₀ <= 0) {
118  /* in ⟹ out */ ψₖ, grad_ψₖ, x̂ₖ, grad_̂ψₖ, work_n, work_m);
119  }
120  // Initial Lipschitz constant provided by the user
121  else {
122  Lₖ = params.Lipschitz.L₀;
123  // Calculate ψ(xₖ), ∇ψ(x₀)
124  ψₖ = calc_ψ_grad_ψ(xₖ, /* in ⟹ out */ grad_ψₖ);
125  }
126  if (not std::isfinite(Lₖ)) {
127  s.status = SolverStatus::NotFinite;
128  return s;
129  }
130  real_t γₖ = params.Lipschitz.Lγ_factor / Lₖ;
131  real_t τ = NaN;
132 
133  // First projected gradient step -------------------------------------------
134 
135  // Calculate x̂₀, p₀ (projected gradient step)
136  calc_x̂(γₖ, xₖ, grad_ψₖ, /* in ⟹ out */ x̂ₖ, pₖ);
137  // Calculate ψ(x̂ₖ) and ŷ(x̂ₖ)
138  real_t ψx̂ₖ = calc_ψ_ŷ(x̂ₖ, /* in ⟹ out */ ŷx̂ₖ);
139  real_t grad_ψₖᵀpₖ = grad_ψₖ.dot(pₖ);
140  real_t pₖᵀpₖ = pₖ.squaredNorm();
141  // Compute forward-backward envelope
142  real_t φₖ = ψₖ + 1 / (2 * γₖ) * pₖᵀpₖ + grad_ψₖᵀpₖ;
143  real_t nmΦₖ = φₖ;
144 
145  // Main PANOC loop
146  // =========================================================================
147  for (unsigned k = 0; k <= params.max_iter; ++k) {
148 
149  // Quadratic upper bound -----------------------------------------------
150  if (k == 0 || params.update_lipschitz_in_linesearch == false) {
151  // Decrease step size until quadratic upper bound is satisfied
152  real_t old_γₖ =
153  descent_lemma(xₖ, ψₖ, grad_ψₖ,
154  /* in ⟹ out */ x̂ₖ, pₖ, ŷx̂ₖ,
155  /* inout */ ψx̂ₖ, pₖᵀpₖ, grad_ψₖᵀpₖ, Lₖ, γₖ);
156  if (γₖ != old_γₖ)
157  φₖ = ψₖ + 1 / (2 * γₖ) * pₖᵀpₖ + grad_ψₖᵀpₖ;
158  }
159  // Calculate ∇ψ(x̂ₖ)
160  calc_grad_ψ_from_ŷ(x̂ₖ, ŷx̂ₖ, /* in ⟹ out */ grad_̂ψₖ);
161 
162  // Check stop condition ------------------------------------------------
163  real_t εₖ = detail::calc_error_stop_crit(params.stop_crit, pₖ, γₖ, xₖ,
164  grad_̂ψₖ, grad_ψₖ, problem.C);
165 
166  // Print progress
167  if (params.print_interval != 0 && k % params.print_interval == 0)
168  print_progress(k, ψₖ, grad_ψₖ, pₖᵀpₖ, γₖ, εₖ);
169  if (progress_cb)
170  progress_cb({k, xₖ, pₖ, pₖᵀpₖ, x̂ₖ, φₖ, ψₖ, grad_ψₖ, ψx̂ₖ, grad_̂ψₖ,
171  Lₖ, γₖ, τ, εₖ, Σ, y, problem, params});
172 
173  auto time_elapsed = std::chrono::steady_clock::now() - start_time;
174  auto stop_status = detail::check_all_stop_conditions(
175  params, time_elapsed, k, stop_signal, ε, εₖ, no_progress);
176  if (stop_status != SolverStatus::Unknown) {
177  // TODO: We could cache g(x) and ẑ, but would that faster?
178  // It saves 1 evaluation of g per ALM iteration, but requires
179  // many extra stores in the inner loops of PANOC.
180  // TODO: move the computation of ẑ and g(x) to ALM?
181  if (stop_status == SolverStatus::Converged ||
182  stop_status == SolverStatus::Interrupted ||
183  always_overwrite_results) {
184  calc_err_z(x̂ₖ, /* in ⟹ out */ err_z);
185  x = std::move(x̂ₖ);
186  y = std::move(ŷx̂ₖ);
187  }
188  s.iterations = k;
189  s.ε = εₖ;
190  s.elapsed_time = duration_cast<microseconds>(time_elapsed);
191  s.status = stop_status;
192  return s;
193  }
194 
195  // Calculate Newton step -----------------------------------------------
196 
197  if (k > 0) { // No L-BFGS estimate on first iteration → no Newton step
198  J.clear();
199  // Find inactive indices J
200  for (vec::Index i = 0; i < n; ++i) {
201  real_t gd = xₖ(i) - γₖ * grad_ψₖ(i);
202  if (gd < problem.C.lowerbound(i)) { // i ∊ J̲ ⊆ K
203  qₖ(i) = pₖ(i); //
204  } else if (problem.C.upperbound(i) < gd) { // i ∊ J̅ ⊆ K
205  qₖ(i) = pₖ(i); //
206  } else { // i ∊ J
207  J.push_back(i);
208  qₖ(i) = 0;
209  }
210  }
211 
212  if (not J.empty()) { // There are inactive indices J
213  if (J.size() == n) { // There are no active indices K
214  qₖ = -grad_ψₖ;
215  } else { // There are active indices K
218  problem, xₖ, y, Σ, grad_ψₖ, qₖ, HqK, work_n,
219  work_n2, work_m);
220  } else {
221  problem.hess_L_prod(xₖ, y, qₖ, HqK);
223  auto &g = work_m;
224  problem.g(xₖ, g);
225  for (vec::Index i = 0; i < m; ++i) {
226  real_t ζ = g(i) + y(i) / Σ(i);
227  bool inactive = problem.D.lowerbound(i) < ζ &&
228  ζ < problem.D.upperbound(i);
229  if (not inactive) {
230  problem.grad_gi(xₖ, i, work_n);
231  auto t = Σ(i) * work_n.dot(qₖ);
232  // TODO: the dot product is more work than
233  // strictly necessary (only over K)
234  for (auto j : J)
235  HqK(j) += work_n(j) * t;
236  }
237  }
238  }
239  }
240 
241  for (auto j : J) // Compute right-hand side of 6.1c
242  qₖ(j) = -grad_ψₖ(j) - HqK(j);
243  }
244 
245  real_t stepsize = params.lbfgs_stepsize ==
246  LBFGSStepSize::BasedOnGradientStepSize
247  ? γₖ
248  : -1;
249  // If all indices are inactive, we can use standard L-BFGS,
250  // if there are active indices, we need the specialized version
251  // that only applies L-BFGS to the inactive indices
252  bool success = lbfgs.apply(qₖ, stepsize, J);
253  // If L-BFGS application failed, qₖ(J) still contains
254  // -∇ψ(x)(J) - HqK(J) or -∇ψ(x)(J), which is not a valid step.
255  // A good alternative is to use H₀ = γI as an L-BFGS estimate.
256  // This seems to be better than just falling back to a projected
257  // gradient step.
258  if (not success) {
259  if (J.size() == n)
260  qₖ *= γₖ;
261  else
262  for (auto j : J)
263  qₖ(j) *= γₖ;
264  }
265  }
266  }
267 
268  // Line search initialization ------------------------------------------
269  τ = 1;
270  real_t σₖγₖ⁻¹pₖᵀpₖ = (1 - γₖ * Lₖ) * pₖᵀpₖ / (2 * γₖ);
271  real_t φₖ₊₁, ψₖ₊₁, ψx̂ₖ₊₁, grad_ψₖ₊₁ᵀpₖ₊₁, pₖ₊₁ᵀpₖ₊₁;
272  real_t Lₖ₊₁, γₖ₊₁;
273  real_t ls_cond;
275  nmΦₖ = k == 0 ? φₖ : w * nmΦₖ + (1 - w) * φₖ;
276  // TODO: make separate parameter
277  real_t margin =
278  (1 + std::abs(nmΦₖ)) * params.quadratic_upperbound_tolerance_factor;
279 
280  // Make sure quasi-Newton step is valid
281  if (k == 0) {
282  τ = 0; // Always use prox step on first iteration
283  } else if (not qₖ.allFinite()) {
284  τ = 0;
285  ++s.lbfgs_failures;
286  } else if (J.empty()) {
287  τ = 0; // All constraints are active, no Newton step possible
288  }
289 
290  // Line search loop ----------------------------------------------------
291  do {
292  Lₖ₊₁ = Lₖ;
293  γₖ₊₁ = γₖ;
294 
295  // Calculate xₖ₊₁
296  if (τ / 2 < params.τ_min) { // line search failed
297  xₖ₊₁.swap(x̂ₖ); // → safe prox step
298  ψₖ₊₁ = ψx̂ₖ;
299  grad_ψₖ₊₁.swap(grad_̂ψₖ);
300  } else { // line search didn't fail (yet)
301  if (τ == 1) // → faster quasi-Newton step
302  xₖ₊₁ = xₖ + qₖ;
303  else
304  xₖ₊₁ = xₖ + (1 - τ) * pₖ + τ * qₖ;
305  // Calculate ψ(xₖ₊₁), ∇ψ(xₖ₊₁)
306  ψₖ₊₁ = calc_ψ_grad_ψ(xₖ₊₁, /* in ⟹ out */ grad_ψₖ₊₁);
307  }
308 
309  // Calculate x̂ₖ₊₁, pₖ₊₁ (projected gradient step in xₖ₊₁)
310  calc_x̂(γₖ₊₁, xₖ₊₁, grad_ψₖ₊₁, /* in ⟹ out */ x̂ₖ₊₁, pₖ₊₁);
311  // Calculate ψ(x̂ₖ₊₁) and ŷ(x̂ₖ₊₁)
312  ψx̂ₖ₊₁ = calc_ψ_ŷ(x̂ₖ₊₁, /* in ⟹ out */ ŷx̂ₖ₊₁);
313 
314  // Quadratic upper bound -------------------------------------------
315  grad_ψₖ₊₁ᵀpₖ₊₁ = grad_ψₖ₊₁.dot(pₖ₊₁);
316  pₖ₊₁ᵀpₖ₊₁ = pₖ₊₁.squaredNorm();
317  real_t pₖ₊₁ᵀpₖ₊₁_ₖ = pₖ₊₁ᵀpₖ₊₁; // prox step with step size γₖ
318 
320  // Decrease step size until quadratic upper bound is satisfied
321  (void)descent_lemma(xₖ₊₁, ψₖ₊₁, grad_ψₖ₊₁,
322  /* in ⟹ out */ x̂ₖ₊₁, pₖ₊₁, ŷx̂ₖ₊₁,
323  /* inout */ ψx̂ₖ₊₁, pₖ₊₁ᵀpₖ₊₁,
324  grad_ψₖ₊₁ᵀpₖ₊₁, Lₖ₊₁, γₖ₊₁);
325  }
326 
327  // Compute forward-backward envelope
328  φₖ₊₁ = ψₖ₊₁ + 1 / (2 * γₖ₊₁) * pₖ₊₁ᵀpₖ₊₁ + grad_ψₖ₊₁ᵀpₖ₊₁;
329  // Compute line search condition
330  ls_cond = φₖ₊₁ - (nmΦₖ - σₖγₖ⁻¹pₖᵀpₖ);
332  ls_cond -= (0.5 / γₖ₊₁ - 0.5 / γₖ) * pₖ₊₁ᵀpₖ₊₁_ₖ;
333 
334  τ /= 2;
335  } while (ls_cond > margin && τ >= params.τ_min);
336 
337  // If τ < τ_min the line search failed and we accepted the prox step
338  if (τ < params.τ_min && k != 0) {
340  τ = 0;
341  }
342  if (k != 0) {
343  s.count_τ += 1;
344  s.sum_τ += τ * 2;
345  s.τ_1_accepted += τ * 2 == 1;
346  // std::cout << J.size() << "," << τ*2 << "\r\n";
347  }
348 
349  // Check if we made any progress
350  if (no_progress > 0 || k % params.max_no_progress == 0)
351  no_progress = xₖ == xₖ₊₁ ? no_progress + 1 : 0;
352 
353  // Update L-BFGS
354  const bool force = true;
355  s.lbfgs_rejected += not lbfgs.update(xₖ, xₖ₊₁, grad_ψₖ, grad_ψₖ₊₁,
356  LBFGS::Sign::Positive, force);
357 
358  // Advance step --------------------------------------------------------
359  Lₖ = Lₖ₊₁;
360  γₖ = γₖ₊₁;
361 
362  ψₖ = ψₖ₊₁;
363  ψx̂ₖ = ψx̂ₖ₊₁;
364  φₖ = φₖ₊₁;
365 
366  xₖ.swap(xₖ₊₁);
367  x̂ₖ.swap(x̂ₖ₊₁);
368  ŷx̂ₖ.swap(ŷx̂ₖ₊₁);
369  pₖ.swap(pₖ₊₁);
370  grad_ψₖ.swap(grad_ψₖ₊₁);
371  grad_ψₖᵀpₖ = grad_ψₖ₊₁ᵀpₖ₊₁;
372  pₖᵀpₖ = pₖ₊₁ᵀpₖ₊₁;
373  }
374  throw std::logic_error("[PANOC] loop error");
375 }
376 
377 } // namespace pa
pa::detail::print_progress
void print_progress(unsigned k, real_t ψₖ, crvec grad_ψₖ, real_t pₖᵀpₖ, real_t γₖ, real_t εₖ)
Definition: standalone/panoc.hpp:139
pa::StructuredPANOCLBFGSSolver::Stats::lbfgs_failures
unsigned lbfgs_failures
Definition: decl/structured-panoc-lbfgs.hpp:90
bicycle-obstacle-avoidance-mpc.t
t
Definition: bicycle-obstacle-avoidance-mpc.py:109
codegen-rosenbrock.w
w
Definition: codegen-rosenbrock.py:13
pa::detail::calc_x̂
void calc_x̂(const ProblemT &prob, real_t γ, crvec x, crvec grad_ψ, rvec x̂, rvec p)
Definition: panoc-helpers.hpp:271
pa::StructuredPANOCLBFGSSolver::operator()
Stats operator()(const Problem &problem, crvec Σ, real_t ε, bool always_overwrite_results, rvec x, rvec y, rvec err_z)
Definition: structured-panoc-lbfgs.hpp:18
structured-panoc-lbfgs.hpp
pa::StructuredPANOCLBFGSParams::L_min
real_t L_min
Minimum Lipschitz constant estimate.
Definition: decl/structured-panoc-lbfgs.hpp:30
pa::SolverStatus::Unknown
@ Unknown
Initial value.
pa::rvec
Eigen::Ref< vec > rvec
Default type for mutable references to vectors.
Definition: vec.hpp:16
pa::StructuredPANOCLBFGSSolver::Stats::status
SolverStatus status
Definition: decl/structured-panoc-lbfgs.hpp:85
pa::StructuredPANOCLBFGSParams::hessian_vec_finited_differences
bool hessian_vec_finited_differences
Definition: decl/structured-panoc-lbfgs.hpp:51
pa::NaN
constexpr real_t NaN
Not a number.
Definition: vec.hpp:28
lbfgs.hpp
pa::StructuredPANOCLBFGSSolver::Stats::elapsed_time
std::chrono::microseconds elapsed_time
Definition: decl/structured-panoc-lbfgs.hpp:87
panocpy.test.y
y
Definition: test.py:76
panocpy.test.err_z
err_z
Definition: test.py:76
panocpy.test.ε
int ε
Definition: test.py:71
pa::vec
realvec vec
Default type for vectors.
Definition: vec.hpp:14
pa::StructuredPANOCLBFGSParams::full_augmented_hessian
bool full_augmented_hessian
Definition: decl/structured-panoc-lbfgs.hpp:52
pa::LipschitzEstimateParams::L₀
real_t L₀
Initial estimate of the Lipschitz constant of ∇ψ(x)
Definition: lipschitz.hpp:9
pa::StructuredPANOCLBFGSParams::nonmonotone_linesearch
real_t nonmonotone_linesearch
Factor used in update for exponentially weighted nonmonotone line search.
Definition: decl/structured-panoc-lbfgs.hpp:35
pa::StructuredPANOCLBFGSParams::print_interval
unsigned print_interval
When to print progress.
Definition: decl/structured-panoc-lbfgs.hpp:43
pa::StructuredPANOCLBFGSParams::quadratic_upperbound_tolerance_factor
real_t quadratic_upperbound_tolerance_factor
Definition: decl/structured-panoc-lbfgs.hpp:45
pa::StructuredPANOCLBFGSSolver::lbfgs
LBFGS lbfgs
Definition: decl/structured-panoc-lbfgs.hpp:128
pa::StructuredPANOCLBFGSSolver::Stats::lbfgs_rejected
unsigned lbfgs_rejected
Definition: decl/structured-panoc-lbfgs.hpp:91
pa::StructuredPANOCLBFGSSolver::Stats
Definition: decl/structured-panoc-lbfgs.hpp:84
main.problem
problem
Definition: main.py:16
pa::StructuredPANOCLBFGSSolver::Stats::count_τ
unsigned count_τ
Definition: decl/structured-panoc-lbfgs.hpp:93
pa::StructuredPANOCLBFGSParams::τ_min
real_t τ_min
Minimum weight factor between Newton step and projected gradient step.
Definition: decl/structured-panoc-lbfgs.hpp:28
pa::StructuredPANOCLBFGSParams::L_max
real_t L_max
Maximum Lipschitz constant estimate.
Definition: decl/structured-panoc-lbfgs.hpp:32
pa
Definition: alm.hpp:10
pa::detail::calc_ψ_grad_ψ
real_t calc_ψ_grad_ψ(const Problem &p, crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m)
Calculate both ψ(x) and its gradient ∇ψ(x).
Definition: panoc-helpers.hpp:113
panocpy.test.x
x
Definition: test.py:40
pa::StructuredPANOCLBFGSParams::lbfgs_stepsize
LBFGSStepSize lbfgs_stepsize
Definition: decl/structured-panoc-lbfgs.hpp:54
pa::StructuredPANOCLBFGSSolver::params
Params params
Definition: decl/structured-panoc-lbfgs.hpp:123
pa::detail::check_all_stop_conditions
SolverStatus check_all_stop_conditions(const ParamsT &params, DurationT time_elapsed, unsigned iteration, const AtomicStopSignal &stop_signal, real_t ε, real_t εₖ, unsigned no_progress)
Check all stop conditions (required tolerance reached, out of time, maximum number of iterations exce...
Definition: panoc-helpers.hpp:469
pa::detail::calc_ψ_ŷ
real_t calc_ψ_ŷ(const Problem &p, crvec x, crvec y, crvec Σ, rvec ŷ)
Calculate both ψ(x) and the vector ŷ that can later be used to compute ∇ψ.
Definition: panoc-helpers.hpp:16
pa::LBFGS::update
bool update(crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁, Sign sign, bool forced=false)
Update the inverse Hessian approximation using the new vectors xₖ₊₁ and pₖ₊₁.
Definition: lbfgs.hpp:32
panoc-helpers.hpp
pa::detail::calc_err_z
void calc_err_z(const Problem &p, crvec x̂, crvec y, crvec Σ, rvec err_z)
Calculate the error between ẑ and g(x).
Definition: panoc-helpers.hpp:212
pa::detail::calc_error_stop_crit
real_t calc_error_stop_crit(PANOCStopCrit crit, crvec pₖ, real_t γ, crvec xₖ, crvec grad_̂ψₖ, crvec grad_ψₖ, const Box &C)
Definition: panoc-helpers.hpp:284
pa::LipschitzEstimateParams::ε
real_t ε
Relative step size for initial finite difference Lipschitz estimate.
Definition: lipschitz.hpp:11
pa::StructuredPANOCLBFGSSolver::Stats::linesearch_failures
unsigned linesearch_failures
Definition: decl/structured-panoc-lbfgs.hpp:89
pa::LBFGS::resize
void resize(size_t n)
Re-allocate storage for a problem with a different size.
Definition: lbfgs.hpp:187
pa::detail::descent_lemma
real_t descent_lemma(const Problem &problem, real_t rounding_tolerance, real_t L_max, crvec xₖ, real_t ψₖ, crvec grad_ψₖ, crvec y, crvec Σ, rvec x̂ₖ, rvec pₖ, rvec ŷx̂ₖ, real_t &ψx̂ₖ, real_t &norm_sq_pₖ, real_t &grad_ψₖᵀpₖ, real_t &Lₖ, real_t &γₖ)
Increase the estimate of the Lipschitz constant of the objective gradient and decrease the step size ...
Definition: panoc-helpers.hpp:327
pa::crvec
Eigen::Ref< const vec > crvec
Default type for immutable references to vectors.
Definition: vec.hpp:18
pa::StructuredPANOCLBFGSParams::Lipschitz
LipschitzEstimateParams Lipschitz
Parameters related to the Lipschitz constant estimate and step size.
Definition: decl/structured-panoc-lbfgs.hpp:22
pa::LipschitzEstimateParams::δ
real_t δ
Minimum step size for initial finite difference Lipschitz estimate.
Definition: lipschitz.hpp:13
pa::StructuredPANOCLBFGSSolver::stop_signal
AtomicStopSignal stop_signal
Definition: decl/structured-panoc-lbfgs.hpp:124
pa::detail::calc_augmented_lagrangian_hessian_prod_fd
void calc_augmented_lagrangian_hessian_prod_fd(const Problem &problem, crvec xₖ, crvec y, crvec Σ, crvec grad_ψ, crvec v, rvec Hv, rvec work_n1, rvec work_n2, rvec work_m)
Compute the Hessian matrix of the augmented Lagrangian function multiplied by the given vector,...
Definition: panoc-helpers.hpp:541
panocpy.test.Σ
int Σ
Definition: test.py:70
pa::StructuredPANOCLBFGSSolver::Stats::ε
real_t ε
Definition: decl/structured-panoc-lbfgs.hpp:86
panocpy.test.m
int m
Definition: test.py:39
pa::StructuredPANOCLBFGSParams::stop_crit
PANOCStopCrit stop_crit
What stopping criterion to use.
Definition: decl/structured-panoc-lbfgs.hpp:37
pa::StructuredPANOCLBFGSSolver::progress_cb
std::function< void(const ProgressInfo &)> progress_cb
Definition: decl/structured-panoc-lbfgs.hpp:125
pa::detail::initial_lipschitz_estimate
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_n1, rvec work_n2, rvec work_n3, rvec work_m)
Estimate the Lipschitz constant of the gradient using finite differences.
Definition: panoc-helpers.hpp:574
panocpy.test.g
g
Definition: test.py:51
pa::StructuredPANOCLBFGSParams::max_no_progress
unsigned max_no_progress
Maximum number of iterations without any progress before giving up.
Definition: decl/structured-panoc-lbfgs.hpp:39
pa::detail::calc_grad_ψ_from_ŷ
void calc_grad_ψ_from_ŷ(const Problem &p, crvec x, crvec ŷ, rvec grad_ψ, rvec work_n)
Calculate ∇ψ(x) using ŷ.
Definition: panoc-helpers.hpp:81
panocpy.test.p
p
Definition: test.py:57
pa::StructuredPANOCLBFGSSolver::Stats::τ_1_accepted
unsigned τ_1_accepted
Definition: decl/structured-panoc-lbfgs.hpp:92
pa::StructuredPANOCLBFGSSolver::Stats::iterations
unsigned iterations
Definition: decl/structured-panoc-lbfgs.hpp:88
pa::LipschitzEstimateParams::Lγ_factor
real_t Lγ_factor
Factor that relates step size γ and Lipschitz constant.
Definition: lipschitz.hpp:15
panocpy.test.n
int n
Definition: test.py:38
pa::LBFGS::apply
bool apply(Vec &&q, real_t γ)
Apply the inverse Hessian approximation to the given vector q.
Definition: lbfgs.hpp:58
pa::real_t
double real_t
Default floating point type.
Definition: vec.hpp:8
pa::LBFGS::Sign::Positive
@ Positive
pa::StructuredPANOCLBFGSParams::max_iter
unsigned max_iter
Maximum number of inner PANOC iterations.
Definition: decl/structured-panoc-lbfgs.hpp:24
pa::StructuredPANOCLBFGSSolver::Stats::sum_τ
real_t sum_τ
Definition: decl/structured-panoc-lbfgs.hpp:94
pa::Problem
Problem description for minimization problems.
Definition: include/panoc-alm/util/problem.hpp:24
pa::StructuredPANOCLBFGSParams::update_lipschitz_in_linesearch
bool update_lipschitz_in_linesearch
Definition: decl/structured-panoc-lbfgs.hpp:48
pa::StructuredPANOCLBFGSParams::alternative_linesearch_cond
bool alternative_linesearch_cond
Definition: decl/structured-panoc-lbfgs.hpp:49