Scippy

SCIP

Solving Constraint Integer Programs

probdata_cyc.c
Go to the documentation of this file.
1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2 /* */
3 /* This file is part of the program and library */
4 /* SCIP --- Solving Constraint Integer Programs */
5 /* */
6 /* Copyright (C) 2002-2020 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not visit scipopt.org. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file probdata_cyc.c
17  * @brief problem data for cycle clustering problem
18  * @author Leon Eifler
19  *
20  * This file implements the problem data for the cycle clustering problem.
21  *
22  * The problem data contains original transition matrix, the scaling parameter that appears in the objective function,
23  * and all variables that appear in the problem.
24  */
25 
26 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
27 #include "probdata_cyc.h"
28 
29 #include "scip/cons_linear.h"
30 #include "scip/cons_logicor.h"
31 #include "scip/cons_quadratic.h"
32 #include "scip/var.h"
33 #include <assert.h>
34 
35 struct SCIP_ProbData
36 {
37  SCIP_VAR**** edgevars; /**< variables for the edges (pairs of states) inside and between consecutive clusters */
38  SCIP_DIGRAPH* edgegraph; /**< digraph which tells us which variables are actually created */
39  SCIP_VAR*** binvars; /**< variable-matrix belonging to the state(bin)-cluster assignment */
40  SCIP_Real** cmatrix; /**< matrix to save the transition matrix */
41  SCIP_Real scale; /**< the weight-factor for the coherence in the objective function */
42  char model_variant; /**< the model that is used. w for weighted objective, e for normal edge-representation */
43  int nbins; /**< the number of states */
44  int ncluster; /**< the number of clusters */
45 };
46 
47 /** Check if the clustering has exactly one state in every cluster. */
49  SCIP* scip, /**< SCIP data structure */
50  SCIP_Real** solclustering, /**< matrix with the clustering */
51  int nbins, /**< the number of bins */
52  int ncluster /**< the number of clusters */
53  )
54 {
55  int i;
56  int j;
57 
58  /* check if the assignment violates paritioning, e.g. because we are in a subscip */
59  for( i = 0; i < nbins; ++i )
60  {
61  SCIP_Real sum = 0.0;
62 
63  for( j = 0; j < ncluster; ++j )
64  {
65  if( !SCIPisIntegral(scip, solclustering[i][j]) )
66  return FALSE;
67  if( !SCIPisZero(scip, solclustering[i][j]) )
68  sum += solclustering[i][j];
69  }
70  if( !SCIPisEQ(scip, sum, 1.0) )
71  return FALSE;
72 
73  }
74 
75  return TRUE;
76 }
77 
78 /** Assign the variables in scip according to the found clustering. */
80  SCIP* scip, /**< SCIP data structure */
81  SCIP_SOL* sol, /**< the SCIP solution */
82  SCIP_Real** clustering, /**< the matrix with the clusterassignment */
83  int nbins, /**< the number of bins */
84  int ncluster /**< the number of cluster */
85  )
86 {
87  SCIP_VAR* var;
88  SCIP_VAR*** binvars;
89  SCIP_VAR**** edgevars;
90  int i;
91  int j;
92  int c;
93 
94  binvars = SCIPcycGetBinvars(scip);
95  edgevars = SCIPcycGetEdgevars(scip);
96 
97  assert(nbins > 0 && ncluster > 0);
98  assert(binvars != NULL);
99  assert(edgevars != NULL);
100 
101  for ( c = 0; c < ncluster; ++c )
102  {
103  /* set values of state-variables */
104  for ( i = 0; i < nbins; ++i )
105  {
106  if( NULL != binvars[i][c] )
107  {
108  if( SCIPvarIsTransformed(binvars[i][c]) )
109  var = binvars[i][c];
110  else
111  var = SCIPvarGetTransVar(binvars[i][c] );
112 
113  /* check if the clusterassignment is feasible for the variable bounds. If not do not assign the variable */
114  if( var != NULL && SCIPisLE(scip, SCIPvarGetLbGlobal(var), clustering[i][c]) &&
115  SCIPisGE(scip, SCIPvarGetUbGlobal(var), clustering[i][c]) &&
117  {
118  SCIP_CALL( SCIPsetSolVal( scip, sol, var, clustering[i][c]) );
119  }
120 
121  assert(SCIPisIntegral(scip, clustering[i][c]));
122  }
123  }
124 
125  /* set the value for the edgevariables for each pair of states */
126  for( i = 0; i < nbins; ++i )
127  {
128  for( j = 0; j < i; ++j )
129  {
130  if( NULL == edgevars[i][j] || NULL == edgevars[j][i])
131  continue;
132 
133  /* check if bins are in same cluster */
134  if( SCIPisEQ(scip, 1.0, clustering[i][c] * clustering[j][c]) )
135  {
136  var = edgevars[i][j][0];
137  if( var != NULL && SCIPisGE(scip, SCIPvarGetUbGlobal(var), clustering[j][c] * clustering[i][c])
139  {
140  SCIP_CALL( SCIPsetSolVal( scip, sol, var, 1.0 ) );
141  }
142  }
143 
144  /* check if bins are in consecutive clusters */
145  else if( SCIPisEQ(scip, 1.0, clustering[i][c] * clustering[j][phi(c, ncluster)]) )
146  {
147  var = edgevars[i][j][1];
148  if( var != NULL && SCIPvarGetStatus(var) != SCIP_VARSTATUS_MULTAGGR &&
149  SCIPisGE(scip, SCIPvarGetUbGlobal(var), clustering[j][phi(c, ncluster)] * clustering[i][c]) )
150  {
151  SCIP_CALL( SCIPsetSolVal( scip, sol, var, 1.0 ) );
152  }
153  }
154 
155  else if( SCIPisEQ(scip, 1.0, clustering[j][c] * clustering[i][phi(c, ncluster)]) )
156  {
157  var = edgevars[j][i][1];
158  if( var != NULL && SCIPvarGetStatus(var) != SCIP_VARSTATUS_MULTAGGR &&
159  SCIPisGE(scip, SCIPvarGetUbGlobal(var), clustering[j][c] * clustering[i][phi(c, ncluster)]) )
160  {
161  SCIP_CALL( SCIPsetSolVal( scip, sol, var, 1.0 ) );
162  }
163  }
164  }
165  }
166  }
167 
168  return SCIP_OKAY;
169 }
170 
171 /** function that returns the successive cluster along the cycle */
172 int phi(
173  int k, /**< the cluster */
174  int ncluster /**< the number of clusters*/
175  )
176 {
177  assert(k < ncluster && k >= 0);
178  assert(ncluster > 0);
179 
180  return (k+1) % ncluster;
181 }
182 
183 /** function that returns the predecessor-cluster along the cycle */
184 int phiinv(
185  int k, /**< the cluster */
186  int ncluster /**< the number of clusters */
187  )
188 {
189  assert(k < ncluster && k >= 0);
190  assert(ncluster > 0);
191 
192  if( k - 1 < 0 )
193  return ncluster - 1;
194  else
195  return k - 1;
196 }
197 
198 /** creates all the variables for the problem. The constraints are added later, depending on the model that is used */
199 static
201  SCIP* scip, /**< SCIP data Structure */
202  SCIP_PROBDATA* probdata /**< the problem data */
203  )
204 {
205  int i;
206  int j;
207  int c;
208  int edgetype;
209  int nbins = probdata->nbins;
210  int ncluster = probdata->ncluster;
211  char varname[SCIP_MAXSTRLEN];
212 
213  /* create variables for bins */
215  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(probdata->binvars), nbins) );
216 
217  for( i = 0; i < nbins; ++i )
218  {
219  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(probdata->binvars[i]), ncluster) ); /*lint !e866*/
220  }
221 
222  for( i = 0; i < nbins; ++i )
223  {
224  for( c = 0; c < ncluster; ++c )
225  {
226  (void)SCIPsnprintf(varname, SCIP_MAXSTRLEN, "x_%d_%d", i, c);
227  SCIP_CALL( SCIPcreateVarBasic(scip, &probdata->binvars[i][c], varname, 0.0, 1.0, 0.0, SCIP_VARTYPE_BINARY) );
228  SCIP_CALL( SCIPaddVar(scip, probdata->binvars[i][c]) );
229  }
230  }
231 
232  /* create variables for the edges in each cluster combination. Index 0 are edges within cluster, 1 edges between
233  * consequtive clusters and 2 edges between non-consequtive clusters
234  */
235  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(probdata->edgevars), nbins) );
236 
237  for( i = 0; i < nbins; ++i )
238  {
239  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &(probdata->edgevars[i]), nbins) ); /*lint !e866*/
240 
241  for( j = 0; j < nbins; ++j )
242  {
243  if( i == j || (SCIPisZero(scip, (probdata->cmatrix[i][j] - probdata->cmatrix[j][i]))
244  && SCIPisZero(scip, (probdata->cmatrix[i][j] + probdata->cmatrix[j][i]) )) )
245  continue;
246 
247  SCIP_CALL( SCIPdigraphAddArc(probdata->edgegraph, i, j, NULL) );
248 
249  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &(probdata->edgevars[i][j]), 3) ); /*lint !e866*/
250 
251  for( edgetype = 0; edgetype < 3; ++edgetype )
252  {
253  if( edgetype == 0 && i < j )
254  continue;
255 
256  (void)SCIPsnprintf(varname, SCIP_MAXSTRLEN, "y_%d_%d_%d", i, j, edgetype);
257  SCIP_CALL( SCIPcreateVarBasic(scip, &probdata->edgevars[i][j][edgetype], varname,
258  0.0, 1.0, 0.0, SCIP_VARTYPE_BINARY) );
259  SCIP_CALL( SCIPaddVar(scip, probdata->edgevars[i][j][edgetype]) );
260  }
261  }
262  }
263 
264  return SCIP_OKAY;
265 }
266 
267 /** create the problem without variable amount of clusters, use simpler non-facet-defining inequalities */
268 static
270  SCIP* scip, /**< SCIP Data Structure */
271  SCIP_PROBDATA* probdata /**< the problem data */
272  )
273 {
274  int i;
275  int j;
276  int c1;
277  int c2;
278  char consname[SCIP_MAXSTRLEN];
279  SCIP_CONS* temp;
280  int nbins = probdata->nbins;
281  int ncluster = probdata->ncluster;
282  SCIP_Real scale;
283 
284  SCIP_CALL( SCIPgetRealParam(scip, "cycleclustering/scale_coherence", &scale) );
285  probdata->scale = scale;
286 
287  /* create constraints */
288 
289  /* create the set-partitioning constraints of the bins */
290  for( i = 0; i < nbins; ++i )
291  {
292  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "setpart_%d", i+1 );
293  SCIP_CALL( SCIPcreateConsSetpart(scip, &temp, consname, 0, NULL, TRUE, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE,
294  FALSE, FALSE, FALSE) );
295 
296  for ( c1 = 0; c1 < ncluster; ++c1 )
297  {
298  SCIP_CALL( SCIPaddCoefSetppc(scip, temp, probdata->binvars[i][c1]) );
299  }
300  SCIP_CALL( SCIPaddCons(scip, temp) );
301  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
302  }
303 
304  /* create constraints for the edge-cut variables */
305  SCIPinfoMessage(scip, NULL, "Using edge-representation with simplified structure. Fixed amount of cluster. \n");
306  for( i = 0; i < nbins; ++i )
307  {
308  for( j = 0; j < i; ++j )
309  {
310  if( probdata->edgevars[i][j] == NULL )
311  continue;
312 
313  /* these edges are not within a cluster */
314  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[i][j][0],
315  (probdata->cmatrix[i][j] + probdata->cmatrix[j][i]) * scale) );
316 
317  /* these are the edges that are between consequtive clusters */
318  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[i][j][1], (probdata->cmatrix[i][j] - probdata->cmatrix[j][i])) );
319  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[j][i][1], (probdata->cmatrix[j][i] - probdata->cmatrix[i][j])) );
320 
321  /* create constraints that determine when the edge-variables have to be non-zero */
322  for( c1 = 0; c1 < ncluster; ++c1 )
323  {
324  /* constraints for edges within clusters */
325  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bins_%d_%d_incluster_%d", i+1, j+1, c1+1 );
326  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
327  TRUE, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE, FALSE, FALSE, FALSE ) );
328 
329  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][0], -1.0) );
330  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][c1], 1.0) );
331  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][c1], 1.0) );
332 
333  SCIP_CALL( SCIPaddCons(scip, temp) );
334  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
335 
336  /* constraints for edges between clusters */
337  for( c2 = 0; c2 < ncluster; ++c2 )
338  {
339  SCIP_VAR* var;
340 
341  if( c2 == c1 )
342  continue;
343 
344  if( c2 == c1 + 1 || ( c2 == 0 && c1 == ncluster -1) )
345  var = probdata->edgevars[i][j][1];
346  else if( c2 == c1 - 1 || ( c1 == 0 && c2 == ncluster -1) )
347  var = probdata->edgevars[j][i][1];
348  else
349  var = probdata->edgevars[i][j][2];
350 
351  /* if two bins are in a different cluster -> the corresponding edge must be cut */
352  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bins_%d_%d_inclusters_%d_%d", i+1, j+1, c1+1, c2+1 );
353  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
354  TRUE, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE, FALSE, FALSE, FALSE ) );
355 
356  SCIP_CALL( SCIPaddCoefLinear(scip, temp, var, -1.0) );
357  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][c1], 1.0) );
358  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][c2], 1.0) );
359 
360  SCIP_CALL( SCIPaddCons(scip, temp) );
361  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
362  }
363  }
364  }
365  }
366 
367  /* only one cluster-pair at the same time can be active for an edge*/
368  for( i = 0; i < nbins; ++i )
369  {
370  for( j = 0; j < i; ++j )
371  {
372  if( NULL == probdata->edgevars[i][j] )
373  continue;
374 
375  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "sumedge_%d_%d", i+1, j+1 );
376  SCIP_CALL( SCIPcreateConsBasicLinear(scip, &temp, consname, 0, NULL, NULL, 1.0, 1.0 ) );
377 
378  for( c1 = 0; c1 < 3; ++c1 )
379  {
380  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][c1], 1.0) );
381  }
382  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][i][1], 1.0) );
383 
384  SCIP_CALL( SCIPaddCons(scip, temp) );
385  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
386  }
387  }
388 
389  /* add constraint that ensures that each cluster is used */
390  for( c1 = 0; c1 < ncluster; ++c1 )
391  {
392  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "cluster_%d_used", c1+1 );
393  SCIP_CALL( SCIPcreateConsBasicLogicor(scip, &temp, consname, 0, NULL) );
394 
395  for( i = 0; i < nbins; ++i )
396  {
397  SCIP_CALL( SCIPaddCoefLogicor(scip, temp, probdata->binvars[i][c1]) );
398  }
399 
400  SCIP_CALL( SCIPaddCons(scip, temp) );
401  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
402  }
403 
404  return SCIP_OKAY;
405 }
406 
407 /** create the problem without variable amount of clusters, using three edge-variables for each pair of states.
408  * This is the tested default version.
409  */
410 static
412  SCIP* scip, /**< SCIP Data Structure */
413  SCIP_PROBDATA* probdata /**< The problem data */
414  )
415 {
416  int i;
417  int j;
418  int c1;
419  int nbins = probdata->nbins;
420  int ncluster = probdata->ncluster;
421  char consname[SCIP_MAXSTRLEN];
422  SCIP_CONS* temp;
423  SCIP_Real scale;
424 
425  SCIP_CALL( SCIPgetRealParam(scip, "cycleclustering/scale_coherence", &scale) );
426  probdata->scale = scale;
427 
428  /* create constraints */
429 
430  /* create the set-partitioning constraints of the bins */
431  for( i = 0; i < nbins; ++i )
432  {
433  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "setpart_%d", i+1);
434  SCIP_CALL( SCIPcreateConsSetpart(scip, &temp, consname, 0, NULL, TRUE, TRUE, TRUE, TRUE, TRUE,
435  FALSE, FALSE, FALSE, FALSE, FALSE) );
436 
437  for ( c1 = 0; c1 < ncluster; ++c1 )
438  {
439  SCIP_CALL( SCIPaddCoefSetppc(scip, temp, probdata->binvars[i][c1]) );
440  }
441 
442  SCIP_CALL( SCIPaddCons(scip, temp) );
443  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
444  }
445 
446  /* create constraints for the edge-variables */
448  "Using edge-representation with simplified structure. No variable amount of cluster. \n");
449 
450  for( i = 0; i < nbins; ++i )
451  {
452  for( j = 0; j < i; ++j )
453  {
454  if( NULL == probdata->edgevars[i][j] )
455  continue;
456 
457  /* the general formulation is needed if there are more than 3 clusters. In the case of three clusters the
458  * formulation is simplified
459  */
460  if( ncluster > 3 )
461  {
462  /* these edges are within a cluster */
463  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[i][j][0],
464  (probdata->cmatrix[i][j] + probdata->cmatrix[j][i]) * scale) );
465 
466  /* these are the edges that are between consequtive clusters */
467  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[i][j][1],
468  (probdata->cmatrix[i][j] - probdata->cmatrix[j][i])) );
469  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[j][i][1],
470  (probdata->cmatrix[j][i] - probdata->cmatrix[i][j])) );
471 
472  /* create constraints that determine when the edge-variables have to be non-zero*/
473  for( c1 = 0; c1 < ncluster; ++c1 )
474  {
475  /* constraints for edges within clusters and between clusters*/
476  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bins_%d_%d_incluster_%d", i, j, c1);
477  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
479 
480  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][0], -1.0) );
481  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][c1], 1.0) );
482  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][c1], 1.0) );
483  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][1], 1.0) );
484  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][phiinv(c1, ncluster)], -1.0) );
485  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][phi(c1, ncluster)], -1.0) );
486 
487  SCIP_CALL( SCIPaddCons(scip, temp) );
488  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
489 
490  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bins_%d_%d_incluster_part2_%d", i, j, c1);
491  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
493 
494  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][0], -1.0) );
495  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][c1], 1.0) );
496  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][c1], 1.0) );
497  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][i][1], 1.0) );
498  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][phi(c1, ncluster)], -1.0) );
499  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][phiinv(c1, ncluster)], -1.0) );
500 
501  SCIP_CALL( SCIPaddCons(scip, temp) );
502  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
503 
504  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bins_%d_%d_incluster_%d_%d", i, j, c1, phi(c1, ncluster));
505  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
507 
508  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][1], -1.0) );
509  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][c1], 1.0) );
510  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][phi(c1, ncluster)], 1.0) );
511  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][0], 1.0) );
512  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][phi(c1, ncluster)], -1.0) );
513  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][c1], -1.0) );
514 
515  SCIP_CALL( SCIPaddCons(scip, temp) );
516  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
517 
518  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bins_%d_%d_incluster_%d_%d", i, j, c1, phiinv(c1, ncluster));
519  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
521 
522  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][i][1], -1.0) );
523  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][c1], 1.0) );
524  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][phiinv(c1,ncluster)], 1.0) );
525  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][0], 1.0) );
526  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][phiinv(c1, ncluster)], -1.0) );
527  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][c1], -1.0) );
528 
529  SCIP_CALL( SCIPaddCons(scip, temp) );
530  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
531  }
532  }
533  /* some variables become obsolete with three clusters */
534  else
535  {
536  /* these are the edges that are between consequtive clusters */
537  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[i][j][1],
538  (probdata->cmatrix[i][j] - probdata->cmatrix[j][i])
539  - (probdata->cmatrix[i][j] + probdata->cmatrix[j][i]) * scale) );
540  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[j][i][1],
541  (probdata->cmatrix[j][i] - probdata->cmatrix[i][j])
542  - (probdata->cmatrix[i][j] + probdata->cmatrix[j][i]) * scale) );
543 
544  SCIP_CALL( SCIPaddOrigObjoffset(scip, (probdata->cmatrix[i][j] + probdata->cmatrix[j][i]) * scale) );
545 
546  /* create constraints that determine when the edge-variables have to be non-zero*/
547  for( c1 = 0; c1 < ncluster; ++c1 )
548  {
549  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bins_%d_%d_incluster_%d", i, j, c1);
550  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
552  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][1], -1.0) );
553 
554  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][i][1], 1.0) );
555  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][c1], 1.0) );
556  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][phi(c1, ncluster)], 1.0) );
557  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][phiinv(c1, ncluster)], -1.0) );
558  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][phiinv(c1, ncluster)], -1.0) );
559 
560  SCIP_CALL( SCIPaddCons(scip, temp) );
561  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
562 
563  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bins_%d_%d_incluster_%d", j, i, c1);
564  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
566 
567  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][i][1], -1.0) );
568  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][1], 1.0) );
569  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][c1], 1.0) );
570  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][phi(c1, ncluster)], 1.0) );
571  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[i][phiinv(c1, ncluster)], -1.0) );
572  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->binvars[j][phiinv(c1, ncluster)], -1.0) );
573 
574  SCIP_CALL( SCIPaddCons(scip, temp) );
575  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
576  }
577  }
578  }
579  }
580 
581  /* only one cluster-pair at the same time can be active for an edge*/
582  for( i = 0; i < nbins; ++i )
583  {
584  for( j = 0; j < i; ++j )
585  {
586  if( NULL == probdata->edgevars[i][j] || (SCIPisZero(scip, (probdata->cmatrix[i][j] - probdata->cmatrix[j][i]))
587  && SCIPisZero(scip, (probdata->cmatrix[i][j] + probdata->cmatrix[j][i]) * scale) ) )
588  continue;
589 
590  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "sumedge_%d_%d", i, j);
591  SCIP_CALL( SCIPcreateConsBasicLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0) );
592 
593  for( c1 = 0; c1 < 2; ++c1 )
594  {
595  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][c1], 1.0) );
596  }
597 
598  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][i][1], 1.0) );
599 
600  SCIP_CALL( SCIPaddCons(scip, temp) );
601  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
602  }
603  }
604 
605  /* add constraint that ensures that each cluster is used */
606  for( c1 = 0; c1 < ncluster; ++c1 )
607  {
608  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "cluster_%d_used", c1);
609  SCIP_CALL( SCIPcreateConsBasicLogicor(scip, &temp, consname, 0, NULL) );
610 
611  for( i = 0; i < nbins; ++i )
612  {
613  SCIP_CALL( SCIPaddCoefLogicor(scip, temp, probdata->binvars[i][c1]) );
614  }
615 
616  SCIP_CALL( SCIPaddCons(scip, temp) );
617  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
618  }
619 
620  return SCIP_OKAY;
621 }
622 
623 /** create the problem without variable amount of clusters, using quadratic formulations. This is inferior to the
624  * simplified variant. Only useful for comparing relaxations.
625  */
626 static
628  SCIP* scip, /**< SCIP Data Structure */
629  SCIP_PROBDATA* probdata /**< The problem data */
630  )
631 {
632  SCIP_VAR** edgevars;
633  SCIP_CONS* temp;
634  SCIP_Real scale;
635  char varname[SCIP_MAXSTRLEN];
636  char consname[SCIP_MAXSTRLEN];
637  int i;
638  int j;
639  int c1;
640  int c;
641  int nbins = probdata->nbins;
642  int ncluster = probdata->ncluster;
643 
644  SCIP_CALL( SCIPgetRealParam(scip, "cycleclustering/scale_coherence", &scale) );
645  probdata->scale = scale;
646  /* create variables for bins */
648 
649  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(probdata->binvars), nbins) );
650 
651  for( i = 0; i < nbins; ++i )
652  {
653  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(probdata->binvars[i]), ncluster) ); /*lint !e866*/
654  }
655 
656  for( i = 0; i < nbins; ++i )
657  {
658  for( c = 0; c < ncluster; ++c )
659  {
660  (void)SCIPsnprintf(varname, SCIP_MAXSTRLEN, "x_%d_%d", i, c);
661  SCIP_CALL( SCIPcreateVarBasic(scip, &probdata->binvars[i][c], varname, 0.0, 1.0, 0.0, SCIP_VARTYPE_BINARY) );
662  SCIP_CALL( SCIPaddVar(scip, probdata->binvars[i][c]) );
663  }
664  }
665 
666  /* create variables for the edges in each cluster combination. Index 0 are edges within cluster, 1 edges between
667  * consequtive clusters and 2 edges between non-consequtive clusters
668  */
669  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &(edgevars), (SCIP_Longint) ncluster * 2) );
670 
671  for( i = 0; i < 2 * ncluster; ++i )
672  {
673  (void) SCIPsnprintf(varname, SCIP_MAXSTRLEN, "f_%d", i);
674 
675  SCIP_CALL( SCIPcreateVarBasic(scip, &edgevars[i], varname, -SCIPinfinity(scip), SCIPinfinity(scip),
676  1.0, SCIP_VARTYPE_CONTINUOUS) );
677  SCIP_CALL( SCIPaddVar(scip, edgevars[i]) );
678  }
679 
680  /* create variables for the edges in each cluster combination. Index 0 are edges within cluster, 1 edges between
681  * consequtive clusters and 2 edges between non-consequtive clusters
682  */
683  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &(probdata->edgevars), nbins) );
684 
685  for( i = 0; i < nbins; ++i )
686  {
687  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &(probdata->edgevars[i]), nbins) ); /*lint !e866*/
688 
689  for( j = 0; j < nbins; ++j )
690  probdata->edgevars[i][j] = NULL;
691  }
692 
693  /*
694  * create constraints
695  */
696 
697  /* create the set-partitioning constraints of the bins */
698  for( i = 0; i < nbins; ++i )
699  {
700  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "setpart_%d", i+1);
701  SCIP_CALL( SCIPcreateConsSetpart(scip, &temp, consname, 0, NULL, TRUE, TRUE, TRUE, TRUE, TRUE,
702  FALSE, FALSE, FALSE, FALSE, FALSE) );
703 
704  for ( c1 = 0; c1 < ncluster; ++c1 )
705  {
706  SCIP_CALL( SCIPaddCoefSetppc(scip, temp, probdata->binvars[i][c1]) );
707  }
708 
709  SCIP_CALL( SCIPaddCons(scip, temp) );
710  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
711  }
712 
713  /* add constraint that ensures that each cluster is used */
714  for( c1 = 0; c1 < ncluster; ++c1 )
715  {
716  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "cluster_%d_used", c1);
717  SCIP_CALL( SCIPcreateConsBasicLogicor(scip, &temp, consname, 0, NULL) );
718 
719  for( i = 0; i < nbins; ++i )
720  {
721  SCIP_CALL( SCIPaddCoefLogicor(scip, temp, probdata->binvars[i][c1]) );
722  }
723 
724  SCIP_CALL( SCIPaddCons(scip, temp) );
725  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
726  }
727 
728  for( c = 0; c < ncluster; ++c)
729  {
730  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "irrev_%d", c);
731  SCIP_CALL( SCIPcreateConsBasicQuadratic(scip, &temp, consname, 0, NULL, NULL, 0, NULL, NULL, NULL,
732  -SCIPinfinity(scip), 0.0) );
733 
734  SCIP_CALL( SCIPaddLinearVarQuadratic(scip, temp, edgevars[c], 1.0) );
735 
736  for( i = 0; i < nbins; ++i )
737  {
738  for( j = 0; j < nbins; ++j )
739  {
740  if( i != j )
741  {
742  SCIP_CALL( SCIPaddBilinTermQuadratic(scip, temp, probdata->binvars[i][c],
743  probdata->binvars[j][phi(c,ncluster)], probdata->cmatrix[j][i] - probdata->cmatrix[i][j]) );
744  }
745  }
746  }
747 
748  SCIP_CALL( SCIPaddCons(scip, temp) );
749  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
750  }
751 
752  for( c = 0; c < ncluster; ++c )
753  {
754  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "coh_%d", c);
755  SCIP_CALL( SCIPcreateConsBasicQuadratic(scip, &temp, consname, 0, NULL, NULL, 0, NULL, NULL, NULL,
756  -SCIPinfinity(scip), 0.0) );
757 
758  SCIP_CALL( SCIPaddLinearVarQuadratic(scip, temp, edgevars[c+ncluster], 1.0) );
759 
760  for( i = 0; i < nbins; ++i)
761  {
762  for( j = 0; j < nbins; ++j )
763  {
764  if( i > j )
765  {
766  SCIP_CALL( SCIPaddBilinTermQuadratic(scip, temp, probdata->binvars[i][c], probdata->binvars[j][c],
767  -scale * (probdata->cmatrix[i][j] + probdata->cmatrix[j][i])) );
768  }
769  }
770  }
771 
772  SCIP_CALL( SCIPaddCons(scip, temp) );
773  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
774  }
775 
776  for (i = 0; i < 2*ncluster; i++ )
777  {
778  SCIP_CALL( SCIPreleaseVar(scip, &edgevars[i]) );
779  }
780 
781  SCIPfreeBlockMemoryArray(scip, &edgevars, (SCIP_Longint) 2 * ncluster);
782 
783  return SCIP_OKAY;
784 }
785 
786 
787 /** create the problem with variable amount of clusters. Very large number of constraints not viable for large scale
788  * problems.
789  */
790 static
792  SCIP* scip, /**< SCIP Data Structure */
793  SCIP_PROBDATA* probdata /**< The problem data */
794  )
795 {
796  SCIP_CONS* temp;
797  SCIP_Real scale;
798  SCIP_Real sign[3][3];
799  char consname[SCIP_MAXSTRLEN];
800  int nbins = probdata->nbins;
801  int i;
802  int j;
803  int k;
804  int l;
805 
806  SCIP_CALL( SCIPgetRealParam(scip, "cycleclustering/scale_coherence", &scale) );
807  probdata->scale = scale;
808 
809  for( i = 0; i < 3; ++i )
810  {
811  for( j = 0; j < 3; ++j )
812  {
813  sign[i][j] = 1.0;
814  }
815  sign[i][i] = -1.0;
816  }
817  /*
818  * create constraints
819  */
820 
821  /* create constraints for the edge-cut variables */
823  "Using edge-representation with simplified structure. No variable amount of cluster. \n");
824 
825  for( i = 0; i < nbins; ++i )
826  {
827  for( j = 0; j < nbins; ++j )
828  {
829  /* set the objective weight for the edge-variables */
830 
831  /* these edges are not within a cluster */
832  if( j < i )
833  {
834  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[i][j][0], (probdata->cmatrix[i][j] + probdata->cmatrix[j][i]) * scale) );
835  /* these are the edges that are between consequtive clusters */
836  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[i][j][1], (probdata->cmatrix[i][j] - probdata->cmatrix[j][i])) );
837  SCIP_CALL( SCIPchgVarObj(scip, probdata->edgevars[j][i][1], (probdata->cmatrix[j][i] - probdata->cmatrix[i][j])) );
838 
839  (void)SCIPsnprintf(consname, SCIP_MAXSTRLEN, "bins_%d_%d", i+1, j+1);
840  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
842 
843  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][0], 1.0) );
844  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][1], 1.0) );
845  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][i][1], 1.0) );
846 
847  SCIP_CALL( SCIPaddCons(scip, temp) );
848  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
849  }
850 
851  for( k = 0; k < nbins; ++k )
852  {
853  if( i == k || i == j || k == j )
854  continue;
855 
856  if( k < j && j < i )
857  {
858  for( l = 0; l < 3; l++ )
859  {
860  (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "tri_%d_%d_%d", i, j, k);
861  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
863 
864  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][j][0], sign[l][0]) );
865  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][k][0], sign[l][1]) );
866  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][k][0], sign[l][2]) );
867 
868  SCIP_CALL( SCIPaddCons(scip, temp) );
869  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
870  }
871  }
872 
873  (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "tri_%d_%d_%d", i, j, k);
874  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
876 
877  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[MAX(i,j)][MIN(i,j)][0], 1.0) );
878  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][k][1], 1.0) );
879  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][k][1], -1.0) );
880 
881  SCIP_CALL( SCIPaddCons(scip, temp) );
882  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
883 
884  (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "tri_%d_%d_%d", i, j, k);
885  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
887 
888  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[MAX(i,j)][MIN(i,j)][0], 1.0) );
889  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[k][j][1], 1.0) );
890  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[k][i][1], -1.0) );
891 
892  SCIP_CALL( SCIPaddCons(scip, temp) );
893  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
894 
895  (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "tri_%d_%d_%d", i, j, k);
896  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
898 
899  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[MAX(i,j)][MIN(i,j)][0], -1.0) );
900  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[j][k][1], 1.0) );
901  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[i][k][1], 1.0) );
902 
903  SCIP_CALL( SCIPaddCons(scip, temp) );
904  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
905 
906  (void) SCIPsnprintf(consname, SCIP_MAXSTRLEN, "tri_%d_%d_%d", i, j, k);
907  SCIP_CALL( SCIPcreateConsLinear(scip, &temp, consname, 0, NULL, NULL, -SCIPinfinity(scip), 1.0,
909 
910  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[MAX(i,j)][MIN(i,j)][0], -1.0) );
911  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[k][j][1], 1.0) );
912  SCIP_CALL( SCIPaddCoefLinear(scip, temp, probdata->edgevars[k][i][1], 1.0) );
913 
914  SCIP_CALL( SCIPaddCons(scip, temp) );
915  SCIP_CALL( SCIPreleaseCons(scip, &temp) );
916  }
917  }
918  }
919 
920  return SCIP_OKAY;
921 }
922 
923 /** Scip callback to transform the problem */
924 static
925 SCIP_DECL_PROBTRANS(probtransCyc)
926 {
927  int i;
928  int j;
929  int c;
930  int edgetype;
931  int nbins = sourcedata->nbins;
932  int ncluster = sourcedata->ncluster;
933 
934  assert(scip != NULL);
935  assert(sourcedata != NULL);
936  assert(targetdata != NULL);
937 
938  /* allocate memory */
939  SCIP_CALL( SCIPallocBlockMemory(scip, targetdata) );
940 
941  (*targetdata)->nbins = sourcedata->nbins;
942  (*targetdata)->ncluster = sourcedata->ncluster;
943  (*targetdata)->model_variant = sourcedata->model_variant;
944  (*targetdata)->scale = sourcedata->scale;
945 
946  /* allocate memory */
947  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &((*targetdata)->cmatrix), nbins) );
948  for( i = 0; i < nbins; ++i )
949  {
950  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &((*targetdata)->cmatrix[i]), nbins) ); /*lint !e866*/
951  }
952  /* copy the matrizes */
953  for ( i = 0; i < nbins; ++i )
954  {
955  for ( j = 0; j < nbins; ++j )
956  {
957  (*targetdata)->cmatrix[i][j] = sourcedata->cmatrix[i][j];
958  }
959  }
960 
961  /* copy the variables */
962  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &((*targetdata)->binvars), nbins) );
963  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &((*targetdata)->edgevars), nbins) );
964 
965  for( i = 0; i < nbins; ++i )
966  {
967  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &((*targetdata)->edgevars[i]), nbins) ); /*lint !e866*/
968  for( j = 0; j < nbins; ++j )
969  {
970  if( sourcedata->edgevars[i][j] == NULL || i == j)
971  continue;
972  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &((*targetdata)->edgevars[i][j]), 3) ); /*lint !e866*/
973  for( edgetype = 0; edgetype < 3; ++edgetype )
974  {
975  if( edgetype == 0 && i < j )
976  continue;
977  if( sourcedata->edgevars[i][j][edgetype] != NULL )
978  {
979  SCIP_CALL( SCIPtransformVar(scip, sourcedata->edgevars[i][j][edgetype],
980  &((*targetdata)->edgevars[i][j][edgetype])) );
981  }
982  else
983  ((*targetdata)->edgevars[i][j][edgetype]) = NULL;
984  }
985  }
986  }
987 
988  for( i = 0; i < nbins; ++i )
989  {
990  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &((*targetdata)->binvars[i]), ncluster) ); /*lint !e866*/
991 
992  for( c = 0; c < ncluster; ++c )
993  {
994  if( sourcedata->binvars[i][c] != NULL )
995  {
996  SCIP_CALL( SCIPtransformVar(scip, sourcedata->binvars[i][c], &((*targetdata)->binvars[i][c])) );
997  }
998  else
999  (*targetdata)->binvars[i][c] = NULL;
1000  }
1001  }
1002 
1003  SCIP_CALL( SCIPcopyDigraph(scip, &((*targetdata)->edgegraph), sourcedata->edgegraph) );
1004 
1005  return SCIP_OKAY;
1006 }
1007 
1008 /** delete-callback method of scip */
1009 static
1010 SCIP_DECL_PROBDELORIG(probdelorigCyc)
1011 {
1012  int c;
1013  int edgetype;
1014  int i;
1015  int j;
1016 
1017  assert(probdata != NULL);
1018  assert(*probdata != NULL);
1019 
1020  /* release all the variables */
1021 
1022  /* binvars */
1023  for ( c = 0; c < (*probdata)->nbins; ++c )
1024  {
1025  for ( i = 0; i < (*probdata)->ncluster; ++i )
1026  {
1027  if( (*probdata)->binvars[c][i] != NULL )
1028  {
1029  SCIP_CALL( SCIPreleaseVar(scip, &((*probdata)->binvars[c][i])) );
1030  }
1031  }
1032  SCIPfreeBlockMemoryArray( scip, &((*probdata)->binvars[c]), (*probdata)->ncluster);
1033  }
1034 
1035  /* cut-edge vars */
1036  for ( i = 0; i < (*probdata)->nbins; ++i )
1037  {
1038  for( j = 0; j < (*probdata)->nbins; ++j )
1039  {
1040  if( (*probdata)->edgevars[i][j] != NULL && j != i )
1041  {
1042  for ( edgetype = 0; edgetype < 3; ++edgetype )
1043  {
1044  if( edgetype == 0 && i < j )
1045  continue;
1046 
1047  if( (*probdata)->edgevars[i][j][edgetype] != NULL )
1048  {
1049  SCIP_CALL( SCIPreleaseVar( scip, &((*probdata)->edgevars[i][j][edgetype])) );
1050  }
1051  }
1052 
1053  SCIPfreeBlockMemoryArray(scip, &((*probdata)->edgevars[i][j]), 3);
1054  }
1055  }
1056 
1057  SCIPfreeBlockMemoryArray(scip, &((*probdata)->edgevars[i]), (*probdata)->nbins);
1058  }
1059 
1060  SCIPfreeBlockMemoryArray(scip, &((*probdata)->edgevars), (*probdata)->nbins);
1061  SCIPfreeBlockMemoryArray(scip, &((*probdata)->binvars), (*probdata)->nbins);
1062 
1063  SCIPdigraphFreeComponents((*probdata)->edgegraph);
1064  SCIPdigraphFree(&((*probdata)->edgegraph));
1065 
1066  for ( i = 0; i < (*probdata)->nbins; ++i )
1067  {
1068  SCIPfreeBlockMemoryArray(scip, &((*probdata)->cmatrix[i]), (*probdata)->nbins);
1069  }
1070  SCIPfreeBlockMemoryArray(scip, &(*probdata)->cmatrix, (*probdata)->nbins);
1071 
1072  SCIPfreeBlockMemory(scip, probdata);
1073 
1074  return SCIP_OKAY;
1075 }
1076 
1077 /** scip-callback to delete the transformed problem */
1078 static
1079 SCIP_DECL_PROBDELTRANS(probdeltransCyc)
1080 {
1081  int c;
1082  int edgetype;
1083  int i;
1084  int j;
1085 
1086  assert(probdata != NULL);
1087  assert(*probdata != NULL);
1088 
1089  /* release all the variables */
1090 
1091  /* binvars */
1092  for ( i = 0; i < (*probdata)->nbins; ++i )
1093  {
1094  for ( c = 0; c < (*probdata)->ncluster; ++c )
1095  {
1096  if( (*probdata)->binvars[i][c] != NULL )
1097  {
1098  SCIP_CALL( SCIPreleaseVar(scip, &((*probdata)->binvars[i][c])) );
1099  }
1100  }
1101  SCIPfreeBlockMemoryArray(scip, &((*probdata)->binvars[i]), (*probdata)->ncluster);
1102  }
1103 
1104  /* cut-edge vars */
1105  for ( i = 0; i < (*probdata)->nbins; ++i )
1106  {
1107  for( j = 0; j < (*probdata)->nbins; ++j )
1108  {
1109  if( (*probdata)->edgevars[i][j] != NULL && j != i )
1110  {
1111  for ( edgetype = 0; edgetype < 3; ++edgetype )
1112  {
1113  if( 0 == edgetype && j > i )
1114  continue;
1115 
1116  if( (*probdata)->edgevars[i][j][edgetype] != NULL )
1117  {
1118  SCIP_CALL( SCIPreleaseVar( scip, &((*probdata)->edgevars[i][j][edgetype])) );
1119  }
1120  }
1121 
1122  SCIPfreeBlockMemoryArray(scip, &((*probdata)->edgevars[i][j]), 3);
1123  }
1124  }
1125 
1126  SCIPfreeBlockMemoryArray(scip, &((*probdata)->edgevars[i]), (*probdata)->nbins);
1127  }
1128 
1129  SCIPfreeBlockMemoryArray(scip, &((*probdata)->edgevars), (*probdata)->nbins);
1130  SCIPfreeBlockMemoryArray(scip, &((*probdata)->binvars), (*probdata)->nbins);
1131 
1132  SCIPdigraphFreeComponents((*probdata)->edgegraph);
1133  SCIPdigraphFree(&((*probdata)->edgegraph));
1134 
1135  for ( i = 0; i < (*probdata)->nbins; ++i )
1136  {
1137  SCIPfreeBlockMemoryArray(scip, &((*probdata)->cmatrix[i]), (*probdata)->nbins);
1138  }
1139  SCIPfreeBlockMemoryArray(scip, &(*probdata)->cmatrix, (*probdata)->nbins);
1140 
1141  SCIPfreeBlockMemory(scip, probdata);
1142 
1143  return SCIP_OKAY;
1144 }
1145 
1146 /** callback method of scip for copying the probdata */
1147 static
1149 {
1150  SCIP_Bool success;
1151  SCIP_VAR* var;
1152  int nbins;
1153  int ncluster;
1154  int edgetype;
1155  int i;
1156  int j;
1157  int c;
1158 
1159  assert(scip != NULL);
1160  assert(sourcescip != NULL);
1161  assert(sourcedata != NULL);
1162  assert(targetdata != NULL);
1163 
1164  /* set up data */
1165  SCIP_CALL( SCIPallocBlockMemory(scip, targetdata) );
1166 
1167  nbins = sourcedata->nbins;
1168  ncluster = sourcedata->ncluster;
1169  success = FALSE;
1170 
1171  (*targetdata)->nbins = nbins;
1172  (*targetdata)->ncluster = ncluster;
1173  (*targetdata)->model_variant = sourcedata->model_variant;
1174  (*targetdata)->scale = sourcedata->scale;
1175 
1176  /* allocate memory */
1177  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &((*targetdata)->cmatrix), nbins) );
1178  for( i = 0; i < nbins; ++i )
1179  {
1180  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &((*targetdata)->cmatrix[i]), nbins) ); /*lint !e866*/
1181  }
1182  /* copy the matrizes */
1183  for ( i = 0; i < nbins; ++i )
1184  {
1185  for ( j = 0; j < nbins; ++j )
1186  {
1187  (*targetdata)->cmatrix[i][j] = sourcedata->cmatrix[i][j];
1188  }
1189  }
1190 
1191  /* copy the variables */
1192  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &((*targetdata)->binvars), nbins) );
1193  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &((*targetdata)->edgevars), nbins) );
1194 
1195  for( i = 0; i < nbins; ++i )
1196  {
1197  SCIP_CALL( SCIPallocClearBlockMemoryArray(scip, &((*targetdata)->edgevars[i]), nbins) ); /*lint !e866*/
1198 
1199  for( j = 0; j < nbins; ++j )
1200  {
1201  if( (sourcedata)->edgevars[i][j] == NULL || j == i )
1202  continue;
1203 
1204  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &((*targetdata)->edgevars[i][j]), 3) ); /*lint !e866*/
1205 
1206  for( edgetype = 0; edgetype < 3; ++edgetype )
1207  {
1208  if( edgetype == 0 && j > i )
1209  continue;
1210 
1211  if( sourcedata->edgevars[i][j][edgetype] != NULL )
1212  {
1213  SCIP_CALL( SCIPgetTransformedVar(sourcescip, sourcedata->edgevars[i][j][edgetype], &var) );
1214 
1215  if( SCIPvarIsActive(var) )
1216  {
1217  SCIP_CALL( SCIPgetVarCopy(sourcescip, scip, var, &((*targetdata)->edgevars[i][j][edgetype]),
1218  varmap, consmap, global, &success) );
1219 
1220  assert(success);
1221  assert((*targetdata)->edgevars[i][j][edgetype] != NULL);
1222 
1223  SCIP_CALL( SCIPcaptureVar(scip, (*targetdata)->edgevars[i][j][edgetype]) );
1224  }
1225  else
1226  (*targetdata)->edgevars[i][j][edgetype] = NULL;
1227  }
1228  else
1229  (*targetdata)->edgevars[i][j][edgetype] = NULL;
1230  }
1231  }
1232  }
1233 
1234  for( i = 0; i < nbins; ++i )
1235  {
1236  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &((*targetdata)->binvars[i]), ncluster) ); /*lint !e866*/
1237 
1238  for( c = 0; c < ncluster; ++c )
1239  {
1240  if( sourcedata->binvars[i][c] != NULL )
1241  {
1242  SCIP_CALL( SCIPgetTransformedVar(sourcescip, sourcedata->binvars[i][c], &var) );
1243 
1244  if( SCIPvarIsActive(var) )
1245  {
1246  SCIP_CALL( SCIPgetVarCopy(sourcescip, scip, var, &((*targetdata)->binvars[i][c]),
1247  varmap, consmap, global, &success) );
1248 
1249  assert(success);
1250  assert((*targetdata)->binvars[i][c] != NULL);
1251 
1252  SCIP_CALL( SCIPcaptureVar(scip, (*targetdata)->binvars[i][c]) );
1253  }
1254  else
1255  (*targetdata)->binvars[i][c] = NULL;
1256  }
1257  else
1258  (*targetdata)->binvars[i][c] = NULL;
1259  }
1260  }
1261 
1262  SCIP_CALL( SCIPcopyDigraph(scip, &((*targetdata)->edgegraph), sourcedata->edgegraph) );
1263 
1264  assert(success);
1265 
1266  *result = SCIP_SUCCESS;
1267 
1268  return SCIP_OKAY;
1269 }
1270 
1271 /**
1272  * Create the probdata for an cyc-clustering problem
1273  */
1275  SCIP* scip, /**< SCIP data structure */
1276  const char* name, /**< problem name */
1277  int nbins, /**< number of bins */
1278  int ncluster, /**< number of cluster */
1279  SCIP_Real** cmatrix /**< The transition matrix */
1280  )
1281 {
1282  SCIP_PROBDATA* probdata = NULL;
1283  int i;
1284  int j;
1285  char model;
1286 
1287  assert(nbins > 0); /* at least one node */
1288  assert(ncluster <= nbins);
1289 
1290  SCIP_CALL( SCIPcreateProbBasic(scip, name) );
1291 
1292  /* Set up the problem */
1293  SCIP_CALL( SCIPallocBlockMemory(scip, &probdata) );
1294 
1295  SCIP_CALL( SCIPcreateDigraph(scip, &(probdata->edgegraph), nbins) );
1296 
1297  /* allocate memory for the matrizes and create them from the edge-arrays */
1298  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(probdata->cmatrix), nbins) );
1299  for ( i = 0; i < nbins; ++i )
1300  {
1301  SCIP_CALL( SCIPallocBlockMemoryArray(scip, &(probdata->cmatrix[i]), nbins) ); /*lint !e866*/
1302  for( j = 0; j < nbins; ++j )
1303  {
1304  probdata->cmatrix[i][j] = cmatrix[i][j];
1305  }
1306  }
1307 
1308  SCIPverbMessage(scip, SCIP_VERBLEVEL_NORMAL, NULL, "Creating problem: %s \n", name);
1309 
1310  /* create variables */
1311  probdata->nbins=nbins;
1312  probdata->ncluster=ncluster;
1313 
1314  SCIP_CALL( SCIPgetCharParam(scip, "cycleclustering/model", &model) );
1315 
1316  /* create constraints depending on model selection */
1317  switch( model )
1318  {
1319  case 's':
1320  SCIP_CALL( createVariables(scip, probdata) );
1321  SCIP_CALL( createProbSimplified(scip, probdata) );
1322  break;
1323  case 't':
1324  SCIP_CALL( createVariables(scip, probdata) );
1325  SCIP_CALL( createProbSimplifiedTest(scip, probdata) );
1326  break;
1327  case 'e':
1328  SCIP_CALL( createVariables(scip, probdata) );
1329  SCIP_CALL( createProbOnlyEdge(scip, probdata) );
1330  break;
1331  case 'q':
1332  SCIP_CALL( createProbQP(scip, probdata) );
1333  break;
1334  default:
1335  SCIPABORT();
1336  break;
1337  }
1338 
1339  /** add callback methods to scip */
1340  SCIP_CALL( SCIPsetProbDelorig(scip, probdelorigCyc) );
1341  SCIP_CALL( SCIPsetProbCopy(scip, probcopyCyc) );
1342  SCIP_CALL( SCIPsetProbData(scip, probdata) );
1343  SCIP_CALL( SCIPsetProbTrans(scip, probtransCyc) );
1344  SCIP_CALL( SCIPsetProbDeltrans(scip, probdeltransCyc) );
1345 
1346  return SCIP_OKAY;
1347 }
1348 
1349 /** Getter methods for the various parts of the probdata */
1350 
1351 
1352 /** Returns the transition matrix*/
1354  SCIP* scip /**< SCIP data structure */
1355  )
1356 {
1357  SCIP_PROBDATA* probdata;
1358 
1359  assert(scip != NULL);
1360 
1361  probdata = SCIPgetProbData(scip);
1362 
1363  assert(probdata != NULL);
1364 
1365  return probdata->cmatrix;
1366 }
1367 
1368 /** Returns the number of states */
1370  SCIP* scip /**< SCIP data structure */
1371  )
1372 {
1373  SCIP_PROBDATA* probdata;
1374 
1375  assert(scip != NULL);
1376 
1377  probdata = SCIPgetProbData(scip);
1378 
1379  assert(probdata != NULL);
1380 
1381  return probdata->nbins;
1382 }
1383 
1384 /** Returns the number of clusters*/
1386  SCIP* scip /**< SCIP data structure */
1387  )
1388 {
1389  SCIP_PROBDATA* probdata;
1390 
1391  assert(scip!= NULL);
1392 
1393  probdata = SCIPgetProbData(scip);
1394 
1395  assert(probdata != NULL);
1396 
1397  return probdata->ncluster;
1398 }
1399 
1400 /** Returns the state-variable-matrix*/
1402  SCIP* scip /**< SCIP data structure */
1403  )
1404 {
1405  SCIP_PROBDATA* probdata;
1406 
1407  assert(scip!= NULL);
1408 
1409  probdata = SCIPgetProbData(scip);
1410 
1411  assert(probdata != NULL);
1412  assert(probdata->binvars != NULL);
1413 
1414  return probdata->binvars;
1415 }
1416 
1417 /** Returns the scaling parameter*/
1419  SCIP* scip /**< SCIP data structure */
1420  )
1421 {
1422  SCIP_PROBDATA* probdata;
1423 
1424  assert(scip!= NULL);
1425 
1426  probdata = SCIPgetProbData(scip);
1427 
1428  assert(probdata != NULL);
1429 
1430  return probdata->scale;
1431 }
1432 
1433 /** Returns the edge variables */
1435  SCIP* scip /**< SCIP data structure */
1436  )
1437 {
1438  SCIP_PROBDATA* probdata;
1439 
1440  assert(scip!= NULL);
1441 
1442  probdata = SCIPgetProbData(scip);
1443 
1444  assert(probdata != NULL);
1445  assert(probdata->edgevars != NULL);
1446 
1447  return probdata->edgevars;
1448 }
1449 
1450 /** return one specific edge variable */
1452  SCIP_VAR**** edgevars, /**< edgevar data structure*/
1453  int state1, /**< first state */
1454  int state2, /**< second state */
1455  int direction /**< direction, 0 = incluster, 1 = forward */
1456  )
1457 {
1458  assert(edgevars != NULL);
1459  assert(edgevars[state1] != NULL);
1460  assert(edgevars[state1][state2] != NULL);
1461  assert(edgevars[state1][state2][direction] != NULL);
1462 
1463  return edgevars[state1][state2][direction];
1464 }
1465 
1466 /** check for an array of states, if all possible edge-combinations exist */
1468  SCIP_VAR**** edgevars, /**< edgevar data structure */
1469  int* states, /**< state array */
1470  int nstates /**< size of state array */
1471  )
1472 {
1473  int i;
1474  int j;
1475 
1476  assert(edgevars != NULL);
1477  assert(states != NULL);
1478 
1479  for( i = 0; i < nstates; ++i )
1480  {
1481  assert(edgevars[states[i]] != NULL);
1482 
1483  for( j = 0; j < nstates; ++j )
1484  {
1485  if( j != i )
1486  {
1487  assert(edgevars[states[j]] != NULL);
1488 
1489  if( edgevars[states[i]][states[j]] == NULL )
1490  return FALSE;
1491  }
1492  }
1493  }
1494 
1495  return TRUE;
1496 }
1497 
1498 /** Returns the edge-graph */
1500  SCIP* scip /**< SCIP data structure */
1501  )
1502 {
1503  SCIP_PROBDATA* probdata;
1504 
1505  assert(scip!= NULL);
1506 
1507  probdata = SCIPgetProbData(scip);
1508 
1509  assert(probdata != NULL);
1510  assert(probdata->edgegraph != NULL);
1511 
1512  return probdata->edgegraph;
1513 }
1514 
1515 
1516 /** print the model-values like coherence in the clusters and transition-probabilities between clusters that are not
1517  * evident from the scip-solution
1518  */
1520  SCIP* scip, /**< SCIP data structure */
1521  SCIP_SOL* sol /**< The solution containg the values */
1522  )
1523 {
1524  SCIP_PROBDATA* probdata;
1525  SCIP_Real value;
1526  SCIP_Real objvalue = 0;
1527  SCIP_Real flow = 0;
1528  SCIP_Real coherence = 0;
1529  int c1;
1530  int c2;
1531  int c3;
1532  int i;
1533  int j;
1534 
1535  assert(scip!= NULL);
1536 
1537  probdata = SCIPgetProbData(scip);
1538 
1539  assert(probdata != NULL);
1540 
1541  SCIPverbMessage(scip, SCIP_VERBLEVEL_NORMAL, NULL, "\nDisplaying aggregated solution data: \n");
1542 
1543  for( c1 = 0; c1 < probdata->ncluster; ++c1 )
1544  {
1545  value = 0;
1546 
1547  for( i = 0; i < probdata->nbins; ++i )
1548  {
1549  for( j = 0; j < probdata->nbins; ++j )
1550  {
1551  if( j == i )
1552  continue;
1553 
1554  value+= probdata->cmatrix[i][j] * SCIPgetSolVal(scip, sol, probdata->binvars[i][c1])
1555  * SCIPgetSolVal(scip, sol, probdata->binvars[j][c1]);
1556  }
1557  }
1558 
1559  SCIPverbMessage(scip, SCIP_VERBLEVEL_NORMAL, NULL, " Coherence in cluster %d : %f \n", c1 + 1, value);
1560  coherence += value;
1561  objvalue += probdata->scale * value;
1562  }
1563 
1564  for( c1 = 0; c1 < probdata->ncluster; ++c1 )
1565  {
1566  for( c2 = 0; c2 < probdata->ncluster; ++c2 )
1567  {
1568  value = 0;
1569 
1570  for( i = 0; i < probdata->nbins; ++i )
1571  {
1572  for( j = 0; j < probdata->nbins; ++j )
1573  {
1574  value+= probdata->cmatrix[i][j] * SCIPgetSolVal(scip, sol, probdata->binvars[i][c1]) *
1575  SCIPgetSolVal(scip, sol, probdata->binvars[j][c2]);
1576  }
1577  }
1578  }
1579 
1580  c3 = (c1 + 1) % probdata->ncluster;
1581  value = 0;
1582 
1583  for( i = 0; i < probdata->nbins; ++i )
1584  {
1585  for( j = 0; j < probdata->nbins; ++j )
1586  {
1587  value+= (probdata->cmatrix[i][j] - probdata->cmatrix[j][i])
1588  * SCIPgetSolVal(scip, sol, probdata->binvars[i][c1]) * SCIPgetSolVal(scip, sol, probdata->binvars[j][c3]);
1589  }
1590  }
1591 
1593  " Net flow from %d to %d : %f \n", c1, phi(c1, probdata->ncluster), value);
1594 
1595  flow += value;
1596  objvalue += value;
1597  }
1598 
1599  SCIPverbMessage(scip, SCIP_VERBLEVEL_NORMAL, NULL, " Total coherence : %f \n", coherence);
1600  SCIPverbMessage(scip, SCIP_VERBLEVEL_NORMAL, NULL, " Total net flow : %f \n", flow);
1601  SCIPverbMessage(scip, SCIP_VERBLEVEL_NORMAL, NULL, " Total objective value : %f \n", objvalue);
1602 
1603  return SCIP_OKAY;
1604 }
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
#define SCIPfreeBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:97
SCIP_RETCODE SCIPcreateConsLinear(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_Real *vals, SCIP_Real lhs, SCIP_Real rhs, SCIP_Bool initial, SCIP_Bool separate, SCIP_Bool enforce, SCIP_Bool check, SCIP_Bool propagate, SCIP_Bool local, SCIP_Bool modifiable, SCIP_Bool dynamic, SCIP_Bool removable, SCIP_Bool stickingatnode)
SCIP_RETCODE SCIPcreateProbCyc(SCIP *scip, const char *name, int nbins, int ncluster, SCIP_Real **cmatrix)
SCIP_EXPORT SCIP_Bool SCIPvarIsTransformed(SCIP_VAR *var)
Definition: var.c:17154
SCIP_Bool SCIPisIntegral(SCIP *scip, SCIP_Real val)
#define SCIPallocBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:80
SCIP_RETCODE SCIPcreateDigraph(SCIP *scip, SCIP_DIGRAPH **digraph, int nnodes)
SCIP_RETCODE SCIPsetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var, SCIP_Real val)
Definition: scip_sol.c:1213
SCIP_Bool SCIPisGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_RETCODE SCIPcreateConsBasicQuadratic(SCIP *scip, SCIP_CONS **cons, const char *name, int nlinvars, SCIP_VAR **linvars, SCIP_Real *lincoefs, int nquadterms, SCIP_VAR **quadvars1, SCIP_VAR **quadvars2, SCIP_Real *quadcoefs, SCIP_Real lhs, SCIP_Real rhs)
SCIP_EXPORT SCIP_VAR * SCIPvarGetTransVar(SCIP_VAR *var)
Definition: var.c:17365
SCIP_VAR * getEdgevar(SCIP_VAR ****edgevars, int state1, int state2, int direction)
#define SCIP_MAXSTRLEN
Definition: def.h:273
SCIP_RETCODE SCIPsetObjsense(SCIP *scip, SCIP_OBJSENSE objsense)
Definition: scip_prob.c:1240
#define SCIPallocClearBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:84
SCIP_RETCODE assignVars(SCIP *scip, SCIP_SOL *sol, SCIP_Real **clustering, int nbins, int ncluster)
Definition: probdata_cyc.c:79
SCIP_RETCODE SCIPcaptureVar(SCIP *scip, SCIP_VAR *var)
Definition: scip_var.c:1218
SCIP_Real SCIPgetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var)
Definition: scip_sol.c:1353
SCIP_RETCODE SCIPaddCoefSetppc(SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var)
Definition: cons_setppc.c:9226
SCIP_RETCODE SCIPaddCoefLinear(SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real val)
SCIP_VAR **** SCIPcycGetEdgevars(SCIP *scip)
SCIP_RETCODE SCIPcreateVarBasic(SCIP *scip, SCIP_VAR **var, const char *name, SCIP_Real lb, SCIP_Real ub, SCIP_Real obj, SCIP_VARTYPE vartype)
Definition: scip_var.c:185
void SCIPverbMessage(SCIP *scip, SCIP_VERBLEVEL msgverblevel, FILE *file, const char *formatstr,...)
Definition: scip_message.c:216
#define FALSE
Definition: def.h:73
SCIP_RETCODE SCIPcreateConsBasicLinear(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_Real *vals, SCIP_Real lhs, SCIP_Real rhs)
#define TRUE
Definition: def.h:72
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:54
SCIP_PROBDATA * SCIPgetProbData(SCIP *scip)
Definition: scip_prob.c:962
static SCIP_DECL_PROBTRANS(probtransCyc)
Definition: probdata_cyc.c:925
void SCIPdigraphFree(SCIP_DIGRAPH **digraph)
Definition: misc.c:7459
SCIP_EXPORT SCIP_Bool SCIPvarIsActive(SCIP_VAR *var)
Definition: var.c:17335
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip_mem.h:95
SCIP_Real SCIPcycGetScale(SCIP *scip)
SCIP_EXPORT SCIP_VARSTATUS SCIPvarGetStatus(SCIP_VAR *var)
Definition: var.c:17131
#define SCIPallocBlockMemory(scip, ptr)
Definition: scip_mem.h:78
SCIP_RETCODE SCIPgetCharParam(SCIP *scip, const char *name, char *value)
Definition: scip_param.c:317
static SCIP_DECL_PROBCOPY(probcopyCyc)
SCIP_Bool SCIPisLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_RETCODE SCIPgetVarCopy(SCIP *sourcescip, SCIP *targetscip, SCIP_VAR *sourcevar, SCIP_VAR **targetvar, SCIP_HASHMAP *varmap, SCIP_HASHMAP *consmap, SCIP_Bool global, SCIP_Bool *success)
Definition: scip_copy.c:697
SCIP_RETCODE SCIPaddOrigObjoffset(SCIP *scip, SCIP_Real addval)
Definition: scip_prob.c:1288
static SCIP_RETCODE createProbOnlyEdge(SCIP *scip, SCIP_PROBDATA *probdata)
Definition: probdata_cyc.c:791
int phi(int k, int ncluster)
Definition: probdata_cyc.c:172
SCIP_RETCODE SCIPsetProbData(SCIP *scip, SCIP_PROBDATA *probdata)
Definition: scip_prob.c:1012
SCIP_RETCODE SCIPaddCoefLogicor(SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var)
Constraint handler for logicor constraints (equivalent to set covering, but algorithms are suited fo...
SCIP_RETCODE SCIPgetTransformedVar(SCIP *scip, SCIP_VAR *var, SCIP_VAR **transvar)
Definition: scip_var.c:1443
static SCIP_RETCODE createProbSimplifiedTest(SCIP *scip, SCIP_PROBDATA *probdata)
Definition: probdata_cyc.c:269
SCIPInterval sign(const SCIPInterval &x)
SCIP_Bool SCIPisZero(SCIP *scip, SCIP_Real val)
constraint handler for quadratic constraints
#define NULL
Definition: lpi_spx1.cpp:155
SCIP_Bool edgesExist(SCIP_VAR ****edgevars, int *states, int nstates)
#define SCIP_CALL(x)
Definition: def.h:364
SCIP_RETCODE SCIPaddBilinTermQuadratic(SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var1, SCIP_VAR *var2, SCIP_Real coef)
int SCIPcycGetNBins(SCIP *scip)
SCIP_RETCODE SCIPcopyDigraph(SCIP *scip, SCIP_DIGRAPH **targetdigraph, SCIP_DIGRAPH *sourcedigraph)
SCIP_RETCODE SCIPcreateConsSetpart(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_Bool initial, SCIP_Bool separate, SCIP_Bool enforce, SCIP_Bool check, SCIP_Bool propagate, SCIP_Bool local, SCIP_Bool modifiable, SCIP_Bool dynamic, SCIP_Bool removable, SCIP_Bool stickingatnode)
Definition: cons_setppc.c:9055
SCIP_RETCODE SCIPsetProbDelorig(SCIP *scip, SCIP_DECL_PROBDELORIG((*probdelorig)))
Definition: scip_prob.c:190
SCIP_RETCODE SCIPgetRealParam(SCIP *scip, const char *name, SCIP_Real *value)
Definition: scip_param.c:298
internal methods for problem variables
int phiinv(int k, int ncluster)
Definition: probdata_cyc.c:184
SCIP_Real SCIPinfinity(SCIP *scip)
#define SCIP_Bool
Definition: def.h:70
static SCIP_RETCODE createProbSimplified(SCIP *scip, SCIP_PROBDATA *probdata)
Definition: probdata_cyc.c:411
SCIP_RETCODE SCIPtransformVar(SCIP *scip, SCIP_VAR *var, SCIP_VAR **transvar)
Definition: scip_var.c:1353
SCIP_EXPORT SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:17672
#define MAX(x, y)
Definition: tclique_def.h:83
problem data for cycle clustering problem
Constraint handler for linear constraints in their most general form, .
SCIP_RETCODE SCIPaddVar(SCIP *scip, SCIP_VAR *var)
Definition: scip_prob.c:1666
static SCIP_DECL_PROBDELTRANS(probdeltransCyc)
SCIP_RETCODE SCIPsetProbTrans(SCIP *scip, SCIP_DECL_PROBTRANS((*probtrans)))
Definition: scip_prob.c:211
SCIP_RETCODE SCIPsetProbCopy(SCIP *scip, SCIP_DECL_PROBCOPY((*probcopy)))
Definition: scip_prob.c:296
SCIP_RETCODE SCIPaddLinearVarQuadratic(SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real coef)
static SCIP_RETCODE createVariables(SCIP *scip, SCIP_PROBDATA *probdata)
Definition: probdata_cyc.c:200
struct SCIP_ProbData SCIP_PROBDATA
Definition: type_prob.h:44
SCIP_RETCODE SCIPchgVarObj(SCIP *scip, SCIP_VAR *var, SCIP_Real newobj)
Definition: scip_var.c:4514
int SCIPcycGetNCluster(SCIP *scip)
SCIP_EXPORT SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:17662
int SCIPsnprintf(char *t, int len, const char *s,...)
Definition: misc.c:10590
SCIP_RETCODE SCIPreleaseVar(SCIP *scip, SCIP_VAR **var)
Definition: scip_var.c:1252
SCIP_DIGRAPH * SCIPcycGetEdgeGraph(SCIP *scip)
#define SCIP_Real
Definition: def.h:163
SCIP_RETCODE SCIPcreateConsBasicLogicor(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars)
static SCIP_RETCODE createProbQP(SCIP *scip, SCIP_PROBDATA *probdata)
Definition: probdata_cyc.c:627
#define SCIP_Longint
Definition: def.h:148
SCIP_RETCODE SCIPcycPrintSolutionValues(SCIP *scip, SCIP_SOL *sol)
SCIP_RETCODE SCIPaddCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip_prob.c:2764
void SCIPinfoMessage(SCIP *scip, FILE *file, const char *formatstr,...)
Definition: scip_message.c:199
SCIP_Bool isPartition(SCIP *scip, SCIP_Real **solclustering, int nbins, int ncluster)
Definition: probdata_cyc.c:48
SCIP_RETCODE SCIPreleaseCons(SCIP *scip, SCIP_CONS **cons)
Definition: scip_cons.c:1110
#define SCIPABORT()
Definition: def.h:336
SCIP_RETCODE SCIPsetProbDeltrans(SCIP *scip, SCIP_DECL_PROBDELTRANS((*probdeltrans)))
Definition: scip_prob.c:232
void SCIPdigraphFreeComponents(SCIP_DIGRAPH *digraph)
Definition: misc.c:8405
SCIP_VAR *** SCIPcycGetBinvars(SCIP *scip)
SCIP_Real ** SCIPcycGetCmatrix(SCIP *scip)
SCIP_RETCODE SCIPcreateProbBasic(SCIP *scip, const char *name)
Definition: scip_prob.c:170
static SCIP_DECL_PROBDELORIG(probdelorigCyc)
SCIP_RETCODE SCIPdigraphAddArc(SCIP_DIGRAPH *digraph, int startnode, int endnode, void *data)
Definition: misc.c:7553