QPALM 1.1.5a1
Proximal Augmented Lagrangian method for Quadratic Programs
Loading...
Searching...
No Matches
qpalm.c
Go to the documentation of this file.
1/**
2 * @file qpalm.c
3 * @author Ben Hermans
4 * @brief QPALM main solver API.
5 * @details This file contains the main functions that can be called by the user.
6 * The user can load the default settings, setup the workspace with data and settings,
7 * run the solver, and cleanup the workspace afterwards.
8 */
9# ifdef __cplusplus
10extern "C" {
11# endif // ifdef __cplusplus
12
13#include "qpalm.h"
14#include "global_opts.h"
15#include "constants.h"
16#include "validate.h"
17#include "lin_alg.h"
18#include "util.h"
19#include "scaling.h"
20#include "linesearch.h"
21#include "termination.h"
22#include "solver_interface.h"
23#include "newton.h"
24#include "nonconvex.h"
25#include "iteration.h"
26
27#include "ladel.h"
28
29/**********************
30* Main API Functions *
31**********************/
32
34{
35 settings->max_iter = MAX_ITER; /* maximum iterations */
36 settings->inner_max_iter = INNER_MAX_ITER; /* maximum iterations per subproblem */
37 settings->eps_abs = (c_float)EPS_ABS; /* absolute convergence tolerance */
38 settings->eps_rel = (c_float)EPS_REL; /* relative convergence tolerance */
39 settings->eps_abs_in = (c_float)EPS_ABS_IN; /* intermediate absolute convergence tolerance */
40 settings->eps_rel_in = (c_float)EPS_REL_IN; /* intermediate relative convergence tolerance */
41 settings->rho = (c_float)RHO; /* tolerance scaling factor */
42 settings->eps_prim_inf = (c_float)EPS_PRIM_INF; /* primal infeasibility tolerance */
43 settings->eps_dual_inf = (c_float)EPS_DUAL_INF; /* dual infeasibility tolerance */
44 settings->theta = (c_float)THETA; /* penalty update criterion parameter */
45 settings->delta = (c_float)DELTA; /* penalty update factor */
46 settings->sigma_max = (c_float)SIGMA_MAX; /* penalty parameter cap */
47 settings->sigma_init = (c_float)SIGMA_INIT; /* initial penalty parameter (guideline) */
48 settings->proximal = PROXIMAL; /* boolean, proximal method of multipliers*/
49 settings->gamma_init = (c_float)GAMMA_INIT; /* proximal penalty parameter */
50 settings->gamma_upd = (c_float)GAMMA_UPD; /* proximal penalty update factor*/
51 settings->gamma_max = (c_float)GAMMA_MAX; /* proximal penalty parameter cap*/
52 settings->scaling = SCALING; /* boolean, scaling */
53 settings->nonconvex = NONCONVEX; /* boolean, nonconvex */
54 settings->warm_start = WARM_START; /* boolean, warm start solver */
55 settings->verbose = VERBOSE; /* boolean, write out progress */
56 settings->print_iter = PRINT_ITER; /* frequency of printing */
57 settings->reset_newton_iter = RESET_NEWTON_ITER; /* frequency of performing a full Cholesky factorization */
58 settings->enable_dual_termination = ENABLE_DUAL_TERMINATION; /* allow for dual termination (useful in branch and bound) */
59 settings->dual_objective_limit = (c_float)DUAL_OBJECTIVE_LIMIT; /* termination value for the dual objective (useful in branch and bound) */
60 settings->time_limit = (c_float)TIME_LIMIT; /* time limit */
61 settings->ordering = ORDERING; /* ordering */
62 settings->factorization_method = FACTORIZATION_METHOD; /* factorization method (kkt or schur) */
63 settings->max_rank_update = MAX_RANK_UPDATE; /* maximum rank of the update */
64 settings->max_rank_update_fraction = (c_float)MAX_RANK_UPDATE_FRACTION; /* maximum rank (relative to n+m) of the update */
65}
66
67
68QPALMWorkspace* qpalm_setup(const QPALMData *data, const QPALMSettings *settings)
69{
70 QPALMWorkspace *work; // Workspace
71
72 // Validate data
73 if (!validate_data(data))
74 {
75 # ifdef QPALM_PRINTING
76 qpalm_eprint("Data validation returned failure");
77 # endif /* ifdef QPALM_PRINTING */
78 return QPALM_NULL;
79 }
80
81 // Validate settings
82 if (!validate_settings(settings))
83 {
84 # ifdef QPALM_PRINTING
85 qpalm_eprint("Settings validation returned failure");
86 # endif /* ifdef QPALM_PRINTING */
87 return QPALM_NULL;
88 }
89
90 // Allocate empty workspace
91 work = qpalm_calloc(1, sizeof(QPALMWorkspace));
92
93 if (!work)
94 {
95 # ifdef QPALM_PRINTING
96 qpalm_eprint("allocating work failure");
97 # endif /* ifdef QPALM_PRINTING */
98 return QPALM_NULL;
99 }
100
101 // Start and allocate directly timer
102 # ifdef QPALM_TIMING
103 work->timer = qpalm_malloc(sizeof(QPALMTimer));
104 qpalm_tic(work->timer);
105 # endif /* ifdef QPALM_TIMING */
106
107 // Copy settings
108 work->settings = copy_settings(settings);
109 work->sqrt_delta = c_sqrt(work->settings->delta);
110 work->gamma = work->settings->gamma_init;
111
112 size_t n = data->n;
113 size_t m = data->m;
114
115 //Initialize the solver for the linear system
116 work->solver = qpalm_calloc(1, sizeof(QPALMSolver));
117 solver_common common, *c;
118 c = &common;
119
120 // Copy problem data into workspace
121 work->data = qpalm_calloc(1, sizeof(QPALMData));
122 work->data->n = data->n;
123 work->data->m = data->m;
124 work->data->bmin = vec_copy(data->bmin, m);
125 work->data->bmax = vec_copy(data->bmax, m);
126 work->data->q = vec_copy(data->q, n);
127 work->data->c = data->c;
128
129 work->data->A = ladel_sparse_allocate_and_copy(data->A);
130 work->data->Q = ladel_sparse_allocate_and_copy(data->Q);
131 ladel_to_upper_diag(work->data->Q);
132
133 // Allocate internal solver variables
134 work->x = qpalm_calloc(n, sizeof(c_float));
135 work->y = qpalm_calloc(m, sizeof(c_float));
136 work->Ax = qpalm_calloc(m, sizeof(c_float));
137 work->Qx = qpalm_calloc(n, sizeof(c_float));
138 work->x_prev = qpalm_calloc(n, sizeof(c_float));
139 work->Aty = qpalm_calloc(n, sizeof(c_float));
140
141 work->x0 = qpalm_calloc(n, sizeof(c_float));
142
143 work->initialized = FALSE;
144
145 // Workspace variables
146 work->temp_m = qpalm_calloc(m, sizeof(c_float));
147 work->temp_n = qpalm_calloc(n, sizeof(c_float));
148 work->sigma = qpalm_calloc(m, sizeof(c_float));
149 work->sigma_inv = qpalm_calloc(m, sizeof(c_float));
150 work->nb_sigma_changed = 0;
151
152 work->z = qpalm_calloc(m, sizeof(c_float));
153 work->Axys = qpalm_calloc(m, sizeof(c_float));
154 work->pri_res = qpalm_calloc(m, sizeof(c_float));
155 work->pri_res_in = qpalm_calloc(m, sizeof(c_float));
156 work->df = qpalm_calloc(n, sizeof(c_float));
157
158 work->xx0 = qpalm_calloc(n, sizeof(c_float));
159 work->dphi = qpalm_calloc(n, sizeof(c_float));
160 work->dphi_prev = qpalm_calloc(n, sizeof(c_float));
161
162 // Linesearch variables
163 work->sqrt_sigma = qpalm_calloc(m, sizeof(c_float));
164 work->delta = qpalm_calloc(m*2, sizeof(c_float));
165 work->alpha = qpalm_calloc(m*2, sizeof(c_float));
166 work->delta2 = qpalm_calloc(m*2, sizeof(c_float));
167 work->delta_alpha = qpalm_calloc(m*2, sizeof(c_float));
168 work->temp_2m = qpalm_calloc(m*2, sizeof(c_float));
169 work->s = qpalm_calloc(m*2, sizeof(array_element));
170 work->index_L = qpalm_calloc(m*2, sizeof(c_int));
171 work->index_P = qpalm_calloc(m*2, sizeof(c_int));
172 work->index_J = qpalm_calloc(m*2, sizeof(c_int));
173
174 // Primal infeasibility variables
175 work->delta_y = qpalm_calloc(m, sizeof(c_float));
176 work->Atdelta_y = qpalm_calloc(n, sizeof(c_float));
177
178 // Dual infeasibility variables
179 work->delta_x = qpalm_calloc(n, sizeof(c_float));
180 work->Qdelta_x = qpalm_calloc(n, sizeof(c_float));
181 work->Adelta_x = qpalm_calloc(m, sizeof(c_float));
182
184 // c = &common;
185
186 // Allocate scaling structure
187 work->scaling = qpalm_malloc(sizeof(QPALMScaling));
188 work->scaling->D = qpalm_calloc(n, sizeof(c_float));
189 work->scaling->Dinv = qpalm_calloc(n, sizeof(c_float));
190 work->scaling->E = qpalm_calloc(m, sizeof(c_float));
191 work->scaling->Einv = qpalm_calloc(m, sizeof(c_float));
192
193 work->solver->E_temp = qpalm_calloc(m, sizeof(c_float));
194 work->E_temp = work->solver->E_temp;
195 work->solver->D_temp = qpalm_calloc(n, sizeof(c_float));
196 work->D_temp = work->solver->D_temp;
197
198 // Solver variables
199 work->solver->active_constraints = qpalm_calloc(m, sizeof(c_int));
200 work->solver->active_constraints_old = qpalm_calloc(m, sizeof(c_int));
202 work->solver->reset_newton = TRUE;
203 work->solver->enter = qpalm_calloc(m, sizeof(c_int));
204 work->solver->leave = qpalm_calloc(m, sizeof(c_int));
205
207 {
208 work->solver->rhs_kkt = qpalm_malloc((n+m)*sizeof(c_float));
209 work->solver->sol_kkt = qpalm_malloc((n+m)*sizeof(c_float));
210 c_int kkt_nzmax = work->data->Q->nzmax + work->data->A->nzmax + m;
211 work->solver->kkt_full = ladel_sparse_alloc(n+m, n+m, kkt_nzmax, UPPER, TRUE, FALSE);
212 work->solver->kkt = ladel_sparse_alloc(n+m, n+m, kkt_nzmax, UPPER, TRUE, TRUE);
213 work->solver->first_row_A = qpalm_malloc(m*sizeof(c_int));
214 work->solver->first_elem_A = qpalm_malloc(m*sizeof(c_float));
215 work->solver->sym = ladel_symbolics_alloc(m+n);
216 }
218 {
219 work->solver->sym = ladel_symbolics_alloc(n);
220 }
221
222 work->solver->neg_dphi = qpalm_calloc(n, sizeof(c_float));
223 work->neg_dphi = work->solver->neg_dphi;
224 work->solver->d = qpalm_calloc(n, sizeof(c_float));
225 work->d = work->solver->d;
226 work->solver->Qd = qpalm_calloc(n, sizeof(c_float));
227 work->Qd = work->solver->Qd;
228 work->solver->Ad = qpalm_calloc(m, sizeof(c_float));
229 work->Ad = work->solver->Ad;
230 work->solver->yh = qpalm_calloc(m, sizeof(c_float));
231 work->yh = work->solver->yh;
232 work->solver->Atyh = qpalm_calloc(n, sizeof(c_float));
233 work->Atyh = work->solver->Atyh;
234 work->solver->At_scale = qpalm_calloc(m, sizeof(c_float));
235
237
239 work->solver->sym_Q = ladel_symbolics_alloc(n);
240
241 // Allocate solution
242 work->solution = qpalm_calloc(1, sizeof(QPALMSolution));
243 work->solution->x = qpalm_calloc(1, n * sizeof(c_float));
244 work->solution->y = qpalm_calloc(1, m * sizeof(c_float));
245
246 // Allocate and initialize information
247 work->info = qpalm_calloc(1, sizeof(QPALMInfo));
249 # ifdef QPALM_TIMING
250 work->info->solve_time = 0.0; // Solve time to zero
251 work->info->run_time = 0.0; // Total run time to zero
252 work->info->setup_time = qpalm_toc(work->timer); // Update timer information
253 # endif /* ifdef QPALM_TIMING */
254
255 return work;
256}
257
258
259void qpalm_warm_start(QPALMWorkspace *work, const c_float *x_warm_start, const c_float *y_warm_start)
260{
261 // If we have previously solved the problem, then reset the setup time
262 if (work->info->status_val != QPALM_UNSOLVED)
263 {
264 #ifdef QPALM_TIMING
265 work->info->setup_time = 0;
266 #endif /* ifdef QPALM_TIMING */
268 }
269 #ifdef QPALM_TIMING
270 qpalm_tic(work->timer);
271 #endif
272
273 size_t n = work->data->n;
274 size_t m = work->data->m;
275
276 if (x_warm_start != NULL)
277 {
278 prea_vec_copy(x_warm_start, work->x, n);
279 }
280 else
281 {
282 qpalm_free(work->x);
283 work->x = NULL;
284 }
285
286 if (y_warm_start != NULL)
287 {
288 prea_vec_copy(y_warm_start, work->y, m);
289 }
290 else
291 {
292 qpalm_free(work->y);
293 work->y = NULL;
294 }
295
296 work->initialized = TRUE;
297
298 #ifdef QPALM_TIMING
299 work->info->setup_time += qpalm_toc(work->timer); // Start timer
300 #endif /* ifdef QPALM_TIMING */
301
302}
303
304static void qpalm_initialize(QPALMWorkspace *work, solver_common **common1, solver_common **common2)
305{
306 // If we have previously solved the problem, then reset the setup time
307 if (work->info->status_val != QPALM_UNSOLVED)
308 {
309 #ifdef QPALM_TIMING
310 work->info->setup_time = 0;
311 #endif /* ifdef QPALM_TIMING */
313 }
314 #ifdef QPALM_TIMING
315 qpalm_tic(work->timer);
316 #endif
317
318 // Print header
319 #ifdef QPALM_PRINTING
320 if (work->settings->verbose)
321 {
322 print_header();
323 }
324 #endif
325
326 size_t n = work->data->n;
327 size_t m = work->data->m;
328 *common1 = ladel_workspace_allocate(n+m);
329 if (work->settings->enable_dual_termination) *common2 = ladel_workspace_allocate(n);
330 else *common2 = *common1;
331 solver_common *c = *common1, *c2 = *common2;
332
333 if (!work->initialized)
334 {
335 qpalm_warm_start(work, NULL, NULL);
336 }
337
338 work->eps_abs_in = work->settings->eps_abs_in;
339 work->eps_rel_in = work->settings->eps_rel_in;
341 work->solver->reset_newton = TRUE;
342 work->gamma = work->settings->gamma_init;
343 work->gamma_maxed = FALSE;
345
346 if (work->x == NULL)
347 {
348 work->x = qpalm_calloc(n, sizeof(c_float));
349 vec_set_scalar(work->x, 0., n);
350 vec_set_scalar(work->x_prev, 0., n);
351 vec_set_scalar(work->x0, 0., n);
352 vec_set_scalar(work->Qx, 0., n);
353 vec_set_scalar(work->Ax, 0., m);
354 work->info->objective = work->data->c;
355 }
356 else
357 {
358 mat_vec(work->data->Q, work->x, work->Qx, c);
359 mat_vec(work->data->A, work->x, work->Ax, c);
360 }
361
362 if (work->y == NULL)
363 {
364 work->y = qpalm_calloc(m, sizeof(c_float));
365 vec_set_scalar(work->y, 0., m);
366 }
367
368 for (size_t i = 0; i < work->data->m; i++)
369 {
370 if (work->data->bmax[i] > QPALM_INFTY) work->data->bmax[i] = QPALM_INFTY;
371 if (work->data->bmin[i] < -QPALM_INFTY) work->data->bmin[i] = -QPALM_INFTY;
372 }
373
374 if (work->settings->scaling)
375 {
376 scale_data(work);
377 }
378
379 //Actions to perform after scaling
380 prea_vec_copy(work->x, work->x0, n);
381 prea_vec_copy(work->x, work->x_prev, n);
382
384 {
385 if (work->solver->At) ladel_sparse_free(work->solver->At);
386 work->solver->At = ladel_transpose(work->data->A, TRUE, c);
387 }
388
389 if (work->settings->nonconvex)
390 {
391 set_settings_nonconvex(work, c);
392 }
393 if (work->settings->proximal)
394 {
395 vec_add_scaled(work->Qx, work->x, work->Qx, 1/work->gamma, n);
396 }
397
398 work->info->objective = compute_objective(work);
399
400 initialize_sigma(work, c);
401
402 //Provide LD factor of Q in case dual_termination is enabled
403 //NB assume Q is positive definite
405 {
406 if (work->solver->LD_Q) ladel_factor_free(work->solver->LD_Q);
407 ladel_factorize(work->data->Q, work->solver->sym_Q, work->settings->ordering, &work->solver->LD_Q, c2);
408 work->info->dual_objective = compute_dual_objective(work, c2);
409 } else
410 {
412 }
413
414 #ifdef QPALM_TIMING
415 work->info->setup_time += qpalm_toc(work->timer); // Start timer
416 #endif /* ifdef QPALM_TIMING */
417}
418
419static void qpalm_termination(QPALMWorkspace *work, solver_common* c, solver_common *c2, c_int iter, c_int iter_out)
420{
421 if (work->info->status_val == QPALM_SOLVED ||
425 {
426 store_solution(work);
427 }
428 else if (work->info->status_val == QPALM_PRIMAL_INFEASIBLE)
429 {
430 if (work->settings->scaling)
431 {
432 vec_self_mult_scalar(work->delta_y, work->scaling->cinv, work->data->m);
433 vec_ew_prod(work->scaling->E, work->delta_y, work->delta_y, work->data->m);
434 }
435 }
436 else if (work->info->status_val == QPALM_DUAL_INFEASIBLE)
437 {
438 if (work->settings->scaling)
439 {
440 vec_ew_prod(work->scaling->D, work->delta_x, work->delta_x, work->data->n);
441 }
442 }
443
444 unscale_data(work);
445
446 work->initialized = FALSE;
447 work->info->iter = iter;
448 work->info->iter_out = iter_out;
449
450 /* Update solve time and run time */
451 #ifdef QPALM_TIMING
452 work->info->solve_time = qpalm_toc(work->timer);
453 work->info->run_time = work->info->setup_time +
454 work->info->solve_time;
455 #endif /* ifdef QPALM_TIMING */
456
457 c = ladel_workspace_free(c);
459 c2 = ladel_workspace_free(c2);
460
461 #ifdef QPALM_PRINTING
462 if (work->settings->verbose)
463 {
464 print_iteration(iter, work);
466 }
467 #endif
468}
469
470static void qpalm_terminate_on_status(QPALMWorkspace *work, solver_common *c, solver_common *c2, c_int iter, c_int iter_out, c_int status_val)
471{
472 update_status(work->info, status_val);
473 qpalm_termination(work, c, c2, iter, iter_out);
474}
475
477{
478 #if defined(QPALM_PRINTING) && defined(_WIN32) && defined(_MSC_VER) && _MSC_VER < 1900
479 unsigned int print_exponent_format = _set_output_format(_TWO_DIGIT_EXPONENT);
480 #endif
481
482 //Initialize ladel workspace, qpalm variables and perform scaling (timing added to setup)
483 solver_common *c, *c2;
484 qpalm_initialize(work, &c, &c2);
485
486 // Start the timer for the solve routine
487 #ifdef QPALM_TIMING
488 qpalm_tic(work->timer);
489 c_float current_time;
490 #endif /* ifdef QPALM_TIMING */
491
492 size_t m = work->data->m;
493 c_int iter;
494 c_int iter_out = 0;
495 c_int prev_iter = 0; /* iteration number at which the previous subproblem finished*/
496 c_float eps_k_abs = work->settings->eps_abs_in;
497 c_float eps_k_rel = work->settings->eps_rel_in;
498 c_int no_change_in_active_constraints = 0;
499
500 for (iter = 0; iter < work->settings->max_iter; iter++)
501 {
502 /* Check whether we passed the time limit */
503 #ifdef QPALM_TIMING
504 current_time = work->info->setup_time + qpalm_toc(work->timer); // Start timer
505 if (current_time > work->settings->time_limit)
506 {
507 qpalm_terminate_on_status(work, c, c2, iter, iter_out, QPALM_TIME_LIMIT_REACHED);
508 return;
509 }
510 #endif /* ifdef QPALM_TIMING */
511
512 /*Perform the iteration */
513 compute_residuals(work, c);
515
516 if (is_solved(work))
517 {
518 qpalm_terminate_on_status(work, c, c2, iter, iter_out, QPALM_SOLVED);
519 return;
520 }
521 else if (is_primal_infeasible(work))
522 {
523 qpalm_terminate_on_status(work, c, c2, iter, iter_out, QPALM_PRIMAL_INFEASIBLE);
524 return;
525 }
526 else if (is_dual_infeasible(work))
527 {
528 qpalm_terminate_on_status(work, c, c2, iter, iter_out, QPALM_DUAL_INFEASIBLE);
529 return;
530 }
531 else if (check_subproblem_termination(work) || (no_change_in_active_constraints == 3))
532 {
533 update_dual_iterate_and_parameters(work, c, iter_out, &eps_k_abs, &eps_k_rel);
534
536 {
539 {
540 qpalm_terminate_on_status(work, c, c2, iter, iter_out, QPALM_DUAL_TERMINATED);
541 return;
542 }
543 }
544
545 no_change_in_active_constraints = 0;
546 iter_out++;
547 prev_iter = iter;
548
549 #ifdef QPALM_PRINTING
550 if (work->settings->verbose && mod(iter, work->settings->print_iter) == 0)
551 {
552 qpalm_print("%4" LADEL_PRIi " | ---------------------------------------------------\n", iter);
553 }
554 #endif
555 }
556 else if (iter == prev_iter + work->settings->inner_max_iter) /*subproblem is hanging so try updating params*/
557 {
558 if (iter_out > 0 && work->info->pri_res_norm > work->eps_pri)
559 {
560 update_sigma(work, c);
561 }
562
563 if(work->settings->proximal)
564 {
565 update_gamma(work);
566 if (!work->settings->nonconvex) prea_vec_copy(work->x, work->x0, work->data->n);
567 }
568
569 prea_vec_copy(work->pri_res, work->pri_res_in, m);
570
571 no_change_in_active_constraints = 0;
572 iter_out++;
573 prev_iter = iter;
574 }
575 else /*primal update*/
576 {
577 if (work->solver->nb_enter + work->solver->nb_leave) no_change_in_active_constraints = 0;
578 else no_change_in_active_constraints++;
579
580 if (mod(iter, work->settings->reset_newton_iter) == 0) work->solver->reset_newton = TRUE;
581 update_primal_iterate(work, c);
582
583 #ifdef QPALM_PRINTING
584 if (work->settings->verbose && mod(iter, work->settings->print_iter) == 0)
585 {
586 work->info->objective = compute_objective(work);
587 print_iteration(iter, work);
588 }
589 #endif
590 }
591 }
592
593 /* If we get here, qpalm has unfortunately hit the maximum number of iterations */
594 qpalm_terminate_on_status(work, c, c2, iter, iter_out, QPALM_MAX_ITER_REACHED);
595 return;
596}
597
599{
600 // If we have previously solved the problem, then reset the setup time
601 if (work->info->status_val != QPALM_UNSOLVED)
602 {
603 #ifdef QPALM_TIMING
604 work->info->setup_time = 0;
605 #endif /* ifdef QPALM_TIMING */
607 }
608 #ifdef QPALM_TIMING
609 qpalm_tic(work->timer); // Start timer
610 #endif /* ifdef QPALM_TIMING */
611
612 // Validate settings
613 if (!validate_settings(settings))
614 {
615 # ifdef QPALM_PRINTING
616 qpalm_eprint("Settings validation returned failure");
617 # endif /* ifdef QPALM_PRINTING */
619 return;
620 }
621
622 // Copy settings
623 qpalm_free(work->settings);
624 work->settings = copy_settings(settings);
625 work->sqrt_delta = c_sqrt(work->settings->delta);
626 # ifdef QPALM_TIMING
627 work->info->setup_time += qpalm_toc(work->timer);
628 # endif /* ifdef QPALM_TIMING */
629}
630
631void qpalm_update_bounds(QPALMWorkspace *work, const c_float *bmin, const c_float *bmax)
632{
633 // If we have previously solved the problem, then reset the setup time
634 if (work->info->status_val != QPALM_UNSOLVED)
635 {
636 #ifdef QPALM_TIMING
637 work->info->setup_time = 0;
638 #endif /* ifdef QPALM_TIMING */
640 }
641 #ifdef QPALM_TIMING
642 qpalm_tic(work->timer); // Start timer
643 #endif /* ifdef QPALM_TIMING */
644
645 // Validate bounds
646 size_t j;
647 size_t m = work->data->m;
648 if (bmin != NULL && bmax != NULL)
649 {
650 for (j = 0; j < m; j++)
651 {
652 if (bmin[j] > bmax[j])
653 {
654 # ifdef QPALM_PRINTING
655 qpalm_eprint("Lower bound at index %d is greater than upper bound: %.4e > %.4e",
656 (int)j, work->data->bmin[j], work->data->bmax[j]);
657 # endif /* ifdef QPALM_PRINTING */
659 return;
660 }
661 }
662 }
663
664 if (bmin != NULL)
665 {
666 prea_vec_copy(bmin, work->data->bmin, m);
667 }
668 if (bmax != NULL)
669 {
670 prea_vec_copy(bmax, work->data->bmax, m);
671 }
672
673 # ifdef QPALM_TIMING
674 work->info->setup_time += qpalm_toc(work->timer);
675 # endif /* ifdef QPALM_TIMING */
676}
677
679{
680 // If we have previously solved the problem, then reset the setup time
681 if (work->info->status_val != QPALM_UNSOLVED)
682 {
683 #ifdef QPALM_TIMING
684 work->info->setup_time = 0;
685 #endif /* ifdef QPALM_TIMING */
687 }
688 #ifdef QPALM_TIMING
689 qpalm_tic(work->timer); // Start timer
690 #endif /* ifdef QPALM_TIMING */
691
692 size_t n = work->data->n;
693 prea_vec_copy(q, work->data->q, n);
694 # ifdef QPALM_TIMING
695 work->info->setup_time += qpalm_toc(work->timer);
696 # endif /* ifdef QPALM_TIMING */
697}
698
699void qpalm_update_Q_A(QPALMWorkspace *work, const c_float *Qx, const c_float *Ax)
700{
702 // If we have previously solved the problem, then reset the setup time
703 if (work->info->status_val != QPALM_UNSOLVED)
704 {
705 #ifdef QPALM_TIMING
706 work->info->setup_time = 0;
707 #endif /* ifdef QPALM_TIMING */
709 }
710 #ifdef QPALM_TIMING
711 qpalm_tic(work->timer); // Start timer
712 #endif /* ifdef QPALM_TIMING */
713
714 ladel_sparse_matrix *Q = work->data->Q, *A = work->data->A;
715 prea_vec_copy(Qx, Q->x, Q->nzmax);
716 prea_vec_copy(Ax, A->x, A->nzmax);
717
718 # ifdef QPALM_TIMING
719 work->info->setup_time += qpalm_toc(work->timer);
720 # endif /* ifdef QPALM_TIMING */
721}
722
724{
725 if (work)
726 {
727 // Free Data
728 if (work->data)
729 {
730 work->data->Q = ladel_sparse_free(work->data->Q);
731
732 work->data->A = ladel_sparse_free(work->data->A);
733
734 if (work->data->q) qpalm_free(work->data->q);
735
736 if (work->data->bmin) qpalm_free(work->data->bmin);
737
738 if (work->data->bmax) qpalm_free(work->data->bmax);
739 qpalm_free(work->data);
740 }
741
742 // Free scaling
743 if (work->scaling->D) qpalm_free(work->scaling->D);
744
745 if (work->scaling->Dinv) qpalm_free(work->scaling->Dinv);
746
747 if (work->scaling->E) qpalm_free(work->scaling->E);
748
749 if (work->scaling->Einv) qpalm_free(work->scaling->Einv);
750
751 qpalm_free(work->scaling);
752
753 // Free other Variables
754 if (work->x) qpalm_free(work->x);
755
756 if (work->y) qpalm_free(work->y);
757
758 if (work->Ax) qpalm_free(work->Ax);
759
760 if (work->Qx) qpalm_free(work->Qx);
761
762 if (work->x_prev) qpalm_free(work->x_prev);
763
764 if (work->Aty) qpalm_free(work->Aty);
765
766 if (work->temp_m) qpalm_free(work->temp_m);
767
768 if (work->temp_n) qpalm_free(work->temp_n);
769
770 if (work->sigma) qpalm_free(work->sigma);
771
772 if (work->sigma_inv) qpalm_free(work->sigma_inv);
773
774 if (work->z) qpalm_free(work->z);
775
776 if (work->Axys) qpalm_free(work->Axys);
777
778 if (work->pri_res) qpalm_free(work->pri_res);
779
780 if (work->pri_res_in) qpalm_free(work->pri_res_in);
781
782 if (work->df) qpalm_free(work->df);
783
784 if (work->x0) qpalm_free(work->x0);
785
786 if (work->xx0) qpalm_free(work->xx0);
787
788 if (work->dphi) qpalm_free(work->dphi);
789
790 if (work->dphi_prev) qpalm_free(work->dphi_prev);
791
792 if (work->sqrt_sigma) qpalm_free(work->sqrt_sigma);
793
794 if (work->delta) qpalm_free(work->delta);
795
796 if (work->alpha) qpalm_free(work->alpha);
797
798 if (work->delta2) qpalm_free(work->delta2);
799
800 if (work->delta_alpha) qpalm_free(work->delta_alpha);
801
802 if (work->temp_2m) qpalm_free(work->temp_2m);
803
804 if (work->s) qpalm_free(work->s);
805
806 if (work->index_L) qpalm_free(work->index_L);
807
808 if (work->index_P) qpalm_free(work->index_P);
809
810 if (work->index_J) qpalm_free(work->index_J);
811
812 if (work->delta_y) qpalm_free(work->delta_y);
813
814 if (work->Atdelta_y) qpalm_free(work->Atdelta_y);
815
816 if (work->delta_x) qpalm_free(work->delta_x);
817
818 if (work->Qdelta_x) qpalm_free(work->Qdelta_x);
819
820 if (work->Adelta_x) qpalm_free(work->Adelta_x);
821
822 // Free Settings
823 if (work->settings) qpalm_free(work->settings);
824
825 //Free chol struct
826 if (work->solver)
827 {
829
831
832 if (work->solver->enter) qpalm_free(work->solver->enter);
833
834 if (work->solver->leave) qpalm_free(work->solver->leave);
835
836
837 work->solver->sol_kkt = ladel_free(work->solver->sol_kkt);
838
839 work->solver->rhs_kkt = ladel_free(work->solver->rhs_kkt);
840
841 work->solver->D_temp = ladel_free(work->solver->D_temp);
842
843 work->solver->E_temp = ladel_free(work->solver->E_temp);
844
845 work->solver->neg_dphi = ladel_free(work->solver->neg_dphi);
846
847 work->solver->d = ladel_free(work->solver->d);
848
849 work->solver->Qd = ladel_free(work->solver->Qd);
850
851 work->solver->Ad = ladel_free(work->solver->Ad);
852
853 work->solver->yh = ladel_free(work->solver->yh);
854
855 work->solver->Atyh = ladel_free(work->solver->Atyh);
856
857 work->solver->LD = ladel_factor_free(work->solver->LD);
858
859 work->solver->LD_Q = ladel_factor_free(work->solver->LD_Q);
860
861 work->solver->sym = ladel_symbolics_free(work->solver->sym);
862
863 work->solver->sym_Q = ladel_symbolics_free(work->solver->sym_Q);
864
865 work->solver->At_scale = ladel_free(work->solver->At_scale);
866
867 work->solver->At_sqrt_sigma = ladel_sparse_free(work->solver->At_sqrt_sigma);
868
869 work->solver->At = ladel_sparse_free(work->solver->At);
870
871 work->solver->kkt = ladel_sparse_free(work->solver->kkt);
872
873 work->solver->kkt_full = ladel_sparse_free(work->solver->kkt_full);
874
875 work->solver->first_row_A = ladel_free(work->solver->first_row_A);
876
877 work->solver->first_elem_A = ladel_free(work->solver->first_elem_A);
878
879
880 qpalm_free(work->solver);
881 }
882
883 // Free solution
884 if (work->solution)
885 {
886 if (work->solution->x) qpalm_free(work->solution->x);
887
888 if (work->solution->y) qpalm_free(work->solution->y);
889 qpalm_free(work->solution);
890 }
891
892 // Free timer
893 # ifdef QPALM_TIMING
894 if (work->timer) qpalm_free(work->timer);
895 # endif /* ifdef QPALM_TIMING */
896
897 // Free information
898 if (work->info) qpalm_free(work->info);
899
900 // Free work
901 qpalm_free(work);
902 }
903
904}
905
906
907# ifdef __cplusplus
908}
909# endif // ifdef __cplusplus
Constants used in QPALM.
#define ORDERING
ordering in the factorization
Definition constants.h:113
#define EPS_REL_IN
default intermediate relative convergence tolerance
Definition constants.h:70
#define SCALING
default number of scaling iterations
Definition constants.h:83
#define NONCONVEX
default use of nonconvex adjustments
Definition constants.h:87
#define QPALM_MAX_ITER_REACHED
status to indicate termination due to reaching the maximum number of iterations
Definition constants.h:32
#define WARM_START
default warm start setting
Definition constants.h:88
#define GAMMA_UPD
default proximal penalty update factor
Definition constants.h:80
#define PRINT_ITER
default frequency of printing
Definition constants.h:90
#define EPS_ABS
default absolute convergence tolerance
Definition constants.h:67
#define QPALM_ERROR
status to indicate an error has occured (this error should automatically be printed)
Definition constants.h:37
#define FACTORIZE_SCHUR
factorize the Schur complement
Definition constants.h:107
#define RHO
default tolerance scaling factor
Definition constants.h:71
#define DELTA
default penalty update factor
Definition constants.h:75
#define VERBOSE
default write out progress setting
Definition constants.h:89
#define INNER_MAX_ITER
default maximum number of iterations per subproblem
Definition constants.h:66
#define QPALM_DUAL_INFEASIBLE
status to indicate the problem is dual infeasible
Definition constants.h:34
#define RESET_NEWTON_ITER
default frequency of performing a full Cholesky factorization
Definition constants.h:92
#define GAMMA_INIT
default initial proximal penalty parameter
Definition constants.h:79
#define SIGMA_INIT
default initial penalty parameter (guideline)
Definition constants.h:77
#define EPS_REL
default relative convergence tolerance
Definition constants.h:68
#define MAX_RANK_UPDATE
maximum rank for the sparse factorization update
Definition constants.h:98
#define PROXIMAL
default use of proximal method of multipliers
Definition constants.h:78
#define QPALM_TIME_LIMIT_REACHED
status to indicate the problem's runtime has exceeded the specified time limit
Definition constants.h:35
#define ENABLE_DUAL_TERMINATION
enable termination after dual objective > something (useful in branch and bound)
Definition constants.h:94
#define FACTORIZATION_METHOD
default method for solving the linear system
Definition constants.h:110
#define QPALM_PRIMAL_INFEASIBLE
status to indicate the problem is primal infeasible
Definition constants.h:33
#define MAX_RANK_UPDATE_FRACTION
maximum rank (relative to n+m) for the factorization update
Definition constants.h:99
#define QPALM_SOLVED
status to indicate the problem is solved to optimality given the specified tolerances
Definition constants.h:30
#define QPALM_DUAL_TERMINATED
status to indicate the problem has a dual objective that is higher than the specified bound
Definition constants.h:31
#define TRUE
Definition constants.h:18
#define FALSE
Definition constants.h:19
#define THETA
default penalty update criterion parameter
Definition constants.h:74
#define QPALM_NULL
NULL, if something goes wrong during setup, the workspace pointer is set to this.
Definition constants.h:53
#define SIGMA_MAX
default penalty cap
Definition constants.h:76
#define FACTORIZE_KKT
factorize the kkt system
Definition constants.h:106
#define TIME_LIMIT
time limit after which the solver aborts
Definition constants.h:96
#define QPALM_UNSOLVED
status to indicate the problem is unsolved.
Definition constants.h:36
#define GAMMA_MAX
default proximal penalty cap
Definition constants.h:81
#define MAX_ITER
default maximum number of iterations
Definition constants.h:65
#define EPS_PRIM_INF
default primal infeasibility tolerance
Definition constants.h:72
#define EPS_ABS_IN
default intermediate absolute convergence tolerance
Definition constants.h:69
#define QPALM_INFTY
infinity, used to indicate one-sided constraints
Definition constants.h:61
#define DUAL_OBJECTIVE_LIMIT
termination value for the dual objective (useful in branch and bound)
Definition constants.h:95
#define EPS_DUAL_INF
default dual infeasibility tolerance
Definition constants.h:73
void qpalm_free(void *ptr)
Definition global_opts.c:16
void * qpalm_malloc(size_t size)
Definition global_opts.c:7
void * qpalm_calloc(size_t num, size_t size)
Definition global_opts.c:3
Custom memory allocation, print and utility functions, and data types for floats and ints.
#define c_sqrt
square root
ladel_int c_int
type for integer numbers
Definition global_opts.h:42
#define qpalm_print
Definition global_opts.h:75
#define mod(a, b)
modulo operation (positive result for all values)
ladel_double c_float
type for floating point numbers
Definition global_opts.h:41
#define qpalm_eprint(...)
Definition global_opts.h:81
void qpalm_solve(QPALMWorkspace *work)
Solve the quadratic program.
Definition qpalm.c:476
QPALMWorkspace * qpalm_setup(const QPALMData *data, const QPALMSettings *settings)
Initialize QPALM solver allocating memory.
Definition qpalm.c:68
void qpalm_update_settings(QPALMWorkspace *work, const QPALMSettings *settings)
Update the settings to the new settings.
Definition qpalm.c:598
void qpalm_set_default_settings(QPALMSettings *settings)
Set default settings from constants.h file.
Definition qpalm.c:33
void qpalm_update_bounds(QPALMWorkspace *work, const c_float *bmin, const c_float *bmax)
Update the lower and upper bounds.
Definition qpalm.c:631
void qpalm_update_Q_A(QPALMWorkspace *work, const c_float *Qx, const c_float *Ax)
Update the matrix entries of Q and A.
Definition qpalm.c:699
void qpalm_warm_start(QPALMWorkspace *work, const c_float *x_warm_start, const c_float *y_warm_start)
Warm start workspace variables x, x_0, x_prev, Ax, Qx, y and sigma.
Definition qpalm.c:259
void qpalm_cleanup(QPALMWorkspace *work)
Cleanup the workspace by deallocating memory.
Definition qpalm.c:723
void qpalm_update_q(QPALMWorkspace *work, const c_float *q)
Update the linear part of the cost.
Definition qpalm.c:678
void initialize_sigma(QPALMWorkspace *work, solver_common *c)
Initialize penalty factors from initial x.
Definition iteration.c:48
c_float compute_dual_objective(QPALMWorkspace *work, solver_common *c)
Compute the (unscaled) dual objective value at the current iterate.
Definition iteration.c:328
void update_dual_iterate_and_parameters(QPALMWorkspace *work, solver_common *c, c_int iter_out, c_float *eps_k_abs, c_float *eps_k_rel)
Definition iteration.c:249
void update_gamma(QPALMWorkspace *work)
Update the proximal penalty.
Definition iteration.c:126
void update_sigma(QPALMWorkspace *work, solver_common *c)
Update the penalty factors.
Definition iteration.c:73
void update_primal_iterate(QPALMWorkspace *work, solver_common *c)
Update the primal iterate.
Definition iteration.c:269
c_float compute_objective(QPALMWorkspace *work)
Compute the (unscaled) primal objective value at the current iterate.
Definition iteration.c:287
void compute_residuals(QPALMWorkspace *work, solver_common *c)
Compute the residuals (in vector form)
Definition iteration.c:22
QPALM main solver routines.
void vec_add_scaled(const c_float *a, const c_float *b, c_float *c, c_float sc, size_t n)
Scaled addition of one vector to another vector, .
Definition lin_alg.c:110
c_float * vec_copy(const c_float *a, size_t n)
Copy vector a into output.
Definition lin_alg.c:11
void vec_set_scalar(c_float *a, c_float sc, size_t n)
Fill float vector with a scalar value.
Definition lin_alg.c:40
void vec_set_scalar_int(c_int *a, c_int sc, size_t n)
Fill int vector with a scalar value.
Definition lin_alg.c:48
void vec_self_mult_scalar(c_float *a, c_float sc, size_t n)
Mulitply vector with a constant scale factor.
Definition lin_alg.c:56
void vec_ew_prod(const c_float *a, const c_float *b, c_float *c, size_t n)
Elementwise product, .
Definition lin_alg.c:92
void prea_vec_copy(const c_float *a, c_float *b, size_t n)
Copy vector a into preallocated vector b.
Definition lin_alg.c:24
Linear algebra with vectors.
Routines to perform exact linesearch.
Functions to calculate the semismooth Newton direction.
void set_settings_nonconvex(QPALMWorkspace *work, solver_common *c)
Set the proximal parameters for nonconvex QPs.
Definition nonconvex.c:331
Routines to deal with nonconvex QPs.
QPALM main solver API.
void scale_data(QPALMWorkspace *work)
Scale problem matrices.
Definition scaling.c:31
void unscale_data(QPALMWorkspace *work)
Unscale the problem data.
Definition scaling.c:99
Problem data scaling during setup.
void mat_vec(solver_sparse *A, solver_dense *x, solver_dense *y, solver_common *c)
Matrix-vector multiplication.
void qpalm_set_factorization_method(QPALMWorkspace *work, solver_common *c)
Choose the linear systems solver method based on the problem data sizes.
Interface and wrapper to matrix/factorization (ladel) functions.
Data structure.
Definition types.h:102
size_t m
number of constraints m
Definition types.h:104
c_float * bmin
dense array for lower bounds (size m)
Definition types.h:109
c_float c
constant part of cost
Definition types.h:108
size_t n
number of variables n
Definition types.h:103
c_float * q
dense array for linear part of cost function (size n)
Definition types.h:107
solver_sparse * A
sparse linear constraints matrix A (size m x n)
Definition types.h:106
c_float * bmax
dense array for upper bounds (size m)
Definition types.h:110
solver_sparse * Q
sparse quadratic part of the cost Q (size n x n)
Definition types.h:105
Solver return information.
Definition types.h:74
c_float run_time
total time (seconds)
Definition types.h:90
c_int iter_out
number of outer iterations (i.e. dual updates)
Definition types.h:76
c_float pri_res_norm
norm of primal residual
Definition types.h:80
c_float dual_objective
dual objective function value (= NaN if enable_dual_termination is false)
Definition types.h:85
c_int iter
number of iterations taken
Definition types.h:75
c_float solve_time
time taken for solve phase (seconds)
Definition types.h:89
c_float objective
objective function value
Definition types.h:84
c_int status_val
status as c_int, defined in constants.h
Definition types.h:78
c_float setup_time
time taken for setup phase (seconds)
Definition types.h:88
Problem scaling matrices stored as vectors.
Definition types.h:61
c_float * Dinv
primal variable rescaling
Definition types.h:63
c_float * E
dual variable scaling
Definition types.h:64
c_float * Einv
dual variable rescaling
Definition types.h:65
c_float cinv
objective rescaling
Definition types.h:67
c_float * D
primal variable scaling
Definition types.h:62
Settings struct.
Definition types.h:117
c_float gamma_upd
proximal penalty update factor
Definition types.h:133
c_float sigma_max
penalty factor cap
Definition types.h:129
c_float eps_abs_in
intermediate absolute convergence tolerance
Definition types.h:122
c_float gamma_max
proximal penalty parameter cap
Definition types.h:134
c_int proximal
boolean, use proximal method of multipliers or not
Definition types.h:131
c_float delta
penalty update factor
Definition types.h:128
c_int warm_start
boolean, warm start
Definition types.h:139
c_float sigma_init
initial penalty parameter (guideline)
Definition types.h:130
c_float eps_dual_inf
dual infeasibility tolerance
Definition types.h:126
c_int reset_newton_iter
frequency of performing a complete Cholesky factorization
Definition types.h:140
c_float dual_objective_limit
termination value for the dual objective (useful in branch and bound)
Definition types.h:142
c_float eps_rel_in
intermediate relative convergence tolerance
Definition types.h:123
c_float rho
tolerance scaling factor
Definition types.h:124
c_float time_limit
time limit
Definition types.h:143
c_float theta
penalty update criterion parameter
Definition types.h:127
c_int max_rank_update
maximum rank for the sparse factorization update
Definition types.h:146
c_int enable_dual_termination
boolean, enable termination based on dual objective (useful in branch and bound)
Definition types.h:141
c_float max_rank_update_fraction
maximum rank (relative to n+m) for the factorization update
Definition types.h:147
c_float gamma_init
initial proximal penalty parameter
Definition types.h:132
c_float eps_prim_inf
primal infeasibility tolerance
Definition types.h:125
c_int inner_max_iter
maximum number of iterations per subproblem
Definition types.h:119
c_float eps_rel
relative convergence tolerance
Definition types.h:121
c_int verbose
boolean, write out progress
Definition types.h:137
c_int max_iter
maximum number of iterations
Definition types.h:118
c_int scaling
scaling iterations, if 0 then scaling is disabled
Definition types.h:135
c_int nonconvex
boolean, indicates whether the QP is nonconvex
Definition types.h:136
c_float eps_abs
absolute convergence tolerance
Definition types.h:120
c_int ordering
ordering method for factorization
Definition types.h:144
c_int factorization_method
factorize KKT or Schur complement
Definition types.h:145
c_int print_iter
frequency of printing
Definition types.h:138
Solution structure.
Definition types.h:43
c_float * y
dual solution
Definition types.h:45
c_float * x
primal solution
Definition types.h:44
Variables for linear system solving.
Definition types.h:155
solver_dense * E_temp
temporary constraints scaling vectors
Definition types.h:166
solver_dense * neg_dphi
-gradient of the Lagrangian
Definition types.h:168
solver_dense * D_temp
temporary primal variable scaling vectors
Definition types.h:167
solver_dense * Ad
A * d.
Definition types.h:172
solver_sparse * At_sqrt_sigma
A' * sqrt(sigma)
Definition types.h:186
c_int * first_row_A
row index of the first element in each column of A'
Definition types.h:160
c_int * active_constraints
index set of active constraints
Definition types.h:178
solver_dense * yh
candidate dual update
Definition types.h:174
c_int reset_newton
boolean, after sigma is updated perform a new factorization
Definition types.h:177
solver_dense * d
primal update step
Definition types.h:171
c_int * leave
index set of leaving constraints
Definition types.h:183
solver_dense * rhs_kkt
[-dphi; zeros(m,1)]
Definition types.h:169
solver_factor * LD
LD factor (part of LDL' factorization)
Definition types.h:162
c_int factorization_method
factorize KKT or Schur complement
Definition types.h:156
c_int * enter
index set of entering constraints
Definition types.h:181
solver_sparse * At
A'.
Definition types.h:159
c_int * active_constraints_old
index set of active constraints in the previous iteration
Definition types.h:179
c_int first_factorization
boolean, indicate we have not factorized previously
Definition types.h:176
solver_dense * At_scale
running vector of sqrt(sigma), used to scale At_sqrt_sigma
Definition types.h:185
solver_dense * Atyh
A' * yh.
Definition types.h:175
solver_dense * Qd
Q * d.
Definition types.h:173
solver_dense * sol_kkt
sol_kkt = kkt \ rhs_kkt
Definition types.h:170
solver_factor * LD_Q
LD factor of Q (useful in computing dual objective)
Definition types.h:164
c_float * first_elem_A
value of the first element in each column of A'
Definition types.h:161
c_int nb_leave
number of leaving constraints
Definition types.h:184
solver_symbolics * sym_Q
symbolics for Q (only used in LADEL)
Definition types.h:165
solver_symbolics * sym
symbolics for kkt (only used in LADEL)
Definition types.h:163
solver_sparse * kkt_full
KKT matrix if all constraints would be active.
Definition types.h:158
solver_sparse * kkt
KKT matrix.
Definition types.h:157
c_int nb_enter
number of entering constraints
Definition types.h:182
QPALM Workspace.
Definition types.h:197
c_float * delta_y
difference of consecutive dual iterates
Definition types.h:279
c_float * x
primal iterate
Definition types.h:204
c_float * y
dual iterate
Definition types.h:205
c_float * yh
candidate dual update
Definition types.h:229
QPALMScaling * scaling
scaling vectors
Definition types.h:306
c_float * delta
linesearch parameter
Definition types.h:252
c_int nb_sigma_changed
number of sigma-components that changed in an outer iteration (relevant for factorization update)
Definition types.h:222
c_int initialized
flag whether the iterates were initialized or not
Definition types.h:210
c_float * D_temp
temporary primal variable scaling vectors
Definition types.h:299
c_float * delta2
delta .* delta
Definition types.h:255
c_float eps_pri
primal tolerance
Definition types.h:268
c_float * Ax
scaled A * x
Definition types.h:206
QPALMInfo * info
solver information
Definition types.h:308
c_float * z
projection of Axys onto the constraint set [bmin, bmax]
Definition types.h:226
c_float * pri_res_in
intermediate primal residual
Definition types.h:228
c_int * index_J
index set J (L xor P)
Definition types.h:260
c_float * temp_m
placeholder for vector of size m
Definition types.h:217
QPALMSolution * solution
problem solution
Definition types.h:307
c_int * index_P
index set P (where delta>0)
Definition types.h:259
c_float * pri_res
primal residual
Definition types.h:227
c_float * sqrt_sigma
elementwise sqrt(sigma)
Definition types.h:248
c_float * sigma
penalty vector
Definition types.h:219
array_element * s
alpha ./ delta
Definition types.h:257
c_float * delta_alpha
delta .* alpha
Definition types.h:256
c_float * x_prev
previous primal iterate
Definition types.h:209
c_float * Adelta_x
A * delta_x.
Definition types.h:290
c_float * neg_dphi
-dphi, required as the rhs in SCHUR
Definition types.h:235
c_float * Axys
Ax + y./sigma.
Definition types.h:225
c_float * Qdelta_x
Q * delta_x.
Definition types.h:289
c_int gamma_maxed
flag to indicate whether gamma has been maximized when the primal residual was low
Definition types.h:224
c_float * delta_x
difference of consecutive primal iterates
Definition types.h:288
QPALMTimer * timer
timer object
Definition types.h:311
c_float * Qx
scaled Q * x
Definition types.h:207
c_float * dphi_prev
previous gradient of the Lagrangian
Definition types.h:236
c_float eps_rel_in
intermediate relative tolerance
Definition types.h:272
QPALMSettings * settings
problem settings
Definition types.h:305
c_float * x0
record of the primal iterate during the last dual update
Definition types.h:232
c_float * sigma_inv
1./sigma
Definition types.h:220
c_float gamma
proximal penalty factor
Definition types.h:223
c_float * dphi
gradient of the Lagrangian
Definition types.h:234
c_float * E_temp
temporary constraints scaling vectors
Definition types.h:300
c_float * df
gradient of the primal objective (+proximal term)
Definition types.h:231
c_float * alpha
linesearch parameter
Definition types.h:253
QPALMSolver * solver
linsys variables
Definition types.h:304
c_float sqrt_delta
sqrt(penalty update factor)
Definition types.h:249
c_float * xx0
x - x0
Definition types.h:233
c_float * Atdelta_y
A' * delta_y.
Definition types.h:280
c_float * temp_2m
placeholder for vector of size 2m
Definition types.h:254
c_float * Ad
A * d.
Definition types.h:247
c_float * Aty
A' * y (useful for saving one mat_tpose_vec)
Definition types.h:208
c_float eps_abs_in
intermediate absolute tolerance
Definition types.h:271
QPALMData * data
problem data to work on (possibly scaled)
Definition types.h:198
c_int * index_L
index set L (where s>0)
Definition types.h:258
c_float * Qd
Q * d.
Definition types.h:246
c_float * temp_n
placeholder for vector of size n
Definition types.h:218
c_float * d
primal update step
Definition types.h:237
c_float * Atyh
A' * yh.
Definition types.h:230
Array to sort in linesearch.
Definition types.h:27
c_int is_primal_infeasible(QPALMWorkspace *work)
Check whether the problem is primal infeasible.
void calculate_residual_norms_and_tolerances(QPALMWorkspace *work)
Calls the routines that compute the norm of the residuals and the tolerances.
Definition termination.c:19
void store_solution(QPALMWorkspace *work)
Helper function to store the (unscaled) solution in the solution struct.
c_int is_solved(QPALMWorkspace *work)
Check whether the primal and dual residual norms are smaller than the primal and dual tolerances.
c_int is_dual_infeasible(QPALMWorkspace *work)
Check whether the problem is dual infeasible.
c_int check_subproblem_termination(QPALMWorkspace *work)
Check whether the subproblem is solved.
Routines to check the termination and infeasibility criteria.
ladel_work solver_common
Definition types.h:18
struct QPALM_TIMER QPALMTimer
QPALM Timer for statistics.
Definition types.h:53
QPALMSettings * copy_settings(const QPALMSettings *settings)
Copy settings creating a new settings structure.
Definition util.c:25
void update_status(QPALMInfo *info, c_int status_val)
Update solver status (value and string).
Definition util.c:62
void print_iteration(c_int iter, QPALMWorkspace *work)
Print information about the current iteration.
Definition util.c:114
void print_header(void)
Print the header with QPALM version number and fields.
Definition util.c:108
void print_final_message(QPALMWorkspace *work)
Print final message as a box with info.
Definition util.c:122
Utility functions.
c_float qpalm_toc(QPALMTimer *t)
Report time in seconds since last call to qpalm_tic.
void qpalm_tic(QPALMTimer *t)
Start timer.
c_int validate_data(const QPALMData *data)
Validate problem data.
Definition validate.c:18
c_int validate_settings(const QPALMSettings *settings)
Validate problem settings.
Definition validate.c:43
Validation of the user provided settings and data.