SCIP-SDP  3.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
cons_sdp.c
Go to the documentation of this file.
1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2 /* */
3 /* This file is part of SCIPSDP - a solving framework for mixed-integer */
4 /* semidefinite programs based on SCIP. */
5 /* */
6 /* Copyright (C) 2011-2013 Discrete Optimization, TU Darmstadt */
7 /* EDOM, FAU Erlangen-Nürnberg */
8 /* 2014-2017 Discrete Optimization, TU Darmstadt */
9 /* */
10 /* */
11 /* This program is free software; you can redistribute it and/or */
12 /* modify it under the terms of the GNU Lesser General Public License */
13 /* as published by the Free Software Foundation; either version 3 */
14 /* of the License, or (at your option) any later version. */
15 /* */
16 /* This program is distributed in the hope that it will be useful, */
17 /* but WITHOUT ANY WARRANTY; without even the implied warranty of */
18 /* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the */
19 /* GNU Lesser General Public License for more details. */
20 /* */
21 /* You should have received a copy of the GNU Lesser General Public License */
22 /* along with this program; if not, write to the Free Software */
23 /* Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301, USA.*/
24 /* */
25 /* */
26 /* Based on SCIP - Solving Constraint Integer Programs */
27 /* Copyright (C) 2002-2017 Zuse Institute Berlin */
28 /* SCIP is distributed under the terms of the SCIP Academic Licence, */
29 /* see file COPYING in the SCIP distribution. */
30 /* */
31 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
32 
44 /* #define SCIP_DEBUG*/
45 /* #define SCIP_MORE_DEBUG *//* shows all cuts added and prints constraint after parsing */
46 /* #define PRINT_HUMAN_READABLE *//* change the output of PRINTCONS to a better readable format (dense instead of sparse), WHICH CAN NO LONGER BE PARSED */
47 
48 #include "cons_sdp.h"
49 
50 #include <assert.h> /*lint !e451*/
51 #include <string.h> /* for NULL, strcmp */
52 #include <ctype.h> /* for isspace */
53 
54 #include "sdpi/lapack.h"
55 
56 #include "scipsdp/SdpVarmapper.h"
57 #include "scipsdp/SdpVarfixer.h"
58 
59 #include "scip/cons_linear.h" /* for SCIPcreateConsLinear */
60 #include "scip/scip.h" /* for SCIPallocBufferArray, etc */
61 
62 #ifdef OMP
63 #include "omp.h" /* for changing the number of threads */
64 #endif
65 
66 /* turn off lint warnings for whole file: */
67 /*lint --e{788,818}*/
68 
69 #define CONSHDLR_NAME "SDP"
70 #define CONSHDLR_DESC "SDP constraints of the form \\sum_{j} A_j y_j - A_0 psd"
71 #define CONSHDLR_SEPAPRIORITY +1000000
72 #define CONSHDLR_ENFOPRIORITY -2000000
73 #define CONSHDLR_CHECKPRIORITY -2000000
74 #define CONSHDLR_SEPAFREQ 1
75 #define CONSHDLR_EAGERFREQ 100
77 #define CONSHDLR_MAXPREROUNDS -1
78 #define CONSHDLR_DELAYSEPA FALSE
79 #define CONSHDLR_NEEDSCONS TRUE
81 #define CONSHDLR_PRESOLTIMING SCIP_PRESOLTIMING_FAST
82 #define PARSE_STARTSIZE 1
83 #define PARSE_SIZEFACTOR 10
84 #define DEFAULT_DIAGGEZEROCUTS TRUE
85 #define DEFAULT_DIAGZEROIMPLCUTS TRUE
86 #ifdef OMP
87 #define DEFAULT_NTHREADS 1
88 #endif
89 
91 struct SCIP_ConsData
92 {
93  int nvars;
94  int nnonz;
95  int blocksize;
96  int* nvarnonz;
97  int** col;
98  int** row;
99  SCIP_Real** val;
100  SCIP_VAR** vars;
101  int constnnonz;
102  int* constcol;
103  int* constrow;
104  SCIP_Real* constval;
105  SCIP_Real maxrhsentry;
106 };
107 
109 struct SCIP_ConshdlrData
110 {
111  int neigveccuts;
112  SCIP_Bool diaggezerocuts;
113  int ndiaggezerocuts;
114  int n1x1blocks;
115  SCIP_Bool diagzeroimplcuts;
116 #ifdef OMP
117  int nthreads;
118 #endif
119 };
120 
122 static
123 SCIP_RETCODE expandSymMatrix(
124  int size,
125  SCIP_Real* symMat,
126  SCIP_Real* fullMat
127  )
128 {
129  int i;
130  int j;
131  int ind = 0;
132 
133  assert( size >= 0 );
134  assert( symMat != NULL );
135  assert( fullMat != NULL );
136 
137  /* traverse the lower triangular part in the order of the indices and copy the values to both lower and upper triangular part */
138  for (i = 0; i < size; i++)
139  {
140  for (j = 0; j <= i; j++)
141  {
142  assert( ind == SCIPconsSdpCompLowerTriangPos(i,j) );
143  fullMat[i*size + j] = symMat[ind]; /*lint !e679*/
144  fullMat[j*size + i] = symMat[ind]; /*lint !e679*/
145  ind++;
146  }
147  }
148 
149  return SCIP_OKAY;
150 }
151 
155 static
156 SCIP_RETCODE computeSdpMatrix(
157  SCIP* scip,
158  SCIP_CONS* cons,
159  SCIP_SOL* y,
160  SCIP_Real* matrix
161  )
162 {
163  SCIP_CONSDATA* consdata;
164  SCIP_Real yval;
165  int blocksize;
166  int nvars;
167  int ind;
168  int i;
169 
170  assert( cons != NULL );
171  assert( matrix != NULL );
172 
173  consdata = SCIPconsGetData(cons);
174  nvars = consdata->nvars;
175  blocksize = consdata->blocksize;
176 
177  /* initialize the matrix with 0 */
178  for (i = 0; i < (blocksize * (blocksize + 1))/2; i++)
179  matrix[i] = 0.0;
180 
181  /* add the non-constant-part */
182  for (i = 0; i < nvars; i++)
183  {
184  yval = SCIPgetSolVal(scip, y, consdata->vars[i]);
185  for (ind = 0; ind < consdata->nvarnonz[i]; ind++)
186  matrix[SCIPconsSdpCompLowerTriangPos(consdata->row[i][ind], consdata->col[i][ind])] += yval * consdata->val[i][ind];
187  }
188 
189  /* substract the constant part */
190  for (ind = 0; ind < consdata->constnnonz; ind++)
191  matrix[SCIPconsSdpCompLowerTriangPos(consdata->constrow[ind], consdata->constcol[ind])] -= consdata->constval[ind];
192 
193  return SCIP_OKAY;
194 }
195 
197 static
198 SCIP_RETCODE multiplyConstraintMatrix(
199  SCIP_CONS* cons,
200  int j,
201  SCIP_Real* v,
202  SCIP_Real* vAv
203  )
204 {
205  SCIP_CONSDATA* consdata;
206  int i;
207 
208  assert( cons != NULL );
209  assert( j >= 0 );
210  assert( vAv != NULL );
211 
212  consdata = SCIPconsGetData(cons);
213  assert( consdata != NULL );
214 
215  assert( j < consdata->nvars );
216 
217  /* initialize the product with 0 */
218  *vAv = 0.0;
219 
220  for (i = 0; i < consdata->nvarnonz[j]; i++)
221  {
222  if (consdata->col[j][i] == consdata->row[j][i])
223  *vAv += v[consdata->col[j][i]] * consdata->val[j][i] * v[consdata->row[j][i]];
224  else
225  {
226  /* Multiply by 2, because the matrix is symmetric and there is one identical contribution each from lower and upper triangular part. */
227  *vAv += 2.0 * v[consdata->col[j][i]] * consdata->val[j][i] * v[consdata->row[j][i]];
228  }
229  }
230 
231  return SCIP_OKAY;
232 }
233 
237 static
238 SCIP_RETCODE setMaxRhsEntry(
239  SCIP_CONS* cons
240  )
241 {
242  SCIP_CONSDATA* consdata;
243  SCIP_Real max;
244  int i;
245 
246  consdata = SCIPconsGetData(cons);
247  assert( consdata != NULL );
248 
249  /* initialize max with zero (this is used if there is no constant-matrix) */
250  max = 0.0;
251 
252  /* iterate over the entries of the constant matrix, updating max if a higher absolute value is found */
253  for (i = 0; i < consdata->constnnonz; i++)
254  {
255  if ( REALABS(consdata->constval[i]) > max )
256  max = REALABS(consdata->constval[i]);
257  }
258 
259  consdata->maxrhsentry = max;
260 
261  return SCIP_OKAY;
262 }
263 
264 
270 static
271 SCIP_RETCODE cutUsingEigenvector(
272  SCIP* scip,
273  SCIP_CONS* cons,
274  SCIP_SOL* sol,
275  SCIP_Real* coeff,
276  SCIP_Real* lhs
277  )
278 {
279  SCIP_CONSDATA* consdata;
280  SCIP_Real* eigenvalues = NULL;
281  SCIP_Real* matrix = NULL;
282  SCIP_Real* fullmatrix = NULL;
283  SCIP_Real* fullconstmatrix = NULL;
284  SCIP_Real* eigenvector = NULL;
285  SCIP_Real* output_vector = NULL;
286  int blocksize;
287  int j;
288 
289  assert( cons != NULL );
290  assert( lhs != NULL );
291 
292  *lhs = 0.0;
293 
294  consdata = SCIPconsGetData(cons);
295  assert( consdata != NULL );
296 
297  blocksize = consdata->blocksize;
298 
299  SCIP_CALL( SCIPallocBufferArray(scip, &eigenvalues, blocksize) );
300  SCIP_CALL( SCIPallocBufferArray(scip, &matrix, (blocksize * (blocksize+1))/2 ) ); /*lint !e647*/
301  SCIP_CALL( SCIPallocBufferArray(scip, &fullmatrix, blocksize * blocksize ) ); /*lint !e647*/
302  SCIP_CALL( SCIPallocBufferArray(scip, &fullconstmatrix, blocksize * blocksize) ); /*lint !e647*/
303  SCIP_CALL( SCIPallocBufferArray(scip, &eigenvector, blocksize) );
304  SCIP_CALL( SCIPallocBufferArray(scip, &output_vector, blocksize) );
305 
306  /* compute the matrix \f$ \sum_j A_j y_j - A_0 \f$ */
307  SCIP_CALL( computeSdpMatrix(scip, cons, sol, matrix) );
308 
309  /* expand it because LAPACK wants the full matrix instead of the lower triangular part */
310  SCIP_CALL( expandSymMatrix(blocksize, matrix, fullmatrix) );
311 
312  SCIP_CALL( SCIPlapackComputeIthEigenvalue(SCIPbuffer(scip), TRUE, blocksize, fullmatrix, 1, eigenvalues, eigenvector) );
313 
314  /* get full constant matrix */
315  SCIP_CALL( SCIPconsSdpGetFullConstMatrix(scip, cons, fullconstmatrix) );
316 
317  /* multiply eigenvector with constant matrix to get lhs (after multiplying again with eigenvector from the left) */
318  SCIP_CALL( SCIPlapackMatrixVectorMult(blocksize, blocksize, fullconstmatrix, eigenvector, output_vector) );
319 
320  for (j = 0; j < blocksize; ++j)
321  *lhs += eigenvector[j] * output_vector[j];
322 
323  /* compute \f$ v^T A_j v \f$ for eigenvector v and each matrix \f$ A_j \f$ to get the coefficients of the LP cut */
324  for (j = 0; j < consdata->nvars; ++j)
325  {
326  SCIP_CALL( multiplyConstraintMatrix(cons, j, eigenvector, &coeff[j]) );
327  }
328 
329  SCIPfreeBufferArray(scip, &output_vector);
330  SCIPfreeBufferArray(scip, &eigenvector);
331  SCIPfreeBufferArray(scip, &fullconstmatrix);
332  SCIPfreeBufferArray(scip, &fullmatrix);
333  SCIPfreeBufferArray(scip, &matrix);
334  SCIPfreeBufferArray(scip, &eigenvalues);
335 
336  return SCIP_OKAY;
337 }
338 
340 SCIP_RETCODE SCIPconsSdpCheckSdpCons(
341  SCIP* scip,
342  SCIP_CONS* cons,
343  SCIP_SOL* sol,
344  SCIP_Bool checkintegrality,
345  SCIP_Bool checklprows,
346  SCIP_Bool printreason,
347  SCIP_RESULT* result
348  )
349 { /*lint --e{715}*/
350  SCIP_CONSDATA* consdata;
351  int blocksize;
352  SCIP_Real check_value;
353  SCIP_Real eigenvalue;
354  SCIP_Real* matrix = NULL;
355  SCIP_Real* fullmatrix = NULL;
356 
357  assert( scip != NULL );
358  assert( cons != NULL );
359  assert( result != NULL );
360  assert( strcmp(SCIPconshdlrGetName(SCIPconsGetHdlr(cons)), "SDP") == 0);
361 
362  consdata = SCIPconsGetData(cons);
363  assert( consdata != NULL );
364  blocksize = consdata->blocksize;
365 
366  SCIP_CALL( SCIPallocBufferArray(scip, &matrix, (blocksize * (blocksize+1)) / 2) ); /*lint !e647*/
367  SCIP_CALL( SCIPallocBufferArray(scip, &fullmatrix, blocksize * blocksize) ); /*lint !e647*/
368  SCIP_CALL( computeSdpMatrix(scip, cons, sol, matrix) );
369  SCIP_CALL( expandSymMatrix(blocksize, matrix, fullmatrix) );
370 
371  SCIP_CALL( SCIPlapackComputeIthEigenvalue(SCIPbuffer(scip), FALSE, blocksize, fullmatrix, 1, &eigenvalue, NULL) );
372 
373  /* This enables checking the second DIMACS Error Norm: err=max{0, -lambda_min(x)/(1+maximumentry of rhs)}, in that case it also needs
374  * to be changed in the sdpi (and implemented there first), when checking feasibility of problems where all variables are fixed */
375 #ifdef DIMACS
376  check_value = (-eigenvalue) / (1.0 + consdata->maxrhsentry);
377 #else
378  check_value = -eigenvalue;
379 #endif
380 
381  if ( SCIPisFeasLE(scip, check_value, 0.0) )
382  *result = SCIP_FEASIBLE;
383  else
384  {
385  *result = SCIP_INFEASIBLE;
386 #ifdef SCIP_DEBUG
387  if ( printreason )
388  {
389  SCIP_CALL( SCIPprintCons(scip, cons, NULL) );
390  SCIP_CALL( SCIPprintSol(scip, sol, NULL, FALSE) );
391  SCIPinfoMessage(scip, NULL, "non psd matrix (eigenvalue %f, dimacs error norm = %f).\n", eigenvalue, check_value);
392  }
393 #endif
394  }
395 
396  if ( sol != NULL )
397  SCIPupdateSolConsViolation(scip, sol, -eigenvalue, (-eigenvalue) / (1.0 + consdata->maxrhsentry));
398 
399  SCIPfreeBufferArray(scip, &fullmatrix);
400  SCIPfreeBufferArray(scip, &matrix);
401 
402  return SCIP_OKAY;
403 }
404 
406 static
407 SCIP_RETCODE separateSol(
408  SCIP* scip,
409  SCIP_CONSHDLR* conshdlr,
410  SCIP_CONS* cons,
411  SCIP_SOL* sol,
412  SCIP_RESULT* result
413  )
414 {
415  char cutname[SCIP_MAXSTRLEN];
416  SCIP_CONSDATA* consdata;
417  SCIP_CONSHDLRDATA* conshdlrdata;
418  SCIP_Real lhs = 0.0;
419  SCIP_Real* coeff = NULL;
420  SCIP_COL** cols;
421  SCIP_Real* vals;
422  SCIP_ROW* row;
423  int nvars;
424  int j;
425  int len;
426 #ifndef NDEBUG
427  int snprintfreturn; /* this is used to assert that the SCIP string concatenation works */
428 #endif
429 
430  assert( cons != NULL );
431 
432  consdata = SCIPconsGetData(cons);
433  assert( consdata != NULL );
434 
435  nvars = consdata->nvars;
436  SCIP_CALL( SCIPallocBufferArray(scip, &coeff, nvars ) );
437 
438  SCIP_CALL( cutUsingEigenvector(scip, cons, sol, coeff, &lhs) );
439 
440  SCIP_CALL( SCIPallocBufferArray(scip, &cols, nvars ) );
441  SCIP_CALL( SCIPallocBufferArray(scip, &vals, nvars ) );
442 
443  len = 0;
444  for (j = 0; j < nvars; ++j)
445  {
446  if ( SCIPisZero(scip, coeff[j]) )
447  continue;
448 
449  cols[len] = SCIPvarGetCol(SCIPvarGetProbvar(consdata->vars[j]));
450  vals[len] = coeff[j];
451  ++len;
452  }
453 
454  conshdlrdata = SCIPconshdlrGetData(conshdlr);
455  assert( conshdlrdata != NULL );
456 
457 #ifndef NDEBUG
458  snprintfreturn = SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "sepa_eig_sdp_%d", ++(conshdlrdata->neigveccuts));
459  assert( snprintfreturn < SCIP_MAXSTRLEN ); /* check whether name fit into string */
460 #else
461  (void) SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "sepa_eig_sdp_%d", ++(conshdlrdata->neigveccuts));
462 #endif
463  SCIP_CALL( SCIPcreateRowCons(scip, &row, conshdlr, cutname , len, cols, vals, lhs, SCIPinfinity(scip), FALSE, FALSE, TRUE) );
464 
465  if ( SCIPisCutEfficacious(scip, sol, row) )
466  {
467  SCIP_Bool infeasible;
468 #ifdef SCIP_MORE_DEBUG
469  SCIPinfoMessage(scip, NULL, "Added cut %s: ", cutname);
470  SCIPinfoMessage(scip, NULL, "%f <= ", lhs);
471  for (j = 0; j < len; j++)
472  SCIPinfoMessage(scip, NULL, "+ (%f)*%s", vals[j], SCIPvarGetName(SCIPcolGetVar(cols[j])));
473  SCIPinfoMessage(scip, NULL, "\n");
474 #endif
475 
476  SCIP_CALL( SCIPaddRow(scip, row, FALSE, &infeasible) );
477  if ( infeasible )
478  *result = SCIP_CUTOFF;
479  else
480  *result = SCIP_SEPARATED;
481  SCIP_CALL( SCIPresetConsAge(scip, cons) );
482  }
483 
484  SCIP_CALL( SCIPreleaseRow(scip, &row) );
485 
486  SCIPfreeBufferArray(scip, &vals);
487  SCIPfreeBufferArray(scip, &cols);
488  SCIPfreeBufferArray(scip, &coeff);
489 
490  return SCIP_OKAY;
491 }
492 
496 static
497 SCIP_RETCODE diagGEzero(
498  SCIP* scip,
499  SCIP_CONS** conss,
500  int nconss,
501  int* naddconss
502  )
503 {
504  char cutname[SCIP_MAXSTRLEN];
505  SCIP_CONSDATA* consdata;
506  SCIP_Real* matrix;
507  SCIP_Real* cons_array;
508  SCIP_Real* lhs_array;
509  SCIP_Real rhs;
510  int blocksize;
511  int nvars;
512  int i;
513  int j;
514  int k;
515  int c;
516 #ifndef NDEBUG
517  int snprintfreturn; /* used to check if sdnprintf worked */
518 #endif
519 
520  for (c = 0; c < nconss; ++c)
521  {
522  SCIP_CONSHDLR* conshdlr;
523 
524  conshdlr = SCIPconsGetHdlr(conss[c]);
525  assert( conshdlr != NULL );
526 #ifndef NDEBUG
527  assert( strcmp(SCIPconshdlrGetName(conshdlr), "SDP") == 0);
528 #endif
529 
530  consdata = SCIPconsGetData(conss[c]);
531  assert( consdata != NULL );
532 
533  blocksize = consdata->blocksize;
534  nvars = consdata->nvars;
535  rhs = SCIPinfinity(scip);
536 
537  SCIP_CALL( SCIPallocBufferArray(scip, &matrix, (blocksize * (blocksize + 1)) / 2) ); /*lint !e647*/
538  SCIP_CALL( SCIPconsSdpGetLowerTriangConstMatrix(scip, conss[c], matrix) );
539 
540  SCIP_CALL( SCIPallocBufferArray(scip, &cons_array, blocksize * nvars) ); /*lint !e647*/
541  SCIP_CALL( SCIPallocBufferArray(scip, &lhs_array, blocksize) );
542 
543  /* the lhs is the (k,k)-entry of the constant matrix */
544  for (k = 0; k < blocksize; ++k)
545  lhs_array[k] = matrix[SCIPconsSdpCompLowerTriangPos(k,k)];
546 
547  /* get the (k,k)-entry of every matrix A_j */
548  for (j = 0; j < nvars; ++j)
549  {
550  for (k = 0; k < blocksize; ++k)
551  {
552  /* initialize these with 0 */
553  cons_array[k * nvars + j] = 0.0; /*lint !e679*/
554  }
555 
556  /* go through the nonzeros of A_j and look for diagonal entries */
557  for (i = 0; i < consdata->nvarnonz[j]; i++)
558  {
559  if ( consdata->col[j][i] == consdata->row[j][i] )
560  cons_array[consdata->col[j][i] * nvars + j] = consdata->val[j][i]; /*lint !e679*/
561  }
562  }
563 
564  /* add the LP-cuts to SCIP */
565  for (k = 0; k < blocksize; ++k)
566  {
567  SCIP_CONS* cons;
568  SCIP_CONSHDLRDATA* conshdlrdata;
569 
570  conshdlrdata = SCIPconshdlrGetData(conshdlr);
571 #ifndef NDEBUG
572  snprintfreturn = SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "diag_ge_zero_%d", ++(conshdlrdata->ndiaggezerocuts));
573  assert( snprintfreturn < SCIP_MAXSTRLEN );
574 #else
575  (void) SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "diag_ge_zero_%d", ++(conshdlrdata->ndiaggezerocuts));
576 #endif
577 
578 #ifdef SCIP_MORE_DEBUG
579  SCIPinfoMessage(scip, NULL, "Added lp-constraint %s: ", cutname);
580  SCIPinfoMessage(scip, NULL, "%f <= ", lhs_array[k]);
581  for (i = 0; i < consdata->nvars; i++)
582  {
583  if ( ! SCIPisZero(scip, cons_array[k * consdata->nvars + i]) )
584  SCIPinfoMessage(scip, NULL, "+ (%f)*%s", cons_array[k * consdata->nvars + i], SCIPvarGetName(consdata->vars[i]));
585  }
586  SCIPinfoMessage(scip, NULL, "\n");
587 #endif
588 
589  SCIP_CALL( SCIPcreateConsLinear(scip, &cons, cutname, consdata->nvars, consdata->vars, cons_array + k * consdata->nvars, lhs_array[k], rhs,
590  TRUE, TRUE, FALSE, FALSE, TRUE, FALSE, FALSE, TRUE, TRUE, FALSE) ); /*lint !e679*/
591 
592  SCIP_CALL( SCIPaddCons(scip, cons) );
593  SCIP_CALL( SCIPreleaseCons(scip, &cons) );
594  ++(*naddconss);
595  }
596 
597  SCIPfreeBufferArray(scip, &lhs_array);
598  SCIPfreeBufferArray(scip, &cons_array);
599  SCIPfreeBufferArray(scip, &matrix);
600  }
601 
602  return SCIP_OKAY;
603 }
604 
612 static
613 SCIP_RETCODE diagZeroImpl(
614  SCIP* scip,
615  SCIP_CONS** conss,
616  int nconss,
617  int* naddconss
618  )
619 {
620  char cutname[SCIP_MAXSTRLEN];
621  /* if entry k is >= 0, gives the number of non-diagonal nonzero-entries in row k of the constant matrix (and the number
622  * of entries in constnonzeroentries), if -1 C_kk =!= 0 (which means we didnot allocate memory for diagvars), if -2
623  * either A_jk =!= 0 for all j with C_jk =!= 0 or A_kk =!= 0 for some continuous variable (which means we did allocate
624  * memory for diagvars but cannot use the cut */
625  int* nconstnonzeroentries;
626  int** constnonzeroentries;
627  SCIP_CONSDATA* consdata;
628  SCIP_CONS* cons;
629  SCIP_VAR** vars;
630  SCIP_Real* vals;
631  int** diagvars;
632  int* ndiagvars;
633  int blocksize;
634  int i;
635  int j;
636  int nvars;
637  int v;
638  int k;
639  int l;
640  SCIP_Bool anycutvalid;
641 #ifndef NDEBUG
642  int snprintfreturn;
643 #endif
644 
645  assert( scip != NULL );
646  assert( conss != NULL );
647  assert( nconss >= 0 );
648  assert( naddconss != NULL );
649 
650  for (i = 0; i < nconss; ++i)
651  {
652  assert( conss[i] != NULL );
653  assert( strcmp(SCIPconshdlrGetName(SCIPconsGetHdlr(conss[i])), "SDP") == 0 );
654 
655  consdata = SCIPconsGetData(conss[i]);
656  assert( consdata != NULL );
657 
658  blocksize = consdata->blocksize;
659  nvars = consdata->nvars;
660  SCIP_CALL( SCIPallocBufferArray(scip, &nconstnonzeroentries, blocksize) );
661  SCIP_CALL( SCIPallocBufferArray(scip, &constnonzeroentries, blocksize) );
662  for (j = 0; j < blocksize; j++)
663  {
664  nconstnonzeroentries[j] = 0;
665  SCIP_CALL( SCIPallocBufferArray(scip, &constnonzeroentries[j], blocksize) );
666  }
667 
668  /* iterate over all nonzeros of the constant matrix and check which diagonal and non-diagonal entries are nonzero */
669  for (j = 0; j < consdata->constnnonz; j++)
670  {
671  /* if it is a nondiagonal-entry we add this row/column to the constnonzeroentries entries unless we already found a
672  * diagonal entry for this row/column */
673  if ( (consdata->constcol[j] != consdata->constrow[j]) )
674  {
675  assert( ! SCIPisZero(scip, consdata->constval[j]) );
676  if ( nconstnonzeroentries[consdata->constcol[j]] >= 0 )
677  {
678  constnonzeroentries[consdata->constcol[j]][nconstnonzeroentries[consdata->constcol[j]]] = consdata->constrow[j];
679  nconstnonzeroentries[consdata->constcol[j]]++;
680  }
681  if ( nconstnonzeroentries[consdata->constrow[j]] >= 0 )
682  {
683  constnonzeroentries[consdata->constrow[j]][nconstnonzeroentries[consdata->constrow[j]]] = consdata->constcol[j];
684  nconstnonzeroentries[consdata->constrow[j]]++;
685  }
686  }
687  /* if we find a diagonal entry in the constant matrix, we remember that we cannot add a cut for this index */
688  else
689  {
690  assert( ! SCIPisZero(scip, consdata->constval[j]) );
691  nconstnonzeroentries[consdata->constcol[j]] = -1;
692  }
693  }
694 
695  /* diagvars[j] is an array with all variables with a diagonal entry (j,j) in the corresponding matrix, if nconstnonzeroentries[j] =!= -1 or NULL otherwise
696  * the outer array goes over all rows to ease the access, but only for those that are really needed memory will be allocated */
697  SCIP_CALL( SCIPallocBufferArray(scip, &diagvars, blocksize) );
698  SCIP_CALL( SCIPallocBufferArray(scip, &ndiagvars, blocksize) );
699  anycutvalid = FALSE;
700  for (j = 0; j < blocksize; ++j)
701  {
702  ndiagvars[j] = 0;
703  if ( nconstnonzeroentries[j] > 0 )
704  {
705  SCIP_CALL( SCIPallocBufferArray(scip, &(diagvars[j]), nvars) );
706  anycutvalid = TRUE;
707  }
708  }
709 
710  /* if no cuts are valid for this block, we free all memory and continue with the next block */
711  if ( ! anycutvalid )
712  {
713  SCIPfreeBufferArray(scip, &ndiagvars);
714  SCIPfreeBufferArray(scip, &diagvars);
715  for (j = blocksize - 1; j >= 0; j--)
716  {
717  SCIPfreeBufferArray(scip, &constnonzeroentries[j]);
718  }
719  SCIPfreeBufferArray(scip, &constnonzeroentries);
720  SCIPfreeBufferArray(scip, &nconstnonzeroentries);
721  continue;
722  }
723 
724  /* find all variables with corresponding diagonal entries for a row with nonzero non-diagonal constant entry, also check for entries
725  * that prevent the cut from being valid */
726  for (v = 0; v < nvars; v++)
727  {
728  for (j = 0; j < consdata->nvarnonz[v]; j++)
729  {
730  /* if it is a diagonal entry for an index that might have a valid cut, we add the variable to the corresponding array if it
731  * is an integer variable and mark the cut invalid otherwise */
732  if ( (consdata->col[v][j] == consdata->row[v][j]) && (nconstnonzeroentries[consdata->col[v][j]] > 0) )
733  {
734  if ( SCIPvarIsIntegral(consdata->vars[v]) )
735  {
736  assert( ! SCIPisEQ(scip, consdata->val[v][j], 0.0) );
737  diagvars[consdata->col[v][j]][ndiagvars[consdata->col[v][j]]] = v;
738  ndiagvars[consdata->col[v][j]]++;
739  }
740  else
741  {
742  nconstnonzeroentries[consdata->col[v][j]] = -2;
743  }
744  }
745  /* If it is a non-diagonal entry, we can no longer use this entry for a cut. If the last entry is removed for a column/row,
746  * mark this column/row invalid (but we still have to free memory later, so we have to set it to -2 instead of 0) */
747  else if ( consdata->col[v][j] != consdata->row[v][j] )
748  {
749  if ( nconstnonzeroentries[consdata->col[v][j]] > 0 )
750  {
751  /* search for the corresponding row-entry in constnonzeroentries */
752  for (k = 0; k < nconstnonzeroentries[consdata->col[v][j]]; k++)
753  {
754  if ( constnonzeroentries[consdata->col[v][j]][k] == consdata->row[v][j] )
755  {
756  /* if there are remaining entries, we shift them back */
757  if ( nconstnonzeroentries[consdata->col[v][j]] > k + 1 )
758  {
759  for (l = k + 1; l < nconstnonzeroentries[consdata->col[v][j]]; l++)
760  constnonzeroentries[consdata->col[v][j]][l - 1] = constnonzeroentries[consdata->col[v][j]][l];
761  }
762  nconstnonzeroentries[consdata->col[v][j]]--;
763  /* if this was the last entry for this index, we mark it invalid */
764  if ( nconstnonzeroentries[consdata->col[v][j]] == 0 )
765  nconstnonzeroentries[consdata->col[v][j]] = -2;
766  break; /* we should not have another entry for this combination of row and column */
767  }
768  }
769  }
770  /* do the same for the row */
771  if ( nconstnonzeroentries[consdata->row[v][j]] > 0 )
772  {
773  /* search for the corresponding row-entry in constnonzeroentries */
774  for (k = 0; k < nconstnonzeroentries[consdata->row[v][j]]; k++)
775  {
776  if ( constnonzeroentries[consdata->row[v][j]][k] == consdata->col[v][j] )
777  {
778  /* if there are remaining entries, we shift them back */
779  if ( nconstnonzeroentries[consdata->row[v][j]] > k + 1 )
780  {
781  for (l = k + 1; l < nconstnonzeroentries[consdata->row[v][j]]; l++)
782  constnonzeroentries[consdata->row[v][j]][l - 1] = constnonzeroentries[consdata->row[v][j]][l];
783  }
784  nconstnonzeroentries[consdata->row[v][j]]--;
785  /* if this was the last entry for this index, we mark it invalid */
786  if ( nconstnonzeroentries[consdata->row[v][j]] == 0 )
787  nconstnonzeroentries[consdata->row[v][j]] = -2;
788  break; /* we should not have another entry for this combination of row and column */
789  }
790  }
791  }
792  }
793  }
794  }
795 
796  for (j = 0; j < blocksize; ++j)
797  {
798  if ( nconstnonzeroentries[j] > 0 )
799  {
800  SCIP_CALL( SCIPallocBufferArray(scip, &vals, ndiagvars[j]) );
801  SCIP_CALL( SCIPallocBufferArray(scip, &vars, ndiagvars[j]) );
802 
803  /* get the corresponding SCIP variables and set all coefficients to 1 */
804  for (v = 0; v < ndiagvars[j]; ++v)
805  {
806  vars[v] = consdata->vars[diagvars[j][v]];
807  vals[v] = 1.0;
808  }
809 #ifndef NDEBUG
810  snprintfreturn = SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "diag_zero_impl_block_%d_row_%d", i, j);
811  assert( snprintfreturn < SCIP_MAXSTRLEN ); /* check whether name fits into string */
812 #else
813  (void) SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "diag_zero_impl_block_%d_row_%d", i, j);
814 #endif
815 
816 #ifdef SCIP_MORE_DEBUG
817  SCIPinfoMessage(scip, NULL, "Added lp-constraint %s: ", cutname);
818  SCIPinfoMessage(scip, NULL, "1 <=");
819  for (l = 0; l < ndiagvars[j]; l++)
820  SCIPinfoMessage(scip, NULL, " + (%f)*%s", vals[l], SCIPvarGetName(vars[l]));
821  SCIPinfoMessage(scip, NULL, "\n");
822 #endif
823 
824  /* add the linear constraint sum_j 1.0 * diagvars[j] >= 1.0 */
825  SCIP_CALL(SCIPcreateConsLinear(scip, &cons, cutname , ndiagvars[j], vars, vals, 1.0, SCIPinfinity(scip), TRUE, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE, TRUE, TRUE, FALSE));
826  SCIP_CALL(SCIPaddCons(scip, cons));
827  SCIP_CALL(SCIPreleaseCons(scip, &cons));
828  (*naddconss)++;
829 
830  SCIPfreeBufferArray(scip, &vars);
831  SCIPfreeBufferArray(scip, &vals);
832  }
833  if ( nconstnonzeroentries[j] == -2 || nconstnonzeroentries[j] > 0 )
834  {
835  SCIPfreeBufferArray(scip, &diagvars[j]);
836  }
837  }
838 
839  SCIPfreeBufferArray(scip, &ndiagvars);
840  SCIPfreeBufferArray(scip, &diagvars);
841  for (j = blocksize - 1; j >= 0; j--)
842  {
843  SCIPfreeBufferArray(scip, &constnonzeroentries[j]);
844  }
845  SCIPfreeBufferArray(scip, &constnonzeroentries);
846  SCIPfreeBufferArray(scip, &nconstnonzeroentries);
847  }
848 
849  return SCIP_OKAY;
850 }
851 
853 static
854 SCIP_RETCODE move_1x1_blocks_to_lp(
855  SCIP* scip,
856  SCIP_CONS** conss,
857  int nconss,
858  int* naddconss,
859  int* ndelconss,
860  SCIP_RESULT* result
861  )
862 {
863  char cutname[SCIP_MAXSTRLEN];
864  SCIP_CONSHDLR* hdlr;
865  SCIP_CONSDATA* consdata;
866  SCIP_CONSHDLRDATA* conshdlrdata;
867  SCIP_CONS* cons;
868  SCIP_VAR** vars;
869  SCIP_Real* coeffs;
870  SCIP_Real rhs;
871  int nnonz;
872  int nvars;
873  int i;
874  int j;
875  int count;
876  int var;
877 #ifndef NDEBUG
878  int snprintfreturn; /* used to assert the return code of snprintf */
879 #endif
880 
881  assert( result != NULL );
882  *result = SCIP_SUCCESS;
883 
884  for (i = 0; i < nconss; ++i)
885  {
886  hdlr = SCIPconsGetHdlr(conss[i]);
887  assert(hdlr != NULL);
888 
889 #ifndef NDEBUG
890  assert( strcmp(SCIPconshdlrGetName(hdlr), "SDP") == 0);
891 #endif
892 
893  consdata = SCIPconsGetData(conss[i]);
894  assert( consdata != NULL );
895 
896  /* if there is a 1x1 SDP-Block */
897  if ( consdata->blocksize == 1 )
898  {
899  nvars = consdata->nvars;
900  nnonz = consdata->nnonz;
901  SCIP_CALL( SCIPallocBufferArray(scip, &vars, nvars) );
902  SCIP_CALL( SCIPallocBufferArray(scip, &coeffs, nnonz) );
903 
904  /* get all lhs-entries */
905  count = 0;
906 
907  for (var = 0; var < nvars; var++)
908  {
909  assert( consdata->nvarnonz[var] <= 1 ); /* in a 1x1 block there may be at most one entry per variable */
910 
911  for (j = 0; j < consdata->nvarnonz[var]; j++)
912  {
913  assert( consdata->col[var][j] == 0 && consdata->row[var][j] == 0 ); /* if the block is size one, all entries should have row and col equal to 0 */
914  coeffs[count] = consdata->val[var][j];
915  vars[count] = consdata->vars[var];
916  count++;
917  }
918  }
919 
920  /* get rhs */
921  assert( consdata->constnnonz <= 1 ); /* the 1x1 constant matrix may only have one entry */
922 
923  rhs = (consdata->constnnonz == 1) ? consdata->constval[0] : 0.0; /* if this one entry is not 0, than this is the rhs, otherwise it's 0 */
924 
925  /* if there is more than one left-hand-side-entry, add a linear constraint, otherwise update the variable bound */
926  if ( count > 1 )
927  {
928  /* add new linear cons */
929  conshdlrdata = SCIPconshdlrGetData(hdlr);
930 #ifndef NDEBUG
931  snprintfreturn = SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "1x1block_%d", ++(conshdlrdata->n1x1blocks));
932  assert( snprintfreturn < SCIP_MAXSTRLEN ); /* check whether name fits into string */
933 #else
934  (void) SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "1x1block_%d", ++(conshdlrdata->n1x1blocks));
935 #endif
936 
937 #ifdef SCIP_MORE_DEBUG
938  SCIPinfoMessage(scip, NULL, "Added lp-constraint %s: ", cutname);
939  for (i = 0; i < consdata->nvars; i++)
940  SCIPinfoMessage(scip, NULL, "+ (%f)*%s", coeffs[i], SCIPvarGetName(vars[i]));
941  SCIPinfoMessage(scip, NULL, " <= %f \n", rhs);
942 #endif
943 
944  SCIP_CALL( SCIPcreateConsLinear(scip, &cons, cutname, consdata->nvars, vars, coeffs, rhs, SCIPinfinity(scip),
945  TRUE, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE, TRUE, TRUE, FALSE) );
946 
947  SCIP_CALL( SCIPaddCons(scip, cons) );
948  SCIP_CALL( SCIPreleaseCons(scip, &cons) );
949 
950  (*naddconss)++;
951  }
952  else
953  {
954  /* we compare this new variable with the current (local) bounds, we don't do an epsilon check here, because 10^(-21)*x >= 10^(-19) for
955  * epsilon = 10^(-20) would be infeasible then, even though it only says x >= 100 */
956  if ( coeffs[0] > 0.0 )
957  {
958  /* this gives a lower bound, if it is bigger than the current one, we need to update it */
959  if ( SCIPisFeasGT(scip, rhs / coeffs[0], SCIPvarGetLbLocal(vars[0])) )
960  {
961  /* check if the changed bound renders the problem infeasible */
962  if( SCIPisFeasGT(scip, rhs / coeffs[0], SCIPvarGetUbLocal(vars[0])) )
963  {
964  SCIPdebugMessage("Problem detected to be infeasible during presolving, 1x1-SDP-constraint %s caused change"
965  "of lower bound for variable %s from %f to %f, which is bigger than upper bound of %f\n",
966  SCIPconsGetName(conss[i]), SCIPvarGetName(vars[0]), SCIPvarGetLbLocal(vars[0]), rhs / coeffs[0],
967  SCIPvarGetUbLocal(vars[0]));
968 
969  *result = SCIP_CUTOFF;
970 
971  /* delete old 1x1 sdpcone */
972  SCIP_CALL( SCIPdelCons(scip, conss[i]) );
973  (*ndelconss)++;
974 
975  SCIPfreeBufferArray(scip, &coeffs);
976  SCIPfreeBufferArray(scip, &vars);
977 
978  return SCIP_OKAY; /* the node is infeasible, we don't care for the other constraints */
979  }
980 
981  SCIPdebugMessage("Changing lower bound of variable %s from %f to %f because of 1x1-SDP-constraint %s!\n",
982  SCIPvarGetName(vars[0]), SCIPvarGetLbLocal(vars[0]), rhs / coeffs[0], SCIPconsGetName(conss[i]));
983  SCIP_CALL( SCIPchgVarLb(scip, vars[0], rhs / coeffs[0]) );
984  }
985  else
986  {
987  SCIPdebugMessage("Deleting 1x1-SDP-constraint %s, new lower bound %f for variable %s no improvement over old bound %f!\n",
988  SCIPconsGetName(conss[i]), rhs / coeffs[0], SCIPvarGetName(vars[0]), SCIPvarGetLbLocal(vars[0]));
989  }
990  }
991  else if ( coeffs[0] < 0.0 )
992  {
993  /* this gives an upper bound, if it is lower than the current one, we need to update it */
994  if (SCIPisFeasLT(scip, rhs / coeffs[0], SCIPvarGetUbLocal(vars[0])))
995  {
996  /* check if the changed bound renders the problem infeasible */
997  if( SCIPisFeasLT(scip, rhs / coeffs[0], SCIPvarGetLbLocal(vars[0])) )
998  {
999  SCIPdebugMessage("Problem detected to be infeasible during presolving, 1x1-SDP-constraint %s caused change"
1000  "of upper bound for variable %s from %f to %f, which is less than lower bound of %f\n",
1001  SCIPconsGetName(conss[i]), SCIPvarGetName(vars[0]), SCIPvarGetUbLocal(vars[0]), rhs / coeffs[0],
1002  SCIPvarGetLbLocal(vars[0]));
1003 
1004  *result = SCIP_CUTOFF;
1005 
1006  /* delete old 1x1 sdpcone */
1007  SCIP_CALL( SCIPdelCons(scip, conss[i]) );
1008  (*ndelconss)++;
1009 
1010  SCIPfreeBufferArray(scip, &coeffs);
1011  SCIPfreeBufferArray(scip, &vars);
1012 
1013  return SCIP_OKAY; /* the node is infeasible, we don't care for the other constraints */
1014  }
1015 
1016  SCIPdebugMessage("Changing upper bound of variable %s from %f to %f because of 1x1-SDP-constraint %s!\n",
1017  SCIPvarGetName(vars[0]), SCIPvarGetUbLocal(vars[0]), -rhs / coeffs[0], SCIPconsGetName(conss[i]));
1018  SCIP_CALL( SCIPchgVarUb(scip, vars[0], rhs / coeffs[0]) );
1019  }
1020  else
1021  {
1022  SCIPdebugMessage("Deleting 1x1-SDP-constraint %s, new upper bound %f for variable %s no improvement over old bound %f!\n",
1023  SCIPconsGetName(conss[i]), rhs / coeffs[0], SCIPvarGetName(vars[0]), SCIPvarGetUbLocal(vars[0]));
1024  }
1025  }
1026  else
1027  {
1028  SCIPdebugMessage("Detected 1x1 SDP-block without any nonzero coefficients \n");
1029  if ( SCIPisFeasGT(scip, rhs, 0.0) )
1030  {
1031  SCIPdebugMessage("Detected infeasibility in 1x1 SDP-block without any nonzero coefficients but with strictly positive rhs\n");
1032  *result = SCIP_CUTOFF;
1033 
1034  /* delete old 1x1 sdpcone */
1035  SCIP_CALL( SCIPdelCons(scip, conss[i]) );
1036  (*ndelconss)++;
1037 
1038  SCIPfreeBufferArray(scip, &coeffs);
1039  SCIPfreeBufferArray(scip, &vars);
1040 
1041  return SCIP_OKAY; /* the node is infeasible, we don't care for the other constraints */
1042  }
1043  }
1044  }
1045 
1046  /* delete old 1x1 sdpcone */
1047  SCIP_CALL( SCIPdelCons(scip, conss[i]) );
1048  (*ndelconss)++;
1049 
1050  SCIPfreeBufferArray(scip, &coeffs);
1051  SCIPfreeBufferArray(scip, &vars);
1052  }
1053  }
1054  return SCIP_OKAY;
1055 }
1056 
1058 static
1059 SCIP_RETCODE multiaggrVar(
1060  SCIP* scip,
1061  SCIP_CONS* cons,
1062  int* v,
1063  SCIP_VAR** aggrvars,
1064  SCIP_Real* scalars,
1065  int naggrvars,
1066  SCIP_Real constant,
1067  int* savedcol,
1068  int* savedrow,
1069  SCIP_Real* savedval,
1070  int* nfixednonz,
1071  int* vararraylength
1072  )
1073 {
1074  int i;
1075  SCIP_CONSDATA* consdata;
1076  int startind;
1077  int aggrind;
1078  int aggrtargetlength;
1079  int globalnvars;
1080  int aggrconsind;
1081  SCIP_Real epsilon;
1082 
1083  assert( scip != NULL );
1084  assert( cons != NULL );
1085  assert( scalars != NULL );
1086  assert( naggrvars > 0 );
1087  assert( savedcol != NULL );
1088  assert( savedrow != NULL );
1089  assert( savedval != NULL );
1090  assert( nfixednonz != NULL );
1091 
1092  consdata = SCIPconsGetData(cons);
1093  assert( consdata != NULL );
1094 
1095  SCIP_CALL( SCIPgetRealParam(scip, "numerics/epsilon", &epsilon) );
1096 
1097  /* save the current nfixednonz-index, all entries starting from here will need to be added to the variables this is aggregated to */
1098  startind = *nfixednonz;
1099 
1100  if ( SCIPisEQ(scip, constant, 0.0) )
1101  {
1102  /* if there is no constant part, we just save the nonzeros to add them to the variables they are aggregated to, we do this to be able to remove
1103  * the nonzero-arrays for this variable to be able to fill it with a newly inserted variable, as copying all variables, if we created an empty
1104  * gap in the variable arrays, will be more time consuming then copying all variables (as we will usually have more variables than nonzeros per
1105  * variable */
1106  for (i = 0; i < consdata->nvarnonz[*v]; i++)
1107  {
1108  savedcol[*nfixednonz] = consdata->col[*v][i];
1109  savedrow[*nfixednonz] = consdata->row[*v][i];
1110  savedval[*nfixednonz] = consdata->val[*v][i];
1111  (*nfixednonz)++;
1112  }
1113  }
1114  else
1115  {
1116  for (i = 0; i < consdata->nvarnonz[*v]; i++)
1117  {
1118  savedcol[*nfixednonz] = consdata->col[*v][i];
1119  savedrow[*nfixednonz] = consdata->row[*v][i];
1120  savedval[*nfixednonz] = consdata->val[*v][i] * constant;
1121  (*nfixednonz)++;
1122  }
1123  }
1124  assert( *nfixednonz - startind == consdata->nvarnonz[*v] );
1125 
1126  /* sort them by nondecreasing row and then col to make the search for already existing entries easier (this is done here, because it
1127  * only needs to be done once and not for each variable this is multiaggregated to), we add startind to the pointers to only start where we started
1128  * inserting, the number of elements added to the saved arrays for this variable is nfixednonz - startind */
1129  SCIPsdpVarfixerSortRowCol(savedrow + startind, savedcol + startind, savedval + startind, *nfixednonz - startind);
1130 
1131  /* free the memory for the entries of the aggregated variable */
1132  SCIPfreeBlockMemoryArray(scip, &(consdata->val[*v]), consdata->nvarnonz[*v]);
1133  SCIPfreeBlockMemoryArray(scip, &(consdata->row[*v]), consdata->nvarnonz[*v]);
1134  SCIPfreeBlockMemoryArray(scip, &(consdata->col[*v]), consdata->nvarnonz[*v]);
1135 
1136  /* fill the empty spot of the (multi-)aggregated variable with the last variable of this constraint (as they don't have to be sorted) */
1137  SCIP_CALL( SCIPreleaseVar(scip, &consdata->vars[*v]) );
1138  consdata->col[*v] = consdata->col[consdata->nvars - 1];
1139  consdata->row[*v] = consdata->row[consdata->nvars - 1];
1140  consdata->val[*v] = consdata->val[consdata->nvars - 1];
1141  consdata->nvarnonz[*v] = consdata->nvarnonz[consdata->nvars - 1];
1142  consdata->vars[*v] = consdata->vars[consdata->nvars - 1];
1143  (consdata->nvars)--;
1144  (*v)--; /* we need to check again if the variable we just shifted to this position also needs to be (multi-)aggregated */
1145 
1146  /* iterate over all variables this was aggregated to and insert the corresponding nonzeros */
1147  for (aggrind = 0; aggrind < naggrvars; aggrind++)
1148  {
1149  /* check if the variable already exists in this block */
1150  aggrconsind = -1;
1151  for (i = 0; i < consdata->nvars; i++)
1152  {
1153  if ( consdata->vars[i] == aggrvars[aggrind] )
1154  {
1155  aggrconsind = i;
1156  break;
1157  }
1158  }
1159 
1160  if ( aggrconsind > -1 )
1161  {
1162  /* if the varialbe to aggregate to is already part of this sdp-constraint, just add the nonzeros of the old variable to it */
1163 
1164  /* resize the arrays to the maximally needed length */
1165  aggrtargetlength = consdata->nvarnonz[aggrconsind] + *nfixednonz - startind;
1166  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->row[aggrconsind]), consdata->nvarnonz[aggrconsind], aggrtargetlength) );
1167  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->col[aggrconsind]), consdata->nvarnonz[aggrconsind], aggrtargetlength) );
1168  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->val[aggrconsind]), consdata->nvarnonz[aggrconsind], aggrtargetlength) );
1169 
1170  if ( SCIPisEQ(scip, constant, 0.0) )
1171  {
1172  /* in this case we saved the original values in savedval, we add startind to the pointers to only add those from
1173  * the current variable, the number of entries is the current position minus the position whre we started */
1174  SCIP_CALL( SCIPsdpVarfixerMergeArrays(SCIPblkmem(scip), epsilon, savedrow + startind, savedcol + startind, savedval + startind,
1175  *nfixednonz - startind, TRUE, scalars[aggrind], consdata->row[aggrconsind], consdata->col[aggrconsind],
1176  consdata->val[aggrconsind], &(consdata->nvarnonz[aggrconsind]), aggrtargetlength) );
1177  }
1178  else
1179  {
1180  /* in this case we saved the original values * constant, so we now have to divide by constant, we add startind to the pointers
1181  * to only add those from the current variable, the number of entries is the current position minus the position whre we started */
1182  SCIP_CALL( SCIPsdpVarfixerMergeArrays(SCIPblkmem(scip), epsilon, savedrow + startind, savedcol + startind, savedval + startind,
1183  *nfixednonz - startind, TRUE, scalars[aggrind] / constant, consdata->row[aggrconsind], consdata->col[aggrconsind],
1184  consdata->val[aggrconsind], &(consdata->nvarnonz[aggrconsind]), aggrtargetlength) );
1185  }
1186 
1187  /* shrink them again if nonzeros could be combined */
1188  assert( consdata->nvarnonz[aggrconsind] <= aggrtargetlength );
1189  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->row[aggrconsind]), aggrtargetlength, consdata->nvarnonz[aggrconsind]) );
1190  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->col[aggrconsind]), aggrtargetlength, consdata->nvarnonz[aggrconsind]) );
1191  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->val[aggrconsind]), aggrtargetlength, consdata->nvarnonz[aggrconsind]) );
1192  }
1193  else
1194  {
1195  /* the variable has to be added to this constraint */
1196 
1197  SCIPdebugMessage("adding variable %s to SDP constraint %s because of (multi-)aggregation\n", SCIPvarGetName(aggrvars[aggrind]), SCIPconsGetName(cons));
1198 
1199  /* check if we have to enlarge the arrays */
1200  if ( consdata->nvars == *vararraylength )
1201  {
1202  globalnvars = SCIPgetNVars(scip);
1203 
1204  /* we don't want to enlarge this by one for every variable added, so we immediately set it to the maximum possible size */
1205  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->col, *vararraylength, globalnvars) );
1206  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->row, *vararraylength, globalnvars) );
1207  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->val, *vararraylength, globalnvars) );
1208  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->nvarnonz, *vararraylength, globalnvars) );
1209  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->vars, *vararraylength, globalnvars) );
1210  *vararraylength = globalnvars;
1211  }
1212 
1213  /* we insert this variable at the last position, as the ordering doesn't matter */
1214  SCIP_CALL( SCIPcaptureVar(scip, aggrvars[aggrind]) );
1215  consdata->vars[consdata->nvars] = aggrvars[aggrind];
1216 
1217  /* as there were no nonzeros thus far, we can just duplicate the saved arrays to get the nonzeros for the new variable */
1218  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(consdata->col[consdata->nvars]), savedcol + startind, *nfixednonz - startind) ); /*lint !e776*/
1219  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(consdata->row[consdata->nvars]), savedrow + startind, *nfixednonz - startind) ); /*lint !e776*/
1220 
1221  /* if scalars[aggrind] = constant, we would multiply with 1, if constant = 0, we didn't divide by constant, so in these cases, we can just
1222  * memcopy the array of nonzero-values */
1223  /* TODO: only checking scalar and constant for feas-equality might lead to big differences, if the nonzeros they are multiplied with are big,
1224  * e.g. scalar = 0, constant = 10^(-6), nonzero = 10^(10) leads to new nonzero of 10^4 instead of 0 */
1225  if ( SCIPisEQ(scip, scalars[aggrind], constant) || SCIPisEQ(scip, constant, 0.0) )
1226  {
1227  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(consdata->val[consdata->nvars]), savedval + startind, *nfixednonz - startind) ); /*lint !e776*/
1228  consdata->nvarnonz[consdata->nvars] = *nfixednonz - startind; /* as there were no nonzeros thus far, the number of nonzeros equals
1229  * the number of nonzeros of the aggregated variable */
1230  }
1231  else /* we have to multiply all entries by scalar before inserting them */
1232  {
1233  SCIP_CALL( SCIPgetRealParam(scip, "numerics/epsilon", &epsilon) );
1234 
1235  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(consdata->val[consdata->nvars]), *nfixednonz - startind) );
1236 
1237  consdata->nvarnonz[consdata->nvars] = 0;
1238 
1239  for (i = 0; i < *nfixednonz - startind; i++)
1240  {
1241  /* if both scalar and savedval are small this might become too small */
1242  if ( (scalars[i] / constant) * savedval[startind + i] >= epsilon ) /*lint !e679*/
1243  {
1244  consdata->val[consdata->nvars][consdata->nvarnonz[consdata->nvars]] = (scalars[i] / constant) * savedval[startind + i]; /*lint !e679*/
1245  consdata->nvarnonz[consdata->nvars]++;
1246  }
1247  }
1248  }
1249 
1250  if ( consdata->nvarnonz[consdata->nvars] > 0 ) /* if scalar and all savedvals were to small */
1251  consdata->nvars++;
1252  }
1253  }
1254 
1255  /* if the constant part is zero, we may delete the nonzero-entries from the saved arrays (by resetting the nfixednonz entry to where
1256  * it started, so that these entries will be overwritten */
1257  if ( SCIPisEQ(scip, constant, 0.0) )
1258  *nfixednonz = startind;
1259 
1260  return SCIP_OKAY;
1261 }
1262 
1263 
1265 static
1266 SCIP_RETCODE fixAndAggrVars(
1267  SCIP* scip,
1268  SCIP_CONS** conss,
1269  int nconss,
1270  SCIP_Bool aggregate
1271  )
1272 {
1273  SCIP_CONSDATA* consdata;
1274  int i;
1275  int nfixednonz;
1276  int* savedcol;
1277  int* savedrow;
1278  SCIP_Real* savedval;
1279  int c;
1280  int v;
1281  int arraylength;
1282  SCIP_VAR* var;
1283  SCIP_VAR** aggrvars;
1284  SCIP_Real scalar;
1285  SCIP_Real* scalars;
1286  int naggrvars;
1287  SCIP_Real constant;
1288  int requiredsize;
1289  int globalnvars;
1290  int vararraylength;
1291  SCIP_Bool negated;
1292  SCIP_Real epsilon;
1293 
1294 
1295  /* loop over all variables once, add all fixed to savedrow/col/val, for all multiaggregated variables, if constant-scalar =!= 0, add
1296  * constant-scalar * entry to savedrow/col/val and call mergeArrays for all aggrvars for savedrow[startindex of this var] and scalar/constant-scalar,
1297  * if constant-scalar = 0, add 1*entry to savedrow/col/val, call mergeArrays for all aggrvars for savedrow[startindex of this var] and scalar and later
1298  * reduze the saved size of savedrow/col/val by the number of nonzeros of the mutliagrregated variable to not add them to the constant part later */
1299 
1300  assert( scip != NULL );
1301  assert( conss != NULL );
1302  assert( nconss >= 0 );
1303 
1304  SCIPdebugMessage("Calling fixAndAggrVars with aggregate = %u\n", aggregate);
1305 
1306  SCIP_CALL( SCIPgetRealParam(scip, "numerics/epsilon", &epsilon) );
1307 
1308  for (c = 0; c < nconss; ++c)
1309  {
1310  assert( conss[c] != NULL );
1311  assert( strcmp(SCIPconshdlrGetName(SCIPconsGetHdlr(conss[c])), "SDP") == 0);
1312 
1313  consdata = SCIPconsGetData(conss[c]);
1314  assert( consdata != NULL );
1315 
1316  /* allocate memory to save nonzeros that need to be fixed */
1317  SCIP_CALL( SCIPallocBufferArray(scip, &savedcol, consdata->nnonz) );
1318  SCIP_CALL( SCIPallocBufferArray(scip, &savedrow, consdata->nnonz) );
1319  SCIP_CALL( SCIPallocBufferArray(scip, &savedval, consdata->nnonz) );
1320 
1321  /* initialize this with zero for each block */
1322  nfixednonz = 0;
1323 
1324  vararraylength = consdata->nvars;
1325  globalnvars = SCIPgetNVars(scip);
1326 
1327  for (v = 0; v < consdata->nvars; v++)/*lint --e{850}*/
1328  {
1329  negated = FALSE;
1330  /* if the variable is negated, get the negation var */
1331  if ( SCIPvarIsBinary(consdata->vars[v]) && SCIPvarIsNegated(consdata->vars[v]) )
1332  {
1333  negated = TRUE;
1334  var = SCIPvarGetNegationVar(consdata->vars[v]);
1335  }
1336  else
1337  var = consdata->vars[v];
1338 
1339  /* check if the variable is fixed in SCIP */
1340  if ( SCIPvarGetStatus(var) == SCIP_VARSTATUS_FIXED || SCIPisEQ(scip, SCIPvarGetLbGlobal(var), SCIPvarGetUbGlobal(var)))
1341  {
1342  assert( SCIPisEQ(scip, SCIPvarGetLbGlobal(var), SCIPvarGetUbGlobal(var)) );
1343 
1344  SCIPdebugMessage("Globally fixing Variable %s to value %f !\n", SCIPvarGetName(var), SCIPvarGetLbGlobal(var));
1345 
1346  if ( ((! negated) && (! SCIPisEQ(scip, SCIPvarGetLbGlobal(var), 0.0))) || (negated && SCIPisEQ(scip, SCIPvarGetLbGlobal(var), 0.0)) )
1347  {
1348  /* the nonzeros are saved to later be inserted into the constant part (this is only done after all nonzeros of fixed variables have been
1349  * assembled, because we need to sort the constant nonzeros and loop over them, which we only want to do once and not once for each fixed
1350  * variable) */
1351  for (i = 0; i < consdata->nvarnonz[v]; i++)
1352  {
1353  savedcol[nfixednonz] = consdata->col[v][i];
1354  savedrow[nfixednonz] = consdata->row[v][i];
1355  /* this is the final value to add, we no longer have to remember from which variable this comes, minus because we have +A_i but -A_0 */
1356  if ( ! negated )
1357  savedval[nfixednonz] = consdata->val[v][i] * SCIPvarGetLbGlobal(var);
1358  else
1359  savedval[nfixednonz] = consdata->val[v][i]; /* if it is the negation of a variable fixed to zero, this variable is fixed to one */
1360 
1361  nfixednonz++;
1362  consdata->nnonz--;
1363  }
1364  }
1365  else
1366  {
1367  /* if the variable is fixed to zero, the nonzeros will just vanish, so we only reduce the number of nonzeros */
1368  consdata->nnonz -= consdata->nvarnonz[v];
1369  }
1370  /* free the memory of the corresponding entries in col/row/val */
1371  SCIPfreeBlockMemoryArrayNull(scip, &(consdata->val[v]), consdata->nvarnonz[v]);
1372  SCIPfreeBlockMemoryArrayNull(scip, &(consdata->row[v]), consdata->nvarnonz[v]);
1373  SCIPfreeBlockMemoryArrayNull(scip, &(consdata->col[v]), consdata->nvarnonz[v]);
1374 
1375  /* as the variables don't need to be sorted, we just put the last variable into the empty spot and decrease sizes by one (at the end) */
1376  SCIP_CALL( SCIPreleaseVar(scip, &(consdata->vars[v])) );
1377  consdata->col[v] = consdata->col[consdata->nvars - 1];
1378  consdata->row[v] = consdata->row[consdata->nvars - 1];
1379  consdata->val[v] = consdata->val[consdata->nvars - 1];
1380  consdata->nvarnonz[v] = consdata->nvarnonz[consdata->nvars - 1];
1381  consdata->vars[v] = consdata->vars[consdata->nvars - 1];
1382  consdata->nvars--;
1383  v--; /* we need to check again if the variable we just shifted to this position also needs to be fixed */
1384  }
1385  else if ( (SCIPvarGetStatus(var) == SCIP_VARSTATUS_AGGREGATED ||
1386  SCIPvarGetStatus(var) == SCIP_VARSTATUS_MULTAGGR)
1387  && aggregate )
1388  {
1389  SCIP_CALL( SCIPallocBufferArray(scip, &aggrvars, globalnvars) );
1390  SCIP_CALL( SCIPallocBufferArray(scip, &scalars, globalnvars) );
1391 
1392  /* this is how they should be initialized before calling SCIPgetProbvarLinearSum */
1393  if (! negated)
1394  {
1395  aggrvars[0] = consdata->vars[v];
1396  naggrvars = 1;
1397  constant = 0.0;
1398  scalars[0] = 1.0;
1399  }
1400  else
1401  {
1402  /* if this variable is the negation of var, than we look for a representation of 1-var */
1403  aggrvars[0] = consdata->vars[v];
1404  naggrvars = 1;
1405  constant = 1.0;
1406  scalars[0] = -1.0;
1407  }
1408 
1409  /* get the variables this var was aggregated to */
1410  SCIP_CALL( SCIPgetProbvarLinearSum(scip, aggrvars, scalars, &naggrvars, globalnvars, &constant, &requiredsize, TRUE) );
1411  assert( requiredsize <= globalnvars ); /* requiredsize is the number of empty spots in aggrvars needed, globalnvars is the number
1412  * of spots we provided */
1413 
1414  /* Debugmessages for the (multi-)aggregation */
1415 #ifdef SCIP_DEBUG
1416  if ( SCIPvarGetStatus(consdata->vars[v]) == SCIP_VARSTATUS_AGGREGATED )
1417  SCIPdebugMessage("aggregating variable %s to ", SCIPvarGetName(var));
1418  else
1419  SCIPdebugMessage("multiaggregating variable %s to ", SCIPvarGetName(var));
1420  for (i = 0; i < naggrvars; i++)
1421  SCIPdebugMessage("+ (%f2) * %s ", scalars[i], SCIPvarGetName(aggrvars[i]));
1422  SCIPdebugMessage("+ (%f2) \n", constant);
1423 #endif
1424 
1425  /* add the nonzeros to the saved-arrays for the constant part, remove the nonzeros for the old variables and add them to the variables this variable
1426  * was (multi-)aggregated to */
1427  SCIP_CALL( multiaggrVar(scip, conss[c], &v, aggrvars, scalars, naggrvars, constant, savedcol, savedrow, savedval, &nfixednonz, &vararraylength) );
1428 
1429  SCIPfreeBufferArray(scip, &aggrvars);
1430  SCIPfreeBufferArray(scip, &scalars);
1431  }
1432  else if ( negated && (SCIPvarGetStatus(var) == SCIP_VARSTATUS_LOOSE || SCIPvarGetStatus(var) == SCIP_VARSTATUS_COLUMN)
1433  && aggregate)
1434  {
1435  /* if var1 is the negation of var2, then this is equivalent to it being aggregated to -var2 + 1 = 1 - var2 */
1436 
1437  SCIPdebugMessage("Changing variable %s to negation of variable %s !\n", SCIPvarGetName(consdata->vars[v]), SCIPvarGetName(var));
1438 
1439  scalar = -1.0;
1440 
1441  SCIP_CALL( multiaggrVar(scip, conss[c], &v, &var, &scalar, 1, 1.0, savedcol, savedrow, savedval, &nfixednonz, &vararraylength) );
1442  }
1443  }
1444 
1445  /* shrink the variable arrays if they were enlarged too much (or more vars were removed than added) */
1446  assert( consdata->nvars <= vararraylength );
1447  if ( consdata->nvars < vararraylength )
1448  {
1449  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->col, vararraylength, consdata->nvars) );
1450  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->row, vararraylength, consdata->nvars) );
1451  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->val, vararraylength, consdata->nvars) );
1452  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->nvarnonz, vararraylength, consdata->nvars) );
1453  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->vars, vararraylength, consdata->nvars) );
1454  }
1455 
1456  /* allocate the maximally needed memory for inserting the fixed variables into the constant part */
1457  arraylength = consdata->constnnonz + nfixednonz;
1458  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->constcol), consdata->constnnonz, arraylength) );
1459  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->constrow), consdata->constnnonz, arraylength) );
1460  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->constval), consdata->constnnonz, arraylength) );
1461 
1462  /* insert the fixed variables into the constant arrays, as we have +A_i but -A_0 we mutliply them by -1 */
1463  SCIP_CALL( SCIPsdpVarfixerMergeArrays(SCIPblkmem(scip), epsilon, savedrow, savedcol, savedval, nfixednonz, FALSE, -1.0, consdata->constrow,
1464  consdata->constcol, consdata->constval, &(consdata->constnnonz), arraylength) );
1465 
1466  assert( consdata->constnnonz <= arraylength ); /* the allocated memory should always be sufficient */
1467 
1468  /* shrink the arrays if nonzeros could be combined */
1469  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->constcol), arraylength, consdata->constnnonz) );
1470  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->constrow), arraylength, consdata->constnnonz) );
1471  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &(consdata->constval), arraylength, consdata->constnnonz) );
1472 
1473  /* free the saved arrays */
1474  SCIPfreeBufferArray(scip, &savedval);
1475  SCIPfreeBufferArray(scip, &savedrow);
1476  SCIPfreeBufferArray(scip, &savedcol);
1477 
1478  /* recompute sdpnnonz */
1479  consdata->nnonz = 0;
1480  for (v = 0; v < consdata->nvars; v++)
1481  consdata->nnonz += consdata->nvarnonz[v];
1482  }
1483 
1484  return SCIP_OKAY;
1485 }
1486 
1488 static
1489 SCIP_RETCODE EnforceConstraint(
1490  SCIP* scip,
1491  SCIP_CONSHDLR* conshdlr,
1492  SCIP_CONS** conss,
1493  int nconss,
1494  SCIP_SOL* sol,
1495  SCIP_RESULT* result
1496  )
1497 {
1498  char cutname[SCIP_MAXSTRLEN];
1499  SCIP_CONSDATA* consdata;
1500  SCIP_CONSHDLRDATA* conshdlrdata;
1501  SCIP_Bool all_feasible = TRUE;
1502  SCIP_Bool separated = FALSE;
1503  SCIP_ROW* row;
1504  SCIP_Bool infeasible;
1505  SCIP_Real lhs;
1506  SCIP_Real* coeff;
1507  SCIP_Real rhs;
1508  int nvars;
1509  int i;
1510  int j;
1511 #ifndef NDEBUG
1512  int snprintfreturn; /* used to check the return code of snprintf */
1513 #endif
1514 
1515  *result = SCIP_FEASIBLE;
1516 
1517  for (i = 0; i < nconss; ++i)
1518  {
1519  consdata = SCIPconsGetData(conss[i]);
1520  SCIP_CALL( SCIPconsSdpCheckSdpCons(scip, conss[i], sol, 0, 0, 0, result) );
1521  if ( *result == SCIP_FEASIBLE )
1522  continue;
1523 
1524  all_feasible = FALSE;
1525 
1526  nvars = consdata->nvars;
1527  lhs = 0.0;
1528  coeff = NULL;
1529 
1530  SCIP_CALL( SCIPallocBufferArray(scip, &coeff, nvars) );
1531  SCIP_CALL( cutUsingEigenvector(scip, conss[i], sol, coeff, &lhs) );
1532 
1533  rhs = SCIPinfinity(scip);
1534  conshdlrdata = SCIPconshdlrGetData(conshdlr);
1535 
1536 #ifndef NDEBUG
1537  snprintfreturn = SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "sepa_eig_sdp_%d", ++(conshdlrdata->neigveccuts));
1538  assert( snprintfreturn < SCIP_MAXSTRLEN ); /* check whether the name fits into the string */
1539 #else
1540  (void) SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "sepa_eig_sdp_%d", ++(conshdlrdata->neigveccuts));
1541 #endif
1542 
1543  SCIP_CALL( SCIPcreateEmptyRowCons(scip, &row, conshdlr, cutname , lhs, rhs, FALSE, FALSE, TRUE) );
1544  SCIP_CALL( SCIPcacheRowExtensions(scip, row) );
1545 
1546  for (j = 0; j < nvars; ++j)
1547  {
1548  SCIP_CALL( SCIPaddVarToRow(scip, row, consdata->vars[j], coeff[j]) );
1549  }
1550 
1551  SCIP_CALL( SCIPflushRowExtensions(scip, row) );
1552 
1553 #ifdef SCIP_MORE_DEBUG
1554  SCIPinfoMessage(scip, NULL, "Added cut %s: ", cutname);
1555  SCIPinfoMessage(scip, NULL, "%f <= ", lhs);
1556  for (j = 0; j < nvars; j++)
1557  SCIPinfoMessage(scip, NULL, "+ (%f)*%s", coeff[j], SCIPvarGetName(consdata->vars[j]));
1558  SCIPinfoMessage(scip, NULL, "\n");
1559 #endif
1560 
1561  SCIP_CALL( SCIPaddRow(scip, row, FALSE, &infeasible) );
1562 
1563  if ( infeasible )
1564  {
1565  *result = SCIP_CUTOFF;
1566 
1567  SCIP_CALL( SCIPreleaseRow(scip, &row) );
1568  SCIPfreeBufferArray(scip, &coeff);
1569 
1570  return SCIP_OKAY;
1571  }
1572  else
1573  {
1574  SCIP_CALL( SCIPaddPoolCut(scip, row) );
1575 
1576  SCIP_CALL( SCIPresetConsAge(scip, conss[i]) );
1577  *result = SCIP_SEPARATED;
1578  separated = TRUE;
1579  }
1580  SCIP_CALL( SCIPreleaseRow(scip, &row) );
1581  SCIPfreeBufferArray(scip, &coeff);
1582  }
1583 
1584  if ( all_feasible )
1585  return SCIP_OKAY;
1586 
1587  if ( separated )
1588  *result = SCIP_SEPARATED;
1589 
1590  return SCIP_OKAY;
1591 }
1592 
1594 static
1595 SCIP_DECL_CONSINITPRE(consInitpreSdp)
1596 {/*lint --e{715}*/
1597  SCIP_CONSHDLRDATA* conshdlrdata;
1598 
1599  assert( conshdlr != NULL );
1600 
1601  conshdlrdata = SCIPconshdlrGetData(conshdlr);
1602  assert( conshdlrdata != NULL );
1603 
1604  conshdlrdata->neigveccuts = 0; /* this is used to give the eigenvector-cuts distinguishable names */
1605  conshdlrdata->ndiaggezerocuts = 0; /* this is used to give the diagGEzero-cuts distinguishable names */
1606  conshdlrdata->n1x1blocks = 0; /* this is used to give the lp constraints resulting from 1x1 sdp-blocks distinguishable names */
1607 
1608  return SCIP_OKAY;
1609 }
1610 
1612 static
1613 SCIP_DECL_CONSLOCK(consLockSdp)
1614 {/*lint --e{715}*/
1615  SCIP_Real* Aj;
1616  SCIP_CONSDATA* consdata;
1617  int blocksize;
1618  int var;
1619  int nvars;
1620  SCIP_Real eigenvalue;
1621 
1622  consdata = SCIPconsGetData(cons);
1623  assert( consdata != NULL );
1624  blocksize = consdata->blocksize;
1625  nvars = consdata->nvars;
1626 
1627  SCIP_CALL( SCIPallocBufferArray(scip, &Aj, blocksize * blocksize) ); /*lint !e647*/
1628 
1629  for (var = 0; var < nvars; var++)
1630  {
1631  SCIP_CALL( SCIPconsSdpGetFullAj(scip, cons, var, Aj) );
1632 
1633  /* compute the smallest eigenvalue */
1634  SCIP_CALL( SCIPlapackComputeIthEigenvalue(SCIPbuffer(scip), FALSE, blocksize, Aj, 1, &eigenvalue, NULL) );
1635  if ( SCIPisNegative(scip, eigenvalue) )
1636  {
1637  /* as the lowest eigenvalue is negative, the matrix is not positive semidefinite, so adding more of it can remove positive
1638  * semidefiniteness of the SDP-matrix */
1639  SCIP_CALL( SCIPaddVarLocks(scip, consdata->vars[var], nlocksneg, nlockspos) );
1640  }
1641 
1642  /* if the smallest eigenvalue is already positive, we don't need to compute the biggest one */
1643  if ( SCIPisPositive(scip, eigenvalue) )
1644  {
1645  /* as an eigenvalue is positive, the matrix is not negative semidefinite, so substracting more of it can remove positive
1646  * semidefiniteness of the SDP-matrix */
1647  SCIP_CALL( SCIPaddVarLocks(scip, consdata->vars[var], nlockspos, nlocksneg) );
1648  }
1649  else
1650  {
1651  /* compute the biggest eigenvalue */
1652  SCIP_CALL( SCIPlapackComputeIthEigenvalue(SCIPbuffer(scip), FALSE, blocksize, Aj, blocksize, &eigenvalue, NULL) );
1653  if ( SCIPisPositive(scip, eigenvalue) )
1654  {
1655  /* as the biggest eigenvalue is positive, the matrix is not negative semidefinite, so substracting more of it can remove positive
1656  * semidefiniteness of the SDP-matrix */
1657  SCIP_CALL( SCIPaddVarLocks(scip, consdata->vars[var], nlockspos, nlocksneg) );
1658  }
1659  }
1660  }
1661 
1662  SCIPfreeBufferArray(scip, &Aj);
1663 
1664  return SCIP_OKAY;
1665 }
1666 
1668 static
1669 SCIP_DECL_CONSEXITPRE(consExitpreSdp)
1670 {/*lint --e{715}*/
1671  assert( scip != NULL );
1672 
1673  if ( conss == NULL )
1674  return SCIP_OKAY;
1675 
1676  SCIP_CALL( fixAndAggrVars(scip, conss, nconss, TRUE) );
1677 
1678  return SCIP_OKAY;
1679 }
1680 
1682 static
1683 SCIP_DECL_CONSPRESOL(consPresolSdp)
1684 {/*lint --e{715}*/
1685  SCIP_CONSHDLRDATA* conshdlrdata;
1686 
1687  assert( conshdlr != NULL );
1688  assert( result != 0 );
1689 
1690  conshdlrdata = SCIPconshdlrGetData(conshdlr);
1691  assert( conshdlrdata != NULL );
1692 
1693  if ( nrounds == 0 )
1694  {
1695  SCIP_CALL( move_1x1_blocks_to_lp(scip, conss, nconss, naddconss, ndelconss, result) );
1696  if ( conshdlrdata->diaggezerocuts )
1697  {
1698  SCIP_CALL( diagGEzero(scip, conss, nconss, naddconss) );
1699  }
1700  if ( conshdlrdata->diagzeroimplcuts )
1701  {
1702  SCIP_CALL( diagZeroImpl(scip, conss, nconss, naddconss) );
1703  }
1704  }
1705 
1706  return SCIP_OKAY;
1707 }
1708 
1710 static
1711 SCIP_DECL_CONSTRANS(consTransSdp)
1712 {/*lint --e{715}*/
1713  SCIP_CONSDATA* sourcedata;
1714  SCIP_CONSDATA* targetdata;
1715 #ifdef OMP
1716  SCIP_CONSHDLRDATA* conshdlrdata;
1717 #endif
1718 #ifndef NDEBUG
1719  int snprintfreturn; /* used to check the return code of snprintf */
1720 #endif
1721  int i;
1722  char transname[SCIP_MAXSTRLEN];
1723 
1724  sourcedata = SCIPconsGetData(sourcecons);
1725  assert( sourcedata != NULL );
1726 
1727  SCIPdebugMessage("Transforming constraint <%s>\n", SCIPconsGetName(sourcecons));
1728 
1729 #ifdef OMP
1730  conshdlrdata = SCIPconshdlrGetData(conshdlr);
1731  SCIPdebugMessage("Setting number of threads to %d via OpenMP in Openblas.\n", conshdlrdata->nthreads);
1732  omp_set_num_threads(conshdlrdata->nthreads);
1733 #endif
1734 
1735  SCIP_CALL( SCIPallocBlockMemory(scip, &targetdata) );
1736 
1737  /* copy some general data */
1738  targetdata->nvars = sourcedata->nvars;
1739  targetdata->nnonz = sourcedata->nnonz;
1740  targetdata->blocksize = sourcedata->blocksize;
1741 
1742  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(targetdata->nvarnonz), sourcedata->nvarnonz, sourcedata->nvars) );
1743 
1744  /* copy the non-constant nonzeros */
1745  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(targetdata->col), sourcedata->nvars) );
1746  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(targetdata->row), sourcedata->nvars) );
1747  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(targetdata->val), sourcedata->nvars) );
1748 
1749  for (i = 0; i < sourcedata->nvars; i++)
1750  {
1751  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(targetdata->col[i]), sourcedata->col[i], sourcedata->nvarnonz[i]) );
1752  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(targetdata->row[i]), sourcedata->row[i], sourcedata->nvarnonz[i]) );
1753  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(targetdata->val[i]), sourcedata->val[i], sourcedata->nvarnonz[i]) );
1754  }
1755  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(targetdata->vars), sourcedata->nvars));
1756 
1757  /* copy & transform the vars array */
1758  for (i = 0; i < sourcedata->nvars; i++)
1759  {
1760  targetdata->vars[i] = SCIPvarGetTransVar(sourcedata->vars[i]);
1761  SCIP_CALL( SCIPcaptureVar(scip, targetdata->vars[i]) );
1762  }
1763 
1764  /* copy the constant nonzeros */
1765  targetdata->constnnonz = sourcedata->constnnonz;
1766 
1767  if ( sourcedata->constnnonz > 0 )
1768  {
1769  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(targetdata->constcol), sourcedata->constcol, sourcedata->constnnonz));
1770  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(targetdata->constrow), sourcedata->constrow, sourcedata->constnnonz));
1771  SCIP_CALL( SCIPduplicateBlockMemoryArray(scip, &(targetdata->constval), sourcedata->constval, sourcedata->constnnonz));
1772  }
1773  else
1774  {
1775  targetdata->constcol = NULL;
1776  targetdata->constrow = NULL;
1777  targetdata->constval = NULL;
1778  }
1779 
1780  /* copy the maxrhsentry */
1781  targetdata->maxrhsentry = sourcedata->maxrhsentry;
1782 
1783  /* name the transformed constraint */
1784 #ifndef NDEBUG
1785  snprintfreturn = SCIPsnprintf(transname, SCIP_MAXSTRLEN, "t_%s", SCIPconsGetName(sourcecons));
1786  assert( snprintfreturn < SCIP_MAXSTRLEN ); /* check whether the name fits into the string */
1787 #else
1788  (void) SCIPsnprintf(transname, SCIP_MAXSTRLEN, "t_%s", SCIPconsGetName(sourcecons));
1789 #endif
1790 
1791  /* create target constraint */
1792  SCIP_CALL( SCIPcreateCons(scip, targetcons, transname, conshdlr, targetdata,
1793  SCIPconsIsInitial(sourcecons), SCIPconsIsSeparated(sourcecons), SCIPconsIsEnforced(sourcecons),
1794  SCIPconsIsChecked(sourcecons), SCIPconsIsPropagated(sourcecons), SCIPconsIsLocal(sourcecons),
1795  SCIPconsIsModifiable(sourcecons), SCIPconsIsDynamic(sourcecons), SCIPconsIsRemovable(sourcecons),
1796  SCIPconsIsStickingAtNode(sourcecons)) );
1797 
1798  return SCIP_OKAY;
1799 }
1800 
1802 static
1803 SCIP_DECL_CONSCHECK(consCheckSdp)
1804 {/*lint --e{715}*/
1805  int i;
1806 
1807  assert( scip != NULL );
1808  assert( result != NULL );
1809 
1810  *result = SCIP_FEASIBLE;
1811 
1812  for (i = 0; i < nconss; ++i)
1813  {
1814  SCIP_CALL( SCIPconsSdpCheckSdpCons(scip, conss[i], sol, checkintegrality, checklprows, printreason, result) );
1815  if ( *result == SCIP_INFEASIBLE )
1816  return SCIP_OKAY;
1817  }
1818 
1819  return SCIP_OKAY;
1820 }
1821 
1826 static
1827 SCIP_DECL_CONSENFOPS(consEnfopsSdp)
1828 {/*lint --e{715}*/
1829  int i;
1830 
1831  assert( scip != NULL );
1832  assert( result != NULL );
1833  assert( conss != NULL );
1834 
1835  *result = SCIP_DIDNOTRUN;
1836 
1837  if ( objinfeasible )
1838  {
1839  SCIPdebugMessage("-> pseudo solution is objective infeasible, return.\n");
1840  return SCIP_OKAY;
1841  }
1842 
1843  for (i = 0; i < nconss; ++i)
1844  {
1845  SCIP_CALL( SCIPconsSdpCheckSdpCons(scip, conss[i], NULL, 0, 0, 0, result) );
1846 
1847  if (*result == SCIP_INFEASIBLE)
1848  {
1849  /* if it is infeasible for one SDP constraint, it is infeasible for the whole problem */
1850  SCIPdebugMessage("-> pseudo solution infeasible for SDP-constraint %s, return.\n", SCIPconsGetName(conss[i]));
1851  return SCIP_OKAY;
1852  }
1853  }
1854 
1855  SCIPdebugMessage("-> pseudo solution feasible for all SDP-constraints.\n");
1856 
1857  return SCIP_OKAY;
1858 }
1859 
1860 
1863 static
1864 SCIP_DECL_CONSENFOLP(consEnfolpSdp)
1865 {/*lint --e{715}*/
1866  assert( scip != NULL );
1867  assert( conshdlr != NULL );
1868  assert( conss != NULL );
1869  assert( result != NULL );
1870 
1871  return EnforceConstraint(scip, conshdlr, conss, nconss, NULL, result);
1872 }
1873 
1876 static
1877 SCIP_DECL_CONSENFORELAX(consEnforelaxSdp)
1878 {/*lint --e{715}*/
1879  assert( scip != NULL );
1880  assert( conshdlr != NULL );
1881  assert( conss != NULL );
1882  assert( result != NULL );
1883 
1884  return EnforceConstraint(scip, conshdlr, conss, nconss, sol, result);
1885 }
1886 
1888 static
1889 SCIP_DECL_CONSSEPASOL(consSepasolSdp)
1890 {/*lint --e{715}*/
1891  int i;
1892 
1893  assert( result != NULL );
1894  *result = SCIP_DIDNOTFIND;
1895 
1896  for (i = 0; i < nusefulconss; ++i)
1897  {
1898  SCIP_CALL( separateSol(scip, conshdlr, conss[i], sol, result) );
1899  }
1900 
1901  return SCIP_OKAY;
1902 }
1903 
1905 static
1906 SCIP_DECL_CONSSEPALP(consSepalpSdp)
1907 {/*lint --e{715}*/
1908  int i;
1909 
1910  assert( result != NULL );
1911  *result = SCIP_DIDNOTFIND;
1912 
1913  for (i = 0; i < nusefulconss; ++i)
1914  {
1915  SCIP_CALL( separateSol(scip, conshdlr, conss[i], NULL, result) );
1916  }
1917 
1918  return SCIP_OKAY;
1919 }
1920 
1922 static
1923 SCIP_DECL_CONSDELETE(consDeleteSdp)
1924 {/*lint --e{715}*/
1925  int i;
1926 
1927  assert( cons != NULL );
1928  assert( consdata != NULL );
1929 
1930  SCIPdebugMessage("deleting SDP constraint <%s>.\n", SCIPconsGetName(cons));
1931 
1932  for (i = 0; i < (*consdata)->nvars; i++)
1933  {
1934  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->val[i], (*consdata)->nvarnonz[i]);
1935  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->row[i], (*consdata)->nvarnonz[i]);
1936  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->col[i], (*consdata)->nvarnonz[i]);
1937  }
1938 
1939  /* release all variables */
1940  for (i = 0; i < (*consdata)->nvars; i++)
1941  {
1942  SCIP_CALL( SCIPreleaseVar(scip, &((*consdata)->vars[i])) );
1943  }
1944 
1945  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->vars, (*consdata)->nvars);
1946  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->constval, (*consdata)->constnnonz);
1947  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->constrow, (*consdata)->constnnonz);
1948  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->constcol, (*consdata)->constnnonz);
1949  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->val, (*consdata)->nvars);
1950  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->row, (*consdata)->nvars);
1951  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->col, (*consdata)->nvars);
1952  SCIPfreeBlockMemoryArrayNull(scip, &(*consdata)->nvarnonz, (*consdata)->nvars);
1953  SCIPfreeBlockMemory(scip, consdata);
1954 
1955  return SCIP_OKAY;
1956 }
1957 
1959 static
1960 SCIP_DECL_CONSFREE(consFreeSdp)
1961 {/*lint --e{715}*/
1962  SCIP_CONSHDLRDATA* conshdlrdata;
1963 
1964  conshdlrdata = SCIPconshdlrGetData(conshdlr);
1965  assert( conshdlrdata != NULL );
1966 
1967  SCIPfreeMemory(scip, &conshdlrdata);
1968  SCIPconshdlrSetData(conshdlr, NULL);
1969 
1970  return SCIP_OKAY;
1971 }
1972 
1974 static
1975 SCIP_DECL_CONSHDLRCOPY(conshdlrCopySdp)
1977  assert(scip != NULL);
1978  assert(conshdlr != NULL);
1979  assert(strcmp(SCIPconshdlrGetName(conshdlr), CONSHDLR_NAME) == 0);
1980 
1981  SCIP_CALL( SCIPincludeConshdlrSdp(scip) );
1982 
1983  *valid = TRUE;
1984 
1985  return SCIP_OKAY;
1986 }
1987 
1989 static
1990 SCIP_DECL_CONSCOPY(consCopySdp)
1991 {/*lint --e{715}*/
1992  SCIP_CONSDATA* sourcedata;
1993  SCIP_Bool success;
1994  SCIP_VAR** targetvars;
1995  SCIP_VAR* var;
1996  int i;
1997 
1998  assert( scip != NULL );
1999  assert( sourcescip != NULL );
2000  assert( sourcecons != NULL );
2001  assert( strcmp(SCIPconshdlrGetName(SCIPconsGetHdlr(sourcecons)), CONSHDLR_NAME) == 0 );
2002  assert( valid != NULL );
2003 
2004  SCIPdebugMessage("Copying SDP constraint <%s>\n", SCIPconsGetName(sourcecons));
2005 
2006  *valid = TRUE;
2007 
2008  /* as we can only map active variables, we have to make sure, that the constraint contains no fixed or (multi-)aggregated vars, after
2009  * exitpresolve (stage 6) this should always be the case, earlier than that we need to call fixAndAggrVars */
2010  if ( SCIPgetStage(sourcescip) <= SCIP_STAGE_EXITPRESOLVE )
2011  {
2012  SCIP_CALL( fixAndAggrVars(sourcescip, &sourcecons, 1, TRUE) );
2013  }
2014 
2015  sourcedata = SCIPconsGetData(sourcecons);
2016  assert( sourcedata != NULL );
2017 
2018  SCIP_CALL( SCIPallocBufferArray(scip, &targetvars, sourcedata->nvars) );
2019 
2020  /* map all variables in the constraint */
2021  for (i = 0; i < sourcedata->nvars; i++)
2022  {
2023  SCIP_CALL( SCIPgetVarCopy(sourcescip, scip, sourcedata->vars[i], &var, varmap, consmap, global, &success) );
2024  if ( success )
2025  targetvars[i] = var;
2026  else
2027  *valid = FALSE;
2028  }
2029 
2030  /* create the new constraint, using an adjusted source name if no new name was given */
2031  if ( name )
2032  {
2033 #ifndef NDEBUG
2034  int snprintfreturn; /* used to check the return code of snprintf */
2035 #endif
2036  char copyname[SCIP_MAXSTRLEN];
2037 
2038  /* name the copied constraint */
2039  #ifndef NDEBUG
2040  snprintfreturn = SCIPsnprintf(copyname, SCIP_MAXSTRLEN, "c_%s", name);
2041  assert( snprintfreturn < SCIP_MAXSTRLEN ); /* check whether the name fits into the string */
2042  #else
2043  (void) SCIPsnprintf(copyname, SCIP_MAXSTRLEN, "c_%s", name);
2044  #endif
2045  SCIP_CALL( SCIPcreateConsSdp( scip, cons, copyname, sourcedata->nvars, sourcedata->nnonz, sourcedata->blocksize, sourcedata->nvarnonz,
2046  sourcedata->col, sourcedata->row, sourcedata->val, targetvars, sourcedata->constnnonz,
2047  sourcedata->constcol, sourcedata->constrow, sourcedata->constval) );
2048  }
2049  else
2050  {
2051 #ifndef NDEBUG
2052  int snprintfreturn; /* used to check the return code of snprintf */
2053 #endif
2054  char copyname[SCIP_MAXSTRLEN];
2055 
2056  /* name the copied constraint */
2057  #ifndef NDEBUG
2058  snprintfreturn = SCIPsnprintf(copyname, SCIP_MAXSTRLEN, "c_%s", SCIPconsGetName(sourcecons));
2059  assert( snprintfreturn < SCIP_MAXSTRLEN ); /* check whether the name fits into the string */
2060  #else
2061  (void) SCIPsnprintf(copyname, SCIP_MAXSTRLEN, "c_%s", SCIPconsGetName(sourcecons));
2062  #endif
2063  SCIP_CALL( SCIPcreateConsSdp( scip, cons, SCIPconsGetName(sourcecons), sourcedata->nvars, sourcedata->nnonz, sourcedata->blocksize,
2064  sourcedata->nvarnonz, sourcedata->col, sourcedata->row, sourcedata->val, targetvars, sourcedata->constnnonz,
2065  sourcedata->constcol, sourcedata->constrow, sourcedata->constval) );
2066  }
2067 
2068  SCIPfreeBufferArray(scip, &targetvars);
2069 
2070  return SCIP_OKAY;
2071 }
2072 
2074 static
2075 SCIP_DECL_CONSPRINT(consPrintSdp)
2076 {/*lint --e{715}*/
2077 #ifdef PRINT_HUMAN_READABLE
2078  SCIP_CONSDATA* consdata;
2079  SCIP_Real* fullmatrix;
2080  int v;
2081  int i;
2082  int j;
2083 
2084  assert( scip != NULL );
2085  assert( cons != NULL );
2086 
2087  consdata = SCIPconsGetData(cons);
2088  assert( consdata != NULL );
2089 
2090  SCIP_CALL( SCIPallocBufferArray(scip, &fullmatrix, consdata->blocksize * consdata->blocksize) );
2091 
2092  /* print the non-constant matrices, for this they first have to be assembled in fullmatrix */
2093  for (v = 0; v < consdata->nvars; v++)
2094  {
2095  /* assemble the matrix */
2096 
2097  /* first initialize it with zero */
2098  for (i = 0; i < consdata->blocksize; i++)
2099  {
2100  for (j = 0; j < consdata->blocksize; j++)
2101  fullmatrix[i * consdata->blocksize + j] = 0.0; /*lint !e679*/
2102  }
2103 
2104  /* then add the nonzeros */
2105  for (i = 0; i < consdata->nvarnonz[v]; i++)
2106  {
2107  fullmatrix[consdata->row[v][i] * consdata->blocksize + consdata->col[v][i]] = consdata->val[v][i]; /* lower triangular entry */ /*lint !e679*/
2108  fullmatrix[consdata->col[v][i] * consdata->blocksize + consdata->row[v][i]] = consdata->val[v][i]; /* upper triangular entry */ /*lint !e679*/
2109  }
2110 
2111  /* print it */
2112  SCIPinfoMessage(scip, file, "+\n");
2113  for (i = 0; i < consdata->blocksize; i++)
2114  {
2115  SCIPinfoMessage(scip, file, "( ");
2116  for (j = 0; j < consdata->blocksize; j++)
2117  SCIPinfoMessage(scip, file, "%f ", fullmatrix[i * consdata->blocksize + j]); /*lint !e679*/
2118  SCIPinfoMessage(scip, file, ")\n");
2119  }
2120  SCIPinfoMessage(scip, file, "* %s\n", SCIPvarGetName(consdata->vars[v]));
2121  }
2122 
2123  /* print the constant-matrix */
2124 
2125  /* assemble the matrix */
2126 
2127  /* first initialize it with zero */
2128  for (i = 0; i < consdata->blocksize; i++)
2129  {
2130  for (j = 0; j < consdata->blocksize; j++)
2131  fullmatrix[i * consdata->blocksize + j] = 0.0; /*lint !e679*/
2132  }
2133 
2134  /* then add the nonzeros */
2135  for (i = 0; i < consdata->constnnonz; i++)
2136  {
2137  fullmatrix[consdata->constrow[i] * consdata->blocksize + consdata->constcol[i]] = consdata->constval[i]; /* lower triangular entry */ /*lint !e679*/
2138  fullmatrix[consdata->constcol[i] * consdata->blocksize + consdata->constrow[i]] = consdata->constval[i]; /* upper triangular entry */ /*lint !e679*/
2139  }
2140 
2141  /* print it */
2142  SCIPinfoMessage(scip, file, "-\n");
2143  for (i = 0; i < consdata->blocksize; i++)
2144  {
2145  SCIPinfoMessage(scip, file, "( ");
2146  for (j = 0; j < consdata->blocksize; j++)
2147  SCIPinfoMessage(scip, file, "%f ", fullmatrix[i * consdata->blocksize + j]); /*lint !e679*/
2148  SCIPinfoMessage(scip, file, ")\n");
2149  }
2150  SCIPinfoMessage(scip, file, ">=0\n");
2151 
2152  SCIPfreeBufferArray(scip, &fullmatrix);
2153 
2154  return SCIP_OKAY;
2155 #else
2156  SCIP_CONSDATA* consdata;
2157  int i;
2158  int v;
2159 
2160  assert( scip != NULL );
2161  assert( cons != NULL );
2162 
2163  consdata = SCIPconsGetData(cons);
2164 
2165  /* print blocksize */
2166  SCIPinfoMessage(scip, file, "%d\n", consdata->blocksize);
2167 
2168  /* print A_0 if it exists */
2169  if ( consdata->constnnonz > 0 )
2170  {
2171  SCIPinfoMessage(scip, file, "A_0: ");
2172 
2173  for (i = 0; i < consdata->constnnonz; i++)
2174  {
2175  SCIPinfoMessage(scip, file, "(%d,%d):%f, ", consdata->constrow[i], consdata->constcol[i], consdata->constval[i]);
2176  }
2177  SCIPinfoMessage(scip, file, "\n");
2178  }
2179 
2180  /* print other matrices */
2181  for (v = 0; v < consdata->nvars; v++)
2182  {
2183  SCIPinfoMessage(scip, file, "<%s>: ", SCIPvarGetName(consdata->vars[v]));
2184  for (i = 0; i < consdata->nvarnonz[v]; i++)
2185  {
2186  SCIPinfoMessage(scip, file, "(%d,%d):%f, ", consdata->row[v][i], consdata->col[v][i], consdata->val[v][i]);
2187  }
2188  /* if this is not the last variable, add a newline */
2189  if (v < consdata->nvars - 1)
2190  {
2191  SCIPinfoMessage(scip, file, "\n");
2192  }
2193  }
2194 
2195  return SCIP_OKAY;
2196 #endif
2197 }
2198 
2200 static
2201 SCIP_DECL_CONSPARSE(consParseSdp)
2202 { /*lint --e{715}*/
2203  SCIP_Bool parsesuccess;
2204  SCIP_CONSDATA* consdata = NULL;
2205  char* pos;
2206  int currentsize;
2207  int nvars;
2208  int i;
2209 
2210  assert( scip != NULL );
2211  assert( str != NULL );
2212 
2213  nvars = SCIPgetNVars(scip);
2214 
2215  assert( success != NULL );
2216  *success = TRUE;
2217 
2218  /* create constraint data */
2219  SCIP_CALL( SCIPallocBlockMemory(scip, &consdata) );
2220  consdata->nvars = 0;
2221  consdata->nnonz = 0;
2222  consdata->constnnonz = 0;
2223  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->nvarnonz, nvars) );
2224  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->col, nvars) );
2225  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->row, nvars) );
2226  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->val, nvars) );
2227  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->vars, nvars));
2228  consdata->constcol = NULL;
2229  consdata->constrow = NULL;
2230  consdata->constval = NULL;
2231 
2232  /* parse the blocksize */
2233  parsesuccess = SCIPstrToIntValue(str, &(consdata->blocksize), &pos);
2234  *success = *success && parsesuccess;
2235 
2236  /* skip whitespace */
2237  while( isspace((unsigned char)*pos) )
2238  pos++;
2239 
2240  /* check if there is a constant part */
2241  if ( pos[0] == 'A' && pos[1] == '_' && pos[2] == '0' )
2242  {
2243  pos += 5; /* we skip "A_0: " */
2244 
2245  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->constcol, PARSE_STARTSIZE) );
2246  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->constrow, PARSE_STARTSIZE) );
2247  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->constval, PARSE_STARTSIZE) );
2248 
2249  currentsize = PARSE_STARTSIZE;
2250 
2251  /* as long as there is another entry for the constant part, parse it */
2252  while (pos[0] == '(')
2253  {
2254  pos++; /* remove the '(' */
2255 
2256  /* check if we need to enlarge the arrays */
2257  if ( consdata->constnnonz == currentsize )
2258  {
2259  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->constcol, currentsize, PARSE_SIZEFACTOR * currentsize) );
2260  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->constrow, currentsize, PARSE_SIZEFACTOR * currentsize) );
2261  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->constval, currentsize, PARSE_SIZEFACTOR * currentsize) );
2262  currentsize *= PARSE_SIZEFACTOR;
2263  }
2264 
2265  parsesuccess = SCIPstrToIntValue(pos, &(consdata->constrow[consdata->constnnonz]), &pos);
2266  *success = *success && parsesuccess;
2267  assert( consdata->constrow[consdata->constnnonz] < consdata->blocksize );
2268  pos++; /* remove the ',' */
2269  parsesuccess = SCIPstrToIntValue(pos, &(consdata->constcol[consdata->constnnonz]), &pos);
2270  *success = *success && parsesuccess;
2271  assert( consdata->constcol[consdata->constnnonz] < consdata->blocksize );
2272  pos += 2; /* remove the "):" */
2273  parsesuccess = SCIPstrToRealValue(pos, &(consdata->constval[consdata->constnnonz]), &pos);
2274  *success = *success && parsesuccess;
2275  pos ++; /* remove the "," */
2276 
2277  /* if we got an entry in the upper triangular part, switch the entries for lower triangular */
2278  if ( consdata->constcol[consdata->constnnonz] > consdata->constrow[consdata->constnnonz] )
2279  {
2280  i = consdata->constcol[consdata->constnnonz];
2281  consdata->constcol[consdata->constnnonz] = consdata->constrow[consdata->constnnonz];
2282  consdata->constrow[consdata->constnnonz] = i;
2283  }
2284 
2285  consdata->constnnonz++;
2286 
2287  /* skip whitespace */
2288  while( isspace((unsigned char)*pos) )
2289  pos++;
2290  }
2291 
2292  /* resize the arrays to their final size */
2293  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->constcol, currentsize, consdata->constnnonz) );
2294  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->constrow, currentsize, consdata->constnnonz) );
2295  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->constval, currentsize, consdata->constnnonz) );
2296  }
2297 
2298  /* skip whitespace */
2299  while( isspace((unsigned char)*pos) )
2300  pos++;
2301 
2302  /* parse the non-constant part */
2303 
2304  /* while there is another variable, parse it */
2305  while (pos[0] == '<')
2306  {
2307  /* add the variable to consdata->vars and create the corresponding nonzero arrays */
2308  SCIP_CALL( SCIPparseVarName(scip, pos, &(consdata->vars[consdata->nvars]), &pos) );
2309  SCIP_CALL( SCIPcaptureVar(scip, consdata->vars[consdata->nvars]) );
2310 
2311  consdata->nvarnonz[consdata->nvars] = 0;
2312  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(consdata->col[consdata->nvars]), PARSE_STARTSIZE));
2313  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(consdata->row[consdata->nvars]), PARSE_STARTSIZE));
2314  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(consdata->val[consdata->nvars]), PARSE_STARTSIZE));
2315  consdata->nvars++;
2316  currentsize = PARSE_STARTSIZE;
2317 
2318  pos += 2; /* remove the ": " */
2319 
2320  /* while there is another entry, parse it */
2321  while (pos[0] == '(')
2322  {
2323  pos++; /* remove the '(' */
2324 
2325  /* check if we need to enlarge the arrays */
2326  if ( consdata->nvarnonz[consdata->nvars - 1] == currentsize )
2327  {
2328  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->col[consdata->nvars - 1], currentsize, PARSE_SIZEFACTOR * currentsize) );
2329  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->row[consdata->nvars - 1], currentsize, PARSE_SIZEFACTOR * currentsize) );
2330  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->val[consdata->nvars - 1], currentsize, PARSE_SIZEFACTOR * currentsize) );
2331  currentsize *= PARSE_SIZEFACTOR;
2332  }
2333 
2334  parsesuccess = SCIPstrToIntValue(pos, &(consdata->row[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]]), &pos);
2335  *success = *success && parsesuccess;
2336  assert( consdata->row[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]] < consdata->blocksize );
2337  pos++; /* remove the ',' */
2338  parsesuccess = SCIPstrToIntValue(pos, &(consdata->col[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]]), &pos);
2339  *success = *success && parsesuccess;
2340  assert( consdata->col[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]] < consdata->blocksize );
2341  pos += 2; /* remove the "):" */
2342  parsesuccess = SCIPstrToRealValue(pos, &(consdata->val[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]]), &pos);
2343  *success = *success && parsesuccess;
2344  pos ++; /* remove the "," */
2345 
2346  /* if we got an entry in the upper triangular part, switch the entries for lower triangular */
2347  if ( consdata->col[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]] >
2348  consdata->row[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]] )
2349  {
2350  i = consdata->col[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]];
2351  consdata->col[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]] =
2352  consdata->row[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]];
2353  consdata->row[consdata->nvars - 1][consdata->nvarnonz[consdata->nvars - 1]] = i;
2354  }
2355 
2356  consdata->nvarnonz[consdata->nvars - 1]++;
2357 
2358  /* skip whitespace */
2359  while( isspace((unsigned char)*pos) )
2360  pos++;
2361  }
2362 
2363  /* resize the arrays to their final size */
2364  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->col[consdata->nvars - 1], currentsize, consdata->nvarnonz[consdata->nvars - 1]) );
2365  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->row[consdata->nvars - 1], currentsize, consdata->nvarnonz[consdata->nvars - 1]) );
2366  SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &consdata->val[consdata->nvars - 1], currentsize, consdata->nvarnonz[consdata->nvars - 1]) );
2367 
2368  /* skip whitespace */
2369  while( isspace((unsigned char)*pos) )
2370  pos++;
2371  }
2372 
2373  /* create the constraint */
2374  SCIP_CALL( SCIPcreateCons(scip, cons, name, conshdlr, consdata, initial, separate, enforce, check, propagate, local, modifiable,
2375  dynamic, removable, stickingatnode) );
2376 
2377  /* compute maximum rhs entry for later use in the DIMACS Error Norm */
2378  SCIP_CALL( setMaxRhsEntry(*cons) );
2379 
2380 #ifdef SCIP_MORE_DEBUG
2381  SCIP_CALL( SCIPprintCons(scip, *cons, NULL) );
2382 #endif
2383 
2384  return SCIP_OKAY;
2385 }
2386 
2388 static
2389 SCIP_DECL_CONSGETVARS(consGetVarsSdp)
2390 {/*lint --e{715}*/
2391  SCIP_CONSDATA* consdata;
2392  int nvars;
2393  int i;
2394 
2395  assert( scip != NULL );
2396  assert( cons != NULL );
2397  assert( vars != NULL );
2398  assert( success != NULL );
2399  assert( varssize >= 0 );
2400 
2401  consdata = SCIPconsGetData(cons);
2402  assert( consdata != NULL );
2403 
2404  nvars = consdata->nvars;
2405 
2406  if ( nvars > varssize )
2407  {
2408  SCIPdebugMessage("consGetVarsIndicator called for array of size %d, needed size %d.\n", varssize, nvars);
2409  *success = FALSE;
2410  return SCIP_OKAY;
2411  }
2412 
2413  for (i = 0; i < nvars; i++)
2414  vars[i] = consdata->vars[i];
2415 
2416  *success = TRUE;
2417 
2418  return SCIP_OKAY;
2419 }
2420 
2422 static
2423 SCIP_DECL_CONSGETNVARS(consGetNVarsSdp)
2424 {/*lint --e{715}*/
2425  SCIP_CONSDATA* consdata;
2426 
2427  assert( scip != NULL );
2428  assert( cons != NULL );
2429  assert( nvars != NULL );
2430  assert( success != NULL );
2431 
2432  consdata = SCIPconsGetData(cons);
2433  assert( consdata != NULL );
2434 
2435  *nvars = consdata->nvars;
2436  *success = TRUE;
2437 
2438  return SCIP_OKAY;
2439 }
2440 
2442 SCIP_RETCODE SCIPincludeConshdlrSdp(
2443  SCIP* scip
2444  )
2445 {
2446  SCIP_CONSHDLR* conshdlr = NULL;
2447  SCIP_CONSHDLRDATA* conshdlrdata = NULL;
2448 
2449  assert( scip != 0 );
2450 
2451  /* allocate memory for the conshdlrdata */
2452  SCIP_CALL( SCIPallocMemory(scip, &conshdlrdata) );
2453 
2454  /* include constraint handler */
2455  SCIP_CALL( SCIPincludeConshdlrBasic(scip, &conshdlr, CONSHDLR_NAME, CONSHDLR_DESC,
2457  consEnfolpSdp, consEnfopsSdp, consCheckSdp, consLockSdp, conshdlrdata) );
2458 
2459  assert( conshdlr != NULL );
2460 
2461  /* set non-fundamental callbacks via specific setter functions */
2462  SCIP_CALL( SCIPsetConshdlrDelete(scip, conshdlr, consDeleteSdp) );
2463  SCIP_CALL( SCIPsetConshdlrFree(scip, conshdlr, consFreeSdp) );
2464  SCIP_CALL( SCIPsetConshdlrCopy(scip, conshdlr, conshdlrCopySdp, consCopySdp) );
2465  SCIP_CALL( SCIPsetConshdlrInitpre(scip, conshdlr,consInitpreSdp) );
2466  SCIP_CALL( SCIPsetConshdlrExitpre(scip, conshdlr, consExitpreSdp) );
2467  SCIP_CALL( SCIPsetConshdlrPresol(scip, conshdlr, consPresolSdp, CONSHDLR_MAXPREROUNDS, CONSHDLR_PRESOLTIMING) );
2468  SCIP_CALL( SCIPsetConshdlrSepa(scip, conshdlr, consSepalpSdp, consSepasolSdp, CONSHDLR_SEPAFREQ,
2470  SCIP_CALL( SCIPsetConshdlrEnforelax(scip, conshdlr, consEnforelaxSdp) );
2471  SCIP_CALL( SCIPsetConshdlrTrans(scip, conshdlr, consTransSdp) );
2472  SCIP_CALL( SCIPsetConshdlrPrint(scip, conshdlr, consPrintSdp) );
2473  SCIP_CALL( SCIPsetConshdlrParse(scip, conshdlr, consParseSdp) );
2474  SCIP_CALL( SCIPsetConshdlrGetVars(scip, conshdlr, consGetVarsSdp) );
2475  SCIP_CALL( SCIPsetConshdlrGetNVars(scip, conshdlr, consGetNVarsSdp) );
2476 
2477  /* add parameter */
2478 #ifdef OMP
2479  SCIP_CALL( SCIPaddIntParam(scip, "constraints/SDP/threads", "number of threads used for OpenBLAS",
2480  &(conshdlrdata->nthreads), TRUE, DEFAULT_NTHREADS, 1, INT_MAX, NULL, NULL) );
2481 #endif
2482  SCIP_CALL( SCIPaddBoolParam(scip, "constraints/SDP/diaggezerocuts",
2483  "Should linear cuts enforcing the non-negativity of diagonal entries of SDP-matrices be added?",
2484  &(conshdlrdata->diaggezerocuts), TRUE, DEFAULT_DIAGGEZEROCUTS, NULL, NULL) );
2485  SCIP_CALL( SCIPaddBoolParam(scip, "constraints/SDP/diagzeroimplcuts",
2486  "Should linear cuts enforcing the implications of diagonal entries of zero in SDP-matrices be added?",
2487  &(conshdlrdata->diagzeroimplcuts), TRUE, DEFAULT_DIAGZEROIMPLCUTS, NULL, NULL) );
2488 
2489  return SCIP_OKAY;
2490 }
2491 
2496  int i,
2497  int j
2498  )
2499 {
2500  assert( j >= 0 );
2501  assert( i >= j );
2502 
2503  return i*(i+1)/2 + j;
2504 }
2505 
2512 SCIP_RETCODE SCIPconsSdpGetData(
2513  SCIP* scip,
2514  SCIP_CONS* cons,
2515  int* nvars,
2516  int* nnonz,
2517  int* blocksize,
2518  int* arraylength,
2519  int* nvarnonz,
2521  int** col,
2522  int** row,
2523  SCIP_Real** val,
2524  SCIP_VAR** vars,
2525  int* constnnonz,
2527  int* constcol,
2528  int* constrow,
2529  SCIP_Real* constval
2530  )
2531 {
2532  SCIP_CONSDATA* consdata;
2533  const char* name;
2534  int i;
2535 
2536  assert( scip != NULL );
2537  assert( cons != NULL );
2538  assert( nvars != NULL );
2539  assert( nnonz != NULL );
2540  assert( blocksize != NULL );
2541  assert( arraylength != NULL );
2542  assert( nvarnonz != NULL );
2543  assert( col != NULL );
2544  assert( row != NULL );
2545  assert( val != NULL );
2546  assert( vars != NULL );
2547  assert( constnnonz != NULL );
2548 
2549  consdata = SCIPconsGetData(cons);
2550  name = SCIPconsGetName(cons);
2551 
2552  assert( consdata->constnnonz == 0 || ( constcol != NULL && constrow != NULL && constval != NULL ) );
2553 
2554  *nvars = consdata->nvars;
2555  *nnonz = consdata->nnonz;
2556  *blocksize = consdata->blocksize;
2557 
2558  for (i = 0; i < consdata->nvars; i++)
2559  vars[i] = consdata->vars[i];
2560 
2561  /* check that the sdp-arrays are long enough to store the information */
2562  if ( *arraylength < consdata->nvars )
2563  {
2564  SCIPdebugMessage("nvarnonz, col, row and val arrays were not long enough to store the information for cons %s, they need to be at least"
2565  "size %d, given was only length %d! \n", name, consdata->nvars, *arraylength);
2566  *arraylength = consdata->nvars;
2567  }
2568  else
2569  {
2570  for (i = 0; i < consdata->nvars; i++)
2571  {
2572  nvarnonz[i] = consdata->nvarnonz[i];
2573  /* set the pointers for each variable */
2574  col[i] = consdata->col[i];
2575  row[i] = consdata->row[i];
2576  val[i] = consdata->val[i];
2577  }
2578  }
2579 
2580  /* set the constant pointers (if a constant part exists) */
2581  if ( consdata->constnnonz > 0 )
2582  {
2583  if ( consdata->constnnonz > *constnnonz )
2584  {
2585  SCIPdebugMessage("The constant nonzeros arrays were not long enough to store the information for cons %s, they need to be at least"
2586  "size %d, given was only length %d! \n", name, consdata->constnnonz, *constnnonz);
2587  }
2588  else
2589  {
2590  for (i = 0; i < consdata->constnnonz; i++)
2591  {
2592  constcol[i] = consdata->constcol[i];
2593  constrow[i] = consdata->constrow[i];
2594  constval[i] = consdata->constval[i];
2595  }
2596  }
2597  }
2598 
2599  *constnnonz = consdata->constnnonz;
2600 
2601  return SCIP_OKAY;
2602 }
2603 
2608 SCIP_RETCODE SCIPconsSdpGetNNonz(
2609  SCIP* scip,
2610  SCIP_CONS* cons,
2611  int* nnonz,
2612  int* constnnonz
2613  )
2614 {
2615  SCIP_CONSDATA* consdata;
2616 
2617  assert( scip != NULL );
2618  assert( cons != NULL );
2619 
2620  consdata = SCIPconsGetData(cons);
2621  assert( consdata != NULL );
2622 
2623  if ( nnonz != NULL )
2624  *nnonz = consdata->nnonz;
2625 
2626  if ( constnnonz != NULL )
2627  *constnnonz = consdata->constnnonz;
2628 
2629  return SCIP_OKAY;
2630 }
2631 
2634  SCIP* scip,
2635  SCIP_CONS* cons
2636  )
2637 {
2638  SCIP_CONSDATA* consdata;
2639 
2640  assert( scip != NULL );
2641  assert( cons != NULL );
2642 
2643  consdata = SCIPconsGetData(cons);
2644  assert( consdata != NULL );
2645 
2646  return consdata->blocksize;
2647 }
2648 
2650 SCIP_RETCODE SCIPconsSdpGetFullAj(
2651  SCIP* scip,
2652  SCIP_CONS* cons,
2653  int j,
2654  SCIP_Real* Aj
2655  )
2656 {
2657  SCIP_CONSDATA* consdata;
2658  int blocksize;
2659  int i;
2660 
2661  assert( scip != NULL );
2662  assert( cons != NULL );
2663  assert( j >= 0 );
2664  assert( Aj != NULL );
2665 
2666  consdata = SCIPconsGetData(cons);
2667  assert( consdata != NULL );
2668  blocksize = consdata->blocksize;
2669 
2670  assert( j < consdata->nvars );
2671 
2672  for (i = 0; i < blocksize * blocksize; i++)
2673  Aj[i] = 0;
2674 
2675  for (i = 0; i < consdata->nvarnonz[j]; i++)
2676  {
2677  Aj[consdata->col[j][i] * blocksize + consdata->row[j][i]] = consdata->val[j][i]; /*lint !e679*/
2678  Aj[consdata->row[j][i] * blocksize + consdata->col[j][i]] = consdata->val[j][i]; /*lint !e679*/
2679  }
2680 
2681  return SCIP_OKAY;
2682 }
2683 
2685 SCIP_RETCODE SCIPconsSdpGetFullConstMatrix(
2686  SCIP* scip,
2687  SCIP_CONS* cons,
2688  SCIP_Real* mat
2689  )
2690 {
2691  SCIP_CONSDATA* consdata;
2692  int blocksize;
2693  int i;
2694  int j;
2695 
2696  assert( scip != NULL );
2697  assert( cons != NULL );
2698  assert( mat != NULL );
2699 
2700  consdata = SCIPconsGetData(cons);
2701  blocksize = consdata->blocksize;
2702 
2703  for (i = 0; i < blocksize; i++)
2704  {
2705  for (j = 0; j < blocksize; j++)
2706  mat[i * blocksize + j] = 0.0; /*lint !e679*/
2707  }
2708 
2709  for (i = 0; i < consdata->constnnonz; i++)
2710  {
2711  mat[consdata->constcol[i] * blocksize + consdata->constrow[i]] = consdata->constval[i]; /*lint !e679*/
2712  mat[consdata->constrow[i] * blocksize + consdata->constcol[i]] = consdata->constval[i]; /*lint !e679*/
2713  }
2714 
2715  return SCIP_OKAY;
2716 }
2717 
2720  SCIP* scip,
2721  SCIP_CONS* cons,
2722  SCIP_Real* mat
2723  )
2724 {
2725  SCIP_CONSDATA* consdata;
2726  int blocksize;
2727  int i;
2728 
2729  assert( scip != NULL );
2730  assert( cons != NULL );
2731  assert( mat != NULL );
2732 
2733  consdata = SCIPconsGetData(cons);
2734  assert( consdata != NULL );
2735 
2736  blocksize = consdata->blocksize;
2737 
2738  /* initialize the matrix with 0 */
2739  for (i = 0; i < (blocksize * (blocksize + 1)) / 2; i++)
2740  mat[i] = 0.0;
2741 
2742  for (i = 0; i < consdata->constnnonz; i++)
2743  mat[SCIPconsSdpCompLowerTriangPos(consdata->constrow[i], consdata->constcol[i])] = consdata->constval[i];
2744 
2745  return SCIP_OKAY;
2746 }
2747 
2758 SCIP_RETCODE SCIPconsSdpGuessInitialPoint(
2759  SCIP* scip,
2760  SCIP_CONS* cons,
2761  SCIP_Real* lambdastar
2762  )
2763 {
2764  SCIP_CONSDATA* consdata;
2765  SCIP_Real sparsity;
2766  SCIP_Real maxinfnorm;
2767  SCIP_Real maxconst;
2768  SCIP_Real mininfnorm;
2769  SCIP_Real maxobj;
2770  SCIP_Real maxbound;
2771  SCIP_Real primalguess;
2772  SCIP_Real dualguess;
2773  SCIP_Real compval;
2774  int blocksize;
2775  int i;
2776  int v;
2777 
2778  assert( scip != NULL );
2779  assert( cons != NULL );
2780  assert( lambdastar != NULL );
2781  assert( strcmp(SCIPconshdlrGetName(SCIPconsGetHdlr(cons)), CONSHDLR_NAME) == 0 );
2782 
2783  consdata = SCIPconsGetData(cons);
2784  assert( consdata != NULL );
2785 
2786  /* If there are no nonzeros, we cannot use the usual formula, since it divides through the number of nonzeros. In this case,
2787  * however, we will not solve an SDP anyways but at most an LP (more likely the problem will be solved in local presolving,
2788  * if all variables are fixed and not only those in the SDP-part), so we just take the default value of SDPA.
2789  */
2790  if ( consdata->nnonz == 0 )
2791  {
2792  *lambdastar = 100.0;
2793  return SCIP_OKAY;
2794  }
2795 
2796  blocksize = consdata->blocksize;
2797 
2798  sparsity = consdata->nnonz / (0.5 * blocksize * (blocksize + 1));
2799 
2800  /* compute the maximum entry of the A_i */
2801  maxinfnorm = 0.0;
2802  mininfnorm = SCIPinfinity(scip);
2803  for (v = 0; v < consdata->nvars; v++)
2804  {
2805  for (i = 0; i < consdata->nvarnonz[v]; i++)
2806  {
2807  if ( SCIPisGT(scip, REALABS(consdata->val[v][i]), maxinfnorm ) )
2808  maxinfnorm = REALABS(consdata->val[v][i]);
2809  if ( SCIPisLT(scip, REALABS(consdata->val[v][i]), mininfnorm) )
2810  mininfnorm = REALABS(consdata->val[v][i]);
2811  }
2812  }
2813  maxconst = 0.0;
2814  for (i = 0; i < consdata->constnnonz; i++)
2815  {
2816  if ( SCIPisGT(scip, REALABS(consdata->constval[i]), maxconst ) )
2817  maxconst = REALABS(consdata->constval[i]);
2818  }
2819 
2820  assert( SCIPisGT(scip, mininfnorm, 0.0) );
2821 
2822  /* compute maximum b_i and bound */
2823  maxobj = 0.0;
2824  maxbound = 0.0;
2825  for (v = 0; v < consdata->nvars; v++)
2826  {
2827  if ( SCIPisGT(scip, REALABS(SCIPvarGetObj(consdata->vars[v])), maxobj) )
2828  maxobj = REALABS(SCIPvarGetObj(consdata->vars[v]));
2829  compval = SCIPisInfinity(scip, REALABS(SCIPvarGetUbGlobal(consdata->vars[v]))) ? 1e+6 : REALABS(SCIPvarGetUbGlobal(consdata->vars[v]));
2830  if ( SCIPisGT(scip, compval, maxbound) )
2831  maxbound = compval;
2832  compval = SCIPisInfinity(scip, REALABS(SCIPvarGetLbGlobal(consdata->vars[v]))) ? 1e+6 : REALABS(SCIPvarGetUbGlobal(consdata->vars[v]));
2833  if ( SCIPisGT(scip, compval, maxbound) )
2834  maxbound = compval;
2835  }
2836 
2837  /* if all variables were unbounded, we set the value to 10^6 */
2838  if ( SCIPisEQ(scip, maxbound, 0.0) )
2839  maxbound = 1E+6;
2840 
2841  /* compute primal and dual guess */
2842  primalguess = maxobj / (sparsity * mininfnorm);
2843  dualguess = sparsity * maxinfnorm * maxbound + maxconst;
2844 
2845  if ( SCIPisGT(scip, primalguess, dualguess) )
2846  *lambdastar = primalguess;
2847  else
2848  *lambdastar = dualguess;
2849 
2850  return SCIP_OKAY;
2851 }
2852 
2854 SCIP_Real SCIPconsSdpGetMaxConstEntry(
2855  SCIP* scip,
2856  SCIP_CONS* cons
2857  )
2858 {
2859  SCIP_CONSDATA* consdata;
2860 
2861  assert( scip != NULL );
2862  assert( cons != NULL );
2863 
2864  consdata = SCIPconsGetData(cons);
2865 
2866  return consdata->maxrhsentry;
2867 }
2868 
2870 SCIP_Real SCIPconsSdpGetMaxSdpCoef(
2871  SCIP* scip,
2872  SCIP_CONS* cons
2873  )
2874 {
2875  SCIP_CONSDATA* consdata;
2876  SCIP_Real maxcoef;
2877  int v;
2878  int i;
2879 
2880  assert( scip != NULL );
2881  assert( cons != NULL );
2882 
2883  consdata = SCIPconsGetData(cons);
2884 
2885  maxcoef = 0.0;
2886 
2887  for (v = 0; v < consdata->nvars; v++)
2888  {
2889  for (i = 0; i < consdata->nvarnonz[v]; i++)
2890  {
2891  if ( SCIPisGT(scip, REALABS(consdata->val[v][i]), maxcoef) )
2892  maxcoef = REALABS(consdata->val[v][i]);
2893  }
2894  }
2895 
2896  return maxcoef;
2897 }
2898 
2905  SCIP* scip,
2906  SCIP_CONS* cons
2907  )
2908 {
2909  SCIP_CONSDATA* consdata;
2910  int v;
2911  int ub;
2912  int denselength;
2913 
2914  assert( cons != NULL );
2915 
2916  consdata = SCIPconsGetData(cons);
2917  assert( consdata != NULL );
2918 
2919  ub = consdata->constnnonz;
2920 
2921  for (v = 0; v < consdata->nvars; v++)
2922  ub += consdata->nvarnonz[v];
2923 
2924  denselength = consdata->blocksize * (consdata->blocksize + 1) / 2;
2925 
2926  return (ub <= denselength ? ub : denselength);
2927 }
2928 
2934  SCIP* scip,
2935  SCIP_CONS* cons,
2936  SCIP_SOL* sol,
2937  int* length,
2939  int* row,
2940  int* col,
2941  SCIP_Real* val
2942  )
2943 {
2944  SCIP_CONSDATA* consdata;
2945  int i;
2946  int v;
2947  int nnonz;
2948  SCIP_Real epsilon;
2949 
2950  assert( scip != NULL );
2951  assert( cons != NULL );
2952  assert( sol != NULL );
2953  assert( length != NULL );
2954  assert( row != NULL );
2955  assert( col != NULL );
2956  assert( val != NULL );
2957 
2958  consdata = SCIPconsGetData(cons);
2959  assert( consdata != NULL );
2960 
2961  /* initialize nnonz/row/col/val with constant arrays */
2962  nnonz = consdata->constnnonz;
2963  if ( *length < nnonz )
2964  {
2965  *length = -1;
2966  SCIPerrorMessage("Arrays not long enough in SCIPconsSdpComputeSparseSdpMatrix, length %d given, need at least %d (probably more)\n",
2967  *length, nnonz);
2968  return SCIP_ERROR;
2969  }
2970  for (i = 0; i < consdata->constnnonz; i++)
2971  {
2972  row[i] = consdata->constrow[i];
2973  col[i] = consdata->constcol[i];
2974  val[i] = -1 * consdata->constval[i];
2975  }
2976 
2977  /* add all variable arrays multiplied by corresponding solution value */
2978  SCIP_CALL( SCIPgetRealParam(scip, "numerics/epsilon", &epsilon) );
2979 
2980  for (v = 0; v < consdata->nvars; v++)
2981  {
2982  SCIP_CALL( SCIPsdpVarfixerMergeArrays(SCIPblkmem(scip), epsilon, consdata->row[v], consdata->col[v], consdata->val[v], consdata->nvarnonz[v],
2983  FALSE, SCIPgetSolVal(scip, sol, consdata->vars[v]), row, col, val, &nnonz, *length) );
2984  if ( nnonz > *length )
2985  {
2986  *length = -1;
2987  SCIPerrorMessage("Arrays not long enough in SCIPconsSdpComputeSparseSdpMatrix, length %d given, need at least %d (probably more)\n",
2988  *length, nnonz);
2989  return SCIP_ERROR;
2990  }
2991  }
2992 
2993  /* update length pointer */
2994  *length = nnonz;
2995 
2996  return SCIP_OKAY;
2997 }
2998 
3000 SCIP_RETCODE SCIPcreateConsSdp(
3001  SCIP* scip,
3002  SCIP_CONS** cons,
3003  const char* name,
3004  int nvars,
3005  int nnonz,
3006  int blocksize,
3007  int* nvarnonz,
3008  int** col,
3009  int** row,
3010  SCIP_Real** val,
3011  SCIP_VAR** vars,
3012  int constnnonz,
3013  int* constcol,
3014  int* constrow,
3015  SCIP_Real* constval
3016  )
3017 {
3018  SCIP_CONSHDLR* conshdlr;
3019  SCIP_CONSDATA* consdata = NULL;
3020  int i;
3021  int j;
3022 
3023  assert( scip != NULL );
3024  assert( cons != NULL );
3025  assert( name != NULL );
3026  assert( nvars >= 0 );
3027  assert( nnonz >= 0 );
3028  assert( blocksize >= 0 );
3029  assert( constnnonz >= 0 );
3030  assert( nvars == 0 || vars != NULL );
3031  assert( nnonz == 0 || (nvarnonz != NULL && col != NULL && row != NULL && val != NULL ));
3032  assert( constnnonz == 0 || (constcol != NULL && constrow != NULL && constval != NULL ));
3033 
3034  conshdlr = SCIPfindConshdlr(scip, "SDP");
3035  if ( conshdlr == NULL )
3036  {
3037  SCIPerrorMessage("SDP constraint handler not found\n");
3038  return SCIP_PLUGINNOTFOUND;
3039  }
3040 
3041  /* create constraint data */
3042  SCIP_CALL( SCIPallocBlockMemory(scip, &consdata) );
3043  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->nvarnonz, nvars) );
3044  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->col, nvars) );
3045  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->row, nvars) );
3046  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->val, nvars) );
3047  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->constcol, constnnonz) );
3048  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->constrow, constnnonz) );
3049  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->constval, constnnonz) );
3050  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->vars, nvars));
3051 
3052  for (i = 0; i < nvars; i++)
3053  {
3054  assert( nvarnonz[i] >= 0 );
3055 
3056  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->col[i], nvarnonz[i]));
3057  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->row[i], nvarnonz[i]));
3058  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &consdata->val[i], nvarnonz[i]));
3059  }
3060 
3061  consdata->nvars = nvars;
3062  consdata->nnonz = nnonz;
3063  consdata->constnnonz = constnnonz;
3064  consdata->blocksize = blocksize;
3065 
3066  for (i = 0; i < nvars; i++)
3067  {
3068  consdata->nvarnonz[i] = nvarnonz[i];
3069 
3070  for (j = 0; j < nvarnonz[i]; j++)
3071  {
3072  assert( col[i][j] >= 0 );
3073  assert( col[i][j] < blocksize );
3074  assert( row[i][j] >= 0 );
3075  assert( row[i][j] < blocksize );
3076 
3077  consdata->col[i][j] = col[i][j];
3078  consdata->row[i][j] = row[i][j];
3079  consdata->val[i][j] = val[i][j];
3080  }
3081  }
3082 
3083  for (i = 0; i < constnnonz; i++)
3084  {
3085  consdata->constcol[i] = constcol[i];
3086  consdata->constrow[i] = constrow[i];
3087  consdata->constval[i] = constval[i];
3088  }
3089 
3090  for (i = 0; i < nvars; i++)
3091  {
3092  consdata->vars[i] = vars[i];
3093  SCIP_CALL( SCIPcaptureVar(scip, consdata->vars[i]) );
3094  }
3095  SCIPdebugMessage("creating cons %s\n", name);
3096  /* create constraint */
3097  SCIP_CALL( SCIPcreateCons(scip, cons, name, conshdlr, consdata, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE, FALSE, FALSE, FALSE, FALSE) );
3098 
3099  /* compute maximum rhs entry for later use in the DIMACS Error Norm */
3100  SCIP_CALL( setMaxRhsEntry(*cons) );
3101 
3102  return SCIP_OKAY;
3103 }
SCIP_Real SCIPconsSdpGetMaxConstEntry(SCIP *scip, SCIP_CONS *cons)
Definition: cons_sdp.c:2855
static SCIP_DECL_CONSTRANS(consTransSdp)
Definition: cons_sdp.c:1712
EXTERN SCIP_RETCODE SCIPlapackMatrixVectorMult(int nrows, int ncols, SCIP_Real *matrix, SCIP_Real *vector, SCIP_Real *result)
Definition: lapack_dsdp.c:321
EXTERN SCIP_RETCODE SCIPlapackComputeIthEigenvalue(BMS_BUFMEM *bufmem, SCIP_Bool geteigenvectors, int n, SCIP_Real *A, int i, SCIP_Real *eigenvalue, SCIP_Real *eigenvector)
static SCIP_RETCODE separateSol(SCIP *scip, SCIP_CONSHDLR *conshdlr, SCIP_CONS *cons, SCIP_SOL *sol, SCIP_RESULT *result)
Definition: cons_sdp.c:408
SCIP_RETCODE SCIPconsSdpGetFullAj(SCIP *scip, SCIP_CONS *cons, int j, SCIP_Real *Aj)
Definition: cons_sdp.c:2651
void SCIPsdpVarfixerSortRowCol(int *row, int *col, SCIP_Real *val, int length)
Definition: SdpVarfixer.c:50
#define CONSHDLR_PRESOLTIMING
Definition: cons_sdp.c:82
#define CONSHDLR_NAME
Definition: cons_sdp.c:69
SCIP_RETCODE SCIPincludeConshdlrSdp(SCIP *scip)
Definition: cons_sdp.c:2443
static SCIP_RETCODE multiplyConstraintMatrix(SCIP_CONS *cons, int j, SCIP_Real *v, SCIP_Real *vAv)
Definition: cons_sdp.c:199
#define PARSE_SIZEFACTOR
Definition: cons_sdp.c:84
static SCIP_DECL_CONSENFORELAX(consEnforelaxSdp)
Definition: cons_sdp.c:1878
#define CONSHDLR_NEEDSCONS
Definition: cons_sdp.c:80
static SCIP_RETCODE diagGEzero(SCIP *scip, SCIP_CONS **conss, int nconss, int *naddconss)
Definition: cons_sdp.c:498
static SCIP_DECL_CONSDELETE(consDeleteSdp)
Definition: cons_sdp.c:1924
#define CONSHDLR_SEPAFREQ
Definition: cons_sdp.c:74
#define CONSHDLR_DELAYSEPA
Definition: cons_sdp.c:79
static SCIP_RETCODE computeSdpMatrix(SCIP *scip, SCIP_CONS *cons, SCIP_SOL *y, SCIP_Real *matrix)
Definition: cons_sdp.c:157
#define CONSHDLR_EAGERFREQ
Definition: cons_sdp.c:75
static SCIP_RETCODE move_1x1_blocks_to_lp(SCIP *scip, SCIP_CONS **conss, int nconss, int *naddconss, int *ndelconss, SCIP_RESULT *result)
Definition: cons_sdp.c:855
#define DEFAULT_DIAGGEZEROCUTS
Definition: cons_sdp.c:85
SCIP_RETCODE SCIPconsSdpComputeSparseSdpMatrix(SCIP *scip, SCIP_CONS *cons, SCIP_SOL *sol, int *length, int *row, int *col, SCIP_Real *val)
Definition: cons_sdp.c:2934
SCIP_RETCODE SCIPconsSdpGetLowerTriangConstMatrix(SCIP *scip, SCIP_CONS *cons, SCIP_Real *mat)
Definition: cons_sdp.c:2720
static SCIP_DECL_CONSSEPALP(consSepalpSdp)
Definition: cons_sdp.c:1907
static SCIP_RETCODE setMaxRhsEntry(SCIP_CONS *cons)
Definition: cons_sdp.c:239
static SCIP_RETCODE multiaggrVar(SCIP *scip, SCIP_CONS *cons, int *v, SCIP_VAR **aggrvars, SCIP_Real *scalars, int naggrvars, SCIP_Real constant, int *savedcol, int *savedrow, SCIP_Real *savedval, int *nfixednonz, int *vararraylength)
Definition: cons_sdp.c:1060
static SCIP_DECL_CONSCHECK(consCheckSdp)
Definition: cons_sdp.c:1804
#define CONSHDLR_DESC
Definition: cons_sdp.c:70
Constraint handler for SDP-constraints.
static SCIP_RETCODE expandSymMatrix(int size, SCIP_Real *symMat, SCIP_Real *fullMat)
Definition: cons_sdp.c:124
maps SCIP variables to SDP indices (the SCIP variables are given SDP indices in the order in which th...
SCIP_RETCODE SCIPconsSdpGetFullConstMatrix(SCIP *scip, SCIP_CONS *cons, SCIP_Real *mat)
Definition: cons_sdp.c:2686
SCIP_RETCODE SCIPcreateConsSdp(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, int nnonz, int blocksize, int *nvarnonz, int **col, int **row, SCIP_Real **val, SCIP_VAR **vars, int constnnonz, int *constcol, int *constrow, SCIP_Real *constval)
Definition: cons_sdp.c:3001
SCIP_Real SCIPconsSdpGetMaxSdpCoef(SCIP *scip, SCIP_CONS *cons)
Definition: cons_sdp.c:2871
SCIP_RETCODE SCIPconsSdpGetData(SCIP *scip, SCIP_CONS *cons, int *nvars, int *nnonz, int *blocksize, int *arraylength, int *nvarnonz, int **col, int **row, SCIP_Real **val, SCIP_VAR **vars, int *constnnonz, int *constcol, int *constrow, SCIP_Real *constval)
Definition: cons_sdp.c:2513
int SCIPconsSdpCompLowerTriangPos(int i, int j)
Definition: cons_sdp.c:2496
#define CONSHDLR_CHECKPRIORITY
Definition: cons_sdp.c:73
static SCIP_DECL_CONSENFOPS(consEnfopsSdp)
Definition: cons_sdp.c:1828
#define DEFAULT_DIAGZEROIMPLCUTS
Definition: cons_sdp.c:86
adds the main functionality to fix/unfix/(multi-)aggregate variables by merging two three-tuple-array...
#define CONSHDLR_ENFOPRIORITY
Definition: cons_sdp.c:72
static SCIP_DECL_CONSFREE(consFreeSdp)
Definition: cons_sdp.c:1961
static SCIP_DECL_CONSPRESOL(consPresolSdp)
Definition: cons_sdp.c:1684
static SCIP_DECL_CONSEXITPRE(consExitpreSdp)
Definition: cons_sdp.c:1670
static SCIP_DECL_CONSPRINT(consPrintSdp)
Definition: cons_sdp.c:2076
static SCIP_DECL_CONSHDLRCOPY(conshdlrCopySdp)
Definition: cons_sdp.c:1976
#define CONSHDLR_MAXPREROUNDS
Definition: cons_sdp.c:78
static SCIP_DECL_CONSINITPRE(consInitpreSdp)
Definition: cons_sdp.c:1596
int SCIPconsSdpComputeUbSparseSdpMatrixLength(SCIP *scip, SCIP_CONS *cons)
Definition: cons_sdp.c:2905
static SCIP_DECL_CONSPARSE(consParseSdp)
Definition: cons_sdp.c:2202
SCIP_RETCODE SCIPconsSdpCheckSdpCons(SCIP *scip, SCIP_CONS *cons, SCIP_SOL *sol, SCIP_Bool checkintegrality, SCIP_Bool checklprows, SCIP_Bool printreason, SCIP_RESULT *result)
Definition: cons_sdp.c:341
static SCIP_RETCODE diagZeroImpl(SCIP *scip, SCIP_CONS **conss, int nconss, int *naddconss)
Definition: cons_sdp.c:614
#define PARSE_STARTSIZE
Definition: cons_sdp.c:83
int SCIPconsSdpGetBlocksize(SCIP *scip, SCIP_CONS *cons)
Definition: cons_sdp.c:2634
static SCIP_DECL_CONSSEPASOL(consSepasolSdp)
Definition: cons_sdp.c:1890
SCIP_RETCODE SCIPsdpVarfixerMergeArrays(BMS_BLKMEM *blkmem, SCIP_Real epsilon, int *originrow, int *origincol, SCIP_Real *originval, int originlength, SCIP_Bool originsorted, SCIP_Real scalar, int *targetrow, int *targetcol, SCIP_Real *targetval, int *targetlength, int targetmemory)
Definition: SdpVarfixer.c:88
static SCIP_RETCODE EnforceConstraint(SCIP *scip, SCIP_CONSHDLR *conshdlr, SCIP_CONS **conss, int nconss, SCIP_SOL *sol, SCIP_RESULT *result)
Definition: cons_sdp.c:1490
static SCIP_DECL_CONSGETVARS(consGetVarsSdp)
Definition: cons_sdp.c:2390
#define CONSHDLR_SEPAPRIORITY
Definition: cons_sdp.c:71
SCIP_RETCODE SCIPconsSdpGetNNonz(SCIP *scip, SCIP_CONS *cons, int *nnonz, int *constnnonz)
Definition: cons_sdp.c:2609
static SCIP_RETCODE fixAndAggrVars(SCIP *scip, SCIP_CONS **conss, int nconss, SCIP_Bool aggregate)
Definition: cons_sdp.c:1267
static SCIP_DECL_CONSLOCK(consLockSdp)
Definition: cons_sdp.c:1614
char name[SCIP_MAXSTRLEN]
static SCIP_DECL_CONSGETNVARS(consGetNVarsSdp)
Definition: cons_sdp.c:2424
static SCIP_DECL_CONSCOPY(consCopySdp)
Definition: cons_sdp.c:1991
interface methods for eigenvector computation and matrix multiplication using different versions of L...
SCIP_RETCODE SCIPconsSdpGuessInitialPoint(SCIP *scip, SCIP_CONS *cons, SCIP_Real *lambdastar)
Definition: cons_sdp.c:2759
static SCIP_RETCODE cutUsingEigenvector(SCIP *scip, SCIP_CONS *cons, SCIP_SOL *sol, SCIP_Real *coeff, SCIP_Real *lhs)
Definition: cons_sdp.c:272
static SCIP_DECL_CONSENFOLP(consEnfolpSdp)
Definition: cons_sdp.c:1865