Nonconvex constrained optimization
Loading...
Searching...
No Matches
gradient-checker.cpp
Go to the documentation of this file.
1#ifdef _WIN32
2#define NOMINMAX
3#include <Windows.h>
4#endif
5
10#include <alpaqa/util/print.hpp>
11#include <guanaqo/demangled-typename.hpp>
12#include <alpaqa-version.h>
13
14#include <Eigen/Sparse>
15
16#ifdef ALPAQA_WITH_EXTERNAL_CASADI
17#include <casadi/config.h>
18#endif
19#ifdef ALPAQA_WITH_IPOPT
20#include <IpoptConfig.h>
21#endif
22
23#include <algorithm>
24#include <filesystem>
25#include <fstream>
26#include <iostream>
27#include <random>
28#include <span>
29#include <stdexcept>
30#include <string>
31#include <string_view>
32#include <tuple>
33#include <type_traits>
34namespace fs = std::filesystem;
35
37
38const auto *docs = R"==(
39problem types:
40 dl: Dynamically loaded problem using the DLProblem class.
41 Specify the name of the registration function using the
42 problem.register option, e.g. problem.register=register_alpaqa_problem.
43 Further options can be passed to the problem using
44 problem.<key>[=<value>].
45 cs: Load a CasADi problem using the CasADiProblem class.
46 If a .tsv file with the same name as the shared library file exists,
47 the bounds and parameters will be loaded from that file. See
48 CasADiProblem::load_numerical_data for more details.
49 The problem parameter can be set using the problem.param option.
50 cu: Load a CUTEst problem using the CUTEstProblem class.
51
52options:
53 --full-print
54 Print the full gradients.
55 --no-hessians
56 Do not check any Hessian matrices.
57 --seed=<seed>
58 Seed for the random number generator.
59)==";
60
61void print_usage(const char *a0) {
62 const auto *opts = " [<problem-type>:][<path>/]<name> [options...]\n";
63 std::cout << "alpaqa gradient-checker (" ALPAQA_VERSION_FULL ")\n\n"
64 " Command-line utility to check problem gradients.\n"
65 " alpaqa is published under the LGPL-3.0.\n"
66 " https://github.com/kul-optec/alpaqa"
67 "\n\n"
68 " Usage: "
69 << a0 << opts << docs << std::endl;
70 std::cout << "Third-party libraries:\n"
71 << " * Eigen " << EIGEN_WORLD_VERSION << '.'
72 << EIGEN_MAJOR_VERSION << '.' << EIGEN_MINOR_VERSION
73 << " (https://gitlab.com/libeigen/eigen) - MPL-2.0\n"
74#ifdef ALPAQA_WITH_EXTERNAL_CASADI
75 << " * CasADi " CASADI_VERSION_STRING
76 " (https://github.com/casadi/casadi) - LGPL-3.0-or-later\n"
77#endif
78#ifdef ALPAQA_WITH_CUTEST
79 << " * CUTEst"
80 " (https://github.com/ralna/CUTEst) - BSD-3-Clause\n"
81#endif
82 << std::endl;
83}
84
85/// Split the string @p s on the first occurrence of @p tok.
86/// Returns ("", s) if tok was not found.
87auto split_once(std::string_view s, char tok = '.') {
88 auto tok_pos = s.find(tok);
89 if (tok_pos == s.npos)
90 return std::make_tuple(std::string_view{}, s);
91 std::string_view key{s.begin(), s.begin() + tok_pos};
92 std::string_view rem{s.begin() + tok_pos + 1, s.end()};
93 return std::make_tuple(key, rem);
94}
95
96auto get_problem_path(const char *const *argv) {
97 bool rel_to_exe = argv[1][0] == '^';
98 std::string_view prob_path_s = argv[1] + static_cast<ptrdiff_t>(rel_to_exe);
99 std::string_view prob_type;
100 std::tie(prob_type, prob_path_s) = split_once(prob_path_s, ':');
101 fs::path prob_path{prob_path_s};
102 if (rel_to_exe)
103 prob_path = fs::canonical(fs::path(argv[0])).parent_path() / prob_path;
104 return std::make_tuple(std::move(prob_path), prob_type);
105}
106
112
113void check_gradients(alpaqa::LoadedProblem &, std::ostream &,
114 const CheckGradientsOpts &);
115
116int main(int argc, const char *argv[]) try {
117#ifdef _WIN32
118 SetConsoleOutputCP(CP_UTF8);
119#endif
120 // Check command line options
121 if (argc < 1)
122 return -1;
123 if (argc == 1)
124 return print_usage(argv[0]), 0;
125 if (argc < 2)
126 return print_usage(argv[0]), -1;
127 std::span args{argv, static_cast<size_t>(argc)};
128 alpaqa::Options opts{argc - 2, argv + 2};
129
130 // Check where to write the output to
131 std::ostream &os = std::cout;
132
133 // Check which problem to load
134 auto [prob_path, prob_type] = get_problem_path(argv);
135
136 // Load problem
137 os << "Loading problem " << prob_path << std::endl;
138 auto problem = load_problem(prob_type, prob_path, opts);
139 os << "Loaded problem " << problem.path.stem().string() << " from "
140 << problem.path << "\nnvar: " << problem.problem.get_num_variables()
141 << "\nncon: " << problem.problem.get_num_constraints()
142 << "\nProvided functions:\n";
143 alpaqa::print_provided_functions(os, problem.problem);
144 os << std::endl;
145
146 // Options
147 auto has_opt = [&opts](std::string_view o) {
148 auto o_it = std::ranges::find(opts.options(), o);
149 if (o_it == opts.options().end())
150 return false;
151 auto index = static_cast<size_t>(o_it - opts.options().begin());
152 ++opts.used()[index];
153 return true;
154 };
155 CheckGradientsOpts cg_opts{
156 .print_full = has_opt("--full-print"),
157 .hessians = !has_opt("--no-hessians"),
158 .scale_perturbations = 1e-2,
159 };
160 set_params(cg_opts.scale_perturbations, "--scale-perturbations", opts);
161
162 // Seed rand
163 auto seed = static_cast<unsigned int>(std::time(nullptr));
164 set_params(seed, "--seed", opts);
165 std::srand(seed);
166
167 // Check options
168 auto used = opts.used();
169 auto unused_opt = std::ranges::find(used, 0);
170 auto unused_idx = static_cast<size_t>(unused_opt - used.begin());
171 if (unused_opt != used.end())
172 throw std::invalid_argument("Unused option: " +
173 std::string(opts.options()[unused_idx]));
174
175 // Check gradients
176 check_gradients(problem, os, cg_opts);
177
178} catch (std::exception &e) {
179 std::cerr << "Error: " << guanaqo::demangled_typename(typeid(e)) << ":\n "
180 << e.what() << std::endl;
181 return -1;
182}
183
185 const alpaqa::TypeErasedProblem<config_t> problem, crvec x, crvec y,
186 rvec grad_L, rvec work_n) {
187 if (y.size() == 0) /* [[unlikely]] */
188 return problem.eval_objective_gradient(x, grad_L);
190 x, y, grad_L, work_n);
191 grad_L += work_n;
192}
193
195 const alpaqa::TypeErasedProblem<config_t> &problem, crvec x, crvec y,
196 crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) -> real_t {
197 if (y.size() == 0) /* [[unlikely]] */
198 return problem.eval_objective_and_gradient(x, grad_ψ);
199
200 auto &ŷ = work_m;
201 // ψ(x) = f(x) + ½ dᵀŷ
202 auto f = problem.eval_objective_and_constraints(x, ŷ);
203 auto dᵀŷ = problem.calc_ŷ_dᵀŷ(ŷ, y, Σ);
204 auto ψ = f + real_t(0.5) * dᵀŷ;
205 // ∇ψ(x) = ∇f(x) + ∇g(x) ŷ
206 try {
207 default_eval_lagrangian_gradient(problem, x, ŷ, grad_ψ, work_n);
208 } catch (alpaqa::not_implemented_error &) {
209 problem.eval_lagrangian_gradient(x, ŷ, grad_ψ, work_n);
210 }
211 return ψ;
212}
213
214auto finite_diff(const std::function<real_t(crvec)> &f, crvec x) {
215 const auto n = x.size();
216 vec grad(n);
217 vec h = vec::Zero(n);
218 const auto ε = 5e-6;
219 const auto δ = 1e-2 * ε;
220 const real_t fx = f(x);
221 for (index_t i = 0; i < n; ++i) {
222 real_t hh = std::abs(x(i)) * ε > δ ? x(i) * ε : δ;
223 h(i) = hh;
224 grad.coeffRef(i) = (f(x + h) - fx) / hh;
225 h(i) = 0;
226 }
227 return std::make_tuple(fx, std::move(grad));
228}
229
230auto finite_diff_hess_sparse(const std::function<void(crvec, rvec)> &grad_L,
231 crvec x) {
232 const auto n = x.size();
233 std::vector<Eigen::Triplet<real_t, index_t>> coo;
234 vec h = vec::Zero(n);
235 const auto ε = 5e-6;
236 const auto δ = 1e-2 * ε;
237 vec grad_x(n), grad_xh(n);
238 grad_L(x, grad_x);
239 for (index_t i = 0; i < n; ++i) {
240 real_t hh = std::abs(x(i)) * ε > δ ? x(i) * ε : δ;
241 h(i) = hh;
242 grad_L(x + h, grad_xh);
243 grad_xh = (grad_xh - grad_x) / hh;
244 for (index_t j = 0; j < n; ++j)
245 if (real_t v = grad_xh(j); v != 0)
246 coo.emplace_back(std::min(j, i), std::max(i, j),
247 v * (i == j ? 1 : 0.5));
248 h(i) = 0;
249 }
250 Eigen::SparseMatrix<real_t, 0, index_t> hess(n, n);
251 hess.setFromTriplets(coo.begin(), coo.end());
252 return hess;
253}
254
255auto finite_diff_hess(const std::function<void(crvec, rvec)> &grad_L, crvec x) {
256 const auto n = x.size();
257 vec h = vec::Zero(n);
258 const auto ε = 5e-6;
259 const auto δ = 1e-2 * ε;
260 vec grad_x(n), grad_xh(n);
261 mat hess(n, n);
262 grad_L(x, grad_x);
263 for (index_t i = 0; i < n; ++i) {
264 real_t hh = std::abs(x(i)) * ε > δ ? x(i) * ε : δ;
265 h(i) = hh;
266 grad_L(x + h, grad_xh);
267 hess.col(i) = (grad_xh - grad_x) / hh;
268 h(i) = 0;
269 }
270 return hess;
271}
272
273void check_gradients(alpaqa::LoadedProblem &lproblem, std::ostream &log,
274 const CheckGradientsOpts &opts) {
275 auto &te_problem = lproblem.problem;
276
277 auto x0 = lproblem.initial_guess_x;
278 auto y0 = lproblem.initial_guess_y;
279 auto sc = opts.scale_perturbations + x0.norm();
280 auto n = te_problem.get_num_variables(),
281 m = te_problem.get_num_constraints();
282
283 std::srand(static_cast<unsigned int>(std::time(nullptr)));
284 vec Σ = 1.5 * vec::Random(m).array() + 2;
285 vec y = y0 + (opts.scale_perturbations + y0.norm()) * vec::Random(m);
286 vec x = x0 + sc * vec::Random(n);
287 vec v = 5e-6 * sc * vec::Random(n);
288 vec gx(m);
289 vec wn(n), wm(m);
290
291 auto print_compare = [&log, &opts](const auto &other, const auto &ref) {
292 auto abs_err =
293 (ref - other).reshaped().template lpNorm<Eigen::Infinity>();
294 auto rel_err =
295 abs_err / ref.reshaped().template lpNorm<Eigen::Infinity>();
296 log << " abs error = " << alpaqa::float_to_str(abs_err) << '\n';
297 log << " rel error = " << alpaqa::float_to_str(rel_err) << '\n';
298 if (opts.print_full) {
299 alpaqa::print_python(log << " (1) = ", other);
300 alpaqa::print_python(log << " (2) = ", ref) << std::endl;
301 }
302 };
303 auto print_compare_scal = [&log, &opts](const auto &fd, const auto &ad) {
304 auto abs_err = std::abs(fd - ad);
305 auto rel_err = abs_err / std::abs(fd);
306 log << " abs error = " << alpaqa::float_to_str(abs_err) << '\n';
307 log << " rel error = " << alpaqa::float_to_str(rel_err) << '\n';
308 if (opts.print_full) {
309 log << " (1) = " << alpaqa::float_to_str(fd) << '\n';
310 log << " (2) = " << alpaqa::float_to_str(ad) << '\n' << std::endl;
311 }
312 };
313
314 auto f = [&](crvec x) { return te_problem.eval_objective(x); };
315 log << "Gradient verification: ∇f(x) (grad_f compared to finite "
316 "differences of f)\n";
317 auto [fx, fd_grad_f] = finite_diff(f, x);
318 vec grad_f(n);
319 te_problem.eval_objective_gradient(x, grad_f);
320 print_compare(grad_f, fd_grad_f);
321
322 if (te_problem.provides_eval_objective_and_gradient()) {
323 log << "Gradient verification: ∇f(x) (f_grad_f compared to grad_f)\n";
324 vec f_grad_f(n);
325 auto f2 = te_problem.eval_objective_and_gradient(x, f_grad_f);
326 print_compare(f_grad_f, grad_f);
327 log << "Function verification: f(x) (f_grad_f compared to f)\n";
328 print_compare_scal(f2, fx);
329 }
330
331 log << "Gradient verification: ∇L(x) (grad_L compared to finite "
332 "differences of f + yᵀg)\n";
333 auto L = [&](crvec x) {
334 te_problem.eval_constraints(x, gx);
335 return te_problem.eval_objective(x) + gx.dot(y);
336 };
337 auto [Lx, fd_grad_L] = finite_diff(L, x);
338 vec grad_L(n);
339 te_problem.eval_lagrangian_gradient(x, y, grad_L, wn);
340 print_compare(grad_L, fd_grad_L);
341
342 log << "Gradient verification: ∇ψ(x) (grad_ψ compared to finite "
343 "differences of ψ)\n";
344 auto ψ = [&](crvec x) {
345 return te_problem.eval_augmented_lagrangian(x, y, Σ, wm);
346 };
347 auto [ψx, fd_grad_ψ] = finite_diff(ψ, x);
348 vec grad_ψ(n);
349 te_problem.eval_augmented_lagrangian_gradient(x, y, Σ, grad_ψ, wn, wm);
350 print_compare(grad_ψ, fd_grad_ψ);
351
352 log << "Gradient verification: ∇ψ(x) (grad_ψ compared to reference "
353 "implementation based on g, ∇f, ∇g)\n";
354 vec grad_ψ_default(n);
356 te_problem, x, y, Σ, grad_ψ_default, wn, wm);
357 print_compare(grad_ψ, grad_ψ_default);
358 log << "Function verification: ψ(x) (ψ compared to reference "
359 "implementation based on f, g)\n";
360 print_compare_scal(ψx, ψ_default);
361
362 if (te_problem.provides_eval_augmented_lagrangian_and_gradient()) {
363 log << "Gradient verification: ∇ψ(x) (grad_ψ compared to ψ_grad_ψ)\n";
364 vec ψ_grad_ψ(n);
365 real_t ψ2 = te_problem.eval_augmented_lagrangian_and_gradient(
366 x, y, Σ, ψ_grad_ψ, wn, wm);
367 print_compare(grad_ψ, ψ_grad_ψ);
368 log << "Function verification: ψ(x) (ψ compared to ψ_grad_ψ)\n";
369 print_compare_scal(ψx, ψ2);
370 }
371
372 if (te_problem.provides_eval_lagrangian_hessian_product()) {
373 log << "Hessian product verification: ∇²L(x) (hess_L_prod compared to "
374 "finite differences of grad_L)\n";
375 vec grad_Lv(n);
376 vec xv = x + v;
377 te_problem.eval_lagrangian_gradient(xv, y, grad_Lv, wn);
378 vec fd_hess_Lv = grad_Lv - grad_L;
379 vec hess_Lv(n);
380 te_problem.eval_lagrangian_hessian_product(x, y, 1, v, hess_Lv);
381 print_compare(hess_Lv, fd_hess_Lv);
382 }
383
384 if (te_problem.provides_eval_augmented_lagrangian_hessian_product()) {
385 log << "Hessian product verification: ∇²ψ(x) (hess_ψ_prod compared to "
386 "finite differences of grad_ψ)\n";
387 vec grad_ψv(n);
388 vec xv = x + v;
389 te_problem.eval_augmented_lagrangian_gradient(xv, y, Σ, grad_ψv, wn,
390 wm);
391 vec fd_hess_ψv = grad_ψv - grad_ψ;
392 vec hess_ψv(n);
393 te_problem.eval_augmented_lagrangian_hessian_product(x, y, Σ, 1, v,
394 hess_ψv);
395 print_compare(hess_ψv, fd_hess_ψv);
396 }
397
398 if (opts.hessians && te_problem.provides_eval_lagrangian_hessian()) {
399 log << "Hessian verification: ∇²L(x) (hess_L compared to finite "
400 "differences of grad_L)\n";
401 auto sp = te_problem.get_lagrangian_hessian_sparsity();
403 auto eval_h = [&](rvec v) {
404 te_problem.eval_lagrangian_hessian(x, y, 1., v);
405 };
406 auto hess_L = cvt.eval(eval_h);
407 mat fd_hess_L = finite_diff_hess(
408 [&](crvec x, rvec g) {
409 te_problem.eval_lagrangian_gradient(x, y, g, wn);
410 },
411 x);
412 print_compare(hess_L, fd_hess_L);
413 }
414
415 if (opts.hessians &&
416 te_problem.provides_eval_augmented_lagrangian_hessian()) {
417 log << "Hessian verification: ∇²ψ(x) (hess_ψ compared to finite "
418 "differences of grad_ψ)\\n";
419 auto sp = te_problem.get_augmented_lagrangian_hessian_sparsity();
421 auto eval_h = [&](rvec v) {
422 te_problem.eval_augmented_lagrangian_hessian(x, y, Σ, 1., v);
423 };
424 auto hess_ψ = cvt.eval(eval_h);
425 mat fd_hess_ψ = finite_diff_hess(
426 [&](crvec x, rvec g) {
427 te_problem.eval_augmented_lagrangian_gradient(x, y, Σ, g, wn,
428 wm);
429 },
430 x);
431 print_compare(hess_ψ, fd_hess_ψ);
432 }
433}
std::span< unsigned > used()
Definition options.hpp:63
std::span< const std::string_view > options() const
Definition options.hpp:60
The main polymorphic minimization problem interface.
void eval_objective_gradient_and_constraints_gradient_product(crvec x, crvec y, rvec grad_f, rvec grad_gxy) const
[Optional] Evaluate both and .
void eval_objective_gradient(crvec x, rvec grad_fx) const
[Required] Function that evaluates the gradient of the cost,
#define USING_ALPAQA_CONFIG(Conf)
Definition config.hpp:77
auto finite_diff(const std::function< real_t(crvec)> &f, crvec x)
auto finite_diff_hess(const std::function< void(crvec, rvec)> &grad_L, crvec x)
void default_eval_lagrangian_gradient(const alpaqa::TypeErasedProblem< config_t > problem, crvec x, crvec y, rvec grad_L, rvec work_n)
const auto * docs
auto get_problem_path(const char *const *argv)
void print_usage(const char *a0)
void check_gradients(alpaqa::LoadedProblem &, std::ostream &, const CheckGradientsOpts &)
auto default_eval_augmented_lagrangian_and_gradient(const alpaqa::TypeErasedProblem< config_t > &problem, crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) -> real_t
auto split_once(std::string_view s, char tok='.')
Split the string s on the first occurrence of tok.
int main(int argc, const char *argv[])
auto finite_diff_hess_sparse(const std::function< void(crvec, rvec)> &grad_L, crvec x)
void print_provided_functions(std::ostream &os, const TypeErasedProblem< Conf > &problem)
vec initial_guess_y
Multipliers g.
EigenConfigd DefaultConfig
Definition config.hpp:31
alpaqa::TypeErasedProblem< config_t > problem
std::ostream & print_python(std::ostream &os, const Eigen::DenseBase< Derived > &M)
Definition print.hpp:13