Scippy

SCIP

Solving Constraint Integer Programs

pricer_rpa.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-2021 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not visit scipopt.org. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file pricer_rpa.c
17  * @brief Ringpacking variable pricer
18  * @author Benjamin Mueller
19  *
20  * This file implements the variable pricer which checks if variables negative reduced cost exist. See
21  * @ref RINGPACKING_PRICER for more details.
22  *
23  * @page RINGPACKING_PRICER Pricing new variables
24  *
25  * The task of the pricer is to search for new variables with negative reduced costs. For this, the following non-linear
26  * program is solved:
27  *
28  * \f[
29  * \begin{equation}
30  * \min_{P \in \mathcal{RP}} \left\{1 - \sum_{t \in \mathcal{T}} \lambda_t P_t\right\},
31  * \end{equation}
32  * \f]
33  *
34  * where \f$\lambda\f$ is given by the dual solution of the restricted master problem. See the
35  * \ref RINGPACKING_PROBLEM "problem description" for more details.
36  *
37  * This problem is very hard, but can be modeled as a weighted circle packing problem for a single rectangle. Therefore,
38  * we first use a simple greedy heuristic to solve the problem. If the heuristic fails, the MINLP is solved with
39  * conventional methods on a new \SCIP instance and a given time limit. If the problem can be solved and the optimal
40  * value is non-negative, the LP relaxation has been solved to optimality and what remains is ensuring integrality of
41  * the solution by the normal \SCIP framework. If, on the other hand, the best solution found by both methods is negative,
42  * we have found an improving pattern, whose corresponding variable needs to be added to the restricted master problem.
43  * It is possible (and not unlikely) that neither method succeeds in finding a pattern with negative solution value. In
44  * that case, we also exit the pricing loop, just as if we had found an optimal solution, and proceed with enforcing
45  * integrality resulting in a feasible primal solution to the whole problem.
46  *
47  * In case that the pricing problem cannot be solved to optimality, we cannot directly deduce a lower bound for the
48  * master problem. However, the following theorem by Farley shows how a valid dual bound can be computed from the
49  * LP solution and the pricing solution, see
50  * <a href="https://doi.org/10.1287/opre.38.5.922">A Note on Bounding a Class of Linear Programming Problems, Including
51  * Cutting Stock Problems</a> for more details.
52  *
53  */
54 
55 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
56 
57 #include <assert.h>
58 #include <string.h>
59 #include <sys/stat.h>
60 #include <sys/types.h>
61 
62 #include "scip/scipdefplugins.h"
63 #include "scip/scip.h"
64 #include "pricer_rpa.h"
65 #include "probdata_rpa.h"
66 #include "pattern.h"
67 
68 /**@name Pricer properties
69  *
70  * @{
71  */
72 
73 #define PRICER_NAME "ringpacking"
74 #define PRICER_DESC "pricer for ringpacking"
75 #define PRICER_PRIORITY 0
76 #define PRICER_DELAY TRUE /* only call pricer if all problem variables have non-negative reduced costs */
77 
78 /* default values of pricing parameters */
79 #define DEFAULT_PRICING_NLPTILIM 1e+20 /**< time limit for each pricing NLP */
80 #define DEFAULT_PRICING_NLPNODELIM 1000L /**< node limit for each pricing NLP */
81 #define DEFAULT_PRICING_HEURTILIM 1e+20 /**< time limit for each heuristic pricing */
82 #define DEFAULT_PRICING_HEURITERLIM 1000 /**< iteration limit for each heuristic pricing */
83 #define DEFAULT_PRICING_TOTALTILIM 1e+20 /**< total time limit for all pricing NLPs and heuristic calls */
84 
85 /**@} */
86 
87 #ifndef M_PI
88 #define M_PI 3.141592653589793238462643
89 #endif
90 
91 /*
92  * Data structures
93  */
94 
95 /** @brief Variable pricer data used in the \ref pricer_ringpacking.c "pricer" */
96 struct SCIP_PricerData
97 {
98  SCIP_RANDNUMGEN* randnumgen; /**< random number generator */
99  SCIP_Real timeleft; /**< time left for solving pricing problems (with NLP or heuristic) */
100 
101  /* parameters */
102  SCIP_Real nlptilim; /**< time limit for pricing NLP */
103  SCIP_Real heurtilim; /**< time limit for pricing heuristic */
104  SCIP_Longint nlpnodelim; /**< node limit for pricing NLP */
105  int heuriterlim; /**< iteration limit for pricing heuristic */
106 };
107 
108 
109 /**@name Local methods
110  *
111  * @{
112  */
113 
114 /** returns an upper bound on the density for n equal circles in a square (holds also for rectangles); this is a result
115  * of Groemer's formula (see, 'Ueber die Einlagerung von Kreisen in einen konvexen Bereich')
116  */
117 static
119 {
120  assert(n > 0);
121  if( n == 1 )
122  return M_PI / 4.0;
123  return (n * M_PI) / SQR(2.0 - SQRT(3.0) + SQRT(7.0 - M_PI + SQRT(3.0)*(2.0*n - 6.0 + M_PI)) );/*lint !e666*/
124 }
125 
126 /** helper function to count how many circles of a given type are needed */
127 static
129  SCIP* scip, /**< SCIP data structure */
130  SCIP_Real rext, /**< external radius */
131  int demand, /**< demand */
132  SCIP_Real width, /**< width of the rectangle */
133  SCIP_Real height, /**< height of the rectangle */
134  SCIP_Real lambda /**< objective coefficient of each circle of the given type */
135  )
136 {
137  SCIP_Real volrect;
138  SCIP_Real volcircle;
139  int ncircles;
140 
141  assert(!SCIPisFeasNegative(scip, lambda));
142 
143  /* objective coefficient of this circle type is zero -> not needed */
144  if( !SCIPisFeasPositive(scip, lambda) )
145  return 0;
146 
147  volrect = width * height;
148  volcircle = M_PI * SQR(rext);
149  assert(volcircle != 0.0);
150 
151  ncircles = (int)SCIPfeasFloor(scip, (getDensityUb(demand) * volrect) / volcircle);
152 
153  /* special cases where too big circles have a minimum distance to each other (in x direction) */
154  if( SCIPisGT(scip, rext, height / 4.0) )
155  {
156  SCIP_Real c = SQRT(4.0 * rext * height - SQR(height));
157  ncircles = (int)MIN(ncircles, 1 + (int)SCIPfloor(scip, (width - 2.0*rext) / c)); /*lint !e666*/
158  }
159  if( SCIPisGT(scip, rext, height / 6.0) && SCIPisLE(scip, rext, height / 4.0) )
160  {
161  SCIP_Real c = MIN(SQRT(3.0*rext*rext + rext * height - height * height / 4.0),
162  SQRT(8.0 * rext * height - height * height - 12.0 * rext * rext)); /*lint !e666*/
163  SCIP_Real k = SCIPfloor(scip, height / (2.0 * rext)) + 1;
164  SCIP_Real l = (width - 2.0 * rext) / c;
165  ncircles = (int)MIN(ncircles, k + l*(k-1) - 1);
166  }
167  assert(ncircles > 0);
168 
169  return MIN(ncircles, demand);
170 }
171 
172 /** adds a variable to the restricted master problem */
173 static
175  SCIP* scip, /**< SCIP data structure */
176  SCIP_PROBDATA* probdata, /**< problem data */
177  int* types, /**< types of elements */
178  SCIP_Real* xs, /**< x coordinate of each element */
179  SCIP_Real* ys, /**< y coordinate of each element */
180  SCIP_Bool* ispacked, /**< checks whether an element has been packed (might be NULL) */
181  int nelems /**< total number of elements */
182  )
183 {
184  SCIP_CONS** conss;
185  SCIP_PATTERN* pattern;
186  SCIP_VAR* var;
187  char name[SCIP_MAXSTRLEN];
188  char strtmp[SCIP_MAXSTRLEN];
189  int i;
190 
191  assert(probdata != NULL);
192  assert(types != NULL);
193  assert(xs != NULL);
194  assert(ys != NULL);
195  assert(nelems > 0);
196 
197  conss = SCIPprobdataGetPatternConss(probdata);
198  assert(conss != NULL);
199 
200  /* create rectangular pattern */
201  SCIP_CALL( SCIPpatternCreateRectangular(scip, &pattern) );
202 
203  /* create variable name */
204  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "r");
205  for( i = 0; i < nelems; ++i )
206  {
207  if( ispacked == NULL || ispacked[i] )
208  {
209  (void) SCIPsnprintf(strtmp, SCIP_MAXSTRLEN, "_%d", types[i]);
210  (void) strcat(name, strtmp);
211  SCIP_CALL( SCIPpatternAddElement(pattern, types[i], xs[i], ys[i]) );
212  }
213  }
214 
215  /* mark pattern to be packable */
217 
218  /* create and add variable */
219  SCIP_CALL( SCIPcreateVarBasic(scip, &var, name, 0.0, SCIPinfinity(scip), 1.0, SCIP_VARTYPE_INTEGER) );
223  SCIP_CALL( SCIPaddPricedVar(scip, var, 1.0) );
224  SCIPdebugMsg(scip, "added variable %s\n", name);
225 
226  /* add coefficients */
227  for( i = 0; i < nelems; ++i )
228  {
229  if( ispacked == NULL || ispacked[i] )
230  {
231  assert(types[i] >= 0 && types[i] < SCIPprobdataGetNTypes(probdata));
232  SCIP_CALL( SCIPaddCoefLinear(scip, conss[types[i]], var, 1.0) );
233  }
234  }
235 
236  /* add pattern and variable to the problem data */
237  SCIP_CALL( SCIPprobdataAddVar(scip, probdata, pattern, var) );
238 
239  /* release memory */
240  SCIP_CALL( SCIPreleaseVar(scip, &var) );
241  SCIPpatternRelease(scip, &pattern);
242 
243  return SCIP_OKAY;
244 }
245 
246 /* extracts and adds a variable with (hopefully) negative reduced costs */
247 static
249  SCIP* scip, /**< SCIP data structure */
250  SCIP_PROBDATA* probdata, /**< problem data */
251  SCIP* subscip, /**< sub-SCIP data structure */
252  SCIP_VAR** x, /**< x variables of sub-SCIP */
253  SCIP_VAR** y, /**< y variables of sub-SCIP */
254  SCIP_VAR** w, /**< w variables of sub-SCIP */
255  int* types, /**< type corresponding to each variable */
256  int nvars, /**< number of variables */
257  SCIP_Bool* success /**< pointer to store if we could add at least one variable with negative reduced costs */
258  )
259 {
260  SCIP_SOL* sol;
261  SCIP_Real* xs;
262  SCIP_Real* ys;
263  int* selectedtypes;
264  int nselected;
265  int i;
266 
267  assert(success != NULL);
268 
269  if( SCIPgetNSols(subscip) == 0 )
270  {
271  *success = FALSE;
272  return SCIP_OKAY;
273  }
274 
275  sol = SCIPgetBestSol(subscip);
276  assert(sol != NULL);
277  SCIPdebugMsg(scip, "found column with reduced cost = %f\n", 1.0 + SCIPgetSolOrigObj(subscip, sol));
278 
279  /* reduced cost is non-negative */
280  if( SCIPisFeasGE(subscip, 1.0 + SCIPgetSolOrigObj(subscip, sol), 0.0) )
281  {
282  *success = FALSE;
283  return SCIP_OKAY;
284  }
285  else
286  *success = TRUE;
287 
288  /* allocate memory */
289  SCIP_CALL( SCIPallocBufferArray(scip, &selectedtypes, nvars) );
290  SCIP_CALL( SCIPallocBufferArray(scip, &xs, nvars) );
291  SCIP_CALL( SCIPallocBufferArray(scip, &ys, nvars) );
292 
293  /* scan which circles have been selected */
294  nselected = 0;
295  for( i = 0; i < nvars; ++i )
296  {
297  SCIP_Real solval = SCIPgetSolVal(subscip, sol, w[i]);
298 
299  if( solval >= 0.5 )
300  {
301  selectedtypes[nselected] = types[i];
302  xs[nselected] = SCIPgetSolVal(subscip, sol, x[i]);
303  ys[nselected] = SCIPgetSolVal(subscip, sol, y[i]);
304  ++nselected;
305  }
306  }
307  assert(nselected > 0); /* otherwise the reduced cost can not be negative */
308 
309  /* add variable to main SCIP */
310  SCIP_CALL( addVariable(scip, probdata, selectedtypes, xs, ys, NULL, nselected) );
311 
312  /* free memory */
313  SCIPfreeBufferArray(scip, &ys);
314  SCIPfreeBufferArray(scip, &xs);
315  SCIPfreeBufferArray(scip, &selectedtypes);
316 
317  return SCIP_OKAY;
318 }
319 
320 /** array to compute the score of each element */
321 static
323  SCIP* scip, /**< SCIP data structure */
324  SCIP_Real* rexts, /**< external radii for each type */
325  SCIP_RANDNUMGEN* randnumgen, /**< random number generator */
326  int* elements, /**< type of each element */
327  int nelements, /**< total number of elements */
328  SCIP_Real* lambdas, /**< dual multipliers for each type */
329  SCIP_Real* scores, /**< array to store the score of each element */
330  int iter, /**< iteration round */
331  int iterlim /**< total iteration limit */
332  )
333 {
334  int i;
335 
336  assert(iter < iterlim);
337 
338  for( i = 0; i < nelements; ++i )
339  {
340  int elemtype = elements[i];
341  SCIP_Real rext = rexts[elemtype];
342 
343  /* use items with largest multipliers first */
344  if( iter == 0 )
345  scores[i] = lambdas[elemtype];
346 
347  /* use largest elements first */
348  else if( iter == 1 )
349  scores[i] = rext;
350 
351  /* use smallest elements first */
352  else if( iter == 2 )
353  scores[i] = -rext;
354 
355  /* use [0,1] * radius^2 */
356  else if( iter <= iterlim * 0.1 )
357  scores[i] = SCIPrandomGetReal(randnumgen, 0.0, 1.0) * rext * rext;
358 
359  /* use [0,1] * radius * lambda */
360  else if( iter <= iterlim * 0.4 )
361  scores[i] = SCIPrandomGetReal(randnumgen, 0.0, 1.0) * rext * lambdas[elemtype];
362 
363  /* use [-1,1] * lambda / radius */
364  else if( iter <= iterlim * 0.8 )
365  scores[i] = SCIPrandomGetReal(randnumgen, -1.0, 1.0) * rext * lambdas[elemtype];
366 
367  /* use a random order */
368  else
369  scores[i] = SCIPrandomGetReal(randnumgen, 0.0, 1.0);
370  }
371 }
372 
373 /** tries to find a column with negative reduced cost by using a greedy packing heuristic */
374 static
376  SCIP* scip, /**< SCIP data structure */
377  SCIP_PROBDATA* probdata, /**< problem data */
378  SCIP_PRICERDATA* pricerdata, /**< pricer data */
379  SCIP_Real* lambdas, /**< dual multipliers for pattern constraints */
380  SCIP_Real timelim, /**< time limit */
381  int iterlim, /**< iteration limit */
382  SCIP_Bool* addedvar /**< pointer to store whether a variable with negative reduced cost has been added */
383  )
384 {
385  SCIP_Real* scores;
386  SCIP_Real* rexts;
387  SCIP_Real* xs;
388  SCIP_Real* ys;
389  SCIP_Bool* ispacked;
390  int* demands;
391  int* elements;
392  SCIP_Real width;
393  SCIP_Real height;
394  SCIP_Real timestart;
395  SCIP_Real bestredcosts;
396  SCIP_Real bestvol;
397  int nelements;
398  int ntypes;
399  int npacked;
400  int niters;
401  int t;
402 
403  assert(pricerdata != NULL);
404  assert(addedvar != NULL);
405 
406  *addedvar = FALSE;
407  niters = 0;
408  timestart = SCIPgetTotalTime(scip);
409  bestredcosts = 0.0;
410  bestvol = 0.0;
411 
412  /* get problem data */
413  rexts = SCIPprobdataGetRexts(probdata);
414  demands = SCIPprobdataGetDemands(probdata);
415  width = SCIPprobdataGetWidth(probdata);
416  height = SCIPprobdataGetHeight(probdata);
417  ntypes = SCIPprobdataGetNTypes(probdata);
418 
419  /* compute the total number of elements that need to be considered */
420  nelements = 0;
421  for( t = 0; t < ntypes; ++t )
422  nelements += getNCircles(scip, rexts[t], demands[t], width, height, lambdas[t]);
423 
424  /* allocate enough memory */
425  SCIP_CALL( SCIPallocBufferArray(scip, &elements, nelements) );
426  SCIP_CALL( SCIPallocBufferArray(scip, &xs, nelements) );
427  SCIP_CALL( SCIPallocBufferArray(scip, &ys, nelements) );
428  SCIP_CALL( SCIPallocBufferArray(scip, &scores, nelements) );
429  SCIP_CALL( SCIPallocBufferArray(scip, &ispacked, nelements) );
430 
431  /* create entry for each element */
432  nelements = 0;
433  for( t = 0; t < ntypes; ++t )
434  {
435  int i;
436  int n;
437 
438  n = getNCircles(scip, rexts[t], demands[t], width, height, lambdas[t]);
439 
440  for( i = 0; i < n; ++i )
441  {
442  elements[nelements] = t;
443  ++nelements;
444  }
445  }
446 
447  /* main loop */
448  while( niters < iterlim
449  && SCIPgetTotalTime(scip) - timestart <= timelim
450  && !SCIPisStopped(scip) )
451  {
452  SCIP_Real redcosts = 1.0;
453  SCIP_Real vol = 0.0;
454  int i;
455 
456  /* compute scores */
457  computeScores(scip, rexts, pricerdata->randnumgen, elements, nelements, lambdas, scores, niters, iterlim);
458 
459  /* sort elements in non-increasing order */
460  SCIPsortDownRealInt(scores, elements, nelements);
461 
462  /* call heuristic */
463  SCIPpackCirclesGreedy(scip, rexts, xs, ys, -1.0, width, height, ispacked, elements, nelements,
464  SCIP_PATTERNTYPE_RECTANGULAR, &npacked, niters);
465 
466  /* compute reduced costs */
467  for( i = 0; i < nelements; ++i )
468  {
469  if( ispacked[i] )
470  {
471  redcosts -= lambdas[elements[i]];
472  vol += rexts[elements[i]];
473  }
474  }
475 
476  /* add pattern if reduced costs are negative */
477  if( SCIPisFeasLT(scip, redcosts, bestredcosts) || (SCIPisGT(scip, vol, bestvol) && SCIPisFeasNegative(scip, redcosts)) )
478  {
479  SCIPdebugMsg(scip, "pricing heuristic found column with reduced costs %g and volume %g after %d iterations\n", redcosts, vol, niters + 1);
480 
481  SCIP_CALL( addVariable(scip, probdata, elements, xs, ys, ispacked, nelements) );
482  *addedvar = TRUE;
483  bestredcosts = MIN(bestredcosts, redcosts);
484  bestvol = MAX(bestvol, vol);
485  }
486 
487  ++niters;
488  }
489 
490  /* free memory */
491  SCIPfreeBufferArray(scip, &ispacked);
492  SCIPfreeBufferArray(scip, &scores);
493  SCIPfreeBufferArray(scip, &ys);
494  SCIPfreeBufferArray(scip, &xs);
495  SCIPfreeBufferArray(scip, &elements);
496 
497  return SCIP_OKAY;
498 }
499 
500 /** auxiliary function for solving the pricing problem exactly */
501 static
503  SCIP* scip, /**< SCIP data structure */
504  SCIP_PROBDATA* probdata, /**< problem data */
505  SCIP_Real* lambdas, /**< dual multipliers for pattern constraints */
506  SCIP_Real timelim, /**< time limit */
507  SCIP_Longint nodelim, /**< node limit */
508  SCIP_Bool* addedvar, /**< pointer to store whether a variable with negative reduced cost has been added */
509  SCIP_STATUS* solstat, /**< pointer to store the solution status */
510  SCIP_Real* dualbound /**< pointer to store the dual bound */
511  )
512 {
513  SCIP* subscip;
514  SCIP_VAR** x;
515  SCIP_VAR** y;
516  SCIP_VAR** w;
517  SCIP_VAR* quadvars1[6];
518  SCIP_VAR* quadvars2[6];
519  SCIP_VAR* linvars[2];
520  SCIP_Real quadcoefs[6];
521  SCIP_Real lincoefs[2];
522  char name[SCIP_MAXSTRLEN];
523  SCIP_CONS* cons;
524  SCIP_Real* rexts;
525  SCIP_Real* vols;
526  int* types;
527  int* demands;
528  SCIP_Real width;
529  SCIP_Real height;
530  int nvars;
531  int ntypes;
532  int pos;
533  int t;
534  int i;
535 
536  assert(probdata != NULL);
537  assert(lambdas != NULL);
538  assert(nodelim >= -1L);
539  assert(addedvar != NULL);
540  assert(solstat != NULL);
541  assert(dualbound != NULL);
542 
543  *addedvar = FALSE;
544  *solstat = SCIP_STATUS_UNKNOWN;
545  *dualbound = -SCIPinfinity(scip);
546 
547  /* no time left */
548  if( timelim <= 0.0 )
549  return SCIP_OKAY;
550 
551  width = SCIPprobdataGetWidth(probdata);
552  height = SCIPprobdataGetHeight(probdata);
553  ntypes = SCIPprobdataGetNTypes(probdata);
554  rexts = SCIPprobdataGetRexts(probdata);
555  demands = SCIPprobdataGetDemands(probdata);
556  assert(ntypes > 0);
557 
558  /* create a sub-SCIP */
559  SCIP_CALL( SCIPcreate(&subscip) );
560  SCIP_CALL( SCIPcreateProbBasic(subscip, "pricing problem") );
562 
563  /* set heuristics to aggressive */
565  SCIP_CALL( SCIPsetIntParam(subscip, "heuristics/mpec/freq", -1) );
566 
567 #ifndef SCIP_DEBUG
568  SCIPsetMessagehdlrQuiet(subscip, TRUE);
569 #endif
570 
571  SCIP_CALL( SCIPsetRealParam(subscip, "limits/time", timelim) );
572  SCIP_CALL( SCIPsetLongintParam(subscip, "limits/nodes", nodelim) );
573 
574  /* count how many variables are needed */
575  nvars = 0;
576  for( t = 0; t < ntypes; ++t )
577  {
578  nvars += getNCircles(scip, rexts[t], demands[t], width, height, lambdas[t]);
579  SCIPdebugMsg(scip, "use %d/%d circles of type %d\n", getNCircles(scip, rexts[t], demands[t], width, height, lambdas[t]), demands[t], t);
580  }
581  assert(nvars > 0);
582 
583  /* allocate enough memory */
584  SCIP_CALL( SCIPallocBlockMemoryArray(subscip, &types, nvars) );
585  SCIP_CALL( SCIPallocBlockMemoryArray(subscip, &vols, nvars) );
586  SCIP_CALL( SCIPallocBlockMemoryArray(subscip, &x, nvars) );
587  SCIP_CALL( SCIPallocBlockMemoryArray(subscip, &y, nvars) );
588  SCIP_CALL( SCIPallocBlockMemoryArray(subscip, &w, nvars) );
589 
590  /* create variables */
591  pos = 0;
592  for( t = 0; t < ntypes; ++t )
593  {
594  SCIP_Real obj = lambdas[t];
595  int ncircles = getNCircles(scip, rexts[t], demands[t], width, height, lambdas[t]);
596  int k;
597 
598  for( k = 0; k < ncircles; ++k )
599  {
600  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "x_%d", pos);
601  SCIP_CALL( SCIPcreateVarBasic(subscip, &x[pos], name, rexts[t], width - rexts[t], 0.0, SCIP_VARTYPE_CONTINUOUS) );
602  SCIP_CALL( SCIPaddVar(subscip, x[pos]) );
603 
604  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "y_%d", pos);
605  SCIP_CALL( SCIPcreateVarBasic(subscip, &y[pos], name, rexts[t], height - rexts[t], 0.0, SCIP_VARTYPE_CONTINUOUS) );
606  SCIP_CALL( SCIPaddVar(subscip, y[pos]) );
607 
608  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "w_%d", pos);
609  SCIP_CALL( SCIPcreateVarBasic(subscip, &w[pos], name, 0.0, 1.0, -obj, SCIP_VARTYPE_BINARY) );
610  SCIP_CALL( SCIPaddVar(subscip, w[pos]) );
611 
612  vols[pos] = SQR(rexts[t]) * M_PI;
613  types[pos] = t;
614  ++pos;
615  }
616  }
617 
618  /* create non-overlapping constraints */
619  for( i = 0; i < nvars; ++i )
620  {
621  int j;
622  int type1;
623 
624  type1 = types[i];
625  assert(type1 >= 0 && type1 < ntypes);
626 
627  for( j = i+1; j < nvars; ++j )
628  {
629  SCIP_Real c;
630  int types2;
631 
632  types2 = types[j];
633  assert(types2 >= 0 && types2 < ntypes);
634 
635  c = (rexts[type1] + rexts[types2]) * (rexts[type1] + rexts[types2]);
636 
637  /* linear part */
638  linvars[0] = w[i]; lincoefs[0] = -c;
639  linvars[1] = w[j]; lincoefs[1] = -c;
640 
641  /* quadratic part */
642  quadvars1[0] = x[i]; quadvars2[0] = x[i]; quadcoefs[0] = 1.0;
643  quadvars1[1] = x[i]; quadvars2[1] = x[j]; quadcoefs[1] = -2.0;
644  quadvars1[2] = x[j]; quadvars2[2] = x[j]; quadcoefs[2] = 1.0;
645  quadvars1[3] = y[i]; quadvars2[3] = y[i]; quadcoefs[3] = 1.0;
646  quadvars1[4] = y[i]; quadvars2[4] = y[j]; quadcoefs[4] = -2.0;
647  quadvars1[5] = y[j]; quadvars2[5] = y[j]; quadcoefs[5] = 1.0;
648 
649  /* create quadratic constraint */
650  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "nonoverlap_%d_%d", i, j);
651  SCIP_CALL( SCIPcreateConsBasicQuadratic(subscip, &cons, name, 2, linvars, lincoefs, 6, quadvars1, quadvars2,
652  quadcoefs, -c, SCIPinfinity(subscip)) );
653 
654  /* add and release constraint */
655  SCIP_CALL( SCIPaddCons(subscip, cons) );
656  SCIP_CALL( SCIPreleaseCons(subscip, &cons) );
657  }
658  }
659 
660  /* w_i >= w_{i+1} if type(i) == type(i+1) */
661  for( i = 0; i < nvars - 1; ++i )
662  {
663  int type1;
664  int type2;
665 
666  type1 = types[i];
667  type2 = types[i+1];
668 
669  if( type1 != type2 )
670  continue;
671 
672  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "cons_w_%d_%d", i, i+1);
673 
674  SCIP_CALL( SCIPcreateConsBasicLinear(subscip, &cons, name, 0, NULL, NULL,
675  0.0, SCIPinfinity(subscip)) );
676  SCIP_CALL( SCIPaddCoefLinear(subscip, cons, w[i], 1.0) );
677  SCIP_CALL( SCIPaddCoefLinear(subscip, cons, w[i+1], -1.0) );
678 
679  /* add and release constraint */
680  SCIP_CALL( SCIPaddCons(subscip, cons) );
681  SCIP_CALL( SCIPreleaseCons(subscip, &cons) );
682  }
683 
684  /* x_i <= x_{i+1} if type(i) == type(i+1) */
685  for( i = 0; i < nvars - 1; ++i )
686  {
687  int type1;
688  int type2;
689 
690  type1 = types[i];
691  type2 = types[i+1];
692 
693  if( type1 != type2 )
694  continue;
695 
696  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "symmetry_%d_%d", i, i+1);
697  SCIP_CALL( SCIPcreateConsBasicLinear(subscip, &cons, name, 0, NULL, NULL,
698  0.0, SCIPinfinity(subscip)) );
699  SCIP_CALL( SCIPaddCoefLinear(subscip, cons, x[i], -1.0) );
700  SCIP_CALL( SCIPaddCoefLinear(subscip, cons, x[i+1], 1.0) );
701 
702  /* add and release constraint */
703  SCIP_CALL( SCIPaddCons(subscip, cons) );
704  SCIP_CALL( SCIPreleaseCons(subscip, &cons) );
705  }
706 
707  /* sum_{i} vol_i w_i <= W*H */
708  (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "volume");
709  SCIP_CALL( SCIPcreateConsBasicLinear(subscip, &cons, name, nvars, w, vols, 0.0, width * height) );
710  SCIP_CALL( SCIPaddCons(subscip, cons) );
711  SCIP_CALL( SCIPreleaseCons(subscip, &cons) );
712 
713  /* solve the pricing problem */
714  SCIPdebugMsg(scip, "----------------------- solve pricing problem -----------------------\n");
715  SCIP_CALL( SCIPsolve(subscip) );
716  SCIPdebugMsg(scip, "---------------------------------------------------------------------\n");
717 
718  /* check solution status */
719  *dualbound = SCIPgetDualbound(subscip);
720  *solstat = SCIPgetStatus(subscip);
721 
722  /* add variable with negative reduced costs */
723  SCIP_CALL( extractVariablesMINLP(scip, probdata, subscip, x, y, w, types, nvars, addedvar) );
724 
725  /* free variables */
726  for( i = 0; i < nvars; ++i )
727  {
728  SCIP_CALL( SCIPreleaseVar(subscip, &w[i]) );
729  SCIP_CALL( SCIPreleaseVar(subscip, &y[i]) );
730  SCIP_CALL( SCIPreleaseVar(subscip, &x[i]) );
731  }
732 
733  /* free memory */
734  SCIPfreeBlockMemoryArray(subscip, &w, nvars);
735  SCIPfreeBlockMemoryArray(subscip, &y, nvars);
736  SCIPfreeBlockMemoryArray(subscip, &x, nvars);
737  SCIPfreeBlockMemoryArray(subscip, &vols, nvars);
738  SCIPfreeBlockMemoryArray(subscip, &types, nvars);
739  SCIP_CALL( SCIPfree(&subscip) );
740 
741  return SCIP_OKAY;
742 }
743 
744 /**@} */
745 
746 /**name Callback methods
747  *
748  * @{
749  */
750 
751 /** destructor of variable pricer to free user data (called when SCIP is exiting) */
752 static
753 SCIP_DECL_PRICERFREE(pricerFreeRingpacking)
754 { /*lint --e{715}*/
755  SCIP_PRICERDATA* pricerdata;
756 
757  pricerdata = SCIPpricerGetData(pricer);
758  assert(pricerdata->randnumgen == NULL);
759 
760  SCIPfreeBlockMemoryNull(scip, &pricerdata);
761 
762  return SCIP_OKAY;
763 }
764 
765 /** initialization method of variable pricer (called after problem was transformed and pricer is active) */
766 static
767 SCIP_DECL_PRICERINIT(pricerInitRingpacking)
768 { /*lint --e{715}*/
769  SCIP_PRICERDATA* pricerdata = SCIPpricerGetData(pricer);
770  assert(pricerdata != NULL);
771  assert(pricerdata->randnumgen == NULL);
772 
773  /* create random number generator */
774  SCIP_CALL( SCIPcreateRandom(scip, &pricerdata->randnumgen, 0, TRUE) );
775 
776  return SCIP_OKAY;
777 }
778 
779 /** deinitialization method of variable pricer (called before transformed problem is freed and pricer is active) */
780 static
781 SCIP_DECL_PRICEREXIT(pricerExitRingpacking)
782 { /*lint --e{715}*/
783  SCIP_PRICERDATA* pricerdata = SCIPpricerGetData(pricer);
784  assert(pricerdata != NULL);
785  assert(pricerdata->randnumgen != NULL);
786 
787  SCIPfreeRandom(scip, &pricerdata->randnumgen);
788 
789  return SCIP_OKAY;
790 }
791 
792 /** reduced cost pricing method of variable pricer for feasible LPs */
793 static
794 SCIP_DECL_PRICERREDCOST(pricerRedcostRingpacking)
795 { /*lint --e{715}*/
796  SCIP_PRICERDATA* pricerdata;
797  SCIP_PROBDATA* probdata;
798  SCIP_CONS** conss;
799  SCIP_Real* lambdas;
800  SCIP_STATUS solstat;
801  SCIP_Real redcostslb;
802  SCIP_Real nlptimelimit;
803  SCIP_Real heurtimelimit;
804  SCIP_Real totaltilim;
805  SCIP_Bool success;
806  int t;
807 
808  pricerdata = SCIPpricerGetData(pricer);
809  assert(pricerdata != NULL);
810  probdata = SCIPgetProbData(scip);
811  assert(probdata != NULL);
812 
813  /* switch to price-and-price algorithm when dual bound has become invalid */
815 
816  /* only run pricer in the root node */
817  if( SCIPgetDepth(scip) > 0 )
818  {
820  *result = SCIP_SUCCESS;
821  return SCIP_OKAY;
822  }
823 
824  conss = SCIPprobdataGetPatternConss(probdata);
825  assert(conss != NULL);
826 
827  /* collect dual multipliers */
828  SCIP_CALL( SCIPallocBufferArray(scip, &lambdas, SCIPprobdataGetNTypes(probdata)) );
829  for( t = 0; t < SCIPprobdataGetNTypes(probdata); ++t )
830  {
831  assert(conss[t] != NULL);
832  assert( !strncmp(SCIPconshdlrGetName(SCIPconsGetHdlr(conss[t])), "linear", 6) );
833 
834  lambdas[t] = SCIPgetDualsolLinear(scip, conss[t]);
835  SCIPdebugMsg(scip, "lambda_%d = %g\n", t, lambdas[t]);
836  }
837 
838  /* collect working limits */
839  SCIP_CALL( SCIPgetRealParam(scip, "limits/time", &totaltilim) );
840 
841  /* solve pricing problem with heuristic */
842  heurtimelimit = MIN(pricerdata->heurtilim, totaltilim - SCIPgetSolvingTime(scip)); /*lint !e666*/
843  pricerdata->timeleft += SCIPgetSolvingTime(scip);
844  SCIP_CALL( solvePricingHeuristic(scip, probdata, pricerdata, lambdas, heurtimelimit, pricerdata->heuriterlim, &success) );
845  pricerdata->timeleft -= SCIPgetSolvingTime(scip);
846 
847  if( success )
848  {
849  *result = SCIP_SUCCESS;
850  }
851  /* solve pricing problem as MINLP if heuristic was not successful and dual bound is still valid */
852  else if ( !SCIPprobdataIsDualboundInvalid(probdata) )
853  {
854  nlptimelimit = MIN3(pricerdata->timeleft, pricerdata->nlptilim, totaltilim - SCIPgetSolvingTime(scip)); /*lint !e666*/
855  pricerdata->timeleft += SCIPgetSolvingTime(scip);
856  SCIP_CALL( solvePricingMINLP(scip, probdata, lambdas, nlptimelimit, pricerdata->nlpnodelim, &success, &solstat,
857  &redcostslb) );
858  pricerdata->timeleft -= SCIPgetSolvingTime(scip);
859  redcostslb += 1.0;
860  SCIPdebugMsg(scip, "result of pricing MINLP: addedvar=%u soltat=%d\n", success, solstat);
861 
862  /* check whether pricing problem could be solved to optimality */
863  if( SCIPisFeasGE(scip, redcostslb, 0.0) )
864  {
865  *lowerbound = SCIPgetLPObjval(scip);
866  SCIPinfoMessage(scip, NULL, "+++++++++++++ LP(master) = ceil(%g) = %g\n", *lowerbound, SCIPfeasCeil(scip, *lowerbound));
867  }
868  else
869  {
870  /* compute Farley's bound */
871  *lowerbound = SCIPgetLPObjval(scip) / (1.0 - redcostslb);
872  SCIPinfoMessage(scip, NULL, "+++++++++++++ Farley's bound = ceil(%g/%g) = %g\n", SCIPgetLPObjval(scip), 1.0 - redcostslb,
873  SCIPfeasCeil(scip, *lowerbound));
874  }
875  *lowerbound = SCIPfeasCeil(scip, *lowerbound);
876 
877  /* updates dual bound that is stored in the problem data */
878  SCIPprobdataUpdateDualbound(scip, probdata, *lowerbound);
879 
880  /* MINLP found an improving column or pricing problem could have been solved to optimality */
881  if( success || solstat == SCIP_STATUS_OPTIMAL || SCIPisFeasGE(scip, redcostslb, 0.0) )
882  *result = SCIP_SUCCESS;
883  }
884 
885  /* free memory */
886  SCIPfreeBufferArray(scip, &lambdas);
887 
888  return SCIP_OKAY;
889 }
890 
891 /** farkas pricing method of variable pricer for infeasible LPs */
892 static
893 SCIP_DECL_PRICERFARKAS(pricerFarkasRingpacking)
894 { /*lint --e{715}*/
895 
896  /* farkas pricing should not happen */
897  SCIPABORT();
898 
899  return SCIP_OKAY;
900 }
901 
902 /**@} */
903 
904 
905 /**@name Interface methods
906  *
907  * @{
908  */
909 
910 /** creates the ringpacking variable pricer and includes it in SCIP */
912  SCIP* scip /**< SCIP data structure */
913  )
914 {
915  SCIP_PRICERDATA* pricerdata;
916  SCIP_PRICER* pricer;
917 
918  /* create ringpacking variable pricer data */
919  SCIP_CALL( SCIPallocBlockMemory(scip, &pricerdata) );
920  BMSclearMemory(pricerdata);
921 
922  /* include variable pricer */
924  pricerRedcostRingpacking, pricerFarkasRingpacking, pricerdata) );
925 
926  SCIP_CALL( SCIPsetPricerFree(scip, pricer, pricerFreeRingpacking) );
927  SCIP_CALL( SCIPsetPricerInit(scip, pricer, pricerInitRingpacking) );
928  SCIP_CALL( SCIPsetPricerExit(scip, pricer, pricerExitRingpacking) );
929 
930  /* variable pricer parameters */
932  "ringpacking/pricing/nlptilim",
933  "time limit for each pricing NLP",
934  &pricerdata->nlptilim, FALSE, DEFAULT_PRICING_NLPTILIM, 0.0, SCIP_REAL_MAX, NULL, NULL) );
935 
937  "ringpacking/pricing/nlpnodelim",
938  "node limit for each pricing NLP",
939  &pricerdata->nlpnodelim, FALSE, DEFAULT_PRICING_NLPNODELIM, 0L, SCIP_LONGINT_MAX, NULL, NULL) );
940 
942  "ringpacking/pricing/heurtilim",
943  "time limit for each heuristic pricing",
944  &pricerdata->heurtilim, FALSE, DEFAULT_PRICING_HEURTILIM, 0.0, SCIP_REAL_MAX, NULL, NULL) );
945 
947  "ringpacking/pricing/heuriterlim",
948  "iteration limit for each heuristic pricing",
949  &pricerdata->heuriterlim, FALSE, DEFAULT_PRICING_HEURITERLIM, 0, INT_MAX, NULL, NULL) );
950 
952  "ringpacking/pricing/totaltilim",
953  "total time limit for all pricing NLPs and heuristic calls",
954  &pricerdata->timeleft, FALSE, DEFAULT_PRICING_TOTALTILIM, 0.0, SCIP_REAL_MAX, NULL, NULL) );
955 
956  return SCIP_OKAY;
957 }
958 
959 /** added problem specific data to pricer and activates pricer */
961  SCIP* scip /**< SCIP data structure */
962  )
963 {
964  SCIP_PRICER* pricer;
965 
966  assert(scip != NULL);
967 
968  pricer = SCIPfindPricer(scip, PRICER_NAME);
969  assert(pricer != NULL);
970 
971  /* activate pricer */
972  SCIP_CALL( SCIPactivatePricer(scip, pricer) );
973 
974  return SCIP_OKAY;
975 }
976 
977 /**@} */
#define SCIPfreeBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:97
static void computeScores(SCIP *scip, SCIP_Real *rexts, SCIP_RANDNUMGEN *randnumgen, int *elements, int nelements, SCIP_Real *lambdas, SCIP_Real *scores, int iter, int iterlim)
Definition: pricer_rpa.c:322
SCIP_RETCODE SCIPsetPricerInit(SCIP *scip, SCIP_PRICER *pricer, SCIP_DECL_PRICERINIT((*pricerinit)))
Definition: scip_pricer.c:214
#define DEFAULT_PRICING_TOTALTILIM
Definition: pricer_rpa.c:83
SCIP_Bool SCIPisStopped(SCIP *scip)
Definition: scip_general.c:687
#define SCIPallocBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:80
SCIP_RETCODE SCIPpatternAddElement(SCIP_PATTERN *pattern, int type, SCIP_Real x, SCIP_Real y)
Definition: pattern.c:173
void SCIPsetMessagehdlrQuiet(SCIP *scip, SCIP_Bool quiet)
Definition: scip_message.c:111
void SCIPfreeRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen)
SCIP_RETCODE SCIPcreateConsBasicQuadratic(SCIP *scip, SCIP_CONS **cons, const char *name, int nlinvars, SCIP_VAR **linvars, SCIP_Real *lincoefs, int nquadterms, SCIP_VAR **quadvars1, SCIP_VAR **quadvars2, SCIP_Real *quadcoefs, SCIP_Real lhs, SCIP_Real rhs)
SCIP_STATUS SCIPgetStatus(SCIP *scip)
Definition: scip_general.c:467
static int getNCircles(SCIP *scip, SCIP_Real rext, int demand, SCIP_Real width, SCIP_Real height, SCIP_Real lambda)
Definition: pricer_rpa.c:128
#define SCIP_MAXSTRLEN
Definition: def.h:279
SCIP_Real SCIPgetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var)
Definition: scip_sol.c:1353
Problem data for ringpacking problem.
SCIP_RETCODE SCIPaddCoefLinear(SCIP *scip, SCIP_CONS *cons, SCIP_VAR *var, SCIP_Real val)
SCIP_EXPORT void SCIPvarMarkDeletable(SCIP_VAR *var)
Definition: var.c:17250
#define PRICER_DESC
Definition: pricer_rpa.c:74
SCIP_RETCODE SCIPpricerRpaActivate(SCIP *scip)
Definition: pricer_rpa.c:960
SCIP_RETCODE SCIPcreateVarBasic(SCIP *scip, SCIP_VAR **var, const char *name, SCIP_Real lb, SCIP_Real ub, SCIP_Real obj, SCIP_VARTYPE vartype)
Definition: scip_var.c:185
SCIP_RETCODE SCIPsetPricerFree(SCIP *scip, SCIP_PRICER *pricer, SCIP_DECL_PRICERFREE((*pricerfree)))
Definition: scip_pricer.c:190
SCIP_Real SCIPgetSolvingTime(SCIP *scip)
Definition: scip_timing.c:360
int * SCIPprobdataGetDemands(SCIP_PROBDATA *probdata)
#define FALSE
Definition: def.h:73
void SCIPprobdataInvalidateDualbound(SCIP *scip, SCIP_PROBDATA *probdata)
static SCIP_RETCODE extractVariablesMINLP(SCIP *scip, SCIP_PROBDATA *probdata, SCIP *subscip, SCIP_VAR **x, SCIP_VAR **y, SCIP_VAR **w, int *types, int nvars, SCIP_Bool *success)
Definition: pricer_rpa.c:248
SCIP_Real SCIPrandomGetReal(SCIP_RANDNUMGEN *randnumgen, SCIP_Real minrandval, SCIP_Real maxrandval)
Definition: misc.c:9981
SCIP_PRICER * SCIPfindPricer(SCIP *scip, const char *name)
Definition: scip_pricer.c:302
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_Bool SCIPisFeasNegative(SCIP *scip, SCIP_Real val)
#define TRUE
Definition: def.h:72
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:54
SCIP_PROBDATA * SCIPgetProbData(SCIP *scip)
Definition: scip_prob.c:962
SCIP_RETCODE SCIPsolve(SCIP *scip)
Definition: scip_solve.c:2555
#define DEFAULT_PRICING_HEURITERLIM
Definition: pricer_rpa.c:82
SCIP_RETCODE SCIPsetHeuristics(SCIP *scip, SCIP_PARAMSETTING paramsetting, SCIP_Bool quiet)
Definition: scip_param.c:922
#define SCIP_LONGINT_MAX
Definition: def.h:149
#define SCIPfreeBufferArray(scip, ptr)
Definition: scip_mem.h:123
#define SCIPallocBlockMemory(scip, ptr)
Definition: scip_mem.h:78
#define SCIPdebugMsg
Definition: scip_message.h:69
SCIP_VAR ** x
Definition: circlepacking.c:54
SCIP_RETCODE SCIPaddPricedVar(SCIP *scip, SCIP_VAR *var, SCIP_Real score)
Definition: scip_prob.c:1731
SCIP_Bool SCIPprobdataIsDualboundInvalid(SCIP_PROBDATA *probdata)
SCIP_RETCODE SCIPcreateRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen, unsigned int initialseed, SCIP_Bool useglobalseed)
static SCIP_DECL_PRICERREDCOST(pricerRedcostRingpacking)
Definition: pricer_rpa.c:794
int SCIPgetNSols(SCIP *scip)
Definition: scip_sol.c:2206
SCIP_Bool SCIPisLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_RETCODE SCIPsetRealParam(SCIP *scip, const char *name, SCIP_Real value)
Definition: scip_param.c:619
SCIP_VAR * w
Definition: circlepacking.c:58
SCIP_RETCODE SCIPcreate(SCIP **scip)
Definition: scip_general.c:283
SCIP_EXPORT SCIP_RETCODE SCIPvarSetInitial(SCIP_VAR *var, SCIP_Bool initial)
Definition: var.c:17104
void SCIPpackCirclesGreedy(SCIP *scip, SCIP_Real *rexts, SCIP_Real *xs, SCIP_Real *ys, SCIP_Real rbounding, SCIP_Real width, SCIP_Real height, SCIP_Bool *ispacked, int *elements, int nelements, SCIP_PATTERNTYPE patterntype, int *npacked, int ncalls)
static SCIP_RETCODE solvePricingHeuristic(SCIP *scip, SCIP_PROBDATA *probdata, SCIP_PRICERDATA *pricerdata, SCIP_Real *lambdas, SCIP_Real timelim, int iterlim, SCIP_Bool *addedvar)
Definition: pricer_rpa.c:375
int SCIPprobdataGetNTypes(SCIP_PROBDATA *probdata)
SCIP_Real SCIPprobdataGetHeight(SCIP_PROBDATA *probdata)
SCIP_EXPORT void SCIPsortDownRealInt(SCIP_Real *realarray, int *intarray, int len)
#define PRICER_NAME
Definition: pricer_rpa.c:73
static SCIP_RETCODE addVariable(SCIP *scip, SCIP_PROBDATA *probdata, int *types, SCIP_Real *xs, SCIP_Real *ys, SCIP_Bool *ispacked, int nelems)
Definition: pricer_rpa.c:174
#define PRICER_PRIORITY
Definition: pricer_rpa.c:75
#define DEFAULT_PRICING_NLPNODELIM
Definition: pricer_rpa.c:80
pattern data for ringpacking problem
SCIP_CONS ** SCIPprobdataGetPatternConss(SCIP_PROBDATA *probdata)
#define NULL
Definition: lpi_spx1.cpp:155
#define SCIP_CALL(x)
Definition: def.h:370
#define SCIPfreeBlockMemoryNull(scip, ptr)
Definition: scip_mem.h:96
SCIP_Real SCIPfeasFloor(SCIP *scip, SCIP_Real val)
SCIP_Real SCIPgetLPObjval(SCIP *scip)
Definition: scip_lp.c:238
SCIP_RETCODE SCIPgetRealParam(SCIP *scip, const char *name, SCIP_Real *value)
Definition: scip_param.c:298
#define SCIPallocBufferArray(scip, ptr, num)
Definition: scip_mem.h:111
SCIP_Real SCIPinfinity(SCIP *scip)
int SCIPgetDepth(SCIP *scip)
Definition: scip_tree.c:638
#define SCIP_Bool
Definition: def.h:70
SCIP_RETCODE SCIPincludeDefaultPlugins(SCIP *scip)
SCIP_Real * SCIPprobdataGetRexts(SCIP_PROBDATA *probdata)
SCIP_Real SCIPgetTotalTime(SCIP *scip)
Definition: scip_timing.c:333
#define DEFAULT_PRICING_NLPTILIM
Definition: pricer_rpa.c:79
enum SCIP_Status SCIP_STATUS
Definition: type_stat.h:58
SCIP_EXPORT SCIP_RETCODE SCIPvarSetRemovable(SCIP_VAR *var, SCIP_Bool removable)
Definition: var.c:17120
SCIP_RETCODE SCIPincludePricerRpa(SCIP *scip)
Definition: pricer_rpa.c:911
static SCIP_DECL_PRICERINIT(pricerInitRingpacking)
Definition: pricer_rpa.c:767
#define MAX(x, y)
Definition: tclique_def.h:83
SCIP_RETCODE SCIPincludePricerBasic(SCIP *scip, SCIP_PRICER **pricerptr, const char *name, const char *desc, int priority, SCIP_Bool delay, SCIP_DECL_PRICERREDCOST((*pricerredcost)), SCIP_DECL_PRICERFARKAS((*pricerfarkas)), SCIP_PRICERDATA *pricerdata)
Definition: scip_pricer.c:118
#define PRICER_DELAY
Definition: pricer_rpa.c:76
SCIP_RETCODE SCIPpatternCreateRectangular(SCIP *scip, SCIP_PATTERN **pattern)
Definition: pattern.c:98
SCIP_Real SCIPprobdataGetWidth(SCIP_PROBDATA *probdata)
SCIP_RETCODE SCIPaddVar(SCIP *scip, SCIP_VAR *var)
Definition: scip_prob.c:1666
SCIP_Bool SCIPisFeasGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
#define BMSclearMemory(ptr)
Definition: memory.h:121
static SCIP_DECL_PRICEREXIT(pricerExitRingpacking)
Definition: pricer_rpa.c:781
SCIP_RETCODE SCIPaddRealParam(SCIP *scip, const char *name, const char *desc, SCIP_Real *valueptr, SCIP_Bool isadvanced, SCIP_Real defaultvalue, SCIP_Real minvalue, SCIP_Real maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:130
#define M_PI
Definition: pricer_rpa.c:88
#define SCIP_REAL_MAX
Definition: def.h:164
SCIP_Real SCIPfloor(SCIP *scip, SCIP_Real val)
SCIP_RETCODE SCIPsetLongintParam(SCIP *scip, const char *name, SCIP_Longint value)
Definition: scip_param.c:561
SCIP_PRICERDATA * SCIPpricerGetData(SCIP_PRICER *pricer)
Definition: pricer.c:501
struct SCIP_ProbData SCIP_PROBDATA
Definition: type_prob.h:44
void SCIPprobdataUpdateDualbound(SCIP *scip, SCIP_PROBDATA *probdata, SCIP_Real dualbound)
void SCIPpatternRelease(SCIP *scip, SCIP_PATTERN **pattern)
Definition: pattern.c:117
SCIP_Bool SCIPisFeasLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
int ncircles
Definition: circlepacking.c:47
int SCIPsnprintf(char *t, int len, const char *s,...)
Definition: misc.c:10604
SCIP_RETCODE SCIPaddIntParam(SCIP *scip, const char *name, const char *desc, int *valueptr, SCIP_Bool isadvanced, int defaultvalue, int minvalue, int maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:74
SCIP_RETCODE SCIPreleaseVar(SCIP *scip, SCIP_VAR **var)
Definition: scip_var.c:1245
#define SCIP_Real
Definition: def.h:163
SCIP_VAR ** y
Definition: circlepacking.c:55
SCIP_CONSHDLR * SCIPconsGetHdlr(SCIP_CONS *cons)
Definition: cons.c:8097
SCIP_Real SCIPgetDualsolLinear(SCIP *scip, SCIP_CONS *cons)
#define DEFAULT_PRICING_HEURTILIM
Definition: pricer_rpa.c:81
SCIP_Bool SCIPisGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_RETCODE SCIPsetPricerExit(SCIP *scip, SCIP_PRICER *pricer, SCIP_DECL_PRICEREXIT((*pricerexit)))
Definition: scip_pricer.c:238
void SCIPpatternSetPackableStatus(SCIP_PATTERN *pattern, SCIP_PACKABLE packable)
Definition: pattern.c:336
static SCIP_DECL_PRICERFREE(pricerFreeRingpacking)
Definition: pricer_rpa.c:753
#define SCIP_Longint
Definition: def.h:148
const char * SCIPconshdlrGetName(SCIP_CONSHDLR *conshdlr)
Definition: cons.c:4167
struct SCIP_PricerData SCIP_PRICERDATA
Definition: type_pricer.h:36
SCIP_RETCODE SCIPaddLongintParam(SCIP *scip, const char *name, const char *desc, SCIP_Longint *valueptr, SCIP_Bool isadvanced, SCIP_Longint defaultvalue, SCIP_Longint minvalue, SCIP_Longint maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:102
SCIP_RETCODE SCIPactivatePricer(SCIP *scip, SCIP_PRICER *pricer)
Definition: scip_pricer.c:375
SCIP_RETCODE SCIPaddCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip_prob.c:2764
SCIP_Real SCIPfeasCeil(SCIP *scip, SCIP_Real val)
static SCIP_RETCODE solvePricingMINLP(SCIP *scip, SCIP_PROBDATA *probdata, SCIP_Real *lambdas, SCIP_Real timelim, SCIP_Longint nodelim, SCIP_Bool *addedvar, SCIP_STATUS *solstat, SCIP_Real *dualbound)
Definition: pricer_rpa.c:502
SCIP_RETCODE SCIPfree(SCIP **scip)
Definition: scip_general.c:315
void SCIPinfoMessage(SCIP *scip, FILE *file, const char *formatstr,...)
Definition: scip_message.c:199
Ringpacking variable pricer.
SCIP_RETCODE SCIPreleaseCons(SCIP *scip, SCIP_CONS **cons)
Definition: scip_cons.c:1110
SCIP_RETCODE SCIPprobdataAddVar(SCIP *scip, SCIP_PROBDATA *probdata, SCIP_PATTERN *pattern, SCIP_VAR *var)
#define SCIPABORT()
Definition: def.h:342
SCIP_Real SCIPgetDualbound(SCIP *scip)
default SCIP plugins
static SCIP_DECL_PRICERFARKAS(pricerFarkasRingpacking)
Definition: pricer_rpa.c:893
SCIP_SOL * SCIPgetBestSol(SCIP *scip)
Definition: scip_sol.c:2305
SCIP_RETCODE SCIPsetIntParam(SCIP *scip, const char *name, int value)
Definition: scip_param.c:503
SCIP_RETCODE SCIPcreateProbBasic(SCIP *scip, const char *name)
Definition: scip_prob.c:170
SCIP callable library.
SCIP_Bool SCIPisFeasPositive(SCIP *scip, SCIP_Real val)
static SCIP_Real getDensityUb(int n)
Definition: pricer_rpa.c:118
SCIP_Real SCIPgetSolOrigObj(SCIP *scip, SCIP_SOL *sol)
Definition: scip_sol.c:1436