|
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
dive.c
Go to the documentation of this file.
24 /* the indicator and SOS1 constraint handlers are included for the diving algorithm SCIPperformGenericDivingAlgorithm() */
56 /* Errors in the LP solver should not kill the overall solving process, if the LP is just needed for a heuristic.
57 * Hence in optimized mode, the return code is caught and a warning is printed, only in debug mode, SCIP will stop.
62 SCIPwarningMessage(scip, "Error while solving LP in %s diving heuristic; LP solve terminated with code <%d>.\n", SCIPdivesetGetName(diveset), retstat);
86 SCIP_Bool* lpcandroundup, /**< array to remember whether the preferred branching direction is upwards */
87 int* nviollpcands, /**< pointer to store the number of LP candidates whose solution value already violates local bounds */
90 SCIP_Bool* infeasible /**< pointer to store whether the diving can be immediately aborted because it is infeasible */
119 /* search for the candidate that maximizes the dive set score function and whose solution value is still feasible */
122 assert(SCIPgetSolVal(scip, worksol, lpcands[c]) == lpcandssol[c]); /*lint !e777 doesn't like comparing floats for equality */
127 SCIP_CALL( SCIPgetDivesetScore(scip, diveset, SCIP_DIVETYPE_INTEGRALITY, lpcands[c], lpcandssol[c], lpcandsfrac[c], &lpcandsscores[c], &lpcandroundup[c]) );
131 /* update the best candidate if it has a higher score and a solution value which does not violate one of the local bounds */
132 if( SCIPisLE(scip, SCIPvarGetLbLocal(lpcands[c]), lpcandssol[c]) && SCIPisGE(scip, SCIPvarGetUbLocal(lpcands[c]), lpcandssol[c]) )
144 /* there is no guarantee that a candidate is found since local bounds might render all solution values infeasible */
162 * This method performs a diving according to the settings defined by the diving settings @p diveset; Contrary to the
163 * name, SCIP enters probing mode (not diving mode) and dives along a path into the tree. Domain propagation
167 * score defined by the @p diveset and whose solution value has not yet been rendered infeasible by propagation,
170 * The algorithm iteratively selects the the next (unfixed) candidate in the list, until either enough domain changes
171 * or the resolve frequency of the LP trigger an LP resolve (and hence, the set of potential candidates changes),
175 * After the set of remaining candidates is empty or the targeted depth is reached, the node LP is
178 * @see heur_guideddiving.c for an example implementation of a dive set controlling the diving algorithm.
180 * @note the node from where the algorithm is called is checked for a basic LP solution. If the solution
181 * is non-basic, e.g., when barrier without crossover is used, the method returns without performing a dive.
183 * @note currently, when multiple diving heuristics call this method and solve an LP at the same node, only the first
270 /* only try to dive, if we are in the correct part of the tree, given by minreldepth and maxreldepth */
274 if( depth < SCIPdivesetGetMinRelDepth(diveset) * maxdepth || depth > SCIPdivesetGetMaxRelDepth(diveset) * maxdepth )
283 maxnlpiterations = (SCIP_Longint)((1.0 + 10*(oldsolsuccess+1.0)/(ncalls+1.0)) * SCIPdivesetGetMaxLPIterQuot(diveset) * nlpiterations);
296 /* these constraint handlers are required for polishing an LP relaxation solution beyond rounding */
300 SCIP_CALL( SCIPgetLPBranchCands(scip, &lpcands, &lpcandssol, &lpcandsfrac, &nlpcands, NULL, NULL) );
320 searchubbound = SCIPgetLowerbound(scip) + ubquot * (SCIPgetCutoffbound(scip) - SCIPgetLowerbound(scip));
325 searchavgbound = SCIPgetLowerbound(scip) + avgquot * (SCIPgetAvgLowerbound(scip) - SCIPgetLowerbound(scip));
349 SCIPdebugMessage("(node %" SCIP_LONGINT_FORMAT ") executing %s heuristic: depth=%d, %d fractionals, dualbound=%g, avgbound=%g, cutoffbound=%g, searchbound=%g\n",
350 SCIPgetNNodes(scip), SCIPheurGetName(heur), SCIPgetDepth(scip), nlpcands, SCIPgetDualbound(scip), SCIPgetAvgDualbound(scip),
354 /* allocate buffer storage for previous candidates and their branching values for pseudo cost updates */
390 * dive as long we are in the given objective, depth and iteration limits and fractional variables exist, but
392 * - if the number of fractional variables decreased at least with 1 variable per 2 dive depths, we continue diving
397 || (SCIPgetProbingDepth(scip) < maxdivedepth && SCIPdivesetGetNLPIterations(diveset) < maxnlpiterations && SCIPgetLPObjval(scip) < searchbound))
448 SCIPdebugMessage("%s found roundable primal solution: obj=%g\n", SCIPdivesetGetName(diveset), SCIPgetSolOrigObj(scip, worksol));
481 /* in case we do not solve LP's at every probing node, we must explicitly store the solution values by unlinking the
482 * solution, otherwise, the working solution may contain wrong entries, if, e.g., a backtrack occurred after an
507 /* select the next diving action by selecting appropriate dive bound changes for the preferred and alternative child */
508 SCIP_CALL( selectNextDiving(scip, diveset, worksol, onlylpbranchcands, SCIPgetProbingDepth(scip) == lastlpdepth,
512 /* if we did not succeed finding an enforcement, the solution is potentially feasible and we break immediately */
531 /* ensure that a new candidate was successfully determined (usually at the end of the previous loop iteration) */
547 SCIPgetDiveBoundChangeData(scip, &bdchgvars, &bdchgdirs, &bdchgvals, &nbdchanges, !backtracked);
573 /* if the variable is already fixed or if the solution value is outside the domain, numerical troubles may have
583 SCIPdebugMessage("Selected variable <%s> already fixed to [%g,%g] (solval: %.9f), diving aborted \n",
591 SCIPdebugMessage("selected variable's <%s> solution value is outside the domain [%g,%g] (solval: %.9f), diving aborted\n",
597 SCIPdebugMessage(" dive %d/%d, LP iter %" SCIP_LONGINT_FORMAT "/%" SCIP_LONGINT_FORMAT ": var <%s>, sol=%g, oldbounds=[%g,%g],",
598 SCIPgetProbingDepth(scip), maxdivedepth, SCIPdivesetGetNLPIterations(diveset), maxnlpiterations,
623 SCIPerrorMessage("Error: Unsupported bound change direction <%d> specified for diving, aborting\n",bdchgdirs[d]);
639 /* add the number of bound changes we applied by ourselves after propagation, otherwise the counter would have been reset */
642 /* resolve the diving LP if the diving resolve frequency is reached or a sufficient number of intermediate bound changes
648 || (onlylpbranchcands && nviollpcands > (int)(SCIPdivesetGetLPResolveDomChgQuot(diveset) * nlpcands))) )
663 SCIPdebugMessage(" *** cutoff detected at level %d - backtracking\n", SCIPgetProbingDepth(scip));
678 /* store candidate for pseudo cost update and choose next candidate only if no cutoff was detected */
681 if( nbdchanges == 1 && (bdchgdir == SCIP_BRANCHDIR_UPWARDS || bdchgdir == SCIP_BRANCHDIR_DOWNWARDS) )
710 SCIP_CALL( selectNextDiving(scip, diveset, worksol, onlylpbranchcands, SCIPgetProbingDepth(scip) == lastlpdepth,
734 assert(cutoff || (SCIPgetLPSolstat(scip) != SCIP_LPSOLSTAT_OBJLIMIT && SCIPgetLPSolstat(scip) != SCIP_LPSOLSTAT_INFEASIBLE &&
735 (SCIPgetLPSolstat(scip) != SCIP_LPSOLSTAT_OPTIMAL || SCIPisLT(scip, SCIPgetLPObjval(scip), SCIPgetCutoffbound(scip)))));
758 /* todo: should the gain be shared with a smaller weight, instead of dividing the gain itself? */
759 /* it may happen that a variable had an integral solution value beforehand, e.g., for indicator variables */
768 SCIPdebugMessage(" -> lpsolstat=%d, objval=%g/%g, nfrac=%d\n", SCIPgetLPSolstat(scip), SCIPgetLPObjval(scip), searchbound, nlpcands);
777 SCIPdebugMessage("%s found primal solution: obj=%g\n", SCIPdivesetGetName(diveset), SCIPgetSolOrigObj(scip, worksol));
790 SCIPupdateDivesetStats(scip, diveset, totalnprobingnodes, totalnbacktracks, SCIPgetNSolsFound(scip) - oldnsolsfound,
793 SCIPdebugMessage("(node %" SCIP_LONGINT_FORMAT ") finished %s heuristic: %d fractionals, dive %d/%d, LP iter %" SCIP_LONGINT_FORMAT "/%" SCIP_LONGINT_FORMAT ", objval=%g/%g, lpsolstat=%d, cutoff=%u\n",
794 SCIPgetNNodes(scip), SCIPdivesetGetName(diveset), nlpcands, SCIPgetProbingDepth(scip), maxdivedepth, SCIPdivesetGetNLPIterations(diveset), maxnlpiterations,
795 SCIPretransformObj(scip, SCIPgetLPObjval(scip)), SCIPretransformObj(scip, searchbound), SCIPgetLPSolstat(scip), cutoff);
SCIP_RETCODE SCIPsolveProbingLP(SCIP *scip, int itlim, SCIP_Bool *lperror, SCIP_Bool *cutoff) Definition: scip.c:32398 SCIP_RETCODE SCIPgetLPBranchCands(SCIP *scip, SCIP_VAR ***lpcands, SCIP_Real **lpcandssol, SCIP_Real **lpcandsfrac, int *nlpcands, int *npriolpcands, int *nfracimplvars) Definition: scip.c:32735 SCIP_RETCODE SCIPaddDiveBoundChange(SCIP *scip, SCIP_VAR *var, SCIP_BRANCHDIR dir, SCIP_Real value, SCIP_Bool preferred) Definition: scip.c:32649 Definition: type_result.h:33 SCIP_CONSHDLR * SCIPfindConshdlr(SCIP *scip, const char *name) Definition: scip.c:5847 Definition: type_result.h:47 void SCIPupdateDivesetLPStats(SCIP *scip, SCIP_DIVESET *diveset, SCIP_Longint niterstoadd) Definition: scip.c:32535 SCIP_RETCODE SCIPmakeIndicatorsFeasible(SCIP *scip, SCIP_CONSHDLR *conshdlr, SCIP_SOL *sol, SCIP_Bool *changed) Definition: cons_indicator.c:7606 SCIP_RETCODE SCIPmakeSOS1sFeasible(SCIP *scip, SCIP_CONSHDLR *conshdlr, SCIP_SOL *sol, SCIP_Bool *changed, SCIP_Bool *success) Definition: cons_sos1.c:9923 Definition: type_result.h:34 SCIP_Bool SCIPdivesetUseBacktrack(SCIP_DIVESET *diveset) Definition: heur.c:522 Definition: struct_scip.h:53 SCIP_RETCODE SCIPbacktrackProbing(SCIP *scip, int probingdepth) Definition: scip.c:31860 void SCIPwarningMessage(SCIP *scip, const char *formatstr,...) Definition: scip.c:1220 SCIP_Bool SCIPisFeasLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2) Definition: scip.c:41528 library methods for diving heuristics Definition: struct_var.h:196 constraint handler for indicator constraints static SCIP_RETCODE solveLP(SCIP *scip, SCIP_DIVESET *diveset, SCIP_Longint maxnlpiterations, SCIP_Bool *lperror, SCIP_Bool *cutoff) Definition: dive.c:33 SCIP_Real SCIPdivesetGetLPResolveDomChgQuot(SCIP_DIVESET *diveset) Definition: heur.c:540 SCIP_RETCODE SCIPpropagateProbing(SCIP *scip, int maxproprounds, SCIP_Bool *cutoff, SCIP_Longint *ndomredsfound) Definition: scip.c:32162 SCIP_RETCODE SCIPupdateVarPseudocost(SCIP *scip, SCIP_VAR *var, SCIP_Real solvaldelta, SCIP_Real objdelta, SCIP_Real weight) Definition: scip.c:23171 SCIP_Real SCIPgetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var) Definition: scip.c:34593 Definition: type_history.h:37 Definition: type_lp.h:37 void SCIPupdateDivesetStats(SCIP *scip, SCIP_DIVESET *diveset, int nprobingnodes, int nbacktracks, SCIP_Longint nsolsfound, SCIP_Longint nbestsolsfound, SCIP_Bool leavewassol) Definition: scip.c:32548 Definition: struct_sol.h:50 SCIP_Real SCIPdivesetGetMaxLPIterQuot(SCIP_DIVESET *diveset) Definition: heur.c:475 Definition: type_result.h:35 Definition: struct_cons.h:116 SCIP_Bool SCIPisFeasEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2) Definition: scip.c:41515 Definition: type_history.h:36 Definition: type_retcode.h:33 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) Definition: scip.c:32511 SCIP_Bool SCIPisFeasGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2) Definition: scip.c:41554 Definition: struct_heur.h:75 SCIP_RETCODE SCIPgetDiveBoundChanges(SCIP *scip, SCIP_DIVESET *diveset, SCIP_SOL *sol, SCIP_Bool *success, SCIP_Bool *infeasible) Definition: scip.c:32588 Definition: type_lp.h:34 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: dive.c:76 void SCIPgetDiveBoundChangeData(SCIP *scip, SCIP_VAR ***variables, SCIP_BRANCHDIR **directions, SCIP_Real **values, int *ndivebdchgs, SCIP_Bool preferred) Definition: scip.c:32675 SCIP_Real SCIPdivesetGetUbQuotNoSol(SCIP_DIVESET *diveset) Definition: heur.c:491 SCIP_Real SCIPdivesetGetMaxRelDepth(SCIP_DIVESET *diveset) Definition: heur.c:349 Definition: struct_heur.h:36 SCIP_Real SCIPdivesetGetMinRelDepth(SCIP_DIVESET *diveset) Definition: heur.c:341 SCIP_Longint SCIPdivesetGetNLPIterations(SCIP_DIVESET *diveset) Definition: heur.c:445 void SCIPdivesetSetWorkSolution(SCIP_DIVESET *diveset, SCIP_SOL *sol) Definition: heur.c:320 Definition: type_history.h:34 SCIP_RETCODE SCIPperformGenericDivingAlgorithm(SCIP *scip, SCIP_DIVESET *diveset, SCIP_SOL *worksol, SCIP_HEUR *heur, SCIP_RESULT *result, SCIP_Bool nodeinfeasible) Definition: dive.c:188 Definition: type_history.h:35 SCIP_RETCODE SCIPchgVarUbProbing(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound) Definition: scip.c:31962 SCIP_Bool SCIPisFeasLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2) Definition: scip.c:41541 SCIP_RETCODE SCIPchgVarLbProbing(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound) Definition: scip.c:31928 constraint handler for SOS type 1 constraints SCIP_Bool SCIPisFeasGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2) Definition: scip.c:41567 SCIP_Real SCIPdivesetGetAvgQuotNoSol(SCIP_DIVESET *diveset) Definition: heur.c:499 int SCIPdivesetGetMaxLPIterOffset(SCIP_DIVESET *diveset) Definition: heur.c:483 Definition: type_lp.h:35 Definition: type_lp.h:33 public methods for primal heuristics Definition: type_retcode.h:43 SCIP_Longint SCIPdivesetGetSolSuccess(SCIP_DIVESET *diveset) Definition: heur.c:357 SCIP_RETCODE SCIPtrySol(SCIP *scip, SCIP_SOL *sol, SCIP_Bool printreason, SCIP_Bool checkbounds, SCIP_Bool checkintegrality, SCIP_Bool checklprows, SCIP_Bool *stored) Definition: scip.c:35827 SCIP_Bool SCIPdivesetUseOnlyLPBranchcands(SCIP_DIVESET *diveset) Definition: heur.c:552 SCIP_RETCODE SCIProundSol(SCIP *scip, SCIP_SOL *sol, SCIP_Bool *success) Definition: scip.c:35519 |