This example shows how to define optimization problems using ordinary C++ functions.
This example shows how to define optimization problems using ordinary C++ functions.
It solves a simple quadratic program of the form:
\[ \begin{aligned}
& \underset{x}{\text{minimize}} && \tfrac12 \tp x Q x \\
& \text{subject to} && Ax \le b \\
\end{aligned} \]
Problem specification
The problem is specified by creating a class (Problem
) that inherits from alpaqa's alpaqa::BoxConstrProblem class. It defines the problem-specific functions for the evaluation of the cost function \( f(x) = \tfrac12 \tp xQx \) (eval_f) and its gradient \( \nabla f(x) = Qx \) (eval_grad_f), as well as the constraint function \( g(x) = Ax \) (eval_g) and the function that evaluates the product of the gradient of the constraints and a given vector \( y \), \( \nabla g(x)\,y = \tp A y \) (eval_grad_g_prod).
If you have more efficient ways to combine evaluations of these functions and gradients, you can specify them as well, see alpaqa::TypeErasedProblem for the full list of supported functions.
The alpaqa::BoxConstrProblem class exposes the two constraint sets, \( x \in C \) and \( g(x) \in D \), and provides the projection functions (eval_prox_grad_step, eval_proj_diff_g and eval_proj_multipliers) for you. The alpaqa::BoxConstrProblem constructor accepts the number of variables \( n = 2 \) and the number of constraints \( m = 1 \) as arguments.
- Note
- Alpaqa uses structural typing for problem definitions. This means that you just have to provide the supported functions with the correct names and arguments, and the library will pick them up automatically, you don't have to inherit from any abstract interfaces or override any virtual functions. See Problem formulations for more information.
- See also
- If you haven't already, be sure to go through the Problem formulations page.
Solver selection
The solver consists of three layers:
- The outer augmented Lagrangian solver (alpaqa::ALMSolver) that handles the general constraints \( g(x) \in D \);
- The inner PANOC solver (alpaqa::PANOCSolver) that is used to solve the ALM subproblems;
- The L-BFGS direction (alpaqa::LBFGSDirection) that provides fast Newton-type directions to speed up PANOC.
You can try out different inner solvers (Inner solvers) and direction providers (Direction providers), you can even write your own.
Solver configuration
Each solver and direction class has a set of parameters, which are collected in a struct. Tuning these parameters can often significantly improve solver performance. You can also use them to limit the run time or the number of iterations, and to enable verbose output from the solvers.
Solver invocation
The solver is invoked by calling the solver object with an instance of the problem and an initial guess for the Lagrange multipliers and the decision variables. The solver will overwrite these guesses with the (approximate) solution, and returns a struct containing solver statistics, the most important of which is the status, which reports whether convergence to the desired tolerance was achieved.
33 real_t eval_f(crvec x)
const {
35 return 0.5 * x.dot(Qx);
38 void eval_grad_f(crvec x, rvec gr)
const { gr.noalias() = Q * x; }
40 void eval_g(crvec x, rvec g)
const { g.noalias() = A * x; }
42 void eval_grad_g_prod(crvec x, crvec y, rvec gr)
const {
44 gr.noalias() = A.transpose() * y;
59 OuterSolver::Params almparam;
60 almparam.tolerance = 1e-8;
61 almparam.dual_tolerance = 1e-8;
62 almparam.penalty_update_factor = 10;
63 almparam.max_iter = 20;
64 almparam.print_interval = 1;
67 InnerSolver::Params panocparam;
68 panocparam.max_iter = 500;
69 panocparam.print_interval = 1;
71 Direction::AcceleratorParams lbfgsparam;
72 lbfgsparam.memory = 2;
77 {panocparam, lbfgsparam},
87 auto stats = solver(counted_problem, x, y);
91 std::cout <<
'\n' << *counted_problem.evaluations <<
'\n';
92 std::cout <<
"status: " << stats.status <<
'\n'
93 <<
"f = " << problem.eval_f(x) <<
'\n'
94 <<
"inner iterations: " << stats.inner.iterations <<
'\n'
95 <<
"outer iterations: " << stats.outer_iterations <<
'\n'
96 <<
"ε = " << stats.ε <<
'\n'
97 <<
"δ = " << stats.δ <<
'\n'
99 << std::chrono::duration<double>{stats.elapsed_time}.count()
101 <<
"x = " << x.transpose() <<
'\n'
102 <<
"y = " << y.transpose() <<
'\n'
103 <<
"avg τ = " << (stats.inner.sum_τ / stats.inner.count_τ) <<
'\n'
104 <<
"L-BFGS rejected = " << stats.inner.lbfgs_rejected <<
'\n'
105 <<
"L-BFGS failures = " << stats.inner.lbfgs_failures <<
'\n'
106 <<
"Line search failures = " << stats.inner.linesearch_failures
int main(int argc, const char *argv[])
Augmented Lagrangian Method solver.
Implements common problem functions for minimization problems with box constraints.
length_t m
Number of constraints, dimension of g(x) and z.
length_t n
Number of decision variables, dimension of x.
#define USING_ALPAQA_CONFIG(Conf)
auto problem_with_counters_ref(Problem &p)
Wraps the given problem into a ProblemWithCounters and keeps track of how many times each function is...
Double-precision double configuration.