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