alpaqa 1.0.0a8
Nonconvex constrained optimization
Loading...
Searching...
No Matches
cutest-loader.cpp
Go to the documentation of this file.
3
4#include <cutest.h>
5#include <dlfcn.h>
6
7#include <cassert>
8#include <functional>
9#include <iostream>
10#include <memory>
11#include <stdexcept>
12#include <string>
13#include <string_view>
14#include <vector>
15
16/*
17CUTEST_cfn : function and constraints values
18CUTEST_cofg : function value and possibly gradient
19CUTEST_cofsg : function value and possibly gradient in sparse format
20CUTEST_ccfg : constraint functions values and possibly gradients
21CUTEST_clfg : Lagrangian function value and possibly gradient
22CUTEST_cgr : constraints gradients and gradient of objective/Lagrangian function
23CUTEST_csgr : constraints gradients and gradient of objective/Lagrangian function
24CUTEST_ccfsg : constraint functions values and possibly their gradients in sparse format
25CUTEST_ccifg : single constraint function value and possibly its gradient
26CUTEST_ccifsg : single constraint function value and possibly gradient in sparse format
27CUTEST_cgrdh : constraints gradients, Hessian of Lagrangian function and gradient of objective/Lagrangian function
28CUTEST_cdh : Hessian of the Lagrangian
29CUTEST_cdhc : Hessian of the constraint part of the Lagrangian
30CUTEST_cshp : sparsity pattern of the Hessian of the Lagrangian function
31CUTEST_csh : Hessian of the Lagrangian, in sparse format
32CUTEST_cshc : Hessian of the constraint part of the Lagrangian, in sparse format
33CUTEST_ceh : sparse Lagrangian Hessian matrix in finite element format
34CUTEST_cifn : problem function value
35CUTEST_cigr : gradient of a problem function
36CUTEST_cisgr : gradient of a problem function in sparse format
37CUTEST_cidh : Hessian of a problem function
38CUTEST_cish : Hessian of an individual problem function, in sparse format
39CUTEST_csgrsh : constraints gradients, sparse Lagrangian Hessian and the gradient of either the objective/Lagrangian in sparse format
40CUTEST_csgreh : constraint gradients, the Lagrangian Hessian in finite element format and the gradient of either the objective/Lagrangian in sparse format
41CUTEST_chprod : matrix-vector product of a vector with the Hessian matrix of the Lagrangian
42CUTEST_cshprod : matrix-vector product of a sparse vector with the Hessian matrix of the Lagrangian
43CUTEST_chcprod : matrix-vector product of a vector with the Hessian matrix of the constraint part of the Lagrangian
44CUTEST_cshcprod : matrix-vector product of a spaarse vector with the Hessian matrix of the constraint part of the Lagrangian
45CUTEST_cjprod : matrix-vector product of a vector with the Jacobian of the constraints, or its transpose
46CUTEST_csjprod : matrix-vector product of a sparse vector with the Jacobian of the constraints, or its transpose
47CUTEST_cchprods : matrix-vector products of a vector with each of the Hessian matrices of the constraint functions
48*/
49
50#ifndef CUTEST_csjp
51#define CUTEST_csjp FUNDERSCORE(cutest_csjp)
52void CUTEST_csjp(integer *status, integer *nnzj, const integer *lj,
53 integer *J_var, integer *J_con);
54#endif
55
56#define STR(x) #x
57#define LOAD_DL_FUNC(f) dlfun<decltype(f)>(STR(f))
58
59using namespace std::string_literals;
60
61namespace {
62void throw_error(std::string s, int code) {
63 throw std::runtime_error(s + " (" + std::to_string(code) + ")");
64}
65void throw_if_error(std::string s, int code) {
66 if (code)
67 throw_error(s, code);
68}
69void log_if_error(std::string s, int code) {
70 if (code)
71 std::cerr << s << " (" << code << ")\n";
72}
73
74std::shared_ptr<void> load_lib(const char *so_filename) {
75 assert(so_filename);
76 ::dlerror();
77 void *h = ::dlopen(so_filename, RTLD_LOCAL | RTLD_NOW);
78 if (auto *err = ::dlerror())
79 throw std::runtime_error(err);
80 return std::shared_ptr<void>{h, &::dlclose};
81}
82} // namespace
83
84namespace alpaqa {
85
87 public:
88 USING_ALPAQA_CONFIG(CUTEstProblem::config_t);
90
91 private:
92 using cleanup_t = std::shared_ptr<void>;
93 template <class F>
94 cleanup_t cleanup(F &&func) {
95 return cleanup_t{nullptr,
96 [func{std::forward<F>(func)}](void *) { func(); }};
97 }
98
99 cleanup_t load_outsdif(const char *outsdif_fname) {
100 integer status;
101 auto fptr_open = LOAD_DL_FUNC(FORTRAN_open);
102 auto fptr_close = LOAD_DL_FUNC(FORTRAN_close);
103 fptr_open(&funit, outsdif_fname, &status);
104 throw_if_error("Failed to open "s + outsdif_fname, status);
105 return cleanup([funit{this->funit}, fptr_close] {
106 integer status;
107 fptr_close(&funit, &status);
108 log_if_error("Failed to close OUTSDIF.d file", status);
109 });
110 }
111
113 auto fptr_cterminate = LOAD_DL_FUNC(CUTEST_cterminate);
114 return cleanup([fptr_cterminate] {
115 integer status;
116 fptr_cterminate(&status);
117 log_if_error("Failed to call cutest_cterminate", status);
118 });
119 }
120
121 public:
122 CUTEstLoader(const char *so_fname, const char *outsdif_fname) {
123 // Open the shared library
124 so_handle = load_lib(so_fname);
125
126 // Open the OUTSDIF.d file
127 cleanup_outsdif = load_outsdif(outsdif_fname);
128
129 // Get the dimensions of the problem
130 integer status;
131 auto fptr_cdimen = LOAD_DL_FUNC(CUTEST_cdimen);
132 fptr_cdimen(&status, &funit, &nvar, &ncon);
133 throw_if_error("Failed to call cutest_cdimen", status);
134 }
135
136 struct ConstrFuncs {
137 decltype(CUTEST_cfn) *cfn;
138 decltype(CUTEST_cofg) *cofg;
139 decltype(CUTEST_ccfg) *ccfg;
140 decltype(CUTEST_clfg) *clfg;
141 decltype(CUTEST_cjprod) *cjprod;
142 decltype(CUTEST_ccifg) *ccifg;
143 decltype(CUTEST_cigr) *cigr;
144 decltype(CUTEST_cdimsj) *cdimsj;
145 decltype(CUTEST_csjp) *csjp;
146 decltype(CUTEST_ccfsg) *ccfsg;
147 decltype(CUTEST_cdh) *cdh;
148 decltype(CUTEST_cdimsh) *cdimsh;
149 decltype(CUTEST_cshp) *cshp;
150 decltype(CUTEST_csh) *csh;
151 decltype(CUTEST_chprod) *chprod;
152 };
153
154 void setup_problem(rvec x0, rvec y0, Box &C, Box &D) {
155 assert(x0.size() == static_cast<length_t>(nvar));
156 assert(C.lowerbound.size() == static_cast<length_t>(nvar));
157 assert(C.upperbound.size() == static_cast<length_t>(nvar));
158 assert(y0.size() == static_cast<length_t>(ncon));
159 assert(D.lowerbound.size() == static_cast<length_t>(ncon));
160 assert(D.upperbound.size() == static_cast<length_t>(ncon));
161 equatn.resize(static_cast<length_t>(ncon));
162 linear.resize(static_cast<length_t>(ncon));
163 integer e_order = 0; // no specific order of equality constraints
164 integer l_order = 0; // no specific order of linear constraints
165 integer v_order = 0; // no specific order of linear variables
166 integer status;
167
168 LOAD_DL_FUNC(CUTEST_csetup)
169 (&status, &funit, &iout, &io_buffer, &nvar, &ncon, x0.data(),
170 C.lowerbound.data(), C.upperbound.data(), y0.data(),
171 D.lowerbound.data(), D.upperbound.data(), equatn.data(), linear.data(),
172 &e_order, &l_order, &v_order);
173 throw_if_error("Failed to call cutest_csetup", status);
175 if (ncon == 0)
176 throw std::runtime_error(
177 "Unconstrained CUTEst problems are currently unsupported");
178 work.resize(std::max(nvar, ncon));
179 std::replace(C.lowerbound.begin(), C.lowerbound.end(), -CUTE_INF,
180 -inf<config_t>);
181 std::replace(C.upperbound.begin(), C.upperbound.end(), +CUTE_INF,
182 +inf<config_t>);
183 std::replace(D.lowerbound.begin(), D.lowerbound.end(), -CUTE_INF,
184 -inf<config_t>);
185 std::replace(D.upperbound.begin(), D.upperbound.end(), +CUTE_INF,
186 +inf<config_t>);
187 funcs = {
188 .cfn = LOAD_DL_FUNC(CUTEST_cfn),
189 .cofg = LOAD_DL_FUNC(CUTEST_cofg),
190 .ccfg = LOAD_DL_FUNC(CUTEST_ccfg),
191 .clfg = LOAD_DL_FUNC(CUTEST_clfg),
192 .cjprod = LOAD_DL_FUNC(CUTEST_cjprod),
193 .ccifg = LOAD_DL_FUNC(CUTEST_ccifg),
194 .cigr = LOAD_DL_FUNC(CUTEST_cigr),
195 .cdimsj = LOAD_DL_FUNC(CUTEST_cdimsj),
196 .csjp = LOAD_DL_FUNC(CUTEST_csjp),
197 .ccfsg = LOAD_DL_FUNC(CUTEST_ccfsg),
198 .cdh = LOAD_DL_FUNC(CUTEST_cdh),
199 .cdimsh = LOAD_DL_FUNC(CUTEST_cdimsh),
200 .cshp = LOAD_DL_FUNC(CUTEST_cshp),
201 .csh = LOAD_DL_FUNC(CUTEST_csh),
202 .chprod = LOAD_DL_FUNC(CUTEST_chprod),
203 };
204 }
205
206 std::string get_name() {
207 std::string name(FSTRING_LEN, '0');
208 integer status;
209 LOAD_DL_FUNC(CUTEST_probname)(&status, name.data());
210 throw_if_error("Failed to call CUTEST_probname", status);
211 name.resize(name.find_last_not_of(' '));
212 return name;
213 }
214
215 integer get_report(doublereal *calls, doublereal *time) {
216 integer status;
217 auto fptr_report = ncon > 0 ? LOAD_DL_FUNC(CUTEST_creport)
218 : LOAD_DL_FUNC(CUTEST_ureport);
219 fptr_report(&status, calls, time);
220 return status;
221 }
222
223 template <class T>
224 T *dlfun(const char *name) {
225 (void)dlerror();
226 auto res = reinterpret_cast<T *>(::dlsym(so_handle.get(), name));
227 if (const char *error = dlerror())
228 throw std::runtime_error(error);
229 return res;
230 }
231
232 // Order of cleanup is important!
233 std::shared_ptr<void> so_handle; ///< dlopen handle to shared library
234 cleanup_t cleanup_outsdif; ///< Responsible for closing the OUTSDIF.d file
235 cleanup_t cutest_terminate; ///< Responsible for calling CUTEST_xterminate
236
237 integer funit = 42; ///< Fortran Unit Number for OUTSDIF.d file
238 integer iout = 6; ///< Fortran Unit Number for standard output
239 integer io_buffer = 11; ///< Fortran Unit Number for internal IO
240
241 integer nvar; ///< Number of decision variabls
242 integer ncon; ///< Number of constraints
243 ConstrFuncs funcs; /// Pointers to loaded problem functions
244
245 using logical_vec = Eigen::VectorX<logical>;
246 logical_vec equatn; ///< whether the constraint is an equality
247 logical_vec linear; ///< whether the constraint is linear
248 mutable vec work; ///< work vector
249};
250
251CUTEstProblem::CUTEstProblem(const char *so_fname, const char *outsdif_fname,
252 bool sparse)
253 : BoxConstrProblem<config_t>{0, 0}, sparse{sparse} {
254 impl = std::make_unique<CUTEstLoader>(so_fname, outsdif_fname);
255 name = impl->get_name();
256 resize(static_cast<length_t>(impl->nvar),
257 static_cast<length_t>(impl->ncon));
258 x0.resize(n);
259 y0.resize(m);
260 impl->setup_problem(x0, y0, C, D);
261}
262
263CUTEstProblem::CUTEstProblem(const CUTEstProblem &) = default;
264CUTEstProblem &CUTEstProblem::operator=(const CUTEstProblem &) = default;
265CUTEstProblem::CUTEstProblem(CUTEstProblem &&) noexcept = default;
266CUTEstProblem &CUTEstProblem::operator=(CUTEstProblem &&) noexcept = default;
267CUTEstProblem::~CUTEstProblem() = default;
268
269CUTEstProblem::Report CUTEstProblem::get_report() const {
270 double calls[7];
271 double time[2];
272 Report r;
273 using stat_t = decltype(r.status);
274 r.status = static_cast<stat_t>(impl->get_report(calls, time));
275 r.name = impl->get_name();
276 r.nvar = impl->nvar;
277 r.ncon = impl->ncon;
278 r.calls.objective = static_cast<unsigned>(calls[0]);
279 r.calls.objective_grad = static_cast<unsigned>(calls[1]);
280 r.calls.objective_hess = static_cast<unsigned>(calls[2]);
281 r.calls.hessian_times_vector = static_cast<unsigned>(calls[3]);
282 r.calls.constraints = impl->ncon > 0 ? static_cast<unsigned>(calls[4]) : 0;
283 r.calls.constraints_grad =
284 impl->ncon > 0 ? static_cast<unsigned>(calls[5]) : 0;
285 r.calls.constraints_hess =
286 impl->ncon > 0 ? static_cast<unsigned>(calls[6]) : 0;
287 r.time_setup = time[0];
288 r.time = time[1];
289 return r;
290}
291
293 assert(x.size() == static_cast<length_t>(impl->nvar));
294 integer status;
295 real_t f;
296 logical grad = FALSE_;
297 impl->funcs.cofg(&status, &impl->nvar, x.data(), &f, nullptr, &grad);
298 throw_if_error("eval_f: CUTEST_cofg", status);
299 return f;
300}
301void CUTEstProblem::eval_grad_f(crvec x, rvec grad_fx) const {
302 assert(x.size() == static_cast<length_t>(impl->nvar));
303 assert(grad_fx.size() == static_cast<length_t>(impl->nvar));
304 integer status;
305 real_t f;
306 logical grad = TRUE_;
307 impl->funcs.cofg(&status, &impl->nvar, x.data(), &f, grad_fx.data(), &grad);
308 throw_if_error("eval_grad_f: CUTEST_cofg", status);
309}
311 assert(x.size() == static_cast<length_t>(impl->nvar));
312 assert(gx.size() == static_cast<length_t>(impl->ncon));
313 integer status;
314 logical jtrans = TRUE_, grad = FALSE_;
315 integer zero = 0;
316 impl->funcs.ccfg(&status, &impl->nvar, &impl->ncon, x.data(), gx.data(),
317 &jtrans, &zero, &zero, nullptr, &grad);
318 throw_if_error("eval_g: CUTEST_ccfg", status);
319}
321 assert(x.size() == static_cast<length_t>(impl->nvar));
322 assert(y.size() == static_cast<length_t>(impl->ncon));
323 assert(grad_gxy.size() == static_cast<length_t>(impl->nvar));
324 integer status;
325 auto lvector = static_cast<integer>(y.size()),
326 lresult = static_cast<integer>(grad_gxy.size());
327 logical gotj = FALSE_, jtrans = TRUE_;
328 impl->funcs.cjprod(&status, &impl->nvar, &impl->ncon, &gotj, &jtrans,
329 x.data(), y.data(), &lvector, grad_gxy.data(), &lresult);
330 throw_if_error("eval_grad_g_prod: CUTEST_cjprod", status);
331}
332
333void CUTEstProblem::eval_jac_g(crvec x, [[maybe_unused]] rindexvec inner_idx,
334 [[maybe_unused]] rindexvec outer_ptr,
335 rvec J_values) const {
336 // Compute the nonzero values
337 if (J_values.size() > 0) {
338 assert(x.size() == static_cast<length_t>(impl->nvar));
339 integer status;
340 // Sparse Jacobian
341 if (sparse) {
342 assert(nnz_J >= 0);
343 assert(J_values.size() == static_cast<length_t>(nnz_J));
344 assert(J_work.size() == static_cast<length_t>(nnz_J));
345 assert(inner_idx.size() == static_cast<length_t>(nnz_J));
346 assert(outer_ptr.size() == static_cast<length_t>(impl->nvar + 1));
347 const integer nnz = nnz_J;
348 logical grad = TRUE_;
349 impl->funcs.ccfsg(&status, &impl->nvar, &impl->ncon, x.data(),
350 impl->work.data(), &nnz_J, &nnz, J_work.data(),
351 J_col.data(), J_row.data(), &grad);
352 throw_if_error("eval_jac_g: CUTEST_ccfsg", status);
353 auto t0 = std::chrono::steady_clock::now();
354 J_values = J_work(J_perm);
355 auto t1 = std::chrono::steady_clock::now();
356 std::cout << "Permutation of J took: "
357 << std::chrono::duration<double>{t1 - t0}.count() * 1e6
358 << " µs\n";
359 }
360 // Dense Jacobian
361 else {
362 assert(J_values.size() == static_cast<length_t>(impl->nvar) *
363 static_cast<length_t>(impl->ncon));
364 integer status;
365 logical jtrans = FALSE_, grad = TRUE_;
366 impl->funcs.ccfg(&status, &impl->nvar, &impl->ncon, x.data(),
367 impl->work.data(), &jtrans, &impl->ncon,
368 &impl->nvar, J_values.data(), &grad);
369 throw_if_error("eval_jac_g: CUTEST_ccfg", status);
370 }
371 }
372 // Compute sparsity pattern without values
373 else {
374 assert(nnz_J >= 0);
375 assert(inner_idx.size() == static_cast<length_t>(nnz_J));
376 assert(outer_ptr.size() == static_cast<length_t>(impl->nvar + 1));
377 integer status;
378 const integer nnz = nnz_J;
379 impl->funcs.csjp(&status, &nnz_J, &nnz, J_col.data(), J_row.data());
380 throw_if_error("eval_jac_g: CUTEST_csjp", status);
381 std::iota(J_perm.begin(), J_perm.end(), index_t{0});
382 util::sort_triplets(J_row, J_col, J_perm);
383 util::convert_triplets_to_ccs<config_t>(J_row, J_col, inner_idx,
384 outer_ptr, 1);
385 }
386}
388 if (!sparse)
389 return -1;
390 if (nnz_J < 0) {
391 integer status;
392 impl->funcs.cdimsj(&status, &nnz_J);
393 throw_if_error("get_jac_g_num_nonzeros: CUTEST_cdimsj", status);
394 nnz_J -= impl->nvar;
395 assert(nnz_J >= 0);
396 J_col.resize(nnz_J);
397 J_row.resize(nnz_J);
398 J_perm.resize(nnz_J);
399 J_work.resize(nnz_J);
400 }
401 return nnz_J;
402}
404 assert(x.size() == static_cast<length_t>(impl->nvar));
405 assert(grad_gi.size() == static_cast<length_t>(impl->nvar));
406 integer status;
407 auto iprob = static_cast<integer>(i + 1);
408 impl->funcs.cigr(&status, &impl->nvar, &iprob, x.data(), grad_gi.data());
409 throw_if_error("eval_grad_gi: CUTEST_cigr", status);
410}
412 rvec Hv) const {
413 assert(x.size() == static_cast<length_t>(impl->nvar));
414 assert(y.size() == static_cast<length_t>(impl->ncon));
415 assert(v.size() == static_cast<length_t>(impl->nvar));
416 assert(Hv.size() == static_cast<length_t>(impl->nvar));
417 const auto *mult = y.data();
418 if (scale != 1) {
419 impl->work = y * (real_t(1) / scale);
420 mult = impl->work.data();
421 }
422 integer status;
423 logical goth = FALSE_;
424 impl->funcs.chprod(&status, &impl->nvar, &impl->ncon, &goth, x.data(), mult,
425 const_cast<real_t *>(v.data()), Hv.data());
426 throw_if_error("eval_hess_L_prod: CUTEST_chprod", status);
427 if (scale != 1)
428 Hv *= scale;
429}
431 rindexvec inner_idx, rindexvec outer_ptr,
432 rvec H_values) const {
433 // Compute the nonzero values
434 if (H_values.size() > 0) {
435 assert(x.size() == static_cast<length_t>(impl->nvar));
436 assert(y.size() == static_cast<length_t>(impl->ncon));
437 const auto *mult = y.data();
438 if (scale != 1) {
439 impl->work = y * (real_t(1) / scale);
440 mult = impl->work.data();
441 }
442 integer status;
443 // Sparse Hessian
444 if (sparse) {
445 assert(nnz_H >= 0);
446 assert(H_values.size() == static_cast<length_t>(nnz_H));
447 assert(H_work.size() == static_cast<length_t>(nnz_H));
448 assert(inner_idx.size() == static_cast<length_t>(nnz_H));
449 assert(outer_ptr.size() == static_cast<length_t>(impl->nvar + 1));
450 const integer nnz = nnz_H;
451 impl->funcs.csh(&status, &impl->nvar, &impl->ncon, x.data(),
452 y.data(), &nnz_H, &nnz, H_work.data(), H_col.data(),
453 H_row.data());
454 throw_if_error("eval_hess_L: CUTEST_csh", status);
455 auto t0 = std::chrono::steady_clock::now();
456 H_values = H_work(H_perm);
457 auto t1 = std::chrono::steady_clock::now();
458 std::cout << "Permutation of H took: "
459 << std::chrono::duration<double>{t1 - t0}.count() * 1e6
460 << " µs\n";
461 }
462 // Dense Hessian
463 else {
464 assert(H_values.size() == static_cast<length_t>(impl->nvar) *
465 static_cast<length_t>(impl->nvar));
466 impl->funcs.cdh(&status, &impl->nvar, &impl->ncon, x.data(), mult,
467 &impl->nvar, H_values.data());
468 throw_if_error("eval_hess_L: CUTEST_cdh", status);
469 }
470 if (scale != 1)
471 H_values *= scale;
472 }
473 // Compute sparsity pattern without values
474 else {
475 assert(nnz_H >= 0);
476 assert(inner_idx.size() == static_cast<length_t>(nnz_H));
477 assert(outer_ptr.size() == static_cast<length_t>(impl->nvar + 1));
478 integer status;
479 const integer nnz = nnz_H;
480 impl->funcs.cshp(&status, &impl->nvar, &nnz_H, &nnz, H_col.data(),
481 H_row.data());
482 throw_if_error("eval_hess_L: CUTEST_cshp", status);
483 std::iota(H_perm.begin(), H_perm.end(), index_t{0});
484 util::sort_triplets(H_row, H_col, H_perm);
485 util::convert_triplets_to_ccs<config_t>(H_row, H_col, inner_idx,
486 outer_ptr, 1);
487 }
488}
490 if (!sparse)
491 return -1;
492 if (nnz_H < 0) {
493 integer status;
494 impl->funcs.cdimsh(&status, &nnz_H);
495 throw_if_error("get_hess_L_num_nonzeros: CUTEST_cdimsh", status);
496 assert(nnz_H >= 0);
497 H_col.resize(nnz_H);
498 H_row.resize(nnz_H);
499 H_perm.resize(nnz_H);
500 H_work.resize(nnz_H);
501 }
502 return nnz_H;
503}
505 assert(x.size() == static_cast<length_t>(impl->nvar));
506 assert(grad_fx.size() == static_cast<length_t>(impl->nvar));
507 integer status;
508 real_t f;
509 logical grad = TRUE_;
510 impl->funcs.cofg(&status, &impl->nvar, x.data(), &f, grad_fx.data(), &grad);
511 throw_if_error("eval_f_grad_f: CUTEST_cofg", status);
512 return f;
513}
515 assert(x.size() == static_cast<length_t>(impl->nvar));
516 assert(g.size() == static_cast<length_t>(impl->ncon));
517 integer status;
518 real_t f;
519 impl->funcs.cfn(&status, &impl->nvar, &impl->ncon, x.data(), &f, g.data());
520 throw_if_error("eval_f_g: CUTEST_cfn", status);
521 return f;
522}
523void CUTEstProblem::eval_grad_L(crvec x, crvec y, rvec grad_L, rvec) const {
524 assert(x.size() == static_cast<length_t>(impl->nvar));
525 assert(y.size() == static_cast<length_t>(impl->ncon));
526 assert(grad_L.size() == static_cast<length_t>(impl->nvar));
527 integer status;
528 real_t L;
529 logical grad = TRUE_;
530 impl->funcs.clfg(&status, &impl->nvar, &impl->ncon, x.data(), y.data(), &L,
531 grad_L.data(), &grad);
532 throw_if_error("eval_f_g: CUTEST_clfg", status);
533}
534
536 using Status = CUTEstProblem::Report::Status;
537 switch (s) {
538 case Status::Success: return "Success";
539 case Status::AllocationError: return "AllocationError";
540 case Status::ArrayBoundError: return "ArrayBoundError";
541 case Status::EvaluationError: return "EvaluationError";
542 default:;
543 }
544 throw std::out_of_range(
545 "invalid value for pa::CUTEstProblem::Report::Status");
546}
547
548std::ostream &operator<<(std::ostream &os, CUTEstProblem::Report::Status s) {
549 return os << enum_name(s);
550}
551
552std::ostream &operator<<(std::ostream &os, const CUTEstProblem::Report &r) {
553 os << "CUTEst problem: " << r.name << "\r\n\n" //
554 << "Number of variables: " << r.nvar << "\r\n" //
555 << "Number of constraints: " << r.ncon << "\r\n\n" //
556 << "Status: " << r.status << " (" << +r.status << ")\r\n\n" //
557 << "Objective function evaluations: " << r.calls.objective
558 << "\r\n"
559 << "Objective function gradient evaluations: "
560 << r.calls.objective_grad << "\r\n"
561 << "Objective function Hessian evaluations: "
562 << r.calls.objective_hess << "\r\n"
563 << "Hessian times vector products: "
564 << r.calls.objective_hess << "\r\n\n";
565 if (r.ncon > 0) {
566 os << "Constraint function evaluations: "
567 << r.calls.constraints << "\r\n"
568 << "Constraint function gradients evaluations: "
569 << r.calls.constraints_grad << "\r\n"
570 << "Constraint function Hessian evaluations: "
571 << r.calls.constraints_hess << "\r\n\n";
572 }
573 return os << "Setup time: " << r.time_setup << "s\r\n"
574 << "Time since setup: " << r.time << "s";
575}
576
577} // namespace alpaqa
Implements common problem functions for minimization problems with box constraints.
cleanup_t cutest_terminate
Responsible for calling CUTEST_xterminate.
decltype(CUTEST_cigr) * cigr
decltype(CUTEST_ccfg) * ccfg
integer ncon
Number of constraints.
decltype(CUTEST_cjprod) * cjprod
std::shared_ptr< void > cleanup_t
cleanup_t cleanup_outsdif
Responsible for closing the OUTSDIF.d file.
decltype(CUTEST_cofg) * cofg
decltype(CUTEST_chprod) * chprod
CUTEstLoader(const char *so_fname, const char *outsdif_fname)
decltype(CUTEST_ccifg) * ccifg
integer iout
Fortran Unit Number for standard output.
integer get_report(doublereal *calls, doublereal *time)
std::string get_name()
integer funit
Fortran Unit Number for OUTSDIF.d file.
decltype(CUTEST_clfg) * clfg
cleanup_t load_outsdif(const char *outsdif_fname)
decltype(CUTEST_cdimsj) * cdimsj
cleanup_t cleanup(F &&func)
integer nvar
Number of decision variabls.
T * dlfun(const char *name)
std::shared_ptr< void > so_handle
dlopen handle to shared library
decltype(CUTEST_cdimsh) * cdimsh
Eigen::VectorX< logical > logical_vec
Pointers to loaded problem functions.
decltype(FUNDERSCORE(cutest_csjp)) * csjp
decltype(CUTEST_ccfsg) * ccfsg
void setup_problem(rvec x0, rvec y0, Box &C, Box &D)
logical_vec equatn
whether the constraint is an equality
decltype(CUTEST_cshp) * cshp
logical_vec linear
whether the constraint is linear
integer io_buffer
Fortran Unit Number for internal IO.
Wrapper for CUTEst problems loaded from an external shared library.
void eval_grad_gi(crvec x, index_t i, rvec grad_gi) const
real_t eval_f_g(crvec x, rvec g) const
CUTEstProblem(const char *so_fname, const char *outsdif_fname, bool sparse=false)
Load a CUTEst problem from the given shared library and OUTSDIF.d file.
void eval_grad_L(crvec x, crvec y, rvec grad_L, rvec work_n) const
Eigen::VectorX< int > J_row
Eigen::VectorX< int > J_col
Eigen::VectorX< int > H_col
void eval_hess_L(crvec x, crvec y, real_t scale, rindexvec inner_idx, rindexvec outer_ptr, rvec H_values) const
util::copyable_unique_ptr< class CUTEstLoader > impl
real_t eval_f_grad_f(crvec x, rvec grad_fx) const
void eval_grad_g_prod(crvec x, crvec y, rvec grad_gxy) const
void eval_hess_L_prod(crvec x, crvec y, real_t scale, crvec v, rvec Hv) const
Eigen::VectorX< int > H_row
void eval_grad_f(crvec x, rvec grad_fx) const
real_t eval_f(crvec x) const
void eval_g(crvec x, rvec gx) const
length_t get_jac_g_num_nonzeros() const
CUTEstProblem & operator=(const CUTEstProblem &)
void eval_jac_g(crvec x, rindexvec inner_idx, rindexvec outer_ptr, rvec J_values) const
length_t get_hess_L_num_nonzeros() const
#define USING_ALPAQA_CONFIG(Conf)
Definition: config.hpp:42
#define CUTEST_csjp
#define LOAD_DL_FUNC(f)
std::ostream & operator<<(std::ostream &os, PANOCStopCrit s)
typename Conf::real_t real_t
Definition: config.hpp:51
typename Conf::rindexvec rindexvec
Definition: config.hpp:65
typename Conf::index_t index_t
Definition: config.hpp:63
typename Conf::length_t length_t
Definition: config.hpp:62
typename Conf::rvec rvec
Definition: config.hpp:55
typename Conf::crvec crvec
Definition: config.hpp:56
typename Conf::vec vec
Definition: config.hpp:52
constexpr const char * enum_name(PANOCStopCrit s)
vec upperbound
Definition: box.hpp:26
vec lowerbound
Definition: box.hpp:25
The report generated by CUTEst.
int ncon
Number of constraints.
double time_setup
CPU time (in seconds) for CUTEST_csetup.
enum alpaqa::CUTEstProblem::Report::Status status
Exit status.
int nvar
Number of independent variables.
Status
Status returned by CUTEst.
struct alpaqa::CUTEstProblem::Report::@0 calls
Function call counters.
std::string name
Name of the problem.
double time
CPU time (in seconds) since the end of CUTEST_csetup.