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