This example shows how to define optimization problems using Fortran routines.
This example shows how to define optimization problems using Fortran routines.
\[ \begin{aligned}
& \underset{x}{\text{minimize}} && \tfrac12 \tp x Q x \\
& \text{subject to} && Ax \le b \\
\end{aligned} \]
1module problem
2
3 use, intrinsic::iso_c_binding, only: c_double, c_ptrdiff_t
4 implicit none (type,external)
5
6 integer(c_ptrdiff_t), parameter :: n = 2
7 integer(c_ptrdiff_t), parameter :: m = 1
8 real(c_double) :: Q(n, n)
9 data q /3, -1, -1, 3/
10 real(c_double) :: A(m, n)
11 data a /2, 1/
12
13contains
14
15
16 pure real(c_double) function problem_eval_f(x) bind(C)
17 real(c_double), intent(in) :: x(n)
18
19 problem_eval_f = 0.5_c_double * dot_product(x, matmul(q, x))
20 end function
21
22
23 pure subroutine problem_eval_grad_f(x, grad_fx) bind(C)
24 real(c_double), intent(in) :: x(n)
25 real(c_double), intent(out) :: grad_fx(n)
26
27 grad_fx = matmul(q, x)
28 end subroutine
29
30
31 pure subroutine problem_eval_g(x, gx) bind(C)
32 real(c_double), intent(in) :: x(n)
33 real(c_double), intent(out) :: gx(m)
34
35 gx = matmul(a, x)
36 end subroutine
37
38
39
40 pure subroutine problem_eval_grad_g_prod(x, y, grad_gxy) bind(C)
41 real(c_double), intent(in) :: x(n)
42 real(c_double), intent(in) :: y(m)
43 real(c_double), intent(out) :: grad_gxy(n)
44
45 grad_gxy = matmul(transpose(a), y)
46 end subroutine
47
48
49 pure integer(c_ptrdiff_t) function problem_get_num_vars() bind(C)
50 problem_get_num_vars = n
51 end function
52
53
54 pure integer(c_ptrdiff_t) function problem_get_num_constr() bind(C)
55 problem_get_num_constr = m
56 end function
57
58end module
1#include <alpaqa/example-util.hpp>
13ptrdiff_t problem_get_num_vars(
void);
14ptrdiff_t problem_get_num_constr(
void);
15double problem_eval_f(
const double *);
16void problem_eval_grad_f(
const double *,
double *);
17void problem_eval_g(
const double *,
double *);
18void problem_eval_grad_g_prod(
const double *,
const double *,
double *);
26 void eval_grad_f(crvec x, rvec fx)
const {
29 void eval_g(crvec x, rvec gx)
const {
problem_eval_g(x.data(), gx.data()); }
30 void eval_grad_g_prod(crvec x, crvec y, rvec grad_gxy)
const {
39 FortranProblem problem{problem_get_num_vars(), problem_get_num_constr()};
44 vec b = vec::Constant(problem.m, -1);
48 problem.D.upperbound = b;
56 OuterSolver::Params almparam;
57 almparam.tolerance = 1e-8;
58 almparam.dual_tolerance = 1e-8;
59 almparam.penalty_update_factor = 10;
60 almparam.max_iter = 20;
61 almparam.print_interval = 1;
64 InnerSolver::Params panocparam;
65 panocparam.max_iter = 500;
66 panocparam.print_interval = 10;
68 Direction::AcceleratorParams lbfgsparam;
69 lbfgsparam.memory = 2;
74 {panocparam, lbfgsparam},
84 auto stats = solver(counted_problem, x, y);
88 std::cout <<
'\n' << *counted_problem.evaluations <<
'\n';
89 std::cout <<
"status: " << stats.status <<
'\n'
90 <<
"solver: " << solver.get_name() <<
'\n'
91 <<
"f = " << problem.eval_f(x) <<
'\n'
92 <<
"inner iterations: " << stats.inner.iterations <<
'\n'
93 <<
"outer iterations: " << stats.outer_iterations <<
'\n'
94 <<
"ε = " << stats.ε <<
'\n'
95 <<
"δ = " << stats.δ <<
'\n'
97 << std::chrono::duration<double>{stats.elapsed_time}.count()
99 <<
"x = " << x.transpose() <<
'\n'
100 <<
"y = " << y.transpose() <<
'\n'
101 <<
"avg τ = " << (stats.inner.sum_τ / stats.inner.count_τ) <<
'\n'
102 <<
"L-BFGS rejected = " << stats.inner.lbfgs_rejected <<
'\n'
103 <<
"L-BFGS failures = " << stats.inner.lbfgs_failures <<
'\n'
104 <<
"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.
#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.