Scippy

    SCIP

    Solving Constraint Integer Programs

    heuristics.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 heuristics.c
    26 * @ingroup OTHER_CFILES
    27 * @brief methods commonly used by primal heuristics
    28 * @author Gregor Hendel
    29 */
    30
    31#include "scip/heuristics.h"
    32#include "scip/cons_linear.h"
    33#include "scip/scip_exact.h"
    34#include "scip/scipdefplugins.h"
    35#include "scip/stat.h"
    36#include "scip/struct_scip.h"
    37
    38#include "scip/pub_heur.h"
    39
    40/* the indicator and SOS1 constraint handlers are included for the diving algorithm SCIPperformGenericDivingAlgorithm() */
    41#include "scip/cons_indicator.h"
    42#include "scip/cons_sos1.h"
    43
    44#define MINLPITER 10000 /**< minimal number of LP iterations allowed in each LP solving call */
    45
    46
    47/** solve probing LP */
    48static
    50 SCIP* scip, /**< SCIP data structure */
    51 SCIP_DIVESET* diveset, /**< diving settings */
    52 SCIP_Longint maxnlpiterations, /**< maximum number of allowed LP iterations */
    53 SCIP_DIVECONTEXT divecontext, /**< context for diving statistics */
    54 SCIP_Bool* lperror, /**< pointer to store if an internal LP error occurred */
    55 SCIP_Bool* cutoff /**< pointer to store whether the LP was infeasible */
    56 )
    57{
    58 int lpiterationlimit;
    59 SCIP_RETCODE retstat;
    60 SCIP_Longint nlpiterations;
    61
    62 assert(lperror != NULL);
    63 assert(cutoff != NULL);
    64
    65 nlpiterations = SCIPgetNLPIterations(scip);
    66
    67 /* allow at least MINLPITER more iterations so as not to run out of LP iterations during this solve */
    68 lpiterationlimit = (int)(maxnlpiterations - SCIPdivesetGetNLPIterations(diveset, divecontext));
    69 lpiterationlimit = MAX(lpiterationlimit, MINLPITER);
    70
    71 retstat = SCIPsolveProbingLP(scip, lpiterationlimit, lperror, cutoff);
    72
    73 /* Errors in the LP solver should not kill the overall solving process, if the LP is just needed for a heuristic.
    74 * Hence in optimized mode, the return code is caught and a warning is printed, only in debug mode, SCIP will stop.
    75 */
    76#ifdef NDEBUG
    77 if( retstat != SCIP_OKAY )
    78 {
    79 SCIPwarningMessage(scip, "Error while solving LP in %s diving heuristic; LP solve terminated with code <%d>.\n", SCIPdivesetGetName(diveset), retstat);
    80 }
    81#else
    82 SCIP_CALL( retstat );
    83#endif
    84
    85 /* update iteration count */
    86 SCIPupdateDivesetLPStats(scip, diveset, SCIPgetNLPIterations(scip) - nlpiterations, divecontext);
    87
    88 return SCIP_OKAY;
    89}
    90
    91/** select the next variable and type of diving */
    92static
    94 SCIP* scip, /**< SCIP data structure */
    95 SCIP_DIVESET* diveset, /**< dive set */
    96 SCIP_SOL* worksol, /**< current working solution */
    97 SCIP_Bool onlylpbranchcands, /**< should only LP branching candidates be considered? */
    98 SCIP_Bool storelpcandscores, /**< should the scores of the LP candidates be updated? */
    99 SCIP_VAR** lpcands, /**< LP branching candidates, or NULL if not needed */
    100 SCIP_Real * lpcandssol, /**< solution values LP branching candidates, or NULL if not needed */
    101 SCIP_Real* lpcandsfrac, /**< fractionalities of LP branching candidates, or NULL if not needed*/
    102 SCIP_Real* lpcandsscores, /**< array with LP branching candidate scores, or NULL */
    103 SCIP_Bool* lpcandroundup, /**< array to remember whether the preferred branching direction is upwards */
    104 int* nviollpcands, /**< pointer to store the number of LP candidates whose solution value already violates local bounds */
    105 int nlpcands, /**< number of current LP cands */
    106 SCIP_Bool* enfosuccess, /**< pointer to store whether a candidate was sucessfully found */
    107 SCIP_Bool* infeasible /**< pointer to store whether the diving can be immediately aborted because it is infeasible */
    108 )
    109{
    110 assert(scip != NULL);
    111 assert(worksol != NULL);
    112 assert(!onlylpbranchcands || lpcandsscores != NULL);
    113 assert(!onlylpbranchcands || lpcandroundup != NULL);
    114 assert(enfosuccess != NULL);
    115 assert(infeasible != NULL);
    116 assert(nviollpcands != NULL);
    117
    118 *nviollpcands = 0;
    119 /* we use diving solution enforcement provided by the constraint handlers */
    120 if( !onlylpbranchcands )
    121 {
    122 SCIP_CALL( SCIPgetDiveBoundChanges(scip, diveset, worksol, enfosuccess, infeasible) );
    123 }
    124 else
    125 {
    126 int c;
    127 int bestcandidx;
    128 SCIP_Real bestscore;
    129 SCIP_Real score;
    130
    131 bestscore = SCIP_REAL_MIN;
    132 bestcandidx = -1;
    133
    135
    136 /* search for the candidate that maximizes the dive set score function and whose solution value is still feasible */
    137 for( c = 0; c < nlpcands; ++c )
    138 {
    139 assert(SCIPgetSolVal(scip, worksol, lpcands[c]) == lpcandssol[c]); /*lint !e777 doesn't like comparing floats for equality */
    140
    141 /* scores are kept in arrays for faster reuse */
    142 if( storelpcandscores )
    143 {
    144 SCIP_CALL( SCIPgetDivesetScore(scip, diveset, SCIP_DIVETYPE_INTEGRALITY, lpcands[c], lpcandssol[c],
    145 lpcandsfrac[c], &lpcandsscores[c], &lpcandroundup[c]) );
    146 }
    147
    148 score = lpcandsscores[c];
    149 /* update the best candidate if it has a higher score and a solution value which does not violate one of the local bounds */
    150 if( SCIPisLE(scip, SCIPvarGetLbLocal(lpcands[c]), lpcandssol[c]) && SCIPisGE(scip, SCIPvarGetUbLocal(lpcands[c]), lpcandssol[c]) )
    151 {
    152 if( score > bestscore )
    153 {
    154 bestcandidx = c;
    155 bestscore = score;
    156 }
    157 }
    158 else
    159 ++(*nviollpcands);
    160 }
    161
    162 /* there is no guarantee that a candidate is found since local bounds might render all solution values infeasible */
    163 *enfosuccess = (bestcandidx >= 0);
    164 if( *enfosuccess )
    165 {
    166 /* if we want to round up the best candidate, it is added as the preferred bound change */
    168 SCIPceil(scip, lpcandssol[bestcandidx]), lpcandroundup[bestcandidx]) );
    170 SCIPfloor(scip, lpcandssol[bestcandidx]), ! lpcandroundup[bestcandidx]) );
    171 }
    172 }
    173 return SCIP_OKAY;
    174}
    175
    176/** return the LP iteration budget suggestion for this dive set */
    177static
    179 SCIP* scip, /**< SCIP data structure */
    180 SCIP_DIVESET* diveset, /**< dive set data structure */
    181 SCIP_DIVECONTEXT divecontext /**< context for diving statistics */
    182 )
    183{
    184 SCIP_Longint iterlimit;
    185
    186 iterlimit = (SCIP_Longint)(SCIPdivesetGetMaxLPIterQuot(diveset) * (SCIPdivesetGetSolSuccess(diveset, divecontext)+1.0)/(SCIPdivesetGetNCalls(diveset, divecontext)+1.0) * SCIPgetNNodeLPIterations(scip));
    187 iterlimit += SCIPdivesetGetMaxLPIterOffset(diveset);
    188 iterlimit -= SCIPdivesetGetNLPIterations(diveset, divecontext);
    189
    190 return iterlimit;
    191}
    192
    193/** performs a diving within the limits of the @p diveset parameters
    194 *
    195 * This method performs a diving according to the settings defined by the diving settings @p diveset; Contrary to the
    196 * name, SCIP enters probing mode (not diving mode) and dives along a path into the tree. Domain propagation
    197 * is applied at every node in the tree, whereas probing LPs might be solved less frequently.
    198 *
    199 * Starting from the current LP solution, the algorithm selects candidates which maximize the
    200 * score defined by the @p diveset and whose solution value has not yet been rendered infeasible by propagation,
    201 * and propagates the bound change on this candidate.
    202 *
    203 * The algorithm iteratively selects the the next (unfixed) candidate in the list, until either enough domain changes
    204 * or the resolve frequency of the LP trigger an LP resolve (and hence, the set of potential candidates changes),
    205 * or the last node is proven to be infeasible. It optionally backtracks and tries the
    206 * other branching direction.
    207 *
    208 * After the set of remaining candidates is empty or the targeted depth is reached, the node LP is
    209 * solved, and the old candidates are replaced by the new LP candidates.
    210 *
    211 * @see heur_guideddiving.c for an example implementation of a dive set controlling the diving algorithm.
    212 *
    213 * @note the node from where the algorithm is called is checked for a basic LP solution. If the solution
    214 * is non-basic, e.g., when barrier without crossover is used, the method returns without performing a dive.
    215 *
    216 * @note currently, when multiple diving heuristics call this method and solve an LP at the same node, only the first
    217 * call will be executed, see SCIPgetLastDiveNode()
    218 *
    219 * @todo generalize method to work correctly with pseudo or external branching/diving candidates
    220 */
    222 SCIP* scip, /**< SCIP data structure */
    223 SCIP_DIVESET* diveset, /**< settings for diving */
    224 SCIP_SOL* worksol, /**< non-NULL working solution */
    225 SCIP_HEUR* heur, /**< the calling primal heuristic */
    226 SCIP_RESULT* result, /**< SCIP result pointer */
    227 SCIP_Bool nodeinfeasible, /**< is the current node known to be infeasible? */
    228 SCIP_Longint iterlim, /**< nonnegative iteration limit for the LP solves, or -1 for dynamic setting */
    229 int nodelimit, /**< nonnegative probing node limit or -1 if no limit should be used */
    230 SCIP_Real lpresolvedomchgquot,/**< percentage of immediate domain changes during probing to trigger LP resolve or -1
    231 * if diveset specific default should be used */
    232 SCIP_DIVECONTEXT divecontext /**< context for diving statistics */
    233 )
    234{
    235 SCIP_CONSHDLR* indconshdlr; /* constraint handler for indicator constraints */
    236 SCIP_CONSHDLR* sos1conshdlr; /* constraint handler for SOS1 constraints */
    237 SCIP_VAR** lpcands;
    238 SCIP_Real* lpcandssol;
    239
    240 SCIP_VAR** previouscands;
    241 SCIP_Real* lpcandsscores;
    242 SCIP_Real* previousvals;
    243 SCIP_Real* lpcandsfrac;
    244 SCIP_Bool* lpcandroundup;
    245 SCIP_Real searchubbound;
    246 SCIP_Real searchavgbound;
    247 SCIP_Real searchbound;
    248 SCIP_Real ubquot;
    249 SCIP_Real avgquot;
    250 SCIP_Longint maxnlpiterations;
    251 SCIP_Longint domreds;
    252 int startndivecands;
    253 int depth;
    254 int maxdepth;
    255 int maxdivedepth;
    256 int totalnbacktracks;
    257 int totalnprobingnodes;
    258 int lastlpdepth;
    259 int previouscandssize;
    260 int lpcandsscoressize;
    261 int nviollpcands;
    262 SCIP_Longint oldnsolsfound;
    263 SCIP_Longint oldnbestsolsfound;
    264 SCIP_Longint oldnconflictsfound;
    265
    266 SCIP_Bool success;
    267 SCIP_Bool leafsol;
    268 SCIP_Bool enfosuccess;
    269 SCIP_Bool lperror;
    270 SCIP_Bool cutoff;
    271 SCIP_Bool backtracked;
    272 SCIP_Bool backtrack;
    273 SCIP_Bool onlylpbranchcands;
    274
    275 int nlpcands;
    276 int lpsolvefreq;
    277
    278 assert(scip != NULL);
    279 assert(heur != NULL);
    280 assert(result != NULL);
    281 assert(worksol != NULL);
    282
    283 *result = SCIP_DELAYED;
    284
    285 /* do not call heuristic in node that was already detected to be infeasible */
    286 if( nodeinfeasible )
    287 return SCIP_OKAY;
    288
    289 /* only call heuristic, if an optimal LP solution is at hand */
    291 return SCIP_OKAY;
    292
    293 /* only call heuristic, if the LP solution is basic (which allows fast resolve in diving) */
    294 if( !SCIPisLPSolBasic(scip) )
    295 return SCIP_OKAY;
    296
    297 /* when heuristic was called by scheduler, be less restrictive */
    298 if( divecontext != SCIP_DIVECONTEXT_SCHEDULER )
    299 {
    300 /* don't dive two times at the same node */
    302 return SCIP_OKAY;
    303
    304 /* only call heuristic, if the LP objective value is smaller than the cutoff bound */
    306 return SCIP_OKAY;
    307 }
    308
    309 *result = SCIP_DIDNOTRUN;
    310
    311 /* only try to dive, if we are in the correct part of the tree, given by minreldepth and maxreldepth */
    312 depth = SCIPgetDepth(scip);
    313 maxdepth = SCIPgetMaxDepth(scip);
    314 maxdepth = MAX(maxdepth, 30);
    315 if( divecontext != SCIP_DIVECONTEXT_SCHEDULER &&
    316 (depth < SCIPdivesetGetMinRelDepth(diveset) * maxdepth || depth > SCIPdivesetGetMaxRelDepth(diveset) * maxdepth) )
    317 {
    318 return SCIP_OKAY;
    319 }
    320
    321 /* calculate the maximal number of LP iterations */
    322 if( iterlim < 0 )
    323 {
    324 maxnlpiterations = SCIPdivesetGetNLPIterations(diveset, divecontext) + getDivesetIterLimit(scip, diveset, divecontext);
    325 }
    326 else
    327 {
    328 maxnlpiterations = SCIPdivesetGetNLPIterations(diveset, divecontext) + iterlim;
    329 }
    330
    331 /* don't try to dive, if we took too many LP iterations during diving */
    332 if( SCIPdivesetGetNLPIterations(diveset, divecontext) >= maxnlpiterations )
    333 return SCIP_OKAY;
    334
    335 /* allow at least a certain number of LP iterations in this dive */
    336 if( SCIPdivesetGetNLPIterations(diveset, divecontext) + MINLPITER > maxnlpiterations )
    337 maxnlpiterations = SCIPdivesetGetNLPIterations(diveset, divecontext) + MINLPITER;
    338
    339 /* these constraint handlers are required for polishing an LP relaxation solution beyond rounding */
    340 indconshdlr = SCIPfindConshdlr(scip, "indicator");
    341 sos1conshdlr = SCIPfindConshdlr(scip, "SOS1");
    342
    343 SCIP_CALL( SCIPgetLPBranchCands(scip, &lpcands, &lpcandssol, &lpcandsfrac, &nlpcands, NULL, NULL) );
    344
    345 onlylpbranchcands = SCIPdivesetUseOnlyLPBranchcands(diveset);
    346 /* don't try to dive, if there are no diving candidates */
    347 if( onlylpbranchcands && nlpcands == 0 )
    348 return SCIP_OKAY;
    349
    350 /* calculate the objective search bound */
    351 if( SCIPgetNSolsFound(scip) == 0 )
    352 {
    353 ubquot = SCIPdivesetGetUbQuotNoSol(diveset);
    354 avgquot = SCIPdivesetGetAvgQuotNoSol(diveset);
    355 }
    356 else
    357 {
    358 ubquot = SCIPdivesetGetUbQuot(diveset);
    359 avgquot = SCIPdivesetGetAvgQuot(diveset);
    360 }
    361
    362 if( ubquot > 0.0 )
    363 searchubbound = SCIPgetLowerbound(scip) + ubquot * (SCIPgetCutoffbound(scip) - SCIPgetLowerbound(scip));
    364 else
    365 searchubbound = SCIPinfinity(scip);
    366
    367 if( avgquot > 0.0 )
    368 searchavgbound = SCIPgetLowerbound(scip) + avgquot * (SCIPgetAvgLowerbound(scip) - SCIPgetLowerbound(scip));
    369 else
    370 searchavgbound = SCIPinfinity(scip);
    371
    372 searchbound = MIN(searchubbound, searchavgbound);
    373
    375 searchbound = SCIPceil(scip, searchbound);
    376
    377 /* calculate the maximal diving depth: 10 * min{number of integer variables, max depth}*/
    379 assert(maxdivedepth >= 0);
    380 if ( sos1conshdlr != NULL )
    381 maxdivedepth += SCIPgetNSOS1Vars(sos1conshdlr);
    382 maxdivedepth = MIN(maxdivedepth, maxdepth);
    383 maxdivedepth *= 10;
    384
    385 /* if lpresolvedomchgquot is not explicitly given (so -1), then use the default for the current diveset */
    386 if( lpresolvedomchgquot < 0 )
    387 lpresolvedomchgquot = SCIPdivesetGetLPResolveDomChgQuot(diveset);
    388 assert(lpresolvedomchgquot > 0.0 && lpresolvedomchgquot <= 1.0);
    389
    390 *result = SCIP_DIDNOTFIND;
    391
    392 /* start probing mode */
    394
    395 /* enables collection of variable statistics during probing */
    397
    398 SCIPdebugMsg(scip, "(node %" SCIP_LONGINT_FORMAT ") executing %s heuristic: depth=%d, %d fractionals, dualbound=%g, avgbound=%g, cutoffbound=%g, searchbound=%g\n",
    401
    402 /* allocate buffer storage for previous candidates and their branching values for pseudo cost updates */
    403 lpsolvefreq = SCIPdivesetGetLPSolveFreq(diveset);
    404 previouscandssize = MAX(1, lpsolvefreq);
    405 SCIP_CALL( SCIPallocBufferArray(scip, &previouscands, previouscandssize) );
    406 SCIP_CALL( SCIPallocBufferArray(scip, &previousvals, previouscandssize) );
    407
    408 /* keep some statistics */
    409 lperror = FALSE;
    410 cutoff = FALSE;
    411 lastlpdepth = -1;
    412 startndivecands = nlpcands;
    413 totalnbacktracks = 0;
    414 totalnprobingnodes = 0;
    415 oldnsolsfound = SCIPgetNSolsFound(scip);
    416 oldnbestsolsfound = SCIPgetNBestSolsFound(scip);
    417 oldnconflictsfound = SCIPgetNConflictConssFound(scip);
    418
    419 /* link the working solution to the dive set */
    420 SCIPdivesetSetWorkSolution(diveset, worksol);
    421
    422 if( onlylpbranchcands )
    423 {
    424 SCIP_CALL( SCIPallocBufferArray(scip, &lpcandsscores, nlpcands) );
    425 SCIP_CALL( SCIPallocBufferArray(scip, &lpcandroundup, nlpcands) );
    426
    427 lpcandsscoressize = nlpcands;
    428 }
    429 else
    430 {
    431 lpcandroundup = NULL;
    432 lpcandsscores = NULL;
    433 lpcandsscoressize = 0;
    434 }
    435
    436 enfosuccess = TRUE;
    437 leafsol = FALSE;
    438
    439 /* LP loop; every time a new LP was solved, conditions are checked
    440 * dive as long we are in the given objective, depth and iteration limits and fractional variables exist, but
    441 * - if possible, we dive at least with the depth 10
    442 * - if the number of fractional variables decreased at least with 1 variable per 2 dive depths, we continue diving
    443 */
    444 while( !lperror && !cutoff && SCIPgetLPSolstat(scip) == SCIP_LPSOLSTAT_OPTIMAL && enfosuccess
    445 && ( nodelimit == -1 || totalnprobingnodes < nodelimit )
    446 && (SCIPgetProbingDepth(scip) < 10
    447 || nlpcands <= startndivecands - SCIPgetProbingDepth(scip) / 2
    448 || (SCIPgetProbingDepth(scip) < maxdivedepth && SCIPdivesetGetNLPIterations(diveset, divecontext) < maxnlpiterations && SCIPgetLPObjval(scip) < searchbound))
    449 && !SCIPisStopped(scip) )
    450 {
    451 SCIP_Real lastlpobjval;
    452 int c;
    453 SCIP_Bool allroundable;
    454 SCIP_Bool infeasible;
    455 int prevcandsinsertidx;
    456
    457 /* remember the last LP depth */
    458 assert(lastlpdepth < SCIPgetProbingDepth(scip));
    459 lastlpdepth = SCIPgetProbingDepth(scip);
    460 domreds = 0;
    461
    462 SCIPdebugMsg(scip, "%s heuristic continues diving at depth %d, %d candidates left\n",
    463 SCIPdivesetGetName(diveset), lastlpdepth, nlpcands);
    464
    465 /* loop over candidates and determine if they are roundable */
    466 allroundable = TRUE;
    467 c = 0;
    468 while( allroundable && c < nlpcands )
    469 {
    470 if( SCIPvarMayRoundDown(lpcands[c]) || SCIPvarMayRoundUp(lpcands[c]) )
    471 allroundable = TRUE;
    472 else
    473 allroundable = FALSE;
    474 ++c;
    475 }
    476
    477 /* if all candidates are roundable, try to round the solution */
    478 if( allroundable )
    479 {
    480 success = FALSE;
    481
    482 /* working solution must be linked to LP solution */
    483 SCIP_CALL( SCIPlinkLPSol(scip, worksol) );
    484 /* create solution from diving LP and try to round it */
    485 SCIP_CALL( SCIProundSol(scip, worksol, &success) );
    486
    487 /* adjust SOS1 constraints */
    488 if( success && sos1conshdlr != NULL )
    489 {
    490 SCIP_Bool changed;
    491 SCIP_CALL( SCIPmakeSOS1sFeasible(scip, sos1conshdlr, worksol, &changed, &success) );
    492 }
    493
    494 /* succesfully rounded solutions are tried for primal feasibility */
    495 if( success )
    496 {
    497 SCIP_Bool changed = FALSE;
    498 SCIPdebugMsg(scip, "%s found roundable primal solution: obj=%g\n", SCIPdivesetGetName(diveset), SCIPgetSolOrigObj(scip, worksol));
    499
    500 /* adjust indicator constraints */
    501 if( indconshdlr != NULL )
    502 {
    503 SCIP_CALL( SCIPmakeIndicatorsFeasible(scip, indconshdlr, worksol, &changed) );
    504 }
    505
    506 success = FALSE;
    507
    508 /* try to add solution to SCIP */
    509 SCIP_CALL( SCIPtrySol(scip, worksol, FALSE, FALSE, FALSE, FALSE, changed, &success) );
    510
    511 /* check, if solution was feasible and good enough */
    512 if( success )
    513 {
    514 SCIPdebugMsg(scip, " -> solution was feasible and good enough\n");
    515 *result = SCIP_FOUNDSOL;
    516 leafsol = (nlpcands == 0);
    517
    518 /* the rounded solution found above led to a cutoff of the node LP solution */
    520 {
    521 cutoff = TRUE;
    522 break;
    523 }
    524 }
    525 }
    526 }
    527
    528 /* working solution must be linked to LP solution */
    530 lastlpobjval = SCIPgetLPObjval(scip);
    531 SCIP_CALL( SCIPlinkLPSol(scip, worksol) );
    532
    533 /* we must explicitly store the solution values by unlinking the solution, otherwise,
    534 * the working solution may contain wrong entries, if, e.g., a backtrack occurred after an
    535 * intermediate LP resolve or the LP was resolved during conflict analysis.
    536 */
    537 SCIP_CALL( SCIPunlinkSol(scip, worksol) );
    538
    539 /* ensure array sizes for the diving on the fractional variables */
    540 if( onlylpbranchcands && nlpcands > lpcandsscoressize )
    541 {
    542 assert(nlpcands > 0);
    543 assert(lpcandsscores != NULL);
    544 assert(lpcandroundup != NULL);
    545
    546 SCIP_CALL( SCIPreallocBufferArray(scip, &lpcandsscores, nlpcands) );
    547 SCIP_CALL( SCIPreallocBufferArray(scip, &lpcandroundup, nlpcands) );
    548
    549 lpcandsscoressize = nlpcands;
    550 }
    551
    552 enfosuccess = FALSE;
    553 /* select the next diving action by selecting appropriate dive bound changes for the preferred and alternative child */
    554 SCIP_CALL( selectNextDiving(scip, diveset, worksol, onlylpbranchcands, SCIPgetProbingDepth(scip) == lastlpdepth,
    555 lpcands, lpcandssol, lpcandsfrac, lpcandsscores, lpcandroundup, &nviollpcands, nlpcands,
    556 &enfosuccess, &infeasible) );
    557
    558 /* if we did not succeed finding an enforcement, the solution is potentially feasible and we break immediately */
    559 if( ! enfosuccess )
    560 break;
    561
    562 prevcandsinsertidx = -1;
    563
    564 /* start propagating candidate variables
    565 * - until the desired targetdepth is reached,
    566 * - or there is no further candidate variable left because of intermediate bound changes,
    567 * - or a cutoff is detected
    568 */
    569 do
    570 {
    571 SCIP_VAR* bdchgvar;
    572 SCIP_Real bdchgvalue;
    573 SCIP_Longint localdomreds = 0;
    574 SCIP_BRANCHDIR bdchgdir;
    575 int nbdchanges;
    576
    577 /* ensure that a new candidate was successfully determined (usually at the end of the previous loop iteration) */
    578 assert(enfosuccess);
    579 bdchgvar = NULL;
    580 bdchgvalue = SCIP_INVALID;
    581 bdchgdir = SCIP_BRANCHDIR_AUTO;
    582
    583 backtracked = FALSE;
    584 do
    585 {
    586 int d;
    587 SCIP_VAR** bdchgvars;
    588 SCIP_BRANCHDIR* bdchgdirs;
    589 SCIP_Real* bdchgvals;
    590
    591 nbdchanges = 0;
    592 /* get the bound change information stored in the dive set */
    593 SCIPgetDiveBoundChangeData(scip, &bdchgvars, &bdchgdirs, &bdchgvals, &nbdchanges, !backtracked);
    594
    595 assert(nbdchanges > 0);
    596 assert(bdchgvars != NULL);
    597 assert(bdchgdirs != NULL);
    598 assert(bdchgvals != NULL);
    599
    600 /* return if we reached the depth limit */
    602 {
    603 SCIPdebugMsg(scip, "depth limit reached, we stop diving immediately.\n");
    604 goto TERMINATE;
    605 }
    606
    607 /* return if we reached the node limit */
    608 if( nodelimit != -1 && totalnprobingnodes >= nodelimit )
    609 {
    610 SCIPdebugMsg(scip, "node limit reached, we stop diving immediately.\n");
    611 goto TERMINATE;
    612 }
    613
    614 /* dive deeper into the tree */
    616 ++totalnprobingnodes;
    617
    618 /* apply all suggested domain changes of the variables */
    619 for( d = 0; d < nbdchanges; ++d )
    620 {
    621 SCIP_Real lblocal;
    622 SCIP_Real ublocal;
    623 SCIP_Bool infeasbdchange;
    624
    625 /* get the next bound change data: variable, direction, and value */
    626 bdchgvar = bdchgvars[d];
    627 bdchgvalue = bdchgvals[d];
    628 bdchgdir = bdchgdirs[d];
    629
    630 assert(bdchgvar != NULL);
    631 lblocal = SCIPvarGetLbLocal(bdchgvar);
    632 ublocal = SCIPvarGetUbLocal(bdchgvar);
    633
    634 SCIPdebugMsg(scip, " dive %d/%d, LP iter %" SCIP_LONGINT_FORMAT "/%" SCIP_LONGINT_FORMAT ": var <%s>, oldbounds=[%g,%g], newbounds=[%g,%g]\n",
    635 SCIPgetProbingDepth(scip), maxdivedepth, SCIPdivesetGetNLPIterations(diveset, divecontext), maxnlpiterations,
    636 SCIPvarGetName(bdchgvar),
    637 lblocal, ublocal,
    638 bdchgdir == SCIP_BRANCHDIR_DOWNWARDS ? lblocal : bdchgvalue,
    639 bdchgdir == SCIP_BRANCHDIR_UPWARDS ? ublocal : bdchgvalue);
    640
    641 infeasbdchange = FALSE;
    642 /* tighten the lower and/or upper bound depending on the bound change type */
    643 switch( bdchgdir )
    644 {
    646 /* test if bound change is possible, otherwise propagation might have deduced the same
    647 * bound already or numerical troubles might have occurred */
    648 if( SCIPisFeasLE(scip, bdchgvalue, lblocal) || SCIPisFeasGT(scip, bdchgvalue, ublocal) )
    649 infeasbdchange = TRUE;
    650 else
    651 {
    652 /* round variable up */
    653 SCIP_CALL( SCIPchgVarLbProbing(scip, bdchgvar, bdchgvalue) );
    654 }
    655 break;
    657 /* test if bound change is possible, otherwise propagation might have deduced the same
    658 * bound already or numerical troubles might have occurred */
    659 if( SCIPisFeasGE(scip, bdchgvalue, ublocal) || SCIPisFeasLT(scip, bdchgvalue, lblocal) )
    660 infeasbdchange = TRUE;
    661 else
    662 {
    663 /* round variable down */
    664 SCIP_CALL( SCIPchgVarUbProbing(scip, bdchgvar, bdchgvalue) );
    665 }
    666 break;
    668 /* test if bound change is possible, otherwise propagation might have deduced the same
    669 * bound already or numerical troubles might have occurred */
    670 if( SCIPisFeasLT(scip, bdchgvalue, lblocal) || SCIPisFeasGT(scip, bdchgvalue, ublocal) ||
    671 (SCIPisFeasEQ(scip, lblocal, ublocal) && nbdchanges < 2) )
    672 infeasbdchange = TRUE;
    673 else
    674 {
    675 /* fix variable to the bound change value */
    676 if( SCIPisFeasLT(scip, lblocal, bdchgvalue) )
    677 {
    678 SCIP_CALL( SCIPchgVarLbProbing(scip, bdchgvar, bdchgvalue) );
    679 }
    680 if( SCIPisFeasGT(scip, ublocal, bdchgvalue) )
    681 {
    682 SCIP_CALL( SCIPchgVarUbProbing(scip, bdchgvar, bdchgvalue) );
    683 }
    684 }
    685 break;
    687 default:
    688 SCIPerrorMessage("Error: Unsupported bound change direction <%d> specified for diving, aborting\n",bdchgdirs[d]);
    689 SCIPABORT();
    690 return SCIP_INVALIDDATA; /*lint !e527*/
    691 }
    692 /* if the variable domain has been shrunk in the meantime, numerical troubles may have
    693 * occured or variable was fixed by propagation while backtracking => Abort diving!
    694 */
    695 if( infeasbdchange )
    696 {
    697 SCIPdebugMsg(scip, "\nSelected variable <%s> domain already [%g,%g] as least as tight as desired bound change, diving aborted \n",
    698 SCIPvarGetName(bdchgvar), SCIPvarGetLbLocal(bdchgvar), SCIPvarGetUbLocal(bdchgvar));
    699 cutoff = TRUE;
    700 break;
    701 }
    702 }
    703 /* break loop immediately if we detected a cutoff */
    704 if( cutoff )
    705 break;
    706
    707 /* apply domain propagation */
    708 localdomreds = 0;
    709 SCIP_CALL( SCIPpropagateProbing(scip, 0, &cutoff, &localdomreds) );
    710
    711 /* add the number of bound changes we applied by ourselves after propagation, otherwise the counter would have been reset */
    712 localdomreds += nbdchanges;
    713
    714 SCIPdebugMsg(scip, "domain reductions: %" SCIP_LONGINT_FORMAT " (total: %" SCIP_LONGINT_FORMAT ")\n",
    715 localdomreds, domreds + localdomreds);
    716
    717 /* resolve the diving LP if the diving resolve frequency is reached or a sufficient number of intermediate bound changes
    718 * was reached
    719 */
    720 if( ! cutoff
    721 && ((lpsolvefreq > 0 && ((SCIPgetProbingDepth(scip) - lastlpdepth) % lpsolvefreq) == 0)
    722 || (domreds + localdomreds > SCIPdivesetGetLPResolveDomChgQuot(diveset) * SCIPgetNVars(scip))
    723 || (onlylpbranchcands && nviollpcands > (int)(SCIPdivesetGetLPResolveDomChgQuot(diveset) * nlpcands))) )
    724 {
    725 SCIP_CALL( solveLP(scip, diveset, maxnlpiterations, divecontext, &lperror, &cutoff) );
    726
    727 /* lp errors lead to early termination */
    728 if( lperror )
    729 {
    730 cutoff = TRUE;
    731 break;
    732 }
    733 }
    734
    735 /* perform backtracking if a cutoff was detected */
    736 if( cutoff && !backtracked && SCIPdivesetUseBacktrack(diveset) )
    737 {
    738 SCIPdebugMsg(scip, " *** cutoff detected at level %d - backtracking\n", SCIPgetProbingDepth(scip));
    740 ++totalnbacktracks;
    741 backtracked = TRUE;
    742 backtrack = TRUE;
    743 cutoff = FALSE;
    744 }
    745 else
    746 backtrack = FALSE;
    747 }
    748 while( backtrack );
    749
    750 /* we add the domain reductions from the last evaluated node */
    751 domreds += localdomreds;
    752
    753 /* store candidate for pseudo cost update and choose next candidate only if no cutoff was detected */
    754 if( ! cutoff )
    755 {
    756 if( nbdchanges == 1 && (bdchgdir == SCIP_BRANCHDIR_UPWARDS || bdchgdir == SCIP_BRANCHDIR_DOWNWARDS) )
    757 {
    758 ++prevcandsinsertidx;
    759 assert(prevcandsinsertidx <= SCIPgetProbingDepth(scip) - lastlpdepth - 1);
    760 assert(SCIPgetProbingDepth(scip) > 0);
    761 assert(bdchgvar != NULL);
    762 assert(bdchgvalue != SCIP_INVALID); /*lint !e777 doesn't like comparing floats for equality */
    763
    764 /* extend array in case of a dynamic, domain change based LP resolve strategy */
    765 if( prevcandsinsertidx >= previouscandssize )
    766 {
    767 previouscandssize *= 2;
    768 SCIP_CALL( SCIPreallocBufferArray(scip, &previouscands, previouscandssize) );
    769 SCIP_CALL( SCIPreallocBufferArray(scip, &previousvals, previouscandssize) );
    770 }
    771 assert(previouscandssize > prevcandsinsertidx);
    772
    773 /* store candidate for pseudo cost update */
    774 previouscands[prevcandsinsertidx] = bdchgvar;
    775 previousvals[prevcandsinsertidx] = bdchgvalue;
    776 }
    777
    778 /* choose next candidate variable and resolve the LP if none is found. */
    780 {
    781 assert(SCIPgetProbingDepth(scip) > lastlpdepth);
    782 enfosuccess = FALSE;
    783
    784 /* select the next diving action */
    785 SCIP_CALL( selectNextDiving(scip, diveset, worksol, onlylpbranchcands, SCIPgetProbingDepth(scip) == lastlpdepth,
    786 lpcands, lpcandssol, lpcandsfrac, lpcandsscores, lpcandroundup, &nviollpcands, nlpcands,
    787 &enfosuccess, &infeasible) );
    788
    789 /* in case of an unsuccesful candidate search, we solve the node LP */
    790 if( !enfosuccess )
    791 {
    792 SCIP_CALL( solveLP(scip, diveset, maxnlpiterations, divecontext, &lperror, &cutoff) );
    793
    794 /* check for an LP error and terminate in this case, cutoffs lead to termination anyway */
    795 if( lperror )
    796 cutoff = TRUE;
    797
    798 /* enfosuccess must be set to TRUE for entering the main LP loop again */
    799 enfosuccess = TRUE;
    800 }
    801 }
    802 }
    803 }
    804 while( !cutoff && SCIPgetLPSolstat(scip) == SCIP_LPSOLSTAT_NOTSOLVED );
    805
    806 assert(cutoff || lperror || SCIPgetLPSolstat(scip) != SCIP_LPSOLSTAT_NOTSOLVED);
    807
    810
    811 /* check new LP candidates and use the LP Objective gain to update pseudo cost information */
    812 if( ! cutoff && SCIPgetLPSolstat(scip) == SCIP_LPSOLSTAT_OPTIMAL )
    813 {
    814 int v;
    815 SCIP_Real gain;
    816
    817 SCIP_CALL( SCIPgetLPBranchCands(scip, &lpcands, &lpcandssol, NULL, &nlpcands, NULL, NULL) );
    818
    819 SCIPdebugMsg(scip, " -> lpsolstat=%d, objval=%g/%g, nlpcands=%d\n", SCIPgetLPSolstat(scip), SCIPgetLPObjval(scip), searchbound, nlpcands);
    820 /* distribute the gain equally over all variables that we rounded since the last LP */
    821 gain = SCIPgetLPObjval(scip) - lastlpobjval;
    822 gain = MAX(gain, 0.0);
    823 gain /= (1.0 * (SCIPgetProbingDepth(scip) - lastlpdepth));
    824
    825 /* loop over previously fixed candidates and share gain improvement */
    826 for( v = 0; v <= prevcandsinsertidx; ++v )
    827 {
    828 SCIP_VAR* cand = previouscands[v];
    829 SCIP_Real val = previousvals[v];
    830 SCIP_Real solval = SCIPgetSolVal(scip, worksol, cand);
    831
    832 /* todo: should the gain be shared with a smaller weight, instead of dividing the gain itself? */
    833 /* it may happen that a variable had an integral solution value beforehand, e.g., for indicator variables */
    834 if( ! SCIPisZero(scip, val - solval) )
    835 {
    836 SCIP_CALL( SCIPupdateVarPseudocost(scip, cand, val - solval, gain, 1.0) );
    837 }
    838 }
    839 }
    840 else
    841 nlpcands = 0;
    842 }
    843
    844 success = FALSE;
    845 /* check if a solution has been found */
    846 if( !enfosuccess && !lperror && !cutoff && SCIPgetLPSolstat(scip) == SCIP_LPSOLSTAT_OPTIMAL )
    847 {
    848 /* create solution from diving LP */
    849 SCIP_CALL( SCIPlinkLPSol(scip, worksol) );
    850 SCIPdebugMsg(scip, "%s found primal solution: obj=%g\n", SCIPdivesetGetName(diveset), SCIPgetSolOrigObj(scip, worksol));
    851
    852 /* try to add solution to SCIP */
    853 SCIP_CALL( SCIPtrySol(scip, worksol, FALSE, FALSE, FALSE, FALSE, FALSE, &success) );
    854
    855 /* check, if solution was feasible and good enough */
    856 if( success )
    857 {
    858 SCIPdebugMsg(scip, " -> solution was feasible and good enough\n");
    859 *result = SCIP_FOUNDSOL;
    860 leafsol = TRUE;
    861 }
    862 }
    863
    864 SCIPupdateDivesetStats(scip, diveset, totalnprobingnodes, totalnbacktracks, SCIPgetNSolsFound(scip) - oldnsolsfound,
    865 SCIPgetNBestSolsFound(scip) - oldnbestsolsfound, SCIPgetNConflictConssFound(scip) - oldnconflictsfound, leafsol, divecontext);
    866
    867 SCIPdebugMsg(scip, "(node %" SCIP_LONGINT_FORMAT ") finished %s diveset (%s heur): %d fractionals, dive %d/%d, LP iter %" SCIP_LONGINT_FORMAT "/%" SCIP_LONGINT_FORMAT ", objval=%g/%g, lpsolstat=%d, cutoff=%u\n",
    868 SCIPgetNNodes(scip), SCIPdivesetGetName(diveset), SCIPheurGetName(heur), nlpcands, SCIPgetProbingDepth(scip), maxdivedepth, SCIPdivesetGetNLPIterations(diveset, divecontext), maxnlpiterations,
    870
    871 TERMINATE:
    872 /* end probing mode */
    874
    876
    877 if( onlylpbranchcands )
    878 {
    879 SCIPfreeBufferArray(scip, &lpcandroundup);
    880 SCIPfreeBufferArray(scip, &lpcandsscores);
    881 }
    882 SCIPfreeBufferArray(scip, &previousvals);
    883 SCIPfreeBufferArray(scip, &previouscands);
    884
    885 return SCIP_OKAY;
    886}
    887
    888
    889/** creates the rows of the subproblem */
    890static
    892 SCIP* scip, /**< original SCIP data structure */
    893 SCIP* subscip, /**< SCIP data structure for the subproblem */
    894 SCIP_HASHMAP* varmap /**< a hashmap to store the mapping of source variables to the corresponding
    895 * target variables */
    896 )
    897{
    898 SCIP_ROW** rows; /* original scip rows */
    899 SCIP_CONS* cons; /* new constraint */
    900 SCIP_VAR** consvars; /* new constraint's variables */
    901 SCIP_COL** cols; /* original row's columns */
    902
    903 SCIP_Real constant; /* constant added to the row */
    904 SCIP_Real lhs; /* left hand side of the row */
    905 SCIP_Real rhs; /* left right side of the row */
    906 SCIP_Real* vals; /* variables' coefficient values of the row */
    907
    908 int nrows;
    909 int nnonz;
    910 int i;
    911 int j;
    912
    913 /* get the rows and their number */
    914 SCIP_CALL( SCIPgetLPRowsData(scip, &rows, &nrows) );
    915
    916 /* copy all rows to linear constraints */
    917 for( i = 0; i < nrows; i++ )
    918 {
    919 /* ignore rows that are only locally valid */
    920 if( SCIProwIsLocal(rows[i]) )
    921 continue;
    922
    923 /* get the row's data */
    924 constant = SCIProwGetConstant(rows[i]);
    925 lhs = SCIProwGetLhs(rows[i]) - constant;
    926 rhs = SCIProwGetRhs(rows[i]) - constant;
    927 vals = SCIProwGetVals(rows[i]);
    928 nnonz = SCIProwGetNNonz(rows[i]);
    929 cols = SCIProwGetCols(rows[i]);
    930
    931 assert(lhs <= rhs);
    932
    933 /* allocate memory array to be filled with the corresponding subproblem variables */
    934 SCIP_CALL( SCIPallocBufferArray(scip, &consvars, nnonz) );
    935 for( j = 0; j < nnonz; j++ )
    936 consvars[j] = (SCIP_VAR*) SCIPhashmapGetImage(varmap, (SCIPcolGetVar(cols[j])));
    937
    938 /* create a new linear constraint and add it to the subproblem */
    939 SCIP_CALL( SCIPcreateConsLinear(subscip, &cons, SCIProwGetName(rows[i]), nnonz, consvars, vals, lhs, rhs,
    941 SCIP_CALL( SCIPaddCons(subscip, cons) );
    942 SCIP_CALL( SCIPreleaseCons(subscip, &cons) );
    943
    944 /* free temporary memory */
    945 SCIPfreeBufferArray(scip, &consvars);
    946 }
    947
    948 return SCIP_OKAY;
    949}
    950
    951
    952/** get a sub-SCIP copy of the transformed problem */
    954 SCIP* sourcescip, /**< source SCIP data structure */
    955 SCIP* subscip, /**< sub-SCIP used by the heuristic */
    956 SCIP_HASHMAP* varmap, /**< a hashmap to store the mapping of source variables to the corresponding
    957 * target variables */
    958 const char* suffix, /**< suffix for the problem name */
    959 SCIP_VAR** fixedvars, /**< source variables whose copies should be fixed in the target SCIP environment, or NULL */
    960 SCIP_Real* fixedvals, /**< array of fixing values for target SCIP variables, or NULL */
    961 int nfixedvars, /**< number of source variables whose copies should be fixed in the target SCIP environment, or NULL */
    962 SCIP_Bool uselprows, /**< should the linear relaxation of the problem defined by LP rows be copied? */
    963 SCIP_Bool copycuts, /**< should cuts be copied (only if uselprows == FALSE) */
    964 SCIP_Bool* success, /**< was the copying successful? */
    965 SCIP_Bool* valid /**< pointer to store whether the copying was valid, or NULL */
    966 )
    967{
    968 assert(sourcescip != NULL);
    969 assert(suffix != NULL);
    970 assert(subscip != NULL);
    971 assert(varmap != NULL);
    972 assert(success != NULL);
    973
    974 if( uselprows )
    975 {
    976 SCIP_Bool msghdlrquiet;
    977 char probname[SCIP_MAXSTRLEN];
    978
    979 /* copy all plugins */
    981
    982 /* store the quiet state of the message handler and explicitly suppress output when copying parameters */
    983 msghdlrquiet = SCIPmessagehdlrIsQuiet(subscip->messagehdlr);
    985 SCIP_CALL( SCIPcopyParamSettings(sourcescip, subscip) );
    986
    987 /* even when solving exactly, sub-SCIP heuristics should be run in floating-point mode, since the exactsol constraint
    988 * handler is in place to perform a final repair step
    989 */
    991
    992 /* restore original quiet state */
    993 SCIPsetMessagehdlrQuiet(subscip, msghdlrquiet);
    994
    995 /* set name to the original problem name and possibly add a suffix */
    996 (void) SCIPsnprintf(probname, SCIP_MAXSTRLEN, "%s_%s", SCIPgetProbName(sourcescip), suffix);
    997
    998 /* create the subproblem */
    999 SCIP_CALL( SCIPcreateProb(subscip, probname, NULL, NULL, NULL, NULL, NULL, NULL, NULL) );
    1000
    1001 /* copy all variables */
    1002 SCIP_CALL( SCIPcopyVars(sourcescip, subscip, varmap, NULL, fixedvars, fixedvals, nfixedvars, TRUE) );
    1003
    1004 /* create linear constraints from LP rows of the source problem */
    1005 SCIP_CALL( createRows(sourcescip, subscip, varmap) );
    1006 }
    1007 else
    1008 {
    1009 SCIP_CALL( SCIPcopyConsCompression(sourcescip, subscip, varmap, NULL, suffix, fixedvars, fixedvals, nfixedvars,
    1010 TRUE, FALSE, FALSE, TRUE, valid) );
    1011
    1012 if( copycuts )
    1013 {
    1014 /* copies all active cuts from cutpool of sourcescip to linear constraints in targetscip */
    1015 SCIP_CALL( SCIPcopyCuts(sourcescip, subscip, varmap, NULL, TRUE, NULL) );
    1016 }
    1017 }
    1018
    1019 /* disable bound limits in subscip since objective might be changed */
    1020 SCIP_CALL( SCIPsetRealParam(subscip, "limits/primal", SCIP_INVALID) );
    1021 SCIP_CALL( SCIPsetRealParam(subscip, "limits/dual", SCIP_INVALID) );
    1022
    1023 *success = TRUE;
    1024
    1025 return SCIP_OKAY;
    1026}
    1027
    1028/** adds a trust region neighborhood constraint to the @p targetscip
    1029 *
    1030 * a trust region constraint measures the deviation from the current incumbent solution \f$x^*\f$ by an auxiliary
    1031 * continuous variable \f$v \geq 0\f$:
    1032 * \f[
    1033 * \sum\limits_{j\in B} |x_j^* - x_j| = v
    1034 * \f]
    1035 * Only binary variables are taken into account. The deviation is penalized in the objective function using
    1036 * a positive \p violpenalty.
    1037 *
    1038 * @note: the trust region constraint creates an auxiliary variable to penalize the deviation from
    1039 * the current incumbent solution. This variable can afterwards be accessed using SCIPfindVar() by its name
    1040 * 'trustregion_violationvar'
    1041 */
    1043 SCIP* sourcescip, /**< the data structure for the main SCIP instance */
    1044 SCIP* targetscip, /**< SCIP data structure of the subproblem */
    1045 SCIP_VAR** subvars, /**< variables of the subproblem, NULL entries are ignored */
    1046 SCIP_Real violpenalty /**< the penalty for violating the trust region */
    1047 )
    1048{
    1049 SCIP_VAR* violvar;
    1050 SCIP_CONS* trustregioncons;
    1051 SCIP_VAR** consvars;
    1052 SCIP_VAR** vars;
    1053 SCIP_SOL* bestsol;
    1054
    1055 int nvars;
    1056 int nbinvars;
    1057 int nconsvars;
    1058 int i;
    1059 SCIP_Real rhs;
    1060 SCIP_Real* consvals;
    1061 char name[SCIP_MAXSTRLEN];
    1062
    1063 /* get the data of the variables and the best solution */
    1064 SCIP_CALL( SCIPgetVarsData(sourcescip, &vars, &nvars, &nbinvars, NULL, NULL, NULL) );
    1065 bestsol = SCIPgetBestSol(sourcescip);
    1066 assert(bestsol != NULL);
    1067 /* otherwise, this neighborhood would not be active in the first place */
    1068 assert(nbinvars > 0);
    1069
    1070 /* memory allocation */
    1071 SCIP_CALL( SCIPallocBufferArray(sourcescip, &consvars, nbinvars + 1) );
    1072 SCIP_CALL( SCIPallocBufferArray(sourcescip, &consvals, nbinvars + 1) );
    1073 nconsvars = 0;
    1074
    1075 /* set initial left and right hand sides of trust region constraint */
    1076 rhs = 0.0;
    1077
    1078 /* create the distance (to incumbent) function of the binary variables */
    1079 for( i = 0; i < nbinvars; i++ )
    1080 {
    1081 SCIP_Real solval;
    1082
    1083 if( subvars[i] == NULL )
    1084 continue;
    1085
    1086 solval = SCIPgetSolVal(sourcescip, bestsol, vars[i]);
    1087 assert( SCIPisFeasIntegral(sourcescip,solval) );
    1088
    1089 /* is variable i part of the binary support of bestsol? */
    1090 if( SCIPisFeasEQ(sourcescip, solval, 1.0) )
    1091 {
    1092 consvals[nconsvars] = -1.0;
    1093 rhs -= 1.0;
    1094 }
    1095 else
    1096 consvals[nconsvars] = 1.0;
    1097 consvars[nconsvars] = subvars[i];
    1098 assert(SCIPvarGetType(consvars[nconsvars]) == SCIP_VARTYPE_BINARY && !SCIPvarIsImpliedIntegral(consvars[nconsvars]));
    1099 ++nconsvars;
    1100 }
    1101
    1102 /* adding the violation variable */
    1103 (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "%s_trustregionviolvar", SCIPgetProbName(sourcescip));
    1104 SCIP_CALL( SCIPcreateVarBasic(targetscip, &violvar, name, 0.0, SCIPinfinity(targetscip), violpenalty, SCIP_VARTYPE_CONTINUOUS) );
    1105 SCIP_CALL( SCIPaddVar(targetscip, violvar) );
    1106 consvars[nconsvars] = violvar;
    1107 consvals[nconsvars] = -1.0;
    1108 ++nconsvars;
    1109
    1110 /* creates trustregion constraint and adds it to subscip */
    1111 (void) SCIPsnprintf(name, SCIP_MAXSTRLEN, "%s_trustregioncons", SCIPgetProbName(sourcescip));
    1112
    1113 SCIP_CALL( SCIPcreateConsLinear(targetscip, &trustregioncons, name, nconsvars, consvars, consvals,
    1114 rhs, rhs, TRUE, TRUE, TRUE, TRUE, TRUE, FALSE, FALSE, TRUE, TRUE, FALSE) );
    1115 SCIP_CALL( SCIPaddCons(targetscip, trustregioncons) );
    1116 SCIP_CALL( SCIPreleaseCons(targetscip, &trustregioncons) );
    1117
    1118 /* releasing the violation variable */
    1119 SCIP_CALL( SCIPreleaseVar(targetscip, &violvar) );
    1120
    1121 /* free local memory */
    1122 SCIPfreeBufferArray(sourcescip, &consvals);
    1123 SCIPfreeBufferArray(sourcescip, &consvars);
    1124
    1125 return SCIP_OKAY;
    1126}
    constraint handler for indicator constraints
    Constraint handler for linear constraints in their most general form, .
    constraint handler for SOS type 1 constraints
    #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_INVALID
    Definition: def.h:178
    #define SCIP_Bool
    Definition: def.h:91
    #define MIN(x, y)
    Definition: def.h:224
    #define SCIP_Real
    Definition: def.h:156
    #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 SCIP_REAL_MIN
    Definition: def.h:159
    #define SCIP_CALL(x)
    Definition: def.h:355
    SCIP_RETCODE SCIPmakeSOS1sFeasible(SCIP *scip, SCIP_CONSHDLR *conshdlr, SCIP_SOL *sol, SCIP_Bool *changed, SCIP_Bool *success)
    Definition: cons_sos1.c:10996
    SCIP_RETCODE SCIPmakeIndicatorsFeasible(SCIP *scip, SCIP_CONSHDLR *conshdlr, SCIP_SOL *sol, SCIP_Bool *changed)
    int SCIPgetNSOS1Vars(SCIP_CONSHDLR *conshdlr)
    Definition: cons_sos1.c:10894
    SCIP_RETCODE SCIPcreateConsLinear(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_Real *vals, SCIP_Real lhs, SCIP_Real rhs, SCIP_Bool initial, SCIP_Bool separate, SCIP_Bool enforce, SCIP_Bool check, SCIP_Bool propagate, SCIP_Bool local, SCIP_Bool modifiable, SCIP_Bool dynamic, SCIP_Bool removable, SCIP_Bool stickingatnode)
    SCIP_RETCODE SCIPcopyVars(SCIP *sourcescip, SCIP *targetscip, SCIP_HASHMAP *varmap, SCIP_HASHMAP *consmap, SCIP_VAR **fixedvars, SCIP_Real *fixedvals, int nfixedvars, SCIP_Bool global)
    Definition: scip_copy.c:1167
    SCIP_RETCODE SCIPcopyConsCompression(SCIP *sourcescip, SCIP *targetscip, SCIP_HASHMAP *varmap, SCIP_HASHMAP *consmap, const char *suffix, SCIP_VAR **fixedvars, SCIP_Real *fixedvals, int nfixedvars, SCIP_Bool global, SCIP_Bool enablepricing, SCIP_Bool threadsafe, SCIP_Bool passmessagehdlr, SCIP_Bool *valid)
    Definition: scip_copy.c:2961
    SCIP_RETCODE SCIPcopyCuts(SCIP *sourcescip, SCIP *targetscip, SCIP_HASHMAP *varmap, SCIP_HASHMAP *consmap, SCIP_Bool global, int *ncutsadded)
    Definition: scip_copy.c:2113
    SCIP_RETCODE SCIPcopyParamSettings(SCIP *sourcescip, SCIP *targetscip)
    Definition: scip_copy.c:2547
    SCIP_Bool SCIPisStopped(SCIP *scip)
    Definition: scip_general.c:759
    SCIP_RETCODE SCIPaddVar(SCIP *scip, SCIP_VAR *var)
    Definition: scip_prob.c:1907
    const char * SCIPgetProbName(SCIP *scip)
    Definition: scip_prob.c:1242
    int SCIPgetNContVars(SCIP *scip)
    Definition: scip_prob.c:2569
    SCIP_RETCODE SCIPgetVarsData(SCIP *scip, SCIP_VAR ***vars, int *nvars, int *nbinvars, int *nintvars, int *nimplvars, int *ncontvars)
    Definition: scip_prob.c:2115
    int SCIPgetNVars(SCIP *scip)
    Definition: scip_prob.c:2246
    SCIP_RETCODE SCIPaddCons(SCIP *scip, SCIP_CONS *cons)
    Definition: scip_prob.c:3274
    int SCIPgetNContImplVars(SCIP *scip)
    Definition: scip_prob.c:2522
    SCIP_RETCODE SCIPcreateProb(SCIP *scip, const char *name, SCIP_DECL_PROBDELORIG((*probdelorig)), SCIP_DECL_PROBTRANS((*probtrans)), SCIP_DECL_PROBDELTRANS((*probdeltrans)), SCIP_DECL_PROBINITSOL((*probinitsol)), SCIP_DECL_PROBEXITSOL((*probexitsol)), SCIP_DECL_PROBCOPY((*probcopy)), SCIP_PROBDATA *probdata)
    Definition: scip_prob.c:119
    SCIP_Bool SCIPisObjIntegral(SCIP *scip)
    Definition: scip_prob.c:1801
    void * SCIPhashmapGetImage(SCIP_HASHMAP *hashmap, void *origin)
    Definition: misc.c:3284
    #define SCIPdebugMsg
    Definition: scip_message.h:78
    void SCIPsetMessagehdlrQuiet(SCIP *scip, SCIP_Bool quiet)
    Definition: scip_message.c:108
    void SCIPwarningMessage(SCIP *scip, const char *formatstr,...)
    Definition: scip_message.c:120
    SCIP_RETCODE SCIPsetRealParam(SCIP *scip, const char *name, SCIP_Real value)
    Definition: scip_param.c:603
    SCIP_RETCODE SCIPgetLPBranchCands(SCIP *scip, SCIP_VAR ***lpcands, SCIP_Real **lpcandssol, SCIP_Real **lpcandsfrac, int *nlpcands, int *npriolpcands, int *nfracimplvars)
    Definition: scip_branch.c:402
    SCIP_VAR * SCIPcolGetVar(SCIP_COL *col)
    Definition: lp.c:17425
    SCIP_CONSHDLR * SCIPfindConshdlr(SCIP *scip, const char *name)
    Definition: scip_cons.c:940
    SCIP_RETCODE SCIPreleaseCons(SCIP *scip, SCIP_CONS **cons)
    Definition: scip_cons.c:1173
    SCIP_Real SCIPdivesetGetMaxLPIterQuot(SCIP_DIVESET *diveset)
    Definition: heur.c:655
    SCIP_Real SCIPdivesetGetUbQuotNoSol(SCIP_DIVESET *diveset)
    Definition: heur.c:671
    SCIP_Real SCIPdivesetGetMaxRelDepth(SCIP_DIVESET *diveset)
    Definition: heur.c:463
    SCIP_Real SCIPdivesetGetAvgQuot(SCIP_DIVESET *diveset)
    Definition: heur.c:694
    SCIP_Bool SCIPdivesetUseBacktrack(SCIP_DIVESET *diveset)
    Definition: heur.c:702
    int SCIPdivesetGetLPSolveFreq(SCIP_DIVESET *diveset)
    Definition: heur.c:710
    SCIP_Real SCIPdivesetGetMinRelDepth(SCIP_DIVESET *diveset)
    Definition: heur.c:455
    SCIP_Real SCIPdivesetGetUbQuot(SCIP_DIVESET *diveset)
    Definition: heur.c:686
    SCIP_Bool SCIPdivesetUseOnlyLPBranchcands(SCIP_DIVESET *diveset)
    Definition: heur.c:743
    SCIP_Longint SCIPdivesetGetNLPIterations(SCIP_DIVESET *diveset, SCIP_DIVECONTEXT divecontext)
    Definition: heur.c:589
    SCIP_Longint SCIPdivesetGetSolSuccess(SCIP_DIVESET *diveset, SCIP_DIVECONTEXT divecontext)
    Definition: heur.c:471
    const char * SCIPdivesetGetName(SCIP_DIVESET *diveset)
    Definition: heur.c:445
    int SCIPdivesetGetNCalls(SCIP_DIVESET *diveset, SCIP_DIVECONTEXT divecontext)
    Definition: heur.c:485
    SCIP_Real SCIPdivesetGetAvgQuotNoSol(SCIP_DIVESET *diveset)
    Definition: heur.c:679
    SCIP_Real SCIPdivesetGetLPResolveDomChgQuot(SCIP_DIVESET *diveset)
    Definition: heur.c:731
    int SCIPdivesetGetMaxLPIterOffset(SCIP_DIVESET *diveset)
    Definition: heur.c:663
    SCIP_RETCODE SCIPenableExactSolving(SCIP *scip, SCIP_Bool enable)
    Definition: scip_exact.c:151
    const char * SCIPheurGetName(SCIP_HEUR *heur)
    Definition: heur.c:1467
    SCIP_Longint SCIPgetLastDivenode(SCIP *scip)
    Definition: scip_lp.c:2710
    SCIP_Bool SCIPhasCurrentNodeLP(SCIP *scip)
    Definition: scip_lp.c:87
    SCIP_RETCODE SCIPgetLPRowsData(SCIP *scip, SCIP_ROW ***rows, int *nrows)
    Definition: scip_lp.c:576
    SCIP_LPSOLSTAT SCIPgetLPSolstat(SCIP *scip)
    Definition: scip_lp.c:174
    SCIP_Real SCIPgetLPObjval(SCIP *scip)
    Definition: scip_lp.c:253
    SCIP_Bool SCIPisLPSolBasic(SCIP *scip)
    Definition: scip_lp.c:673
    #define SCIPallocBufferArray(scip, ptr, num)
    Definition: scip_mem.h:124
    #define SCIPreallocBufferArray(scip, ptr, num)
    Definition: scip_mem.h:128
    #define SCIPfreeBufferArray(scip, ptr)
    Definition: scip_mem.h:136
    void SCIPclearDiveBoundChanges(SCIP *scip)
    int SCIPgetProbingDepth(SCIP *scip)
    Definition: scip_probing.c:199
    void SCIPupdateDivesetStats(SCIP *scip, SCIP_DIVESET *diveset, int nprobingnodes, int nbacktracks, SCIP_Longint nsolsfound, SCIP_Longint nbestsolsfound, SCIP_Longint nconflictsfound, SCIP_Bool leavewassol, SCIP_DIVECONTEXT divecontext)
    void SCIPgetDiveBoundChangeData(SCIP *scip, SCIP_VAR ***variables, SCIP_BRANCHDIR **directions, SCIP_Real **values, int *ndivebdchgs, SCIP_Bool preferred)
    SCIP_RETCODE SCIPgetDiveBoundChanges(SCIP *scip, SCIP_DIVESET *diveset, SCIP_SOL *sol, SCIP_Bool *success, SCIP_Bool *infeasible)
    SCIP_RETCODE SCIPgetDivesetScore(SCIP *scip, SCIP_DIVESET *diveset, SCIP_DIVETYPE divetype, SCIP_VAR *divecand, SCIP_Real divecandsol, SCIP_Real divecandfrac, SCIP_Real *candscore, SCIP_Bool *roundup)
    SCIP_RETCODE SCIPchgVarUbProbing(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound)
    Definition: scip_probing.c:346
    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
    void SCIPupdateDivesetLPStats(SCIP *scip, SCIP_DIVESET *diveset, SCIP_Longint niterstoadd, SCIP_DIVECONTEXT divecontext)
    SCIP_RETCODE SCIPstartProbing(SCIP *scip)
    Definition: scip_probing.c:120
    SCIP_RETCODE SCIPaddDiveBoundChange(SCIP *scip, SCIP_VAR *var, SCIP_BRANCHDIR dir, SCIP_Real value, SCIP_Bool preferred)
    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 SCIPendProbing(SCIP *scip)
    Definition: scip_probing.c:261
    SCIP_Real SCIProwGetLhs(SCIP_ROW *row)
    Definition: lp.c:17686
    int SCIProwGetNNonz(SCIP_ROW *row)
    Definition: lp.c:17607
    SCIP_COL ** SCIProwGetCols(SCIP_ROW *row)
    Definition: lp.c:17632
    SCIP_Real SCIProwGetRhs(SCIP_ROW *row)
    Definition: lp.c:17696
    SCIP_Bool SCIProwIsLocal(SCIP_ROW *row)
    Definition: lp.c:17795
    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 SCIPunlinkSol(SCIP *scip, SCIP_SOL *sol)
    Definition: scip_sol.c:1506
    SCIP_RETCODE SCIProundSol(SCIP *scip, SCIP_SOL *sol, SCIP_Bool *success)
    Definition: scip_sol.c:3123
    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_Real SCIPgetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var)
    Definition: scip_sol.c:1765
    SCIP_Real SCIPretransformObj(SCIP *scip, SCIP_Real obj)
    Definition: scip_sol.c:2132
    SCIP_Longint SCIPgetNSolsFound(SCIP *scip)
    int SCIPgetMaxDepth(SCIP *scip)
    SCIP_Longint SCIPgetNNodes(SCIP *scip)
    SCIP_Real SCIPgetDualbound(SCIP *scip)
    SCIP_Real SCIPgetLowerbound(SCIP *scip)
    SCIP_Real SCIPgetAvgLowerbound(SCIP *scip)
    SCIP_Longint SCIPgetNNodeLPIterations(SCIP *scip)
    SCIP_Real SCIPgetAvgDualbound(SCIP *scip)
    SCIP_Longint SCIPgetNBestSolsFound(SCIP *scip)
    SCIP_Real SCIPgetCutoffbound(SCIP *scip)
    SCIP_Longint SCIPgetNConflictConssFound(SCIP *scip)
    SCIP_Longint SCIPgetNLPIterations(SCIP *scip)
    SCIP_RETCODE SCIPperformGenericDivingAlgorithm(SCIP *scip, SCIP_DIVESET *diveset, SCIP_SOL *worksol, SCIP_HEUR *heur, SCIP_RESULT *result, SCIP_Bool nodeinfeasible, SCIP_Longint iterlim, int nodelimit, SCIP_Real lpresolvedomchgquot, SCIP_DIVECONTEXT divecontext)
    Definition: heuristics.c:221
    SCIP_RETCODE SCIPcopyLargeNeighborhoodSearch(SCIP *sourcescip, SCIP *subscip, SCIP_HASHMAP *varmap, const char *suffix, SCIP_VAR **fixedvars, SCIP_Real *fixedvals, int nfixedvars, SCIP_Bool uselprows, SCIP_Bool copycuts, SCIP_Bool *success, SCIP_Bool *valid)
    Definition: heuristics.c:953
    SCIP_RETCODE SCIPaddTrustregionNeighborhoodConstraint(SCIP *sourcescip, SCIP *targetscip, SCIP_VAR **subvars, SCIP_Real violpenalty)
    Definition: heuristics.c:1042
    SCIP_Bool SCIPisFeasGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
    SCIP_Real SCIPinfinity(SCIP *scip)
    SCIP_Bool SCIPisGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
    SCIP_Bool SCIPisFeasEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
    SCIP_Bool SCIPisLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
    SCIP_Real SCIPfloor(SCIP *scip, SCIP_Real val)
    SCIP_Bool SCIPisFeasLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
    SCIP_Bool SCIPisFeasLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
    SCIP_Bool SCIPisFeasIntegral(SCIP *scip, SCIP_Real val)
    SCIP_Real SCIPceil(SCIP *scip, SCIP_Real val)
    SCIP_Bool SCIPisFeasGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
    SCIP_Bool SCIPisZero(SCIP *scip, SCIP_Real val)
    SCIP_Bool SCIPisLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
    int SCIPgetDepth(SCIP *scip)
    Definition: scip_tree.c:672
    SCIP_Bool SCIPvarMayRoundUp(SCIP_VAR *var)
    Definition: var.c:4484
    SCIP_Bool SCIPvarIsImpliedIntegral(SCIP_VAR *var)
    Definition: var.c:23498
    SCIP_Real SCIPvarGetUbLocal(SCIP_VAR *var)
    Definition: var.c:24268
    SCIP_Bool SCIPvarMayRoundDown(SCIP_VAR *var)
    Definition: var.c:4473
    SCIP_VARTYPE SCIPvarGetType(SCIP_VAR *var)
    Definition: var.c:23453
    const char * SCIPvarGetName(SCIP_VAR *var)
    Definition: var.c:23267
    SCIP_RETCODE SCIPreleaseVar(SCIP *scip, SCIP_VAR **var)
    Definition: scip_var.c:1887
    SCIP_Real SCIPvarGetLbLocal(SCIP_VAR *var)
    Definition: var.c:24234
    SCIP_RETCODE SCIPupdateVarPseudocost(SCIP *scip, SCIP_VAR *var, SCIP_Real solvaldelta, SCIP_Real objdelta, SCIP_Real weight)
    Definition: scip_var.c:11122
    SCIP_RETCODE SCIPcreateVarBasic(SCIP *scip, SCIP_VAR **var, const char *name, SCIP_Real lb, SCIP_Real ub, SCIP_Real obj, SCIP_VARTYPE vartype)
    Definition: scip_var.c:184
    void SCIPenableVarHistory(SCIP *scip)
    Definition: scip_var.c:11083
    int SCIPsnprintf(char *t, int len, const char *s,...)
    Definition: misc.c:10827
    void SCIPdivesetSetWorkSolution(SCIP_DIVESET *diveset, SCIP_SOL *sol)
    Definition: heur.c:434
    #define MINLPITER
    Definition: heuristics.c:44
    static SCIP_RETCODE solveLP(SCIP *scip, SCIP_DIVESET *diveset, SCIP_Longint maxnlpiterations, SCIP_DIVECONTEXT divecontext, SCIP_Bool *lperror, SCIP_Bool *cutoff)
    Definition: heuristics.c:49
    static SCIP_RETCODE createRows(SCIP *scip, SCIP *subscip, SCIP_HASHMAP *varmap)
    Definition: heuristics.c:891
    static SCIP_RETCODE selectNextDiving(SCIP *scip, SCIP_DIVESET *diveset, SCIP_SOL *worksol, SCIP_Bool onlylpbranchcands, SCIP_Bool storelpcandscores, SCIP_VAR **lpcands, SCIP_Real *lpcandssol, SCIP_Real *lpcandsfrac, SCIP_Real *lpcandsscores, SCIP_Bool *lpcandroundup, int *nviollpcands, int nlpcands, SCIP_Bool *enfosuccess, SCIP_Bool *infeasible)
    Definition: heuristics.c:93
    static SCIP_Longint getDivesetIterLimit(SCIP *scip, SCIP_DIVESET *diveset, SCIP_DIVECONTEXT divecontext)
    Definition: heuristics.c:178
    methods commonly used by primal heuristics
    SCIP_Bool SCIPmessagehdlrIsQuiet(SCIP_MESSAGEHDLR *messagehdlr)
    Definition: message.c:910
    public methods for primal heuristics
    #define SCIPerrorMessage
    Definition: pub_message.h:64
    public methods for exact solving
    SCIP_RETCODE SCIPincludeDefaultPlugins(SCIP *scip)
    default SCIP plugins
    internal methods for problem statistics
    SCIP_MESSAGEHDLR * messagehdlr
    Definition: struct_scip.h:78
    SCIP main data structure.
    enum SCIP_DiveContext SCIP_DIVECONTEXT
    Definition: type_heur.h:73
    #define SCIP_DIVETYPE_INTEGRALITY
    Definition: type_heur.h:60
    @ SCIP_DIVECONTEXT_SCHEDULER
    Definition: type_heur.h:71
    @ SCIP_BRANCHDIR_DOWNWARDS
    Definition: type_history.h:43
    @ SCIP_BRANCHDIR_FIXED
    Definition: type_history.h:45
    @ SCIP_BRANCHDIR_AUTO
    Definition: type_history.h:46
    @ SCIP_BRANCHDIR_UPWARDS
    Definition: type_history.h:44
    enum SCIP_BranchDir SCIP_BRANCHDIR
    Definition: type_history.h:48
    @ SCIP_LPSOLSTAT_NOTSOLVED
    Definition: type_lp.h:43
    @ SCIP_LPSOLSTAT_OPTIMAL
    Definition: type_lp.h:44
    @ SCIP_LPSOLSTAT_INFEASIBLE
    Definition: type_lp.h:45
    @ SCIP_LPSOLSTAT_OBJLIMIT
    Definition: type_lp.h:47
    @ SCIP_DIDNOTRUN
    Definition: type_result.h:42
    @ SCIP_DELAYED
    Definition: type_result.h:43
    @ SCIP_DIDNOTFIND
    Definition: type_result.h:44
    @ SCIP_FOUNDSOL
    Definition: type_result.h:56
    enum SCIP_Result SCIP_RESULT
    Definition: type_result.h:61
    @ SCIP_INVALIDDATA
    Definition: type_retcode.h:52
    @ SCIP_OKAY
    Definition: type_retcode.h:42
    enum SCIP_Retcode SCIP_RETCODE
    Definition: type_retcode.h:63
    @ SCIP_VARTYPE_CONTINUOUS
    Definition: type_var.h:71
    @ SCIP_VARTYPE_BINARY
    Definition: type_var.h:64