C++ API Reference

namespace pa

Typedefs

using real_t = double

Default floating point type.

using realvec = Eigen::Matrix<real_t, Eigen::Dynamic, 1>

Default type for floating point vectors.

using realmat = Eigen::Matrix<real_t, Eigen::Dynamic, Eigen::Dynamic>

Default type for floating point matrices.

using vec = realvec

Default type for vectors.

using rvec = Eigen::Ref<vec>

Default type for mutable references to vectors.

using crvec = Eigen::Ref<const vec>

Default type for immutable references to vectors.

using mat = realmat

Default type for matrices.

using rmat = Eigen::Ref<mat>

Default type for mutable references to matrices.

using crmat = Eigen::Ref<const mat>

Default type for immutable references to matrices.

Enums

enum LBFGSStepSize

Which method to use to select the L-BFGS step size.

Values:

enumerator BasedOnGradientStepSize
enumerator BasedOnCurvature
enum PANOCStopCrit

Values:

enumerator ApproxKKT

Find a ε-approximate KKT point.

enumerator ProjGradNorm

∞-norm of the projected gradient with step size γ

enumerator ProjGradUnitNorm

∞-norm of the projected gradient with unit step size

enumerator FPRNorm

∞-norm of fixed point residual

enum SolverStatus

Exit status of a numerical solver such as ALM or PANOC.

Values:

enumerator Unknown

Initial value.

enumerator Converged

Converged and reached given tolerance.

enumerator MaxTime

Maximum allowed execution time exceeded.

enumerator MaxIter

Maximum number of iterations exceeded.

enumerator NotFinite

Intermediate results were infinite or not-a-number.

enumerator NoProgress

No progress was made in the last iteration.

enumerator Interrupted

Solver was interrupted by the user.

Functions

inline const char *enum_name(PANOCStopCrit s)
inline std::ostream &operator<<(std::ostream &os, PANOCStopCrit s)
inline InnerStatsAccumulator<PANOCStats> &operator+=(InnerStatsAccumulator<PANOCStats> &acc, const PANOCStats &s)
inline InnerStatsAccumulator<SecondOrderPANOCSolver::Stats> &operator+=(InnerStatsAccumulator<SecondOrderPANOCSolver::Stats> &acc, const SecondOrderPANOCSolver::Stats &s)
inline InnerStatsAccumulator<StructuredPANOCLBFGSSolver::Stats> &operator+=(InnerStatsAccumulator<StructuredPANOCLBFGSSolver::Stats> &acc, const StructuredPANOCLBFGSSolver::Stats &s)
inline void minimize_update_anderson(LimitedMemoryQR &qr, rmat G, crvec rₖ, crvec rₖ₋₁, crvec gₖ, rvec γ_LS, rvec xₖ_aa)

Solve one step of Anderson acceleration to find a fixed point of a function g(x):

\( g(x^\star) - x^\star = 0 \)

Updates the QR factorization of \( \mathcal{R}_k = QR \), solves the least squares problem to find \( \gamma_\text{LS} \), computes the next iterate \( x_{k+1} \), and stores the current function value \( g_k \) in the matrix \( G \), which is used as a circular buffer.

\[\begin{split} \begin{aligned} \mathcal{R}_k &= \begin{pmatrix} \Delta r_{k-m_k} & \dots & \Delta r_{k-1} \end{pmatrix} \in \mathbb{R}^{n\times m_k} \\ \Delta r_i &= r_{i+1} - r_i \\ r_i &= g_i - x_i \\ g_i &= g(x_i) \\ \DeclareMathOperator*{\argmin}{argmin} \gamma_\text{LS} &= \argmin_\gamma \left\| \mathcal{R}_k \gamma - r_k \right\|_2 \\ \alpha_i &= \begin{cases} \gamma_\text{LS}[0] & i = 0 \\ \gamma_\text{LS}[i] - \gamma_\text{LS}[i-1] & 0 < i < m_k \\ 1 - \gamma_\text{LS}[m_k - 1] & i = m_k \end{cases} \\ x_{k+1} &= \sum_{i=0}^{m_k} \alpha_i\,g_i \end{aligned} \end{split}\]

inline InnerStatsAccumulator<GAAPGASolver::Stats> &operator+=(InnerStatsAccumulator<GAAPGASolver::Stats> &acc, const GAAPGASolver::Stats &s)
inline InnerStatsAccumulator<LBFGSBStats> &operator+=(InnerStatsAccumulator<LBFGSBStats> &acc, const LBFGSBStats &s)
inline InnerStatsAccumulator<PGASolver::Stats> &operator+=(InnerStatsAccumulator<PGASolver::Stats> &acc, const PGASolver::Stats &s)
std::function<pa::Problem::f_sig> load_CasADi_objective(const std::string &so_name, const std::string &fun_name = "f")

Load an objective function generated by CasADi.

std::function<pa::Problem::grad_f_sig> load_CasADi_gradient_objective(const std::string &so_name, const std::string &fun_name = "grad_f")

Load the gradient of an objective function generated by CasADi.

std::function<pa::Problem::g_sig> load_CasADi_constraints(const std::string &so_name, const std::string &fun_name = "g")

Load a constraint function generated by CasADi.

std::function<pa::Problem::grad_g_prod_sig> load_CasADi_gradient_constraints_prod(const std::string &so_name, const std::string &fun_name = "grad_g")

Load the gradient-vector product of a constraint function generated by CasADi.

std::function<pa::Problem::hess_L_sig> load_CasADi_hessian_lagrangian(const std::string &so_name, const std::string &fun_name = "hess_L")

Load the Hessian of a Lagrangian function generated by CasADi.

std::function<pa::Problem::hess_L_prod_sig> load_CasADi_hessian_lagrangian_prod(const std::string &so_name, const std::string &fun_name = "hess_L_prod")

Load the Hessian-vector product of a Lagrangian function generated by CasADi.

pa::Problem load_CasADi_problem(const std::string &filename, unsigned n, unsigned m, bool second_order = false)

Load a problem generated by CasADi (without parameters).

The file should contain functions with the names f, grad_f, g and grad_g. These functions evaluate the objective function, its gradient, the constraints, and the constraint gradient times a vector respecitvely. If second_order is true, additional functions hess_L and hess_L_prod should be provided to evaluate the Hessian of the Lagrangian and Hessian-vector products.

Parameters
  • filename – Filename of the shared library to load the functions from.

  • n – Number of decision variables ( \( x \in \mathbb{R}^n \)).

  • m – Number of general constraints ( \( g(x) \in \mathbb{R}^m \)).

  • second_order – Load the additional functions required for second-order PANOC.

pa::ProblemWithParam load_CasADi_problem_with_param(const std::string &filename, unsigned n, unsigned m, bool second_order = false)

Load a problem generated by CasADi (with parameters).

The file should contain functions with the names f, grad_f, g and grad_g. These functions evaluate the objective function, its gradient, the constraints, and the constraint gradient times a vector respecitvely. If second_order is true, additional functions hess_L and hess_L_prod should be provided to evaluate the Hessian of the Lagrangian and Hessian-vector products.

Parameters
  • filename – Filename of the shared library to load the functions from.

  • n – Number of decision variables ( \( x \in \mathbb{R}^n \)).

  • m – Number of general constraints ( \( g(x) \in \mathbb{R}^m \)).

  • second_order – Load the additional functions required for second-order PANOC.

pa::ProblemFull load_CasADi_problem_full(const char *filename, unsigned n, unsigned m1, unsigned m2, bool second_order = false)

Load a problem generated by CasADi (without parameters).

The file should contain functions with the names f, grad_f, g1, grad_g1, g2 and grad_g2. These functions evaluate the objective function, its gradient, the constraints, and the constraint gradient times a vector respecitvely. If second_order is true, additional functions hess_L and hess_L_prod should be provided to evaluate the Hessian of the Lagrangian and Hessian-vector products.

Parameters
  • filename – Filename of the shared library to load the functions from.

  • n – Number of decision variables ( \( x \in \mathbb{R}^n \)).

  • m1 – Number of ALM constraints ( \( g1(x) \in \mathbb{R}^m1 \)).

  • m2 – Number of quadratic penalty constraints ( \( g2(x) \in \mathbb{R}^m2 \)).

  • second_order – Load the additional functions required for second-order PANOC.

pa::ProblemFullWithParam load_CasADi_problem_full_with_param(const char *filename, unsigned n, unsigned m1, unsigned m2, bool second_order = false)

Load a problem generated by CasADi (with parameters).

The file should contain functions with the names f, grad_f, g1, grad_g1, g2 and grad_g2. These functions evaluate the objective function, its gradient, the constraints, and the constraint gradient times a vector respecitvely. If second_order is true, additional functions hess_L and hess_L_prod should be provided to evaluate the Hessian of the Lagrangian and Hessian-vector products.

Parameters
  • filename – Filename of the shared library to load the functions from.

  • n – Number of decision variables ( \( x \in \mathbb{R}^n \)).

  • m1 – Number of ALM constraints ( \( g1(x) \in \mathbb{R}^m1 \)).

  • m2 – Number of quadratic penalty constraints ( \( g2(x) \in \mathbb{R}^m2 \)).

  • second_order – Load the additional functions required for second-order PANOC.

template<class ObjFunT, class ObjGradFunT, class DirectionT>
inline PANOCStats panoc_impl(ObjFunT &ψ, ObjGradFunT &grad_ψ, const Box &C, rvec x, real_t ε, const PANOCParams &params, vec_allocator &alloc, DirectionT &direction_provider)
template<class DirectionProviderT = LBFGS, class ObjFunT, class ObjGradFunT>
PANOCStats panoc(ObjFunT &ψ, ObjGradFunT &grad_ψ, const Box &C, rvec x, real_t ε, const PANOCParams &params, PANOCDirection<DirectionProviderT> direction, vec_allocator &alloc)
template<class DirectionProviderT = LBFGS, class ObjFunT, class ObjGradFunT>
PANOCStats panoc(ObjFunT &ψ, ObjGradFunT &grad_ψ, const Box &C, rvec x, real_t ε, const PANOCParams &params, PANOCDirection<DirectionProviderT> direction)
template<class Vec>
inline auto project(const Vec &v, const Box &box)

Project a vector onto a box.

\[ \Pi_C(v) \]

Parameters
  • v – [in] The vector to project

  • box – [in] The box to project onto

template<class Vec>
inline auto projecting_difference(const Vec &v, const Box &box)

Get the difference between the given vector and its projection.

\[ v - \Pi_C(v) \]

Warning

Beware catastrophic cancellation!

Parameters
  • v – [in] The vector to project

  • box – [in] The box to project onto

inline real_t dist_squared(crvec v, const Box &box)

Get the distance squared between the given vector and its projection.

\[ \left\| v - \Pi_C(v) \right\|_2^2 \]

Warning

Beware catastrophic cancellation!

Parameters
  • v – [in] The vector to project

  • box – [in] The box to project onto

inline real_t dist_squared(crvec v, const Box &box, crvec Σ)

Get the distance squared between the given vector and its projection in the Σ norm.

\[ \left\| v - \Pi_C(v) \right\|_\Sigma^2 = \left(v - \Pi_C(v)\right)^\top \Sigma \left(v - \Pi_C(v)\right) \]

Warning

Beware catastrophic cancellation!

Parameters
  • v – [in] The vector to project

  • box – [in] The box to project onto

  • Σ – [in] Diagonal matrix defining norm

inline EvalCounter &operator+=(EvalCounter &a, EvalCounter b)
inline EvalCounter operator+(EvalCounter a, EvalCounter b)
inline EvalCounterFull &operator+=(EvalCounterFull &a, EvalCounterFull b)
inline EvalCounterFull operator+(EvalCounterFull a, EvalCounterFull b)
const char *enum_name(SolverStatus s)
std::ostream &operator<<(std::ostream &os, SolverStatus s)

Variables

constexpr size_t panoc_min_alloc_size = 10
constexpr real_t inf = std::numeric_limits<real_t>::infinity()

\( \infty \)

constexpr real_t NaN = std::numeric_limits<real_t>::quiet_NaN()

Not a number.

struct ALMParams
#include <panoc-alm/decl/alm.hpp>

Parameters for the Augmented Lagrangian solver.

Public Members

real_t ε = 1e-5

Primal tolerance.

real_t δ = 1e-5

Dual tolerance.

real_t Δ = 10

Factor used in updating the penalty parameters.

real_t Δ_lower = 0.8

Factor to reduce ALMParams::Δ when inner convergence fails.

real_t Σ₀ = 1

Initial penalty parameter.

When set to zero (which is the default), it is computed automatically, based on the constraint violation in the starting point and the parameter ALMParams::σ₀.

real_t σ₀ = 20

Initial penalty parameter factor.

Active if ALMParams::Σ₀ is set to zero.

real_t Σ₀_lower = 0.6

Factor to reduce the initial penalty factor by if convergence fails in in the first iteration.

real_t ε₀ = 1

Initial primal tolerance.

real_t ε₀_increase = 1.1

Factor to increase the initial primal tolerance if convergence fails in the first iteration.

real_t ρ = 1e-1

Update factor for primal tolerance.

real_t ρ_increase = 2

Factor to increase the primal tolerance update factor by if convergence fails.

real_t θ = 0.1

Error tolerance for penalty increase.

real_t M = 1e9

Lagrange multiplier bound.

real_t Σ_max = 1e9

Maximum penalty factor.

real_t Σ_min = 1e-9

Minimum penalty factor (used during initialization).

unsigned int max_iter = 100

Maximum number of outer ALM iterations.

std::chrono::microseconds max_time = std::chrono::minutes(5)

Maximum duration.

unsigned max_num_initial_retries = 20

How many times can the initial penalty ALMParams::Σ₀ or ALMParams::σ₀ and the initial primal tolerance ALMParams::ε₀ be reduced.

unsigned max_num_retries = 20

How many times can the penalty update factor ALMParams::Δ and the primal tolerance factor ALMParams::ρ be reduced.

unsigned max_total_num_retries = 40

Combined limit for ALMParams::max_num_initial_retries and ALMParams::max_num_retries.

unsigned print_interval = 0

When to print progress.

If set to zero, nothing will be printed. If set to N != 0, progress is printed every N iterations.

bool preconditioning = true

Apply preconditioning to the problem, based on the gradients in the starting point.

bool single_penalty_factor = false

Use one penalty factor for all m constraints.

template<class InnerSolverT = PANOCSolver<>>
class ALMSolver
#include <panoc-alm/decl/alm.hpp>

Augmented Lagrangian Method solver.

Public Types

using Params = ALMParams
using InnerSolver = InnerSolverT

Public Functions

inline ALMSolver(Params params, InnerSolver &&inner_solver)
inline ALMSolver(Params params, const InnerSolver &inner_solver)
Stats operator()(const Problem &problem, rvec y, rvec x)
inline std::string get_name() const
inline void stop()

Abort the computation and return the result so far.

Can be called from other threads or signal handlers.

inline const Params &get_params() const

Public Members

InnerSolver inner_solver
struct Stats
#include <panoc-alm/decl/alm.hpp>

Public Members

unsigned outer_iterations = 0
std::chrono::microseconds elapsed_time
unsigned initial_penalty_reduced = 0
unsigned penalty_reduced = 0
unsigned inner_convergence_failures = 0
real_t ε = inf
real_t δ = inf
real_t norm_penalty = 0
SolverStatus status = SolverStatus::Unknown
InnerStatsAccumulator<typename InnerSolver::Stats> inner
template<class InnerSolverT = PANOCSolverFull<>>
class ALMSolverFull
#include <panoc-alm/decl/alm.hpp>

Public Types

using Params = ALMParams
using InnerSolver = InnerSolverT

Public Functions

inline ALMSolverFull(Params params, InnerSolver &&inner_solver)
inline ALMSolverFull(Params params, const InnerSolver &inner_solver)
Stats operator()(const ProblemFull &problem, rvec y, rvec x)
inline std::string get_name() const
inline void stop()

Abort the computation and return the result so far.

Can be called from other threads or signal handlers.

inline const Params &get_params() const

Public Members

InnerSolver inner_solver
struct Stats
#include <panoc-alm/decl/alm.hpp>

Public Functions

inline Stats(int m2)

Public Members

unsigned outer_iterations = 0
std::chrono::microseconds elapsed_time
unsigned initial_penalty_reduced = 0
unsigned penalty_reduced = 0
unsigned inner_convergence_failures = 0
real_t ε = inf
real_t δ₁ = inf
real_t δ₂ = inf
real_t norm_penalty₁ = 0
real_t norm_penalty₂ = 0
vec penalty₂
SolverStatus status = SolverStatus::Unknown
InnerStatsAccumulator<typename InnerSolver::Stats> inner
class AndersonAccel
#include <panoc-alm/inner/directions/anderson-acceleration.hpp>

Anderson Acceleration.

Public Functions

inline void resize(size_t n, size_t m)
inline void initialize(crvec g₀, crvec r₀)
inline void compute(crvec gₖ, crvec rₖ, rvec xₖ_aa)
inline std::string get_name() const
inline void reinitialize(real_t γₖ, real_t old_γₖ)
inline void reset()
class AtomicStopSignal
#include <panoc-alm/util/atomic_stop_signal.hpp>

Public Functions

AtomicStopSignal() = default
inline AtomicStopSignal(const AtomicStopSignal&)
AtomicStopSignal &operator=(const AtomicStopSignal&) = delete
inline AtomicStopSignal(AtomicStopSignal&&)
inline AtomicStopSignal &operator=(AtomicStopSignal&&)
inline void stop()
inline bool stop_requested() const
struct Box
#include <panoc-alm/util/box.hpp>

Public Members

vec upperbound
vec lowerbound
template<class IndexT = size_t>
struct CircularIndexIterator
#include <panoc-alm/util/ringbuffer.hpp>

Public Types

using Index = IndexT
using Indices = CircularIndices<Index>
using value_type = Indices
using reference = value_type
using difference_type = std::ptrdiff_t
using pointer = void
using iterator_category = std::input_iterator_tag

Public Functions

inline CircularIndexIterator()
inline CircularIndexIterator(Indices i, Index max)
inline reference operator*() const
inline CircularIndexIterator &operator++()
inline CircularIndexIterator &operator--()
inline CircularIndexIterator operator++(int)
inline CircularIndexIterator operator--(int)

Public Members

Indices i
Index max
template<class IndexT>
bool operator==(CircularIndexIterator<IndexT> a, CircularIndexIterator<IndexT> b)

Note

Only valid for two indices in the same range.

template<class IndexT>
bool operator!=(CircularIndexIterator<IndexT> a, CircularIndexIterator<IndexT> b)

Note

Only valid for two indices in the same range.

template<class IndexT = size_t>
struct CircularIndices
#include <panoc-alm/util/ringbuffer.hpp>

Public Types

using Index = IndexT

Public Functions

inline CircularIndices(Index zerobased, Index circular)

Public Members

Index zerobased
Index circular
template<class IndexT>
bool operator==(CircularIndices<IndexT> a, CircularIndices<IndexT> b)

Note

Only valid for two indices in the same range.

template<class IndexT>
bool operator!=(CircularIndices<IndexT> a, CircularIndices<IndexT> b)

Note

Only valid for two indices in the same range.

template<class IndexT>
class CircularRange
#include <panoc-alm/util/ringbuffer.hpp>

Public Types

using Index = IndexT
using Indices = CircularIndices<Index>
using const_iterator = CircularIndexIterator<Index>
using iterator = const_iterator
using const_reverse_iterator = ReverseCircularIndexIterator<Index>
using reverse_iterator = const_reverse_iterator

Public Functions

inline CircularRange(Index size, Index idx1, Index idx2, Index max)
inline iterator begin() const
inline iterator end() const
inline const_iterator cbegin() const
inline const_iterator cend() const
inline reverse_iterator rbegin() const
inline reverse_iterator rend() const
inline const_reverse_iterator crbegin() const
inline const_reverse_iterator crend() const
struct EvalCounter
#include <panoc-alm/util/problem.hpp>

Public Functions

inline void reset()

Public Members

unsigned f = 0
unsigned grad_f = 0
unsigned g = 0
unsigned grad_g_prod = 0
unsigned grad_gi = 0
unsigned hess_L_prod = 0
unsigned hess_L = 0
struct EvalCounterFull
#include <panoc-alm/util/problem.hpp>

Public Functions

inline void reset()

Public Members

unsigned f = 0
unsigned grad_f = 0
unsigned g = 0
unsigned grad_g_prod = 0
unsigned grad_gi = 0
unsigned hess_L_prod = 0
unsigned hess_L = 0
struct GAAPGAParams
#include <panoc-alm/inner/guarded-aa-pga.hpp>

Public Members

LipschitzEstimateParams Lipschitz

Parameters related to the Lipschitz constant estimate and step size.

unsigned limitedqr_mem = 10

Length of the history to keep in the limited-memory QR algorithm.

unsigned max_iter = 100

Maximum number of inner iterations.

std::chrono::microseconds max_time = std::chrono::minutes(5)

Maximum duration.

real_t L_min = 1e-5

Minimum Lipschitz constant estimate.

real_t L_max = 1e9

Maximum Lipschitz constant estimate.

PANOCStopCrit stop_crit = PANOCStopCrit::ApproxKKT

What stopping criterion to use.

unsigned print_interval = 0

When to print progress.

If set to zero, nothing will be printed. If set to N != 0, progress is printed every N iterations.

real_t quadratic_upperbound_tolerance_factor = 10 * std::numeric_limits<real_t>::epsilon()
unsigned max_no_progress = 10

Maximum number of iterations without any progress before giving up.

bool full_flush_on_γ_change = true
struct GAAPGAProgressInfo
#include <panoc-alm/inner/guarded-aa-pga.hpp>

Public Members

unsigned k
crvec x
crvec p
real_t norm_sq_p
crvec x_hat
real_t ψ
crvec grad_ψ
real_t ψ_hat
crvec grad_ψ_hat
real_t L
real_t γ
real_t ε
crvec Σ
crvec y
const Problem &problem
const GAAPGAParams &params
class GAAPGASolver
#include <panoc-alm/inner/guarded-aa-pga.hpp>

Guarded Anderson Accelerated Proximal Gradient Algorithm.

Vien V. Mai and Mikael Johansson, Anderson Acceleration of Proximal Gradient Methods. https://arxiv.org/abs/1910.08590v2

Public Types

using Params = GAAPGAParams
using ProgressInfo = GAAPGAProgressInfo

Public Functions

inline GAAPGASolver(const Params &params)
inline Stats operator()(const Problem &problem, crvec Σ, real_t ε, bool always_overwrite_results, rvec x, rvec λ, rvec err_z)
inline GAAPGASolver &set_progress_callback(std::function<void(const ProgressInfo&)> cb)
inline std::string get_name() const
inline void stop()
inline const Params &get_params() const
struct Stats
#include <panoc-alm/inner/guarded-aa-pga.hpp>

Public Members

SolverStatus status = SolverStatus::Unknown
real_t ε = inf
std::chrono::microseconds elapsed_time
unsigned iterations = 0
unsigned accelerated_steps_accepted = 0
template<class InnerSolverStats>
struct InnerStatsAccumulator
template<> Stats >
#include <panoc-alm/inner/guarded-aa-pga.hpp>

Public Members

std::chrono::microseconds elapsed_time
unsigned iterations = 0
unsigned accelerated_steps_accepted = 0
template<>
struct InnerStatsAccumulator<LBFGSBStats>
#include <panoc-alm/inner/lbfgspp.hpp>

Public Members

std::chrono::microseconds elapsed_time
unsigned iterations = 0
template<>
struct InnerStatsAccumulator<PANOCStats>
#include <panoc-alm/inner/decl/panoc.hpp>

Public Members

std::chrono::microseconds elapsed_time
unsigned iterations = 0
unsigned linesearch_failures = 0
unsigned lbfgs_failures = 0
unsigned lbfgs_rejected = 0
unsigned τ_1_accepted = 0
unsigned count_τ = 0
real_t sum_τ = 0
template<> Stats >
#include <panoc-alm/inner/pga.hpp>

Public Members

std::chrono::microseconds elapsed_time
unsigned iterations = 0
template<> Stats >
#include <panoc-alm/inner/decl/second-order-panoc.hpp>

Public Members

std::chrono::microseconds elapsed_time
unsigned iterations = 0
unsigned newton_failures = 0
unsigned linesearch_failures = 0
unsigned τ_1_accepted = 0
unsigned count_τ = 0
real_t sum_τ = 0
template<> Stats >
#include <panoc-alm/inner/decl/structured-panoc-lbfgs.hpp>

Public Members

std::chrono::microseconds elapsed_time
unsigned iterations = 0
unsigned linesearch_failures = 0
unsigned lbfgs_failures = 0
unsigned lbfgs_rejected = 0
unsigned τ_1_accepted = 0
unsigned count_τ = 0
real_t sum_τ = 0
class LBFGS
#include <panoc-alm/inner/directions/decl/lbfgs.hpp>

Limited memory Broyden–Fletcher–Goldfarb–Shanno (L-BFGS) algorithm.

Public Types

enum Sign

The sign of the vectors \( p \) passed to the LBFGS::update method.

Values:

enumerator Positive

\( p \sim \nabla \psi(x) \)

enumerator Negative

\( p \sim -\nabla \psi(x) \)

using Params = LBFGSParams

Public Functions

inline LBFGS(Params params)
inline LBFGS(Params params, size_t n)
inline 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ₖ₊₁.

template<class Vec>
bool apply(Vec &&q, real_t γ)

Apply the inverse Hessian approximation to the given vector q.

template<class Vec, class IndexVec>
bool apply(Vec &&q, real_t γ, const IndexVec &J)

Apply the inverse Hessian approximation to the given vector q, applying only the columns and rows of the Hessian in the index set J.

inline void reset()

Throw away the approximation and all previous vectors s and y.

inline void resize(size_t n)

Re-allocate storage for a problem with a different size.

Causes a reset.

inline void scale_y(real_t factor)

Scale the stored y vectors by the given factor.

inline std::string get_name() const
inline const Params &get_params() const
inline size_t n() const

Get the size of the s and y vectors in the buffer.

inline size_t history() const

Get the number of previous vectors s and y stored in the buffer.

inline size_t succ(size_t i) const

Get the next index in the circular buffer of previous s and y vectors.

inline auto s(size_t i)
inline auto s(size_t i) const
inline auto y(size_t i)
inline auto y(size_t i) const
inline real_t &ρ(size_t i)
inline const real_t &ρ(size_t i) const
inline real_t &α(size_t i)
inline const real_t &α(size_t i) const

Public Static Functions

static inline bool update_valid(LBFGSParams params, real_t yᵀs, real_t sᵀs, real_t pᵀp)

Check if the new vectors s and y allow for a valid BFGS update that preserves the positive definiteness of the Hessian approximation.

template<template<class> class LineSearchT = LBFGSpp::LineSearchBacktracking>
class LBFGSBSolver
#include <panoc-alm/inner/lbfgspp.hpp>

Box-constrained LBFGS solver for ALM.

Public Types

using Params = LBFGSpp::LBFGSBParam<real_t>
using LineSearch = LineSearchT<real_t>
using Stats = LBFGSBStats

Public Functions

inline LBFGSBSolver(Params params)
inline Stats operator()(const Problem &problem, crvec Σ, real_t ε, bool always_overwrite_results, rvec x, rvec y, rvec err_z)
inline std::string get_name() const
inline const Params &get_params() const
inline void stop()
struct LBFGSBStats
#include <panoc-alm/inner/lbfgspp.hpp>

Public Members

SolverStatus status = SolverStatus::Unknown
real_t ε = inf
std::chrono::microseconds elapsed_time
unsigned iterations = 0
struct LBFGSParams
#include <panoc-alm/inner/directions/decl/lbfgs.hpp>

Parameters for the LBFGS and SpecializedLBFGS classes.

Public Members

unsigned memory = 10

Length of the history to keep.

struct pa::LBFGSParams::[anonymous] cbfgs

Parameters in the cautious BFGS update condition.

\[ \frac{y^\top s}{s^\top s} \ge \epsilon \| g \|^\alpha \]
See

https://epubs.siam.org/doi/10.1137/S1052623499354242

bool rescale_when_γ_changes = false
LBFGSParams.cbfgs

Public Members

real_t α
real_t ϵ
template<template<class> class LineSearchT = LBFGSpp::LineSearchBacktracking>
class LBFGSSolver
#include <panoc-alm/inner/lbfgspp.hpp>

Unconstrained LBFGS solver for ALM.

Public Types

using Params = LBFGSpp::LBFGSParam<real_t>
using LineSearch = LineSearchT<real_t>

Public Functions

inline LBFGSSolver(Params params)
inline Stats operator()(const Problem &problem, crvec Σ, real_t ε, bool always_overwrite_results, rvec x, rvec y, rvec err_z)
inline std::string get_name() const
inline const Params &get_params() const
inline void stop()
struct Stats
#include <panoc-alm/inner/lbfgspp.hpp>

Public Members

unsigned iterations = 0
real_t ε = inf
std::chrono::microseconds elapsed_time
SolverStatus status = SolverStatus::Unknown
unsigned linesearch_failures = 0
unsigned lbfgs_failures = 0
unsigned lbfgs_rejected = 0
class LimitedMemoryQR
#include <panoc-alm/inner/detail/limited-memory-qr.hpp>

Incremental QR factorization using modified Gram-Schmidt with reorthogonalization.

Computes A = QR while allowing efficient removal of the first column of A or adding new columns at the end of A.

Public Types

template<class Derived>
using solve_ret_t = std::conditional_t<Eigen::internal::traits<Derived>::ColsAtCompileTime == 1, vec, mat>

Public Functions

LimitedMemoryQR() = default
inline LimitedMemoryQR(size_t n, size_t m)

The maximum dimensions of Q are n×m and the maximum dimensions of R are m×m.

Parameters
  • n – The size of the vectors, the number of rows of A.

  • m – The maximum number of columns of A.

inline size_t n() const
inline size_t m() const
template<class VecV>
inline void add_column(const VecV &v)

Add the given column to the right.

inline void remove_column()

Remove the leftmost column.

template<class VecB, class VecX>
inline void solve_col(const VecB &b, VecX &x) const

Solve the least squares problem Ax = b.

template<class MatB, class MatX>
inline void solve(const MatB &B, MatX &X) const

Solve the least squares problem AX = B.

template<class Derived>
inline solve_ret_t<Derived> solve(const Eigen::DenseBase<Derived> &B)

Solve the least squares problem AX = B.

inline const mat &get_raw_Q() const

Get the full, raw storage for the orthogonal matrix Q.

inline const mat &get_raw_R() const

Get the full, raw storage for the upper triangular matrix R.

The columns of this matrix are permuted because it’s stored as a circular buffer for efficiently appending columns to the end and popping columns from the front.

inline mat get_full_R() const

Get the full storage for the upper triangular matrix R but with the columns in the correct order.

Note

Meant for tests only, creates a permuted copy.

inline mat get_R() const

Get the matrix R such that Q times R is the original matrix.

Note

Meant for tests only, creates a permuted copy.

inline mat get_Q() const

Get the matrix Q such that Q times R is the original matrix.

Note

Meant for tests only, creates a copy.

inline void scale_R(real_t scal)

Multiply the matrix R by a scalar.

inline size_t get_reorth_count() const

Get the number of MGS reorthogonalizations.

inline void clear_reorth_count()

Reset the number of MGS reorthogonalizations.

inline void reset()

Reset all indices, clearing the Q and R matrices.

inline void resize(size_t n, size_t m)

Re-allocate storage for a problem with a different size.

Causes a reset.

inline size_t num_columns() const

Get the number of columns that are currently stored.

inline size_t ring_head() const

Get the head index of the circular buffer (points to the oldest element).

inline size_t ring_tail() const

Get the tail index of the circular buffer (points to one past the most recent element).

inline size_t ring_next(size_t i) const

Get the next index in the circular buffer.

inline size_t ring_prev(size_t i) const

Get the previous index in the circular buffer.

inline CircularRange<size_t> ring_iter() const

Get iterators in the circular buffer.

inline ReverseCircularRange<size_t> ring_reverse_iter() const

Get reverse iterators in the circular buffer.

struct LipschitzEstimateParams
#include <panoc-alm/util/lipschitz.hpp>

Public Members

real_t L₀ = 0

Initial estimate of the Lipschitz constant of ∇ψ(x)

real_t ε = 1e-6

Relative step size for initial finite difference Lipschitz estimate.

real_t δ = 1e-12

Minimum step size for initial finite difference Lipschitz estimate.

real_t Lγ_factor = 0.95

Factor that relates step size γ and Lipschitz constant.

template<class DirectionProviderT>
struct PANOCDirection
#include <panoc-alm/inner/directions/decl/panoc-direction-update.hpp>

Public Static Functions

static void initialize(DirectionProviderT &dp, crvec x₀, crvec x̂₀, crvec p₀, crvec grad₀) = delete
static bool update(DirectionProviderT &dp, crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁, crvec gradₖ₊₁, const Box &C, real_t γₖ₊₁) = delete
static bool apply(DirectionProviderT &dp, crvec xₖ, crvec x̂ₖ, crvec pₖ, real_t γ, rvec qₖ) = delete

Apply the direction estimation in the current point.

Parameters
  • dp[in] Direction provider (e.g. LBFGS, Anderson Acceleration).

  • xₖ[in] Current iterate.

  • x̂ₖ[in] Result of proximal gradient step in current iterate.

  • pₖ[in] Proximal gradient step between x̂ₖ and xₖ.

  • γ[in] H₀ = γI for L-BFGS

  • qₖ[out] Resulting step.

static void changed_γ(DirectionProviderT &dp, real_t γₖ, real_t old_γₖ) = delete
template<>
struct PANOCDirection<AndersonAccel>
#include <panoc-alm/inner/directions/anderson-acceleration.hpp>

Public Static Functions

static inline void initialize(AndersonAccel &aa, crvec x₀, crvec x̂₀, crvec p₀, crvec grad₀)
static inline bool update(AndersonAccel &aa, crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁, crvec grad_new, const Box &C, real_t γ_new)
static inline bool apply(AndersonAccel &aa, crvec xₖ, crvec x̂ₖ, crvec pₖ, real_t γ, rvec qₖ)
static inline void changed_γ(AndersonAccel &aa, real_t γₖ, real_t old_γₖ)
template<>
struct PANOCDirection<LBFGS>
#include <panoc-alm/inner/directions/decl/lbfgs.hpp>

Public Functions

inline PANOCDirection(const LBFGSParams &params)
inline PANOCDirection(const LBFGS &lbfgs)
inline PANOCDirection(LBFGS &&lbfgs)
inline void initialize(crvec x₀, crvec x̂₀, crvec p₀, crvec grad₀)
inline bool update(crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁, crvec grad_new, const Box &C, real_t γ_new)
inline bool apply(crvec xₖ, crvec x̂ₖ, crvec pₖ, real_t γ, rvec qₖ)
inline void changed_γ(real_t γₖ, real_t old_γₖ)
inline void reset()
inline std::string get_name() const
inline LBFGSParams get_params() const

Public Members

LBFGS lbfgs
template<>
struct PANOCDirection<SpecializedLBFGS>
#include <panoc-alm/inner/directions/specialized-lbfgs.hpp>

Public Static Functions

static inline void initialize(SpecializedLBFGS &lbfgs, crvec x₀, crvec x̂₀, crvec p₀, crvec grad₀)
static inline bool update(SpecializedLBFGS &lbfgs, crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁, crvec gradₖ₊₁, const Box &C, real_t γ)
static inline bool apply(SpecializedLBFGS &lbfgs, crvec xₖ, crvec x̂ₖ, crvec pₖ, real_t γ, rvec qₖ)
static inline void changed_γ(SpecializedLBFGS &lbfgs, real_t γₖ, real_t old_γₖ)
struct PANOCFullProgressInfo
#include <panoc-alm/inner/decl/panoc.hpp>

Public Members

unsigned k
crvec x
crvec p
real_t norm_sq_p
crvec x_hat
real_t φγ
real_t ψ
crvec grad_ψ
real_t ψ_hat
crvec grad_ψ_hat
real_t L
real_t γ
real_t τ
real_t ε
crvec Σ1
crvec Σ2
crvec y
const ProblemFull &problem
const PANOCParams &params
struct PANOCParams
#include <panoc-alm/inner/decl/panoc.hpp>

Tuning parameters for the PANOC algorithm.

Public Members

LipschitzEstimateParams Lipschitz

Parameters related to the Lipschitz constant estimate and step size.

unsigned max_iter = 100

Maximum number of inner PANOC iterations.

std::chrono::microseconds max_time = std::chrono::minutes(5)

Maximum duration.

real_t τ_min = 1. / 256

Minimum weight factor between Newton step and projected gradient step.

real_t L_min = 1e-5

Minimum Lipschitz constant estimate.

real_t L_max = 1e9

Maximum Lipschitz constant estimate.

PANOCStopCrit stop_crit = PANOCStopCrit::ApproxKKT

What stopping criterion to use.

unsigned max_no_progress = 10

Maximum number of iterations without any progress before giving up.

unsigned print_interval = 0

When to print progress.

If set to zero, nothing will be printed. If set to N != 0, progress is printed every N iterations.

real_t quadratic_upperbound_tolerance_factor = 10 * std::numeric_limits<real_t>::epsilon()
bool update_lipschitz_in_linesearch = true
bool alternative_linesearch_cond = false
LBFGSStepSize lbfgs_stepsize = LBFGSStepSize::BasedOnCurvature
struct PANOCProgressInfo
#include <panoc-alm/inner/decl/panoc.hpp>

Public Members

unsigned k
crvec x
crvec p
real_t norm_sq_p
crvec x_hat
real_t φγ
real_t ψ
crvec grad_ψ
real_t ψ_hat
crvec grad_ψ_hat
real_t L
real_t γ
real_t τ
real_t ε
crvec Σ
crvec y
const Problem &problem
const PANOCParams &params
template<class DirectionProviderT>
class PANOCSolver
#include <panoc-alm/inner/decl/panoc.hpp>

PANOC solver for ALM.

Public Types

using Params = PANOCParams
using DirectionProvider = DirectionProviderT
using Stats = PANOCStats
using ProgressInfo = PANOCProgressInfo

Public Functions

inline PANOCSolver(Params params, PANOCDirection<DirectionProvider> &&direction_provider)
inline PANOCSolver(Params params, const PANOCDirection<DirectionProvider> &direction_provider)
Stats operator()(const Problem &problem, crvec Σ, real_t ε, bool always_overwrite_results, rvec x, rvec y, rvec err_z, std::chrono::microseconds time_remaining = std::chrono::microseconds(0))
Parameters
  • problem – [in] Problem description

  • Σ – [in] Constraint weights \( \Sigma \)

  • ε – [in] Tolerance \( \varepsilon \)

  • always_overwrite_results – [in] Overwrite x, y and err_z even if not converged

  • x – [inout] Decision variable \( x \)

  • y – [inout] Lagrange multipliers \( y \)

  • err_z – [out] Slack variable error \( g(x) - z \)

  • time_remaining – [in] Time remaining

inline PANOCSolver &set_progress_callback(std::function<void(const ProgressInfo&)> cb)
std::string get_name() const
inline void stop()
inline const Params &get_params() const

Public Members

PANOCDirection<DirectionProvider> direction_provider
template<class DirectionProviderT>
class PANOCSolverFull
#include <panoc-alm/inner/decl/panoc.hpp>

Public Types

using Params = PANOCParams
using DirectionProvider = DirectionProviderT
using Stats = PANOCStats
using ProgressInfo = PANOCFullProgressInfo

Public Functions

inline PANOCSolverFull(Params params, PANOCDirection<DirectionProvider> &&direction_provider)
inline PANOCSolverFull(Params params, const PANOCDirection<DirectionProvider> &direction_provider)
Stats operator()(const ProblemFull &problem, crvec Σ1, crvec Σ2, real_t ε, bool always_overwrite_results, rvec x, rvec y, rvec err_z1, rvec err_z2, std::chrono::microseconds time_remaining = std::chrono::microseconds(0))
Parameters
  • problem – [in] Problem description

  • Σ1 – [in] Constraint weights \( \Sigma \)

  • Σ2 – [in] Constraint weights \( \Sigma \)

  • ε – [in] Tolerance \( \varepsilon \)

  • always_overwrite_results – [in] Overwrite x, y and err_z even if not converged

  • x – [inout] Decision variable \( x \)

  • y – [inout] Lagrange multipliers \( y \)

  • err_z1 – [out] Slack variable error \( g(x) - z \)

  • err_z2 – [out] Slack variable error \( g(x) - z \)

  • time_remaining – [in] Time remaining

inline PANOCSolverFull &set_progress_callback(std::function<void(const ProgressInfo&)> cb)
std::string get_name() const
inline void stop()
inline const Params &get_params() const

Public Members

PANOCDirection<DirectionProvider> direction_provider
struct PANOCStats
#include <panoc-alm/inner/decl/panoc.hpp>

Public Members

SolverStatus status = SolverStatus::Unknown
real_t ε = inf
std::chrono::microseconds elapsed_time
unsigned iterations = 0
unsigned linesearch_failures = 0
unsigned lbfgs_failures = 0
unsigned lbfgs_rejected = 0
unsigned τ_1_accepted = 0
unsigned count_τ = 0
real_t sum_τ = 0
struct PGAParams
#include <panoc-alm/inner/pga.hpp>

Public Members

LipschitzEstimateParams Lipschitz

Parameters related to the Lipschitz constant estimate and step size.

unsigned max_iter = 100

Maximum number of inner iterations.

std::chrono::microseconds max_time = std::chrono::minutes(5)

Maximum duration.

real_t L_min = 1e-5

Minimum Lipschitz constant estimate.

real_t L_max = 1e9

Maximum Lipschitz constant estimate.

PANOCStopCrit stop_crit = PANOCStopCrit::ApproxKKT

What stop criterion to use.

unsigned print_interval = 0

When to print progress.

If set to zero, nothing will be printed. If set to N != 0, progress is printed every N iterations.

real_t quadratic_upperbound_tolerance_factor = 10 * std::numeric_limits<real_t>::epsilon()
struct PGAProgressInfo
#include <panoc-alm/inner/pga.hpp>

Public Members

unsigned k
crvec x
crvec p
real_t norm_sq_p
crvec x_hat
real_t ψ
crvec grad_ψ
real_t ψ_hat
crvec grad_ψ_hat
real_t L
real_t γ
real_t ε
crvec Σ
crvec y
const Problem &problem
const PGAParams &params
class PGASolver
#include <panoc-alm/inner/pga.hpp>

Standard Proximal Gradient Algorithm without any bells and whistles.

Public Types

using Params = PGAParams
using ProgressInfo = PGAProgressInfo

Public Functions

inline PGASolver(const Params &params)
inline Stats operator()(const Problem &problem, crvec Σ, real_t ε, bool always_overwrite_results, rvec x, rvec λ, rvec err_z)
inline PGASolver &set_progress_callback(std::function<void(const ProgressInfo&)> cb)
inline std::string get_name() const
inline void stop()
inline const Params &get_params() const
struct Stats
#include <panoc-alm/inner/pga.hpp>

Public Members

SolverStatus status = SolverStatus::Unknown
real_t ε = inf
std::chrono::microseconds elapsed_time
unsigned iterations = 0
class Problem
#include <panoc-alm/util/problem.hpp>

Problem description for minimization problems.

\[\begin{split} \begin{aligned} & \underset{x}{\text{minimize}} & & f(x) &&&& f : \mathbb{R}^n \rightarrow \mathbb{R} \\ & \text{subject to} & & \underline{x} \le x \le \overline{x} \\ &&& \underline{z} \le g(x) \le \overline{z} &&&& g : \mathbb{R}^n \rightarrow \mathbb{R}^m \end{aligned} \end{split}\]

Subclassed by ProblemOnlyD, ProblemWithCounters, ProblemWithParam

Public Types

using f_sig = real_t(crvec x)

Signature of the function that evaluates the cost \( f(x) \).

Parameters

x[in] Decision variable \( x \in \mathbb{R}^n \)

using grad_f_sig = void(crvec x, rvec grad_fx)

Signature of the function that evaluates the gradient of the cost function \( \nabla f(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • grad_fx[out] Gradient of cost function \( \nabla f(x) \in \mathbb{R}^n \)

using g_sig = void(crvec x, rvec gx)

Signature of the function that evaluates the constraints \( g(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • gx[out] Value of the constraints \( g(x) \in \mathbb{R}^m \)

using grad_g_prod_sig = void(crvec x, crvec y, rvec grad_gxy)

Signature of the function that evaluates the gradient of the constraints times a vector \( \nabla g(x)\ y \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Vector \( y \in \mathbb{R}^m \) to multiply the gradient by

  • grad_gxy[out] Gradient of the constraints \( \nabla g(x)\ y \in \mathbb{R}^n \)

using grad_gi_sig = void(crvec x, unsigned i, rvec grad_gi)

Signature of the function that evaluates the gradient of one specific constraints \( \nabla g_i(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • i[in] Which constraint \( 0 \le i \lt m \)

  • grad_gi[out] Gradient of the constraint \( \nabla g_i(x) \mathbb{R}^n \)

using hess_L_prod_sig = void(crvec x, crvec y, crvec v, rvec Hv)

Signature of the function that evaluates the Hessian of the Lagrangian multiplied by a vector \( \nabla_{xx}^2L(x, y)\ v \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • v[in] Vector to multiply by \( v \in \mathbb{R}^n \)

  • Hv[out] Hessian-vector product \( \nabla_{xx}^2 L(x, y)\ v \in \mathbb{R}^{n} \)

using hess_L_sig = void(crvec x, crvec y, rmat H)

Signature of the function that evaluates the Hessian of the Lagrangian \( \nabla_{xx}^2L(x, y) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • H[out] Hessian \( \nabla_{xx}^2 L(x, y) \in \mathbb{R}^{n\times n} \)

Public Functions

Problem() = default
inline Problem(unsigned int n, unsigned int m)
inline Problem(unsigned n, unsigned int m, Box C, Box D, std::function<f_sig> f, std::function<grad_f_sig> grad_f, std::function<g_sig> g, std::function<grad_g_prod_sig> grad_g_prod, std::function<grad_gi_sig> grad_gi, std::function<hess_L_prod_sig> hess_L_prod, std::function<hess_L_sig> hess_L)

Public Members

unsigned int n

Number of decision variables, dimension of x.

unsigned int m

Number of constraints, dimension of g(x) and z.

Box C

Constraints of the decision variables, \( x \in C \).

Box D

Other constraints, \( g(x) \in D \).

std::function<f_sig> f

Cost function \( f(x) \).

std::function<grad_f_sig> grad_f

Gradient of the cost function \( \nabla f(x) \).

std::function<g_sig> g

Constraint function \( g(x) \).

std::function<grad_g_prod_sig> grad_g_prod

Gradient of the constraint function times vector \( \nabla g(x)\ y \).

std::function<grad_gi_sig> grad_gi

Gradient of a specific constraint \( \nabla g_i(x) \).

std::function<hess_L_prod_sig> hess_L_prod

Hessian of the Lagrangian function times vector \( \nabla_{xx}^2 L(x, y)\ v \).

std::function<hess_L_sig> hess_L

Hessian of the Lagrangian function \( \nabla_{xx}^2 L(x, y) \).

class ProblemFull
#include <panoc-alm/util/problem.hpp>

Problem description for minimization problems.

\[\begin{split} \begin{aligned} & \underset{x}{\text{minimize}} & & f(x) &&&& f : \mathbb{R}^n \rightarrow \mathbb{R} \\ & \text{subject to} & & \underline{x} \le x \le \overline{x} \\ &&& \underline{z} \le g(x) \le \overline{z} &&&& g : \mathbb{R}^n \rightarrow \mathbb{R}^m \end{aligned} \end{split}\]

Subclassed by ProblemFullWithParam

Public Types

using f_sig = real_t(crvec x)

Signature of the function that evaluates the cost \( f(x) \).

Parameters

x[in] Decision variable \( x \in \mathbb{R}^n \)

using grad_f_sig = void(crvec x, rvec grad_fx)

Signature of the function that evaluates the gradient of the cost function \( \nabla f(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • grad_fx[out] Gradient of cost function \( \nabla f(x) \in \mathbb{R}^n \)

using g1_sig = void(crvec x, rvec g1x)

Signature of the function that evaluates the ALM constraints \( g1(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • g1x[out] Value of the constraints \( g1(x) \in \mathbb{R}^m \)

using grad_g1_prod_sig = void(crvec x, crvec y, rvec grad_g1xy)

Signature of the function that evaluates the gradient of the ALM constraints times a vector \( \nabla g1(x)\ y \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Vector \( y \in \mathbb{R}^m1 \) to multiply the gradient by

  • grad_g1xy[out] Gradient of the constraints \( \nabla g1(x)\ y \in \mathbb{R}^n \)

using grad_g1i_sig = void(crvec x, unsigned i, rvec grad_g1i)

Signature of the function that evaluates the gradient of one specific ALM constraint \( \nabla g1_i(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • i[in] Which constraint \( 0 \le i \lt m1 \)

  • grad_g1i[out] Gradient of the constraint \( \nabla g1_i(x) \mathbb{R}^n \)

using g2_sig = void(crvec x, rvec g2x)

Signature of the function that evaluates the quadratic penalty constraints \( g2(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • g2x[out] Value of the constraints \( g2(x) \in \mathbb{R}^m \)

using grad_g2_prod_sig = void(crvec x, crvec y, rvec grad_g2xy)

Signature of the function that evaluates the gradient of the quadratic penalty constraints times a vector \( \nabla g2(x)\ y \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Vector \( y \in \mathbb{R}^m \) to multiply the gradient by

  • grad_g2xy[out] Gradient of the constraints \( \nabla g2(x)\ y \in \mathbb{R}^n \)

using grad_g2i_sig = void(crvec x, unsigned i, rvec grad_g2i)

Signature of the function that evaluates the gradient of one specific quadratic penalty constraints \( \nabla g2_i(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • i[in] Which constraint \( 0 \le i \lt m2 \)

  • grad_g2i[out] Gradient of the constraint \( \nabla g2_i(x) \mathbb{R}^n \)

using hess_L_prod_sig = void(crvec x, crvec y, crvec v, rvec Hv)

Signature of the function that evaluates the Hessian of the Lagrangian multiplied by a vector \( \nabla_{xx}^2L(x, y)\ v \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • v[in] Vector to multiply by \( v \in \mathbb{R}^n \)

  • Hv[out] Hessian-vector product \( \nabla_{xx}^2 L(x, y)\ v \in \mathbb{R}^{n} \)

using hess_L_sig = void(crvec x, crvec y, rmat H)

Signature of the function that evaluates the Hessian of the Lagrangian \( \nabla_{xx}^2L(x, y) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • H[out] Hessian \( \nabla_{xx}^2 L(x, y) \in \mathbb{R}^{n\times n} \)

Public Functions

ProblemFull() = default
inline ProblemFull(unsigned int n, unsigned int m1, unsigned int m2)
inline ProblemFull(unsigned n, unsigned int m1, unsigned int m2, Box C, Box D1, Box D2, std::function<f_sig> f, std::function<grad_f_sig> grad_f, std::function<g1_sig> g1, std::function<grad_g1_prod_sig> grad_g1_prod, std::function<grad_g1i_sig> grad_g1i, std::function<g2_sig> g2, std::function<grad_g2_prod_sig> grad_g2_prod, std::function<grad_g2i_sig> grad_g2i, std::function<hess_L_prod_sig> hess_L_prod, std::function<hess_L_sig> hess_L)

Public Members

unsigned int n

Number of decision variables, dimension of x.

unsigned int m1

Number of ALM constraints, dimension of g1(x) and z1.

unsigned int m2

Number of quadratic penalty constraints, dimension of g2(x) and z2.

Box C

Constraints of the decision variables, \( x \in C \).

Box D1

ALM constraints, \( g1(x) \in D1 \).

Box D2

Quadratic penalty constraints, \( g2(x) \in D2 \).

std::function<f_sig> f

Cost function \( f(x) \).

std::function<grad_f_sig> grad_f

Gradient of the cost function \( \nabla f(x) \).

std::function<g1_sig> g1

Constraint function \( g1(x) \).

std::function<grad_g1_prod_sig> grad_g1_prod

Gradient of the constraint function times vector \( \nabla g1(x)\ y \).

std::function<grad_g1i_sig> grad_g1i

Gradient of a specific constraint \( \nabla g1_i(x) \).

std::function<g2_sig> g2

Constraint function \( g2(x) \).

std::function<grad_g2_prod_sig> grad_g2_prod

Gradient of the constraint function times vector \( \nabla g2(x)\ y \).

std::function<grad_g2i_sig> grad_g2i

Gradient of a specific constraint \( \nabla g2_i(x) \).

std::function<hess_L_prod_sig> hess_L_prod

Hessian of the Lagrangian function times vector \( \nabla_{xx}^2 L(x, y)\ v \).

std::function<hess_L_sig> hess_L

Hessian of the Lagrangian function \( \nabla_{xx}^2 L(x, y) \).

class ProblemFullWithParam : public ProblemFull
#include <panoc-alm/util/problem.hpp>

Public Types

using f_sig = real_t(crvec x)

Signature of the function that evaluates the cost \( f(x) \).

Parameters

x[in] Decision variable \( x \in \mathbb{R}^n \)

using grad_f_sig = void(crvec x, rvec grad_fx)

Signature of the function that evaluates the gradient of the cost function \( \nabla f(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • grad_fx[out] Gradient of cost function \( \nabla f(x) \in \mathbb{R}^n \)

using g1_sig = void(crvec x, rvec g1x)

Signature of the function that evaluates the ALM constraints \( g1(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • g1x[out] Value of the constraints \( g1(x) \in \mathbb{R}^m \)

using grad_g1_prod_sig = void(crvec x, crvec y, rvec grad_g1xy)

Signature of the function that evaluates the gradient of the ALM constraints times a vector \( \nabla g1(x)\ y \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Vector \( y \in \mathbb{R}^m1 \) to multiply the gradient by

  • grad_g1xy[out] Gradient of the constraints \( \nabla g1(x)\ y \in \mathbb{R}^n \)

using grad_g1i_sig = void(crvec x, unsigned i, rvec grad_g1i)

Signature of the function that evaluates the gradient of one specific ALM constraint \( \nabla g1_i(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • i[in] Which constraint \( 0 \le i \lt m1 \)

  • grad_g1i[out] Gradient of the constraint \( \nabla g1_i(x) \mathbb{R}^n \)

using g2_sig = void(crvec x, rvec g2x)

Signature of the function that evaluates the quadratic penalty constraints \( g2(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • g2x[out] Value of the constraints \( g2(x) \in \mathbb{R}^m \)

using grad_g2_prod_sig = void(crvec x, crvec y, rvec grad_g2xy)

Signature of the function that evaluates the gradient of the quadratic penalty constraints times a vector \( \nabla g2(x)\ y \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Vector \( y \in \mathbb{R}^m \) to multiply the gradient by

  • grad_g2xy[out] Gradient of the constraints \( \nabla g2(x)\ y \in \mathbb{R}^n \)

using grad_g2i_sig = void(crvec x, unsigned i, rvec grad_g2i)

Signature of the function that evaluates the gradient of one specific quadratic penalty constraints \( \nabla g2_i(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • i[in] Which constraint \( 0 \le i \lt m2 \)

  • grad_g2i[out] Gradient of the constraint \( \nabla g2_i(x) \mathbb{R}^n \)

using hess_L_prod_sig = void(crvec x, crvec y, crvec v, rvec Hv)

Signature of the function that evaluates the Hessian of the Lagrangian multiplied by a vector \( \nabla_{xx}^2L(x, y)\ v \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • v[in] Vector to multiply by \( v \in \mathbb{R}^n \)

  • Hv[out] Hessian-vector product \( \nabla_{xx}^2 L(x, y)\ v \in \mathbb{R}^{n} \)

using hess_L_sig = void(crvec x, crvec y, rmat H)

Signature of the function that evaluates the Hessian of the Lagrangian \( \nabla_{xx}^2L(x, y) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • H[out] Hessian \( \nabla_{xx}^2 L(x, y) \in \mathbb{R}^{n\times n} \)

Public Functions

inline void set_param(pa::crvec p)
inline void set_param(pa::vec &&p)
inline pa::vec &get_param()
inline const pa::vec &get_param() const
inline std::shared_ptr<pa::vec> get_param_ptr() const

Public Members

unsigned int n

Number of decision variables, dimension of x.

unsigned int m1

Number of ALM constraints, dimension of g1(x) and z1.

unsigned int m2

Number of quadratic penalty constraints, dimension of g2(x) and z2.

Box C

Constraints of the decision variables, \( x \in C \).

Box D1

ALM constraints, \( g1(x) \in D1 \).

Box D2

Quadratic penalty constraints, \( g2(x) \in D2 \).

std::function<f_sig> f

Cost function \( f(x) \).

std::function<grad_f_sig> grad_f

Gradient of the cost function \( \nabla f(x) \).

std::function<g1_sig> g1

Constraint function \( g1(x) \).

std::function<grad_g1_prod_sig> grad_g1_prod

Gradient of the constraint function times vector \( \nabla g1(x)\ y \).

std::function<grad_g1i_sig> grad_g1i

Gradient of a specific constraint \( \nabla g1_i(x) \).

std::function<g2_sig> g2

Constraint function \( g2(x) \).

std::function<grad_g2_prod_sig> grad_g2_prod

Gradient of the constraint function times vector \( \nabla g2(x)\ y \).

std::function<grad_g2i_sig> grad_g2i

Gradient of a specific constraint \( \nabla g2_i(x) \).

std::function<hess_L_prod_sig> hess_L_prod

Hessian of the Lagrangian function times vector \( \nabla_{xx}^2 L(x, y)\ v \).

std::function<hess_L_sig> hess_L

Hessian of the Lagrangian function \( \nabla_{xx}^2 L(x, y) \).

class ProblemOnlyD : public Problem
#include <panoc-alm/util/problem.hpp>

Moves the state constraints in the set C to the set D, resulting in an unconstraint inner problem.

The new constraints function g becomes the concatenation of the original g function and the identity function. The new set D is the cartesian product of the original D × C.

Public Types

using f_sig = real_t(crvec x)

Signature of the function that evaluates the cost \( f(x) \).

Parameters

x[in] Decision variable \( x \in \mathbb{R}^n \)

using grad_f_sig = void(crvec x, rvec grad_fx)

Signature of the function that evaluates the gradient of the cost function \( \nabla f(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • grad_fx[out] Gradient of cost function \( \nabla f(x) \in \mathbb{R}^n \)

using g_sig = void(crvec x, rvec gx)

Signature of the function that evaluates the constraints \( g(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • gx[out] Value of the constraints \( g(x) \in \mathbb{R}^m \)

using grad_g_prod_sig = void(crvec x, crvec y, rvec grad_gxy)

Signature of the function that evaluates the gradient of the constraints times a vector \( \nabla g(x)\ y \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Vector \( y \in \mathbb{R}^m \) to multiply the gradient by

  • grad_gxy[out] Gradient of the constraints \( \nabla g(x)\ y \in \mathbb{R}^n \)

using grad_gi_sig = void(crvec x, unsigned i, rvec grad_gi)

Signature of the function that evaluates the gradient of one specific constraints \( \nabla g_i(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • i[in] Which constraint \( 0 \le i \lt m \)

  • grad_gi[out] Gradient of the constraint \( \nabla g_i(x) \mathbb{R}^n \)

using hess_L_prod_sig = void(crvec x, crvec y, crvec v, rvec Hv)

Signature of the function that evaluates the Hessian of the Lagrangian multiplied by a vector \( \nabla_{xx}^2L(x, y)\ v \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • v[in] Vector to multiply by \( v \in \mathbb{R}^n \)

  • Hv[out] Hessian-vector product \( \nabla_{xx}^2 L(x, y)\ v \in \mathbb{R}^{n} \)

using hess_L_sig = void(crvec x, crvec y, rmat H)

Signature of the function that evaluates the Hessian of the Lagrangian \( \nabla_{xx}^2L(x, y) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • H[out] Hessian \( \nabla_{xx}^2 L(x, y) \in \mathbb{R}^{n\times n} \)

Public Functions

inline ProblemOnlyD(Problem &&p)
inline ProblemOnlyD(const Problem &p)

Public Members

unsigned int n

Number of decision variables, dimension of x.

unsigned int m

Number of constraints, dimension of g(x) and z.

Box C

Constraints of the decision variables, \( x \in C \).

Box D

Other constraints, \( g(x) \in D \).

std::function<f_sig> f

Cost function \( f(x) \).

std::function<grad_f_sig> grad_f

Gradient of the cost function \( \nabla f(x) \).

std::function<g_sig> g

Constraint function \( g(x) \).

std::function<grad_g_prod_sig> grad_g_prod

Gradient of the constraint function times vector \( \nabla g(x)\ y \).

std::function<grad_gi_sig> grad_gi

Gradient of a specific constraint \( \nabla g_i(x) \).

std::function<hess_L_prod_sig> hess_L_prod

Hessian of the Lagrangian function times vector \( \nabla_{xx}^2 L(x, y)\ v \).

std::function<hess_L_sig> hess_L

Hessian of the Lagrangian function \( \nabla_{xx}^2 L(x, y) \).

class ProblemWithCounters : public Problem
#include <panoc-alm/util/problem.hpp>

Public Types

using f_sig = real_t(crvec x)

Signature of the function that evaluates the cost \( f(x) \).

Parameters

x[in] Decision variable \( x \in \mathbb{R}^n \)

using grad_f_sig = void(crvec x, rvec grad_fx)

Signature of the function that evaluates the gradient of the cost function \( \nabla f(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • grad_fx[out] Gradient of cost function \( \nabla f(x) \in \mathbb{R}^n \)

using g_sig = void(crvec x, rvec gx)

Signature of the function that evaluates the constraints \( g(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • gx[out] Value of the constraints \( g(x) \in \mathbb{R}^m \)

using grad_g_prod_sig = void(crvec x, crvec y, rvec grad_gxy)

Signature of the function that evaluates the gradient of the constraints times a vector \( \nabla g(x)\ y \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Vector \( y \in \mathbb{R}^m \) to multiply the gradient by

  • grad_gxy[out] Gradient of the constraints \( \nabla g(x)\ y \in \mathbb{R}^n \)

using grad_gi_sig = void(crvec x, unsigned i, rvec grad_gi)

Signature of the function that evaluates the gradient of one specific constraints \( \nabla g_i(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • i[in] Which constraint \( 0 \le i \lt m \)

  • grad_gi[out] Gradient of the constraint \( \nabla g_i(x) \mathbb{R}^n \)

using hess_L_prod_sig = void(crvec x, crvec y, crvec v, rvec Hv)

Signature of the function that evaluates the Hessian of the Lagrangian multiplied by a vector \( \nabla_{xx}^2L(x, y)\ v \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • v[in] Vector to multiply by \( v \in \mathbb{R}^n \)

  • Hv[out] Hessian-vector product \( \nabla_{xx}^2 L(x, y)\ v \in \mathbb{R}^{n} \)

using hess_L_sig = void(crvec x, crvec y, rmat H)

Signature of the function that evaluates the Hessian of the Lagrangian \( \nabla_{xx}^2L(x, y) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • H[out] Hessian \( \nabla_{xx}^2 L(x, y) \in \mathbb{R}^{n\times n} \)

Public Functions

inline ProblemWithCounters(Problem &&p)
inline ProblemWithCounters(const Problem &p)
ProblemWithCounters() = delete
ProblemWithCounters(const ProblemWithCounters&) = delete
ProblemWithCounters(ProblemWithCounters&&) = delete
ProblemWithCounters &operator=(const ProblemWithCounters&) = delete
ProblemWithCounters &operator=(ProblemWithCounters&&) = delete

Public Members

EvalCounter evaluations
unsigned int n

Number of decision variables, dimension of x.

unsigned int m

Number of constraints, dimension of g(x) and z.

Box C

Constraints of the decision variables, \( x \in C \).

Box D

Other constraints, \( g(x) \in D \).

std::function<f_sig> f

Cost function \( f(x) \).

std::function<grad_f_sig> grad_f

Gradient of the cost function \( \nabla f(x) \).

std::function<g_sig> g

Constraint function \( g(x) \).

std::function<grad_g_prod_sig> grad_g_prod

Gradient of the constraint function times vector \( \nabla g(x)\ y \).

std::function<grad_gi_sig> grad_gi

Gradient of a specific constraint \( \nabla g_i(x) \).

std::function<hess_L_prod_sig> hess_L_prod

Hessian of the Lagrangian function times vector \( \nabla_{xx}^2 L(x, y)\ v \).

std::function<hess_L_sig> hess_L

Hessian of the Lagrangian function \( \nabla_{xx}^2 L(x, y) \).

class ProblemWithParam : public Problem
#include <panoc-alm/util/problem.hpp>

Public Types

using f_sig = real_t(crvec x)

Signature of the function that evaluates the cost \( f(x) \).

Parameters

x[in] Decision variable \( x \in \mathbb{R}^n \)

using grad_f_sig = void(crvec x, rvec grad_fx)

Signature of the function that evaluates the gradient of the cost function \( \nabla f(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • grad_fx[out] Gradient of cost function \( \nabla f(x) \in \mathbb{R}^n \)

using g_sig = void(crvec x, rvec gx)

Signature of the function that evaluates the constraints \( g(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • gx[out] Value of the constraints \( g(x) \in \mathbb{R}^m \)

using grad_g_prod_sig = void(crvec x, crvec y, rvec grad_gxy)

Signature of the function that evaluates the gradient of the constraints times a vector \( \nabla g(x)\ y \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Vector \( y \in \mathbb{R}^m \) to multiply the gradient by

  • grad_gxy[out] Gradient of the constraints \( \nabla g(x)\ y \in \mathbb{R}^n \)

using grad_gi_sig = void(crvec x, unsigned i, rvec grad_gi)

Signature of the function that evaluates the gradient of one specific constraints \( \nabla g_i(x) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • i[in] Which constraint \( 0 \le i \lt m \)

  • grad_gi[out] Gradient of the constraint \( \nabla g_i(x) \mathbb{R}^n \)

using hess_L_prod_sig = void(crvec x, crvec y, crvec v, rvec Hv)

Signature of the function that evaluates the Hessian of the Lagrangian multiplied by a vector \( \nabla_{xx}^2L(x, y)\ v \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • v[in] Vector to multiply by \( v \in \mathbb{R}^n \)

  • Hv[out] Hessian-vector product \( \nabla_{xx}^2 L(x, y)\ v \in \mathbb{R}^{n} \)

using hess_L_sig = void(crvec x, crvec y, rmat H)

Signature of the function that evaluates the Hessian of the Lagrangian \( \nabla_{xx}^2L(x, y) \).

Parameters
  • x[in] Decision variable \( x \in \mathbb{R}^n \)

  • y[in] Lagrange multipliers \( y \in \mathbb{R}^m \)

  • H[out] Hessian \( \nabla_{xx}^2 L(x, y) \in \mathbb{R}^{n\times n} \)

Public Functions

inline void set_param(pa::crvec p)
inline void set_param(pa::vec &&p)
inline pa::vec &get_param()
inline const pa::vec &get_param() const
inline std::shared_ptr<pa::vec> get_param_ptr() const

Public Members

unsigned int n

Number of decision variables, dimension of x.

unsigned int m

Number of constraints, dimension of g(x) and z.

Box C

Constraints of the decision variables, \( x \in C \).

Box D

Other constraints, \( g(x) \in D \).

std::function<f_sig> f

Cost function \( f(x) \).

std::function<grad_f_sig> grad_f

Gradient of the cost function \( \nabla f(x) \).

std::function<g_sig> g

Constraint function \( g(x) \).

std::function<grad_g_prod_sig> grad_g_prod

Gradient of the constraint function times vector \( \nabla g(x)\ y \).

std::function<grad_gi_sig> grad_gi

Gradient of a specific constraint \( \nabla g_i(x) \).

std::function<hess_L_prod_sig> hess_L_prod

Hessian of the Lagrangian function times vector \( \nabla_{xx}^2 L(x, y)\ v \).

std::function<hess_L_sig> hess_L

Hessian of the Lagrangian function \( \nabla_{xx}^2 L(x, y) \).

template<class IndexT = size_t>
struct ReverseCircularIndexIterator
#include <panoc-alm/util/ringbuffer.hpp>

Public Types

using ForwardIterator = CircularIndexIterator<IndexT>
using Index = typename ForwardIterator::Index
using Indices = typename ForwardIterator::Indices
using value_type = Indices
using reference = value_type
using difference_type = std::ptrdiff_t
using pointer = void
using iterator_category = std::input_iterator_tag

Public Functions

inline ReverseCircularIndexIterator()
inline ReverseCircularIndexIterator(Indices i, Index max)
inline ReverseCircularIndexIterator(ForwardIterator forwardit)
inline reference operator*() const
inline ReverseCircularIndexIterator &operator++()
inline ReverseCircularIndexIterator &operator--()
inline ReverseCircularIndexIterator operator++(int)
inline ReverseCircularIndexIterator operator--(int)

Public Members

ForwardIterator forwardit
template<class IndexT>
bool operator==(ReverseCircularIndexIterator<IndexT> a, ReverseCircularIndexIterator<IndexT> b)

Note

Only valid for two indices in the same range.

template<class IndexT>
bool operator!=(ReverseCircularIndexIterator<IndexT> a, ReverseCircularIndexIterator<IndexT> b)

Note

Only valid for two indices in the same range.

template<class IndexT>
class ReverseCircularRange
#include <panoc-alm/util/ringbuffer.hpp>

Public Types

using ForwardRange = CircularRange<IndexT>
using Index = typename ForwardRange::Index
using Indices = typename ForwardRange::Indices
using const_iterator = typename ForwardRange::const_reverse_iterator
using iterator = typename ForwardRange::reverse_iterator
using const_reverse_iterator = typename ForwardRange::const_iterator
using reverse_iterator = typename ForwardRange::iterator

Public Functions

inline ReverseCircularRange(const ForwardRange &forwardrange)
inline ReverseCircularRange(Index size, Index idx1, Index idx2, Index max)
inline iterator begin() const
inline iterator end() const
inline const_iterator cbegin() const
inline const_iterator cend() const
inline reverse_iterator rbegin() const
inline reverse_iterator rend() const
inline const_reverse_iterator crbegin() const
inline const_reverse_iterator crend() const
struct SecondOrderPANOCParams
#include <panoc-alm/inner/decl/second-order-panoc.hpp>

Tuning parameters for the second order PANOC algorithm.

Public Members

LipschitzEstimateParams Lipschitz

Parameters related to the Lipschitz constant estimate and step size.

unsigned max_iter = 100

Maximum number of inner PANOC iterations.

std::chrono::microseconds max_time = std::chrono::minutes(5)

Maximum duration.

real_t τ_min = 1. / 256

Minimum weight factor between Newton step and projected gradient step.

real_t L_min = 1e-5

Minimum Lipschitz constant estimate.

real_t L_max = 1e9

Maximum Lipschitz constant estimate.

PANOCStopCrit stop_crit = PANOCStopCrit::ApproxKKT

What stopping criterion to use.

unsigned max_no_progress = 10

Maximum number of iterations without any progress before giving up.

unsigned print_interval = 0

When to print progress.

If set to zero, nothing will be printed. If set to N != 0, progress is printed every N iterations.

real_t quadratic_upperbound_tolerance_factor = 10 * std::numeric_limits<real_t>::epsilon()
bool update_lipschitz_in_linesearch = true
bool alternative_linesearch_cond = false
class SecondOrderPANOCSolver
#include <panoc-alm/inner/decl/second-order-panoc.hpp>

Second order PANOC solver for ALM.

Public Types

using Params = SecondOrderPANOCParams

Public Functions

inline SecondOrderPANOCSolver(Params params)
inline Stats operator()(const Problem &problem, crvec Σ, real_t ε, bool always_overwrite_results, rvec x, rvec y, rvec err_z)
Parameters
  • problem – [in] Problem description

  • Σ – [in] Constraint weights \( \Sigma \)

  • ε – [in] Tolerance \( \varepsilon \)

  • always_overwrite_results – [in] Overwrite x, y and err_z even if not converged

  • x – [inout] Decision variable \( x \)

  • y – [inout] Lagrange multipliers \( y \)

  • err_z – [out] Slack variable error \( g(x) - z \)

inline SecondOrderPANOCSolver &set_progress_callback(std::function<void(const ProgressInfo&)> cb)
inline std::string get_name() const
inline void stop()
inline const Params &get_params() const
struct ProgressInfo
#include <panoc-alm/inner/decl/second-order-panoc.hpp>

Public Members

unsigned k
crvec x
crvec p
real_t norm_sq_p
crvec x_hat
real_t ψ
crvec grad_ψ
real_t ψ_hat
crvec grad_ψ_hat
real_t L
real_t γ
real_t ε
crvec Σ
crvec y
const Problem &problem
const Params &params
struct Stats
#include <panoc-alm/inner/decl/second-order-panoc.hpp>

Public Members

SolverStatus status = SolverStatus::Unknown
real_t ε = inf
std::chrono::microseconds elapsed_time
unsigned iterations = 0
unsigned newton_failures = 0
unsigned linesearch_failures = 0
unsigned τ_1_accepted = 0
unsigned count_τ = 0
real_t sum_τ = 0
class SpecializedLBFGS
#include <panoc-alm/inner/directions/decl/specialized-lbfgs.hpp>

Limited memory Broyden–Fletcher–Goldfarb–Shanno (L-BFGS) algorithm that can handle updates of the γ parameter.

Public Types

using Params = LBFGSParams

Public Functions

inline SpecializedLBFGS(Params params)
inline SpecializedLBFGS(Params params, size_t n, size_t history)
inline bool standard_update(crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁, crvec gradₖ₊₁)

Standard L-BFGS update without changing the step size γ.

inline bool full_update(crvec xₖ, crvec xₖ₊₁, crvec pₖ_old_γ, crvec pₖ₊₁, crvec gradₖ₊₁, const Box &C, real_t γ)

L-BFGS update when changing the step size γ, recomputing everything.

inline bool update(crvec xₖ, crvec xₖ₊₁, crvec pₖ, crvec pₖ₊₁, crvec gradₖ₊₁, const Box &C, real_t γ)

Update the inverse Hessian approximation using the new vectors xₖ₊₁ and pₖ₊₁.

template<class Vec>
void apply(Vec &&q)

Apply the inverse Hessian approximation to the given vector q.

inline void initialize(crvec x₀, crvec grad₀)

Initialize with the starting point x₀ and the gradient in that point.

inline void reset()

Throw away the approximation and all previous vectors s and y.

inline void resize(size_t n, size_t history)

Re-allocate storage for a problem with a different size.

inline std::string get_name() const
inline size_t n() const

Get the size of the s, y, x and g vectors in the buffer.

inline size_t history() const

Get the number of previous vectors s, y, x and g stored in the buffer.

inline size_t succ(size_t i) const

Get the next index in the circular buffer of previous s, y, x and g vectors.

inline size_t pred(size_t i) const

Get the previous index in the circular buffer of previous s, y, x and g vectors.

inline auto s(size_t i)
inline auto s(size_t i) const
inline auto y(size_t i)
inline auto y(size_t i) const
inline auto x(size_t i)
inline auto x(size_t i) const
inline auto g(size_t i)
inline auto g(size_t i) const
inline auto p()
inline auto p() const
inline auto w()
inline auto w() const
inline real_t &ρ(size_t i)
inline const real_t &ρ(size_t i) const
inline real_t &α(size_t i)
inline const real_t &α(size_t i) const
struct StructuredPANOCLBFGSParams
#include <panoc-alm/inner/decl/structured-panoc-lbfgs.hpp>

Tuning parameters for the second order PANOC algorithm.

Public Members

LipschitzEstimateParams Lipschitz

Parameters related to the Lipschitz constant estimate and step size.

unsigned max_iter = 100

Maximum number of inner PANOC iterations.

std::chrono::microseconds max_time = std::chrono::minutes(5)

Maximum duration.

real_t τ_min = 1. / 256

Minimum weight factor between Newton step and projected gradient step.

real_t L_min = 1e-5

Minimum Lipschitz constant estimate.

real_t L_max = 1e9

Maximum Lipschitz constant estimate.

real_t nonmonotone_linesearch = 0

Factor used in update for exponentially weighted nonmonotone line search.

Zero means monotone line search.

PANOCStopCrit stop_crit = PANOCStopCrit::ApproxKKT

What stopping criterion to use.

unsigned max_no_progress = 10

Maximum number of iterations without any progress before giving up.

unsigned print_interval = 0

When to print progress.

If set to zero, nothing will be printed. If set to N != 0, progress is printed every N iterations.

real_t quadratic_upperbound_tolerance_factor = 10 * std::numeric_limits<real_t>::epsilon()
bool update_lipschitz_in_linesearch = true
bool alternative_linesearch_cond = false
bool hessian_vec_finited_differences = true
bool full_augmented_hessian = true
LBFGSStepSize lbfgs_stepsize = LBFGSStepSize::BasedOnCurvature
struct StructuredPANOCLBFGSProgressInfo
#include <panoc-alm/inner/decl/structured-panoc-lbfgs.hpp>

Public Members

unsigned k
crvec x
crvec p
real_t norm_sq_p
crvec x_hat
real_t φγ
real_t ψ
crvec grad_ψ
real_t ψ_hat
crvec grad_ψ_hat
real_t L
real_t γ
real_t τ
real_t ε
crvec Σ
crvec y
const Problem &problem
const StructuredPANOCLBFGSParams &params
class StructuredPANOCLBFGSSolver
#include <panoc-alm/inner/decl/structured-panoc-lbfgs.hpp>

Second order PANOC solver for ALM.

Public Types

using Params = StructuredPANOCLBFGSParams
using ProgressInfo = StructuredPANOCLBFGSProgressInfo

Public Functions

inline StructuredPANOCLBFGSSolver(Params params, LBFGSParams lbfgsparams)
inline Stats operator()(const Problem &problem, crvec Σ, real_t ε, bool always_overwrite_results, rvec x, rvec y, rvec err_z)
Parameters
  • problem – [in] Problem description

  • Σ – [in] Constraint weights \( \Sigma \)

  • ε – [in] Tolerance \( \varepsilon \)

  • always_overwrite_results – [in] Overwrite x, y and err_z even if not converged

  • x – [inout] Decision variable \( x \)

  • y – [inout] Lagrange multipliers \( y \)

  • err_z – [out] Slack variable error \( g(x) - z \)

inline StructuredPANOCLBFGSSolver &set_progress_callback(std::function<void(const ProgressInfo&)> cb)
inline std::string get_name() const
inline void stop()
inline const Params &get_params() const

Public Members

LBFGS lbfgs
struct Stats
#include <panoc-alm/inner/decl/structured-panoc-lbfgs.hpp>

Public Members

SolverStatus status = SolverStatus::Unknown
real_t ε = inf
std::chrono::microseconds elapsed_time
unsigned iterations = 0
unsigned linesearch_failures = 0
unsigned lbfgs_failures = 0
unsigned lbfgs_rejected = 0
unsigned τ_1_accepted = 0
unsigned count_τ = 0
real_t sum_τ = 0
class vec_allocator
#include <panoc-alm/util/alloc.hpp>

Public Functions

inline vec_allocator(size_t num_vec, Eigen::Index n)
vec_allocator(const vec_allocator&) = delete
vec_allocator(vec_allocator&&) = delete
vec_allocator &operator=(const vec_allocator&) = delete
vec_allocator &operator=(vec_allocator&&) = delete
inline auto alloc()
inline alloc_raii_wrapper alloc_raii()
inline void free(rvec v)
template<class ...Vecs>
inline void free(rvec first, Vecs&&... vecs)
inline size_t size() const
inline size_t used_space() const
inline Eigen::Index vector_size() const
inline size_t highwater() const
struct alloc_raii_wrapper
#include <panoc-alm/util/alloc.hpp>

Public Types

using mvec = Eigen::Map<vec>

Public Functions

inline alloc_raii_wrapper(real_t *dptr, Eigen::Index n, vec_allocator *alloc)
inline alloc_raii_wrapper(mvec &&v, vec_allocator *alloc)
inline ~alloc_raii_wrapper()
alloc_raii_wrapper(const alloc_raii_wrapper&) = delete
inline alloc_raii_wrapper(alloc_raii_wrapper &&o)
alloc_raii_wrapper &operator=(const alloc_raii_wrapper&) = delete
inline alloc_raii_wrapper &operator=(alloc_raii_wrapper &&o)
inline alloc_raii_wrapper &operator=(crvec v)
inline operator crvec() const
inline operator rvec()

Public Members

mvec v
vec_allocator *alloc
namespace detail

Functions

inline void project_y(rvec y, crvec z_lb, crvec z_ub, real_t M)
inline void update_penalty_weights(const ALMParams &params, real_t Δ, bool first_iter, rvec e, rvec old_e, real_t norm_e, real_t old_norm_e, crvec Σ_old, rvec Σ)
inline void initialize_penalty(const Problem &p, const ALMParams &params, crvec x0, rvec Σ)
inline void initialize_penalty(const ProblemFull &p, const ALMParams &params, crvec x0, rvec Σ1, rvec Σ2)
inline void apply_preconditioning(const Problem &problem, Problem &prec_problem, crvec x, real_t &prec_f, vec &prec_g)
inline 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 ∇ψ.

\[ \psi(x^k) = f(x^k) + \frac{1}{2} \text{dist}_\Sigma^2\left(g(x^k) + \Sigma^{-1}y,\;D\right) \]
\[ \hat{y} \]

Parameters
  • p – [in] Problem description

  • x – [in] Decision variable \( x \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ – [in] Penalty weights \( \Sigma \)

  • ŷ – [out] \( \hat{y} \)

inline real_t calc_ψ_ŷ(const ProblemFull &p, crvec x, crvec y, crvec Σ1, crvec Σ2, rvec ŷ1, rvec ŷ2)

Calculate both ψ(x) and the vector ŷ that can later be used to compute ∇ψ.

\[ \psi(x^k) = f(x^k) + \frac{1}{2} \text{dist}_\Sigma^2\left(g(x^k) + \Sigma^{-1}y,\;D\right) \]
\[ \hat{y} \]

Parameters
  • p – [in] Problem description

  • x – [in] Decision variable \( x \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ1 – [in] Penalty weights \( \Sigma_1 \)

  • Σ2 – [in] Penalty weights \( \Sigma_2 \)

  • ŷ1 – [out] \( \hat{y}_1 \)

  • ŷ2 – [out] \( \hat{y}_2 \)

inline void calc_grad_ψ_from_ŷ(const Problem &p, crvec x, crvec ŷ, rvec grad_ψ, rvec work_n)

Calculate ∇ψ(x) using ŷ.

Parameters
  • p – [in] Problem description

  • x – [in] Decision variable \( x \)

  • ŷ – [in] \( \hat{y} \)

  • grad_ψ – [out] \( \nabla \psi(x) \)

  • work_n – Dimension n

inline void calc_grad_ψ_from_ŷ(const ProblemFull &p, crvec x, crvec ŷ1, crvec ŷ2, rvec grad_ψ, rvec work_n)

Calculate ∇ψ(x) using ŷ.

Parameters
  • p – [in] Problem description

  • x – [in] Decision variable \( x \)

  • ŷ1 – [in] \( \hat{y} \)

  • ŷ2 – [in] \( \hat{y} \)

  • grad_ψ – [out] \( \nabla \psi(x) \)

  • work_n – Dimension n

inline 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).

\[ \psi(x^k) = f(x^k) + \frac{1}{2} \text{dist}_\Sigma^2\left(g(x^k) + \Sigma^{-1}y,\;D\right) \]
\[ \nabla \psi(x) = \nabla f(x) + \nabla g(x)\ \hat{y}(x) \]

Parameters
  • p – [in] Problem description

  • x – [in] Decision variable \( x \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ – [in] Penalty weights \( \Sigma \)

  • grad_ψ – [out] \( \nabla \psi(x) \)

  • work_n – Dimension n

  • work_m – Dimension m

inline real_t calc_ψ_grad_ψ(const ProblemFull &p, crvec x, crvec y, crvec Σ1, crvec Σ2, rvec grad_ψ, rvec work_n, rvec work_m1, rvec work_m2)

Calculate both ψ(x) and its gradient ∇ψ(x).

\[ \psi(x^k) = f(x^k) + \frac{1}{2} \text{dist}_\Sigma^2\left(g(x^k) + \Sigma^{-1}y,\;D\right) \]
\[ \nabla \psi(x) = \nabla f(x) + \nabla g(x)\ \hat{y}(x) \]

Parameters
  • p – [in] Problem description

  • x – [in] Decision variable \( x \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ1 – [in] Penalty weights \( \Sigma_1 \)

  • Σ2 – [in] Penalty weights \( \Sigma_2 \)

  • grad_ψ – [out] \( \nabla \psi(x) \)

  • work_n – Dimension n

  • work_m1 – Dimension m1

  • work_m2 – Dimension m2

inline void calc_grad_ψ(const Problem &p, crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m)

Calculate the gradient ∇ψ(x).

\[ \nabla \psi(x) = \nabla f(x) + \nabla g(x)\ \hat{y}(x) \]

Parameters
  • p – [in] Problem description

  • x – [in] Decision variable \( x \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ – [in] Penalty weights \( \Sigma \)

  • grad_ψ – [out] \( \nabla \psi(x) \)

  • work_n – Dimension n

  • work_m – Dimension m

inline void calc_grad_ψ(const ProblemFull &p, crvec x, crvec y, crvec Σ1, crvec Σ2, rvec grad_ψ, rvec work_n, rvec work_m1, rvec work_m2)

Calculate the gradient ∇ψ(x).

\[ \nabla \psi(x) = \nabla f(x) + \nabla g(x)\ \hat{y}(x) \]

Parameters
  • p – [in] Problem description

  • x – [in] Decision variable \( x \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ1 – [in] Penalty weights \( \Sigma_1 \)

  • Σ2 – [in] Penalty weights \( \Sigma_2 \)

  • grad_ψ – [out] \( \nabla \psi(x) \)

  • work_n – Dimension n

  • work_m1 – Dimension m_1

  • work_m2 – Dimension m_2

inline void calc_err_z (const Problem &p, crvec x̂, crvec y, crvec Σ, rvec err_z)

Calculate the error between ẑ and g(x).

\[ \hat{z}^k = \Pi_D\left(g(x^k) + \Sigma^{-1}y\right) \]

Parameters
  • p – [in] Problem description

  • – [in] Decision variable \( \hat{x} \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ – [in] Penalty weights \( \Sigma \)

  • err_z – [out] \( g(\hat{x}) - \hat{z} \)

inline void calc_err_z (const ProblemFull &p, crvec x̂, crvec y, crvec Σ1, rvec err_z1, rvec err_z2)

Calculate the error between ẑ and g(x).

\[ \hat{z}^k = \Pi_D\left(g(x^k) + \Sigma^{-1}y\right) \]

Parameters
  • p – [in] Problem description

  • – [in] Decision variable \( \hat{x} \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ1 – [in] Penalty weights \( \Sigma_1 \)

  • err_z1 – [out] \( g_1(\hat{x}) - \hat{z}_1 \)

  • err_z2 – [out] \( g_2(\hat{x}) - \hat{z}_2 \)

inline auto projected_gradient_step(const Box &C, real_t γ, crvec x, crvec grad_ψ)

Projected gradient step.

\[\begin{split} \begin{aligned} \hat{x}^k &= T_{\gamma^k}\left(x^k\right) \\ &= \Pi_C\left(x^k - \gamma^k \nabla \psi(x^k)\right) \\ p^k &= \hat{x}^k - x^k \\ \end{aligned} \end{split}\]

Parameters
  • C – [in] Set to project onto

  • γ – [in] Step size

  • x – [in] Decision variable \( x \)

  • grad_ψ – [in] \( \nabla \psi(x^k) \)

template<typename ProblemT> inline void calc_x̂ (const ProblemT &prob, real_t γ, crvec x, crvec grad_ψ, rvec x̂, rvec p)
Parameters
  • prob – [in] Problem description

  • γ – [in] Step size

  • x – [in] Decision variable \( x \)

  • grad_ψ – [in] \( \nabla \psi(x^k) \)

  • – [out] \( \hat{x}^k = T_{\gamma^k}(x^k) \)

  • p – [out] \( \hat{x}^k - x^k \)

inline real_t calc_error_stop_crit(PANOCStopCrit crit, crvec pₖ, real_t γ, crvec xₖ, crvec grad_̂ψₖ, crvec grad_ψₖ, const Box &C)

\[ \left\| \gamma_k^{-1} (x^k - \hat x^k) + \nabla \psi(\hat x^k) - \nabla \psi(x^k) \right\|_\infty \]

Parameters
  • crit – [in] What stoppint criterion to use

  • pₖ – [in] Projected gradient step \( \hat x^k - x^k \)

  • γ – [in] Step size

  • xₖ – [in] Current iterate

  • grad_̂ψₖ – [in] Gradient in \( \hat x^k \)

  • grad_ψₖ – [in] Gradient in \( x^k \)

  • C – [in] Feasible set \( C \)

inline 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 until quadratic upper bound or descent lemma is satisfied:

\[ \psi(x + p) \le \psi(x) + \nabla\psi(x)^\top p + \frac{L}{2} \|p\|^2 \]
The projected gradient iterate \( \hat x^k \) and step \( p^k \) are updated with the new step size \( \gamma^k \), and so are other quantities that depend on them, such as \( \nabla\psi(x^k)^\top p^k \) and \( \|p\|^2 \). The intermediate vector \( \hat y(x^k) \) can be used to compute the gradient \( \nabla\psi(\hat x^k) \) or to update the Lagrange multipliers.

Parameters
  • problem – [in] Problem description

  • rounding_tolerance – [in] Tolerance used to ignore rounding errors when the function \( \psi(x) \) is relatively flat or the step size is very small, which could cause \( \psi(x^k) < \psi(\hat x^k) \), which is mathematically impossible but could occur in finite precision floating point arithmetic.

  • L_max – [in] Maximum allowed Lipschitz constant estimate (prevents infinite loop if function or derivatives are discontinuous, and keeps step size bounded away from zero).

  • xₖ – [in] Current iterate \( x^k \)

  • ψₖ – [in] Objective function \( \psi(x^k) \)

  • grad_ψₖ – [in] Gradient of objective \( \nabla\psi(x^k) \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ – [in] Penalty weights \( \Sigma \)

  • x̂ₖ – [out] Projected gradient iterate \( \hat x^k \)

  • pₖ – [out] Projected gradient step \( p^k \)

  • ŷx̂ₖ – [out] Intermediate vector \( \hat y(x^k) \)

  • ψx̂ₖ – [inout] Objective function \( \psi(\hat x^k) \)

  • norm_sq_pₖ – [inout] Squared norm of the step \( \left\| p^k \right\|^2 \)

  • grad_ψₖᵀpₖ – [inout] Gradient of objective times step \( \nabla\psi(x^k)^\top p^k\)

  • Lₖ – [inout] Lipschitz constant estimate \( L_{\nabla\psi}^k \)

  • γₖ – [inout] Step size \( \gamma^k \)

Returns

The original step size, before it was reduced by this function.

inline real_t descent_lemma(const ProblemFull &problem, real_t rounding_tolerance, real_t L_max, crvec xₖ, real_t ψₖ, crvec grad_ψₖ, crvec y, crvec Σ1, crvec Σ2, rvec x̂ₖ, rvec pₖ, rvec ŷx̂ₖ1, rvec ŷx̂ₖ2, 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 until quadratic upper bound or descent lemma is satisfied:

\[ \psi(x + p) \le \psi(x) + \nabla\psi(x)^\top p + \frac{L}{2} \|p\|^2 \]
The projected gradient iterate \( \hat x^k \) and step \( p^k \) are updated with the new step size \( \gamma^k \), and so are other quantities that depend on them, such as \( \nabla\psi(x^k)^\top p^k \) and \( \|p\|^2 \). The intermediate vector \( \hat y(x^k) \) can be used to compute the gradient \( \nabla\psi(\hat x^k) \) or to update the Lagrange multipliers.

Parameters
  • problem – [in] Problem description

  • rounding_tolerance – [in] Tolerance used to ignore rounding errors when the function \( \psi(x) \) is relatively flat or the step size is very small, which could cause \( \psi(x^k) < \psi(\hat x^k) \), which is mathematically impossible but could occur in finite precision floating point arithmetic.

  • L_max – [in] Maximum allowed Lipschitz constant estimate (prevents infinite loop if function or derivatives are discontinuous, and keeps step size bounded away from zero).

  • xₖ – [in] Current iterate \( x^k \)

  • ψₖ – [in] Objective function \( \psi(x^k) \)

  • grad_ψₖ – [in] Gradient of objective \( \nabla\psi(x^k) \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ1 – [in] Penalty weights \( \Sigma_1 \)

  • Σ2 – [in] Penalty weights \( \Sigma_2 \)

  • x̂ₖ – [out] Projected gradient iterate \( \hat x^k \)

  • pₖ – [out] Projected gradient step \( p^k \)

  • ŷx̂ₖ1 – [out] Intermediate vector \( \hat y(x^k)1 \)

  • ŷx̂ₖ2 – [out] Intermediate vector \( \hat y(x^k)2 \)

  • ψx̂ₖ – [inout] Objective function \( \psi(\hat x^k) \)

  • norm_sq_pₖ – [inout] Squared norm of the step \( \left\| p^k \right\|^2 \)

  • grad_ψₖᵀpₖ – [inout] Gradient of objective times step \( \nabla\psi(x^k)^\top p^k\)

  • Lₖ – [inout] Lipschitz constant estimate \( L_{\nabla\psi}^k \)

  • γₖ – [inout] Step size \( \gamma^k \)

Returns

The original step size, before it was reduced by this function.

template<class ParamsT, class DurationT>
inline 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 exceeded, interrupted by user, infinite iterate, no progress made)

Parameters
  • params – [in] Parameters including `max_iter`, `max_time` and `max_no_progress`

  • time_elapsed – [in] Time elapsed since the start of the algorithm

  • iteration – [in] The current iteration number

  • stop_signal – [in] A stop signal for the user to interrupt the algorithm

  • ε – [in] Desired primal tolerance

  • εₖ – [in] Tolerance of the current iterate

  • no_progress – [in] The number of successive iterations no progress was made

inline void calc_augmented_lagrangian_hessian(const Problem &problem, crvec xₖ, crvec ŷxₖ, crvec y, crvec Σ, rvec g, mat &H, rvec work_n)

Compute the Hessian matrix of the augmented Lagrangian function.

\[ \nabla^2_{xx} L_\Sigma(x, y) = \Big. \nabla_{xx}^2 L(x, y) \Big|_{\big(x,\, \hat y(x, y)\big)} + \sum_{i\in\mathcal{I}} \Sigma_i\,\nabla g_i(x) \nabla g_i(x)^\top \]

Parameters
  • problem – [in] Problem description

  • xₖ – [in] Current iterate \( x^k \)

  • ŷxₖ – [in] Intermediate vector \( \hat y(x^k) \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ – [in] Penalty weights \( \Sigma \)

  • g – [out] The constraint values \( g(x^k) \)

  • H – [out] Hessian matrix \( H(x, y) \)

  • work_n – Dimension n

inline 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, using finite differences.

\[ \nabla^2_{xx} L_\Sigma(x, y)\, v \approx \frac{\nabla_x L_\Sigma(x+hv, y) - \nabla_x L_\Sigma(x, y)}{h} \]

Parameters
  • problem – [in] Problem description

  • xₖ – [in] Current iterate \( x^k \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ – [in] Penalty weights \( \Sigma \)

  • grad_ψ – [in] Gradient \( \nabla \psi(x^k) \)

  • v – [in] Vector to multiply with the Hessian

  • Hv – [out] Hessian-vector product

  • work_n1 – Dimension n

  • work_n2 – Dimension n

  • work_m – Dimension m

inline 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 \( \nabla \psi \) using finite differences.

Parameters
  • problem – [in] Problem description

  • xₖ – [in] Current iterate \( x^k \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ – [in] Penalty weights \( \Sigma \)

  • ε – [in] Finite difference step size relative to xₖ

  • δ – [in] Minimum absolute finite difference step size

  • L_min – [in] Minimum allowed Lipschitz estimate.

  • L_max – [in] Maximum allowed Lipschitz estimate.

  • ψ – [out] \( \psi(x^k) \)

  • grad_ψ – [out] Gradient \( \nabla \psi(x^k) \)

  • work_n1 – Dimension n

  • work_n2 – Dimension n

  • work_n3 – Dimension n

  • work_m – Dimension m

inline 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, rvec grad_ψ, rvec work_n1, rvec work_n2, rvec work_n3, rvec work_m)

Estimate the Lipschitz constant of the gradient \( \nabla \psi \) using finite differences.

Parameters
  • problem – [in] Problem description

  • xₖ – [in] Current iterate \( x^k \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ – [in] Penalty weights \( \Sigma \)

  • ε – [in] Finite difference step size relative to xₖ

  • δ – [in] Minimum absolute finite difference step size

  • L_min – [in] Minimum allowed Lipschitz estimate.

  • L_max – [in] Maximum allowed Lipschitz estimate.

  • grad_ψ – [out] Gradient \( \nabla \psi(x^k) \)

  • work_n1 – Dimension n

  • work_n2 – Dimension n

  • work_n3 – Dimension n

  • work_m – Dimension m

inline real_t initial_lipschitz_estimate(const ProblemFull &problem, crvec xₖ, crvec y, crvec Σ1, crvec Σ2, real_t ε, real_t δ, real_t &ψ, rvec grad_ψ, rvec work_n1, rvec work_n2, rvec work_n3, rvec work_m1, rvec work_m2)
Parameters
  • problem – [in] Problem description

  • xₖ – [in] Current iterate \( x^k \)

  • y – [in] Lagrange multipliers \( y \)

  • Σ1 – [in] Penalty weights \( \Sigma \)

  • Σ2 – [in] Penalty weights \( \Sigma \)

  • ε – [in] Finite difference step size relative to xₖ

  • δ – [in] Minimum absolute finite difference step size

  • ψ – [out] \( \psi(x^k) \)

  • grad_ψ – [out] Gradient \( \nabla \psi(x^k) \)

  • work_n1 – Dimension n

  • work_n2 – Dimension n

  • work_n3 – Dimension n

  • work_m1 – Dimension m1

  • work_m2 – Dimension m2

template<class ObjFunT, class ObjGradFunT>
real_t initial_lipschitz_estimate(const ObjFunT &ψ, const ObjGradFunT &grad_ψ, crvec xₖ, real_t ε, real_t δ, real_t L_min, real_t L_max, real_t &ψₖ, rvec grad_ψₖ, vec_allocator &alloc)

Estimate the Lipschitz constant of the gradient \( \nabla \psi \) using finite differences.

Parameters
  • ψ – [in] Objective

  • grad_ψ – [in] Gradient of

  • xₖ – [in] Current iterate \( x^k \)

  • ε – [in] Finite difference step size relative to xₖ

  • δ – [in] Minimum absolute finite difference step size

  • L_min – [in] Minimum allowed Lipschitz estimate.

  • L_max – [in] Maximum allowed Lipschitz estimate.

  • ψₖ – [out] \( \psi(x^k) \)

  • grad_ψₖ – [out] Gradient \( \nabla \psi(x^k) \)

  • alloc – [in]

inline void calc_x̂ (real_t γ, crvec x, const Box &C, crvec grad_ψ, rvec x̂, rvec p)
Parameters
  • γ – [in] Step size

  • x – [in] Decision variable \( x \)

  • C – [in] Box constraints \( C \)

  • grad_ψ – [in] \( \nabla \psi(x^k) \)

  • – [out] \( \hat{x}^k = T_{\gamma^k}(x^k) \)

  • p – [out] \( \hat{x}^k - x^k \)

template<class ObjFunT>
real_t descent_lemma(const ObjFunT &ψ, const Box &C, real_t rounding_tolerance, real_t L_max, crvec xₖ, real_t ψₖ, crvec grad_ψₖ, rvec x̂ₖ, rvec pₖ, 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 until quadratic upper bound or descent lemma is satisfied:

\[ \psi(x + p) \le \psi(x) + \nabla\psi(x)^\top p + \frac{L}{2} \|p\|^2 \]
The projected gradient iterate \( \hat x^k \) and step \( p^k \) are updated with the new step size \( \gamma^k \), and so are other quantities that depend on them, such as \( \nabla\psi(x^k)^\top p^k \) and \( \|p\|^2 \). The intermediate vector \( \hat y(x^k) \) can be used to compute the gradient \( \nabla\psi(\hat x^k) \) or to update the Lagrange multipliers.

Parameters
  • ψ – [in] Objective.

  • C – [in] Box constraints.

  • rounding_tolerance – [in] Tolerance used to ignore rounding errors when the function \( \psi(x) \) is relatively flat or the step size is very small, which could cause \( \psi(x^k) < \psi(\hat x^k) \), which is mathematically impossible but could occur in finite precision floating point arithmetic.

  • L_max – [in] Maximum allowed Lipschitz constant estimate (prevents infinite loop if function or derivatives are discontinuous, and keeps step size bounded away from zero).

  • xₖ – [in] Current iterate \( x^k \)

  • ψₖ – [in] Objective function \( \psi(x^k) \)

  • grad_ψₖ – [in] Gradient of objective \( \nabla\psi(x^k) \)

  • x̂ₖ – [out] Projected gradient iterate \( \hat x^k \)

  • pₖ – [out] Projected gradient step \( p^k \)

  • ψx̂ₖ – [inout] Objective function \( \psi(\hat x^k) \)

  • norm_sq_pₖ – [inout] Squared norm of the step \( \left\| p^k \right\|^2 \)

  • grad_ψₖᵀpₖ – [inout] Gradient of objective times step \( \nabla\psi(x^k)^\top p^k\)

  • Lₖ – [inout] Lipschitz constant estimate \( L_{\nabla\psi}^k \)

  • γₖ – [inout] Step size \( \gamma^k \)

Returns

The original step size, before it was reduced by this function.

inline void print_progress(unsigned k, real_t ψₖ, crvec grad_ψₖ, real_t pₖᵀpₖ, real_t γₖ, real_t εₖ)
namespace problems

Functions

Problem himmelblau_problem()
Problem riskaverse_mpc_problem()
struct RiskaverseProblem

Public Types

using Diag = Eigen::DiagonalMatrix<real_t, Eigen::Dynamic, Eigen::Dynamic>

Public Functions

template<class VecX>
inline auto u(VecX &&x) const
template<class VecX>
inline auto s(VecX &&x) const
template<class VecX>
inline auto y(VecX &&x) const
inline Box get_C() const
inline Box get_D() const
inline RiskaverseProblem()
inline auto mpc_dynamics(crvec x, crvec u) const
inline real_t f(crvec ux) const
inline void grad_f(crvec ux, rvec grad_f) const
inline void g(crvec ux, rvec g_u) const
inline void grad_g(crvec ux, crvec v, rvec grad_u_v) const

Public Members

unsigned nu = 2
unsigned nx = 4
unsigned ns = 2
unsigned ny = 3
unsigned n = nu + ns + ny
unsigned m = 3
real_t Ts = 0.05
mat A
mat B
mat Q
mat R
vec x0
namespace vec_util

Functions

template<class V, class M>
auto norm_squared_weighted(V &&v, M &&Σ)

Get the Σ norm squared of a given vector, with Σ a diagonal matrix.

Returns

\( \langle v, \Sigma v \rangle \)

template<class Vec>
real_t norm_inf(const Vec &v)

Get the maximum or infinity-norm of the given vector.

Returns

\( \left\|v\right\|_\infty \)