Problems#

group grp_Problems

Classes for defining optimization problems.

Functions

template<class Problem>
auto problem_with_counters(Problem &&p)#

Wraps the given problem into a ProblemWithCounters and keeps track of how many times each function is called, and how long these calls took.

The wrapper has its own copy of the given problem. Making copies of the wrapper also copies the underlying problem, but does not copy the evaluation counters, all copies share the same counters.

template<class Problem>
auto problem_with_counters_ref(Problem &p)#

Wraps the given problem into a ProblemWithCounters and keeps track of how many times each function is called, and how long these calls took.

The wrapper keeps only a reference to the given problem, it is the responsibility of the caller to make sure that the wrapper does not outlive the original problem. Making copies of the wrapper does not copy the evaluation counters, all copies share the same counters.

template<Config Conf>
void print_provided_functions(std::ostream &os, const TypeErasedProblem<Conf> &problem)#
template<Config Conf>
class BoxConstrProblem#
#include <alpaqa/problem/box-constr-problem.hpp>

Implements common problem functions for minimization problems with box constraints.

Meant to be used as a base class for custom problem implementations. Supports optional \( \ell_1 \)-regularization.

Subclassed by CasADiProblem< Conf >, FunctionalProblem< Conf >

Public Types

using Box = alpaqa::Box<config_t>#

Public Functions

inline BoxConstrProblem(length_t n, length_t m)#

Create a problem with inactive boxes \( (-\infty, +\infty) \), with no \( \ell_1 \)-regularization, and all general constraints handled using ALM.

Parameters:
  • n – Number of decision variables

  • m – Number of constraints

inline BoxConstrProblem(Box C, Box D, vec l1_reg = vec(0), index_t penalty_alm_split = 0)#
inline void resize(length_t n, length_t m)#

Change the dimensions of the problem (number of decision variables and number of constaints).

Destructive: resizes and/or resets the members C, D, l1_reg and penalty_alm_split.

BoxConstrProblem(const BoxConstrProblem&) = default#
BoxConstrProblem &operator=(const BoxConstrProblem&) = default#
BoxConstrProblem(BoxConstrProblem&&) noexcept = default#
BoxConstrProblem &operator=(BoxConstrProblem&&) noexcept = default#
inline length_t get_n() const#

Number of decision variables, n.

inline length_t get_m() const#

Number of constraints, m.

inline real_t eval_prox_grad_step(real_t γ, crvec x, crvec grad_ψ, rvec x_hat, rvec p) const#

inline void eval_proj_diff_g(crvec z, rvec p) const#

inline void eval_proj_multipliers(rvec y, real_t M) const#

inline const Box &get_box_C() const#

inline const Box &get_box_D() const#

inline index_t eval_inactive_indices_res_lna(real_t γ, crvec x, crvec grad_ψ, rindexvec J) const#

inline void check() const#

Public Members

length_t n#

Number of decision variables, dimension of x.

length_t m#

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

Box C = {this->n}#

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

Box D = {this->m}#

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

vec l1_reg = {}#

\( \ell_1 \) (1-norm) regularization parameter.

Possible dimensions are: \( 0 \) (no regularization), \( 1 \) (a single scalar factor), or \( n \) (a different factor for each variable).

index_t penalty_alm_split = 0#

Components of the constraint function with indices below this number are handled using a quadratic penalty method rather than using an augmented Lagrangian method.

Specifically, the Lagrange multipliers for these components (which determine the shifts in ALM) are kept at zero.

Public Static Functions

static inline real_t eval_proj_grad_step_box(const Box &C, real_t γ, crvec x, crvec grad_ψ, rvec x_hat, rvec p)#

Projected gradient step for rectangular box C.

\[\begin{split} \begin{aligned} \hat x &= \Pi_C(x - \gamma\nabla\psi(x)) \\ p &= \hat x - x \\ &= \max(\underline x - x, \;\min(-\gamma\nabla\psi(x), \overline x - x) \end{aligned} \end{split}\]

static inline void eval_prox_grad_step_box_l1_impl(const Box &C, const auto &λ, real_t γ, crvec x, crvec grad_ψ, rvec x_hat, rvec p)#

Proximal gradient step for rectangular box C with ℓ₁-regularization.

\[\begin{split} \begin{aligned} h(x) &= \|x\|_1 + \delta_C(x) \\ \hat x &= \prox_{\gamma h}(x - \gamma\nabla\psi(x)) \\ &= -\max\big( x - \overline x, \;\min\big( x - \underline x, \;\min\big( \gamma(\nabla\psi(x) + \lambda), \;\max\big( \gamma(\nabla\psi(x) - \lambda), x \big) \big) \big) \big) \end{aligned} \end{split}\]

static inline real_t eval_prox_grad_step_box_l1(const Box &C, const auto &λ, real_t γ, crvec x, crvec grad_ψ, rvec x_hat, rvec p)#

Proximal gradient step for rectangular box C with ℓ₁-regularization.

\[\begin{split} \begin{aligned} h(x) &= \|x\|_1 + \delta_C(x) \\ \hat x &= \prox_{\gamma h}(x - \gamma\nabla\psi(x)) \\ &= -\max\big( x - \overline x, \;\min\big( x - \underline x, \;\min\big( \gamma(\nabla\psi(x) + \lambda), \;\max\big( \gamma(\nabla\psi(x) - \lambda), x \big) \big) \big) \big) \end{aligned} \end{split}\]

static inline real_t eval_prox_grad_step_box_l1_scal(const Box &C, real_t λ, real_t γ, crvec x, crvec grad_ψ, rvec x_hat, rvec p)#

Proximal gradient step for rectangular box C with ℓ₁-regularization.

\[\begin{split} \begin{aligned} h(x) &= \|x\|_1 + \delta_C(x) \\ \hat x &= \prox_{\gamma h}(x - \gamma\nabla\psi(x)) \\ &= -\max\big( x - \overline x, \;\min\big( x - \underline x, \;\min\big( \gamma(\nabla\psi(x) + \lambda), \;\max\big( \gamma(\nabla\psi(x) - \lambda), x \big) \big) \big) \big) \end{aligned} \end{split}\]

static inline void eval_proj_multipliers_box(const Box &D, rvec y, real_t M, index_t penalty_alm_split)#
template<Config Conf = DefaultConfig>
class FunctionalProblem : public alpaqa::BoxConstrProblem<DefaultConfig>#
#include <alpaqa/problem/functional-problem.hpp>

Problem class that allows specifying the basic functions as C++ std::functions.

Public Functions

inline real_t eval_f(crvec x) const#
inline void eval_grad_f(crvec x, rvec grad_fx) const#
inline void eval_g(crvec x, rvec gx) const#
inline void eval_grad_g_prod(crvec x, crvec y, rvec grad_gxy) const#
inline void eval_grad_gi(crvec x, index_t i, rvec grad_gix) const#
inline void eval_hess_L_prod(crvec x, crvec y, real_t scale, crvec v, rvec Hv) const#
inline void eval_hess_ψ_prod(crvec x, crvec y, crvec Σ, real_t scale, crvec v, rvec Hv) const#
inline void eval_jac_g(crvec x, rvec J_values) const#
inline void eval_hess_L(crvec x, crvec y, real_t scale, rvec H_values) const#
inline void eval_hess_ψ(crvec x, crvec y, crvec Σ, real_t scale, rvec H_values) const#
inline bool provides_eval_grad_gi() const#

inline bool provides_eval_jac_g() const#

inline bool provides_eval_hess_L_prod() const#

inline bool provides_eval_hess_L() const#

inline bool provides_eval_hess_ψ_prod() const#

inline bool provides_eval_hess_ψ() const#

FunctionalProblem(const FunctionalProblem&) = default#
FunctionalProblem &operator=(const FunctionalProblem&) = default#
FunctionalProblem(FunctionalProblem&&) noexcept = default#
FunctionalProblem &operator=(FunctionalProblem&&) noexcept = default#

Public Members

std::function<real_t(crvec)> f#
std::function<void(crvec, rvec)> grad_f#
std::function<void(crvec, rvec)> g#
std::function<void(crvec, crvec, rvec)> grad_g_prod#
std::function<void(crvec, index_t, rvec)> grad_gi#
std::function<void(crvec, rmat)> jac_g#
std::function<void(crvec, crvec, real_t, crvec, rvec)> hess_L_prod#
std::function<void(crvec, crvec, real_t, rmat)> hess_L#
std::function<void(crvec, crvec, crvec, real_t, crvec, rvec)> hess_ψ_prod#
std::function<void(crvec, crvec, crvec, real_t, rmat)> hess_ψ#
template<Config Conf = DefaultConfig, class Allocator = std::allocator<std::byte>>
class TypeErasedControlProblem : public alpaqa::util::TypeErased<ControlProblemVTable<DefaultConfig>, std::allocator<std::byte>>#
#include <alpaqa/problem/ocproblem.hpp>

Nonlinear optimal control problem with finite horizon \( N \).

\[\begin{split} \newcommand\U{U} \newcommand\D{D} \newcommand\nnu{{n_u}} \newcommand\nnx{{n_x}} \newcommand\nny{{n_y}} \newcommand\xinit{x_\text{init}} \begin{equation}\label{eq:OCP} \tag{OCP}\hspace{-0.8em} \begin{aligned} &\minimize_{u,x} && \sum_{k=0}^{N-1} \ell_k\big(h_k(x^k, u^k)\big) + \ell_N\big(h_N(x^N)\big)\hspace{-0.8em} \\ &\subjto && u^k \in \U \\ &&& c_k(x^k) \in \D \\ &&& c_N(x^N) \in \D_N \\ &&& x^0 = \xinit \\ &&& x^{k+1} = f(x^k, u^k) \quad\quad (0 \le k \lt N) \end{aligned} \end{equation} \end{split}\]

The function \( f : \R^\nnx \times \R^\nnu \to \R^\nnx \) models the discrete-time, nonlinear dynamics of the system, which starts from an initial state \( \xinit \). The functions \( h_k : \R^\nnx \times \R^\nnu \to \R^{n_h} \) for \( 0 \le k \lt N \) and \( h_N : \R^\nnx \to \R^{n_h^N} \) can be used to represent the (possibly time-varying) output mapping of the system, and the convex functions \( \ell_k : \R^{n_h} \to \R \) and \( \ell_N : \R^{n_h^N} \to \R \) define the stage costs and the terminal cost respectively. Stage constraints and terminal constraints are represented by the functions \( c_k : \R^{n_x} \to \R^{n_c} \) and \( c_N : \R^{n_x} \to \R^{n_c^N} \), and the boxes \( D \) and \( D_N \).

Additional functions for computing Gauss-Newton approximations of the cost Hessian are included as well:

\[\begin{split} \begin{aligned} q^k &\defeq \tp{\jac_{h_k}^x\!(\barxuk)} \nabla \ell_k(\hhbar^k) \\ r^k &\defeq \tp{\jac_{h_k}^u\!(\barxuk)} \nabla \ell_k(\hhbar^k) \\ \Lambda_k &\defeq \partial^2 \ell_k(\hhbar^k) \\ Q_k &\defeq \tp{\jac_{h_k}^x\!(\barxuk)} \Lambda_k\, \jac_{h_k}^x\!(\barxuk) \\ S_k &\defeq \tp{\jac_{h_k}^u\!(\barxuk)} \Lambda_k\, \jac_{h_k}^x\!(\barxuk) \\ R_k &\defeq \tp{\jac_{h_k}^u\!(\barxuk)} \Lambda_k\, \jac_{h_k}^u\!(\barxuk). \\ \end{aligned} \end{split}\]
See [3] for more details.

Problem dimensions

inline length_t get_N() const#

Horizon length.

inline length_t get_nu() const#

Number of inputs.

inline length_t get_nx() const#

Number of states.

inline length_t get_nh() const#

Number of outputs.

inline length_t get_nh_N() const#
inline length_t get_nc() const#

Number of constraints.

inline length_t get_nc_N() const#
inline Dim get_dim() const#

All dimensions.

inline length_t get_n() const#

Total number of variables.

inline length_t get_m() const#

Total number of constraints.

Projections onto constraint sets

inline void eval_proj_diff_g(crvec z, rvec e) const#

[Required] Function that evaluates the difference between the given point \( z \) and its projection onto the constraint set \( D \).

Note

z and e can refer to the same vector.

Parameters:
  • z[in] Slack variable, \( z \in \R^m \)

  • e[out] The difference relative to its projection, \( e = z - \Pi_D(z) \in \R^m \)

inline void eval_proj_multipliers(rvec y, real_t M) const#

[Required] Function that projects the Lagrange multipliers for ALM.

Parameters:
  • y[inout] Multipliers, \( y \leftarrow \Pi_Y(y) \in \R^m \)

  • M[in] The radius/size of the set \( Y \). See ALMParams::max_multiplier.

Constraint sets

inline void get_U(Box &U) const#

Input box constraints \( U \).

inline void get_D(Box &D) const#

Stage box constraints \( D \).

inline void get_D_N(Box &D) const#

Terminal box constraints \( D_N \).

Dynamics and initial state

inline void get_x_init(rvec x_init) const#

Initial state \( x_\text{init} \).

inline void eval_f(index_t timestep, crvec x, crvec u, rvec fxu) const#

Discrete-time dynamics \( x^{k+1} = f_k(x^k, u^k) \).

inline void eval_jac_f(index_t timestep, crvec x, crvec u, rmat J_fxu) const#

Jacobian of discrete-time dynamics \( \jac_f(x^k, u^k) \).

inline void eval_grad_f_prod(index_t timestep, crvec x, crvec u, crvec p, rvec grad_fxu_p) const#

Gradient-vector product of discrete-time dynamics \( \nabla f(x^k, u^k)\,p \).

Output mapping

inline void eval_h(index_t timestep, crvec x, crvec u, rvec h) const#

Stage output mapping \( h_k(x^k, u^k) \).

inline void eval_h_N(crvec x, rvec h) const#

Terminal output mapping \( h_N(x^N) \).

Stage and terminal cost

inline real_t eval_l(index_t timestep, crvec h) const#

Stage cost \( \ell_k(\hbar^k) \).

inline real_t eval_l_N(crvec h) const#

Terminal cost \( \ell_N(\hbar^N) \).

Gauss-Newton approximations

inline void eval_qr(index_t timestep, crvec xu, crvec h, rvec qr) const#

Cost gradients w.r.t.

states and inputs \( q^k = \tp{\jac_{h_k}^x\!(\barxuk)} \nabla \ell_k(\hbar^k) \) and \( r^k = \tp{\jac_{h_k}^u\!(\barxuk)} \nabla \ell_k(\hbar^k) \).

inline void eval_q_N(crvec x, crvec h, rvec q) const#

Terminal cost gradient w.r.t.

states \( q^N = \tp{\jac_{h_N}(\bar x^N)} \nabla \ell_k(\hbar^N) \).

inline void eval_add_Q(index_t timestep, crvec xu, crvec h, rmat Q) const#

Cost Hessian w.r.t.

states \( Q_k = \tp{\jac_{h_k}^x\!(\barxuk)} \partial^2\ell_k(\hbar^k)\, \jac_{h_k}^x\!(\barxuk) \), added to the given matrix Q. \( Q \leftarrow Q + Q_k \).

inline void eval_add_Q_N(crvec x, crvec h, rmat Q) const#

Terminal cost Hessian w.r.t.

states \( Q_N = \tp{\jac_{h_N}(\bar x^N)} \partial^2\ell_N(\hbar^N)\, \jac_{h_N}(\bar x^N) \), added to the given matrix Q. \( Q \leftarrow Q + Q_N \).

inline void eval_add_R_masked(index_t timestep, crvec xu, crvec h, crindexvec mask, rmat R, rvec work) const#

Cost Hessian w.r.t.

inputs \( R_k = \tp{\jac_{h_k}^u\!(\barxuk)} \partial^2\ell_k(\hbar^k)\, \jac_{h_k}^u\!(\barxuk) \), keeping only rows and columns in the mask \( \mathcal J \), added to the given matrix R. \( R \leftarrow R + R_k[\mathcal J, \mathcal J] \). The size of work should be get_R_work_size().

inline void eval_add_S_masked(index_t timestep, crvec xu, crvec h, crindexvec mask, rmat S, rvec work) const#

Cost Hessian w.r.t.

inputs and states \( S_k = \tp{\jac_{h_k}^u\!(\barxuk)} \partial^2\ell_k(\hbar^k)\, \jac_{h_k}^x\!(\barxuk) \), keeping only rows in the mask \( \mathcal J \), added to the given matrix S. \( S \leftarrow S + S_k[\mathcal J, \cdot] \). The size of work should be get_S_work_size().

inline void eval_add_R_prod_masked(index_t timestep, crvec xu, crvec h, crindexvec mask_J, crindexvec mask_K, crvec v, rvec out, rvec work) const#

\( out \leftarrow out + R[\mathcal J, \mathcal K]\,v[\mathcal K] \).

Work should contain the contents written to it by a prior call to eval_add_R_masked() in the same point.

inline void eval_add_S_prod_masked(index_t timestep, crvec xu, crvec h, crindexvec mask_K, crvec v, rvec out, rvec work) const#

\( out \leftarrow out + \tp{S[\mathcal K, \cdot]}\, v[\mathcal K] \).

Work should contain the contents written to it by a prior call to eval_add_S_masked() in the same point.

inline length_t get_R_work_size() const#

Size of the workspace required by eval_add_R_masked() and eval_add_R_prod_masked().

inline length_t get_S_work_size() const#

Size of the workspace required by eval_add_S_masked() and eval_add_S_prod_masked().

Constraints

inline void eval_constr(index_t timestep, crvec x, rvec c) const#

Stage constraints \( c_k(x^k) \).

inline void eval_constr_N(crvec x, rvec c) const#

Terminal constraints \( c_N(x^N) \).

inline void eval_grad_constr_prod(index_t timestep, crvec x, crvec p, rvec grad_cx_p) const#

Gradient-vector product of stage constraints \( \nabla c_k(x^k)\, p \).

inline void eval_grad_constr_prod_N(crvec x, crvec p, rvec grad_cx_p) const#

Gradient-vector product of terminal constraints \( \nabla c_N(x^N)\, p \).

inline void eval_add_gn_hess_constr(index_t timestep, crvec x, crvec M, rmat out) const#

Gauss-Newton Hessian of stage constraints \( \tp{\jac_{c_k}}(x^k)\, \operatorname{diag}(M)\; \jac_{c_k}(x^k) \).

inline void eval_add_gn_hess_constr_N(crvec x, crvec M, rmat out) const#

Gauss-Newton Hessian of terminal constraints \( \tp{\jac_{c_N}}(x^N)\, \operatorname{diag}(M)\; \jac_{c_N}(x^N) \).

Checks

inline void check() const#

Check that the problem formulation is well-defined, the dimensions match, etc.

Throws an exception if this is not the case.

Querying specialized implementations

inline bool provides_get_D() const#
inline bool provides_get_D_N() const#
inline bool provides_eval_h() const#
inline bool provides_eval_h_N() const#
inline bool provides_eval_add_Q_N() const#
inline bool provides_eval_add_R_prod_masked() const#
inline bool provides_eval_add_S_prod_masked() const#
inline bool provides_get_R_work_size() const#
inline bool provides_get_S_work_size() const#
inline bool provides_eval_constr() const#
inline bool provides_eval_constr_N() const#
inline bool provides_eval_grad_constr_prod() const#
inline bool provides_eval_grad_constr_prod_N() const#
inline bool provides_eval_add_gn_hess_constr() const#
inline bool provides_eval_add_gn_hess_constr_N() const#

Public Types

using VTable = ControlProblemVTable<config_t>#
using allocator_type = Allocator#
using Box = typename VTable::Box#
using Dim = OCPDim<config_t>#
using TypeErased = util::TypeErased<VTable, allocator_type>#

Public Functions

TypeErased() noexcept(noexcept(allocator_type())) = default

Default constructor.

template<class Alloc>
inline TypeErased(std::allocator_arg_t, const Alloc &alloc)#

Default constructor (allocator aware).

inline TypeErased(const TypeErased &other)

Copy constructor.

inline TypeErased(const TypeErased &other, allocator_type alloc)

Copy constructor (allocator aware).

inline TypeErased(TypeErased &&other) noexcept

Move constructor.

inline TypeErased(TypeErased &&other, const allocator_type &alloc) noexcept

Move constructor (allocator aware).

template<class T, class Alloc>
inline explicit TypeErased(std::allocator_arg_t, const Alloc &alloc, T &&d)#

Main constructor that type-erases the given argument.

template<class T, class Alloc, class ...Args>
inline explicit TypeErased(std::allocator_arg_t, const Alloc &alloc, te_in_place_t<T>, Args&&... args)#

Main constructor that type-erases the object constructed from the given argument.

template<class T>
inline explicit TypeErased(T &&d)#

Requirement prevents this constructor from taking precedence over the copy and move constructors.

template<class T, class ...Args>
inline explicit TypeErased(te_in_place_t<T>, Args&&... args)#

Main constructor that type-erases the object constructed from the given argument.

Public Static Functions

template<class T, class ...Args>
static inline TypeErasedControlProblem make(Args&&... args)#
template<class Problem>
struct ProblemWithCounters#
#include <alpaqa/problem/problem-with-counters.hpp>

Problem wrapper that keeps track of the number of evaluations and the run time of each function.

You probably want to use problem_with_counters or problem_with_counters_ref instead of instantiating this class directly.

Note

The evaluation counters are stored using a std::shared_pointers, which means that different copies of a ProblemWithCounters instance all share the same counters. To opt out of this behavior, you can use the decouple_evaluations function.

Public Types

using Box = typename TypeErasedProblem<config_t>::Box#
using Sparsity = sparsity::Sparsity<config_t>#

Public Functions

inline void eval_proj_diff_g(crvec z, rvec e) const#
inline void eval_proj_multipliers(rvec y, real_t M) const#
inline real_t eval_prox_grad_step(real_t γ, crvec x, crvec grad_ψ, rvec x_hat, rvec p) const#
inline index_t eval_inactive_indices_res_lna(real_t γ, crvec x, crvec grad_ψ, rindexvec J) const#
inline real_t eval_f(crvec x) const#
inline void eval_grad_f(crvec x, rvec grad_fx) const#
inline void eval_g(crvec x, rvec gx) const#
inline void eval_grad_g_prod(crvec x, crvec y, rvec grad_gxy) const#
inline void eval_grad_gi(crvec x, index_t i, rvec grad_gi) const#
inline void eval_jac_g(crvec x, rvec J_values) const#
inline Sparsity get_jac_g_sparsity() const#
inline void eval_hess_L_prod(crvec x, crvec y, real_t scale, crvec v, rvec Hv) const#
inline void eval_hess_L(crvec x, crvec y, real_t scale, rvec H_values) const#
inline Sparsity get_hess_L_sparsity() const#
inline void eval_hess_ψ_prod(crvec x, crvec y, crvec Σ, real_t scale, crvec v, rvec Hv) const#
inline void eval_hess_ψ(crvec x, crvec y, crvec Σ, real_t scale, rvec H_values) const#
inline Sparsity get_hess_ψ_sparsity() const#
inline real_t eval_f_grad_f(crvec x, rvec grad_fx) const#
inline real_t eval_f_g(crvec x, rvec g) const#
inline void eval_grad_f_grad_g_prod(crvec x, crvec y, rvec grad_f, rvec grad_gxy) const#
inline void eval_grad_L(crvec x, crvec y, rvec grad_L, rvec work_n) const#
inline real_t eval_ψ(crvec x, crvec y, crvec Σ, rvec ŷ) const#
inline void eval_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const#
inline real_t eval_ψ_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const#
inline const Box &get_box_C() const#
inline const Box &get_box_D() const#
inline void check() const#
inline bool provides_eval_grad_gi() const#
inline bool provides_eval_inactive_indices_res_lna() const#
inline bool provides_eval_jac_g() const#
inline bool provides_get_jac_g_sparsity() const#
inline bool provides_eval_hess_L_prod() const#
inline bool provides_eval_hess_L() const#
inline bool provides_get_hess_L_sparsity() const#
inline bool provides_eval_hess_ψ_prod() const#
inline bool provides_eval_hess_ψ() const#
inline bool provides_get_hess_ψ_sparsity() const#
inline bool provides_eval_f_grad_f() const#
inline bool provides_eval_f_g() const#
inline bool provides_eval_grad_f_grad_g_prod() const#
inline bool provides_eval_grad_L() const#
inline bool provides_eval_ψ() const#
inline bool provides_eval_grad_ψ() const#
inline bool provides_eval_ψ_grad_ψ() const#
inline bool provides_get_box_C() const#
inline bool provides_get_box_D() const#
inline bool provides_check() const#
inline length_t get_n() const#
inline length_t get_m() const#
ProblemWithCounters() = default#
template<class P>
inline explicit ProblemWithCounters(P &&problem)#
template<class ...Args>
inline explicit ProblemWithCounters(std::in_place_t, Args&&... args)#
inline void reset_evaluations()#

Reset all evaluation counters and timers to zero.

Affects all instances that share the same evaluations. If you only want to reset the counters of this instance, use decouple_evaluations first.

inline void decouple_evaluations()#

Give this instance its own evaluation counters and timers, decoupling it from any other instances they might have previously been shared with.

The evaluation counters and timers are preserved (a copy is made).

Public Members

std::shared_ptr<EvalCounter> evaluations = std::make_shared<EvalCounter>()#
Problem problem#
template<Config Conf = DefaultConfig, class Allocator = std::allocator<std::byte>>
class TypeErasedProblem : public alpaqa::util::TypeErased<ProblemVTable<DefaultConfig>, std::allocator<std::byte>>#
#include <alpaqa/problem/type-erased-problem.hpp>

The main polymorphic minimization problem interface.

This class wraps the actual problem implementation class, filling in the missing member functions with sensible defaults, and providing a uniform interface that is used by the solvers.

The problem implementations do not inherit from an abstract base class. Instead, structural typing is used. The ProblemVTable constructor uses reflection to discover which member functions are provided by the problem implementation. See Problem formulations for more information, and C++/CustomCppProblem/main.cpp for an example.

Problem dimensions

length_t get_n() const#

[Required] Number of decision variables.

length_t get_m() const#

[Required] Number of constraints.

Required cost and constraint functions

real_t eval_f(crvec x) const#

[Required] Function that evaluates the cost, \( f(x) \)

Parameters:

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

void eval_grad_f(crvec x, rvec grad_fx) const#

[Required] Function that evaluates the gradient of the cost, \( \nabla f(x) \)

Parameters:
  • x[in] Decision variable \( x \in \R^n \)

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

void eval_g(crvec x, rvec gx) const#

[Required] Function that evaluates the constraints, \( g(x) \)

Parameters:
  • x[in] Decision variable \( x \in \R^n \)

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

void eval_grad_g_prod(crvec x, crvec y, rvec grad_gxy) const#

[Required] Function that evaluates the gradient of the constraints times a vector, \( \nabla g(x)\,y = \tp{\jac_g(x)}y \)

Parameters:
  • x[in] Decision variable \( x \in \R^n \)

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

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

Projections onto constraint sets and proximal mappings

void eval_proj_diff_g(crvec z, rvec e) const#

[Required] Function that evaluates the difference between the given point \( z \) and its projection onto the constraint set \( D \).

Note

z and e can refer to the same vector.

Parameters:
  • z[in] Slack variable, \( z \in \R^m \)

  • e[out] The difference relative to its projection, \( e = z - \Pi_D(z) \in \R^m \)

void eval_proj_multipliers(rvec y, real_t M) const#

[Required] Function that projects the Lagrange multipliers for ALM.

Parameters:
  • y[inout] Multipliers, \( y \leftarrow \Pi_Y(y) \in \R^m \)

  • M[in] The radius/size of the set \( Y \). See ALMParams::max_multiplier.

real_t eval_prox_grad_step(real_t γ, crvec x, crvec grad_ψ, rvec x_hat, rvec p) const#

[Required] Function that computes a proximal gradient step.

Note

The vector \( p \) is often used in stopping criteria, so its numerical accuracy is more important than that of \( \hat x \).

Parameters:
  • γ[in] Step size, \( \gamma \in \R_{>0} \)

  • x[in] Decision variable \( x \in \R^n \)

  • grad_ψ[in] Gradient of the subproblem cost, \( \nabla\psi(x) \in \R^n \)

  • [out] Next proximal gradient iterate, \( \hat x = T_\gamma(x) = \prox_{\gamma h}(x - \gamma\nabla\psi(x)) \in \R^n \)

  • p[out] The proximal gradient step, \( p = \hat x - x \in \R^n \)

Returns:

The nonsmooth function evaluated at x̂, \( h(\hat x) \).

index_t eval_inactive_indices_res_lna(real_t γ, crvec x, crvec grad_ψ, rindexvec J) const#

[Optional] Function that computes the inactive indices \( \mathcal J(x) \) for the evaluation of the linear Newton approximation of the residual, as in [4].

For example, in the case of box constraints, we have

\[ \mathcal J(x) \defeq \defset{i \in \N_{[0, n-1]}}{\underline x_i \lt x_i - \gamma\nabla_{\!x_i}\psi(x) \lt \overline x_i}. \]

Parameters:
  • γ[in] Step size, \( \gamma \in \R_{>0} \)

  • x[in] Decision variable \( x \in \R^n \)

  • grad_ψ[in] Gradient of the subproblem cost, \( \nabla\psi(x) \in \R^n \)

  • J[out] The indices of the components of \( x \) that are in the index set \( \mathcal J(x) \). In ascending order, at most n.

Returns:

The number of inactive constraints, \( \# \mathcal J(x) \).

Constraint sets

const Box &get_box_C() const#

[Optional] Get the rectangular constraint set of the decision variables, \( x \in C \).

const Box &get_box_D() const#

[Optional] Get the rectangular constraint set of the general constraint function, \( g(x) \in D \).

Functions for second-order solvers

void eval_jac_g(crvec x, rvec J_values) const#

[Optional] Function that evaluates the nonzero values of the Jacobian matrix of the constraints, \( \jac_g(x) \)

Required for second-order solvers only.

Parameters:
  • x[in] Decision variable \( x \in \R^n \)

  • J_values[out] Nonzero values of the Jacobian \( \jac_g(x) \in \R^{m\times n} \)

Sparsity get_jac_g_sparsity() const#

[Optional] Function that returns (a view of) the sparsity pattern of the Jacobian of the constraints.

Required for second-order solvers only.

void eval_grad_gi(crvec x, index_t i, rvec grad_gi) const#

[Optional] Function that evaluates the gradient of one specific constraint, \( \nabla g_i(x) \)

Required for second-order solvers only.

Parameters:
  • x[in] Decision variable \( x \in \R^n \)

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

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

void eval_hess_L_prod(crvec x, crvec y, real_t scale, crvec v, rvec Hv) const#

[Optional] Function that evaluates the Hessian of the Lagrangian multiplied by a vector, \( \nabla_{xx}^2L(x, y)\,v \)

Required for second-order solvers only.

Parameters:
  • x[in] Decision variable \( x \in \R^n \)

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

  • scale[in] Scale factor for the cost function.

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

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

void eval_hess_L(crvec x, crvec y, real_t scale, rvec H_values) const#

[Optional] Function that evaluates the nonzero values of the Hessian of the Lagrangian, \( \nabla_{xx}^2L(x, y) \)

Required for second-order solvers only.

Parameters:
  • x[in] Decision variable \( x \in \R^n \)

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

  • scale[in] Scale factor for the cost function.

  • H_values[out] Nonzero values of the Hessian \( \nabla_{xx}^2 L(x, y) \in \R^{n\times n} \).

Sparsity get_hess_L_sparsity() const#

[Optional] Function that returns (a view of) the sparsity pattern of the Hessian of the Lagrangian.

Required for second-order solvers only.

void eval_hess_ψ_prod(crvec x, crvec y, crvec Σ, real_t scale, crvec v, rvec Hv) const#

[Optional] Function that evaluates the Hessian of the augmented Lagrangian multiplied by a vector, \( \nabla_{xx}^2L_\Sigma(x, y)\,v \)

Required for second-order solvers only.

Parameters:
  • x[in] Decision variable \( x \in \R^n \)

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

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

  • scale[in] Scale factor for the cost function.

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

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

void eval_hess_ψ(crvec x, crvec y, crvec Σ, real_t scale, rvec H_values) const#

[Optional] Function that evaluates the nonzero values of the Hessian of the augmented Lagrangian, \( \nabla_{xx}^2L_\Sigma(x, y) \)

Required for second-order solvers only.

Parameters:
  • x[in] Decision variable \( x \in \R^n \)

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

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

  • scale[in] Scale factor for the cost function.

  • H_values[out] Nonzero values of the Hessian \( \nabla_{xx}^2 L_\Sigma(x, y) \in \R^{n\times n} \)

Sparsity get_hess_ψ_sparsity() const#

[Optional] Function that returns (a view of) the sparsity pattern of the Hessian of the augmented Lagrangian.

Required for second-order solvers only.

Combined evaluations

real_t eval_f_grad_f(crvec x, rvec grad_fx) const#

[Optional] Evaluate both \( f(x) \) and its gradient, \( \nabla f(x) \).

Default implementation:

ProblemVTable::default_eval_f_grad_f

real_t eval_f_g(crvec x, rvec g) const#

[Optional] Evaluate both \( f(x) \) and \( g(x) \).

Default implementation:

ProblemVTable::default_eval_f_g

void eval_grad_f_grad_g_prod(crvec x, crvec y, rvec grad_f, rvec grad_gxy) const#

[Optional] Evaluate both \( \nabla f(x) \) and \( \nabla g(x)\,y \).

Default implementation:

ProblemVTable::default_eval_grad_f_grad_g_prod

void eval_grad_L(crvec x, crvec y, rvec grad_L, rvec work_n) const#

[Optional] Evaluate the gradient of the Lagrangian \( \nabla_x L(x, y) = \nabla f(x) + \nabla g(x)\,y \)

Default implementation:

ProblemVTable::default_eval_grad_L

Augmented Lagrangian

real_t eval_ψ(crvec x, crvec y, crvec Σ, rvec ŷ) const#

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

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

ProblemVTable::default_eval_ψ

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

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

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

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

void eval_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const#

[Optional] Calculate the gradient ∇ψ(x).

\[ \nabla \psi(x) = \nabla f(x) + \nabla g(x)\,\hat y(x) \]
Default implementation:

ProblemVTable::default_eval_grad_ψ

Parameters:
  • 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 \)

real_t eval_ψ_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const#

[Optional] Calculate both ψ(x) and its gradient ∇ψ(x).

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

ProblemVTable::default_eval_ψ_grad_ψ

Parameters:
  • 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 \)

Checks

void check() const#

[Optional] Check that the problem formulation is well-defined, the dimensions match, etc.

Throws an exception if this is not the case.

Querying specialized implementations

inline bool provides_eval_inactive_indices_res_lna() const#

Returns true if the problem provides an implementation of eval_inactive_indices_res_lna.

inline bool provides_eval_jac_g() const#

Returns true if the problem provides an implementation of eval_jac_g.

inline bool provides_get_jac_g_sparsity() const#

Returns true if the problem provides an implementation of get_jac_g_sparsity.

inline bool provides_eval_grad_gi() const#

Returns true if the problem provides an implementation of eval_grad_gi.

inline bool provides_eval_hess_L_prod() const#

Returns true if the problem provides an implementation of eval_hess_L_prod.

inline bool provides_eval_hess_L() const#

Returns true if the problem provides an implementation of eval_hess_L.

inline bool provides_get_hess_L_sparsity() const#

Returns true if the problem provides an implementation of get_hess_L_sparsity.

inline bool provides_eval_hess_ψ_prod() const#

Returns true if the problem provides an implementation of eval_hess_ψ_prod.

inline bool provides_eval_hess_ψ() const#

Returns true if the problem provides an implementation of eval_hess_ψ.

inline bool provides_get_hess_ψ_sparsity() const#

Returns true if the problem provides an implementation of get_hess_ψ_sparsity.

inline bool provides_eval_f_grad_f() const#

Returns true if the problem provides a specialized implementation of eval_f_grad_f, false if it uses the default implementation.

inline bool provides_eval_f_g() const#

Returns true if the problem provides a specialized implementation of eval_f_g, false if it uses the default implementation.

inline bool provides_eval_grad_f_grad_g_prod() const#

Returns true if the problem provides a specialized implementation of eval_grad_f_grad_g_prod, false if it uses the default implementation.

inline bool provides_eval_grad_L() const#

Returns true if the problem provides a specialized implementation of eval_grad_L, false if it uses the default implementation.

inline bool provides_eval_ψ() const#

Returns true if the problem provides a specialized implementation of eval_ψ, false if it uses the default implementation.

inline bool provides_eval_grad_ψ() const#

Returns true if the problem provides a specialized implementation of eval_grad_ψ, false if it uses the default implementation.

inline bool provides_eval_ψ_grad_ψ() const#

Returns true if the problem provides a specialized implementation of eval_ψ_grad_ψ, false if it uses the default implementation.

inline bool provides_get_box_C() const#

Returns true if the problem provides an implementation of get_box_C.

inline bool provides_get_box_D() const#

Returns true if the problem provides an implementation of get_box_D.

inline bool provides_check() const#

Returns true if the problem provides an implementation of check.

Helpers

real_t calc_ŷ_dᵀŷ(rvec g_ŷ, crvec y, crvec Σ) const#

Given g(x), compute the intermediate results ŷ and dᵀŷ that can later be used to compute ψ(x) and ∇ψ(x).

Computes the result using the following algorithm:

\[\begin{split} \begin{aligned} \zeta &= g(x) + \Sigma^{-1} y \\[] d &= \zeta - \Pi_D(\zeta) = \operatorname{eval\_proj\_diff\_g}(\zeta, \zeta) \\[] \hat y &= \Sigma d \\[] \end{aligned} \end{split}\]

See also

page_math

Parameters:
  • g_ŷ[inout] Input \( g(x) \), outputs \( \hat y \)

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

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

Returns:

The inner product \( d^\top \hat y \)

Public Types

using Sparsity = alpaqa::Sparsity<config_t>#
using Box = alpaqa::Box<config_t>#
using VTable = ProblemVTable<config_t>#
using allocator_type = Allocator#
using TypeErased = util::TypeErased<VTable, allocator_type>#

Public Functions

TypeErased() noexcept(noexcept(allocator_type())) = default

Default constructor.

template<class Alloc>
inline TypeErased(std::allocator_arg_t, const Alloc &alloc)#

Default constructor (allocator aware).

inline TypeErased(const TypeErased &other)

Copy constructor.

inline TypeErased(const TypeErased &other, allocator_type alloc)

Copy constructor (allocator aware).

inline TypeErased(TypeErased &&other) noexcept

Move constructor.

inline TypeErased(TypeErased &&other, const allocator_type &alloc) noexcept

Move constructor (allocator aware).

template<class T, class Alloc>
inline explicit TypeErased(std::allocator_arg_t, const Alloc &alloc, T &&d)#

Main constructor that type-erases the given argument.

template<class T, class Alloc, class ...Args>
inline explicit TypeErased(std::allocator_arg_t, const Alloc &alloc, te_in_place_t<T>, Args&&... args)#

Main constructor that type-erases the object constructed from the given argument.

template<class T>
inline explicit TypeErased(T &&d)#

Requirement prevents this constructor from taking precedence over the copy and move constructors.

template<class T, class ...Args>
inline explicit TypeErased(te_in_place_t<T>, Args&&... args)#

Main constructor that type-erases the object constructed from the given argument.

Public Static Functions

template<class T, class ...Args>
static inline TypeErasedProblem make(Args&&... args)#
template<Config Conf>
class UnconstrProblem#
#include <alpaqa/problem/unconstr-problem.hpp>

Implements common problem functions for minimization problems without constraints.

Meant to be used as a base class for custom problem implementations.

Public Functions

inline UnconstrProblem(length_t n)#
Parameters:

n – Number of decision variables

inline void resize(length_t n)#

Change the number of decision variables.

UnconstrProblem(const UnconstrProblem&) = default#
UnconstrProblem &operator=(const UnconstrProblem&) = default#
UnconstrProblem(UnconstrProblem&&) noexcept = default#
UnconstrProblem &operator=(UnconstrProblem&&) noexcept = default#
inline length_t get_n() const#

Number of decision variables, n.

inline length_t get_m() const#

Number of constraints (always zero)

inline void eval_g(crvec, rvec) const#

No-op, no constraints.

inline void eval_grad_g_prod(crvec, crvec, rvec grad) const#

Constraint gradient is always zero.

inline void eval_jac_g(crvec, rvec) const#

Constraint Jacobian is always empty.

inline void eval_grad_gi(crvec, index_t, rvec grad_gi) const#

Constraint gradient is always zero.

inline real_t eval_prox_grad_step(real_t γ, crvec x, crvec grad_ψ, rvec x_hat, rvec p) const#

No proximal mapping, just a forward (gradient) step.

inline void eval_proj_diff_g(crvec, rvec) const#

inline void eval_proj_multipliers(rvec, real_t) const#

inline index_t eval_inactive_indices_res_lna(real_t, crvec, crvec, rindexvec J) const#

Public Members

length_t n#

Number of decision variables, dimension of x.

template<Config Conf = EigenConfigd>
class CasADiProblem : public alpaqa::BoxConstrProblem<EigenConfigd>#
#include <alpaqa/casadi/CasADiProblem.hpp>

Problem definition for a CasADi problem, loaded from a DLL.

Public Types

using Sparsity = alpaqa::Sparsity<config_t>#

Public Functions

CasADiProblem(const std::string &filename)#

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 respectively. For second order solvers, additional functions hess_L, hess_ψ, hess_L_prod and hess_ψ_prod can be provided to evaluate the Hessian of the (augmented) Lagrangian and Hessian-vector products.

Parameters:

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

Throws:

std::invalid_argument – The dimensions of the loaded functions do not match.

CasADiProblem(const SerializedCasADiFunctions &functions)#

Create a problem from a collection of serialized CasADi functions.

~CasADiProblem()#
CasADiProblem(const CasADiProblem&)#
CasADiProblem &operator=(const CasADiProblem&)#
CasADiProblem(CasADiProblem&&) noexcept#
CasADiProblem &operator=(CasADiProblem&&) noexcept#
void load_numerical_data(const std::filesystem::path &filepath, char sep = ',')#

Load the numerical problem data (bounds and parameters) from a CSV file.

The file should contain 7 rows, with the following contents:

  1. C lower bound [n]

  2. C upper bound [n]

  3. D lower bound [m]

  4. D upper bound [m]

  5. param [p]

  6. l1_reg [0, 1 or n]

  7. penalty_alm_split [1]

Line endings are encoded using a single line feed (\n), and the column separator can be specified using the sep argument.

real_t eval_f(crvec x) const#
void eval_grad_f(crvec x, rvec grad_fx) const#
real_t eval_f_grad_f(crvec x, rvec grad_fx) const#
void eval_g(crvec x, rvec g) const#
void eval_grad_g_prod(crvec x, crvec y, rvec grad_gxy) const#
void eval_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const#
real_t eval_ψ_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const#
void eval_grad_L(crvec x, crvec y, rvec grad_L, rvec work_n) const#
real_t eval_ψ(crvec x, crvec y, crvec Σ, rvec ŷ) const#
void eval_grad_gi(crvec x, index_t i, rvec grad_i) const#
Sparsity get_jac_g_sparsity() const#
void eval_jac_g(crvec x, rvec J_values) const#
void eval_hess_L_prod(crvec x, crvec y, real_t scale, crvec v, rvec Hv) const#
Sparsity get_hess_L_sparsity() const#
void eval_hess_L(crvec x, crvec y, real_t scale, rvec H_values) const#
void eval_hess_ψ_prod(crvec x, crvec y, crvec Σ, real_t scale, crvec v, rvec Hv) const#
Sparsity get_hess_ψ_sparsity() const#
void eval_hess_ψ(crvec x, crvec y, crvec Σ, real_t scale, rvec H_values) const#
bool provides_eval_grad_gi() const#

bool provides_eval_jac_g() const#

bool provides_eval_hess_L_prod() const#

bool provides_eval_hess_L() const#

bool provides_eval_hess_ψ_prod() const#

bool provides_eval_hess_ψ() const#

Public Members

vec param#
class CUTEstProblem : public BoxConstrProblem<alpaqa::EigenConfigd>#
#include <alpaqa/cutest/cutest-loader.hpp>

Wrapper for CUTEst problems loaded from an external shared library.

Public Types

using Sparsity = alpaqa::Sparsity<config_t>#

Public Functions

CUTEstProblem(const char *so_fname, const char *outsdif_fname = nullptr, bool sparse = false)#

Load a CUTEst problem from the given shared library and OUTSDIF.d file.

CUTEstProblem(const CUTEstProblem&)#
CUTEstProblem &operator=(const CUTEstProblem&)#
CUTEstProblem(CUTEstProblem&&) noexcept#
CUTEstProblem &operator=(CUTEstProblem&&) noexcept#
~CUTEstProblem()#
Report get_report() const#
std::ostream &format_report(std::ostream &os, const Report &r) const#
inline std::ostream &format_report(std::ostream &os) const#
real_t eval_f(crvec x) const#
void eval_grad_f(crvec x, rvec grad_fx) const#
void eval_g(crvec x, rvec gx) const#
void eval_grad_g_prod(crvec x, crvec y, rvec grad_gxy) const#
void eval_jac_g(crvec x, rvec J_values) const#
Sparsity get_jac_g_sparsity() const#
void eval_grad_gi(crvec x, index_t i, rvec grad_gi) const#
void eval_hess_L_prod(crvec x, crvec y, real_t scale, crvec v, rvec Hv) const#
void eval_hess_L(crvec x, crvec y, real_t scale, rvec H_values) const#
Sparsity get_hess_L_sparsity() const#
void eval_hess_ψ_prod(crvec x, crvec y, crvec Σ, real_t scale, crvec v, rvec Hv) const#
real_t eval_f_grad_f(crvec x, rvec grad_fx) const#
real_t eval_f_g(crvec x, rvec g) const#
void eval_grad_L(crvec x, crvec y, rvec grad_L, rvec work_n) const#

Public Members

std::string name = "<UNKNOWN>"#

Problem name.

vec x0#

Initial value of decision variables.

vec y0#

Initial value of Lagrange multipliers.

struct Report#
#include <alpaqa/cutest/cutest-loader.hpp>

The report generated by CUTEst.

Public Members

Calls calls#

Function call counters.

double time_setup = 0#

CPU time (in seconds) for CUTEST_csetup.

double time = 0#

CPU time (in seconds) since the end of CUTEST_csetup.

struct Calls#
#include <alpaqa/cutest/cutest-loader.hpp>

Function call counters.

Note

Note that hessian_times_vector, constraints and constraints_grad may account for codes which allow the evaluation of a selection of constraints only and may thus be much smaller than the number of constraints times the number of iterations.

Public Members

unsigned objective = 0#

Number of calls to the objective function.

unsigned objective_grad = 0#

Number of calls to the objective gradient.

unsigned objective_hess = 0#

Number of calls to the objective Hessian.

unsigned hessian_times_vector = 0#

Number of Hessian times vector products.

unsigned constraints = 0#

Number of calls to the constraint functions.

unsigned constraints_grad = 0#

Number of calls to the constraint gradients.

unsigned constraints_hess = 0#

Number of calls to the constraint Hessians.

class DLProblem : public BoxConstrProblem<DefaultConfig>#
#include <alpaqa/dl/dl-problem.hpp>

Class that loads a problem using dlopen.

The shared library should export a C function with the name function_name that accepts a void pointer with user data, and returns a struct of type alpaqa_problem_register_t that contains all data to represent the problem, as well as function pointers for all required operations. See C++/DLProblem/main.cpp and problems/sparse-logistic-regression.cpp for examples.

See also

alpaqa_problem_functions_t

See also

alpaqa_problem_register_t

Note

Copies are shallow, they all share the same problem instance, take that into account when using multiple threads.

Public Types

using Sparsity = alpaqa::Sparsity<config_t>#
using instance_t = ExtraFuncs::instance_t#

Public Functions

DLProblem(const std::string &so_filename, const std::string &function_name = "register_alpaqa_problem", void *user_param = nullptr)#

Load a problem from a shared library.

Parameters:
  • so_filename – Filename of the shared library to load.

  • function_name – Name of the problem registration function. Should have signature alpaqa_problem_register_t(void *).

  • user_param – Pointer to custom user data to pass to the registration function.

real_t eval_prox_grad_step(real_t γ, crvec x, crvec grad_ψ, rvec x_hat, rvec p) const#
real_t eval_f(crvec x) const#
void eval_grad_f(crvec x, rvec grad_fx) const#
void eval_g(crvec x, rvec gx) const#
void eval_grad_g_prod(crvec x, crvec y, rvec grad_gxy) const#
void eval_jac_g(crvec x, rvec J_values) const#
Sparsity get_jac_g_sparsity() const#
void eval_grad_gi(crvec x, index_t i, rvec grad_gi) const#
void eval_hess_L_prod(crvec x, crvec y, real_t scale, crvec v, rvec Hv) const#
void eval_hess_L(crvec x, crvec y, real_t scale, rvec H_values) const#
Sparsity get_hess_L_sparsity() const#
void eval_hess_ψ_prod(crvec x, crvec y, crvec Σ, real_t scale, crvec v, rvec Hv) const#
void eval_hess_ψ(crvec x, crvec y, crvec Σ, real_t scale, rvec H_values) const#
Sparsity get_hess_ψ_sparsity() const#
real_t eval_f_grad_f(crvec x, rvec grad_fx) const#
real_t eval_f_g(crvec x, rvec g) const#
void eval_grad_f_grad_g_prod(crvec x, crvec y, rvec grad_f, rvec grad_gxy) const#
void eval_grad_L(crvec x, crvec y, rvec grad_L, rvec work_n) const#
real_t eval_ψ(crvec x, crvec y, crvec Σ, rvec ŷ) const#
void eval_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const#
real_t eval_ψ_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const#
bool provides_eval_f() const#
bool provides_eval_grad_f() const#
bool provides_eval_g() const#
bool provides_eval_grad_g_prod() const#
bool provides_eval_jac_g() const#
bool provides_get_jac_g_sparsity() const#
bool provides_eval_grad_gi() const#
bool provides_eval_hess_L_prod() const#
bool provides_eval_hess_L() const#
bool provides_get_hess_L_sparsity() const#
bool provides_eval_hess_ψ_prod() const#
bool provides_eval_hess_ψ() const#
bool provides_get_hess_ψ_sparsity() const#
bool provides_eval_f_grad_f() const#
bool provides_eval_f_g() const#
bool provides_eval_grad_f_grad_g_prod() const#
bool provides_eval_grad_L() const#
bool provides_eval_ψ() const#
bool provides_eval_grad_ψ() const#
bool provides_eval_ψ_grad_ψ() const#
bool provides_get_box_C() const#
template<class Signature, class ...Args>
inline decltype(auto) call_extra_func(const std::string &name, Args&&... args) const#
template<class Signature, class ...Args>
inline decltype(auto) call_extra_func(const std::string &name, Args&&... args)#
class DLControlProblem#
#include <alpaqa/dl/dl-problem.hpp>

Class that loads an optimal control problem using dlopen.

The shared library should export a C function with the name function_name that accepts a void pointer with user data, and returns a struct of type alpaqa_control_problem_register_t that contains all data to represent the problem, as well as function pointers for all required operations.

Note

Copies are shallow, they all share the same problem instance, take that into account when using multiple threads.

Public Types

using Box = alpaqa::Box<config_t>#
using instance_t = ExtraFuncs::instance_t#

Public Functions

DLControlProblem(const std::string &so_filename, const std::string &function_name = "register_alpaqa_control_problem", void *user_param = nullptr)#

Load a problem from a shared library.

Parameters:
  • so_filename – Filename of the shared library to load.

  • function_name – Name of the problem registration function. Should have signature alpaqa_control_problem_register_t(void *).

  • user_param – Pointer to custom user data to pass to the registration function.

inline length_t get_N() const#
inline length_t get_nx() const#
inline length_t get_nu() const#
inline length_t get_nh() const#
inline length_t get_nh_N() const#
inline length_t get_nc() const#
inline length_t get_nc_N() const#
inline void check() const#
void get_U(Box &U) const#
void get_D(Box &D) const#
void get_D_N(Box &D) const#
void get_x_init(rvec x_init) const#
void eval_f(index_t timestep, crvec x, crvec u, rvec fxu) const#
void eval_jac_f(index_t timestep, crvec x, crvec u, rmat J_fxu) const#
void eval_grad_f_prod(index_t timestep, crvec x, crvec u, crvec p, rvec grad_fxu_p) const#
void eval_h(index_t timestep, crvec x, crvec u, rvec h) const#
void eval_h_N(crvec x, rvec h) const#
real_t eval_l(index_t timestep, crvec h) const#
real_t eval_l_N(crvec h) const#
void eval_qr(index_t timestep, crvec xu, crvec h, rvec qr) const#
void eval_q_N(crvec x, crvec h, rvec q) const#
void eval_add_Q(index_t timestep, crvec xu, crvec h, rmat Q) const#
void eval_add_Q_N(crvec x, crvec h, rmat Q) const#
void eval_add_R_masked(index_t timestep, crvec xu, crvec h, crindexvec mask, rmat R, rvec work) const#
void eval_add_S_masked(index_t timestep, crvec xu, crvec h, crindexvec mask, rmat S, rvec work) const#
void eval_add_R_prod_masked(index_t timestep, crvec xu, crvec h, crindexvec mask_J, crindexvec mask_K, crvec v, rvec out, rvec work) const#
void eval_add_S_prod_masked(index_t timestep, crvec xu, crvec h, crindexvec mask_K, crvec v, rvec out, rvec work) const#
length_t get_R_work_size() const#
length_t get_S_work_size() const#
void eval_constr(index_t timestep, crvec x, rvec c) const#
void eval_constr_N(crvec x, rvec c) const#
void eval_grad_constr_prod(index_t timestep, crvec x, crvec p, rvec grad_cx_p) const#
void eval_grad_constr_prod_N(crvec x, crvec p, rvec grad_cx_p) const#
void eval_add_gn_hess_constr(index_t timestep, crvec x, crvec M, rmat out) const#
void eval_add_gn_hess_constr_N(crvec x, crvec M, rmat out) const#
bool provides_get_D() const#
bool provides_get_D_N() const#
bool provides_eval_add_Q_N() const#
bool provides_eval_add_R_prod_masked() const#
bool provides_eval_add_S_prod_masked() const#
bool provides_get_R_work_size() const#
bool provides_get_S_work_size() const#
bool provides_eval_constr() const#
bool provides_eval_constr_N() const#
bool provides_eval_grad_constr_prod() const#
bool provides_eval_grad_constr_prod_N() const#
bool provides_eval_add_gn_hess_constr() const#
bool provides_eval_add_gn_hess_constr_N() const#
template<class Signature, class ...Args>
inline decltype(auto) call_extra_func(const std::string &name, Args&&... args) const#
template<class Signature, class ...Args>
inline decltype(auto) call_extra_func(const std::string &name, Args&&... args)#