QPALM 1.0.0
Proximal Augmented Lagrangian method for Quadratic Programs
types.h
Go to the documentation of this file.
1/**
2 * @file types.h
3 * @author Ben Hermans
4 * @brief Internal data structures used in QPALM.
5 * @details See data structures tab.
6 */
7
8#ifndef QPALM_TYPES_H
9# define QPALM_TYPES_H
10
11# ifdef __cplusplus
12extern "C" {
13# endif
14
15#include "global_opts.h"
16
17#include "ladel.h"
18typedef ladel_work solver_common;
19typedef ladel_sparse_matrix solver_sparse;
20typedef ladel_double solver_dense;
21typedef ladel_factor solver_factor;
22typedef ladel_symbolics solver_symbolics;
23
24/**
25 * Array to sort in linesearch
26 */
27typedef struct array_element {
28 c_float x; ///< value of the element
29 size_t i; ///< index
31
32
33/// @addtogroup solver-grp
34/// @{
35
36/******************
37* Internal types *
38******************/
39
40/**
41 * Solution structure
42 */
43typedef struct {
44 c_float *x; ///< primal solution
45 c_float *y; ///< dual solution
47
48/// @}
49
50/**
51 * QPALM Timer for statistics
52 */
53typedef struct QPALM_TIMER QPALMTimer;
54
55/// @addtogroup solver-grp
56/// @{
57
58/**
59 * Problem scaling matrices stored as vectors
60 */
61typedef struct {
62 c_float *D; ///< primal variable scaling
63 c_float *Dinv; ///< primal variable rescaling
64 c_float *E; ///< dual variable scaling
65 c_float *Einv; ///< dual variable rescaling
66 c_float c; ///< objective scaling
67 c_float cinv; ///< objective rescaling
69
70
71/**
72 * Solver return information
73 */
74typedef struct {
75 c_int iter; ///< number of iterations taken
76 c_int iter_out; ///< number of outer iterations (i.e. dual updates)
77 char status[32]; ///< status string, e.g. 'solved'
78 c_int status_val; ///< status as c_int, defined in constants.h
79
80 c_float pri_res_norm; ///< norm of primal residual
81 c_float dua_res_norm; ///< norm of dual residual
82 c_float dua2_res_norm; ///< norm of intermediate dual residual (minus proximal term)
83
84 c_float objective; ///< objective function value
85 c_float dual_objective; ///< dual objective function value (= NaN if enable_dual_termination is false)
86
87 #ifdef PROFILING
88 c_float setup_time; ///< time taken for setup phase (seconds)
89 c_float solve_time; ///< time taken for solve phase (seconds)
90 c_float run_time; ///< total time (seconds)
91 #endif
92
93} QPALMInfo;
94
95/**********************************
96* Main structures and Data Types *
97**********************************/
98
99/**
100 * Data structure
101 */
102typedef struct {
103 size_t n; ///< number of variables n
104 size_t m; ///< number of constraints m
105 solver_sparse *Q; ///< sparse quadratic part of the cost Q (size n x n)
106 solver_sparse *A; ///< sparse linear constraints matrix A (size m x n)
107 c_float *q; ///< dense array for linear part of cost function (size n)
108 c_float c; ///< constant part of cost
109 c_float *bmin; ///< dense array for lower bounds (size m)
110 c_float *bmax; ///< dense array for upper bounds (size m)
111} QPALMData;
112
113
114/**
115 * Settings struct
116 */
117typedef struct {
118 c_int max_iter; ///< maximum number of iterations @details @note Assumption: @f$>0@f$
119 c_int inner_max_iter; ///< maximum number of iterations per subproblem @details @note Assumption: @f$>0@f$
120 c_float eps_abs; ///< absolute convergence tolerance @details @note Assumption: @f$>=0@f$, either eps_abs or eps_rel must be @f$>0@f$
121 c_float eps_rel; ///< relative convergence tolerance @details @note Assumption: @f$>=0@f$, either eps_abs or eps_rel must be @f$>0@f$
122 c_float eps_abs_in; ///< intermediate absolute convergence tolerance @details @note Assumption: @f$>=0@f$, either eps_abs_in or eps_rel_in must be @f$>0@f$
123 c_float eps_rel_in; ///< intermediate relative convergence tolerance @details @note Assumption: @f$>=0@f$, either eps_abs_in or eps_rel_in must be @f$>0@f$
124 c_float rho; ///< tolerance scaling factor @details @note Assumption: @f$0<\rho<1@f$
125 c_float eps_prim_inf; ///< primal infeasibility tolerance @details @note Assumption: @f$>=0@f$
126 c_float eps_dual_inf; ///< dual infeasibility tolerance @details @note Assumption: @f$>=0@f$
127 c_float theta; ///< penalty update criterion parameter @details @note Assumption: @f$<=1@f$
128 c_float delta; ///< penalty update factor @details @note Assumption: @f$>1@f$
129 c_float sigma_max; ///< penalty factor cap @details @note Assumption: @f$>0@f$
130 c_float sigma_init; ///< initial penalty parameter (guideline) @note Assumption: @f$>0@f$
131 c_int proximal; ///< boolean, use proximal method of multipliers or not @details @note Assumption: @f$\in \{0,1\}@f$
132 c_float gamma_init; ///< initial proximal penalty parameter @details @note Assumption: @f$>0@f$
133 c_float gamma_upd; ///< proximal penalty update factor @details @note Assumption: @f$>=1@f$
134 c_float gamma_max; ///< proximal penalty parameter cap @details @note Assumption: @f$>=\gamma_\textrm{init}@f$
135 c_int scaling; ///< scaling iterations, if 0 then scaling is disabled @details @note Assumption: @f$>=0@f$
136 c_int nonconvex; ///< boolean, indicates whether the QP is nonconvex @details @note Assumption: @f$\in \{0,1\}@f$
137 c_int verbose; ///< boolean, write out progress @details @note Assumption:@f$\in \{0,1\}@f$
138 c_int print_iter; ///< frequency of printing @details @note Assumption: @f$>0@f$
139 c_int warm_start; ///< boolean, warm start @details @note Assumption: @f$\in \{0,1\}@f$
140 c_int reset_newton_iter; ///< frequency of performing a complete Cholesky factorization @details @note Assumption: @f$>0@f$
141 c_int enable_dual_termination; ///< boolean, enable termination based on dual objective (useful in branch and bound) @details @note Assumption: @f$\in \{0,1\}@f$
142 c_float dual_objective_limit; ///< termination value for the dual objective (useful in branch and bound) @details @note Assumption: none
143 c_float time_limit; ///< time limit @details @note Assumption: @f$>0@f$
144 c_int ordering; ///< ordering method for factorization
145 c_int factorization_method; ///< factorize KKT or Schur complement
146 c_int max_rank_update; ///< maximum rank for the sparse factorization update
147 c_float max_rank_update_fraction; ///< maximum rank (relative to n+m) for the factorization update
149
150/// @}
151
152/**
153 * Variables for linear system solving
154 */
155typedef struct {
156 c_int factorization_method; ///< factorize KKT or Schur complement
157 solver_sparse *kkt; ///< KKT matrix
158 solver_sparse *kkt_full; ///< KKT matrix if all constraints would be active
160 c_int *first_row_A; ///< row index of the first element in each column of A'
161 c_float *first_elem_A; ///< value of the first element in each column of A'
162 solver_factor *LD; ///< LD factor (part of LDL' factorization)
163 solver_symbolics *sym; ///< symbolics for kkt (only used in LADEL)
164 solver_factor *LD_Q; ///< LD factor of Q (useful in computing dual objective)
165 solver_symbolics *sym_Q; ///< symbolics for Q (only used in LADEL)
166 solver_dense *E_temp; ///< temporary constraints scaling vectors
167 solver_dense *D_temp; ///< temporary primal variable scaling vectors
168 solver_dense *neg_dphi; ///< -gradient of the Lagrangian
169 solver_dense *rhs_kkt; ///< [-dphi; zeros(m,1)]
170 solver_dense *sol_kkt; ///< sol_kkt = kkt \\ rhs_kkt
171 solver_dense *d; ///< primal update step
172 solver_dense *Ad; ///< A * d
173 solver_dense *Qd; ///< Q * d
174 solver_dense *yh; ///< candidate dual update
175 solver_dense *Atyh; ///< A' * yh
176 c_int first_factorization; ///< boolean, indicate we have not factorized previously
177 c_int reset_newton; ///< boolean, after sigma is updated perform a new factorization
178 c_int *active_constraints; ///< index set of active constraints
179 c_int *active_constraints_old; ///< index set of active constraints in the previous iteration
180 c_int nb_active_constraints; ///< number of active constraints
181 c_int *enter; ///< index set of entering constraints
182 c_int nb_enter; ///< number of entering constraints
183 c_int *leave; ///< index set of leaving constraints
184 c_int nb_leave; ///< number of leaving constraints
185 solver_dense *At_scale; ///< running vector of sqrt(sigma), used to scale At_sqrt_sigma
186 solver_sparse *At_sqrt_sigma; ///< A' * sqrt(sigma)
188
189
190/**
191 * QPALM Workspace
192 *
193 * The workspace is the main data structure and is given as a pointer to (almost) all QPALM functions.
194 * It contains pointers to the settings, the data, return info, solution variables and intermediate
195 * workspace variables.
196 */
197typedef struct {
198 QPALMData *data; ///< problem data to work on (possibly scaled)
199
200 /**
201 * @name Iterates
202 * @{
203 */
204 c_float *x; ///< primal iterate
205 c_float *y; ///< dual iterate
206 c_float *Ax; ///< scaled A * x
207 c_float *Qx; ///< scaled Q * x
208 c_float *Aty; ///< A' * y (useful for saving one mat_tpose_vec)
209 c_float *x_prev; ///< previous primal iterate
210 c_int initialized; ///< flag whether the iterates were initialized or not
211 /** @} */
212
213 /**
214 * @name Workspace variables
215 * @{
216 */
217 c_float *temp_m; ///< placeholder for vector of size m
218 c_float *temp_n; ///< placeholder for vector of size n
219 c_float *sigma; ///< penalty vector
220 c_float *sigma_inv; ///< 1./sigma
221 c_float sqrt_sigma_max; ///< sqrt(sigma_max)
222 c_int nb_sigma_changed; ///< number of sigma-components that changed in an outer iteration (relevant for factorization update)
223 c_float gamma; ///< proximal penalty factor
224 c_int gamma_maxed; ///< flag to indicate whether gamma has been maximized when the primal residual was low
225 c_float *Axys; ///< Ax + y./sigma
226 c_float *z; ///< projection of Axys onto the constraint set [bmin, bmax]
227 c_float *pri_res; ///< primal residual
228 c_float *pri_res_in; ///< intermediate primal residual
229 c_float *yh; ///< candidate dual update
230 c_float *Atyh; ///< A' * yh
231 c_float *df; ///< gradient of the primal objective (+proximal term)
232 c_float *x0; ///< record of the primal iterate during the last dual update
233 c_float *xx0; ///< x - x0
234 c_float *dphi; ///< gradient of the Lagrangian
235 c_float *neg_dphi; ///< -dphi, required as the rhs in SCHUR
236 c_float *dphi_prev; ///< previous gradient of the Lagrangian
237 c_float *d; ///< primal update step
238
239 /** @} */
240
241 /**
242 * @name Linesearch variables
243 * @{
244 */
245 c_float tau; ///< stepsize
246 c_float *Qd; ///< Q * d
247 c_float *Ad; ///< A * d
248 c_float *sqrt_sigma; ///< elementwise sqrt(sigma)
249 c_float sqrt_delta; ///< sqrt(penalty update factor)
250 c_float eta; ///< linesearch parameter
251 c_float beta; ///< linesearch parameter
252 c_float *delta; ///< linesearch parameter
253 c_float *alpha; ///< linesearch parameter
254 c_float *temp_2m; ///< placeholder for vector of size 2m
255 c_float *delta2; ///< delta .* delta
256 c_float *delta_alpha; ///< delta .* alpha
257 array_element *s; ///< alpha ./ delta
258 c_int *index_L; ///< index set L (where s>0)
259 c_int *index_P; ///< index set P (where delta>0)
260 c_int *index_J; ///< index set J (L xor P)
261
262 /** @} */
263
264 /**
265 * @name Termination criteria variables
266 * @{
267 */
268 c_float eps_pri; ///< primal tolerance
269 c_float eps_dua; ///< dual tolerance
270 c_float eps_dua_in; ///< intermediate dual tolerance
271 c_float eps_abs_in; ///< intermediate absolute tolerance
272 c_float eps_rel_in; ///< intermediate relative tolerance
273 /** @} */
274
275 /**
276 * @name Primal infeasibility variables
277 * @{
278 */
279 c_float *delta_y; ///< difference of consecutive dual iterates
280 c_float *Atdelta_y; ///< A' * delta_y
281
282 /** @} */
283
284 /**
285 * @name Dual infeasibility variables
286 * @{
287 */
288 c_float *delta_x; ///< difference of consecutive primal iterates
289 c_float *Qdelta_x; ///< Q * delta_x
290 c_float *Adelta_x; ///< A * delta_x
291
292 /** @} */
293
294 /**
295 * @name Temporary vectors used in scaling
296 * @{
297 */
298
299 c_float *D_temp; ///< temporary primal variable scaling vectors
300 c_float *E_temp; ///< temporary constraints scaling vectors
301
302
303 /** @} */
304 QPALMSolver *solver; ///< linsys variables
305 QPALMSettings *settings; ///< problem settings
306 QPALMScaling *scaling; ///< scaling vectors
307 QPALMSolution *solution; ///< problem solution
308 QPALMInfo *info; ///< solver information
309
310 # ifdef PROFILING
311 QPALMTimer *timer; ///< timer object
312 # endif // ifdef PROFILING
313
315
316
317# ifdef __cplusplus
318}
319# endif
320
321#endif // ifndef QPALM_TYPES_H
Custom memory allocation, print and utility functions, and data types for floats and ints.
ladel_int c_int
type for integer numbers
Definition: global_opts.h:42
ladel_double c_float
type for floating point numbers
Definition: global_opts.h:41
::QPALMInfo QPALMInfo
Definition: qpalm.cpp:23
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 dua2_res_norm
norm of intermediate dual residual (minus proximal term)
Definition: types.h:82
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
c_float dua_res_norm
norm of dual residual
Definition: types.h:81
Problem scaling matrices stored as vectors.
Definition: types.h:61
c_float * Dinv
primal variable rescaling
Definition: types.h:63
c_float c
objective scaling
Definition: types.h:66
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
c_int nb_active_constraints
number of active constraints
Definition: types.h:180
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 eta
linesearch parameter
Definition: types.h:250
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 tau
stepsize
Definition: types.h:245
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 eps_dua
dual tolerance
Definition: types.h:269
c_float beta
linesearch parameter
Definition: types.h:251
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 sqrt_sigma_max
sqrt(sigma_max)
Definition: types.h:221
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 eps_dua_in
intermediate dual tolerance
Definition: types.h:270
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
size_t i
index
Definition: types.h:29
c_float x
value of the element
Definition: types.h:28
ladel_work solver_common
Definition: types.h:18
ladel_double solver_dense
Definition: types.h:20
ladel_symbolics solver_symbolics
Definition: types.h:22
ladel_sparse_matrix solver_sparse
Definition: types.h:19
struct QPALM_TIMER QPALMTimer
QPALM Timer for statistics.
Definition: types.h:53
ladel_factor solver_factor
Definition: types.h:21
struct array_element array_element
Array to sort in linesearch.