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