Scippy

SCIP

Solving Constraint Integer Programs

heur_cycgreedy.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_cycgreedy.c
26 * @brief Greedy primal heuristic. States are assigned to clusters iteratively. At each iteration all possible
27 * assignments are computed and the one with the best change in objective value is selected.
28 * @author Leon Eifler
29 */
30
31#include "heur_cycgreedy.h"
32
33#include <assert.h>
34#include <string.h>
35#include <time.h>
36#include <stdlib.h>
37#include "scip/misc.h"
38#include "probdata_cyc.h"
39#include "scip/cons_and.h"
40
41#define HEUR_NAME "cycgreedy"
42#define HEUR_DESC "primal heuristic template"
43#define HEUR_DISPCHAR 'h'
44#define HEUR_PRIORITY 536870911
45#define HEUR_FREQ 1
46#define HEUR_FREQOFS 0
47#define HEUR_MAXDEPTH -1
48#define HEUR_TIMING SCIP_HEURTIMING_BEFORENODE
49#define HEUR_USESSUBSCIP FALSE /**< does the heuristic use a secondary SCIP instance? */
50
51/** primal heuristic data */
52struct SCIP_HeurData
53{
54 int lasteffectrootdepth;/**< index of the last solution for which oneopt was performed */
55 SCIP_Bool local; /**< the heuristic only computes assignments until any improvement is found */
56};
57
58/** calculate the current objective value for a q-matrix */
59static
61 SCIP* scip, /**< SCIP data structure */
62 SCIP_Real** qmatrix, /**< the irreversibility matrix*/
63 SCIP_Real scale, /**< the scaling parameter in the objective function */
64 int ncluster /**< the number of cluster*/
65 )
66{
67 SCIP_Real objective = 0.0;
68 int c;
69 int c2;
70
71 for( c = 0; c < ncluster; ++c )
72 {
73 c2 = ( c + 1 ) % ncluster;
74 objective += qmatrix[c][c2] - qmatrix[c2][c];
75 objective += scale * qmatrix[c][c];
76 }
77
78 /* if we have no transitions at all then irreversibility should be set to 0 */
79 return objective;
80}
81
82/** initialize the q-matrix from a given (possibly incomplete) clusterassignment */
83static
85 SCIP_Real** clusterassignment, /**< the matrix containing the (incomplete) clusterassignment */
86 SCIP_Real** qmatrix, /**< the returned matrix with the irreversibility between two clusters */
87 SCIP_Real** cmatrix, /**< the transition-matrix containg the probability-data */
88 int nbins, /**< the number of bins */
89 int ncluster /**< the number of possible clusters */
90 )
91{
92 int i;
93 int j;
94 int k;
95 int l;
96
97 for( k = 0; k < ncluster; ++k )
98 {
99 for( l = 0; l < ncluster; ++l )
100 {
101 qmatrix[k][l] = 0;
102
103 for( i = 0; i < nbins; ++i )
104 {
105 for( j = 0; j < nbins; ++j )
106 {
107 /* as -1 and 0 are both interpreted as 0, this check is necessary. Compute x_ik*x_jl*c_ij */
108 if( clusterassignment[i][k] < 1 || clusterassignment[j][l] < 1 )
109 continue;
110
111 qmatrix[k][l] += cmatrix[i][j];
112 }
113 }
114 }
115 }
116}
117
118/** update the irreversibility matrix, after the clusterassignment[newcluster][newbin] was either set
119 * from 0 to 1 or from 1 to 0
120 */
121static
123 SCIP_Real** clusterassignment, /**< the matrix containing the (incomplete) clusterassignment */
124 SCIP_Real** qmatrix, /**< the returned matrix with the irreversibility between two clusters */
125 SCIP_Real** cmatrix, /**< the transition-matrix containg the probability-data */
126 int newbin, /**< the bin to be added to the assignment */
127 int newcluster, /**< the bluster in which the bin was changed */
128 int nbins, /**< the number of bins */
129 int ncluster /**< the number of clusters */
130 )
131{
132 int bin;
133 int cluster;
134
135 for( cluster = 0; cluster < ncluster; ++cluster )
136 {
137 for( bin = 0; bin < nbins; ++bin )
138 {
139 /* multiplier is 1 if clusterassignment is 1, and 0 if it is 0 (set to 0) or -1 (unassigned) */
140 int temp = 0;
141 if( clusterassignment[bin][cluster] == 1 )
142 temp = 1;
143
144 if( cluster != newcluster )
145 {
146 qmatrix[newcluster][cluster] += temp * cmatrix[newbin][bin];
147 qmatrix[cluster][newcluster] += temp * cmatrix[bin][newbin];
148 }
149 else
150 {
151 if( bin == newbin )
152 qmatrix[newcluster][newcluster] += cmatrix[newbin][bin];
153 else
154 qmatrix[newcluster][newcluster] += (cmatrix[newbin][bin] + cmatrix[bin][newbin]) * temp;
155 }
156 }
157 }
158}
159
160/** get the temporary objective value bound after newbin would be added to newcluster
161 * but dont not change anything with the clustering
162 */
163static
165 SCIP* scip, /**< SCIP data structure */
166 SCIP_Real** qmatrix, /**< the irreversibility matrix */
167 SCIP_Real** cmatrix, /**< the transition matrix */
168 SCIP_Real** clusterassignment, /**< the clusterassignment */
169 int newbin, /**< the bin that would be added to cluster */
170 int newcluster, /**< the cluster the bin would be added to */
171 int nbins, /**< the number of bins */
172 int ncluster /**< the number of cluster */
173 )
174{
175 SCIP_Real obj;
176 SCIP_Real temp;
177 int i;
178
179 obj = getObjective(scip, qmatrix, SCIPcycGetScale(scip), ncluster);
180
181 /* the coh in cluster changes as well as the flow to the next and the previous cluster */
182 for( i = 0; i < nbins; ++i )
183 {
184 temp = (clusterassignment[i][phiinv(newcluster, ncluster)] < 1 ? 0 : 1);
185 obj += (cmatrix[i][newbin] - cmatrix[newbin][i]) * temp;
186 temp = (clusterassignment[i][phi(newcluster, ncluster)] < 1 ? 0 : 1);
187 obj -= (cmatrix[i][newbin] - cmatrix[newbin][i]) * temp;
188 temp = (clusterassignment[i][newcluster] < 1 ? 0 : 1);
189 obj += (cmatrix[i][newbin] + cmatrix[newbin][i]) * temp;
190 }
191
192 return obj;
193}
194
195/* find and assign the next unassigned bin to an appropriate cluster */
196static
198 SCIP* scip, /**< SCIP data structure */
199 SCIP_Bool localheur, /**< should the heuristic only compute local optimal assignment */
200 SCIP_Real** clusterassignment, /**< the matrix with the Clusterassignment */
201 SCIP_Real** cmatrix, /**< the transition matrix */
202 SCIP_Real** qmatrix, /**< the irreversibility matrix */
203 SCIP_Bool* isassigned, /**< TRUE, if the bin i was already assigned to a cluster*/
204 int nbins, /**< the number of bins*/
205 int ncluster, /**< the number of cluster*/
206 int* amountassigned, /**< the total amount of bins already assigned*/
207 int* binsincluster, /**< the number of bins currently in a cluster*/
208 SCIP_Real* objective /**< the objective */
209 )
210{
211 SCIP_Real* binobjective;
212 SCIP_Bool** clusterispossible;
213 int* bestcluster;
214 SCIP_Real tempobj;
216 int i;
217 int c;
218 int c1;
219 int c2;
220 int save = -1;
221 int ind = -1;
222
223 /* allocate memory */
224 SCIP_CALL( SCIPallocClearBufferArray(scip, &binobjective, nbins) );
225 SCIP_CALL( SCIPallocClearBufferArray(scip, &bestcluster, nbins) );
226 SCIP_CALL( SCIPallocClearBufferArray(scip, &clusterispossible, nbins) );
227
228 for( i = 0; i < nbins; ++i )
229 {
230 SCIP_CALL( SCIPallocClearBufferArray(scip, &clusterispossible[i], ncluster) ); /*lint !e866*/
231 }
232
233 /* make ceratin that each cluster is non-empty*/
234 for( c = 0; c < ncluster; ++c )
235 {
236 tempobj = 0;
237
238 if( binsincluster[c] == 0 )
239 {
240 for( i = 0; i < nbins; ++i )
241 {
242 /* if already assigned do nothing */
243 if( isassigned[i] )
244 continue;
245
246 /* check if assigning this state is better than the previous best state */
247 binobjective[i] = getTempObj(scip, qmatrix, cmatrix, clusterassignment, i, c, nbins, ncluster);
248
249 if( binobjective[i] > tempobj )
250 {
251 save = i;
252 tempobj = binobjective[i];
253 }
254
255 /* ensure that a state is assigned */
256 if( save == -1 )
257 save = i;
258 }
259
260 /* assign the found state to the cluster */
261 for( c1 = 0; c1 < ncluster; ++c1 )
262 {
263 clusterassignment[save][c1] = 0;
264 }
265
266 clusterassignment[save][c] = 1;
267 binsincluster[c]++;
268
269 assert(binsincluster[c] == 1);
270
271 isassigned[save] = TRUE;
272 *amountassigned += 1;
273
274 /* update the q-matrix */
275 updateIrrevMat(clusterassignment, qmatrix, cmatrix, save, c, nbins, ncluster);
276 }
277 }
278
279 /*phase 2: iteratively assign states such that at each iteration the highest objective improvement is achieved */
280 for( i = 0; i < nbins; ++i )
281 {
282 bestcluster[i] = 0;
283 binobjective[i] = -SCIPinfinity(scip);
284 }
285
286 for( i = 0; i < nbins; ++i )
287 {
288 if( isassigned[i] )
289 continue;
290
291 /* check which clusters the bin can be assigned to. -1 means unassigned, 0 means fixed to 0. */
292 for( c1 = 0; c1 < ncluster; ++c1 )
293 {
294 /* if assignment to i would violate abs-var assignment then set clusterpossible to FALSE */
295 if( 0 != clusterassignment[i][c1] )
296 clusterispossible[i][c1] = TRUE;
297 else
298 clusterispossible[i][c1] = FALSE;
299 }
300
301 /* calculate the irrevbound for all possible clusterassignments */
302 for( c2 = 0; c2 < ncluster; ++c2 )
303 {
304 if( !clusterispossible[i][c2] || clusterassignment[i][c2] == 0 )
305 continue;
306
307 /* temporarily assign i to c2 */
308 save = (int) clusterassignment[i][c2];
309 clusterassignment[i][c2] = 1;
310
311 /* save the best possible irrevbound for each bin */
312 tempobj = getTempObj(scip, qmatrix, cmatrix, clusterassignment, i, c2, nbins, ncluster);
313
314 /* check if this is an improvement compared to the best known assignment */
315 if( SCIPisGT(scip, tempobj, binobjective[i]) )
316 {
317 binobjective[i] = tempobj;
318 bestcluster[i] = c2;
319 }
320
321 clusterassignment[i][c2] = save;
322 }
323
324 /* if localheur is true, then the heuristic assigns a state as soon as any improvement is found */
325 if( localheur && SCIPisGT(scip, binobjective[i], *objective) )
326 break;
327 }
328
329 /* take the bin with the highest increase in irrev-bound */
330 for( i = 0; i < nbins; ++i )
331 {
332 if( SCIPisLT(scip, max, binobjective[i]) )
333 {
334 max = binobjective[i];
335 ind = i;
336 }
337 }
338
339 assert(!isassigned[ind] && ind > -1 && ind < nbins);
340
341 /* assign this bin to the found cluster */
342 for( c1 = 0; c1 < ncluster; ++c1 )
343 {
344 clusterassignment[ind][c1] = 0;
345 }
346
347 clusterassignment[ind][bestcluster[ind]] = 1;
348 binsincluster[bestcluster[ind]]++;
349 *amountassigned += 1;
350 isassigned[ind] = TRUE;
351
352 /* update the Irreversibility matrix */
353 updateIrrevMat(clusterassignment, qmatrix, cmatrix, ind, bestcluster[ind], nbins, ncluster);
354 *objective = getObjective(scip, qmatrix, SCIPcycGetScale(scip), ncluster);
355
356 /* free the allocated memory */
357 for( i = 0; i < nbins; ++i )
358 {
359 SCIPfreeBufferArray(scip, &(clusterispossible[i]));
360 }
361 SCIPfreeBufferArray(scip, &clusterispossible);
362 SCIPfreeBufferArray(scip, &bestcluster);
363 SCIPfreeBufferArray(scip, &binobjective);
364
365 return SCIP_OKAY;
366}
367
368/*
369 * Callback methods of primal heuristic
370 */
371
372/** copy method for primal heuristic plugins (called when SCIP copies plugins) */
373static
374SCIP_DECL_HEURCOPY(heurCopyCycGreedy)
375{ /*lint --e{715}*/
376 assert(scip != NULL);
377 assert(heur != NULL);
378 assert(strcmp(SCIPheurGetName(heur), HEUR_NAME) == 0);
379
380 /* call inclusion method of primal heuristic */
382
383 return SCIP_OKAY;
384}
385
386/** destructor of primal heuristic to free user data (called when SCIP is exiting) */
387static
388SCIP_DECL_HEURFREE(heurFreeCycGreedy)
389{ /*lint --e{715}*/
390 SCIP_HEURDATA* heurdata;
391
392 assert(heur != NULL);
393 assert(strcmp(SCIPheurGetName(heur), HEUR_NAME) == 0);
394 assert(scip != NULL);
395
396 /* free heuristic data */
397 heurdata = SCIPheurGetData(heur);
398
399 assert(heurdata != NULL);
400
401 SCIPfreeMemory(scip, &heurdata);
402 SCIPheurSetData(heur, NULL);
403
404 return SCIP_OKAY;
405}
406
407/** solving process deinitialization method of primal heuristic (called before branch and bound process data is freed) */
408static
409SCIP_DECL_HEUREXITSOL(heurExitsolCycGreedy)
410{ /*lint --e{715}*/
411 assert(heur != NULL);
412 assert(strcmp(SCIPheurGetName(heur), HEUR_NAME) == 0);
413
414 /* reset the timing mask to its default value */
416
417 return SCIP_OKAY;
418}
419
420/** initialization method of primal heuristic (called after problem was transformed) */
421static
422SCIP_DECL_HEURINIT(heurInitCycGreedy)
423{ /*lint --e{715}*/
424 SCIP_HEURDATA* heurdata;
425
426 assert(heur != NULL);
427 assert(scip != NULL);
428
429 /* get heuristic data */
430 heurdata = SCIPheurGetData(heur);
431 assert(heurdata != NULL);
432
433 /* initialize last solution index */
434 heurdata->lasteffectrootdepth = -1;
435
436 return SCIP_OKAY;
437}
438
439/** execution method of primal heuristic */
440static
441SCIP_DECL_HEUREXEC(heurExecCycGreedy)
442{ /*lint --e{715}*/
443 SCIP_Real** cmatrix; /* the transition matrixx */
444 SCIP_Real** qmatrix; /* the low-dimensional transition matrix between clusters */
445 SCIP_VAR*** binvars; /* SCIP variables */
446 SCIP_Real** clustering; /* matrix for the assignment of the binary variables */
447 int* binsincluster; /* amount of bins in a given cluster */
448 SCIP_Bool* isassigned; /* TRUE if a bin has already bin assigned to a cluster */
449 SCIP_HEURDATA* heurdata; /* the heurdata */
450 SCIP_SOL* sol; /* pointer to solution */
451 SCIP_Bool possible = TRUE; /* can the heuristic be run */
452 SCIP_Bool feasible = FALSE; /* is the solution feasible */
453 SCIP_Real obj = 0.0; /* objective value */
454 int amountassigned; /* total amount of bins assigned */
455 int nbins; /* number of bins */
456 int ncluster; /* number of cluster */
457 int i; /* running indices */
458 int j;
459 int c;
460
461 *result = SCIP_DIDNOTRUN;
462 amountassigned = 0;
463
464 /* for now: do not use heurisitc if weighted objective is used */
465 heurdata = SCIPheurGetData(heur);
466 if( SCIPgetEffectiveRootDepth(scip) == heurdata->lasteffectrootdepth )
467 return SCIP_OKAY;
468
469 heurdata->lasteffectrootdepth = SCIPgetEffectiveRootDepth(scip);
470
471 /* get the problem data from scip */
472 cmatrix = SCIPcycGetCmatrix(scip);
473 nbins = SCIPcycGetNBins(scip);
474 ncluster = SCIPcycGetNCluster(scip);
475 binvars = SCIPcycGetBinvars(scip);
476
477 assert(nbins > 0 && ncluster > 0);
478
479 /* allocate memory for the assignment */
480 SCIP_CALL( SCIPallocClearBufferArray(scip, &clustering, nbins) );
481 SCIP_CALL( SCIPallocClearBufferArray(scip, &binsincluster, ncluster) );
482 SCIP_CALL( SCIPallocClearBufferArray(scip, &qmatrix, ncluster) );
483 SCIP_CALL( SCIPallocClearBufferArray(scip, &isassigned, nbins) );
484
485 for ( i = 0; i < nbins; ++i )
486 {
487 if( i < ncluster )
488 {
489 SCIP_CALL( SCIPallocClearBufferArray(scip, &qmatrix[i], ncluster) ); /*lint !e866*/
490 }
491
492 SCIP_CALL( SCIPallocClearBufferArray(scip, &clustering[i], ncluster) ); /*lint !e866*/
493
494 for( j = 0; j < ncluster; ++j )
495 {
496 /* unassigned is set to -1 so we can differentiate unassigned and fixed in the branch and bound tree */
497 clustering[i][j] = -1;
498 }
499 }
500
501 /* get the already fixed bin-variables from scip. An assignment of -1 one means unassigned.
502 * 0 is fixed to 0, 1 is fixed to 1
503 */
504 for( i = 0; i < nbins; ++i )
505 {
506 for( j = 0; j < ncluster; ++j )
507 {
508 if( NULL == binvars[i][j] )
509 {
510 possible = FALSE;
511 break;
512 }
513
514 /* if the bounds determine a fixed binary variable, then fix the variable in the clusterassignment */
515 if( SCIPisEQ(scip, SCIPvarGetLbGlobal(binvars[i][j]), SCIPvarGetUbGlobal(binvars[i][j])) )
516 {
517 clustering[i][j] = SCIPvarGetLbGlobal(binvars[i][j]);
518
519 if( SCIPisEQ(scip, 1.0, clustering[i][j]) )
520 {
521 binsincluster[j]++;
522 isassigned[i] = TRUE;
523 amountassigned += 1;
524
525 for( c = 0; c < ncluster; ++c )
526 {
527 if( clustering[i][c] == -1 )
528 clustering[i][c] = 0;
529 }
530 }
531 }
532 }
533 }
534
535 /* check if the assignment violates paritioning, e.g. because we are in a subscip */
536 for( i = 0; i < nbins; ++i )
537 {
538 int amountzeros = 0;
539 int sum = 0;
540
541 for( j = 0; j < ncluster; ++j )
542 {
543 if( 0 == clustering[i][j] )
544 amountzeros++;
545 if( 1 == clustering[i][j] )
546 sum++;
547 }
548
549 if( ncluster == amountzeros || sum > 1 )
550 possible = FALSE;
551 }
552
553 if( amountassigned < nbins && possible )
554 {
555 /* initialize the qmatrix and the lower irreversibility bound */
556 computeIrrevMat(clustering, qmatrix, cmatrix, nbins, ncluster);
557 obj = getObjective(scip, qmatrix, SCIPcycGetScale(scip), ncluster);
558
559 /* assign bins iteratively until all bins are assigned */
560 while( amountassigned < nbins )
561 {
562 SCIP_CALL( assignNextBin(scip, heurdata->local, clustering, cmatrix, qmatrix,
563 isassigned, nbins, ncluster, &amountassigned, binsincluster, &obj ) );
564 }
565
566 /* assert that the assignment is valid in the sense that it is a partition of the bins.
567 * Feasibility is not checked in this method
568 */
569 assert(isPartition(scip,clustering, nbins, ncluster));
570
571 /* update the qmatrix */
572 computeIrrevMat(clustering, qmatrix, cmatrix, nbins, ncluster);
573
574 /* set the variables the problem to the found clustering and test feasibility */
575 SCIP_CALL( SCIPcreateSol(scip, &sol, heur) );
576 SCIP_CALL( assignVars( scip, sol, clustering, nbins, ncluster) );
577 SCIP_CALL( SCIPtrySolFree(scip, &sol, FALSE, TRUE, TRUE, TRUE, TRUE, &feasible) );
578 }
579
580 if( feasible )
581 *result = SCIP_FOUNDSOL;
582 else
583 *result = SCIP_DIDNOTFIND;
584
585 /* free allocated memory */
586 for ( i = 0; i < nbins; ++i )
587 {
588 SCIPfreeBufferArray(scip, &clustering[i]);
589
590 if( i < ncluster )
591 SCIPfreeBufferArray(scip, &qmatrix[i]);
592 }
593
594 SCIPfreeBufferArray(scip, &isassigned);
595 SCIPfreeBufferArray(scip, &qmatrix);
596 SCIPfreeBufferArray(scip, &binsincluster);
597 SCIPfreeBufferArray(scip, &clustering);
598
599 return SCIP_OKAY;
600}
601
602/*
603 * * primal heuristic specific interface methods
604 */
605
606/** creates the CycGreedy - primal heuristic and includes it in SCIP */
608 SCIP* scip /**< SCIP data structure */
609 )
610{
611 SCIP_HEURDATA* heurdata;
612 SCIP_HEUR* heur;
613
614 /* create greedy primal heuristic data */
615 SCIP_CALL( SCIPallocMemory(scip, &heurdata) );
616
617 /* include primal heuristic */
618
621 HEUR_MAXDEPTH, HEUR_TIMING, HEUR_USESSUBSCIP, heurExecCycGreedy, heurdata) );
622
623 assert(heur != NULL);
624
625 /* set non fundamental callbacks via setter functions */
626 SCIP_CALL( SCIPsetHeurCopy(scip, heur, heurCopyCycGreedy) );
627 SCIP_CALL( SCIPsetHeurFree(scip, heur, heurFreeCycGreedy) );
628 SCIP_CALL( SCIPsetHeurExitsol(scip, heur, heurExitsolCycGreedy) );
629 SCIP_CALL( SCIPsetHeurInit(scip, heur, heurInitCycGreedy) );
630
632 "localheur", "If set to true, heuristic assigns bins as soon as any improvement is found",
633 &heurdata->local, FALSE, TRUE, NULL, NULL) );
634
635 return SCIP_OKAY;
636}
Constraint handler for AND constraints, .
#define NULL
Definition: def.h:267
#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 SCIPaddBoolParam(SCIP *scip, const char *name, const char *desc, SCIP_Bool *valueptr, SCIP_Bool isadvanced, SCIP_Bool defaultvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:57
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_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 SCIPallocClearBufferArray(scip, ptr, num)
Definition: scip_mem.h:126
#define SCIPallocMemory(scip, ptr)
Definition: scip_mem.h:60
#define SCIPfreeBufferArray(scip, ptr)
Definition: scip_mem.h:136
#define SCIPfreeMemory(scip, ptr)
Definition: scip_mem.h:78
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 SCIPinfinity(SCIP *scip)
SCIP_Bool SCIPisGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPisLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
int SCIPgetEffectiveRootDepth(SCIP *scip)
Definition: scip_tree.c:127
SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:18088
SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:18078
static SCIP_Real getObjective(SCIP *scip, SCIP_Real **qmatrix, SCIP_Real scale, int ncluster)
static SCIP_DECL_HEURINIT(heurInitCycGreedy)
#define HEUR_TIMING
SCIP_RETCODE SCIPincludeHeurCycGreedy(SCIP *scip)
static SCIP_DECL_HEUREXEC(heurExecCycGreedy)
#define HEUR_FREQOFS
static SCIP_DECL_HEUREXITSOL(heurExitsolCycGreedy)
#define HEUR_DESC
#define HEUR_DISPCHAR
#define HEUR_MAXDEPTH
#define HEUR_PRIORITY
#define HEUR_NAME
static SCIP_DECL_HEURCOPY(heurCopyCycGreedy)
static void computeIrrevMat(SCIP_Real **clusterassignment, SCIP_Real **qmatrix, SCIP_Real **cmatrix, int nbins, int ncluster)
static SCIP_DECL_HEURFREE(heurFreeCycGreedy)
#define HEUR_FREQ
#define HEUR_USESSUBSCIP
static void updateIrrevMat(SCIP_Real **clusterassignment, SCIP_Real **qmatrix, SCIP_Real **cmatrix, int newbin, int newcluster, int nbins, int ncluster)
static SCIP_RETCODE assignNextBin(SCIP *scip, SCIP_Bool localheur, SCIP_Real **clusterassignment, SCIP_Real **cmatrix, SCIP_Real **qmatrix, SCIP_Bool *isassigned, int nbins, int ncluster, int *amountassigned, int *binsincluster, SCIP_Real *objective)
static SCIP_Real getTempObj(SCIP *scip, SCIP_Real **qmatrix, SCIP_Real **cmatrix, SCIP_Real **clusterassignment, int newbin, int newcluster, int nbins, int ncluster)
Greedy primal heuristic. States are assigned to clusters iteratively. At each iteration all possible ...
internal miscellaneous methods
SCIP_RETCODE assignVars(SCIP *scip, SCIP_SOL *sol, SCIP_Real **clustering, int nbins, int ncluster)
Definition: probdata_cyc.c:88
int SCIPcycGetNBins(SCIP *scip)
int phiinv(int k, int ncluster)
Definition: probdata_cyc.c:193
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
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
@ SCIP_OKAY
Definition: type_retcode.h:42
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:63