Scippy

SCIP

Solving Constraint Integer Programs

presol_qpkktref.c
Go to the documentation of this file.
1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2 /* */
3 /* This file is part of the program and library */
4 /* SCIP --- Solving Constraint Integer Programs */
5 /* */
6 /* Copyright (C) 2002-2017 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not email to scip@zib.de. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file presol_qpkktref.c
17  * @brief qpkktref presolver
18  * @author Tobias Fischer
19  *
20  * This presolver tries to add the KKT conditions as additional (redundant) constraints to the (mixed-binary) quadratic
21  * program
22  * \f[
23  * \begin{array}{ll}
24  * \min & x^T Q x + c^T x + d \\
25  * & A x \leq b, \\
26  * & x \in \{0, 1\}^{p} \times R^{n-p}.
27  * \end{array}
28  * \f]
29  *
30  * We first check if the structure of the program is like (QP), see the documentation of the function
31  * checkConsQuadraticProblem().
32  *
33  * If the problem is known to be bounded (all variables have finite lower and upper bounds), then we add the KKT
34  * conditions. For a continuous QPs the KKT conditions have the form
35  * \f[
36  * \begin{array}{ll}
37  * Q x + c + A^T \mu = 0,\\
38  * Ax \leq b,\\
39  * \mu_i \cdot (Ax - b)_i = 0, & i \in \{1, \dots, m\},\\
40  * \mu \geq 0.
41  * \end{array}
42  * \f]
43  * where \f$\mu\f$ are the Lagrangian variables. Each of the complementarity constraints \f$\mu_i \cdot (Ax - b)_i = 0\f$
44  * is enforced via an SOS1 constraint for \f$\mu_i\f$ and an additional slack variable \f$s_i = (Ax - b)_i\f$.
45  *
46  * For mixed-binary QPs, the KKT-like conditions are
47  * \f[
48  * \begin{array}{ll}
49  * Q x + c + A^T \mu + I_J \lambda = 0,\\
50  * Ax \leq b,\\
51  * x_j \in \{0,1\} & j \in J,\\
52  * (1 - x_j) \cdot z_j = 0 & j \in J,\\
53  * x_j \cdot (z_j - \lambda_j) = 0 & j \in J,\\
54  * \mu_i \cdot (Ax - b)_i = 0 & i \in \{1, \dots, m\},\\
55  * \mu \geq 0,
56  * \end{array}
57  * \f]
58  * where \f$J = \{1,\dots, p\}\f$, \f$\mu\f$ and \f$\lambda\f$ are the Lagrangian variables, and \f$I_J\f$ is the
59  * submatrix of the \f$n\times n\f$ identity matrix with columns indexed by \f$J\f$. For the derivation of the KKT-like
60  * conditions, see
61  *
62  * Branch-And-Cut for Complementarity and Cardinality Constrained Linear Programs,@n
63  * Tobias Fischer, PhD Thesis (2016)
64  *
65  * Algorithmically:
66  *
67  * - we handle the quadratic term variables of the quadratic constraint like in the method
68  * presolveAddKKTQuadQuadraticTerms()
69  * - we handle the bilinear term variables of the quadratic constraint like in the method presolveAddKKTQuadBilinearTerms()
70  * - we handle the linear term variables of the quadratic constraint like in the method presolveAddKKTQuadLinearTerms()
71  * - we handle linear constraints in the method presolveAddKKTLinearConss()
72  * - we handle aggregated variables in the method presolveAddKKTAggregatedVars()
73  *
74  * we have a hashmap from each variable to the index of the dual constraint in the KKT conditions.
75  */
76 
77 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
78 
79 #include <stdio.h>
80 #include <assert.h>
81 #include <string.h>
82 
83 #include "scip/presol_qpkktref.h"
84 #include "scip/cons_quadratic.h"
85 #include "scip/cons_linear.h"
86 #include "scip/cons_knapsack.h"
87 #include "scip/cons_logicor.h"
88 #include "scip/cons_setppc.h"
89 #include "scip/cons_varbound.h"
90 #include "scip/cons_sos1.h"
91 
92 
93 #define PRESOL_NAME "qpkktref"
94 #define PRESOL_DESC "adds KKT conditions to (mixed-binary) quadratic programs"
95 #define PRESOL_PRIORITY -1 /**< priority of the presolver (>= 0: before, < 0: after constraint handlers);
96  * combined with propagators */
97 #define PRESOL_MAXROUNDS -1 /**< maximal number of presolving rounds the presolver participates in (-1: no
98  * limit) */
99 #define PRESOL_TIMING SCIP_PRESOLTIMING_MEDIUM /* timing of the presolver (fast, medium, or exhaustive) */
100 
102 /*
103  * Data structures
104  */
105 
106 /** presolver data */
107 struct SCIP_PresolData
108 {
109  SCIP_Bool addkktbinary; /**< if TRUE then allow binary variables for KKT update */
110  SCIP_Bool updatequadbounded; /**< if TRUE then only apply the update to QPs with bounded variables; if
111  * the variables are not bounded then a finite optimal solution might not
112  * exist and the KKT conditions would then be invalid */
113  SCIP_Bool updatequadindef; /**< if TRUE then apply quadratic constraint update even if the quadratic
114  * constraint matrix is known to be indefinite */
115 };
116 
117 
118 /*
119  * Local methods
120  */
121 
122 /** for a linear constraint \f$a^T x \leq b\f$, create the complementarity constraint \f$\mu \cdot s = 0\f$, where
123  * \f$s = b - a^T x\f$ and \f$\mu\f$ is the dual variable associated to the constraint \f$a^T x \leq b\f$
124  */
125 static
127  SCIP* scip, /**< SCIP pointer */
128  const char* namepart, /**< name of linear constraint */
129  SCIP_VAR** vars, /**< variables of linear constraint */
130  SCIP_Real* vals, /**< coefficients of variables in linear constraint */
131  SCIP_Real lhs, /**< left hand side of linear constraint */
132  SCIP_Real rhs, /**< right hand side of linear constraint */
133  int nvars, /**< number of variables of linear constraint */
134  SCIP_VAR* dualvar, /**< dual variable associated to linear constraint */
135  SCIP_Bool takelhs, /**< whether to consider the lhs or the rhs of the constraint */
136  int* naddconss /**< buffer to increase with number of created additional constraints */
137  )
138 {
139  char name[SCIP_MAXSTRLEN];
140  SCIP_CONS* KKTlincons;
141  SCIP_CONS* sos1cons;
142  SCIP_VAR* slack;
143  SCIP_Real slackcoef;
144  SCIP_Real eqval;
145 
146  assert( scip != NULL );
147  assert( namepart != NULL );
148  assert( vars != NULL );
149  assert( vals != NULL );
150  assert( dualvar != NULL );
151  assert( ! takelhs || ! SCIPisInfinity(scip, -lhs) );
152  assert( takelhs || ! SCIPisInfinity(scip, rhs) );
153  assert( naddconss != NULL );
154 
155  if( takelhs )
156  {
157  eqval = lhs;
158  slackcoef = -1.0;
159  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "slack_lhs_%s", namepart);
160  }
161  else
162  {
163  eqval = rhs;
164  slackcoef = 1.0;
165  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "slack_rhs_%s", namepart);
166  }
167 
168  /* create slack variable */
169  SCIP_CALL( SCIPcreateVarBasic(scip, &slack, name, 0.0, SCIPinfinity(scip), 0.0, SCIP_VARTYPE_CONTINUOUS) );
170 
171  /* add skack variable */
172  SCIP_CALL( SCIPaddVar(scip, slack) );
173 
174  /* create a new linear constraint */
175  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKTlin_%s_%d", namepart, takelhs);
176  SCIP_CALL( SCIPcreateConsBasicLinear(scip, &KKTlincons, name, nvars, vars, vals, eqval, eqval) );
177 
178  /* add slack variable to linear constraint */
179  SCIP_CALL( SCIPaddCoefLinear(scip, KKTlincons, slack, slackcoef) );
180 
181  /* create SOS1 (complementarity) constraint involving dual variable of linear constraint and slack variable */
182  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKTsos1_lin_%s_%d", namepart, takelhs);
183  SCIP_CALL( SCIPcreateConsBasicSOS1(scip, &sos1cons, name, 0, NULL, NULL) );
184 
185  /* add slack and dual variable to SOS1 constraint */
186  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons, slack, 1.0) );
187  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons, dualvar, 2.0) );
188 
189  /* add/release constraints */
190  SCIP_CALL( SCIPaddCons(scip, sos1cons) );
191  SCIP_CALL( SCIPaddCons(scip, KKTlincons) );
192  SCIP_CALL( SCIPreleaseCons(scip, &sos1cons) );
193  SCIP_CALL( SCIPreleaseCons(scip, &KKTlincons) );
194  *naddconss = *naddconss + 2;
195 
196  /* release slack variable */
197  SCIP_CALL( SCIPreleaseVar(scip, &slack) );
198 
199  return SCIP_OKAY;
200 }
201 
202 /** create complementarity constraints of KKT conditions associated to bounds of variables
203  * - for an upper bound constraint \f$x_i \leq u_i\f$, create the complementarity constraint \f$\mu_i \cdot s_i = 0\f$,
204  * where \f$s_i = u_i - x_i\f$ and \f$\mu_i\f$ is the dual variable of the upper bound constraint
205  * - for a lower bound constraint \f$x_i \geq l_i\f$, create the complementarity constraint \f$\lambda_i \cdot w_i = 0\f$,
206  * where \f$w_i = x_i - l_i\f$
207  * and \f$\lambda_i\f$ is the dual variable of the lower bound constraint
208  */
209 static
211  SCIP* scip, /**< SCIP pointer */
212  SCIP_VAR* var, /**< variable */
213  SCIP_VAR* dualvar, /**< dual variable associated to bound of variable */
214  SCIP_Bool takelb, /**< whether to consider the lower or upper bound of variable */
215  int* naddconss /**< buffer to increase with number of created additional constraints */
216  )
217 {
218  char name[SCIP_MAXSTRLEN];
219  SCIP_CONS* KKTlincons;
220  SCIP_CONS* sos1cons;
221  SCIP_VAR* slack;
222  SCIP_Real slackcoef;
223  SCIP_Real eqval;
224 
225  assert( scip != NULL );
226  assert( var != NULL );
227  assert( dualvar != NULL );
228  assert( ! takelb || ! SCIPisInfinity(scip, -SCIPvarGetLbGlobal(var)) );
229  assert( takelb || ! SCIPisInfinity(scip, SCIPvarGetUbGlobal(var)) );
230  assert( naddconss != NULL );
231 
232  if( takelb )
233  {
234  eqval = SCIPvarGetLbGlobal(var);
235  slackcoef = -1.0;
236  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "slack_lb_%s", SCIPvarGetName(var));
237  }
238  else
239  {
240  eqval = SCIPvarGetUbGlobal(var);
241  slackcoef = 1.0;
242  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "slack_ub_%s", SCIPvarGetName(var));
243  }
244 
245  /* create complementarity constraint; if bound is nonzero, we additionally need to introduce a slack variable */
246  if( SCIPisFeasZero(scip, eqval) && SCIPvarGetStatus(var) != SCIP_VARSTATUS_MULTAGGR )
247  {
248  /* create SOS1 (complementarity) constraint involving dual variable of linear constraint and slack variable */
249  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKTsos1_bound%s_%d", SCIPvarGetName(var), takelb);
250  SCIP_CALL( SCIPcreateConsBasicSOS1(scip, &sos1cons, name, 0, NULL, NULL) );
251 
252  /* add slack and dual variable to SOS1 constraint */
253  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons, var, 1.0) );
254  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons, dualvar, 2.0) );
255 
256  /* add/release constraint */
257  SCIP_CALL( SCIPaddCons(scip, sos1cons) );
258  SCIP_CALL( SCIPreleaseCons(scip, &sos1cons) );
259  ++(*naddconss);
260  }
261  else
262  {
263  /* create slack variable */
264  SCIP_CALL( SCIPcreateVarBasic(scip, &slack, name, 0.0, SCIPinfinity(scip), 0.0, SCIP_VARTYPE_CONTINUOUS) );
265 
266  /* add skack variable */
267  SCIP_CALL( SCIPaddVar(scip, slack) );
268 
269  /* create a new linear constraint */
270  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKT_bound%s_%d", SCIPvarGetName(var), takelb);
271  SCIP_CALL( SCIPcreateConsBasicLinear(scip, &KKTlincons, name, 0, NULL, NULL, eqval, eqval) );
272 
273  /* add slack variable to linear constraint */
274  SCIP_CALL( SCIPaddCoefLinear(scip, KKTlincons, var, 1.0) );
275  SCIP_CALL( SCIPaddCoefLinear(scip, KKTlincons, slack, slackcoef) );
276 
277  /* create SOS1 (complementarity) constraint involving dual variable of linear constraint and slack variable */
278  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKTsos1_bound%s_%d", SCIPvarGetName(var), takelb);
279  SCIP_CALL( SCIPcreateConsBasicSOS1(scip, &sos1cons, name, 0, NULL, NULL) );
280 
281  /* add slack and dual variable to SOS1 constraint */
282  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons, slack, 1.0) );
283  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons, dualvar, 2.0) );
284 
285  /* add/release constraints */
286  SCIP_CALL( SCIPaddCons(scip, sos1cons) );
287  SCIP_CALL( SCIPaddCons(scip, KKTlincons) );
288  SCIP_CALL( SCIPreleaseCons(scip, &sos1cons) );
289  SCIP_CALL( SCIPreleaseCons(scip, &KKTlincons) );
290  *naddconss = *naddconss + 2;
291 
292  /* release slack variable */
293  SCIP_CALL( SCIPreleaseVar(scip, &slack) );
294  }
295 
296  return SCIP_OKAY;
297 }
298 
299 /** create the complementarity constraints of the KKT-like conditions associated to a binary variable \f$x_i\f$;
300  * these are \f$(1 - x_i) \cdot z_i = 0\f$ and \f$x_i \cdot (z_i - \lambda_i) = 0\f$, where \f$z_i\f$ and
301  * \f$\lambda_i\f$ are dual variables
302  */
303 static
305  SCIP* scip, /**< SCIP pointer */
306  SCIP_VAR* var, /**< variable */
307  SCIP_VAR* dualbin1, /**< first dual variable associated to binary variable */
308  SCIP_VAR* dualbin2, /**< second dual variable associated to binary variable */
309  int* naddconss /**< buffer to increase with number of created additional constraints */
310  )
311 {
312  char name[SCIP_MAXSTRLEN];
313  SCIP_CONS* conslinbin1;
314  SCIP_CONS* conslinbin2;
315  SCIP_CONS* sos1cons1;
316  SCIP_CONS* sos1cons2;
317  SCIP_VAR* slackbin1;
318  SCIP_VAR* slackbin2;
319 
320  assert( scip != NULL );
321  assert( var != NULL );
322  assert( dualbin1 != NULL );
323  assert( dualbin2 != NULL );
324  assert( naddconss != NULL );
325 
326  /* create first slack variable associated to binary constraint; domain [-inf, inf] */
327  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "dual_%s_slackbin1", SCIPvarGetName(var));
328  SCIP_CALL( SCIPcreateVarBasic(scip, &slackbin1, name, -SCIPinfinity(scip), SCIPinfinity(scip), 0.0,
330  SCIP_CALL( SCIPaddVar(scip, slackbin1) );
331  assert( slackbin1 != NULL );
332 
333  /* create a new linear constraint: dualbin1 - dualbin2 = slackbin */
334  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKTBinary1_%s", SCIPvarGetName(var));
335  SCIP_CALL( SCIPcreateConsBasicLinear(scip, &conslinbin1, name, 0, NULL, NULL, 0.0, 0.0) );
336  SCIP_CALL( SCIPaddCoefLinear(scip, conslinbin1, dualbin1, 1.0) );
337  SCIP_CALL( SCIPaddCoefLinear(scip, conslinbin1, dualbin2, -1.0) );
338  SCIP_CALL( SCIPaddCoefLinear(scip, conslinbin1, slackbin1, -1.0) );
339  SCIP_CALL( SCIPaddCons(scip, conslinbin1) );
340  SCIP_CALL( SCIPreleaseCons(scip, &conslinbin1) );
341  ++(*naddconss);
342 
343  /* create SOS1 (complementarity) constraint involving binary variable and slack variable */
344  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKTsos1_bin1%s", SCIPvarGetName(var));
345  SCIP_CALL( SCIPcreateConsBasicSOS1(scip, &sos1cons1, name, 0, NULL, NULL) );
346 
347  /* add slack and dual variable to SOS1 constraint */
348  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons1, var, 1.0) );
349  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons1, slackbin1, 2.0) );
350 
351  /* add/release constraint */
352  SCIP_CALL( SCIPaddCons(scip, sos1cons1) );
353  SCIP_CALL( SCIPreleaseCons(scip, &sos1cons1) );
354  ++(*naddconss);
355 
356  /* release slack variable */
357  SCIP_CALL( SCIPreleaseVar(scip, &slackbin1) );
358 
359 
360  /* create second slack variable associated to binary constraint; domain [0, inf] */
361  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "dual_%s_slackbin2", SCIPvarGetName(var));
362  SCIP_CALL( SCIPcreateVarBasic(scip, &slackbin2, name, 0.0, SCIPinfinity(scip), 0.0, SCIP_VARTYPE_CONTINUOUS) );
363  SCIP_CALL( SCIPaddVar(scip, slackbin2) );
364  assert( slackbin2 != NULL );
365 
366  /* create a new linear constraint: 1.0 - var = slackbin2 */
367  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKTBinary2_%s", SCIPvarGetName(var));
368  SCIP_CALL( SCIPcreateConsBasicLinear(scip, &conslinbin2, name, 0, NULL, NULL, 1.0, 1.0) );
369  SCIP_CALL( SCIPaddCoefLinear(scip, conslinbin2, var, 1.0) );
370  SCIP_CALL( SCIPaddCoefLinear(scip, conslinbin2, slackbin2, 1.0) );
371  SCIP_CALL( SCIPaddCons(scip, conslinbin2) );
372  SCIP_CALL( SCIPreleaseCons(scip, &conslinbin2) );
373  ++(*naddconss);
374 
375  /* create SOS1 (complementarity) constraint involving first dual variable and slack variable */
376  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKTsos1_bin2%s", SCIPvarGetName(var));
377  SCIP_CALL( SCIPcreateConsBasicSOS1(scip, &sos1cons2, name, 0, NULL, NULL) );
378 
379  /* add slack and dual variable to SOS1 constraint */
380  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons2, dualbin1, 1.0) );
381  SCIP_CALL( SCIPaddVarSOS1(scip, sos1cons2, slackbin2, 2.0) );
382 
383  /* add/release constraint */
384  SCIP_CALL( SCIPaddCons(scip, sos1cons2) );
385  SCIP_CALL( SCIPreleaseCons(scip, &sos1cons2) );
386  ++(*naddconss);
387 
388  /* release slack variable */
389  SCIP_CALL( SCIPreleaseVar(scip, &slackbin2) );
390 
391  return SCIP_OKAY;
392 }
393 
394 /** create/get dual constraint of KKT conditions associated to primal variable @n@n
395  * if variable does not already exist in hashmap then
396  * 1. create dual constraint for variable
397  * 2. create a dual variable \f$\mu_i\f$ for the upper bound constraint \f$x_i \leq u_i\f$
398  * 3. create a dual variable \f$\lambda_i\f$ for the lower bound constraint \f$x_i \geq l_i\f$
399  * 4. create the complementarity constraint \f$\mu_i \cdot s_i = 0\f$, where \f$s_i = u_i - x_i\f$
400  * 5. create the complementarity constraint \f$\lambda_i \cdot w_i = 0\f$, where \f$w_i = x_i - l_i\f$
401  * 6. add objective coefficients of dual variables
402  * 7. the treatment of binary variables needs special care see the documentation of createKKTComplementarityBinary()
403  *
404  * if variable exists in hasmap then the dual constraint associated to the variable has already been created and is returned
405  */
406 static
408  SCIP* scip, /**< SCIP pointer */
409  SCIP_CONS* objcons, /**< objective constraint */
410  SCIP_VAR* var, /**< variable */
411  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
412  SCIP_CONS** dualconss, /**< array with dual constraints */
413  int* ndualconss, /**< pointer to store number of dual constraints */
414  SCIP_CONS** dualcons, /**< dual constraint associated to variable */
415  int* naddconss /**< buffer to increase with number of created additional constraints */
416  )
417 {
418  SCIP_VAR* dualub = NULL; /* dual variable associated to upper bound constraint */
419  SCIP_VAR* duallb = NULL; /* dual variable associated to lower bound constraint */
420  SCIP_VAR* dualbin1 = NULL; /* first dual variable associated to binary variable */
421  SCIP_VAR* dualbin2 = NULL; /* second dual variable associated to binary variable */
422 
423  assert( scip != NULL );
424  assert( objcons != NULL );
425  assert( var != NULL );
426  assert( varhash != NULL );
427  assert( dualconss != NULL );
428  assert( ndualconss != NULL );
429  assert( naddconss != NULL );
430 
431  /* if variable exists in hashmap */
432  if( SCIPhashmapExists(varhash, var) )
433  {
434  int ind;
435  ind = (int) (size_t) SCIPhashmapGetImage(varhash, var);
436  *dualcons = dualconss[ind];
437  }
438  else
439  {
440  char name[SCIP_MAXSTRLEN];
441  SCIP_Real lb;
442  SCIP_Real ub;
443 
444  lb = SCIPvarGetLbGlobal(var);
445  ub = SCIPvarGetUbGlobal(var);
446 
447  /* create dual variables corresponding to the bounds of the variables; binary variables have to be treated in a
448  * different way */
449  if( SCIPvarIsBinary(var) )
450  {
451  /* create first dual variable associated to binary constraint; the domain of dualbin is [-inf,inf]; the objective
452  * coefficient is -0.5 */
453  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "dual_%s_bin1", SCIPvarGetName(var));
454  SCIP_CALL( SCIPcreateVarBasic(scip, &dualbin1, name, -SCIPinfinity(scip), SCIPinfinity(scip), 0.0,
456  SCIP_CALL( SCIPaddVar(scip, dualbin1) );
457  assert( dualbin1 != NULL );
458  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, dualbin1, -0.5) );
459 
460  /* create second variable associated to binary constraint; the domain of dualbin2 is [-inf,inf]; the objective
461  * coefficient is zero */
462  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "dual_%s_bin2", SCIPvarGetName(var));
463  SCIP_CALL( SCIPcreateVarBasic(scip, &dualbin2, name, -SCIPinfinity(scip), SCIPinfinity(scip), 0.0,
465  SCIP_CALL( SCIPaddVar(scip, dualbin2) );
466  assert( dualbin2 != NULL );
467  }
468  else
469  {
470  if( ! SCIPisInfinity(scip, -lb) )
471  {
472  /* create dual variable associated to lower bound; the domain of duallb is [0,inf] */
473  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "dual_%s_lb", SCIPvarGetName(var));
474  SCIP_CALL( SCIPcreateVarBasic(scip, &duallb, name, 0.0, SCIPinfinity(scip), 0.0, SCIP_VARTYPE_CONTINUOUS) );
475  SCIP_CALL( SCIPaddVar(scip, duallb) );
476  assert( duallb != NULL );
477  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, duallb, 0.5 * lb) );
478  }
479 
480  if( ! SCIPisInfinity(scip, ub) )
481  {
482  /* create dual variable associated to upper bound; the domain of dualub is [0,inf] */
483  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "dual_%s_ub", SCIPvarGetName(var));
484  SCIP_CALL( SCIPcreateVarBasic(scip, &dualub, name, 0.0, SCIPinfinity(scip), 0.0, SCIP_VARTYPE_CONTINUOUS) );
485  SCIP_CALL( SCIPaddVar(scip, dualub) );
486  assert( dualub != NULL );
487  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, dualub, -0.5 * ub) );
488  }
489  }
490 
491  /* add variable in map */
492  SCIP_CALL( SCIPhashmapInsert(varhash, var, (void*) (size_t) *ndualconss) );/*lint !e571*/
493  assert( *ndualconss == (int) (size_t) SCIPhashmapGetImage(varhash, var) );
494 
495  /* create a new linear constraint */
496  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "KKTref_%s", SCIPvarGetName(var));
497  SCIP_CALL( SCIPcreateConsBasicLinear(scip, dualcons, name, 0, NULL, NULL, 0.0, 0.0) );
498 
499  /* add dual constraint to array for later use */
500  dualconss[(*ndualconss)++] = *dualcons;
501 
502 
503  /* add dual variables to dual constraints and create complementarity constraints; binary variables have to be
504  * treated in a different way */
505  if( SCIPvarIsBinary(var) )
506  {
507  /* add coefficient of second dual variable corresponding to binary variable */
508  SCIP_CALL( SCIPaddCoefLinear(scip, *dualcons, dualbin2, 1.0) );
509 
510  /* create complementarity constraints */
511  SCIP_CALL( createKKTComplementarityBinary(scip, var, dualbin1, dualbin2, naddconss) );
512 
513  SCIP_CALL( SCIPreleaseVar(scip, &dualbin1) );
514  SCIP_CALL( SCIPreleaseVar(scip, &dualbin2) );
515  }
516  else
517  {
518  if( duallb != NULL )
519  {
520  /* add dual variable corresponding to lower bound of variable */
521  SCIP_CALL( SCIPaddCoefLinear(scip, *dualcons, duallb, -1.0) );
522 
523  /* create complementarity constraint between slack variable of lower bound constraint and dual variable of
524  * lower bound */
525  SCIP_CALL( createKKTComplementarityBounds(scip, var, duallb, TRUE, naddconss) );
526 
527  SCIP_CALL( SCIPreleaseVar(scip, &duallb) );
528  }
529 
530  if( dualub != NULL )
531  {
532  /* add dual variable corresponding to upper bound of variable */
533  SCIP_CALL( SCIPaddCoefLinear(scip, *dualcons, dualub, 1.0) );
534 
535  /* create complementarity constraint between slack variable of upper bound constraint and dual variable of
536  * upper bound */
537  SCIP_CALL( createKKTComplementarityBounds(scip, var, dualub, FALSE, naddconss) );
538 
539  SCIP_CALL( SCIPreleaseVar(scip, &dualub) );
540  }
541  }
542  }
543  assert( *dualcons != NULL );
544 
545  return SCIP_OKAY;
546 }
547 
548 /** handle (a single) linear constraint for quadratic constraint update
549  * 1. create the dual constraints (i.e., the two rows of \f$Q x + c + A^T \mu = 0\f$) associated to the variables of the
550  * linear constraint, if not done already
551  * 2. create the dual variables and the complementarity constraints for the lower and upper bound constraints of the
552  * variables of the linear constraint, if not done already
553  * 3. create the dual variable \f$\mu_i\f$ associated to this linear constraint
554  * 4. create the complementarity constraint \f$\mu_i \cdot (Ax - b)_i = 0\f$ associated to this linear constraint
555  * 5. add objective coefficients of dual variables
556  *
557  * for steps 1 and 2 see the documentation of createKKTDualCons() for further information.@n
558  * for step 4 see the documentation of the function createKKTComplementarityLinear() for further information.
559  */
560 static
562  SCIP* scip, /**< SCIP pointer */
563  SCIP_CONS* objcons, /**< objective constraint */
564  const char* namepart, /**< name of linear constraint */
565  SCIP_VAR** vars, /**< variables of linear constraint */
566  SCIP_Real* vals, /**< coefficients of variables in linear constraint */
567  SCIP_Real lhs, /**< left hand side of linear constraint */
568  SCIP_Real rhs, /**< right hand side of linear constraint */
569  int nvars, /**< number of variables of linear constraint */
570  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
571  SCIP_CONS** dualconss, /**< array with dual constraints */
572  int* ndualconss, /**< pointer to store number of dual constraints */
573  int* naddconss /**< buffer to increase with number of created additional constraints */
574  )
575 {
576  int i;
577 
578  assert( scip != NULL );
579  assert( objcons != NULL );
580  assert( namepart != NULL );
581  assert( varhash != NULL );
582  assert( dualconss != NULL );
583  assert( ndualconss != NULL );
584  assert( vars != NULL );
585  assert( vals != NULL );
586  assert( namepart != NULL );
587  assert( naddconss != NULL );
588 
589  /* differ between left hand side and right hand side case (i=0 -> lhs; i=1 -> rhs) */
590  for( i = 0; i < 2; ++i )
591  {
592  char name[SCIP_MAXSTRLEN];
593  SCIP_VAR* duallin = NULL;
594  int j;
595 
596  /* skip one iteration if lhs equals rhs */
597  if( i == 0 && SCIPisFeasEQ(scip, lhs, rhs) )
598  continue;
599 
600  /* create dual variable corresponding to linear constraint */
601  if( i == 0 )
602  {
603  assert( ! SCIPisFeasEQ(scip, lhs, rhs) );
604 
605  if( SCIPisInfinity(scip, -lhs) )
606  continue;
607 
608  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "dual_%s_lhs", namepart);
609  SCIP_CALL( SCIPcreateVarBasic(scip, &duallin, name, 0.0, SCIPinfinity(scip), 0.0, SCIP_VARTYPE_CONTINUOUS) );
610  SCIP_CALL( SCIPaddVar(scip, duallin) );
611  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, duallin, 0.5 * lhs) );
612 
613  /* create complementarity constraint between dual variable and slack variable of linear constraint */
614  SCIP_CALL( createKKTComplementarityLinear(scip, namepart, vars, vals, lhs, rhs, nvars, duallin, TRUE,
615  naddconss) );
616  }
617  else
618  {
619  if( SCIPisInfinity(scip, rhs) )
620  continue;
621 
622  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "dual_%s_rhs", namepart);
623  if( SCIPisFeasEQ(scip, lhs, rhs) )
624  {
625  SCIP_CALL( SCIPcreateVarBasic(scip, &duallin, name, -SCIPinfinity(scip), SCIPinfinity(scip), 0.0,
627  SCIP_CALL( SCIPaddVar(scip, duallin) );
628  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, duallin, -0.5 * rhs) );
629  }
630  else
631  {
632  SCIP_CALL( SCIPcreateVarBasic(scip, &duallin, name, 0.0, SCIPinfinity(scip), 0.0, SCIP_VARTYPE_CONTINUOUS) );
633  SCIP_CALL( SCIPaddVar(scip, duallin) );
634  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, duallin, -0.5 * rhs) );
635 
636  /* create complementarity constraint between dual variable and slack variable of linear constraint */
637  SCIP_CALL( createKKTComplementarityLinear(scip, namepart, vars, vals, lhs, rhs, nvars, duallin, FALSE,
638  naddconss) );
639  }
640  }
641  assert( duallin != NULL );
642 
643 
644  /* loop through variables of linear constraint */
645  for( j = 0; j < nvars; ++j )
646  {
647  SCIP_CONS* dualcons = NULL; /* dual constraint associated to variable */
648  SCIP_VAR* var;
649 
650  var = vars[j];
651 
652  /* create/get dual constraint associated to variable;
653  * if variable does not already exist in hashmap then create dual variables for its bounds */
654  SCIP_CALL( createKKTDualCons(scip, objcons, var, varhash, dualconss, ndualconss, &dualcons, naddconss) );
655  assert( dualcons != NULL );
656 
657  /* add dual variable corresponding to linear constraint */
658  if( i == 0 )
659  {
660  SCIP_CALL( SCIPaddCoefLinear(scip, dualcons, duallin, -vals[j]) );
661  }
662  else
663  {
664  SCIP_CALL( SCIPaddCoefLinear(scip, dualcons, duallin, vals[j]) );
665  }
666  }
667 
668  /* release dual variable */
669  SCIP_CALL( SCIPreleaseVar(scip, &duallin) );
670  }
671 
672  return SCIP_OKAY;
673 }
674 
675 /** handle linear constraints for quadratic constraint update, see the documentation of the function
676  * presolveAddKKTLinearCons() for an explanation
677  */
678 static
680  SCIP* scip, /**< SCIP pointer */
681  SCIP_CONS* objcons, /**< objective constraint */
682  SCIP_CONS** savelinconss, /**< copy of array with linear constraints */
683  int nlinconss, /**< number of linear constraints */
684  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
685  SCIP_CONS** dualconss, /**< array with dual constraints */
686  int* ndualconss, /**< pointer to store number of dual constraints */
687  int* naddconss, /**< buffer to increase with number of created additional constraints */
688  int* ndelconss /**< buffer to increase with number of deleted constraints */
689  )
690 {
691  int c;
692 
693  assert( scip != NULL );
694  assert( objcons != NULL );
695  assert( varhash != NULL );
696  assert( dualconss != NULL );
697  assert( ndualconss != NULL );
698  assert( naddconss != NULL );
699  assert( ndelconss != NULL );
700 
701  /* loop through linear constraints */
702  for( c = 0; c < nlinconss; ++c )
703  {
704  SCIP_CONS* lincons;
705  SCIP_VAR** vars;
706  SCIP_Real* vals;
707  SCIP_Real lhs;
708  SCIP_Real rhs;
709  int nvars;
710 
711  /* get data of constraint */
712  lincons = savelinconss[c];
713  assert( lincons != NULL );
714  lhs = SCIPgetLhsLinear(scip, lincons);
715  rhs = SCIPgetRhsLinear(scip, lincons);
716  nvars = SCIPgetNVarsLinear(scip, lincons);
717  vars = SCIPgetVarsLinear(scip, lincons);
718  vals = SCIPgetValsLinear(scip, lincons);
719 
720  /* handle linear constraint for quadratic constraint update */
721  SCIP_CALL( presolveAddKKTLinearCons(scip, objcons, SCIPconsGetName(lincons),
722  vars, vals, lhs, rhs, nvars, varhash, dualconss, ndualconss, naddconss) );
723  }
724 
725  /* remove linear constraints if lhs != rhs, since they are now redundant; their feasibility is already expressed
726  * by s >= 0, where s is the new slack variable that we introduced for these linear constraints */
727  for( c = nlinconss-1; c >= 0; --c )
728  {
729  SCIP_CONS* lincons;
730 
731  lincons = savelinconss[c];
732  assert( savelinconss[c] != NULL );
733 
734  if( ! SCIPisFeasEQ(scip, SCIPgetLhsLinear(scip, lincons), SCIPgetRhsLinear(scip, lincons)) )
735  {
736  SCIP_CALL( SCIPdelCons(scip, savelinconss[c]) );
737  ++(*ndelconss);
738  }
739  }
740 
741  return SCIP_OKAY;
742 }
743 
744 /** handle knapsack constraints for quadratic constraint update, see the documentation of the function
745  * presolveAddKKTLinearCons() for an explanation
746  */
747 static
749  SCIP* scip, /**< SCIP pointer */
750  SCIP_CONS* objcons, /**< objective constraint */
751  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
752  SCIP_CONS** dualconss, /**< array with dual constraints */
753  int* ndualconss, /**< pointer to store number of dual constraints */
754  int* naddconss, /**< buffer to increase with number of created additional constraints */
755  int* ndelconss /**< buffer to increase with number of deleted constraints */
756  )
757 {
758  SCIP_CONSHDLR* conshdlr;
759  SCIP_CONS** conss;
760  int nconss;
761  int c;
762 
763  assert( scip != NULL );
764  assert( objcons != NULL );
765  assert( varhash != NULL );
766  assert( dualconss != NULL );
767  assert( ndualconss != NULL );
768  assert( naddconss != NULL );
769  assert( ndelconss != NULL );
770 
771  conshdlr = SCIPfindConshdlr(scip, "knapsack");
772  if( conshdlr == NULL )
773  return SCIP_OKAY;
774 
775  nconss = SCIPconshdlrGetNConss(conshdlr);
776  conss = SCIPconshdlrGetConss(conshdlr);
777 
778  /* loop through knapsack constraints */
779  for( c = 0; c < nconss; ++c )
780  {
781  SCIP_CONS* cons;
782  SCIP_VAR** vars;
783  SCIP_Longint* weights;
784  SCIP_Real* vals;
785  SCIP_Real lhs;
786  SCIP_Real rhs;
787  int nvars;
788  int v;
789 
790  /* get data of constraint */
791  cons = conss[c];
792  assert( cons != NULL );
793  lhs = -SCIPinfinity(scip);
794  rhs = (SCIP_Real) SCIPgetCapacityKnapsack(scip, cons);
795  nvars = SCIPgetNVarsKnapsack(scip, cons);
796  vars = SCIPgetVarsKnapsack(scip, cons);
797  weights = SCIPgetWeightsKnapsack(scip, cons);
798 
799  /* set coefficients of variables */
800  SCIP_CALL( SCIPallocBufferArray(scip, &vals, nvars) );
801  for( v = 0; v < nvars; ++v )
802  vals[v] = (SCIP_Real) weights[v];
803 
804  /* handle linear constraint for quadratic constraint update */
805  SCIP_CALL( presolveAddKKTLinearCons(scip, objcons, SCIPconsGetName(cons),
806  vars, vals, lhs, rhs, nvars, varhash, dualconss, ndualconss, naddconss) );
807 
808  /* free buffer array */
809  SCIPfreeBufferArray(scip, &vals);
810  }
811 
812  /* remove knapsack constraints, since they are now redundant; their feasibility is already expressed
813  * by s >= 0, where s is the new slack variable that we introduced for these linear constraints */
814  for( c = nconss-1; c >= 0; --c )
815  {
816  assert( conss[c] != NULL );
817  SCIP_CALL( SCIPdelCons(scip, conss[c]) );
818  ++(*ndelconss);
819  }
820 
821  return SCIP_OKAY;
822 }
823 
824 /** handle set packing constraints for quadratic constraint update, see the documentation of the function
825  * presolveAddKKTLinearCons() for an explanation
826  */
827 static
829  SCIP* scip, /**< SCIP pointer */
830  SCIP_CONS* objcons, /**< objective constraint */
831  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
832  SCIP_CONS** dualconss, /**< array with dual constraints */
833  int* ndualconss, /**< pointer to store number of dual constraints */
834  int* naddconss, /**< buffer to increase with number of created additional constraints */
835  int* ndelconss /**< buffer to increase with number of deleted constraints */
836  )
837 {
838  SCIP_CONSHDLR* conshdlr;
839  SCIP_CONS** conss;
840  int nconss;
841  int c;
842 
843  assert( scip != NULL );
844  assert( objcons != NULL );
845  assert( varhash != NULL );
846  assert( dualconss != NULL );
847  assert( ndualconss != NULL );
848  assert( naddconss != NULL );
849  assert( ndelconss != NULL );
850 
851  conshdlr = SCIPfindConshdlr(scip, "setppc");
852  if( conshdlr == NULL )
853  return SCIP_OKAY;
854 
855  nconss = SCIPconshdlrGetNConss(conshdlr);
856  conss = SCIPconshdlrGetConss(conshdlr);
857 
858  /* loop through linear constraints */
859  for( c = 0; c < nconss; ++c )
860  {
861  SCIP_SETPPCTYPE type;
862  SCIP_CONS* cons;
863  SCIP_VAR** vars;
864  SCIP_Real* vals;
865  SCIP_Real lhs;
866  SCIP_Real rhs;
867  int nvars;
868  int v;
869 
870  /* get data of constraint */
871  cons = conss[c];
872  assert( cons != NULL );
873 
874  /* get setppc type */
875  type = SCIPgetTypeSetppc(scip, cons);
876  lhs = -SCIPinfinity(scip);
877  rhs = SCIPinfinity(scip);
878  switch( type )
879  {
881  lhs = 1.0;
882  rhs = 1.0;
883  break;
885  rhs = 1.0;
886  break;
888  lhs = 1.0;
889  break;
890  default:
891  SCIPerrorMessage("unknown setppc type\n");
892  return SCIP_INVALIDDATA;
893  }
894 
895  nvars = SCIPgetNVarsSetppc(scip, cons);
896  vars = SCIPgetVarsSetppc(scip, cons);
897 
898  /* set coefficients of variables */
899  SCIP_CALL( SCIPallocBufferArray(scip, &vals, nvars) );
900  for( v = 0; v < nvars; ++v )
901  vals[v] = 1.0;
902 
903  /* handle linear constraint for quadratic constraint update */
904  SCIP_CALL( presolveAddKKTLinearCons(scip, objcons, SCIPconsGetName(cons),
905  vars, vals, lhs, rhs, nvars, varhash, dualconss, ndualconss, naddconss) );
906 
907  /* free buffer array */
908  SCIPfreeBufferArray(scip, &vals);
909  }
910 
911  /* remove set packing constraints if lhs != rhs, since they are now redundant; their feasibility is already expressed
912  * by s >= 0, where s is the new slack variable that we introduced for these linear constraints */
913  for( c = nconss-1; c >= 0; --c )
914  {
915  assert( conss[c] != NULL );
916 
917  if( SCIPgetTypeSetppc(scip, conss[c]) != SCIP_SETPPCTYPE_PARTITIONING )
918  {
919  assert( SCIPgetTypeSetppc(scip, conss[c]) == SCIP_SETPPCTYPE_PACKING
920  || SCIPgetTypeSetppc(scip, conss[c]) == SCIP_SETPPCTYPE_COVERING );
921 
922  SCIP_CALL( SCIPdelCons(scip, conss[c]) );
923  ++(*ndelconss);
924  }
925  }
926 
927  return SCIP_OKAY;
928 }
929 
930 /** handle varbound constraints for quadratic constraint update, see the documentation of the function
931  * presolveAddKKTLinearCons() for an explanation
932  */
933 static
935  SCIP* scip, /**< SCIP pointer */
936  SCIP_CONS* objcons, /**< objective constraint */
937  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
938  SCIP_CONS** dualconss, /**< array with dual constraints */
939  int* ndualconss, /**< pointer to store number of dual constraints */
940  int* naddconss, /**< buffer to increase with number of created additional constraints */
941  int* ndelconss /**< buffer to increase with number of deleted constraints */
942  )
943 {
944  SCIP_CONSHDLR* conshdlr;
945  SCIP_CONS** conss;
946  int nconss;
947  int c;
948 
949  assert( scip != NULL );
950  assert( objcons != NULL );
951  assert( varhash != NULL );
952  assert( dualconss != NULL );
953  assert( ndualconss != NULL );
954  assert( naddconss != NULL );
955  assert( ndelconss != NULL );
956 
957  conshdlr = SCIPfindConshdlr(scip, "varbound");
958  if( conshdlr == NULL )
959  return SCIP_OKAY;
960 
961  nconss = SCIPconshdlrGetNConss(conshdlr);
962  conss = SCIPconshdlrGetConss(conshdlr);
963 
964  /* loop through linear constraints */
965  for( c = 0; c < nconss; ++c )
966  {
967  SCIP_CONS* cons;
968  SCIP_VAR** vars;
969  SCIP_Real* vals;
970  SCIP_Real lhs;
971  SCIP_Real rhs;
972  int nvars;
973 
974  /* allocate buffer arrays */
975  SCIP_CALL( SCIPallocBufferArray(scip, &vars, 2) );
976  SCIP_CALL( SCIPallocBufferArray(scip, &vals, 2) );
977 
978  /* get data of constraint */
979  cons = conss[c];
980  assert( cons != NULL );
981 
982  lhs = SCIPgetLhsVarbound(scip, cons);
983  rhs = SCIPgetRhsVarbound(scip, cons);
984  vars[0] = SCIPgetVarVarbound(scip, cons);
985  vars[1] = SCIPgetVbdvarVarbound(scip, cons);
986  vals[0] = 1.0;
987  vals[1] = SCIPgetVbdcoefVarbound(scip, cons);
988  nvars = 2;
989 
990  /* handle linear constraint for quadratic constraint update */
991  SCIP_CALL( presolveAddKKTLinearCons(scip, objcons, SCIPconsGetName(cons),
992  vars, vals, lhs, rhs, nvars, varhash, dualconss, ndualconss, naddconss) );
993 
994  /* free buffer array */
995  SCIPfreeBufferArray(scip, &vals);
996  SCIPfreeBufferArray(scip, &vars);
997  }
998 
999  /* remove varbound constraints if lhs != rhs, since they are now redundant; their feasibility is already expressed
1000  * by s >= 0, where s is the new slack variable that we introduced for these linear constraints */
1001  for( c = nconss-1; c >= 0; --c )
1002  {
1003  SCIP_CONS* cons;
1004 
1005  cons = conss[c];
1006  assert( cons != NULL );
1007 
1008  if( ! SCIPisFeasEQ(scip, SCIPgetLhsVarbound(scip, cons), SCIPgetRhsVarbound(scip, cons)) )
1009  {
1010  SCIP_CALL( SCIPdelCons(scip, cons) );
1011  ++(*ndelconss);
1012  }
1013  }
1014 
1015  return SCIP_OKAY;
1016 }
1017 
1018 /** handle logicor constraints for quadratic constraint update, see the documentation of the function
1019  * presolveAddKKTLinearCons() for an explanation
1020  */
1021 static
1023  SCIP* scip, /**< SCIP pointer */
1024  SCIP_CONS* objcons, /**< objective constraint */
1025  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
1026  SCIP_CONS** dualconss, /**< array with dual constraints */
1027  int* ndualconss, /**< pointer to store number of dual constraints */
1028  int* naddconss, /**< buffer to increase with number of created additional constraints */
1029  int* ndelconss /**< buffer to increase with number of deleted constraints */
1030  )
1031 {
1032  SCIP_CONSHDLR* conshdlr;
1033  SCIP_CONS** conss;
1034  int nconss;
1035  int c;
1036 
1037  assert( scip != NULL );
1038  assert( objcons != NULL );
1039  assert( varhash != NULL );
1040  assert( dualconss != NULL );
1041  assert( ndualconss != NULL );
1042  assert( naddconss != NULL );
1043  assert( ndelconss != NULL );
1044 
1045  conshdlr = SCIPfindConshdlr(scip, "logicor");
1046  if( conshdlr == NULL )
1047  return SCIP_OKAY;
1048 
1049  nconss = SCIPconshdlrGetNConss(conshdlr);
1050  conss = SCIPconshdlrGetConss(conshdlr);
1051 
1052  /* loop through linear constraints */
1053  for( c = 0; c < nconss; ++c )
1054  {
1055  SCIP_CONS* cons;
1056  SCIP_VAR** vars;
1057  SCIP_Real* vals;
1058  SCIP_Real lhs;
1059  SCIP_Real rhs;
1060  int nvars;
1061  int v;
1062 
1063  /* get data of constraint */
1064  cons = conss[c];
1065  assert( cons != NULL );
1066 
1067  /* get setppc type */
1068  lhs = 1.0;
1069  rhs = SCIPinfinity(scip);
1070 
1071  nvars = SCIPgetNVarsLogicor(scip, cons);
1072  vars = SCIPgetVarsLogicor(scip, cons);
1073 
1074  /* set coefficients of variables */
1075  SCIP_CALL( SCIPallocBufferArray(scip, &vals, nvars) );
1076  for( v = 0; v < nvars; ++v )
1077  vals[v] = 1.0;
1078 
1079  /* handle linear constraint for quadratic constraint update */
1080  SCIP_CALL( presolveAddKKTLinearCons(scip, objcons, SCIPconsGetName(cons),
1081  vars, vals, lhs, rhs, nvars, varhash, dualconss, ndualconss, naddconss) );
1082 
1083  /* free buffer array */
1084  SCIPfreeBufferArray(scip, &vals);
1085  }
1086 
1087  /* remove logicor constraints, since they are now redundant; their feasibility is already expressed
1088  * by s >= 0, where s is the new slack variable that we introduced for these linear constraints */
1089  for( c = nconss-1; c >= 0; --c )
1090  {
1091  assert( conss[c] != NULL );
1092 
1093  SCIP_CALL( SCIPdelCons(scip, conss[c]) );
1094  ++(*ndelconss);
1095  }
1096 
1097  return SCIP_OKAY;
1098 }
1099 
1100 /** handle aggregated variables for quadratic constraint update @n
1101  * we apply the function presolveAddKKTLinearCons() to the aggregation constraint, see the documentation of this
1102  * function for further information
1103  */
1104 static
1106  SCIP* scip, /**< SCIP pointer */
1107  SCIP_CONS* objcons, /**< objective constraint */
1108  SCIP_VAR** agrvars, /**< aggregated variables */
1109  int nagrvars, /**< number of aggregated variables */
1110  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
1111  SCIP_CONS** dualconss, /**< array with dual constraints */
1112  int* ndualconss, /**< pointer to store number of dual constraints */
1113  int* naddconss /**< buffer to increase with number of created additional constraints */
1114  )
1115 {
1116  int v;
1117 
1118  assert( scip != NULL );
1119  assert( objcons != NULL );
1120  assert( agrvars != NULL );
1121  assert( varhash != NULL );
1122  assert( dualconss != NULL );
1123  assert( ndualconss != NULL );
1124  assert( naddconss != NULL );
1125 
1126  /* loop through variables */
1127  for( v = 0; v < nagrvars; ++v )
1128  {
1129  SCIP_VAR* var;
1130  SCIP_VAR** vars = NULL;
1131  SCIP_Real* vals = NULL;
1132  SCIP_Real lhs;
1133  SCIP_Real rhs;
1134  int nvars;
1135 
1136  var = agrvars[v];
1137 
1139  {
1140  SCIP_Real constant;
1141 
1142  SCIP_CALL( SCIPallocBufferArray(scip, &vars, 2) );
1143  SCIP_CALL( SCIPallocBufferArray(scip, &vals, 2) );
1144 
1145  /* get aggregation variable */
1146  constant = SCIPvarGetAggrConstant(var);
1147  vars[0] = SCIPvarGetAggrVar(var);
1148  vals[0] = SCIPvarGetAggrScalar(var);
1149  vars[1] = var;
1150  vals[1] = -1.0;
1151  lhs = -constant;
1152  rhs = -constant;
1153  nvars = 2;
1154  }
1155  else if( SCIPvarGetStatus(var) == SCIP_VARSTATUS_MULTAGGR )
1156  {
1157  SCIP_Real* scalars;
1158  SCIP_VAR** multvars;
1159  SCIP_Real constant;
1160  int nmultvars;
1161  int nbuffer;
1162  int j;
1163 
1164  nmultvars = SCIPvarGetMultaggrNVars(var);
1165  nbuffer = nmultvars+1;
1166 
1167  SCIP_CALL( SCIPallocBufferArray(scip, &vars, nbuffer) );
1168  SCIP_CALL( SCIPallocBufferArray(scip, &vals, nbuffer) );
1169 
1170  /* get aggregation variables */
1171  multvars = SCIPvarGetMultaggrVars(var);
1172  scalars = SCIPvarGetMultaggrScalars(var);
1173  constant = SCIPvarGetMultaggrConstant(var);
1174 
1175  /* add multi-aggregated variables to array */
1176  for( j = 0; j < nmultvars; ++j )
1177  {
1178  vars[j] = multvars[j];
1179  vals[j] = scalars[j];
1180  }
1181 
1182  /* add new variable to array */
1183  vars[nmultvars] = var;
1184  vals[nmultvars] = -1.0;
1185  lhs = -constant;
1186  rhs = -constant;
1187  nvars = nmultvars + 1;
1188  }
1189  else if( SCIPvarGetStatus(var) == SCIP_VARSTATUS_NEGATED )
1190  {
1191  SCIP_VAR* negvar;
1192  SCIP_Real negconst;
1193 
1194  /* get negation variable and negation offset */
1195  negvar = SCIPvarGetNegationVar(var);
1196  negconst = SCIPvarGetNegationConstant(var);
1197 
1198  SCIP_CALL( SCIPallocBufferArray(scip, &vars, 2) );
1199  SCIP_CALL( SCIPallocBufferArray(scip, &vals, 2) );
1200 
1201  vars[0] = negvar;
1202  vars[1] = var;
1203  vals[0] = 1.0;
1204  vals[1] = 1.0;
1205  lhs = negconst;
1206  rhs = negconst;
1207  nvars = 2;
1208  }
1209  else if( SCIPvarGetStatus(var) == SCIP_VARSTATUS_FIXED )
1210  {
1211  SCIP_Real lb;
1212  SCIP_Real ub;
1213 
1214  lb = SCIPvarGetLbGlobal(var);
1215  ub = SCIPvarGetUbGlobal(var);
1216  assert( SCIPisFeasEQ(scip, lb, ub) );
1217 
1218  if( SCIPisFeasZero(scip, lb) && SCIPisFeasZero(scip, ub) )
1219  continue;
1220  else
1221  {
1222  SCIP_CALL( SCIPallocBufferArray(scip, &vars, 1) );
1223  SCIP_CALL( SCIPallocBufferArray(scip, &vals, 1) );
1224 
1225  vars[0] = var;
1226  vals[0] = 1.0;
1227  lhs = lb;
1228  rhs = lb;
1229  nvars = 1;
1230  }
1231  }
1232  else
1233  {
1234  SCIPerrorMessage("unexpected variable status\n");
1235  return SCIP_ERROR;
1236  }
1237 
1238  if( nvars > 0 )
1239  {
1240  /* handle aggregation constraint for quadratic constraint update */
1241  SCIP_CALL( presolveAddKKTLinearCons(scip, objcons, SCIPvarGetName(var),
1242  vars, vals, lhs, rhs, nvars, varhash, dualconss, ndualconss, naddconss) );
1243  }
1244 
1245  SCIPfreeBufferArrayNull(scip, &vars);
1246  SCIPfreeBufferArrayNull(scip, &vals);
1247  }
1248 
1249  return SCIP_OKAY;
1250 }
1251 
1252 /** handle bilinear terms of quadratic constraint for quadratic constraint update
1253  *
1254  * For the two variables of each bilinear term
1255  * 1. create the dual constraints (i.e., the two rows of \f$Q x + c + A^T \mu = 0\f$) associated to these variables, if not
1256  * done already
1257  * 2. create the dual variables and the complementarity constraints for the lower and upper bound constraints of the two
1258  * variables of the bilinear term, if not done already
1259  * 3. add the coefficient \f$Q_{ij}\f$ of the bilinear term to the dual constraint
1260  *
1261  * for steps 1 and 2 see the documentation of createKKTDualCons() for further information.
1262  **/
1263 static
1265  SCIP* scip, /**< SCIP pointer */
1266  SCIP_CONS* objcons, /**< objective constraint */
1267  SCIP_BILINTERM* bilinterms, /**< bilinear terms of quadratic constraint */
1268  int nbilinterms, /**< number of bilinear terms of quadratic constraint */
1269  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
1270  SCIP_Real scale, /**< scale factor of quadratic constraint */
1271  SCIP_CONS** dualconss, /**< array with dual constraints */
1272  int* ndualconss, /**< pointer to store number of dual constraints */
1273  int* naddconss /**< buffer to increase with number of created additional constraints */
1274  )
1275 {
1276  int j;
1277 
1278  assert( scip != NULL );
1279  assert( objcons != NULL );
1280  assert( varhash != NULL );
1281  assert( dualconss != NULL );
1282  assert( ndualconss != NULL );
1283  assert( naddconss != NULL );
1284 
1285  /* return if there are no bilinear terms */
1286  if( bilinterms == NULL )
1287  return SCIP_OKAY;
1288 
1289  /* loop through bilinear terms of quadratic constraint */
1290  for( j = 0; j < nbilinterms; ++j )
1291  {
1292  int i;
1293 
1294  /* quadratic matrix has to be symmetric; therefore, split bilinear terms into two parts */
1295  for( i = 0; i < 2; ++i )
1296  {
1297  SCIP_CONS* dualcons = NULL; /* dual constraint associated to variable */
1298  SCIP_VAR* bilvar1;
1299  SCIP_VAR* bilvar2;
1300 
1301  if( i == 0 )
1302  {
1303  bilvar1 = bilinterms[j].var1;
1304  bilvar2 = bilinterms[j].var2;
1305  }
1306  else
1307  {
1308  bilvar1 = bilinterms[j].var2;
1309  bilvar2 = bilinterms[j].var1;
1310  }
1311 
1312  /* create/get dual constraint associated to variable 'bilvar1';
1313  * if variable does not already exist in hashmap then create dual variables for its bounds */
1314  SCIP_CALL( createKKTDualCons(scip, objcons, bilvar1, varhash, dualconss, ndualconss, &dualcons, naddconss) );
1315  assert( dualcons != NULL );
1316 
1317  /* add variable to dual constraint */
1318  assert( ! SCIPisFeasZero(scip, scale) );
1319  SCIP_CALL( SCIPaddCoefLinear(scip, dualcons, bilvar2, bilinterms[j].coef / scale) );
1320  }
1321  }
1322 
1323  return SCIP_OKAY;
1324 }
1325 
1326 /** handle quadratic terms of quadratic constraint for quadratic constraint update
1327  *
1328  * For each quadratic term variable
1329  * 1. create the dual constraint (i.e., a row of \f$Q x + c + A^T \mu = 0\f$) associated to this variable, if not done
1330  * already
1331  * 2. create the dual variables and the complementarity constraints for the lower and upper bound constraints of this
1332  * variable, if not done already
1333  * 3. add the coefficient \f$Q_{ii}\f$ of this variable to the dual constraint
1334  *
1335  * for steps 1 and 2 see the documentation of createKKTDualCons() for further information.
1336  **/
1337 static
1339  SCIP* scip, /**< SCIP pointer */
1340  SCIP_CONS* objcons, /**< objective constraint */
1341  SCIP_QUADVARTERM* quadterms, /**< quadratic terms of quadratic constraint */
1342  int nquadterms, /**< number of quadratic terms of quadratic constraint */
1343  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
1344  SCIP_Real scale, /**< scale factor of quadratic constraint */
1345  SCIP_CONS** dualconss, /**< array with dual constraints */
1346  int* ndualconss, /**< pointer to store number of dual constraints */
1347  int* naddconss /**< buffer to increase with number of created additional constraints */
1348  )
1349 {
1350  int j;
1351 
1352  assert( scip != NULL );
1353  assert( objcons != NULL );
1354  assert( varhash != NULL );
1355  assert( dualconss != NULL );
1356  assert( ndualconss != NULL );
1357  assert( naddconss != NULL );
1358 
1359  /* return if there are no quadratic terms */
1360  if( quadterms == NULL )
1361  return SCIP_OKAY;
1362 
1363  /* loop through quadratic terms */
1364  for( j = 0; j < nquadterms; ++j )
1365  {
1366  SCIP_CONS* dualcons = NULL; /* dual constraint associated to variable */
1367  SCIP_VAR* quadvar;
1368 
1369  quadvar = quadterms[j].var;
1370 
1371  /* create/get dual constraint associated to variable 'bilvar1';
1372  * if variable does not already exist in hashmap then create dual variables for its bounds */
1373  SCIP_CALL( createKKTDualCons(scip, objcons, quadvar, varhash, dualconss, ndualconss, &dualcons, naddconss) );
1374  assert( dualcons != NULL );
1375 
1376  /* add variable to dual constraint */
1377  assert( ! SCIPisFeasZero(scip, scale) );
1378  SCIP_CALL( SCIPaddCoefLinear(scip, dualcons, quadvar, (SCIP_Real) quadterms[j].sqrcoef * 2 / scale) );
1379  }
1380 
1381  return SCIP_OKAY;
1382 }
1383 
1384 /** handle linear terms of quadratic constraint for quadratic constraint update
1385  *
1386  * For each linear term variable
1387  * 1. create the dual constraint (i.e., a row of \f$Q x + c + A^T \mu = 0\f$) associated to this variable, if not done
1388  * already
1389  * 2. create the dual variables and the complementarity constraints for the lower and upper bound constraints of this
1390  * variable, if not done already
1391  * 3. add the right hand side \f$-c_i\f$ to the dual constraint
1392  * 4. add \f$c_i\f$ to the objective constraint \f$1/2 ( c^T x + b^T \mu) = t\f$, where t is the objective variable
1393  *
1394  * for steps 1 and 2 see the documentation of createKKTDualCons() for further information.
1395  **/
1396 static
1398  SCIP* scip, /**< SCIP pointer */
1399  SCIP_CONS* objcons, /**< objective constraint */
1400  SCIP_VAR** lintermvars, /**< linear terms of quadratic constraint */
1401  SCIP_Real* lintermcoefs, /**< coefficients of linear terms of quadratic constraint */
1402  int nlintermvars, /**< number of linear terms of quadratic constraints */
1403  SCIP_QUADVARTERM* quadterms, /**< quadratic terms of quadratic constraint */
1404  int nquadterms, /**< number of quadratic terms of quadratic constraint */
1405  SCIP_HASHMAP* varhash, /**< hash map from variable to index of linear constraint */
1406  SCIP_VAR* objvar, /**< variable of objective function */
1407  SCIP_Real scale, /**< scale factor of quadratic constraint */
1408  SCIP_CONS** dualconss, /**< array with dual constraints */
1409  int* ndualconss, /**< pointer to store number of dual constraints */
1410  int* naddconss /**< buffer to increase with number of created additional constraints */
1411  )
1412 {
1413  int j;
1414 
1415  assert( scip != NULL );
1416  assert( objcons != NULL );
1417  assert( lintermcoefs != NULL );
1418  assert( varhash != NULL );
1419  assert( objvar != NULL );
1420  assert( dualconss != NULL );
1421  assert( ndualconss != NULL );
1422  assert( naddconss != NULL );
1423 
1424  /* loop through linear terms of quadratic constraint */
1425  if( lintermvars != NULL )
1426  {
1427  assert( lintermcoefs != NULL );
1428  for( j = 0; j < nlintermvars; ++j )
1429  {
1430  SCIP_VAR* var;
1431 
1432  var = lintermvars[j];
1433 
1434  if( var != objvar )
1435  {
1436  SCIP_CONS* dualcons = NULL; /* dual constraint associated to variable */
1437  SCIP_Real coef;
1438 
1439  /* create/get dual constraint associated to variable;
1440  * if variable does not already exist in hashmap then create dual variables for its bounds */
1441  SCIP_CALL( createKKTDualCons(scip, objcons, var, varhash, dualconss, ndualconss, &dualcons, naddconss) );
1442  assert( dualcons != NULL );
1443 
1444  /* get coefficient of variable in quadratic constraint */
1445  coef = lintermcoefs[j];
1446 
1447  /* change lhs and rhs of dual constraint */
1448  assert( ! SCIPisFeasZero(scip, scale) );
1449  SCIP_CALL( SCIPchgLhsLinear(scip, dualcons, SCIPgetLhsLinear(scip, dualcons) - coef / scale) );
1450  SCIP_CALL( SCIPchgRhsLinear(scip, dualcons, SCIPgetRhsLinear(scip, dualcons) - coef / scale) );
1451 
1452  /* add variable to objective constraint */
1453  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, var, coef / (scale * 2)) );
1454  }
1455  }
1456  }
1457 
1458  /* loop through linear terms that are part of a quadratic term of quadratic constraint */
1459  if( quadterms != NULL )
1460  {
1461  for( j = 0; j < nquadterms; ++j )
1462  {
1463  SCIP_CONS* dualcons;
1464  SCIP_Real coef;
1465  SCIP_VAR* var;
1466  int ind;
1467 
1468  var = quadterms[j].var;
1469  coef = quadterms[j].lincoef;
1470  assert( var != objvar );
1471 
1472  /* get dual constraint associated to variable (has already been created in function
1473  * presolveAddKKTQuadQuadraticTerms() */
1474  assert( SCIPhashmapExists(varhash, var) );
1475  ind = (int) (size_t) SCIPhashmapGetImage(varhash, var);
1476  dualcons = dualconss[ind];
1477  assert( dualcons != NULL );
1478 
1479  /* change lhs and rhs of dual constraint */
1480  assert( ! SCIPisFeasZero(scip, scale) );
1481  SCIP_CALL( SCIPchgLhsLinear(scip, dualcons, SCIPgetLhsLinear(scip, dualcons) -coef / scale) );
1482  SCIP_CALL( SCIPchgRhsLinear(scip, dualcons, SCIPgetRhsLinear(scip, dualcons) -coef / scale) );
1483 
1484  /* add variable to objective constraint */
1485  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, var, coef / (scale * 2)) );
1486  }
1487  }
1488 
1489  return SCIP_OKAY;
1490 }
1491 
1492 /** checks for a given constraint whether it is the objective function of a (mixed-binary) quadratic program
1493  * \f[
1494  * \begin{array}{ll}
1495  * \min & z \\
1496  * s.t. & x^T Q x + c^T x + d <= z \\
1497  * & A x \leq b, \\
1498  * & x \in \{0, 1\}^{p} \times R^{n-p},
1499  * \end{array}
1500  * \f]
1501  * which is equivalent to
1502  * \f[
1503  * \begin{array}{ll}
1504  * \min & x^T Q x + c^T x + d \\
1505  * s.t. & A x \leq b, \\
1506  * & x \in \{0, 1\}^{p} \times R^{n-p}.
1507  * \end{array}
1508  * \f]
1509  *
1510  *
1511  * We check whether
1512  * 1. there is a single quadratic constraint that can be written as \f$x^T Q x + c^T x + d \leq t\f$
1513  * 2. all other constraints are linear
1514  * 3. all integer variables are binary if allowbinary = TRUE, or all variables are continuous if allowbinary = FALSE
1515  * 4. t is the only variable in the objective and doesn't appear in any other constraint
1516  */
1517 static
1519  SCIP* scip, /**< SCIP data structure */
1520  SCIP_CONSHDLR* quadconshdlr, /**< constraint handler data structure */
1521  SCIP_CONS* cons, /**< quadratic constraint */
1522  SCIP_Bool allowbinary, /**< if TRUE then allow binary variables in the problem, if FALSE then all
1523  * variables have to be continuous */
1524  SCIP_VAR** objvar, /**< pointer to store the objective variable @p z */
1525  SCIP_Real* scale, /**< pointer to store the value by which we have to scale the quadratic
1526  * constraint such that the objective variable @p z has coefficient -1 */
1527  SCIP_Real* objrhs, /**< pointer to store the right hand side @p -d of the objective constraint */
1528  SCIP_Bool* isqp /**< pointer to store whether the problem is a (mixed-binary) QP */
1529  )
1530 {
1531  SCIP_CONSHDLR* conshdlr;
1532  SCIP_VAR** lintermvars;
1533  SCIP_Real* lintermcoefs;
1534  int nconss = 0;
1535  SCIP_Real quadlhs;
1536  SCIP_Real quadrhs;
1537  SCIP_Real coef;
1538  SCIP_Real obj;
1539  int mayincrease;
1540  int maydecrease;
1541  int objind = -1;
1542 
1543  *objrhs = 0.0;
1544  *scale = 0.0;
1545  *isqp = FALSE;
1546  *objvar = NULL;
1547 
1548  /* desired structure: there exists only one variable with nonzero objective value; this is the objective variable 'z' */
1549  if( SCIPgetNObjVars(scip) != 1 )
1550  return SCIP_OKAY;
1551 
1552  /* desired structure: all integer variables are binary; if the parameter 'allowbinary' is set to FALSE, then all
1553  * variables have to be continuous */
1554  if( SCIPgetNIntVars(scip) > 0 || ( ! allowbinary && SCIPgetNBinVars(scip) > 0 ) )
1555  return SCIP_OKAY;
1556 
1557  /* desired structure: there exists only one quadratic constraint */
1558  if( SCIPconshdlrGetNConss(quadconshdlr) != 1 )
1559  return SCIP_OKAY;
1560 
1561  /* desired structure: the constraint has to take one of the three forms
1562  * i) x^T Q x + c^T x <= d
1563  * ii) x^T Q x + c^T x >= d
1564  * iii) x^T Q x + c^T x == d
1565  * the case a <= x^T Q x + c^T x <= d with 'a' and 'd' finite and a != d is not allowed.
1566  */
1567  quadlhs = SCIPgetLhsQuadratic(scip, cons);
1568  quadrhs = SCIPgetRhsQuadratic(scip, cons);
1569  if( ! SCIPisFeasEQ(scip, quadlhs, quadrhs) && ! SCIPisInfinity(scip, -quadlhs) && ! SCIPisInfinity(scip, quadrhs) )
1570  return SCIP_OKAY;
1571 
1572  /* get number of linear constraints (including special cases of linear constraints) */
1573  conshdlr = SCIPfindConshdlr(scip, "linear");
1574  if( conshdlr != NULL )
1575  nconss += SCIPconshdlrGetNConss(conshdlr);
1576 
1577  conshdlr = SCIPfindConshdlr(scip, "setppc");
1578  if( conshdlr != NULL )
1579  nconss += SCIPconshdlrGetNConss(conshdlr);
1580 
1581  conshdlr = SCIPfindConshdlr(scip, "knapsack");
1582  if( conshdlr != NULL )
1583  nconss += SCIPconshdlrGetNConss(conshdlr);
1584 
1585  conshdlr = SCIPfindConshdlr(scip, "varbound");
1586  if( conshdlr != NULL )
1587  nconss += SCIPconshdlrGetNConss(conshdlr);
1588 
1589  conshdlr = SCIPfindConshdlr(scip, "logicor");
1590  if( conshdlr != NULL )
1591  nconss += SCIPconshdlrGetNConss(conshdlr);
1592 
1593  /* desired structure: all the nonquadratic constraints are linear constraints */
1594  if( nconss != SCIPgetNConss(scip) - 1 )
1595  return SCIP_OKAY;
1596 
1597  /* get variables that are in the linear term of the quadratic constraint */
1598  lintermvars = SCIPgetLinearVarsQuadratic(scip, cons);
1599  lintermcoefs = SCIPgetCoefsLinearVarsQuadratic(scip, cons);
1600 
1601  /* compute the objective shift of the QP. Note that
1602  *
1603  * min z s.t. x^T Q x + c^T x <= d + z
1604  * Ax <= b
1605  *
1606  * is equivalent to
1607  *
1608  * min x^T Q x + c^T x - d s.t. Ax <= b
1609  *
1610  * Here, -d is the objective shift. We define b to be the right hand side of the objective constraint.
1611  */
1612  if( ! SCIPisInfinity(scip, -quadlhs) )
1613  *objrhs = quadlhs;
1614  else
1615  *objrhs = quadrhs;
1616  assert( ! SCIPisInfinity(scip, REALABS(*objrhs)) );
1617 
1618  /* search for the objective variable 'objvar' in the linear term of quadratic constraint (it is already known that
1619  * at most one variable has a nonzero objective value); additionally, check the sign of the objective variable */
1620  maydecrease = SCIPgetLinvarMayDecreaseQuadratic(scip, cons);
1621  mayincrease = SCIPgetLinvarMayIncreaseQuadratic(scip, cons);
1622  if( maydecrease < 0 && mayincrease < 0 )
1623  return SCIP_OKAY;
1624  else if( maydecrease >= 0 )
1625  {
1626  objind = maydecrease;
1627 
1628  /* if both mayincrease and maydecrese are nonnegative, then check objective coefficient */
1629  if( mayincrease >= 0 && SCIPisFeasZero(scip, SCIPvarGetObj(lintermvars[maydecrease])) )
1630  objind = mayincrease;
1631  }
1632  else
1633  objind = mayincrease;
1634  assert( objind < SCIPgetNLinearVarsQuadratic(scip, cons) );
1635 
1636  *objvar = lintermvars[objind];
1637  coef = lintermcoefs[objind];
1638  obj = SCIPvarGetObj(*objvar);
1639 
1640  /* check sign of coefficient */
1641  if( SCIPisFeasPositive(scip, obj)
1642  && ( ( SCIPisFeasNegative(scip, coef) && SCIPisFeasEQ(scip, quadrhs, *objrhs) )
1643  || ( SCIPisFeasPositive(scip, coef) && SCIPisFeasEQ(scip, quadlhs, *objrhs) )
1644  )
1645  )
1646  *scale = -1.0/coef; /* value by which we have to scale the quadratic constraint such that the objective variable
1647  * has coefficient -1 */
1648  else if( SCIPisFeasNegative(scip, obj)
1649  && ( ( SCIPisFeasNegative(scip, coef) && SCIPisFeasEQ(scip, quadlhs, *objrhs) )
1650  || ( SCIPisFeasPositive(scip, coef) && SCIPisFeasEQ(scip, quadrhs, *objrhs) )
1651  )
1652  )
1653  *scale = 1.0/coef; /* value by which we have to scale the quadratic constraint such that the objective variable
1654  * has coefficient 1 */
1655  else
1656  return SCIP_OKAY;
1657  assert( *objvar != NULL && ! SCIPisFeasZero(scip, SCIPvarGetObj(*objvar)) );
1658  assert( ! SCIPisFeasZero(scip, *scale) );
1659 
1660  /* scale the right hand side of the objective constraint */
1661  *objrhs = (*objrhs)/(*scale); /*lint !e414*/
1662 
1663  /* check whether 'objvar' is part of a linear constraint; if this is true then return
1664  * whether 'objvar' is part of a linear constraint can be deduced from the variable locks */
1665  if( SCIPisFeasEQ(scip, quadlhs, quadrhs) )
1666  {
1667  if( SCIPvarGetNLocksDown(*objvar) != 1 || SCIPvarGetNLocksUp(*objvar) != 1 )
1668  return SCIP_OKAY;
1669  }
1670  else
1671  {
1672  assert( SCIPisInfinity(scip, -quadlhs) || SCIPisInfinity(scip, quadrhs) );
1673 
1674  if( ( SCIPvarGetNLocksDown(*objvar) != 1 || SCIPvarGetNLocksUp(*objvar) != 0 )
1675  && ( SCIPvarGetNLocksDown(*objvar) != 0 || SCIPvarGetNLocksUp(*objvar) != 1 ) )
1676  return SCIP_OKAY;
1677  }
1678 
1679  *isqp = TRUE;
1680 
1681  return SCIP_OKAY;
1682 }
1683 
1684 /*
1685  * Callback methods of presolver
1686  */
1687 
1688 /** copy method for constraint handler plugins (called when SCIP copies plugins) */
1689 static
1690 SCIP_DECL_PRESOLCOPY(presolCopyQPKKTref)
1691 { /*lint --e{715}*/
1692  assert(scip != NULL);
1693  assert(presol != NULL);
1694  assert(strcmp(SCIPpresolGetName(presol), PRESOL_NAME) == 0);
1695 
1696  /* call inclusion method of presolver */
1698 
1699  return SCIP_OKAY;
1700 }
1701 
1702 
1703 /** destructor of presolver to free user data (called when SCIP is exiting) */
1704 static
1705 SCIP_DECL_PRESOLFREE(presolFreeQPKKTref)
1706 { /*lint --e{715}*/
1707  SCIP_PRESOLDATA* presoldata;
1708 
1709  /* free presolver data */
1710  presoldata = SCIPpresolGetData(presol);
1711  assert(presoldata != NULL);
1712 
1713  SCIPfreeBlockMemory(scip, &presoldata);
1714  SCIPpresolSetData(presol, NULL);
1715 
1716  return SCIP_OKAY;
1717 }
1718 
1719 
1720 /** execution method of presolver */
1721 static
1722 SCIP_DECL_PRESOLEXEC(presolExecQPKKTref)
1723 { /*lint --e{715}*/
1724  SCIP_PRESOLDATA* presoldata;
1725  SCIP_CONSHDLR* linconshdlr;
1726  SCIP_CONSHDLR* quadconshdlr;
1727  SCIP_CONS** conss;
1728  SCIP_CONS* cons;
1729 
1730  SCIP_QUADVARTERM* quadterms;
1731  SCIP_BILINTERM* bilinterms;
1732  int nquadterms;
1733  int nbilinterms;
1734 
1735  SCIP_VAR** lintermvars;
1736  SCIP_Real* lintermcoefs;
1737  int nlintermvars;
1738 
1739  SCIP_CONS** savelinconss = NULL;
1740  SCIP_CONS** linconss = NULL;
1741  int nlinconss = 0;
1742 
1743  SCIP_HASHMAP* varhash; /* hash map from variable to index of dual constraint */
1744  SCIP_CONS** dualconss; /* constraints associated to the Lagrangean function */
1745  int ndualconss = 0;
1746 
1747  SCIP_CONS* objcons;
1748  SCIP_VAR* objvar;
1749  SCIP_Real scale;
1750  SCIP_Real objrhs;
1751  SCIP_Bool isqp;
1752  int j;
1753 
1754  assert( scip != NULL );
1755  assert( naddconss != NULL );
1756  assert( ndelconss != NULL );
1757 
1758  /* desired structure: there exists only one quadratic constraint */
1759  quadconshdlr = SCIPfindConshdlr(scip,"quadratic");
1760  if( quadconshdlr == NULL || SCIPconshdlrGetNConss(quadconshdlr) != 1 )
1761  return SCIP_OKAY;
1762 
1763  /* get quadratic constraint */
1764  conss = SCIPconshdlrGetConss(quadconshdlr);
1765  cons = conss[0];
1766  assert( cons != NULL );
1767 
1768  SCIPdebugMsg(scip, "tries to add the KKT conditions for constraint <%s>\n", SCIPconsGetName(cons));
1769 
1770  /* get presolver data */
1771  presoldata = SCIPpresolGetData(presol);
1772  assert(presoldata != NULL);
1773 
1774  /* desired structure: matrix associated to quadratic constraint is indefinite;
1775  * otherwise, the problem usually can be solved faster by standard methods. */
1776  SCIP_CALL( SCIPcheckCurvatureQuadratic(scip, cons) );
1777  if( ! presoldata->updatequadindef && ( SCIPisConvexQuadratic(scip, cons) || SCIPisConcaveQuadratic(scip, cons) ) )
1778  {
1779  SCIPdebugMsg(scip, "quadratic constraint update failed, since matrix associated to quadratic constraint <%s> is not \
1780  indefinite.\n", SCIPconsGetName(cons) );
1781  return SCIP_OKAY;
1782  }
1783 
1784  /* first, check whether the problem is equivalent to
1785  *
1786  * min z
1787  * s.t. x^T Q x + c^T x <= b + z
1788  * x \in \{0, 1\}^{p} \times R^{n-p}.
1789  *
1790  */
1791  SCIP_CALL( checkConsQuadraticProblem(scip, quadconshdlr, cons, presoldata->addkktbinary, &objvar, &scale, &objrhs, &isqp) );
1792  if( ! isqp )
1793  return SCIP_OKAY;
1794  assert( objvar != NULL );
1795 
1796  /* get constraint handler data of linear constraints */
1797  linconshdlr = SCIPfindConshdlr(scip, "linear");
1798 
1799  /* get linear constraints and number of linear constraints */
1800  if( linconshdlr != NULL )
1801  {
1802  nlinconss = SCIPconshdlrGetNConss(linconshdlr);
1803  linconss = SCIPconshdlrGetConss(linconshdlr);
1804  }
1805 
1806  /* get variables that are in the linear term of the quadratic constraint */
1807  nlintermvars = SCIPgetNLinearVarsQuadratic(scip, cons);
1808  lintermvars = SCIPgetLinearVarsQuadratic(scip, cons);
1809  lintermcoefs = SCIPgetCoefsLinearVarsQuadratic(scip, cons);
1810 
1811  /* get bilinear terms */
1812  bilinterms = SCIPgetBilinTermsQuadratic(scip, cons);
1813  nbilinterms = SCIPgetNBilinTermsQuadratic(scip, cons);
1814 
1815  /* get quadratic terms */
1816  quadterms = SCIPgetQuadVarTermsQuadratic(scip, cons);
1817  nquadterms = SCIPgetNQuadVarTermsQuadratic(scip, cons);
1818 
1819  /* the update is only valid if a finite optimal solution of the problem exists,
1820  * since only finite optimal solutions satisfy the KKT conditions;
1821  * we check whether all variables have finite bounds, otherwise we return */
1822  if( presoldata->updatequadbounded )
1823  {
1824  /* check linear term variables */
1825  for( j = 0; j < nlintermvars; ++j )
1826  {
1827  SCIP_VAR* var;
1828 
1829  var = lintermvars[j];
1830  if( var != objvar &&
1831  ( SCIPisInfinity(scip, -SCIPvarGetLbGlobal(var)) || SCIPisInfinity(scip, SCIPvarGetUbGlobal(var)) )
1832  )
1833  {
1834  SCIPdebugMsg(scip, "failed adding the KKT conditions, since not all variables to quadratic constraint <%s> are \
1835  bounded.\n", SCIPconsGetName(cons) );
1836  return SCIP_OKAY;
1837  }
1838  }
1839 
1840  /* check linear term variables */
1841  for( j = 0; j < nbilinterms; ++j )
1842  {
1843  SCIP_VAR* bilvar1;
1844  SCIP_VAR* bilvar2;
1845 
1846  bilvar1 = bilinterms[j].var1;
1847  bilvar2 = bilinterms[j].var2;
1848  if( ( bilvar1 != objvar && ( SCIPisInfinity(scip, -SCIPvarGetLbGlobal(bilvar1)) || SCIPisInfinity(scip, SCIPvarGetUbGlobal(bilvar1)) ) )
1849  || ( bilvar2 != objvar && ( SCIPisInfinity(scip, -SCIPvarGetLbGlobal(bilvar2)) || SCIPisInfinity(scip, SCIPvarGetUbGlobal(bilvar2)) ) ) )
1850  {
1851  SCIPdebugMsg(scip, "failed adding the KKT conditions, since not all variables to quadratic constraint <%s> \
1852  are bounded.\n", SCIPconsGetName(cons) );
1853  return SCIP_OKAY;
1854  }
1855  }
1856 
1857  /* check quadratic term variables */
1858  for( j = 0; j < nquadterms; ++j )
1859  {
1860  SCIP_VAR* var;
1861 
1862  var = quadterms[j].var;
1863  if( var != objvar
1864  && ( SCIPisInfinity(scip, -SCIPvarGetLbGlobal(var)) || SCIPisInfinity(scip, SCIPvarGetUbGlobal(var)) )
1865  )
1866  {
1867  SCIPdebugMsg(scip, "failed adding the KKT conditions, since not all variables to quadratic constraint <%s> \
1868  are bounded.\n", SCIPconsGetName(cons) );
1869  return SCIP_OKAY;
1870  }
1871  }
1872  }
1873 
1874 
1875  /* add KKT constraints */
1876 
1877  /* set up hash map */
1878  SCIP_CALL( SCIPhashmapCreate(&varhash, SCIPblkmem(scip), SCIPgetNVars(scip) + SCIPgetNFixedVars(scip)) );
1879 
1880  /* allocate buffer array */
1881  SCIP_CALL( SCIPallocBufferArray(scip, &dualconss, 2 * SCIPgetNVars(scip) + 2 * SCIPgetNFixedVars(scip)) ); /*lint !e647*/
1882 
1883  /* duplicate linconss for later use, since in the following, we create new linear constraints */
1884  if( linconss != NULL )
1885  {
1886  SCIP_CALL( SCIPduplicateBufferArray(scip, &savelinconss, linconss, nlinconss) );
1887  }
1888 
1889  /* create new objective constraint */
1890  SCIP_CALL( SCIPcreateConsBasicLinear(scip, &objcons, "objcons", 0, NULL, NULL, objrhs, objrhs) );
1891  if( SCIPisFeasNegative(scip, SCIPvarGetObj(objvar)) )
1892  {
1893  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, objvar, 1.0) );
1894  }
1895  else
1896  {
1897  SCIP_CALL( SCIPaddCoefLinear(scip, objcons, objvar, -1.0) );
1898  }
1899 
1900  /* handle linear constraints */
1901  if( savelinconss != NULL )
1902  {
1903  SCIP_CALL( presolveAddKKTLinearConss(scip, objcons, savelinconss, nlinconss, varhash, dualconss, &ndualconss,
1904  naddconss, ndelconss) );
1905  }
1906 
1907  /* handle set packing constraints */
1908  SCIP_CALL( presolveAddKKTSetppcConss(scip, objcons, varhash, dualconss, &ndualconss, naddconss, ndelconss) );
1909 
1910  /* handle knapsack constraints */
1911  SCIP_CALL( presolveAddKKTKnapsackConss(scip, objcons, varhash, dualconss, &ndualconss, naddconss, ndelconss) );
1912 
1913  /* handle varbound constraints */
1914  SCIP_CALL( presolveAddKKTVarboundConss(scip, objcons, varhash, dualconss, &ndualconss, naddconss, ndelconss) );
1915 
1916  /* handle logicor constraints */
1917  SCIP_CALL( presolveAddKKTLogicorConss(scip, objcons, varhash, dualconss, &ndualconss, naddconss, ndelconss) );
1918 
1919  /* handle linear constraints associated to aggregations of variables */
1920  if( SCIPgetNFixedVars(scip) > 0 )
1921  {
1923  varhash, dualconss, &ndualconss, naddconss) );
1924  }
1925 
1926  /* handle bilinear terms of quadratic constraint */
1927  SCIP_CALL( presolveAddKKTQuadBilinearTerms(scip, objcons, bilinterms, nbilinterms, varhash, scale, dualconss,
1928  &ndualconss, naddconss) );
1929 
1930  /* handle quadratic terms of quadratic constraint */
1931  SCIP_CALL( presolveAddKKTQuadQuadraticTerms(scip, objcons, quadterms, nquadterms, varhash, scale, dualconss,
1932  &ndualconss, naddconss) );
1933 
1934  /* handle linear terms of quadratic constraint */
1935  SCIP_CALL( presolveAddKKTQuadLinearTerms(scip, objcons, lintermvars, lintermcoefs, nlintermvars, quadterms, nquadterms,
1936  varhash, objvar, scale, dualconss, &ndualconss, naddconss) );
1937 
1938 
1939  /* add/release objective constraint */
1940  SCIP_CALL( SCIPaddCons(scip, objcons) );
1941  SCIP_CALL( SCIPreleaseCons(scip, &objcons) );
1942  ++(*naddconss);
1943 
1944  /* add/release dual constraints associated to the KKT conditions */
1945  for( j = 0; j < ndualconss; ++j )
1946  {
1947  SCIP_CALL( SCIPaddCons(scip, dualconss[j]) );
1948  SCIP_CALL( SCIPreleaseCons(scip, &dualconss[j]) );
1949  }
1950  *naddconss = *naddconss + ndualconss;
1951 
1952  /* free buffer array */
1953  SCIPfreeBufferArrayNull(scip, &savelinconss);
1954  SCIPfreeBufferArray(scip, &dualconss);
1955 
1956  /* free hash map */
1957  SCIPhashmapFree(&varhash);
1958 
1959  if( SCIPgetNBinVars(scip) > 0 )
1960  SCIPdebugMsg(scip, "added the KKT conditions to the mixed-binary quadratic program\n");
1961  else
1962  SCIPdebugMsg(scip, "added the KKT conditions to the quadratic program\n");
1963 
1964  /*SCIP_CALL( SCIPwriteTransProblem(scip, "trafoQP.lp", NULL, FALSE ) );*/
1965 
1966  return SCIP_OKAY;
1967 }
1968 
1969 
1970 /*
1971  * presolver specific interface methods
1972  */
1973 
1974 /** creates the QP KKT reformulation presolver and includes it in SCIP */
1976  SCIP* scip /**< SCIP data structure */
1977  )
1978 {
1979  SCIP_PRESOLDATA* presoldata;
1980  SCIP_PRESOL* presol= NULL;
1981 
1982  /* alloc presolve data object */
1983  SCIP_CALL( SCIPallocBlockMemory(scip, &presoldata) );
1984 
1985  /* include presolver */
1987  PRESOL_TIMING, presolExecQPKKTref, presoldata) );
1988  assert(presol != NULL);
1989 
1990  /* set non fundamental callbacks via setter functions */
1991  SCIP_CALL( SCIPsetPresolCopy(scip, presol, presolCopyQPKKTref) );
1992  SCIP_CALL( SCIPsetPresolFree(scip, presol, presolFreeQPKKTref) );
1993 
1994  /* add qpkktref presolver parameters */
1995  SCIP_CALL( SCIPaddBoolParam(scip, "presolving/" PRESOL_NAME "/addkktbinary",
1996  "if TRUE then allow binary variables for KKT update",
1997  &presoldata->addkktbinary, TRUE, FALSE, NULL, NULL) );
1998 
1999  SCIP_CALL( SCIPaddBoolParam(scip, "presolving/" PRESOL_NAME "/updatequadbounded",
2000  "if TRUE then only apply the update to QPs with bounded variables; if the variables are not bounded then a \
2001  finite optimal solution might not exist and the KKT conditions would then be invalid",
2002  &presoldata->updatequadbounded, TRUE, TRUE, NULL, NULL) );
2003 
2004  SCIP_CALL( SCIPaddBoolParam(scip, "presolving/" PRESOL_NAME "/updatequadindef",
2005  "if TRUE then apply quadratic constraint update even if the quadratic constraint matrix is known to be indefinite",
2006  &presoldata->updatequadindef, TRUE, FALSE, NULL, NULL) );
2007 
2008  return SCIP_OKAY;
2009 }
SCIP_VAR ** SCIPgetLinearVarsQuadratic(SCIP *scip, SCIP_CONS *cons)
static SCIP_RETCODE createKKTComplementarityLinear(SCIP *scip, const char *namepart, SCIP_VAR **vars, SCIP_Real *vals, SCIP_Real lhs, SCIP_Real rhs, int nvars, SCIP_VAR *dualvar, SCIP_Bool takelhs, int *naddconss)
SCIP_RETCODE SCIPincludePresolBasic(SCIP *scip, SCIP_PRESOL **presolptr, const char *name, const char *desc, int priority, int maxrounds, SCIP_PRESOLTIMING timing, SCIP_DECL_PRESOLEXEC((*presolexec)), SCIP_PRESOLDATA *presoldata)
Definition: scip.c:6854
int SCIPgetNIntVars(SCIP *scip)
Definition: scip.c:11721
SCIP_Bool SCIPisFeasZero(SCIP *scip, SCIP_Real val)
Definition: scip.c:46151
struct SCIP_PresolData SCIP_PRESOLDATA
Definition: type_presol.h:37
static SCIP_RETCODE presolveAddKKTQuadBilinearTerms(SCIP *scip, SCIP_CONS *objcons, SCIP_BILINTERM *bilinterms, int nbilinterms, SCIP_HASHMAP *varhash, SCIP_Real scale, SCIP_CONS **dualconss, int *ndualconss, int *naddconss)
SCIP_VAR * var2
SCIP_RETCODE SCIPsetPresolFree(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLFREE((*presolfree)))
Definition: scip.c:6905
SCIP_Bool SCIPisFeasEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip.c:46086
static SCIP_DECL_PRESOLCOPY(presolCopyQPKKTref)
SCIP_Real * SCIPvarGetMultaggrScalars(SCIP_VAR *var)
Definition: var.c:16958
SCIP_RETCODE SCIPcreateConsBasicLinear(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_Real *vals, SCIP_Real lhs, SCIP_Real rhs)
Constraint handler for variable bound constraints .
SCIP_CONSHDLR * SCIPfindConshdlr(SCIP *scip, const char *name)
Definition: scip.c:6541
static SCIP_RETCODE presolveAddKKTLogicorConss(SCIP *scip, SCIP_CONS *objcons, SCIP_HASHMAP *varhash, SCIP_CONS **dualconss, int *ndualconss, int *naddconss, int *ndelconss)
int SCIPgetNVarsSetppc(SCIP *scip, SCIP_CONS *cons)
Definition: cons_setppc.c:9172
SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:17166
int SCIPgetNVarsLogicor(SCIP *scip, SCIP_CONS *cons)
static SCIP_RETCODE presolveAddKKTKnapsackConss(SCIP *scip, SCIP_CONS *objcons, SCIP_HASHMAP *varhash, SCIP_CONS **dualconss, int *ndualconss, int *naddconss, int *ndelconss)
#define SCIP_MAXSTRLEN
Definition: def.h:215
SCIP_VAR * var1
SCIP_RETCODE SCIPdelCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip.c:12481
SCIP_VAR ** SCIPvarGetMultaggrVars(SCIP_VAR *var)
Definition: var.c:16946
#define PRESOL_DESC
SCIP_RETCODE SCIPreleaseVar(SCIP *scip, SCIP_VAR **var)
Definition: scip.c:18384
SCIP_Bool SCIPvarIsBinary(SCIP_VAR *var)
Definition: var.c:16732
SCIP_Bool SCIPisFeasNegative(SCIP *scip, SCIP_Real val)
Definition: scip.c:46175
static SCIP_RETCODE createKKTComplementarityBounds(SCIP *scip, SCIP_VAR *var, SCIP_VAR *dualvar, SCIP_Bool takelb, int *naddconss)
SCIP_CONS ** SCIPconshdlrGetConss(SCIP_CONSHDLR *conshdlr)
Definition: cons.c:4485
#define FALSE
Definition: def.h:64
SCIP_RETCODE SCIPhashmapCreate(SCIP_HASHMAP **hashmap, BMS_BLKMEM *blkmem, int mapsize)
Definition: misc.c:2765
SCIP_Real SCIPinfinity(SCIP *scip)
Definition: scip.c:45816
int SCIPsnprintf(char *t, int len, const char *s,...)
Definition: misc.c:9340
#define TRUE
Definition: def.h:63
#define PRESOL_NAME
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:53
SCIP_RETCODE SCIPincludePresolQPKKTref(SCIP *scip)
SCIP_Real SCIPvarGetNegationConstant(SCIP_VAR *var)
Definition: var.c:17003
SCIP_PRESOLDATA * SCIPpresolGetData(SCIP_PRESOL *presol)
Definition: presol.c:466
SCIP_RETCODE SCIPcreateVarBasic(SCIP *scip, SCIP_VAR **var, const char *name, SCIP_Real lb, SCIP_Real ub, SCIP_Real obj, SCIP_VARTYPE vartype)
Definition: scip.c:17317
static SCIP_RETCODE presolveAddKKTLinearCons(SCIP *scip, SCIP_CONS *objcons, const char *namepart, SCIP_VAR **vars, SCIP_Real *vals, SCIP_Real lhs, SCIP_Real rhs, int nvars, SCIP_HASHMAP *varhash, SCIP_CONS **dualconss, int *ndualconss, int *naddconss)
SCIP_Real SCIPvarGetAggrScalar(SCIP_VAR *var)
Definition: var.c:16912
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip.h:21907
SCIP_VAR ** SCIPgetVarsKnapsack(SCIP *scip, SCIP_CONS *cons)
#define PRESOL_PRIORITY
#define SCIPduplicateBufferArray(scip, ptr, source, num)
Definition: scip.h:21933
void * SCIPhashmapGetImage(SCIP_HASHMAP *hashmap, void *origin)
Definition: misc.c:2903
static SCIP_RETCODE createKKTDualCons(SCIP *scip, SCIP_CONS *objcons, SCIP_VAR *var, SCIP_HASHMAP *varhash, SCIP_CONS **dualconss, int *ndualconss, SCIP_CONS **dualcons, int *naddconss)
#define SCIPfreeBufferArray(scip, ptr)
Definition: scip.h:21937
Constraint handler for the set partitioning / packing / covering constraints .
#define SCIPallocBlockMemory(scip, ptr)
Definition: scip.h:21890
SCIP_VAR * SCIPvarGetNegationVar(SCIP_VAR *var)
Definition: var.c:16992
#define SCIPdebugMsg
Definition: scip.h:451
SCIP_Real SCIPgetRhsLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_Real SCIPgetLhsQuadratic(SCIP *scip, SCIP_CONS *cons)
SCIP_RETCODE SCIPaddCoefLinear(SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real val)
int SCIPgetNQuadVarTermsQuadratic(SCIP *scip, SCIP_CONS *cons)
SCIP_Bool SCIPhashmapExists(SCIP_HASHMAP *hashmap, void *origin)
Definition: misc.c:2997
int SCIPgetNBilinTermsQuadratic(SCIP *scip, SCIP_CONS *cons)
int SCIPgetNFixedVars(SCIP *scip)
Definition: scip.c:11948
SCIP_VAR ** SCIPgetFixedVars(SCIP *scip)
Definition: scip.c:11905
SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:17176
static SCIP_DECL_PRESOLEXEC(presolExecQPKKTref)
SCIP_Bool SCIPisConcaveQuadratic(SCIP *scip, SCIP_CONS *cons)
Constraint handler for knapsack constraints of the form , x binary and .
static SCIP_RETCODE presolveAddKKTQuadLinearTerms(SCIP *scip, SCIP_CONS *objcons, SCIP_VAR **lintermvars, SCIP_Real *lintermcoefs, int nlintermvars, SCIP_QUADVARTERM *quadterms, int nquadterms, SCIP_HASHMAP *varhash, SCIP_VAR *objvar, SCIP_Real scale, SCIP_CONS **dualconss, int *ndualconss, int *naddconss)
SCIP_VAR * SCIPgetVarVarbound(SCIP *scip, SCIP_CONS *cons)
#define SCIPerrorMessage
Definition: pub_message.h:45
SCIP_RETCODE SCIPaddCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip.c:12410
static SCIP_RETCODE presolveAddKKTAggregatedVars(SCIP *scip, SCIP_CONS *objcons, SCIP_VAR **agrvars, int nagrvars, SCIP_HASHMAP *varhash, SCIP_CONS **dualconss, int *ndualconss, int *naddconss)
SCIP_Real SCIPgetRhsVarbound(SCIP *scip, SCIP_CONS *cons)
SCIP_RETCODE SCIPaddVarSOS1(SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real weight)
Definition: cons_sos1.c:10459
Constraint handler for logicor constraints (equivalent to set covering, but algorithms are suited fo...
static SCIP_RETCODE presolveAddKKTQuadQuadraticTerms(SCIP *scip, SCIP_CONS *objcons, SCIP_QUADVARTERM *quadterms, int nquadterms, SCIP_HASHMAP *varhash, SCIP_Real scale, SCIP_CONS **dualconss, int *ndualconss, int *naddconss)
void SCIPpresolSetData(SCIP_PRESOL *presol, SCIP_PRESOLDATA *presoldata)
Definition: presol.c:476
SCIP_Real * SCIPgetCoefsLinearVarsQuadratic(SCIP *scip, SCIP_CONS *cons)
SCIP_Real SCIPgetVbdcoefVarbound(SCIP *scip, SCIP_CONS *cons)
#define SCIPfreeBufferArrayNull(scip, ptr)
Definition: scip.h:21938
BMS_BLKMEM * SCIPblkmem(SCIP *scip)
Definition: scip.c:45519
int SCIPgetLinvarMayDecreaseQuadratic(SCIP *scip, SCIP_CONS *cons)
const char * SCIPconsGetName(SCIP_CONS *cons)
Definition: cons.c:7881
SCIP_Real SCIPvarGetAggrConstant(SCIP_VAR *var)
Definition: var.c:16923
SCIP_VAR ** SCIPgetVarsLogicor(SCIP *scip, SCIP_CONS *cons)
const char * SCIPvarGetName(SCIP_VAR *var)
Definition: var.c:16552
void SCIPhashmapFree(SCIP_HASHMAP **hashmap)
Definition: misc.c:2798
constraint handler for quadratic constraints
#define NULL
Definition: lpi_spx1.cpp:137
#define REALABS(x)
Definition: def.h:159
SCIP_VAR * SCIPgetVbdvarVarbound(SCIP *scip, SCIP_CONS *cons)
#define SCIP_CALL(x)
Definition: def.h:306
SCIP_Real SCIPvarGetMultaggrConstant(SCIP_VAR *var)
Definition: var.c:16970
static SCIP_RETCODE presolveAddKKTLinearConss(SCIP *scip, SCIP_CONS *objcons, SCIP_CONS **savelinconss, int nlinconss, SCIP_HASHMAP *varhash, SCIP_CONS **dualconss, int *ndualconss, int *naddconss, int *ndelconss)
SCIP_BILINTERM * SCIPgetBilinTermsQuadratic(SCIP *scip, SCIP_CONS *cons)
int SCIPconshdlrGetNConss(SCIP_CONSHDLR *conshdlr)
Definition: cons.c:4515
SCIP_Longint SCIPgetCapacityKnapsack(SCIP *scip, SCIP_CONS *cons)
SCIP_Real SCIPgetRhsQuadratic(SCIP *scip, SCIP_CONS *cons)
#define SCIPallocBufferArray(scip, ptr, num)
Definition: scip.h:21925
#define SCIP_Bool
Definition: def.h:61
static SCIP_RETCODE checkConsQuadraticProblem(SCIP *scip, SCIP_CONSHDLR *quadconshdlr, SCIP_CONS *cons, SCIP_Bool allowbinary, SCIP_VAR **objvar, SCIP_Real *scale, SCIP_Real *objrhs, SCIP_Bool *isqp)
static SCIP_RETCODE createKKTComplementarityBinary(SCIP *scip, SCIP_VAR *var, SCIP_VAR *dualbin1, SCIP_VAR *dualbin2, int *naddconss)
SCIP_SETPPCTYPE SCIPgetTypeSetppc(SCIP *scip, SCIP_CONS *cons)
Definition: cons_setppc.c:9214
int SCIPvarGetNLocksUp(SCIP_VAR *var)
Definition: var.c:3217
static SCIP_DECL_PRESOLFREE(presolFreeQPKKTref)
const char * SCIPpresolGetName(SCIP_PRESOL *presol)
Definition: presol.c:553
SCIP_Real SCIPvarGetObj(SCIP_VAR *var)
Definition: var.c:17014
SCIP_VAR * SCIPvarGetAggrVar(SCIP_VAR *var)
Definition: var.c:16901
static SCIP_RETCODE presolveAddKKTSetppcConss(SCIP *scip, SCIP_CONS *objcons, SCIP_HASHMAP *varhash, SCIP_CONS **dualconss, int *ndualconss, int *naddconss, int *ndelconss)
Constraint handler for linear constraints in their most general form, .
int SCIPgetNObjVars(SCIP *scip)
Definition: scip.c:11859
SCIP_RETCODE SCIPchgRhsLinear(SCIP *scip, SCIP_CONS *cons, SCIP_Real rhs)
int SCIPvarGetMultaggrNVars(SCIP_VAR *var)
Definition: var.c:16934
SCIP_Bool SCIPisInfinity(SCIP *scip, SCIP_Real val)
Definition: scip.c:45827
int SCIPgetNBinVars(SCIP *scip)
Definition: scip.c:11676
SCIP_RETCODE SCIPcreateConsBasicSOS1(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_Real *weights)
Definition: cons_sos1.c:10443
SCIP_VAR ** SCIPgetVarsSetppc(SCIP *scip, SCIP_CONS *cons)
Definition: cons_setppc.c:9193
int SCIPgetNVars(SCIP *scip)
Definition: scip.c:11631
SCIP_RETCODE SCIPcheckCurvatureQuadratic(SCIP *scip, SCIP_CONS *cons)
int SCIPgetLinvarMayIncreaseQuadratic(SCIP *scip, SCIP_CONS *cons)
int SCIPgetNLinearVarsQuadratic(SCIP *scip, SCIP_CONS *cons)
int SCIPvarGetNLocksDown(SCIP_VAR *var)
Definition: var.c:3162
enum SCIP_SetppcType SCIP_SETPPCTYPE
Definition: cons_setppc.h:77
static const SCIP_Real scalars[]
Definition: lp.c:5563
SCIP_RETCODE SCIPaddVar(SCIP *scip, SCIP_VAR *var)
Definition: scip.c:11311
SCIP_VAR ** SCIPgetVarsLinear(SCIP *scip, SCIP_CONS *cons)
#define PRESOL_MAXROUNDS
int SCIPgetNConss(SCIP *scip)
Definition: scip.c:12679
SCIP_RETCODE SCIPreleaseCons(SCIP *scip, SCIP_CONS **cons)
Definition: scip.c:27323
SCIP_RETCODE SCIPsetPresolCopy(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLCOPY((*presolcopy)))
Definition: scip.c:6889
SCIP_Bool SCIPisFeasPositive(SCIP *scip, SCIP_Real val)
Definition: scip.c:46163
SCIP_VARSTATUS SCIPvarGetStatus(SCIP_VAR *var)
Definition: var.c:16671
#define SCIP_Real
Definition: def.h:135
constraint handler for SOS type 1 constraints
int SCIPgetNVarsKnapsack(SCIP *scip, SCIP_CONS *cons)
#define SCIP_Longint
Definition: def.h:120
SCIP_Real * SCIPgetValsLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_RETCODE SCIPchgLhsLinear(SCIP *scip, SCIP_CONS *cons, SCIP_Real lhs)
SCIP_RETCODE SCIPhashmapInsert(SCIP_HASHMAP *hashmap, void *origin, void *image)
Definition: misc.c:2846
SCIP_QUADVARTERM * SCIPgetQuadVarTermsQuadratic(SCIP *scip, SCIP_CONS *cons)
#define PRESOL_TIMING
SCIP_Longint * SCIPgetWeightsKnapsack(SCIP *scip, SCIP_CONS *cons)
SCIP_Real SCIPgetLhsVarbound(SCIP *scip, SCIP_CONS *cons)
int SCIPgetNVarsLinear(SCIP *scip, SCIP_CONS *cons)
static SCIP_RETCODE presolveAddKKTVarboundConss(SCIP *scip, SCIP_CONS *objcons, SCIP_HASHMAP *varhash, SCIP_CONS **dualconss, int *ndualconss, int *naddconss, int *ndelconss)
SCIP_Real SCIPgetLhsLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_Bool SCIPisConvexQuadratic(SCIP *scip, SCIP_CONS *cons)
SCIP_RETCODE SCIPaddBoolParam(SCIP *scip, const char *name, const char *desc, SCIP_Bool *valueptr, SCIP_Bool isadvanced, SCIP_Bool defaultvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip.c:4176
qpkktref presolver