Scippy

SCIP

Solving Constraint Integer Programs

heur_shiftandpropagate.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-2014 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 email to scip@zib.de. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file heur_shiftandpropagate.c
17  * @brief shiftandpropagate primal heuristic
18  * @author Timo Berthold
19  * @author Gregor Hendel
20  */
21 
22 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
23 
24 #include <assert.h>
25 #include <string.h>
27 
28 #define HEUR_NAME "shiftandpropagate"
29 #define HEUR_DESC "Pre-root heuristic to expand an auxiliary branch-and-bound tree and apply propagation techniques"
30 #define HEUR_DISPCHAR 'T'
31 #define HEUR_PRIORITY 1000
32 #define HEUR_FREQ 0
33 #define HEUR_FREQOFS 0
34 #define HEUR_MAXDEPTH -1
35 #define HEUR_TIMING SCIP_HEURTIMING_BEFORENODE
36 #define HEUR_USESSUBSCIP FALSE /**< does the heuristic use a secondary SCIP instance? */
37 
38 #define DEFAULT_WEIGHT_INEQUALITY 1 /**< the heuristic row weight for inequalities */
39 #define DEFAULT_WEIGHT_EQUALITY 3 /**< the heuristic row weight for equations */
40 #define DEFAULT_RELAX TRUE /**< Should continuous variables be relaxed from the problem? */
41 #define DEFAULT_PROBING TRUE /**< Is propagation of solution values enabled? */
42 #define DEFAULT_ONLYWITHOUTSOL TRUE /**< Should heuristic only be executed if no primal solution was found, yet? */
43 #define DEFAULT_NPROPROUNDS 10 /**< The default number of propagation rounds for each propagation used */
44 #define DEFAULT_PROPBREAKER 65000 /**< fixed maximum number of propagations */
45 #define DEFAULT_CUTOFFBREAKER 15 /**< fixed maximum number of allowed cutoffs before the heuristic stops */
46 #define DEFAULT_RANDSEED 3141598 /**< the default random seed for random number generation */
47 #define DEFAULT_SORTKEY 'v' /**< the default key for variable sorting */
48 #define DEFAULT_SORTVARS TRUE /**< should variables be processed in sorted order? */
49 #define DEFAULT_COLLECTSTATS TRUE /**< should variable statistics be collected during probing? */
50 #define DEFAULT_STOPAFTERFEASIBLE TRUE /**< Should the heuristic stop calculating optimal shift values when no more rows are violated? */
51 #define DEFAULT_PREFERBINARIES TRUE /**< Should binary variables be shifted first? */
52 #define SORTKEYS "nrtuv"/**< options sorting key: (n)orms down, norms (u)p, (v)iolated rows decreasing,
53  * viola(t)ed rows increasing, or (r)andom */
54 #define DEFAULT_NOZEROFIXING FALSE /**< should variables with a zero shifting value be delayed instead of being fixed? */
55 #define DEFAULT_FIXBINLOCKS TRUE /**< should binary variables with no locks in one direction be fixed to that direction? */
56 #define DEFAULT_NORMALIZE TRUE /**< should coefficients and left/right hand sides be normalized by max row coeff? */
57 #define DEFAULT_UPDATEWEIGHTS FALSE /**< should row weight be increased every time the row is violated? */
58 #define DEFAULT_IMPLISCONTINUOUS TRUE /**< should implicit integer variables be treated as continuous variables? */
59 
60 #define EVENTHDLR_NAME "eventhdlrshiftandpropagate"
61 #define EVENTHDLR_DESC "event handler to catch bound changes"
62 
63 /*
64  * Data structures
65  */
66 
67 /** primal heuristic data */
68 struct SCIP_HeurData
69 {
70  SCIP_COL** lpcols; /**< stores lp columns with discrete variables before cont. variables */
71  int* rowweights; /**< row weight storage */
72  SCIP_Bool relax; /**< should continuous variables be relaxed from the problem */
73  SCIP_Bool probing; /**< should probing be executed? */
74  SCIP_Bool onlywithoutsol; /**< Should heuristic only be executed if no primal solution was found, yet? */
75  int nlpcols; /**< the number of lp columns */
76  int nproprounds; /**< The default number of propagation rounds for each propagation used */
77  int cutoffbreaker; /**< the number of cutoffs before heuristic execution is stopped, or -1 for no
78  * limit */
79  SCIP_EVENTHDLR* eventhdlr; /**< event handler to register and process variable bound changes */
80 
81  unsigned int randseed; /**< seed for random number generation */
82  char sortkey; /**< the key by which variables are sorted */
83  SCIP_Bool sortvars; /**< should variables be processed in sorted order? */
84  SCIP_Bool collectstats; /**< should variable statistics be collected during probing? */
85  SCIP_Bool stopafterfeasible; /**< Should the heuristic stop calculating optimal shift values when no
86  * more rows are violated? */
87  SCIP_Bool preferbinaries; /**< Should binary variables be shifted first? */
88  SCIP_Bool nozerofixing; /**< should variables with a zero shifting value be delayed instead of being fixed? */
89  SCIP_Bool fixbinlocks; /**< should binary variables with no locks in one direction be fixed to that direction? */
90  SCIP_Bool normalize; /**< should coefficients and left/right hand sides be normalized by max row coeff? */
91  SCIP_Bool updateweights; /**< should row weight be increased every time the row is violated? */
92  SCIP_Bool impliscontinuous; /**< should implicit integer variables be treated as continuous variables? */
94  SCIP_LPSOLSTAT lpsolstat; /**< the probing status after probing */
95  SCIP_Longint ntotaldomredsfound; /**< the total number of domain reductions during heuristic */
96  SCIP_Longint nlpiters; /**< number of LP iterations which the heuristic needed */
97  int nremainingviols; /**< the number of remaining violations */
98  int nprobings; /**< how many probings has the heuristic executed? */
99  int ncutoffs; /**< has the probing node been cutoff? */
100  int nredundantrows; /**< how many rows were redundant after relaxation? */
101  )
102 };
103 
104 /** status of a variable in heuristic transformation */
105 enum TransformStatus
106 {
107  TRANSFORMSTATUS_NONE = 0, /**< variable has not been transformed yet */
108  TRANSFORMSTATUS_LB = 1, /**< variable has been shifted by using lower bound (x-lb) */
109  TRANSFORMSTATUS_NEG = 2, /**< variable has been negated by using upper bound (ub-x) */
110  TRANSFORMSTATUS_FREE = 3 /**< variable does not have to be shifted */
111 };
112 typedef enum TransformStatus TRANSFORMSTATUS;
114 /** information about the matrix after its heuristic transformation */
115 struct ConstraintMatrix
116 {
117  SCIP_Real* rowmatvals; /**< matrix coefficients row by row */
118  int* rowmatind; /**< the indices of the corresponding variables */
119  int* rowmatbegin; /**< the starting indices of each row */
120  SCIP_Real* colmatvals; /**< matrix coefficients column by column */
121  int* colmatind; /**< the indices of the corresponding rows for each coefficient */
122  int* colmatbegin; /**< the starting indices of each column */
123  TRANSFORMSTATUS* transformstatus; /**< information about transform status of every discrete variable */
124  SCIP_Real* lhs; /**< left hand side vector after normalization */
125  SCIP_Real* rhs; /**< right hand side vector after normalization */
126  SCIP_Real* colnorms; /**< vector norms of all discrete problem variables after normalization */
127  SCIP_Real* upperbounds; /**< the upper bounds of every non-continuous variable after transformation*/
128  SCIP_Real* transformshiftvals; /**< values by which original discrete variable bounds were shifted */
129  int nnonzs; /**< number of nonzero column entries */
130  int nrows; /**< number of rows of matrix */
131  int ncols; /**< the number of columns in matrix (including continuous vars) */
132  int ndiscvars; /**< number of discrete problem variables */
133  SCIP_Bool normalized; /**< indicates if the matrix data has already been normalized */
134 };
135 typedef struct ConstraintMatrix CONSTRAINTMATRIX;
137 struct SCIP_EventhdlrData
138 {
139  CONSTRAINTMATRIX* matrix; /**< the constraint matrix of the heuristic */
140  SCIP_HEURDATA* heurdata; /**< heuristic data */
141  int* violatedrows; /**< all currently violated LP rows */
142  int* violatedrowpos; /**< position in violatedrows array for every row */
143  int* nviolatedrows; /**< pointer to the total number of currently violated rows */
144 };
145 
146 struct SCIP_EventData
147 {
148  int colpos; /**< column position of the event-related variable */
149 };
150 /*
151  * Local methods
152  */
153 
154 /** returns whether a given variable is counted as discrete, depending on the parameter impliscontinuous */
155 static
157  SCIP_VAR* var, /**< variable to check for discreteness */
158  SCIP_Bool impliscontinuous /**< should implicit integer variables be counted as continuous? */
159  )
160 {
161  return SCIPvarIsIntegral(var) && (SCIPvarGetType(var) != SCIP_VARTYPE_IMPLINT || !impliscontinuous);
162 }
163 
164 /** returns whether a given column is counted as discrete, depending on the parameter impliscontinuous */
165 static
167  SCIP_COL* col, /**< column to check for discreteness */
168  SCIP_Bool impliscontinuous /**< should implicit integer variables be counted as continuous? */
169  )
170 {
171  return SCIPcolIsIntegral(col) && (!impliscontinuous || SCIPvarGetType(SCIPcolGetVar(col)) != SCIP_VARTYPE_IMPLINT);
172 }
173 
174 /** returns nonzero values and corresponding columns of given row */
175 static
176 void getRowData(
177  CONSTRAINTMATRIX* matrix, /**< constraint matrix object */
178  int rowindex, /**< index of the desired row */
179  SCIP_Real** valpointer, /**< pointer to store the nonzero coefficients of the row */
180  SCIP_Real* lhs, /**< lhs of the row */
181  SCIP_Real* rhs, /**< rhs of the row */
182  int** indexpointer, /**< pointer to store column indices which belong to the nonzeros */
183  int* nrowvals /**< pointer to store number of nonzeros in the desired row (or NULL) */
184  )
185 {
186  int arrayposition;
187 
188  assert(matrix != NULL);
189  assert(0 <= rowindex && rowindex < matrix->nrows);
190 
191  arrayposition = matrix->rowmatbegin[rowindex];
192 
193  if ( nrowvals != NULL )
194  {
195  if( rowindex == matrix->nrows - 1 )
196  *nrowvals = matrix->nnonzs - arrayposition;
197  else
198  *nrowvals = matrix->rowmatbegin[rowindex + 1] - arrayposition; /*lint !e679*/
199  }
200 
201  if( valpointer != NULL )
202  *valpointer = &(matrix->rowmatvals[arrayposition]);
203  if( indexpointer != NULL )
204  *indexpointer = &(matrix->rowmatind[arrayposition]);
205 
206  if( lhs != NULL )
207  *lhs = matrix->lhs[rowindex];
208 
209  if( rhs != NULL )
210  *rhs = matrix->rhs[rowindex];
211 }
212 
213 /** returns nonzero values and corresponding rows of given column */
214 static
215 void getColumnData(
216  CONSTRAINTMATRIX* matrix, /**< constraint matrix object */
217  int colindex, /**< the index of the desired column */
218  SCIP_Real** valpointer, /**< pointer to store the nonzero coefficients of the column */
219  int** indexpointer, /**< pointer to store row indices which belong to the nonzeros */
220  int* ncolvals /**< pointer to store number of nonzeros in the desired column */
221  )
222 {
223  int arrayposition;
224 
225  assert(matrix != NULL);
226  assert(0 <= colindex && colindex < matrix->ncols);
227 
228  arrayposition = matrix->colmatbegin[colindex];
229 
230  if( ncolvals != NULL )
231  {
232  if( colindex == matrix->ncols - 1 )
233  *ncolvals = matrix->nnonzs - arrayposition;
234  else
235  *ncolvals = matrix->colmatbegin[colindex + 1] - arrayposition; /*lint !e679*/
236  }
237  if( valpointer != NULL )
238  *valpointer = &(matrix->colmatvals[arrayposition]);
239 
240  if( indexpointer != NULL )
241  *indexpointer = &(matrix->colmatind[arrayposition]);
242 }
243 
244 /** relaxes a continuous variable from all its rows, which has influence
245  * on both the left and right hand side of the constraint.
246  */
247 static
248 void relaxVar(
249  SCIP* scip, /**< current scip instance */
250  SCIP_VAR* var, /**< variable which is relaxed from the problem */
251  CONSTRAINTMATRIX* matrix, /**< constraint matrix object */
252  SCIP_Bool normalize /**< should coefficients and be normalized by rows maximum norms? */
253  )
254 {
255  SCIP_ROW** colrows;
256  SCIP_COL* varcol;
257  SCIP_Real* colvals;
258  SCIP_Real ub;
259  SCIP_Real lb;
260  int ncolvals;
261  int r;
262 
263  assert(var != NULL);
264  assert(SCIPvarGetStatus(var) == SCIP_VARSTATUS_COLUMN);
265 
266  varcol = SCIPvarGetCol(var);
267  assert(varcol != NULL);
268 
269  /* get nonzero values and corresponding rows of variable */
270  colvals = SCIPcolGetVals(varcol);
271  ncolvals = SCIPcolGetNLPNonz(varcol);
272  colrows = SCIPcolGetRows(varcol);
273 
274  ub = SCIPvarGetUbGlobal(var);
275  lb = SCIPvarGetLbGlobal(var);
276 
277  assert(colvals != NULL || ncolvals == 0);
278 
279  SCIPdebugMessage("Relaxing variable <%s> with lb <%g> and ub <%g>\n",
280  SCIPvarGetName(var), lb, ub);
281 
282  assert(matrix->normalized);
283  /* relax variable from all its constraints */
284  for( r = 0; r < ncolvals; ++r )
285  {
286  SCIP_ROW* colrow;
287  SCIP_Real lhs;
288  SCIP_Real rhs;
289  SCIP_Real lhsvarbound;
290  SCIP_Real rhsvarbound;
291  SCIP_Real rowabs;
292  SCIP_Real colval;
293  int rowindex;
294 
295  colrow = colrows[r];
296  rowindex = SCIProwGetLPPos(colrow);
297 
298  if( rowindex == -1 )
299  break;
300 
301  rowabs = SCIPgetRowMaxCoef(scip, colrow);
302  assert(colvals != NULL); /* to please flexelint */
303  colval = colvals[r];
304  if( normalize && SCIPisFeasGT(scip, rowabs, 0.0) )
305  colval /= rowabs;
306 
307  assert(0 <= rowindex && rowindex < matrix->nrows);
308  getRowData(matrix, rowindex, NULL, &lhs, &rhs, NULL, NULL);
309  /* variables bound influence the lhs and rhs of current row depending on the sign
310  * of the variables coefficient.
311  */
312  if( SCIPisFeasPositive(scip, colval) )
313  {
314  lhsvarbound = ub;
315  rhsvarbound = lb;
316  }
317  else if( SCIPisFeasNegative(scip, colval) )
318  {
319  lhsvarbound = lb;
320  rhsvarbound = ub;
321  }
322  else
323  continue;
324 
325  /* relax variable from the current row */
326  if( !SCIPisInfinity(scip, -matrix->lhs[rowindex]) && !SCIPisInfinity(scip, ABS(lhsvarbound)) )
327  matrix->lhs[rowindex] -= colval * lhsvarbound;
328  else
329  matrix->lhs[rowindex] = -SCIPinfinity(scip);
330 
331  if( !SCIPisInfinity(scip, matrix->rhs[rowindex]) && !SCIPisInfinity(scip, ABS(rhsvarbound)) )
332  matrix->rhs[rowindex] -= colval * rhsvarbound;
333  else
334  matrix->rhs[rowindex] = SCIPinfinity(scip);
335 
336  SCIPdebugMessage("Row <%s> changed:Coefficient <%g>, LHS <%g> --> <%g>, RHS <%g> --> <%g>\n",
337  SCIProwGetName(colrow), colval, lhs, matrix->lhs[rowindex], rhs, matrix->rhs[rowindex]);
338  }
339 }
340 
341 /** transforms bounds of a given variable s.t. its lower bound equals zero afterwards.
342  * If the variable already has lower bound zero, the variable is not transformed,
343  * if not, the variable's bounds are changed w.r.t. the smaller absolute value of its
344  * bounds in order to avoid numerical inaccuracies. If both lower and upper bound
345  * of the variable differ from infinity, there are two cases. If |lb| <= |ub|,
346  * the bounds are shifted by -lb, else a new variable ub - x replaces x.
347  * The transformation is memorized by the transform status of the variable s.t.
348  * retransformation is possible.
349  */
350 static
351 void transformVariable(
352  SCIP* scip, /**< current scip instance */
353  CONSTRAINTMATRIX* matrix, /**< constraint matrix object */
354  SCIP_HEURDATA* heurdata, /**< heuristic data */
355  int colpos /**< position of variable column in matrix */
356  )
357 {
358  SCIP_COL* col;
359  SCIP_VAR* var;
360  SCIP_Real lb;
361  SCIP_Real ub;
362 
363  SCIP_Bool negatecoeffs; /* do the row coefficients need to be negated? */
364  SCIP_Real deltashift; /* difference from previous transformation */
365 
366  assert(matrix != NULL);
367  assert(0 <= colpos && colpos < heurdata->nlpcols);
368  col = heurdata->lpcols[colpos];
369  assert(col != NULL);
370  assert(SCIPcolIsInLP(col));
371 
372  var = SCIPcolGetVar(col);
373  assert(var != NULL);
374  assert(SCIPvarIsIntegral(var));
375  lb = SCIPvarGetLbLocal(var);
376  ub = SCIPvarGetUbLocal(var);
377 
378  deltashift = 0.0;
379  negatecoeffs = FALSE;
380  /* if both lower and upper bound are -infinity and infinity, resp., this is reflected by a free transform status.
381  * If the lower bound is already zero, this is reflected by identity transform status. In both cases, none of the
382  * corresponding rows needs to be modified.
383  */
384  if( SCIPisInfinity(scip, -lb) && SCIPisInfinity(scip, ub) )
385  {
386  if( matrix->transformstatus[colpos] == TRANSFORMSTATUS_NEG )
387  negatecoeffs = TRUE;
388 
389  deltashift = matrix->transformshiftvals[colpos];
390  matrix->transformshiftvals[colpos] = 0.0;
391  matrix->transformstatus[colpos] = TRANSFORMSTATUS_FREE;
392  }
393  else if( SCIPisFeasLE(scip, ABS(lb), ABS(ub)) )
394  {
395  assert(!SCIPisInfinity(scip, lb));
396  matrix->transformstatus[colpos] = TRANSFORMSTATUS_LB;
397  deltashift = lb;
398  matrix->transformshiftvals[colpos] = lb;
399  }
400  else
401  {
402  assert(!SCIPisInfinity(scip, ub));
403  if( matrix->transformstatus[colpos] != TRANSFORMSTATUS_NEG )
404  negatecoeffs = TRUE;
405  matrix->transformstatus[colpos] = TRANSFORMSTATUS_NEG;
406  deltashift = ub;
407  matrix->transformshiftvals[colpos] = ub;
408  }
409 
410  /* determine the upper bound for this variable in heuristic transformation (lower bound is implicit; always 0) */
411  if( !SCIPisInfinity(scip, ub) && !SCIPisInfinity(scip, lb) )
412  matrix->upperbounds[colpos] = ub - lb;
413  else
414  matrix->upperbounds[colpos] = SCIPinfinity(scip);
415 
416  /* a real transformation is necessary. The variable x is either shifted by -lb or
417  * replaced by ub - x, depending on the smaller absolute of lb and ub.
418  */
419  if( !SCIPisFeasZero(scip, deltashift) || negatecoeffs )
420  {
421  SCIP_Real* vals;
422  int* rows;
423 
424  int nrows;
425  int i;
426 
427  assert(!SCIPisInfinity(scip, deltashift));
428 
429  /* get nonzero values and corresponding rows of column */
430  getColumnData(matrix, colpos, &vals, &rows, &nrows);
431  assert(nrows == 0 ||(vals != NULL && rows != NULL));
432 
433  /* go through rows and modify its lhs, rhs and the variable coefficient, if necessary */
434  for( i = 0; i < nrows; ++i )
435  {
436  assert(rows[i] >= 0);
437  assert(rows[i] < matrix->nrows);
438 
439  if( !SCIPisInfinity(scip, -(matrix->lhs[rows[i]])) )
440  matrix->lhs[rows[i]] -= (vals[i]) * deltashift;
441 
442  if( !SCIPisInfinity(scip, matrix->rhs[rows[i]]) )
443  matrix->rhs[rows[i]] -= (vals[i]) * deltashift;
444 
445  if( negatecoeffs )
446  (vals[i]) = -(vals[i]);
447 
448  assert(SCIPisFeasLE(scip, matrix->lhs[rows[i]], matrix->rhs[rows[i]]));
449  }
450  }
451  SCIPdebugMessage("Variable <%s> at colpos %d transformed. LB <%g> --> <%g>, UB <%g> --> <%g>\n",
452  SCIPvarGetName(var), colpos, lb, 0.0, ub, matrix->upperbounds[colpos]);
453 }
454 
455 /** initializes copy of the original coefficient matrix and applies heuristic specific adjustments: normalizing row
456  * vectors, transforming variable domains such that lower bound is zero, and relaxing continuous variables.
457  */
458 static
460  SCIP* scip, /**< current scip instance */
461  CONSTRAINTMATRIX* matrix, /**< constraint matrix object to be initialized */
462  SCIP_HEURDATA* heurdata, /**< heuristic data */
463  SCIP_Bool normalize, /**< should coefficients and be normalized by rows maximum norms? */
464  int* nmaxrows, /**< maximum number of rows a variable appears in */
465  SCIP_Bool relax, /**< should continuous variables be relaxed from the problem? */
466  SCIP_Bool* initialized, /**< was the initialization successful? */
467  SCIP_Bool* infeasible /**< is the problem infeasible? */
468  )
469 {
470  SCIP_ROW** lprows;
471  SCIP_COL** lpcols;
472  SCIP_Bool impliscontinuous;
473  int i;
474  int j;
475  int currentpointer;
476 
477  int nrows;
478  int ncols;
479 
480  assert(scip != NULL);
481  assert(matrix != NULL);
482  assert(initialized!= NULL);
483  assert(infeasible != NULL);
484  assert(nmaxrows != NULL);
485 
486  SCIPdebugMessage("entering Matrix Initialization method of SHIFTANDPROPAGATE heuristic!\n");
487 
488  /* get LP row data; column data is already initialized in heurdata */
489  SCIP_CALL( SCIPgetLPRowsData(scip, &lprows, &nrows) );
490  lpcols = heurdata->lpcols;
491  ncols = heurdata->nlpcols;
492 
493  matrix->nrows = nrows;
494  matrix->nnonzs = 0;
495  matrix->normalized = FALSE;
496  matrix->ndiscvars = 0;
497  *nmaxrows = 0;
498  impliscontinuous = heurdata->impliscontinuous;
499 
500  /* count the number of nonzeros of the LP constraint matrix */
501  for( j = 0; j < ncols; ++j )
502  {
503  assert(lpcols[j] != NULL);
504  assert(SCIPcolGetLPPos(lpcols[j]) >= 0);
505 
506  if( colIsDiscrete(lpcols[j], impliscontinuous) )
507  {
508  matrix->nnonzs += SCIPcolGetNLPNonz(lpcols[j]);
509  ++matrix->ndiscvars;
510  }
511  }
512 
513  matrix->ncols = matrix->ndiscvars;
514 
515  if( matrix->nnonzs == 0 )
516  {
517  SCIPdebugMessage("No matrix entries - Terminating initialization of matrix.\n");
518 
519  *initialized = FALSE;
520 
521  return SCIP_OKAY;
522  }
523 
524  /* allocate memory for the members of heuristic matrix */
525  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->rowmatvals, matrix->nnonzs) );
526  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->rowmatind, matrix->nnonzs) );
527  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->colmatvals, matrix->nnonzs) );
528  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->colmatind, matrix->nnonzs) );
529  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->rowmatbegin, nrows) );
530  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->colmatbegin, matrix->ncols) );
531  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->lhs, matrix->nrows) );
532  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->rhs, matrix->nrows) );
533  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->colnorms, matrix->ncols) );
534  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->transformstatus, matrix->ndiscvars) );
535  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->upperbounds, matrix->ndiscvars) );
536  SCIP_CALL( SCIPallocMemoryArray(scip, &matrix->transformshiftvals, matrix->ndiscvars) );
537 
538  /* set transform status of variables */
539  for( j = 0; j < matrix->ndiscvars; ++j )
540  matrix->transformstatus[j] = TRANSFORMSTATUS_NONE;
541 
542  currentpointer = 0;
543  *infeasible = FALSE;
544 
545  /* initialize the rows vector of the heuristic matrix together with its corresponding
546  * lhs, rhs.
547  */
548  for( i = 0; i < nrows; ++i )
549  {
550  SCIP_COL** cols;
551  SCIP_ROW* row;
552  SCIP_Real* rowvals;
553  SCIP_Real constant;
554  SCIP_Real maxval;
555  int nrowlpnonz;
556 
557  /* get LP row information */
558  row = lprows[i];
559  rowvals = SCIProwGetVals(row);
560  nrowlpnonz = SCIProwGetNLPNonz(row);
561  maxval = SCIPgetRowMaxCoef(scip, row);
562  cols = SCIProwGetCols(row);
563  constant = SCIProwGetConstant(row);
564 
565  SCIPdebugMessage(" %s : lhs=%g, rhs=%g, maxval=%g \n", SCIProwGetName(row), matrix->lhs[i], matrix->rhs[i], maxval);
566  SCIPdebug( SCIP_CALL( SCIPprintRow(scip, row, NULL) ) );
567  assert(!SCIPisInfinity(scip, constant));
568 
569  matrix->rowmatbegin[i] = currentpointer;
570 
571  /* modify the lhs and rhs w.r.t to the rows constant and normalize by 1-norm, i.e divide the lhs and rhs by the
572  * maximum absolute value of the row
573  */
574  if( !SCIPisInfinity(scip, -SCIProwGetLhs(row)) )
575  matrix->lhs[i] = (SCIProwGetLhs(row) - constant);
576  else
577  matrix->lhs[i] = -SCIPinfinity(scip);
578 
579  if( !SCIPisInfinity(scip, SCIProwGetRhs(row)) )
580  matrix->rhs[i] = (SCIProwGetRhs(row) - constant);
581  else
582  matrix->rhs[i] = SCIPinfinity(scip);
583 
584  /* make sure that maxval is larger than zero before normalization.
585  * Maxval may be zero if the constraint contains no variables but is modifiable, hence not redundant
586  */
587  if( normalize && !SCIPisFeasZero(scip, maxval) )
588  {
589  if( !SCIPisInfinity(scip, -matrix->lhs[i]) )
590  matrix->lhs[i] /= maxval;
591  if( !SCIPisInfinity(scip, matrix->rhs[i]) )
592  matrix->rhs[i] /= maxval;
593  }
594 
595 
596  /* in case of empty rows with a 0 < lhs <= 0.0 or 0.0 <= rhs < 0 we deduce the infeasibility of the problem */
597  if( nrowlpnonz == 0 && (SCIPisFeasPositive(scip, matrix->lhs[i]) || SCIPisFeasNegative(scip, matrix->rhs[i])) )
598  {
599  *infeasible = TRUE;
600  SCIPdebugMessage(" Matrix initialization stopped because of row infeasibility! \n");
601  break;
602  }
603 
604  /* row coefficients are normalized and copied to heuristic matrix */
605  for( j = 0; j < nrowlpnonz; ++j )
606  {
607  if( !colIsDiscrete(cols[j], impliscontinuous) )
608  continue;
609  assert(SCIPcolGetLPPos(cols[j]) >= 0);
610  assert(currentpointer < matrix->nnonzs);
611 
612  matrix->rowmatvals[currentpointer] = rowvals[j];
613  if( normalize && SCIPisFeasGT(scip, maxval, 0.0) )
614  matrix->rowmatvals[currentpointer] /= maxval;
615 
616  matrix->rowmatind[currentpointer] = SCIPcolGetLPPos(cols[j]);
617 
618  ++currentpointer;
619  }
620  }
621 
622  matrix->normalized = TRUE;
623 
624  if( *infeasible )
625  return SCIP_OKAY;
626 
627  assert(currentpointer == matrix->nnonzs);
628 
629  currentpointer = 0;
630 
631  /* copy the nonzero coefficient data column by column to heuristic matrix */
632  for( j = 0; j < matrix->ncols; ++j )
633  {
634  SCIP_COL* currentcol;
635  SCIP_ROW** rows;
636  SCIP_Real* colvals;
637  int ncolnonz;
638 
639 
640  assert(SCIPcolGetLPPos(lpcols[j]) >= 0);
641 
642  currentcol = lpcols[j];
643  assert(colIsDiscrete(currentcol, impliscontinuous));
644 
645  colvals = SCIPcolGetVals(currentcol);
646  rows = SCIPcolGetRows(currentcol);
647  ncolnonz = SCIPcolGetNLPNonz(currentcol);
648  matrix->colnorms[j] = ncolnonz;
649 
650  *nmaxrows = MAX(*nmaxrows, ncolnonz);
651 
652  /* loop over all rows with nonzero coefficients in the column, transform them and add them to the heuristic matrix */
653  matrix->colmatbegin[j] = currentpointer;
654 
655  for( i = 0; i < ncolnonz; ++i )
656  {
657  SCIP_Real maxval;
658 
659  assert(rows[i] != NULL);
660  assert(0 <= SCIProwGetLPPos(rows[i]));
661  assert(SCIProwGetLPPos(rows[i]) < nrows);
662  assert(currentpointer < matrix->nnonzs);
663 
664  /* rows are normalized by maximum norm */
665  maxval = SCIPgetRowMaxCoef(scip, rows[i]);
666 
667  assert(maxval > 0);
668 
669  matrix->colmatvals[currentpointer] = colvals[i];
670  if( normalize && SCIPisFeasGT(scip, maxval, 0.0) )
671  matrix->colmatvals[currentpointer] /= maxval;
672 
673  matrix->colmatind[currentpointer] = SCIProwGetLPPos(rows[i]);
674 
675  /* update the column norm */
676  matrix->colnorms[j] += ABS(matrix->colmatvals[currentpointer]);
677 
678  ++currentpointer;
679  }
680  }
681  assert(currentpointer == matrix->nnonzs);
682 
683  /* each variable is either transformed, if it supposed to be integral, or relaxed */
684  for( j = 0; j < (relax ? ncols : matrix->ndiscvars); ++j )
685  {
686  SCIP_COL* col;
687 
688  col = lpcols[j];
689  if( colIsDiscrete(col, impliscontinuous) )
690  {
691  matrix->transformshiftvals[j] = 0.0;
692  transformVariable(scip, matrix, heurdata, j);
693  }
694  else
695  {
696  SCIP_VAR* var;
697  var = SCIPcolGetVar(col);
698  assert(!varIsDiscrete(var, impliscontinuous));
699  relaxVar(scip, var, matrix, normalize);
700  }
701  }
702  *initialized = TRUE;
703 
704  SCIPdebugMessage("Matrix initialized for %d discrete variables with %d cols, %d rows and %d nonzero entries\n",
705  matrix->ndiscvars, matrix->ncols, matrix->nrows, matrix->nnonzs);
706  return SCIP_OKAY;
707 }
708 
709 /** frees all members of the heuristic matrix */
710 static
711 void freeMatrix(
712  SCIP* scip, /**< current SCIP instance */
713  CONSTRAINTMATRIX** matrix /**< constraint matrix object */
714  )
715 {
716  assert(scip != NULL);
717  assert(matrix != NULL);
718 
719  /* all fields are only allocated, if problem is not empty */
720  if( (*matrix)->nnonzs > 0 )
721  {
722  assert((*matrix) != NULL);
723  assert((*matrix)->rowmatbegin != NULL);
724  assert((*matrix)->rowmatvals != NULL);
725  assert((*matrix)->rowmatind != NULL);
726  assert((*matrix)->colmatbegin != NULL);
727  assert((*matrix)->colmatvals!= NULL);
728  assert((*matrix)->colmatind != NULL);
729  assert((*matrix)->lhs != NULL);
730  assert((*matrix)->rhs != NULL);
731  assert((*matrix)->transformstatus != NULL);
732  assert((*matrix)->transformshiftvals != NULL);
733 
734  /* free all fields */
735  SCIPfreeMemoryArray(scip, &((*matrix)->rowmatbegin));
736  SCIPfreeMemoryArray(scip, &((*matrix)->rowmatvals));
737  SCIPfreeMemoryArray(scip, &((*matrix)->rowmatind));
738  SCIPfreeMemoryArray(scip, &((*matrix)->colmatvals));
739  SCIPfreeMemoryArray(scip, &((*matrix)->colmatind));
740  SCIPfreeMemoryArray(scip, &((*matrix)->colmatbegin));
741  SCIPfreeMemoryArray(scip, &((*matrix)->lhs));
742  SCIPfreeMemoryArray(scip, &((*matrix)->rhs));
743  SCIPfreeMemoryArray(scip, &((*matrix)->colnorms));
744  SCIPfreeMemoryArray(scip, &((*matrix)->transformstatus));
745  SCIPfreeMemoryArray(scip, &((*matrix)->upperbounds));
746  SCIPfreeMemoryArray(scip, &((*matrix)->transformshiftvals));
747 
748  (*matrix)->nrows = 0;
749  (*matrix)->ncols = 0;
750  }
751 
752  /* free matrix */
753  SCIPfreeBuffer(scip, matrix);
754 }
755 
756 /** collects the necessary information about row violations for the zero-solution. That is,
757  * all solution values in heuristic transformation are zero.
758  */
759 static
760 void checkViolations(
761  SCIP* scip, /**< current scip instance */
762  CONSTRAINTMATRIX* matrix, /**< constraint matrix object */
763  int* violatedrows, /**< violated rows */
764  int* violatedrowpos, /**< row positions of violated rows */
765  int* nviolatedrows, /**< pointer to store the number of violated rows */
766  int* nredundantrows, /**< pointer to store the number of redundant rows */
767  int* violatedvarrows /**< reference to number of violated rows for every variable, or NULL */
768  )
769 {
770  SCIP_Real* rhs;
771  SCIP_Real* lhs;
772  int nrows;
773  int i;
774 
775  assert(matrix != NULL);
776  assert(violatedrows != NULL);
777  assert(violatedrowpos != NULL);
778  assert(nviolatedrows != NULL);
779 
780  /* get RHS, LHS and number of the problem rows */
781  rhs = matrix->rhs;
782  lhs = matrix->lhs;
783  nrows = matrix->nrows;
784 
785  SCIPdebugMessage("Entering violation check for %d rows! \n", nrows);
786  *nviolatedrows = 0;
787  if( nredundantrows != NULL )
788  *nredundantrows = 0;
789 
790  /* loop over rows and check if it is violated */
791  for( i = 0; i < nrows; ++i )
792  {
793  /* check, if zero solution violates this row */
794  if( SCIPisFeasLT(scip, rhs[i], 0.0) || SCIPisFeasGT(scip, lhs[i], 0.0) )
795  {
796  violatedrows[*nviolatedrows] = i;
797  (violatedrowpos)[i] = *nviolatedrows;
798  ++(*nviolatedrows);
799 
800  /* if needed, increase the counter for violated rows for every variable column of this row */
801  if( violatedvarrows != NULL )
802  {
803  int* rowcols;
804  int nrowcols;
805  int j;
806 
807  rowcols = NULL;
808  nrowcols = 0;
809  getRowData(matrix, i, NULL, NULL, NULL, &rowcols, &nrowcols);
810  assert(nrowcols == 0 || rowcols != NULL);
811 
812  for( j = 0; j < nrowcols; ++j )
813  {
814  assert(0 <= rowcols[j] && rowcols[j] < matrix->ndiscvars);
815  ++violatedvarrows[rowcols[j]];
816  }
817  }
818  }
819  else
820  violatedrowpos[i] = -1;
821 
822  assert((violatedrowpos[i] == -1 && SCIPisFeasGE(scip, rhs[i], 0.0) && SCIPisFeasGE(scip, -lhs[i], 0.0))
823  || (violatedrowpos[i] >= 0 &&(SCIPisFeasLT(scip, rhs[i], 0.0) || SCIPisFeasLT(scip, -lhs[i], 0.0))));
824 
825  if( SCIPisInfinity(scip, rhs[i]) && SCIPisInfinity(scip, -lhs[i]) && nredundantrows != NULL)
826  ++(*nredundantrows);
827  }
828 }
829 
830 /** retransforms solution values of variables according to their transformation status */
831 static
833  SCIP* scip, /**< current scip instance */
834  CONSTRAINTMATRIX* matrix, /**< constraint matrix object */
835  SCIP_VAR* var, /**< variable whose solution value has to be retransformed */
836  int varindex, /**< permutation of variable indices according to sorting */
837  SCIP_Real solvalue /**< solution value of the variable */
838  )
839 {
840  TRANSFORMSTATUS status;
841 
842  assert(matrix != NULL);
843  assert(var != NULL);
844 
845  status = matrix->transformstatus[varindex];
846  assert(status != TRANSFORMSTATUS_NONE);
847 
848  /* check if original variable has different bounds and transform solution value correspondingly */
849  if( status == TRANSFORMSTATUS_LB )
850  {
851  assert(!SCIPisInfinity(scip, -SCIPvarGetLbLocal(var)));
852 
853  return solvalue + matrix->transformshiftvals[varindex];
854  }
855  else if( status == TRANSFORMSTATUS_NEG )
856  {
857  assert(!SCIPisInfinity(scip, SCIPvarGetUbLocal(var)));
858  return matrix->transformshiftvals[varindex] - solvalue;
859  }
860  return solvalue;
861 }
862 
863 /** determines the best shifting value of a variable
864  * @todo if there is already an incumbent solution, try considering the objective cutoff as additional constraint */
865 static
867  SCIP* scip, /**< current scip instance */
868  CONSTRAINTMATRIX* matrix, /**< constraint matrix object */
869  int varindex, /**< index of variable which should be shifted */
870  int direction, /**< the direction for this variable */
871  int* rowweights, /**< weighting of rows for best shift calculation */
872  SCIP_Real* steps, /**< buffer array to store the individual steps for individual rows */
873  int* violationchange, /**< buffer array to store the individual change of feasibility of row */
874  SCIP_Real* beststep, /**< pointer to store optimal shifting step */
875  int* rowviolations /**< pointer to store new weighted sum of row violations, i.e, v - f */
876  )
877 {
878  SCIP_Real* vals;
879  int* rows;
880 
881  SCIP_Real slacksurplus;
882  SCIP_Real upperbound;
883 
884  int nrows;
885  int sum;
886  int i;
887 
888  SCIP_Bool allzero;
889 
890  assert(beststep != NULL);
891  assert(rowviolations != NULL);
892  assert(rowweights != NULL);
893  assert(steps != NULL);
894  assert(violationchange != NULL);
895  assert(direction == 1 || direction == -1);
896 
897  upperbound = matrix->upperbounds[varindex];
898 
899  /* get nonzero values and corresponding rows of variable */
900  getColumnData(matrix, varindex, &vals, &rows, &nrows);
901 
902  /* loop over rows and calculate, which is the minimum shift to make this row feasible
903  * or the minimum shift to violate this row
904  */
905  allzero = TRUE;
906  slacksurplus = 0.0;
907  for( i = 0; i < nrows; ++i )
908  {
909  SCIP_Real lhs;
910  SCIP_Real rhs;
911  SCIP_Real val;
912  int rowpos;
913  SCIP_Bool rowisviolated;
914  int rowweight;
915 
916  /* get the row data */
917  rowpos = rows[i];
918  assert(rowpos >= 0);
919  lhs = matrix->lhs[rowpos];
920  rhs = matrix->rhs[rowpos];
921  rowweight = rowweights[rowpos];
922  val = direction * vals[i];
923 
924  /* determine if current row is violated or not */
925  rowisviolated =(SCIPisFeasLT(scip, rhs, 0.0) || SCIPisFeasLT(scip, -lhs, 0.0));
926 
927  /* for a feasible row, determine the minimum integer value within the bounds of the variable by which it has to be
928  * shifted to make row infeasible.
929  */
930  if( !rowisviolated )
931  {
932  SCIP_Real maxfeasshift;
933 
934  maxfeasshift = SCIPinfinity(scip);
935 
936  /* feasibility can only be violated if the variable has a lock in the corresponding direction,
937  * i.e. a positive coefficient for a "<="-constraint, a negative coefficient for a ">="-constraint.
938  */
939  if( SCIPisFeasGT(scip, val, 0.0) && !SCIPisInfinity(scip, rhs) )
940  maxfeasshift = SCIPfeasFloor(scip, rhs/val);
941  else if( SCIPisFeasLT(scip, val, 0.0) && !SCIPisInfinity(scip, -lhs) )
942  maxfeasshift = SCIPfeasFloor(scip, lhs/val);
943 
944  /* if the variable has no lock in the current row, it can still help to increase the slack of this row;
945  * we measure slack increase for shifting by one
946  */
947  if( SCIPisFeasGT(scip, val, 0.0) && SCIPisInfinity(scip, rhs) )
948  slacksurplus += val;
949  if( SCIPisFeasLT(scip, val, 0.0) && SCIPisInfinity(scip, -lhs) )
950  slacksurplus -= val;
951 
952  /* check if the least violating shift lies within variable bounds and set corresponding array values */
953  if( SCIPisFeasLE(scip, maxfeasshift + 1.0, upperbound) )
954  {
955  steps[i] = maxfeasshift + 1.0;
956  violationchange[i] = rowweight;
957  allzero = FALSE;
958  }
959  else
960  {
961  steps[i] = upperbound;
962  violationchange[i] = 0;
963  }
964  }
965  /* for a violated row, determine the minimum integral value within the bounds of the variable by which it has to be
966  * shifted to make row feasible.
967  */
968  else
969  {
970  SCIP_Real minfeasshift;
971 
972  minfeasshift = SCIPinfinity(scip);
973 
974  /* if coefficient has the right sign to make row feasible, determine the minimum integer to shift variable
975  * to obtain feasibility
976  */
977  if( SCIPisFeasLT(scip, -lhs, 0.0) && SCIPisFeasGT(scip, val, 0.0) )
978  minfeasshift = SCIPfeasCeil(scip, lhs/val);
979  else if( SCIPisFeasLT(scip, rhs,0.0) && SCIPisFeasLT(scip, val, 0.0) )
980  minfeasshift = SCIPfeasCeil(scip, rhs/val);
981 
982  /* check if the minimum feasibility recovery shift lies within variable bounds and set corresponding array
983  * values
984  */
985  if( !SCIPisInfinity(scip, minfeasshift) && SCIPisFeasLE(scip, minfeasshift, upperbound) )
986  {
987  steps[i] = minfeasshift;
988  violationchange[i] = -rowweight;
989  allzero = FALSE;
990  }
991  else
992  {
993  steps[i] = upperbound;
994  violationchange[i] = 0;
995  }
996  }
997  }
998 
999  /* in case that the variable cannot affect the feasibility of any row, in particular it cannot violate
1000  * a single row, but we can add slack to already feasible rows, we will do this
1001  */
1002  if( allzero )
1003  {
1004  *beststep = SCIPisFeasGT(scip, slacksurplus, 0.0) ? direction * upperbound : 0.0;
1005  return SCIP_OKAY;
1006  }
1007 
1008  /* sorts rows by increasing value of steps */
1009  SCIPsortRealInt(steps, violationchange, nrows);
1010 
1011  *beststep = 0.0;
1012  *rowviolations = 0;
1013  sum = 0;
1014 
1015  /* best shifting step is calculated by summing up the violation changes for each relevant step and
1016  * taking the one which leads to the minimum sum. This sum measures the balance of feasibility recovering and
1017  * violating changes which will be obtained by shifting the variable by this step
1018  * note, the sums for smaller steps have to be taken into account for all bigger steps, i.e., the sums can be
1019  * computed iteratively
1020  */
1021  for( i = 0; i < nrows && !SCIPisInfinity(scip, steps[i]); ++i )
1022  {
1023  sum += violationchange[i];
1024 
1025  /* if we reached the last entry for the current step value, we have finished computing its sum and
1026  * update the step defining the minimum sum
1027  */
1028  if( (i == nrows-1 || steps[i+1] > steps[i]) && sum < *rowviolations ) /*lint !e679*/
1029  {
1030  *rowviolations = sum;
1031  *beststep = direction * steps[i];
1032  }
1033  }
1034  assert(*rowviolations <= 0);
1035  assert(!SCIPisInfinity(scip, *beststep));
1036 
1037  return SCIP_OKAY;
1038 }
1039 
1040 /** updates the information about a row whenever violation status changes */
1041 static
1042 void updateViolations(
1043  SCIP* scip, /**< current SCIP instance */
1044  CONSTRAINTMATRIX* matrix, /**< constraint matrix object */
1045  int rowindex, /**< index of the row */
1046  int* violatedrows, /**< contains all violated rows */
1047  int* violatedrowpos, /**< positions of rows in the violatedrows array */
1048  int* nviolatedrows, /**< pointer to update total number of violated rows */
1049  int* rowweights, /**< row weight storage */
1050  SCIP_Bool updateweights /**< should row weight be increased every time the row is violated? */
1051  )
1052 {
1053  assert(matrix != NULL);
1054  assert(violatedrows != NULL);
1055  assert(violatedrowpos != NULL);
1056  assert(nviolatedrows != NULL);
1057 
1058  /* row is now violated. Enqueue it in the set of violated rows. */
1059  if( SCIPisFeasLT(scip, -(matrix->lhs[rowindex]), 0.0) || SCIPisFeasLT(scip, matrix->rhs[rowindex], 0.0) )
1060  {
1061  assert(violatedrowpos[rowindex] == -1);
1062  assert(*nviolatedrows < matrix->nrows);
1063 
1064  violatedrows[*nviolatedrows] = rowindex;
1065  violatedrowpos[rowindex] = *nviolatedrows;
1066  ++(*nviolatedrows);
1067  if( updateweights )
1068  ++rowweights[rowindex];
1069  }
1070  /* row is now feasible. Remove it from the set of violated rows. */
1071  else
1072  {
1073  assert(violatedrowpos[rowindex] > -1);
1074 
1075  /* swap the row with last violated row */
1076  if( violatedrowpos[rowindex] != *nviolatedrows - 1 )
1077  {
1078  assert(*nviolatedrows - 1 >= 0);
1079  violatedrows[violatedrowpos[rowindex]] = violatedrows[*nviolatedrows - 1];
1080  violatedrowpos[violatedrows[*nviolatedrows - 1]] = violatedrowpos[rowindex];
1081  }
1082 
1083  /* unlink the row from its position in the array and decrease number of violated rows */
1084  violatedrowpos[rowindex] = -1;
1085  --(*nviolatedrows);
1086  }
1087 }
1088 
1089 /** updates transformation of a given variable by taking into account current local bounds. if the bounds have changed
1090  * since last update, updating the heuristic specific upper bound of the variable, its current transformed solution value
1091  * and all affected rows is necessary.
1092  */
1093 static
1095  SCIP* scip, /**< current scip */
1096  CONSTRAINTMATRIX* matrix, /**< constraint matrix object */
1097  SCIP_HEURDATA* heurdata, /**< heuristic data */
1098  int varindex, /**< index of variable in matrix */
1099  SCIP_Real* transformshiftval, /**< value, by which the variable has been shifted during transformation */
1100  SCIP_Real lb, /**< local lower bound of the variable */
1101  SCIP_Real ub, /**< local upper bound of the variable */
1102  int* violatedrows, /**< violated rows */
1103  int* violatedrowpos, /**< violated row positions */
1104  int* nviolatedrows /**< pointer to store number of violated rows */
1105  )
1106 {
1107  TRANSFORMSTATUS status;
1108  SCIP_Real deltashift;
1109 
1110  assert(scip != NULL);
1111  assert(matrix != NULL);
1112  assert(0 <= varindex && varindex < matrix->ndiscvars);
1113 
1114  /* deltashift is the difference between the old and new transformation value. */
1115  deltashift = 0.0;
1116  status = matrix->transformstatus[varindex];
1117 
1118  SCIPdebugMessage(" Variable <%d> [%g,%g], status %d(%g), ub %g \n", varindex, lb, ub, status,
1119  matrix->transformshiftvals[varindex], matrix->upperbounds[varindex]);
1120 
1121  /* depending on the variable status, deltashift is calculated differently. */
1122  if( status == TRANSFORMSTATUS_LB )
1123  {
1124  if( SCIPisInfinity(scip, -lb) )
1125  transformVariable(scip, matrix, heurdata, varindex);
1126  else
1127  {
1128  deltashift = lb - (*transformshiftval);
1129  *transformshiftval = lb;
1130  if( !SCIPisInfinity(scip, ub) )
1131  matrix->upperbounds[varindex] = ub - lb;
1132  else
1133  matrix->upperbounds[varindex] = SCIPinfinity(scip);
1134  }
1135  }
1136 
1137  if( status == TRANSFORMSTATUS_NEG )
1138  {
1139 
1140  if( SCIPisInfinity(scip, ub) )
1141  transformVariable(scip, matrix, heurdata, varindex);
1142  else
1143  {
1144  deltashift = (*transformshiftval) - ub;
1145  *transformshiftval = ub;
1146 
1147  if( !SCIPisInfinity(scip, -lb) )
1148  matrix->upperbounds[varindex] = ub - lb;
1149  }
1150  }
1151 
1152  if( status == TRANSFORMSTATUS_FREE )
1153  {
1154  /* in case of a free transform status, if one of the bounds has become finite, we want
1155  * to transform this variable to a variable with a lowerbound or a negated transform status */
1156  if( !SCIPisInfinity(scip, -lb) || !SCIPisInfinity(scip, ub) )
1157  {
1158  transformVariable(scip, matrix, heurdata, varindex);
1159 
1160  /* violations have to be rechecked for all rows
1161  * @todo : change this and only update violations of rows in which this variable
1162  * appears
1163  */
1164  checkViolations(scip, matrix, violatedrows, violatedrowpos, nviolatedrows, NULL, NULL);
1165 
1166  assert(matrix->transformstatus[varindex] == TRANSFORMSTATUS_LB || TRANSFORMSTATUS_NEG);
1167  assert(SCIPisFeasLE(scip, ABS(lb), ABS(ub)) || matrix->transformstatus[varindex] == TRANSFORMSTATUS_NEG);
1168  }
1169  }
1170 
1171  /* if the bound, by which the variable was shifted, has changed, deltashift is different from zero, which requires
1172  * an update of all affected rows
1173  */
1174  if( !SCIPisFeasZero(scip, deltashift) )
1175  {
1176  int i;
1177  int* rows;
1178  SCIP_Real* vals;
1179  int nrows;
1180 
1181  /* get nonzero values and corresponding rows of variable */
1182  getColumnData(matrix, varindex, &vals, &rows, &nrows);
1183 
1184  /* go through rows, update the rows w.r.t. the influence of the changed transformation of the variable */
1185  for( i = 0; i < nrows; ++i )
1186  {
1187  SCIP_Bool updaterow;
1188  SCIP_Bool leftviolation;
1189  SCIP_Bool rightviolation;
1190 
1191  SCIPdebugMessage(" update slacks of row<%d>: coefficient <%g>, %g <= 0 <= %g \n",
1192  rows[i], vals[i], matrix->lhs[rows[i]], matrix->rhs[rows[i]]);
1193 
1194  /* the row has to be updated if either lhs or rhs changes its sign. */
1195  leftviolation = SCIPisFeasLT(scip, -(matrix->lhs[rows[i]]), 0.0);
1196 
1197  if( !SCIPisInfinity(scip, -(matrix->lhs[rows[i]])) )
1198  matrix->lhs[rows[i]] -= (vals[i]) * deltashift;
1199 
1200  updaterow = leftviolation != SCIPisFeasLT(scip, -(matrix->lhs[rows[i]]), 0.0);
1201 
1202  rightviolation = SCIPisFeasLT(scip,(matrix->rhs[rows[i]]), 0.0);
1203  if( !SCIPisInfinity(scip, matrix->rhs[rows[i]]) )
1204  matrix->rhs[rows[i]] -= (vals[i]) * deltashift;
1205 
1206  updaterow = updaterow != (rightviolation != SCIPisFeasLT(scip,(matrix->rhs[rows[i]]), 0.0));
1207 
1208  /* update the row violation */
1209  if( updaterow )
1210  updateViolations(scip, matrix, rows[i], violatedrows, violatedrowpos, nviolatedrows, heurdata->rowweights, heurdata->updateweights);
1211 
1212  SCIPdebugMessage(" --> %g <= 0 <= %g %s\n",
1213  matrix->lhs[rows[i]], matrix->rhs[rows[i]], updaterow ? ": row violation updated " : "");
1214  }
1215  }
1216  SCIPdebugMessage(" Variable <%d> [%g,%g], status %d(%g), ub %g \n", varindex, lb, ub, status,
1217  matrix->transformshiftvals[varindex], matrix->upperbounds[varindex]);
1218 }
1219 
1220 /** comparison method for columns; binary < integer < implicit < continuous variables */
1221 static
1222 SCIP_DECL_SORTPTRCOMP(heurSortColsShiftandpropagate)
1224  SCIP_COL* col1;
1225  SCIP_COL* col2;
1226  SCIP_VAR* var1;
1227  SCIP_VAR* var2;
1228  SCIP_VARTYPE vartype1;
1229  SCIP_VARTYPE vartype2;
1230  int key1;
1231  int key2;
1232 
1233  col1 = (SCIP_COL*)elem1;
1234  col2 = (SCIP_COL*)elem2;
1235  var1 = SCIPcolGetVar(col1);
1236  var2 = SCIPcolGetVar(col2);
1237  assert(var1 != NULL && var2 != NULL);
1238 
1239  vartype1 = SCIPvarGetType(var1);
1240  vartype2 = SCIPvarGetType(var2);
1241 
1242  switch (vartype1)
1243  {
1244  case SCIP_VARTYPE_BINARY:
1245  key1 = 1;
1246  break;
1247  case SCIP_VARTYPE_INTEGER:
1248  key1 = 2;
1249  break;
1250  case SCIP_VARTYPE_IMPLINT:
1251  key1 = 3;
1252  break;
1254  key1 = 4;
1255  break;
1256  default:
1257  key1 = -1;
1258  SCIPerrorMessage("unknown variable type\n");
1259  SCIPABORT();
1260  break;
1261  }
1262  switch (vartype2)
1263  {
1264  case SCIP_VARTYPE_BINARY:
1265  key2 = 1;
1266  break;
1267  case SCIP_VARTYPE_INTEGER:
1268  key2 = 2;
1269  break;
1270  case SCIP_VARTYPE_IMPLINT:
1271  key2 = 3;
1272  break;
1274  key2 = 4;
1275  break;
1276  default:
1277  key2 = -1;
1278  SCIPerrorMessage("unknown variable type\n");
1279  SCIPABORT();
1280  break;
1281  }
1282  return key1 - key2;
1283 }
1284 
1285 /*
1286  * Callback methods of primal heuristic
1287  */
1288 
1289 /** deinitialization method of primal heuristic(called before transformed problem is freed) */
1290 static
1291 SCIP_DECL_HEUREXIT(heurExitShiftandpropagate)
1292 { /*lint --e{715}*/
1293  /* if statistic mode is enabled, statistics are printed to console */
1294  SCIPstatistic(
1295  SCIP_HEURDATA* heurdata;
1296 
1297  heurdata = SCIPheurGetData(heur);
1298 
1299  assert(heurdata != NULL);
1300 
1302  " DETAILS : %d violations left, %d probing status, %d redundant rows\n",
1303  heurdata->nremainingviols,
1304  heurdata->lpsolstat,
1305  heurdata->nredundantrows);
1307  " SHIFTANDPROPAGATE PROBING : %d probings, %lld domain reductions, ncutoffs: %d , LP iterations: %lld \n ",
1308  heurdata->nprobings,
1309  heurdata->ntotaldomredsfound,
1310  heurdata->ncutoffs,
1311  heurdata->nlpiters);
1312  );
1313 
1314  return SCIP_OKAY;
1315 }
1316 
1317 /** initialization method of primal heuristic(called after problem was transformed). We only need this method for
1318  * statistic mode of heuristic.
1319  */
1320 static
1321 SCIP_DECL_HEURINIT(heurInitShiftandpropagate)
1322 { /*lint --e{715}*/
1323 
1324  SCIP_HEURDATA* heurdata;
1325 
1326  heurdata = SCIPheurGetData(heur);
1327 
1328  assert(heurdata != NULL);
1329 
1330  heurdata->randseed = DEFAULT_RANDSEED;
1331 
1332  SCIPstatistic(
1333  heurdata->lpsolstat = SCIP_LPSOLSTAT_NOTSOLVED;
1334  heurdata->nremainingviols = 0;
1335  heurdata->nprobings = 0;
1336  heurdata->ntotaldomredsfound = 0;
1337  heurdata->ncutoffs = 0;
1338  heurdata->nlpiters = 0;
1339  heurdata->nredundantrows = 0;
1340  )
1341  return SCIP_OKAY;
1342 }
1343 
1344 /** destructor of primal heuristic to free user data(called when SCIP is exiting) */
1345 static
1346 SCIP_DECL_HEURFREE(heurFreeShiftandpropagate)
1347 { /*lint --e{715}*/
1348  SCIP_HEURDATA* heurdata;
1349  SCIP_EVENTHDLR* eventhdlr;
1350  SCIP_EVENTHDLRDATA* eventhdlrdata;
1351 
1352  heurdata = SCIPheurGetData(heur);
1353  eventhdlr = heurdata->eventhdlr;
1354  assert(eventhdlr != NULL);
1355  eventhdlrdata = SCIPeventhdlrGetData(eventhdlr);
1356 
1357  SCIPfreeMemory(scip, &eventhdlrdata);
1358 
1359  /* free heuristic data */
1360  if( heurdata != NULL )
1361  SCIPfreeMemory(scip, &heurdata);
1362 
1363  SCIPheurSetData(heur, NULL);
1364 
1365  return SCIP_OKAY;
1366 }
1367 
1368 
1369 /** copy method for primal heuristic plugins(called when SCIP copies plugins) */
1370 static
1371 SCIP_DECL_HEURCOPY(heurCopyShiftandpropagate)
1372 { /*lint --e{715}*/
1373  assert(scip != NULL);
1374  assert(heur != NULL);
1375  assert(strcmp(SCIPheurGetName(heur), HEUR_NAME) == 0);
1376 
1377  /* call inclusion method of primal heuristic */
1379 
1380  return SCIP_OKAY;
1381 }
1382 
1383 /** execution method of primal heuristic */
1384 static
1385 SCIP_DECL_HEUREXEC(heurExecShiftandpropagate)
1386 { /*lint --e{715}*/
1387  SCIP_HEURDATA* heurdata; /* heuristic data */
1388  SCIP_EVENTHDLR* eventhdlr; /* shiftandpropagate event handler */
1389  SCIP_EVENTHDLRDATA* eventhdlrdata; /* event handler data */
1390  SCIP_EVENTDATA** eventdatas; /* event data for every variable */
1391 
1392  CONSTRAINTMATRIX* matrix; /* constraint matrix object */
1393  SCIP_COL** lpcols; /* lp columns */
1394  SCIP_SOL* sol; /* solution pointer */
1395  SCIP_Real* colnorms; /* contains Euclidean norms of column vectors */
1396 
1397  SCIP_Real* steps; /* buffer arrays for best shift selection in main loop */
1398  int* violationchange;
1399 
1400  int* violatedrows; /* the violated rows */
1401  int* violatedrowpos; /* the array position of a violated row, or -1 */
1402  int* permutation; /* reflects the position of the variables after sorting */
1403  int* violatedvarrows; /* number of violated rows for each variable */
1404  int nlpcols; /* number of lp columns */
1405  int nviolatedrows; /* number of violated rows */
1406  int ndiscvars; /* number of non-continuous variables of the problem */
1407  int lastindexofsusp; /* last variable which has been swapped due to a cutoff */
1408  int nbinvars; /* number of binary variables */
1409  int nintvars; /* number of integer variables */
1410  int i;
1411  int r;
1412  int v;
1413  int c;
1414  int ncutoffs; /* counts the number of cutoffs for this execution */
1415  int nprobings; /* counts the number of probings */
1416  int nredundantrows; /* the number of redundant rows */
1417  int nlprows; /* the number LP rows */
1418  int nmaxrows; /* maximum number of LP rows of a variable */
1419 
1420  SCIP_Bool initialized; /* has the matrix been initialized? */
1421  SCIP_Bool cutoff; /* has current probing node been cutoff? */
1422  SCIP_Bool probing; /* should probing be applied or not? */
1423  SCIP_Bool infeasible; /* FALSE as long as currently infeasible rows have variables left */
1424  SCIP_Bool impliscontinuous;
1425 
1426  heurdata = SCIPheurGetData(heur);
1427  assert(heurdata != NULL);
1428 
1429  eventhdlr = heurdata->eventhdlr;
1430  assert(eventhdlr != NULL);
1431 
1432  eventhdlrdata = SCIPeventhdlrGetData(eventhdlr);
1433  assert(eventhdlrdata != NULL);
1434 
1435  *result = SCIP_DIDNOTRUN;
1436  SCIPdebugMessage("entering execution method of shift and propagate heuristic\n");
1437 
1438  /* heuristic is obsolete if there are only continuous variables */
1439  if( SCIPgetNVars(scip) - SCIPgetNContVars(scip) == 0 )
1440  return SCIP_OKAY;
1441 
1442  /* stop execution method if there is already a primarily feasible solution at hand */
1443  if( SCIPgetBestSol(scip) != NULL && heurdata->onlywithoutsol )
1444  return SCIP_OKAY;
1445 
1446  /* stop if there is no LP available */
1447  if ( ! SCIPhasCurrentNodeLP(scip) )
1448  return SCIP_OKAY;
1449 
1450  if( !SCIPisLPConstructed(scip) )
1451  {
1452  SCIP_Bool nodecutoff;
1453 
1454  nodecutoff = FALSE;
1455  /* @note this call can have the side effect that variables are created */
1456  SCIP_CALL( SCIPconstructLP(scip, &nodecutoff) );
1457  SCIP_CALL( SCIPflushLP(scip) );
1458  }
1459 
1460  SCIPstatistic( heurdata->nlpiters = SCIPgetNLPIterations(scip) );
1461 
1462  nlprows = SCIPgetNLPRows(scip);
1463 
1464  SCIP_CALL( SCIPgetLPColsData(scip, &lpcols, &nlpcols) );
1465  assert(nlpcols == 0 || lpcols != NULL);
1466 
1467  /* we need an LP */
1468  if( nlprows == 0 || nlpcols == 0 )
1469  return SCIP_OKAY;
1470 
1471 
1472  *result = SCIP_DIDNOTFIND;
1473  initialized = FALSE;
1474 
1475  /* allocate lp column array */
1476  SCIP_CALL( SCIPallocBufferArray(scip, &heurdata->lpcols, nlpcols) );
1477  heurdata->nlpcols = nlpcols;
1478 
1479  impliscontinuous = heurdata->impliscontinuous;
1480 
1481 #ifndef NDEBUG
1482  BMSclearMemoryArray(heurdata->lpcols, nlpcols);
1483 #endif
1484 
1485  BMScopyMemoryArray(heurdata->lpcols, lpcols, nlpcols);
1486 
1487  SCIPsortPtr((void**)heurdata->lpcols, heurSortColsShiftandpropagate, nlpcols);
1488 
1489  /* we have to collect the number of different variable types before we start probing since during probing variable
1490  * can be created (e.g., cons_xor.c)
1491  */
1492  ndiscvars = 0;
1493  nbinvars = 0;
1494  nintvars = 0;
1495  for( c = 0; c < nlpcols; ++c )
1496  {
1497  SCIP_COL* col;
1498  SCIP_VAR* colvar;
1499 
1500  col = heurdata->lpcols[c];
1501  assert(col != NULL);
1502  colvar = SCIPcolGetVar(col);
1503  assert(colvar != NULL);
1504 
1505  if( varIsDiscrete(colvar, impliscontinuous) )
1506  ++ndiscvars;
1507  if( SCIPvarGetType(colvar) == SCIP_VARTYPE_BINARY )
1508  ++nbinvars;
1509  else if( SCIPvarGetType(colvar) == SCIP_VARTYPE_INTEGER )
1510  ++nintvars;
1511  }
1512  assert(nbinvars + nintvars <= ndiscvars);
1513 
1514  /* start probing mode */
1515  SCIP_CALL( SCIPstartProbing(scip) );
1516 
1517  /* enables collection of variable statistics during probing */
1518  if( heurdata->collectstats )
1519  SCIPenableVarHistory(scip);
1520  else
1521  SCIPdisableVarHistory(scip);
1522 
1523  SCIP_CALL( SCIPnewProbingNode(scip) );
1524  ncutoffs = 0;
1525  nprobings = 0;
1526  nmaxrows = 0;
1527  infeasible = FALSE;
1528 
1529  /* initialize heuristic matrix and working solution */
1530  SCIP_CALL( SCIPallocBuffer(scip, &matrix) );
1531  SCIP_CALL( initMatrix(scip, matrix, heurdata, heurdata->normalize, &nmaxrows, heurdata->relax, &initialized, &infeasible) );
1532 
1533  /* could not initialize matrix */
1534  if( !initialized || infeasible )
1535  {
1536  SCIPdebugMessage(" MATRIX not initialized -> Execution of heuristic stopped! \n");
1537  goto TERMINATE;
1538  }
1539 
1540  /* the number of discrete LP column variables can be less than the actual number of variables, if, e.g., there
1541  * are nonlinearities in the problem. The heuristic execution can be terminated in that case.
1542  */
1543  if( matrix->ndiscvars < ndiscvars )
1544  {
1545  SCIPdebugMessage(" Not all discrete variables are in the current LP. Shiftandpropagate execution terminated\n");
1546  goto TERMINATE;
1547  }
1548 
1549  assert(nmaxrows > 0);
1550 
1551  eventhdlrdata->matrix = matrix;
1552  eventhdlrdata->heurdata = heurdata;
1553 
1554  SCIP_CALL( SCIPcreateSol(scip, &sol, heur) );
1555  SCIPsolSetHeur(sol, heur);
1556 
1557  /* allocate arrays for execution method */
1558  SCIP_CALL( SCIPallocBufferArray(scip, &permutation, ndiscvars) );
1559  SCIP_CALL( SCIPallocBufferArray(scip, &heurdata->rowweights, matrix->nrows) );
1560 
1561  /* allocate necessary memory for best shift search */
1562  SCIP_CALL( SCIPallocBufferArray(scip, &steps, nmaxrows) );
1563  SCIP_CALL( SCIPallocBufferArray(scip, &violationchange, nmaxrows) );
1564 
1565  /* allocate arrays to store information about infeasible rows */
1566  SCIP_CALL( SCIPallocBufferArray(scip, &violatedrows, matrix->nrows) );
1567  SCIP_CALL( SCIPallocBufferArray(scip, &violatedrowpos, matrix->nrows) );
1568 
1569  eventhdlrdata->violatedrows = violatedrows;
1570  eventhdlrdata->violatedrowpos = violatedrowpos;
1571  eventhdlrdata->nviolatedrows = &nviolatedrows;
1572 
1573 
1574 
1575  /* initialize arrays. Before sorting, permutation is the identity permutation */
1576  for( i = 0; i < ndiscvars; ++i )
1577  permutation[i] = i;
1578 
1579  /* initialize row weights */
1580  for( r = 0; r < matrix->nrows; ++r )
1581  {
1582  if( !SCIPisInfinity(scip, -(matrix->lhs[r])) && !SCIPisInfinity(scip, matrix->rhs[r]) )
1583  heurdata->rowweights[r] = DEFAULT_WEIGHT_EQUALITY;
1584  else
1585  heurdata->rowweights[r] = DEFAULT_WEIGHT_INEQUALITY;
1586 
1587  }
1588  colnorms = matrix->colnorms;
1589 
1590  assert(nbinvars >= 0);
1591  assert(nintvars >= 0);
1592 
1593  /* allocate memory for violatedvarrows array only if variable ordering relies on it */
1594  if( heurdata->sortvars && (heurdata->sortkey == 't' || heurdata->sortkey == 'v') )
1595  {
1596  SCIP_CALL( SCIPallocBufferArray(scip, &violatedvarrows, ndiscvars) );
1597  BMSclearMemoryArray(violatedvarrows, ndiscvars);
1598  }
1599  else
1600  violatedvarrows = NULL;
1601 
1602  /* check rows for infeasibility */
1603  nredundantrows = 0;
1604  checkViolations(scip, matrix, violatedrows, violatedrowpos, &nviolatedrows, &nredundantrows, violatedvarrows);
1605 
1606 
1607  /* sort variables w.r.t. the sorting key parameter. Sorting is indirect, all matrix column data
1608  * stays in place, but permutation array gives access to the sorted order of variables
1609  */
1610  if( heurdata->sortvars )
1611  {
1612  switch (heurdata->sortkey)
1613  {
1614  case 'n':
1615  /* variable ordering w.r.t. column norms nonincreasing */
1616  if( heurdata->preferbinaries )
1617  {
1618  if( nbinvars > 0 )
1619  SCIPsortDownRealInt(colnorms, permutation, nbinvars);
1620  if( nbinvars < ndiscvars )
1621  SCIPsortDownRealInt(&colnorms[nbinvars], &permutation[nbinvars], ndiscvars - nbinvars);
1622  }
1623  else
1624  {
1625  SCIPsortDownRealInt(colnorms, permutation, ndiscvars);
1626  }
1627  SCIPdebugMessage("Variables sorted down w.r.t their normalized columns!\n");
1628  break;
1629  case 'u':
1630  /* variable ordering w.r.t. column norms nondecreasing */
1631  if( heurdata->preferbinaries )
1632  {
1633  if( nbinvars > 0 )
1634  SCIPsortRealInt(colnorms, permutation, nbinvars);
1635  if( nbinvars < ndiscvars )
1636  SCIPsortRealInt(&colnorms[nbinvars], &permutation[nbinvars], ndiscvars - nbinvars);
1637  }
1638  else
1639  {
1640  SCIPsortRealInt(colnorms, permutation, ndiscvars);
1641  }
1642  SCIPdebugMessage("Variables sorted w.r.t their normalized columns!\n");
1643  break;
1644  case 'v':
1645  /* variable ordering w.r.t. nonincreasing number of violated rows */
1646  assert(violatedvarrows != NULL);
1647  if( heurdata->preferbinaries )
1648  {
1649  if( nbinvars > 0 )
1650  SCIPsortDownIntInt(violatedvarrows, permutation, nbinvars);
1651  if( nbinvars < ndiscvars )
1652  SCIPsortDownIntInt(&violatedvarrows[nbinvars], &permutation[nbinvars], ndiscvars - nbinvars);
1653  }
1654  else
1655  {
1656  SCIPsortDownIntInt(violatedvarrows, permutation, ndiscvars);
1657  }
1658 
1659  SCIPdebugMessage("Variables sorted down w.r.t their number of currently infeasible rows!\n");
1660  break;
1661  case 't':
1662  /* variable ordering w.r.t. nondecreasing number of violated rows */
1663  assert(violatedvarrows != NULL);
1664  if( heurdata->preferbinaries )
1665  {
1666  if( nbinvars > 0 )
1667  SCIPsortIntInt(violatedvarrows, permutation, nbinvars);
1668  if( nbinvars < ndiscvars )
1669  SCIPsortIntInt(&violatedvarrows[nbinvars], &permutation[nbinvars], ndiscvars - nbinvars);
1670  }
1671  else
1672  {
1673  SCIPsortIntInt(violatedvarrows, permutation, ndiscvars);
1674  }
1675 
1676  SCIPdebugMessage("Variables sorted (upwards) w.r.t their number of currently infeasible rows!\n");
1677  break;
1678  case 'r':
1679  /* random sorting */
1680  if( heurdata->preferbinaries )
1681  {
1682  if( nbinvars > 0 )
1683  SCIPpermuteIntArray(permutation, 0, nbinvars - 1, &heurdata->randseed);
1684  if( nbinvars < ndiscvars )
1685  SCIPpermuteIntArray(&permutation[nbinvars], nbinvars - 1, ndiscvars - nbinvars - 1, &heurdata->randseed);
1686  }
1687  else
1688  {
1689  SCIPpermuteIntArray(permutation, 0, ndiscvars - 1, &heurdata->randseed);
1690  }
1691  SCIPdebugMessage("Variables permuted randomly!\n");
1692  break;
1693  default:
1694  SCIPdebugMessage("No variable permutation applied\n");
1695  break;
1696  }
1697  }
1698 
1699  SCIP_CALL( SCIPallocBufferArray(scip, &eventdatas, matrix->ndiscvars) );
1700  BMSclearMemoryArray(eventdatas, matrix->ndiscvars);
1701 
1702  /* initialize variable events to catch bound changes during propagation */
1703  for( c = 0; c < matrix->ndiscvars; ++c )
1704  {
1705  SCIP_VAR* var;
1706 
1707  var = SCIPcolGetVar(heurdata->lpcols[c]);
1708  assert(var != NULL);
1709  assert(SCIPvarIsIntegral(var));
1710  assert(eventdatas[c] == NULL);
1711 
1712  SCIP_CALL( SCIPallocBuffer(scip, &(eventdatas[c])) ); /*lint !e866*/
1713 
1714  eventdatas[c]->colpos = c;
1715 
1716  SCIP_CALL( SCIPcatchVarEvent(scip, var, SCIP_EVENTTYPE_BOUNDCHANGED, eventhdlr, eventdatas[c], NULL) );
1717  }
1718 
1719  cutoff = FALSE;
1720  lastindexofsusp = -1;
1721  probing = heurdata->probing;
1722  infeasible = FALSE;
1723 
1724  SCIPdebugMessage("SHIFT_AND_PROPAGATE heuristic starts main loop with %d violations and %d remaining variables!\n",
1725  nviolatedrows, ndiscvars);
1726 
1727  assert(matrix->ndiscvars == ndiscvars);
1728 
1729  /* loop over variables, shift them according to shifting criteria and try to reduce the global infeasibility */
1730  for( c = 0; c < ndiscvars; ++c )
1731  {
1732  SCIP_VAR* var;
1733  SCIP_Longint ndomredsfound;
1734  SCIP_Real optimalshiftvalue;
1735  SCIP_Real origsolval;
1736  SCIP_Real lb;
1737  SCIP_Real ub;
1738  TRANSFORMSTATUS status;
1739  int nviolations;
1740  int permutedvarindex;
1741  SCIP_Bool marksuspicious;
1742 
1743  permutedvarindex = permutation[c];
1744  optimalshiftvalue = 0.0;
1745  nviolations = 0;
1746  var = SCIPcolGetVar(heurdata->lpcols[permutedvarindex]);
1747  lb = SCIPvarGetLbLocal(var);
1748  ub = SCIPvarGetUbLocal(var);
1749  assert(SCIPcolGetLPPos(SCIPvarGetCol(var)) >= 0);
1750  assert(SCIPvarIsIntegral(var));
1751 
1752  /* check whether we hit some limit, e.g. the time limit, in between
1753  * since the check itself consumes some time, we only do it every tenth iteration
1754  */
1755  if( c % 10 == 0 && SCIPisStopped(scip) )
1756  goto TERMINATE2;
1757 
1758  /* if propagation is enabled, check if propagation has changed the variables bounds
1759  * and update the transformed upper bound correspondingly
1760  */
1761  if( heurdata->probing )
1762  updateTransformation(scip, matrix, heurdata, permutedvarindex, &(matrix->transformshiftvals[permutedvarindex]),
1763  lb, ub, violatedrows, violatedrowpos, &nviolatedrows);
1764 
1765  status = matrix->transformstatus[permutedvarindex];
1766 
1767  SCIPdebugMessage("Variable %s with local bounds [%g,%g], status <%d>, matrix bound <%g>\n",
1768  SCIPvarGetName(var), lb, ub, status, matrix->upperbounds[permutedvarindex]);
1769 
1770  /* ignore variable if propagation fixed it(lb and ub will be zero) */
1771  if( SCIPisFeasZero(scip, matrix->upperbounds[permutedvarindex]) )
1772  {
1773  assert(!SCIPisInfinity(scip, ub));
1774  assert(SCIPisFeasEQ(scip, lb, ub));
1775 
1776  SCIP_CALL( SCIPsetSolVal(scip, sol, var, ub) );
1777 
1778  continue;
1779  }
1780 
1781  marksuspicious = FALSE;
1782 
1783  /* check whether the variable is binary and has no locks in one direction, so that we want to fix it to the
1784  * respective bound (only enabled by parameter)
1785  */
1786  if( heurdata->fixbinlocks && SCIPvarIsBinary(var) && (SCIPvarGetNLocksUp(var) == 0 || SCIPvarGetNLocksDown(var) == 0) )
1787  {
1788  if( SCIPvarGetNLocksUp(var) == 0 )
1789  origsolval = SCIPvarGetUbLocal(var);
1790  else
1791  {
1792  assert(SCIPvarGetNLocksDown(var) == 0);
1793  origsolval = SCIPvarGetLbLocal(var);
1794  }
1795  }
1796  else
1797  {
1798  /* only apply the computationally expensive best shift selection, if there is a violated row left */
1799  if( !heurdata->stopafterfeasible || nviolatedrows > 0 )
1800  {
1801  /* compute optimal shift value for variable */
1802  SCIP_CALL( getOptimalShiftingValue(scip, matrix, permutedvarindex, 1, heurdata->rowweights, steps, violationchange,
1803  &optimalshiftvalue, &nviolations) );
1804  assert(SCIPisFeasGE(scip, optimalshiftvalue, 0.0));
1805 
1806  /* Variables with FREE transform have to be dealt with twice */
1807  if( matrix->transformstatus[permutedvarindex] == TRANSFORMSTATUS_FREE )
1808  {
1809  SCIP_Real downshiftvalue;
1810  int ndownviolations;
1811 
1812  downshiftvalue = 0.0;
1813  ndownviolations = 0;
1814  SCIP_CALL( getOptimalShiftingValue(scip, matrix, permutedvarindex, -1, heurdata->rowweights, steps, violationchange,
1815  &downshiftvalue, &ndownviolations) );
1816 
1817  assert(SCIPisLE(scip, downshiftvalue, 0.0));
1818 
1819  /* compare to positive direction and select the direction which makes more rows feasible */
1820  if( ndownviolations < nviolations )
1821  {
1822  optimalshiftvalue = downshiftvalue;
1823  }
1824  }
1825  }
1826  else
1827  optimalshiftvalue = 0.0;
1828 
1829  /* if zero optimal shift values are forbidden by the user parameter, delay the variable by marking it suspicious */
1830  if( heurdata->nozerofixing && nviolations > 0 && SCIPisFeasZero(scip, optimalshiftvalue) )
1831  marksuspicious = TRUE;
1832 
1833  /* retransform the solution value from the heuristic transformation space */
1834  assert(varIsDiscrete(var, impliscontinuous));
1835  origsolval = retransformVariable(scip, matrix, var, permutedvarindex, optimalshiftvalue);
1836  }
1837  assert(SCIPisFeasGE(scip, origsolval, lb) && SCIPisFeasLE(scip, origsolval, ub));
1838 
1839  /* check if propagation should still be performed */
1840  if( nprobings > DEFAULT_PROPBREAKER )
1841  probing = FALSE;
1842 
1843  /* if propagation is enabled, fix the variable to the new solution value and propagate the fixation
1844  * (to fix other variables and to find out early whether solution is already infeasible)
1845  */
1846  if( !marksuspicious && probing )
1847  {
1848  SCIP_CALL( SCIPnewProbingNode(scip) );
1849  SCIP_CALL( SCIPfixVarProbing(scip, var, origsolval) );
1850  ndomredsfound = 0;
1851 
1852  SCIPdebugMessage(" Shift %g(%g originally) is optimal, propagate solution\n", optimalshiftvalue, origsolval);
1853  SCIP_CALL( SCIPpropagateProbing(scip, heurdata->nproprounds, &cutoff, &ndomredsfound) );
1854 
1855  ++nprobings;
1856  SCIPstatistic( heurdata->ntotaldomredsfound += ndomredsfound );
1857  SCIPdebugMessage("Propagation finished! <%lld> domain reductions %s, <%d> probing depth\n", ndomredsfound, cutoff ? "CUTOFF" : "",
1858  SCIPgetProbingDepth(scip));
1859  }
1860  assert(!cutoff || probing);
1861 
1862  /* propagation led to an empty domain, hence we backtrack and postpone the variable */
1863  if( cutoff )
1864  {
1865  assert(probing);
1866 
1867  ++ncutoffs;
1868 
1869  /* only continue heuristic if number of cutoffs occured so far is reasonably small */
1870  if( heurdata->cutoffbreaker >= 0 && ncutoffs >= heurdata->cutoffbreaker )
1871  break;
1872 
1873  cutoff = FALSE;
1874 
1875  /* backtrack to the parent of the current node */
1876  assert(SCIPgetProbingDepth(scip) >= 1);
1878  marksuspicious = TRUE;
1879 
1880  /* if the variable were to be set to one of its bounds, repropagate by tightening this bound by 1.0
1881  * into the direction of the other bound, if possible */
1882  if( SCIPisFeasEQ(scip, SCIPvarGetLbLocal(var), origsolval) )
1883  {
1884  assert(SCIPisFeasGE(scip, SCIPvarGetUbLocal(var), origsolval + 1.0));
1885 
1886  ndomredsfound = 0;
1887  SCIP_CALL( SCIPnewProbingNode(scip) );
1888  SCIP_CALL( SCIPchgVarLbProbing(scip, var, origsolval + 1.0) );
1889  SCIP_CALL( SCIPpropagateProbing(scip, heurdata->nproprounds, &cutoff, &ndomredsfound) );
1890 
1891  SCIPstatistic( heurdata->ntotaldomredsfound += ndomredsfound );
1892 
1893  /* if the tightened bound again leads to a cutoff, both subproblems are proven infeasible and the heuristic
1894  * can be stopped */
1895  if( cutoff )
1896  break;
1897  }
1898  else if( SCIPisFeasEQ(scip, SCIPvarGetUbLocal(var), origsolval) )
1899  {
1900  assert(SCIPisFeasLE(scip, SCIPvarGetLbLocal(var), origsolval - 1.0));
1901 
1902  ndomredsfound = 0;
1903 
1904  SCIP_CALL( SCIPnewProbingNode(scip) );
1905  SCIP_CALL( SCIPchgVarUbProbing(scip, var, origsolval - 1.0) );
1906  SCIP_CALL( SCIPpropagateProbing(scip, heurdata->nproprounds, &cutoff, &ndomredsfound) );
1907 
1908  SCIPstatistic( heurdata->ntotaldomredsfound += ndomredsfound );
1909 
1910  /* if the tightened bound again leads to a cutoff, both subproblems are proven infeasible and the heuristic
1911  * can be stopped */
1912  if( cutoff )
1913  break;
1914  }
1915  }
1916 
1917  if( marksuspicious )
1918  {
1919  /* mark the variable as suspicious */
1920  assert(permutedvarindex == permutation[c]);
1921 
1922  ++lastindexofsusp;
1923  assert(lastindexofsusp >= 0 && lastindexofsusp <= c);
1924 
1925  permutation[c] = permutation[lastindexofsusp];
1926  permutation[lastindexofsusp] = permutedvarindex;
1927 
1928  SCIPdebugMessage(" Suspicious variable! Postponed from pos <%d> to position <%d>\n", c, lastindexofsusp);
1929  }
1930  else
1931  {
1932  SCIPdebugMessage("Variable <%d><%s> successfully shifted by value <%g>!\n", permutedvarindex,
1933  SCIPvarGetName(var), optimalshiftvalue);
1934 
1935  /* update solution */
1936  SCIP_CALL( SCIPsetSolVal(scip, sol, var, origsolval) );
1937 
1938  /* only to ensure that some assertions can be made later on */
1939  if( !probing )
1940  {
1941  SCIP_CALL( SCIPfixVarProbing(scip, var, origsolval) );
1942  }
1943  }
1944  }
1945  SCIPdebugMessage("Heuristic finished with %d remaining violations and %d remaining variables!\n",
1946  nviolatedrows, lastindexofsusp + 1);
1947 
1948  /* if constructed solution might be feasible, go through the queue of suspicious variables and set the solution
1949  * values
1950  */
1951  if( nviolatedrows == 0 && !cutoff )
1952  {
1953  SCIP_Bool stored;
1954 
1955  for( v = 0; v <= lastindexofsusp; ++v )
1956  {
1957  SCIP_VAR* var;
1958  SCIP_Real origsolval;
1959  int permutedvarindex;
1960 
1961  /* get the column position of the variable */
1962  permutedvarindex = permutation[v];
1963  var = SCIPcolGetVar(heurdata->lpcols[permutedvarindex]);
1964  assert(varIsDiscrete(var, impliscontinuous));
1965 
1966  /* update the transformation of the variable, since the bound might have changed after the last update. */
1967  if( heurdata->probing )
1968  updateTransformation(scip, matrix, heurdata, permutedvarindex, &(matrix->transformshiftvals[permutedvarindex]),
1969  SCIPvarGetLbLocal(var), SCIPvarGetUbLocal(var), violatedrows, violatedrowpos, &nviolatedrows);
1970 
1971  /* retransform the solution value from the heuristic transformed space, set the solution value accordingly */
1972  assert(varIsDiscrete(var, impliscontinuous));
1973  origsolval = retransformVariable(scip, matrix, var, permutedvarindex, 0.0);
1974  assert(SCIPisFeasGE(scip, origsolval, SCIPvarGetLbLocal(var))
1975  && SCIPisFeasLE(scip, origsolval, SCIPvarGetUbLocal(var)));
1976  SCIP_CALL( SCIPsetSolVal(scip, sol, var, origsolval) );
1977  SCIP_CALL( SCIPfixVarProbing(scip, var, origsolval) ); /* only to ensure that some assertions can be made later */
1978 
1979  SCIPdebugMessage(" Remaining variable <%s> set to <%g>; %d Violations\n", SCIPvarGetName(var), origsolval,
1980  nviolatedrows);
1981  }
1982  /* Fixing of remaining variables led to infeasibility */
1983  if( nviolatedrows > 0 )
1984  goto TERMINATE2;
1985 
1986  stored = TRUE;
1987  /* if the constructed solution might still be extendable to a feasible solution, try this by
1988  * solving the remaining LP
1989  */
1990  if( nlpcols != matrix->ndiscvars )
1991  {
1992  /* case that remaining LP has to be solved */
1993  SCIP_Bool lperror;
1994 
1995 #ifndef NDEBUG
1996  {
1997  SCIP_VAR** vars;
1998 
1999  vars = SCIPgetVars(scip);
2000  assert(vars != NULL);
2001  /* ensure that all discrete variables in the remaining LP are fixed */
2002  for( v = 0; v < ndiscvars; ++v )
2003  {
2004  if( SCIPvarIsInLP(vars[v]) )
2005  assert(SCIPisFeasEQ(scip, SCIPvarGetLbLocal(vars[v]), SCIPvarGetUbLocal(vars[v])));
2006 
2007  }
2008  }
2009 #endif
2010 
2011  SCIPdebugMessage(" -> old LP iterations: %"SCIP_LONGINT_FORMAT"\n", SCIPgetNLPIterations(scip));
2012 
2013  /* solve LP;
2014  * errors in the LP solver should not kill the overall solving process, if the LP is just needed for a heuristic.
2015  * hence in optimized mode, the return code is caught and a warning is printed, only in debug mode, SCIP will stop.
2016  */
2017 #ifdef NDEBUG
2018  {
2019  SCIP_RETCODE retstat;
2020  retstat = SCIPsolveProbingLP(scip, -1, &lperror, NULL);
2021  if( retstat != SCIP_OKAY )
2022  {
2023  SCIPwarningMessage(scip, "Error while solving LP in SHIFTANDPROPAGATE heuristic; LP solve terminated with code <%d>\n",
2024  retstat);
2025  }
2026  }
2027 #else
2028  SCIP_CALL( SCIPsolveProbingLP(scip, -1, &lperror, NULL) );
2029 #endif
2030 
2031  SCIPdebugMessage(" -> new LP iterations: %"SCIP_LONGINT_FORMAT"\n", SCIPgetNLPIterations(scip));
2032  SCIPdebugMessage(" -> error=%u, status=%d\n", lperror, SCIPgetLPSolstat(scip));
2033 
2034  /* check if this is a feasible solution */
2035  if( !lperror && SCIPgetLPSolstat(scip) == SCIP_LPSOLSTAT_OPTIMAL )
2036  {
2037  /* copy the current LP solution to the working solution */
2038  SCIP_CALL( SCIPlinkLPSol(scip, sol) );
2039  }
2040  else
2041  stored = FALSE;
2042 
2043  SCIPstatistic( heurdata->lpsolstat = SCIPgetLPSolstat(scip) );
2044  }
2045  /* check solution for feasibility, and add it to solution store if possible.
2046  * Neither integrality nor feasibility of LP rows have to be checked, because they
2047  * are guaranteed by the heuristic at this stage.
2048  */
2049  if( stored )
2050  {
2051 #ifndef NDEBUG
2052  SCIP_CALL( SCIPtrySol(scip, sol, FALSE, TRUE, TRUE, TRUE, &stored) );
2053 #else
2054  /* @todo: maybe bounds don't need to be checked, in this case put an assert concerning stored ?????????? */
2055  SCIP_CALL( SCIPtrySol(scip, sol, FALSE, TRUE, FALSE, FALSE, &stored) );
2056 #endif
2057  if( stored )
2058  {
2059  SCIPdebugMessage("found feasible shifted solution:\n");
2060  SCIPdebug( SCIP_CALL( SCIPprintSol(scip, sol, NULL, FALSE) ) );
2061  *result = SCIP_FOUNDSOL;
2062  SCIPstatisticMessage(" Shiftandpropagate solution value: %16.9g \n", SCIPgetSolOrigObj(scip, sol));
2063  }
2064  }
2065  }
2066  else
2067  {
2068  SCIPdebugMessage("Solution constructed by heuristic is already known to be infeasible\n");
2069  }
2070 
2071  SCIPstatistic(
2072  heurdata->nremainingviols = nviolatedrows;
2073  heurdata->nredundantrows = nredundantrows;
2074  );
2075 
2076  TERMINATE2:
2077  /* free allocated memory in reverse order of allocation */
2078  for( c = matrix->ndiscvars - 1; c >= 0; --c )
2079  {
2080  SCIP_VAR* var;
2081 
2082  var = SCIPcolGetVar(heurdata->lpcols[c]);
2083  assert(var != NULL);
2084  assert(eventdatas[c] != NULL);
2085 
2086  SCIP_CALL( SCIPdropVarEvent(scip, var, SCIP_EVENTTYPE_BOUNDCHANGED, eventhdlr, eventdatas[c], -1) );
2087  SCIPfreeBuffer(scip, &(eventdatas[c]));
2088  }
2089  SCIPfreeBufferArray(scip, &eventdatas);
2090 
2091  if( violatedvarrows != NULL )
2092  {
2093  assert(heurdata->sortkey == 'v' || heurdata->sortkey == 't');
2094  SCIPfreeBufferArray(scip, &violatedvarrows);
2095  }
2096  /* free all allocated memory */
2097  SCIPfreeBufferArray(scip, &violatedrowpos);
2098  SCIPfreeBufferArray(scip, &violatedrows);
2099  SCIPfreeBufferArray(scip, &violationchange);
2100  SCIPfreeBufferArray(scip, &steps);
2101  SCIPfreeBufferArray(scip, &heurdata->rowweights);
2102  SCIPfreeBufferArray(scip, &permutation);
2103  SCIP_CALL( SCIPfreeSol(scip, &sol) );
2104 
2105  eventhdlrdata->nviolatedrows = NULL;
2106  eventhdlrdata->violatedrowpos = NULL;
2107  eventhdlrdata->violatedrows = NULL;
2108 
2109  TERMINATE:
2110  /* terminate probing mode and free the remaining memory */
2111  SCIPstatistic(
2112  heurdata->ncutoffs += ncutoffs;
2113  heurdata->nprobings += nprobings;
2114  heurdata->nlpiters = SCIPgetNLPIterations(scip) - heurdata->nlpiters;
2115  );
2116 
2117  SCIP_CALL( SCIPendProbing(scip) );
2118  SCIPfreeBufferArray(scip, &heurdata->lpcols);
2119  freeMatrix(scip, &matrix);
2120  eventhdlrdata->matrix = NULL;
2121 
2122  return SCIP_OKAY;
2123 }
2124 
2125 /** event handler execution method for the heuristic which catches all
2126  * events in which a lower or upper bound were tightened */
2127 static
2128 SCIP_DECL_EVENTEXEC(eventExecShiftandpropagate)
2129 { /*lint --e{715}*/
2130  SCIP_EVENTHDLRDATA* eventhdlrdata;
2131  SCIP_VAR* var;
2132  SCIP_COL* col;
2133  SCIP_Real lb;
2134  SCIP_Real ub;
2135  int colpos;
2136  CONSTRAINTMATRIX* matrix;
2137  SCIP_HEURDATA* heurdata;
2138 
2139  assert(scip != NULL);
2140  assert(eventhdlr != NULL);
2141  assert(strcmp(EVENTHDLR_NAME, SCIPeventhdlrGetName(eventhdlr)) == 0);
2142 
2143  eventhdlrdata = SCIPeventhdlrGetData(eventhdlr);
2144  assert(eventhdlrdata != NULL);
2145 
2146  matrix = eventhdlrdata->matrix;
2147 
2148  heurdata = eventhdlrdata->heurdata;
2149  assert(heurdata != NULL && heurdata->lpcols != NULL);
2150 
2151  colpos = eventdata->colpos;
2152 
2153  assert(0 <= colpos && colpos < matrix->ndiscvars);
2154 
2155  col = heurdata->lpcols[colpos];
2156  var = SCIPcolGetVar(col);
2157 
2158  lb = SCIPvarGetLbLocal(var);
2159  ub = SCIPvarGetUbLocal(var);
2160 
2161  updateTransformation(scip, matrix, eventhdlrdata->heurdata, colpos, &(matrix->transformshiftvals[colpos]),
2162  lb, ub, eventhdlrdata->violatedrows, eventhdlrdata->violatedrowpos, eventhdlrdata->nviolatedrows);
2163 
2164  return SCIP_OKAY;
2165 }
2166 
2167 /*
2168  * primal heuristic specific interface methods
2169  */
2170 
2171 /** creates the shiftandpropagate primal heuristic and includes it in SCIP */
2173  SCIP* scip /**< SCIP data structure */
2174  )
2175 {
2176  SCIP_HEURDATA* heurdata;
2177  SCIP_HEUR* heur;
2178  SCIP_EVENTHDLRDATA* eventhandlerdata;
2179  SCIP_EVENTHDLR* eventhdlr;
2180 
2181 
2182  SCIP_CALL( SCIPallocMemory(scip, &eventhandlerdata) );
2183  eventhandlerdata->matrix = NULL;
2184 
2185  eventhdlr = NULL;
2187  eventExecShiftandpropagate, eventhandlerdata) );
2188  assert(eventhdlr != NULL);
2189 
2190  /* create Shiftandpropagate primal heuristic data */
2191  SCIP_CALL( SCIPallocMemory(scip, &heurdata) );
2192  heurdata->rowweights = NULL;
2193  heurdata->nlpcols = 0;
2194  heurdata->eventhdlr = eventhdlr;
2195 
2196  /* include primal heuristic */
2197  SCIP_CALL( SCIPincludeHeurBasic(scip, &heur,
2199  HEUR_MAXDEPTH, HEUR_TIMING, HEUR_USESSUBSCIP, heurExecShiftandpropagate, heurdata) );
2200 
2201  assert(heur != NULL);
2202 
2203  /* set non-NULL pointers to callback methods */
2204  SCIP_CALL( SCIPsetHeurCopy(scip, heur, heurCopyShiftandpropagate) );
2205  SCIP_CALL( SCIPsetHeurFree(scip, heur, heurFreeShiftandpropagate) );
2206  SCIP_CALL( SCIPsetHeurInit(scip, heur, heurInitShiftandpropagate) );
2207  SCIP_CALL( SCIPsetHeurExit(scip, heur, heurExitShiftandpropagate) );
2208 
2209 
2210  /* add shiftandpropagate primal heuristic parameters */
2211  SCIP_CALL( SCIPaddIntParam(scip, "heuristics/"HEUR_NAME"/nproprounds", "The number of propagation rounds used for each propagation",
2212  &heurdata->nproprounds, TRUE, DEFAULT_NPROPROUNDS, -1, 1000, NULL, NULL) );
2213  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/relax", "Should continuous variables be relaxed?",
2214  &heurdata->relax, TRUE, DEFAULT_RELAX, NULL, NULL) );
2215  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/probing", "Should domains be reduced by probing?",
2216  &heurdata->probing, TRUE, DEFAULT_PROBING, NULL, NULL) );
2217  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/onlywithoutsol", "Should heuristic only be executed if no primal solution was found, yet?",
2218  &heurdata->onlywithoutsol, TRUE, DEFAULT_ONLYWITHOUTSOL, NULL, NULL) );
2219  SCIP_CALL( SCIPaddIntParam(scip, "heuristics/"HEUR_NAME"/cutoffbreaker", "The number of cutoffs before heuristic stops",
2220  &heurdata->cutoffbreaker, TRUE, DEFAULT_CUTOFFBREAKER, -1, 1000000, NULL, NULL) );
2221  SCIP_CALL( SCIPaddCharParam(scip, "heuristics/"HEUR_NAME"/sortkey", "the key for variable sorting: (n)orms down, norms (u)p, (v)iolations down, viola(t)ions up, or (r)andom",
2222  &heurdata->sortkey, TRUE, DEFAULT_SORTKEY, SORTKEYS, NULL, NULL) );
2223  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/sortvars", "Should variables be sorted for the heuristic?",
2224  &heurdata->sortvars, TRUE, DEFAULT_SORTVARS, NULL, NULL));
2225  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/"HEUR_NAME"/collectstats", "should variable statistics be collected during probing?",
2226  &heurdata->collectstats, TRUE, DEFAULT_COLLECTSTATS, NULL, NULL) );
2227  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/stopafterfeasible", "Should the heuristic stop calculating optimal shift values when no more rows are violated?",
2228  &heurdata->stopafterfeasible, TRUE, DEFAULT_STOPAFTERFEASIBLE, NULL, NULL) );
2229  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/preferbinaries", "Should binary variables be shifted first?",
2230  &heurdata->preferbinaries, TRUE, DEFAULT_PREFERBINARIES, NULL, NULL) );
2231  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/nozerofixing", "should variables with a zero shifting value be delayed instead of being fixed?",
2232  &heurdata->nozerofixing, TRUE, DEFAULT_NOZEROFIXING, NULL, NULL) );
2233  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/fixbinlocks", "should binary variables with no locks in one direction be fixed to that direction?",
2234  &heurdata->fixbinlocks, TRUE, DEFAULT_FIXBINLOCKS, NULL, NULL) );
2235  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/normalize", "should coefficients and left/right hand sides be normalized by max row coeff?",
2236  &heurdata->normalize, TRUE, DEFAULT_NORMALIZE, NULL, NULL) );
2237  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/updateweights", "should row weight be increased every time the row is violated?",
2238  &heurdata->updateweights, TRUE, DEFAULT_UPDATEWEIGHTS, NULL, NULL) );
2239  SCIP_CALL( SCIPaddBoolParam(scip, "heuristics/shiftandpropagate/impliscontinuous", "should implicit integer variables be treated as continuous variables?",
2240  &heurdata->impliscontinuous, TRUE, DEFAULT_IMPLISCONTINUOUS, NULL, NULL) );
2241  return SCIP_OKAY;
2242 }
void SCIPsortRealInt(SCIP_Real *realarray, int *intarray, int len)
SCIP_RETCODE SCIPsolveProbingLP(SCIP *scip, int itlim, SCIP_Bool *lperror, SCIP_Bool *cutoff)
Definition: scip.c:29936
int SCIPgetNVars(SCIP *scip)
Definition: scip.c:10071
const char * SCIPvarGetName(SCIP_VAR *var)
Definition: var.c:15788
#define DEFAULT_FIXBINLOCKS
preroot heuristic that alternatingly fixes variables and propagates domains
SCIP_Bool SCIPisInfinity(SCIP *scip, SCIP_Real val)
Definition: scip.c:38360
void SCIPdisableVarHistory(SCIP *scip)
Definition: scip.c:21143
static void relaxVar(SCIP *scip, SCIP_VAR *var, CONSTRAINTMATRIX *matrix, SCIP_Bool normalize)
const char * SCIPheurGetName(SCIP_HEUR *heur)
Definition: heur.c:591
#define HEUR_USESSUBSCIP
SCIP_RETCODE SCIPsetHeurCopy(SCIP *scip, SCIP_HEUR *heur, SCIP_DECL_HEURCOPY((*heurcopy)))
Definition: scip.c:7014
SCIP_RETCODE SCIPbacktrackProbing(SCIP *scip, int probingdepth)
Definition: scip.c:29573
#define DEFAULT_SORTKEY
SCIP_Bool SCIPvarIsBinary(SCIP_VAR *var)
Definition: var.c:15968
void SCIPwarningMessage(SCIP *scip, const char *formatstr,...)
Definition: scip.c:1206
#define HEUR_DESC
SCIP_Bool SCIPisFeasLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip.c:38667
SCIP_Longint SCIPgetNLPIterations(SCIP *scip)
Definition: scip.c:34147
SCIP_VAR ** SCIPgetVars(SCIP *scip)
Definition: scip.c:10026
static SCIP_DECL_HEUREXEC(heurExecShiftandpropagate)
void SCIPenableVarHistory(SCIP *scip)
Definition: scip.c:21124
#define NULL
Definition: lpi_spx.cpp:129
SCIP_Real SCIPvarGetLbLocal(SCIP_VAR *var)
Definition: var.c:16426
#define HEUR_NAME
SCIP_Bool SCIPisStopped(SCIP *scip)
Definition: scip.c:1111
SCIP_Real SCIPfeasFloor(SCIP *scip, SCIP_Real val)
Definition: scip.c:38803
SCIP_Real SCIProwGetLhs(SCIP_ROW *row)
Definition: lp.c:18637
SCIP_COL ** SCIProwGetCols(SCIP_ROW *row)
Definition: lp.c:18583
SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:16380
SCIP_RETCODE SCIPsetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var, SCIP_Real val)
Definition: scip.c:31639
struct SCIP_EventhdlrData SCIP_EVENTHDLRDATA
Definition: type_event.h:129
SCIP_Bool SCIPcolIsInLP(SCIP_COL *col)
Definition: lp.c:18470
SCIP_RETCODE SCIPincludeHeurBasic(SCIP *scip, SCIP_HEUR **heur, const char *name, const char *desc, char dispchar, int priority, int freq, int freqofs, int maxdepth, unsigned int timingmask, SCIP_Bool usessubscip, SCIP_DECL_HEUREXEC((*heurexec)), SCIP_HEURDATA *heurdata)
Definition: scip.c:6969
int SCIProwGetNLPNonz(SCIP_ROW *row)
Definition: lp.c:18572
SCIP_RETCODE SCIPpropagateProbing(SCIP *scip, int maxproprounds, SCIP_Bool *cutoff, SCIP_Longint *ndomredsfound)
Definition: scip.c:29766
#define FALSE
Definition: def.h:52
#define DEFAULT_CUTOFFBREAKER
static SCIP_DECL_EVENTEXEC(eventExecShiftandpropagate)
#define DEFAULT_RELAX
static void freeMatrix(SCIP *scip, CONSTRAINTMATRIX **matrix)
SCIP_RETCODE SCIPincludeEventhdlrBasic(SCIP *scip, SCIP_EVENTHDLR **eventhdlrptr, const char *name, const char *desc, SCIP_DECL_EVENTEXEC((*eventexec)), SCIP_EVENTHDLRDATA *eventhdlrdata)
Definition: scip.c:7209
#define TRUE
Definition: def.h:51
#define SCIPdebug(x)
Definition: pub_message.h:74
static void checkViolations(SCIP *scip, CONSTRAINTMATRIX *matrix, int *violatedrows, int *violatedrowpos, int *nviolatedrows, int *nredundantrows, int *violatedvarrows)
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:53
#define SCIPstatisticMessage
Definition: pub_message.h:104
#define HEUR_DISPCHAR
SCIP_Bool SCIPisFeasZero(SCIP *scip, SCIP_Real val)
Definition: scip.c:38743
#define DEFAULT_SORTVARS
#define DEFAULT_UPDATEWEIGHTS
struct SCIP_HeurData SCIP_HEURDATA
Definition: type_heur.h:44
void SCIPsortDownIntInt(int *intarray1, int *intarray2, int len)
#define SCIPfreeBuffer(scip, ptr)
Definition: scip.h:19219
SCIP_LPSOLSTAT SCIPgetLPSolstat(SCIP *scip)
Definition: scip.c:24136
#define SCIPdebugMessage
Definition: pub_message.h:77
#define SCIPallocBufferArray(scip, ptr, num)
Definition: scip.h:19214
SCIP_Bool SCIPvarIsInLP(SCIP_VAR *var)
Definition: var.c:16104
#define DEFAULT_WEIGHT_INEQUALITY
void SCIPsortDownRealInt(SCIP_Real *realarray, int *intarray, int len)
SCIP_Real SCIProwGetConstant(SCIP_ROW *row)
Definition: lp.c:18603
int SCIPgetNContVars(SCIP *scip)
Definition: scip.c:10251
enum SCIP_LPSolStat SCIP_LPSOLSTAT
Definition: type_lp.h:42
#define DEFAULT_RANDSEED
#define SCIP_EVENTTYPE_BOUNDCHANGED
Definition: type_event.h:99
#define HEUR_TIMING
SCIP_Real SCIPgetRowMaxCoef(SCIP *scip, SCIP_ROW *row)
Definition: scip.c:25739
SCIP_Real SCIPfeasCeil(SCIP *scip, SCIP_Real val)
Definition: scip.c:38815
#define DEFAULT_ONLYWITHOUTSOL
#define HEUR_PRIORITY
SCIP_RETCODE SCIPgetLPColsData(SCIP *scip, SCIP_COL ***cols, int *ncols)
Definition: scip.c:24369
SCIP_Real SCIPgetSolOrigObj(SCIP *scip, SCIP_SOL *sol)
Definition: scip.c:31855
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.c:3388
SCIP_VARSTATUS SCIPvarGetStatus(SCIP_VAR *var)
Definition: var.c:15907
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.c:3414
SCIP_Bool SCIPvarIsIntegral(SCIP_VAR *var)
Definition: var.c:15979
static SCIP_DECL_HEUREXIT(heurExitShiftandpropagate)
SCIP_RETCODE SCIPcatchVarEvent(SCIP *scip, SCIP_VAR *var, SCIP_EVENTTYPE eventtype, SCIP_EVENTHDLR *eventhdlr, SCIP_EVENTDATA *eventdata, int *filterpos)
Definition: scip.c:33378
SCIP_Bool SCIPisLPConstructed(SCIP *scip)
Definition: scip.c:24069
int SCIPcolGetLPPos(SCIP_COL *col)
Definition: lp.c:18448
static void updateTransformation(SCIP *scip, CONSTRAINTMATRIX *matrix, SCIP_HEURDATA *heurdata, int varindex, SCIP_Real *transformshiftval, SCIP_Real lb, SCIP_Real ub, int *violatedrows, int *violatedrowpos, int *nviolatedrows)
#define DEFAULT_NORMALIZE
SCIP_RETCODE SCIPconstructLP(SCIP *scip, SCIP_Bool *cutoff)
Definition: scip.c:24092
SCIP_ROW ** SCIPcolGetRows(SCIP_COL *col)
Definition: lp.c:18506
SCIP_Bool SCIPisFeasPositive(SCIP *scip, SCIP_Real val)
Definition: scip.c:38755
static SCIP_Bool varIsDiscrete(SCIP_VAR *var, SCIP_Bool impliscontinuous)
#define SCIPallocMemory(scip, ptr)
Definition: scip.h:19159
SCIP_Bool SCIPcolIsIntegral(SCIP_COL *col)
Definition: lp.c:18427
#define SCIPallocBuffer(scip, ptr)
Definition: scip.h:19213
#define DEFAULT_STOPAFTERFEASIBLE
#define SCIPerrorMessage
Definition: pub_message.h:45
#define SCIPfreeBufferArray(scip, ptr)
Definition: scip.h:19221
SCIP_Real SCIProwGetRhs(SCIP_ROW *row)
Definition: lp.c:18647
void SCIPsortIntInt(int *intarray1, int *intarray2, int len)
static void transformVariable(SCIP *scip, CONSTRAINTMATRIX *matrix, SCIP_HEURDATA *heurdata, int colpos)
static SCIP_Bool colIsDiscrete(SCIP_COL *col, SCIP_Bool impliscontinuous)
SCIP_Bool SCIPisFeasNegative(SCIP *scip, SCIP_Real val)
Definition: scip.c:38767
int SCIPgetProbingDepth(SCIP *scip)
Definition: scip.c:29546
SCIP_Bool SCIPisLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip.c:38292
SCIP_RETCODE SCIPcreateSol(SCIP *scip, SCIP_SOL **sol, SCIP_HEUR *heur)
Definition: scip.c:30856
void SCIPsolSetHeur(SCIP_SOL *sol, SCIP_HEUR *heur)
Definition: sol.c:2233
#define SORTKEYS
SCIP_Bool SCIPisFeasEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip.c:38648
struct SCIP_EventData SCIP_EVENTDATA
Definition: type_event.h:146
#define DEFAULT_PREFERBINARIES
SCIP_COL * SCIPvarGetCol(SCIP_VAR *var)
Definition: var.c:16093
SCIP_RETCODE SCIPgetLPRowsData(SCIP *scip, SCIP_ROW ***rows, int *nrows)
Definition: scip.c:24447
int SCIPvarGetNLocksUp(SCIP_VAR *var)
Definition: var.c:3214
SCIP_Real SCIPinfinity(SCIP *scip)
Definition: scip.c:38349
void SCIPsortPtr(void **ptrarray, SCIP_DECL_SORTPTRCOMP((*ptrcomp)), int len)
#define SCIP_CALL(x)
Definition: def.h:258
SCIP_EVENTHDLRDATA * SCIPeventhdlrGetData(SCIP_EVENTHDLR *eventhdlr)
Definition: event.c:288
SCIP_Bool SCIPhasCurrentNodeLP(SCIP *scip)
Definition: scip.c:24051
#define DEFAULT_NOZEROFIXING
int SCIPcolGetNLPNonz(SCIP_COL *col)
Definition: lp.c:18495
static SCIP_DECL_HEURINIT(heurInitShiftandpropagate)
const char * SCIProwGetName(SCIP_ROW *row)
Definition: lp.c:18696
SCIP_Bool SCIPisFeasGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip.c:38705
void SCIPheurSetData(SCIP_HEUR *heur, SCIP_HEURDATA *heurdata)
Definition: heur.c:512
struct ConstraintMatrix CONSTRAINTMATRIX
#define EVENTHDLR_NAME
SCIP_Real SCIPvarGetUbLocal(SCIP_VAR *var)
Definition: var.c:16436
#define SCIP_Bool
Definition: def.h:49
#define HEUR_MAXDEPTH
#define HEUR_FREQOFS
SCIP_RETCODE SCIPincludeHeurShiftandpropagate(SCIP *scip)
#define DEFAULT_NPROPROUNDS
#define MAX(x, y)
Definition: tclique_def.h:75
enum TransformStatus TRANSFORMSTATUS
SCIP_RETCODE SCIPfixVarProbing(SCIP *scip, SCIP_VAR *var, SCIP_Real fixedval)
Definition: scip.c:29709
#define DEFAULT_IMPLISCONTINUOUS
static SCIP_RETCODE initMatrix(SCIP *scip, CONSTRAINTMATRIX *matrix, SCIP_HEURDATA *heurdata, SCIP_Bool normalize, int *nmaxrows, SCIP_Bool relax, SCIP_Bool *initialized, SCIP_Bool *infeasible)
#define SCIPfreeMemoryArray(scip, ptr)
Definition: scip.h:19178
SCIP_RETCODE SCIPstartProbing(SCIP *scip)
Definition: scip.c:29476
#define BMScopyMemoryArray(ptr, source, num)
Definition: memory.h:102
#define DEFAULT_WEIGHT_EQUALITY
SCIP_RETCODE SCIPendProbing(SCIP *scip)
Definition: scip.c:29607
#define SCIPallocMemoryArray(scip, ptr, num)
Definition: scip.h:19161
#define HEUR_FREQ
static void getColumnData(CONSTRAINTMATRIX *matrix, int colindex, SCIP_Real **valpointer, int **indexpointer, int *ncolvals)
SCIP_RETCODE SCIPlinkLPSol(SCIP *scip, SCIP_SOL *sol)
Definition: scip.c:31444
SCIP_RETCODE SCIPsetHeurExit(SCIP *scip, SCIP_HEUR *heur, SCIP_DECL_HEUREXIT((*heurexit)))
Definition: scip.c:7062
SCIP_RETCODE SCIPsetHeurInit(SCIP *scip, SCIP_HEUR *heur, SCIP_DECL_HEURINIT((*heurinit)))
Definition: scip.c:7046
SCIP_RETCODE SCIPdropVarEvent(SCIP *scip, SCIP_VAR *var, SCIP_EVENTTYPE eventtype, SCIP_EVENTHDLR *eventhdlr, SCIP_EVENTDATA *eventdata, int filterpos)
Definition: scip.c:33424
SCIP_VARTYPE SCIPvarGetType(SCIP_VAR *var)
Definition: var.c:15953
#define SCIPfreeMemory(scip, ptr)
Definition: scip.h:19176
SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:16370
static SCIP_DECL_HEURCOPY(heurCopyShiftandpropagate)
static SCIP_Real retransformVariable(SCIP *scip, CONSTRAINTMATRIX *matrix, SCIP_VAR *var, int varindex, SCIP_Real solvalue)
SCIP_Real * SCIProwGetVals(SCIP_ROW *row)
Definition: lp.c:18593
const char * SCIPeventhdlrGetName(SCIP_EVENTHDLR *eventhdlr)
Definition: event.c:278
static void getRowData(CONSTRAINTMATRIX *matrix, int rowindex, SCIP_Real **valpointer, SCIP_Real *lhs, SCIP_Real *rhs, int **indexpointer, int *nrowvals)
SCIP_RETCODE SCIPchgVarUbProbing(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound)
Definition: scip.c:29674
SCIP_Bool SCIPisFeasLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip.c:38686
SCIP_RETCODE SCIPsetHeurFree(SCIP *scip, SCIP_HEUR *heur, SCIP_DECL_HEURFREE((*heurfree)))
Definition: scip.c:7030
#define SCIPstatistic(x)
Definition: pub_message.h:101
#define SCIP_Real
Definition: def.h:123
SCIP_RETCODE SCIPnewProbingNode(SCIP *scip)
Definition: scip.c:29513
SCIP_RETCODE SCIPaddCharParam(SCIP *scip, const char *name, const char *desc, char *valueptr, SCIP_Bool isadvanced, char defaultvalue, const char *allowedvalues, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip.c:3498
SCIP_RETCODE SCIPchgVarLbProbing(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound)
Definition: scip.c:29640
SCIP_RETCODE SCIPflushLP(SCIP *scip)
Definition: scip.c:24116
void SCIPpermuteIntArray(int *array, int begin, int end, unsigned int *randseed)
Definition: misc.c:7367
SCIP_RETCODE SCIPprintRow(SCIP *scip, SCIP_ROW *row, FILE *file)
Definition: scip.c:26010
static SCIP_DECL_SORTPTRCOMP(heurSortColsShiftandpropagate)
SCIP_Bool SCIPisFeasGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip.c:38724
#define SCIP_Longint
Definition: def.h:107
int SCIProwGetLPPos(SCIP_ROW *row)
Definition: lp.c:18826
#define EVENTHDLR_DESC
#define DEFAULT_PROBING
SCIP_HEURDATA * SCIPheurGetData(SCIP_HEUR *heur)
Definition: heur.c:502
#define DEFAULT_PROPBREAKER
enum SCIP_Vartype SCIP_VARTYPE
Definition: type_var.h:58
#define DEFAULT_COLLECTSTATS
#define BMSclearMemoryArray(ptr, num)
Definition: memory.h:98
SCIP_VAR * SCIPcolGetVar(SCIP_COL *col)
Definition: lp.c:18407
SCIP_SOL * SCIPgetBestSol(SCIP *scip)
Definition: scip.c:32531
SCIP_Real * SCIPcolGetVals(SCIP_COL *col)
Definition: lp.c:18516
#define SCIPABORT()
Definition: def.h:230
static SCIP_DECL_HEURFREE(heurFreeShiftandpropagate)
int SCIPgetNLPRows(SCIP *scip)
Definition: scip.c:24503
int SCIPvarGetNLocksDown(SCIP_VAR *var)
Definition: var.c:3159
SCIP_RETCODE SCIPtrySol(SCIP *scip, SCIP_SOL *sol, SCIP_Bool printreason, SCIP_Bool checkbounds, SCIP_Bool checkintegrality, SCIP_Bool checklprows, SCIP_Bool *stored)
Definition: scip.c:32978
static SCIP_RETCODE getOptimalShiftingValue(SCIP *scip, CONSTRAINTMATRIX *matrix, int varindex, int direction, int *rowweights, SCIP_Real *steps, int *violationchange, SCIP_Real *beststep, int *rowviolations)
static void updateViolations(SCIP *scip, CONSTRAINTMATRIX *matrix, int rowindex, int *violatedrows, int *violatedrowpos, int *nviolatedrows, int *rowweights, SCIP_Bool updateweights)
SCIP_RETCODE SCIPprintSol(SCIP *scip, SCIP_SOL *sol, FILE *file, SCIP_Bool printzeros)
Definition: scip.c:32161
SCIP_RETCODE SCIPfreeSol(SCIP *scip, SCIP_SOL **sol)
Definition: scip.c:31403