Scippy

SCIP

Solving Constraint Integer Programs

heur_cyckerlin.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 heur_cyckerlin.c
26 * @brief improvement heuristic that exchanges binary variables between clusters.
27 * Similar to the famous kernighan/lin heuristic for graph partitioning
28 * @author Leon Eifler
29 */
30
31/*---+---- 1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
32
33#include "heur_cyckerlin.h"
34
35#include <assert.h>
36#include <string.h>
37
38#include "probdata_cyc.h"
39#include "scip/pub_misc.h"
40
41#define HEUR_NAME "cyckerlin"
42#define HEUR_DESC "switch heuristic that tries to improve solution by trading states betweeen clusters"
43#define HEUR_DISPCHAR '@'
44#define HEUR_PRIORITY 500
45#define HEUR_FREQ 10
46#define HEUR_FREQOFS 0
47#define HEUR_MAXDEPTH -1
48#define MAXPERMUTATIONS 5
49#define DEFAULT_RANDSEED 177 /**< random seed */
50#define HEUR_TIMING SCIP_HEURTIMING_BEFORENODE
51#define HEUR_USESSUBSCIP FALSE /**< does the heuristic use a secondary SCIP instance? */
52
53struct SCIP_HeurData
54{
55 SCIP_SOL** candidates;
56 int ncandidates;
57 int candlength;
58};
59
60/*
61 * Local methods
62 */
63
64/** external method that adds a solution to the list of candidate-solutions that should be improved */
66 SCIP* scip, /**< SCIP data structure */
67 SCIP_SOL* sol /**< the given solution */
68 )
69{
70 SCIP_HEUR* heur;
71 SCIP_HEURDATA* heurdata;
72
73 heur = SCIPfindHeur(scip, "cyckerlin");
74
75 if( heur == NULL )
76 return SCIP_OKAY;
77
78 heurdata = SCIPheurGetData(heur);
79
80 assert(heurdata != NULL);
81 assert(sol != NULL);
82
83 /* realloc candidate array, if necessary */
84 if( heurdata->candlength - 1 <= heurdata->ncandidates )
85 {
86 SCIP_CALL( SCIPreallocMemoryArray(scip, &(heurdata->candidates), (SCIP_Longint) heurdata->candlength * 2) );
87 heurdata->candlength *= 2;
88 }
89
90 heurdata->candidates[heurdata->ncandidates] = sol;
91 heurdata->ncandidates++;
92
93 return SCIP_OKAY;
94}
95
96
97/** get the bin-var assignment from scip and save it as a matrix */
98static
100 SCIP* scip, /**< SCIP data structure */
101 SCIP_SOL* bestsol, /**< the solution */
102 SCIP_Real** solclustering, /**< matrix to save the bin-vars*/
103 SCIP_Bool** binfixed, /**< matrix to save if a bin is fixed in scip */
104 int* clusterofbin, /**< array containing the cluster of each bin */
105 int* nbinsincluster /**< number of bins in each cluster */
106 )
107{
108 SCIP_VAR*** binvars;
109 int nbins;
110 int ncluster;
111 int i;
112 int c;
113
114 assert(bestsol != NULL);
115
116 nbins = SCIPcycGetNBins(scip);
117 ncluster = SCIPcycGetNCluster(scip);
118 binvars = SCIPcycGetBinvars(scip);
119
120 assert(nbins > 0 && ncluster > 0 && nbins > ncluster);
121 assert(binvars != NULL);
122
123 /* get the bin-variable values from the solution */
124 for( i = 0; i < nbins; ++i )
125 {
126 for( c = 0; c < ncluster; ++c )
127 {
128 binfixed[i][c] = FALSE;
129
130 if( binvars[i][c] != NULL)
131 {
132 if( (SCIPisFeasEQ(scip, SCIPvarGetUbGlobal(binvars[i][c]), SCIPvarGetLbGlobal(binvars[i][c]))) )
133 binfixed[i][c] = TRUE;
134
135 solclustering[i][c] = SCIPgetSolVal(scip, bestsol, binvars[i][c]);
136
137 if( SCIPisFeasEQ(scip, solclustering[i][c], 1.0) )
138 {
139 clusterofbin[i] = c;
140 nbinsincluster[c]++;
141 }
142 }
143 }
144 }
145
146 return SCIP_OKAY;
147}
148
149/** Set a bin to a new cluster, update the qmatrix. */
150static
152 SCIP_Real** solclustering, /**< the matrix with the clustering of the bins */
153 SCIP_Real** cmatrix, /**< the transition matrix*/
154 SCIP_Real** qmatrix, /**< the matrix containing the transition probabilities between clusters*/
155 int newbin, /**< the bin to be changed*/
156 int newcluster, /**< the cluster where the bin is changed*/
157 SCIP_Bool setone, /**< TRUE if the assignment is switched from 0 to 1, FALSE if 1 to 0*/
158 int nbins, /**< the number of bins*/
159 int ncluster /**< the number of clusters*/
160 )
161{
162 int bin;
163 int cluster;
164 SCIP_Real sign = setone ? 1.0 : -1.0;
165
166 for( cluster = 0; cluster < ncluster; ++cluster )
167 {
168 for( bin = 0; bin < nbins; ++bin )
169 {
170 if( cluster != newcluster )
171 {
172 qmatrix[newcluster][cluster] += sign * solclustering[bin][cluster] * cmatrix[newbin][bin];
173 qmatrix[cluster][newcluster] += sign * solclustering[bin][cluster] * cmatrix[bin][newbin];
174 }
175 else
176 {
177 if( bin != newbin )
178 {
179 qmatrix[newcluster][newcluster] += sign * (cmatrix[newbin][bin] + cmatrix[bin][newbin])
180 * solclustering[bin][cluster];
181 }
182 }
183 }
184 }
185
186 solclustering[newbin][newcluster] = (sign + 1.0) / 2.0;
187}
188
189/** initialize the q-matrix from a given (possibly incomplete) clusterassignment */
190static
192 SCIP_Real** clustering, /**< the matrix containing the clusterassignment */
193 SCIP_Real** qmatrix, /**< the returned matrix the transition probability between clusters */
194 SCIP_Real** cmatrix, /**< the transition-matrix containg the probability-data */
195 int nbins, /**< the number of bins */
196 int ncluster /**< the number of possible clusters */
197 )
198{
199 int i;
200 int j;
201 int k;
202 int l;
203
204 for( k = 0; k < ncluster; ++k )
205 {
206 for( l = 0; l < ncluster; ++l )
207 {
208 qmatrix[k][l] = 0;
209 for( i = 0; i < nbins; ++i )
210 {
211 for( j = 0; j < nbins; ++j )
212 {
213 /* as -1 and 0 are both interpreted as 0, this check is necessary. Compute x_ik*x_jl*c_ij */
214 if( i == j )
215 continue;
216
217 qmatrix[k][l] += cmatrix[i][j] * clustering[i][k] * clustering[j][l];
218 }
219 }
220 }
221 }
222}
223
224/** calculate the current objective value for a q-matrix */
225static
227 SCIP* scip, /**< SCIP data structure */
228 SCIP_Real** qmatrix, /**< the irreversibility matrix*/
229 SCIP_Real scale, /**< the scaling parameter in the objective function */
230 int ncluster /**< the number of cluster*/
231 )
232{
233 SCIP_Real objective = 0.0;
234 int c;
235 int c2;
236
237 for( c = 0; c < ncluster; ++c )
238 {
239 c2 = ( c + 1 ) % ncluster;
240 objective += ( qmatrix[c][c2] - qmatrix[c2][c]);
241 objective += scale * qmatrix[c][c];
242 }
243
244 /* if we have no transitions at all then irreversibility should be set to 0 */
245 return objective;
246}
247
248/** exchange another bin to a different cluster. No bin may be changed twice */
249static
251 SCIP* scip, /**< SCIP data structure */
252 SCIP_Real** cmatrix, /**< the transition matrix */
253 SCIP_Real** qmatrix, /**< the irreversibility matrix */
254 SCIP_Real** clustering, /**< the clusterassignement */
255 SCIP_Bool** binfixed, /**< array containing information about fixedbins */
256 SCIP_Bool* binprocessed, /**< has the bin already been switched? */
257 int* clusterofbin, /**< contains the cluster each bin is in at the moment */
258 int* nbinsincluster, /**< number of bins in each cluster */
259 int* switchedbin, /**< the bins swithced in each iteration */
260 int* switchedcluster, /**< the cluster to witch the bin was assigned in each iteration */
261 SCIP_Real* switchbound, /**< the objective achieved in each iteration */
262 SCIP_Real* maxbound, /**< the best objective value so far */
263 int* bestlength, /**< the amount of switches with the best objective value so far */
264 int iteration /**< which iteration are we in */
265 )
266{
267 SCIP_Real irrevchg;
268 SCIP_Real cohchg;
269 SCIP_Real maxboundlocal;
270 SCIP_Real scale;
271 SCIP_Real oldobjective;
272 int bin;
273 int k;
274 int i;
275 int l;
276 int indkmin;
277 int indlmin;
278 int maxbin;
279 int maxcluster;
280 int nbins = SCIPcycGetNBins(scip);
281 int ncluster = SCIPcycGetNCluster(scip);
282
283 scale = SCIPcycGetScale(scip);
284 maxboundlocal = -SCIPinfinity(scip);
285 oldobjective = getObjective(scip, qmatrix, scale, ncluster);
286 maxbin = -1;
287 maxcluster = -1;
288
289 assert(isPartition(scip, clustering, nbins, ncluster));
290
291 for( bin = 0; bin < nbins; ++bin )
292 {
293 if( binprocessed[bin] || nbinsincluster[clusterofbin[bin]] == 1 )
294 continue;
295 k = clusterofbin[bin];
296
297 assert(SCIPisFeasEQ(scip, clustering[bin][k], 1.0));
298
299 /* calculate the irreversibility and coherence after bin was moved from k to l */
300 for( l = 0; l < ncluster; ++l )
301 {
302 if( binfixed[bin][k] || binfixed[bin][l] )
303 continue;
304
305 if( k == l )
306 continue;
307
308 if( k == 0 )
309 indkmin = ncluster - 1;
310 else
311 indkmin = k - 1;
312
313 if( l == 0 )
314 indlmin = ncluster - 1;
315 else
316 indlmin = l - 1;
317
318 irrevchg = 0;
319 cohchg = 0;
320
321 assert(SCIPisZero(scip, clustering[bin][l]));
322
323 for( i = 0; i < nbins; ++i )
324 {
325 /*irreversibility change */
326 irrevchg -= clustering[i][indkmin] * (cmatrix[i][bin] - cmatrix[bin][i]);
327 irrevchg -= clustering[i][(k+1) % ncluster] * (cmatrix[bin][i] - cmatrix[i][bin]);
328
329 clustering[bin][k] = 0;
330 clustering[bin][l] = 1;
331
332 irrevchg += clustering[i][indlmin] * (cmatrix[i][bin] - cmatrix[bin][i]);
333 irrevchg += clustering[i][(l+1) % ncluster] * (cmatrix[bin][i] - cmatrix[i][bin]);
334
335 clustering[bin][k] = 1;
336 clustering[bin][l] = 0;
337
338 /*coherence change */
339 if( i != bin )
340 {
341 cohchg -= clustering[i][k] * (cmatrix[bin][i]+ cmatrix[i][bin]);
342 cohchg += clustering[i][l] * (cmatrix[bin][i]+ cmatrix[i][bin]);
343 }
344 }
345
346 if( oldobjective + irrevchg + scale * cohchg > maxboundlocal )
347 {
348 maxboundlocal = oldobjective + irrevchg + scale * cohchg;
349 maxbin = bin;
350 maxcluster = l;
351 }
352 }
353 }
354
355 if( maxbin == -1 )
356 return FALSE;
357
358 assert(maxbin >= 0 && maxcluster >= 0);
359 assert(maxbin < nbins && maxcluster < ncluster);\
360
361 /* assign the exchange and update all saving-structures */
362 setBinToCluster(clustering, cmatrix, qmatrix, maxbin, clusterofbin[maxbin], FALSE, nbins, ncluster);
363 setBinToCluster(clustering, cmatrix, qmatrix, maxbin, maxcluster, TRUE, nbins, ncluster);
364
365 nbinsincluster[clusterofbin[maxbin]]--;
366 nbinsincluster[maxcluster]++;
367
368 clusterofbin[maxbin] = maxcluster;
369 binprocessed[maxbin] = TRUE;
370
371 switchedbin[iteration] = maxbin;
372 switchedcluster[iteration] = maxcluster;
373 switchbound[iteration] = getObjective(scip, qmatrix, scale, ncluster);
374
375 if( switchbound[iteration] > *maxbound )
376 {
377 *maxbound = switchbound[iteration];
378 *bestlength = iteration;
379 }
380
381 return TRUE;
382}
383
384/** Create a solution in scip from the clustering */
385static
387 SCIP* scip, /**< SCIP data structure */
388 SCIP_HEUR* heur, /**< heuristic pointer */
389 SCIP_Real** cmatrix, /**< the transition matrix */
390 SCIP_Real** qmatrix, /**< the projected transition matrix using the clustering */
391 SCIP_Bool** binfixed, /**< matrix that tells which bin-variables cannot be changed */
392 SCIP_Real** startclustering, /**< the start-assignment */
393 SCIP_RESULT* result, /**< result pointer */
394 int nbins, /**< the number of states */
395 int ncluster /**< the number of clusters */
396 )
397{
398 SCIP_SOL* bestsol;
399 SCIP_SOL* worksol;
400 SCIP_Real** clustering;
401 SCIP_Real** solclustering;
402 SCIP_Bool* binprocessed;
403 SCIP_Real max;
404 SCIP_Real objective;
405 SCIP_Real* switchbound;
406 SCIP_Real maxbound;
407 SCIP_Bool heurpossible = TRUE;
408 SCIP_Bool feasible;
409 int c;
410 int i;
411 int* nbinsincluster; /*lint !e771*/
412 int* switchedbin;
413 int* switchedcluster;
414 int* clusterofbin;
415 int bestlength;
416 int nrswitches;
417
418 /* allocate memory */
419 SCIP_CALL( SCIPallocBufferArray(scip, &binprocessed, nbins) );
420 SCIP_CALL( SCIPallocBufferArray(scip, &switchbound, nbins) );
421 SCIP_CALL( SCIPallocBufferArray(scip, &nbinsincluster, ncluster) );
422 SCIP_CALL( SCIPallocBufferArray(scip, &switchedbin, nbins) );
423 SCIP_CALL( SCIPallocBufferArray(scip, &switchedcluster, nbins) );
424 SCIP_CALL( SCIPallocBufferArray(scip, &clusterofbin, nbins) );
425 SCIP_CALL( SCIPallocClearBufferArray(scip, &clustering, nbins) );
426 SCIP_CALL( SCIPallocClearBufferArray(scip, &solclustering, nbins) );
427
428 for( i = 0; i < nbins; ++i )
429 {
430 SCIP_CALL( SCIPallocClearBufferArray(scip, &clustering[i], ncluster) ); /*lint !e866*/
431 SCIP_CALL( SCIPallocClearBufferArray(scip, &solclustering[i], ncluster) ); /*lint !e866*/
432 }
433
434 /* copy the solution so that we may change it and still keep the original one*/
435 for( c = 0; c < ncluster; ++c )
436 {
437 nbinsincluster[c] = 0;
438 }
439
440 for( i = 0; i < nbins; ++i )
441 {
442 for( c = 0; c < ncluster; ++c )
443 {
444 solclustering[i][c] = startclustering[i][c];
445
446 if( SCIPisFeasEQ(scip, startclustering[i][c], 1.0) )
447 {
448 clusterofbin[i] = c;
449 nbinsincluster[c]++; /*lint !e771*/
450 }
451 }
452
453 binprocessed[i] = FALSE;
454 }
455
456 bestsol = SCIPgetBestSol(scip);
457
458 while( heurpossible )
459 {
460 /* we run the heuristic until we cannot find any more improvement */
461 for( c = 0; c < ncluster; ++c )
462 {
463 nbinsincluster[c] = 0;
464 }
465
466 for( i = 0; i < nbins; ++i )
467 {
468 for( c = 0; c < ncluster; ++c )
469 {
470 clustering[i][c] = solclustering[i][c];
471
472 if( SCIPisFeasEQ(scip, solclustering[i][c], 1.0) )
473 {
474 clusterofbin[i] = c;
475 nbinsincluster[c]++;
476 }
477 }
478
479 binprocessed[i] = FALSE;
480 }
481
482 bestlength = -1;
483 nrswitches = nbins;
484
485 /* initialize qmatrix */
486 computeIrrevMat(solclustering, qmatrix, cmatrix, nbins, ncluster);
487 maxbound = SCIPgetSolOrigObj(scip, bestsol);
488
489 /* main part of the heuristic. States are switched until every state has been exchanged exactly once */
490 for( i = 0; i < nrswitches; ++i )
491 {
492 if( !switchNext(scip, cmatrix, qmatrix, clustering, binfixed, binprocessed, clusterofbin,
493 nbinsincluster, switchedbin, switchedcluster, switchbound, &maxbound, &bestlength, i) )
494 {
495 nrswitches = i;
496 break;
497 }
498 }
499
500 /* select the clustering with the best objective and reconstruct it from the start clustering */
501 for( i = 0; i <= bestlength; ++i )
502 {
503 for( c = 0; c < ncluster; ++c )
504 {
505 solclustering[switchedbin[i]][c] = 0; /*lint !e771*/
506 }
507
508 solclustering[switchedbin[i]][switchedcluster[i]] = 1; /*lint !e771*/
509 clusterofbin[switchedbin[i]] = switchedcluster[i];
510 }
511
512 computeIrrevMat(solclustering, qmatrix, cmatrix, nbins, ncluster);
513 max = getObjective(scip, qmatrix, SCIPcycGetScale(scip), ncluster);
514 objective = SCIPgetSolOrigObj(scip, bestsol);
515 feasible = FALSE;
516
517 /* if the solution is an improvement we add it to scip */
518 if( max > objective )
519 {
520 SCIP_CALL( SCIPcreateSol(scip, &worksol, heur) );
521
522 assert(isPartition(scip, solclustering, nbins, ncluster));
523
524 SCIP_CALL( assignVars(scip, worksol, solclustering, nbins, ncluster) );
525 SCIP_CALL( SCIPtrySolFree(scip, &worksol, FALSE, TRUE, TRUE, TRUE, TRUE, &feasible) );
526 }
527
528 if( feasible )
529 {
530 *result = SCIP_FOUNDSOL;
531 objective = max;
532 }
533 else
534 {
535 *result = SCIP_DIDNOTFIND;
536 heurpossible = FALSE;
537 }
538 }
539
540 /* free memory */
541 for( i = 0; i < nbins; ++i )
542 {
543 SCIPfreeBufferArray(scip, &solclustering[i]);
544 SCIPfreeBufferArray(scip, &clustering[i]);
545 }
546
547 SCIPfreeBufferArray(scip, &solclustering);
548 SCIPfreeBufferArray(scip, &clustering);
549 SCIPfreeBufferArray(scip, &clusterofbin);
550 SCIPfreeBufferArray(scip, &switchedcluster);
551 SCIPfreeBufferArray(scip, &switchedbin);
552 SCIPfreeBufferArray(scip, &nbinsincluster);
553 SCIPfreeBufferArray(scip, &switchbound);
554 SCIPfreeBufferArray(scip, &binprocessed);
555
556 return SCIP_OKAY;
557}
558
559/** method that randomly creates a different solution from a given solution. From each cluster, half the states are
560 * randomly selected and added to the next cluster. */
561static
563 SCIP* scip, /**< SCIP data structure */
564 SCIP_Real** startclustering, /**< the solution to be permuted */
565 SCIP_RANDNUMGEN* rnd, /**< a random number generator */
566 int nbins, /**< the number of states */
567 int ncluster /**< the number of clusters */
568 )
569{
570 int i;
571 int t;
572 int c;
573 int rndcluster;
574 int pushed;
575 int* binsincluster;
576 int **bins;
577
578 SCIP_CALL( SCIPallocBufferArray(scip, &binsincluster, ncluster) );
579 SCIP_CALL( SCIPallocBufferArray(scip, &bins, ncluster) );
580
581 for( t = 0; t < ncluster; ++t )
582 {
583 binsincluster[t] = 0;
584
585 for( i = 0; i < nbins; ++i )
586 {
587 if( SCIPisPositive(scip, startclustering[i][t]) )
588 binsincluster[t]++;
589 }
590
591 SCIP_CALL( SCIPallocClearBufferArray(scip, &bins[t], binsincluster[t]) ); /*lint !e866*/
592
593 c = 0;
594
595 for( i = 0; i < nbins; ++i )
596 {
597 if( SCIPisPositive(scip, startclustering[i][t]) )
598 {
599 bins[t][c] = i;
600 c++;
601 }
602 }
603 }
604
605 for( t = 0; t < ncluster; ++t )
606 {
607 pushed = 0;
608
609 while(pushed < binsincluster[t] / 2) /*lint !e771*/
610 {
611 rndcluster = bins[t][SCIPrandomGetInt(rnd, 0, binsincluster[t] - 1)];
612
613 if( rndcluster == nbins -1 )
614 continue;
615 if( SCIPisZero(scip, startclustering[rndcluster][t]) )
616 continue;
617
618 startclustering[rndcluster][t] = 0;
619 startclustering[rndcluster][phi(t,ncluster)] = 1;
620 pushed++;
621 }
622
623 SCIPfreeBufferArray(scip, &bins[t]);
624 }
625
627 SCIPfreeBufferArray(scip, &binsincluster);
628
629 return SCIP_OKAY;
630}
631
632/** executes the exchange heuristic for a given solution */
633static
635 SCIP* scip, /**< SCIP data structure */
636 SCIP_HEUR* heur, /**< heuristic pointer */
637 SCIP_SOL* sol, /**< given solution */
638 SCIP_RESULT* result /**< result pointer */
639 )
640{
641 SCIP_Real** startclustering;
642 SCIP_Bool** binfixed;
643 SCIP_Real** cmatrix;
644 SCIP_Real** qmatrix;
645 SCIP_RANDNUMGEN* rnd;
646 int* clusterofbin;
647 int* nbinsincluster;
648 int nbins;
649 int ncluster;
650 int i;
651 int c;
652
653 /* get problem variables */
654 nbins = SCIPcycGetNBins(scip);
655 ncluster = SCIPcycGetNCluster(scip);
656 cmatrix = SCIPcycGetCmatrix(scip);
658
659 assert(nbins >= 0);
660 assert(ncluster >= 0);
661
662 /* allocate Memory */
663 SCIP_CALL( SCIPallocClearBufferArray(scip, &startclustering, nbins) );
664 SCIP_CALL( SCIPallocClearBufferArray(scip, &clusterofbin, nbins) );
665 SCIP_CALL( SCIPallocClearBufferArray(scip, &binfixed, nbins) );
666 SCIP_CALL( SCIPallocClearBufferArray(scip, &qmatrix, ncluster) );
667 SCIP_CALL( SCIPallocClearBufferArray(scip, &nbinsincluster, ncluster) );
668
669 for( c = 0; c < ncluster; ++c )
670 {
671 SCIP_CALL( SCIPallocBufferArray(scip, &qmatrix[c], ncluster) ); /*lint !e866*/
672 }
673 for( i = 0; i < nbins; ++i )
674 {
675 SCIP_CALL( SCIPallocClearBufferArray(scip, &startclustering[i], ncluster) ); /*lint !e866*/
676 SCIP_CALL( SCIPallocClearBufferArray(scip, &binfixed[i], ncluster) ); /*lint !e866*/
677 }
678
679 /* get the solution values from scip */
680 SCIP_CALL( getSolutionValues(scip, sol, startclustering, binfixed, clusterofbin, nbinsincluster) );
681
682 if( isPartition(scip, startclustering, nbins, ncluster) )
683 {
684 SCIP_CALL( createSwitchSolution(scip, heur, cmatrix, qmatrix, binfixed, startclustering, result, nbins, ncluster) );
685 for( i = 0; i < MAXPERMUTATIONS; ++i )
686 {
687 SCIP_CALL( permuteStartSolution(scip, startclustering, rnd, nbins, ncluster) );
688
689 assert(isPartition(scip, startclustering, nbins, ncluster));
690
691 SCIP_CALL( createSwitchSolution(scip, heur, cmatrix, qmatrix, binfixed, startclustering, result, nbins, ncluster) );
692 }
693 }
694
695 /* free all data-structures */
696 for( i = 0; i < nbins; ++i )
697 {
698 SCIPfreeBufferArray(scip, &binfixed[i]);
699 SCIPfreeBufferArray(scip, &startclustering[i]);
700 }
701
702 for( c = 0; c < ncluster; ++c )
703 {
704 SCIPfreeBufferArray(scip, &qmatrix[c]);
705 }
706
707 SCIPfreeBufferArray(scip, &nbinsincluster);
708 SCIPfreeBufferArray(scip, &qmatrix);
709 SCIPfreeBufferArray(scip, &binfixed);
710 SCIPfreeBufferArray(scip, &clusterofbin);
711 SCIPfreeBufferArray(scip, &startclustering);
712
713 SCIPfreeRandom(scip, &rnd);
714
715 return SCIP_OKAY;
716}
717
718/*
719 * Callback methods of primal heuristic
720 */
721
722/** copy method for primal heuristic plugins (called when SCIP copies plugins) */
723static
724SCIP_DECL_HEURCOPY(heurCopyCyckerlin)
725{ /*lint --e{715}*/
726 assert(scip != NULL);
727 assert(heur != NULL);
728
729 assert(strcmp(SCIPheurGetName(heur), HEUR_NAME) == 0);
730
731 /* call inclusion method of primal heuristic */
733
734 return SCIP_OKAY;
735}
736
737/** destructor of primal heuristic to free user data (called when SCIP is exiting) */
738static
739SCIP_DECL_HEURFREE(heurFreeCyckerlin)
740{ /*lint --e{715}*/
741 SCIP_HEURDATA* heurdata;
742
743 assert(heur != NULL);
744 assert(strcmp(SCIPheurGetName(heur), HEUR_NAME) == 0);
745 assert(scip != NULL);
746
747 /* get heuristic data */
748 heurdata = SCIPheurGetData(heur);
749
750 assert(heurdata != NULL);
751
752 SCIPfreeMemoryArray(scip, &(heurdata->candidates));
753
754 /* free heuristic data */
755 SCIPfreeBlockMemory(scip, &heurdata);
756 SCIPheurSetData(heur, NULL);
757
758 return SCIP_OKAY;
759}
760
761/** solving process deinitialization method of primal heuristic (called before branch and bound process data is freed) */
762static
763SCIP_DECL_HEUREXITSOL(heurExitsolCyckerlin)
764{ /*lint --e{715}*/
765 SCIP_HEURDATA* heurdata;
766
767 assert(heur != NULL);
768 assert(scip != NULL);
769 assert(strcmp(SCIPheurGetName(heur), HEUR_NAME) == 0);
770
771 /* reset the timing mask to its default value */
773
774 /* get heuristic data */
775 heurdata = SCIPheurGetData(heur);
776
777 assert(heurdata != NULL);
778
779 heurdata->candlength = 0;
780
781 return SCIP_OKAY;
782}
783
784/** initialization method of primal heuristic (called after problem was transformed) */
785static
786SCIP_DECL_HEURINIT(heurInitCyckerlin)
787{ /*lint --e{715}*/
788 SCIP_HEURDATA* heurdata;
789
790 assert(heur != NULL);
791 assert(scip != NULL);
792
793 /* get heuristic's data */
794 heurdata = SCIPheurGetData(heur);
795
796 assert(heurdata != NULL);
797
798 heurdata->ncandidates = 0;
799 heurdata->candlength = 10;
800
801 SCIP_CALL( SCIPallocMemoryArray(scip, &(heurdata->candidates), heurdata->candlength) );
802
803 return SCIP_OKAY;
804}
805
806
807/** execution method of primal heuristic */
808static
809SCIP_DECL_HEUREXEC(heurExecCyckerlin)
810{ /*lint --e{715}*/
811 SCIP_Real objective;
812 SCIP_HEURDATA* heurdata;
813 int i;
814
815 assert(heur != NULL);
816 assert(scip != NULL);
817 assert(result != NULL);
818
819 *result = SCIP_DIDNOTRUN;
820 heurdata = SCIPheurGetData(heur);
821
822 assert(NULL != heurdata);
823
824 /* reset the timing mask to its default value (at the root node it could be different) */
825 if( SCIPgetNNodes(scip) > 1 )
827 for( i = 0; i < heurdata->ncandidates; i++ )
828 {
829 objective = SCIPgetSolOrigObj(scip, heurdata->candidates[i]);
830
831 if( !SCIPisZero(scip, objective) )
832 {
833 SCIP_CALL( runCyckerlin(scip, heur, heurdata->candidates[i], result) );
834 }
835
836 heurdata->candidates[i] = NULL;
837 }
838
839 heurdata->ncandidates = 0;
840 return SCIP_OKAY;
841}
842
843/*
844 * primal heuristic specific interface methods
845 */
846
847/** creates the oneopt primal heuristic and includes it in SCIP */
849 SCIP* scip /**< SCIP data structure */
850 )
851{
852 SCIP_HEUR* heur;
853 SCIP_HEURDATA* heurdata;
854
855 SCIP_CALL( SCIPallocBlockMemory(scip, &heurdata) );
856
857 /* include primal heuristic */
860 HEUR_MAXDEPTH, HEUR_TIMING, HEUR_USESSUBSCIP, heurExecCyckerlin, heurdata) );
861
862 assert(heur != NULL);
863
864 /* set non-NULL pointers to callback methods */
865 SCIP_CALL( SCIPsetHeurCopy(scip, heur, heurCopyCyckerlin) );
866 SCIP_CALL( SCIPsetHeurFree(scip, heur, heurFreeCyckerlin) );
867 SCIP_CALL( SCIPsetHeurExitsol(scip, heur, heurExitsolCyckerlin) );
868 SCIP_CALL( SCIPsetHeurInit(scip, heur, heurInitCyckerlin) );
869
870 return SCIP_OKAY;
871}
#define NULL
Definition: def.h:267
#define SCIP_Longint
Definition: def.h:158
#define SCIP_Bool
Definition: def.h:91
#define SCIP_Real
Definition: def.h:173
#define TRUE
Definition: def.h:93
#define FALSE
Definition: def.h:94
#define SCIP_CALL(x)
Definition: def.h:374
SCIP_RETCODE SCIPsetHeurExitsol(SCIP *scip, SCIP_HEUR *heur, SCIP_DECL_HEUREXITSOL((*heurexitsol)))
Definition: scip_heur.c:242
SCIP_RETCODE SCIPsetHeurCopy(SCIP *scip, SCIP_HEUR *heur, SCIP_DECL_HEURCOPY((*heurcopy)))
Definition: scip_heur.c:162
SCIP_HEURDATA * SCIPheurGetData(SCIP_HEUR *heur)
Definition: heur.c:1364
SCIP_RETCODE SCIPincludeHeurBasic(SCIP *scip, SCIP_HEUR **heur, const char *name, const char *desc, char dispchar, int priority, int freq, int freqofs, int maxdepth, SCIP_HEURTIMING timingmask, SCIP_Bool usessubscip, SCIP_DECL_HEUREXEC((*heurexec)), SCIP_HEURDATA *heurdata)
Definition: scip_heur.c:117
SCIP_RETCODE SCIPsetHeurFree(SCIP *scip, SCIP_HEUR *heur, SCIP_DECL_HEURFREE((*heurfree)))
Definition: scip_heur.c:178
void SCIPheurSetTimingmask(SCIP_HEUR *heur, SCIP_HEURTIMING timingmask)
Definition: heur.c:1493
SCIP_HEUR * SCIPfindHeur(SCIP *scip, const char *name)
Definition: scip_heur.c:258
SCIP_RETCODE SCIPsetHeurInit(SCIP *scip, SCIP_HEUR *heur, SCIP_DECL_HEURINIT((*heurinit)))
Definition: scip_heur.c:194
const char * SCIPheurGetName(SCIP_HEUR *heur)
Definition: heur.c:1453
void SCIPheurSetData(SCIP_HEUR *heur, SCIP_HEURDATA *heurdata)
Definition: heur.c:1374
#define SCIPreallocMemoryArray(scip, ptr, newnum)
Definition: scip_mem.h:70
#define SCIPallocMemoryArray(scip, ptr, num)
Definition: scip_mem.h:64
#define SCIPallocClearBufferArray(scip, ptr, num)
Definition: scip_mem.h:126
#define SCIPallocBufferArray(scip, ptr, num)
Definition: scip_mem.h:124
#define SCIPfreeBufferArray(scip, ptr)
Definition: scip_mem.h:136
#define SCIPfreeMemoryArray(scip, ptr)
Definition: scip_mem.h:80
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip_mem.h:108
#define SCIPallocBlockMemory(scip, ptr)
Definition: scip_mem.h:89
SCIP_SOL * SCIPgetBestSol(SCIP *scip)
Definition: scip_sol.c:2169
SCIP_RETCODE SCIPcreateSol(SCIP *scip, SCIP_SOL **sol, SCIP_HEUR *heur)
Definition: scip_sol.c:184
SCIP_RETCODE SCIPtrySolFree(SCIP *scip, SCIP_SOL **sol, SCIP_Bool printreason, SCIP_Bool completely, SCIP_Bool checkbounds, SCIP_Bool checkintegrality, SCIP_Bool checklprows, SCIP_Bool *stored)
Definition: scip_sol.c:3050
SCIP_Real SCIPgetSolOrigObj(SCIP *scip, SCIP_SOL *sol)
Definition: scip_sol.c:1300
SCIP_Real SCIPgetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var)
Definition: scip_sol.c:1217
SCIP_Longint SCIPgetNNodes(SCIP *scip)
SCIP_Real SCIPinfinity(SCIP *scip)
SCIP_Bool SCIPisFeasEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPisPositive(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisZero(SCIP *scip, SCIP_Real val)
SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:18088
SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:18078
void SCIPfreeRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen)
int SCIPrandomGetInt(SCIP_RANDNUMGEN *randnumgen, int minrandval, int maxrandval)
Definition: misc.c:10108
SCIP_RETCODE SCIPcreateRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen, unsigned int initialseed, SCIP_Bool useglobalseed)
static SCIP_Real getObjective(SCIP *scip, SCIP_Real **qmatrix, SCIP_Real scale, int ncluster)
static void setBinToCluster(SCIP_Real **solclustering, SCIP_Real **cmatrix, SCIP_Real **qmatrix, int newbin, int newcluster, SCIP_Bool setone, int nbins, int ncluster)
static SCIP_DECL_HEUREXITSOL(heurExitsolCyckerlin)
#define HEUR_TIMING
static void computeIrrevMat(SCIP_Real **clustering, SCIP_Real **qmatrix, SCIP_Real **cmatrix, int nbins, int ncluster)
#define HEUR_FREQOFS
#define HEUR_DESC
static SCIP_RETCODE createSwitchSolution(SCIP *scip, SCIP_HEUR *heur, SCIP_Real **cmatrix, SCIP_Real **qmatrix, SCIP_Bool **binfixed, SCIP_Real **startclustering, SCIP_RESULT *result, int nbins, int ncluster)
SCIP_RETCODE SCIPincludeHeurCycKerlin(SCIP *scip)
static SCIP_RETCODE runCyckerlin(SCIP *scip, SCIP_HEUR *heur, SCIP_SOL *sol, SCIP_RESULT *result)
#define HEUR_DISPCHAR
#define HEUR_MAXDEPTH
#define HEUR_PRIORITY
static SCIP_RETCODE getSolutionValues(SCIP *scip, SCIP_SOL *bestsol, SCIP_Real **solclustering, SCIP_Bool **binfixed, int *clusterofbin, int *nbinsincluster)
#define HEUR_NAME
static SCIP_DECL_HEUREXEC(heurExecCyckerlin)
#define DEFAULT_RANDSEED
static SCIP_DECL_HEURCOPY(heurCopyCyckerlin)
#define MAXPERMUTATIONS
#define HEUR_FREQ
#define HEUR_USESSUBSCIP
static SCIP_RETCODE permuteStartSolution(SCIP *scip, SCIP_Real **startclustering, SCIP_RANDNUMGEN *rnd, int nbins, int ncluster)
static SCIP_Bool switchNext(SCIP *scip, SCIP_Real **cmatrix, SCIP_Real **qmatrix, SCIP_Real **clustering, SCIP_Bool **binfixed, SCIP_Bool *binprocessed, int *clusterofbin, int *nbinsincluster, int *switchedbin, int *switchedcluster, SCIP_Real *switchbound, SCIP_Real *maxbound, int *bestlength, int iteration)
SCIP_RETCODE addCandSolCyckerlin(SCIP *scip, SCIP_SOL *sol)
static SCIP_DECL_HEURFREE(heurFreeCyckerlin)
static SCIP_DECL_HEURINIT(heurInitCyckerlin)
Improvement heuristic that trades bin-variables between clusters.
SCIP_RETCODE assignVars(SCIP *scip, SCIP_SOL *sol, SCIP_Real **clustering, int nbins, int ncluster)
Definition: probdata_cyc.c:88
int SCIPcycGetNBins(SCIP *scip)
SCIP_Real SCIPcycGetScale(SCIP *scip)
int SCIPcycGetNCluster(SCIP *scip)
SCIP_VAR *** SCIPcycGetBinvars(SCIP *scip)
SCIP_Real ** SCIPcycGetCmatrix(SCIP *scip)
SCIP_Bool isPartition(SCIP *scip, SCIP_Real **solclustering, int nbins, int ncluster)
Definition: probdata_cyc.c:57
problem data for cycle clustering problem
public data structures and miscellaneous methods
static SCIP_Real phi(SCIP *scip, SCIP_Real val, SCIP_Real lb, SCIP_Real ub)
Definition: sepa_eccuts.c:844
struct SCIP_HeurData SCIP_HEURDATA
Definition: type_heur.h:77
@ SCIP_DIDNOTRUN
Definition: type_result.h:42
@ SCIP_DIDNOTFIND
Definition: type_result.h:44
@ SCIP_FOUNDSOL
Definition: type_result.h:56
enum SCIP_Result SCIP_RESULT
Definition: type_result.h:61
@ SCIP_OKAY
Definition: type_retcode.h:42
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:63