Scippy

SCIP

Solving Constraint Integer Programs

branch_relpscost.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-2024 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 branch_relpscost.c
26 * @ingroup DEFPLUGINS_BRANCH
27 * @brief reliable pseudo costs branching rule
28 * @author Tobias Achterberg
29 * @author Timo Berthold
30 * @author Gerald Gamrath
31 * @author Marc Pfetsch
32 */
33
34/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
35
38#include "scip/treemodel.h"
39#include "scip/cons_and.h"
40#include "scip/pub_branch.h"
41#include "scip/pub_cons.h"
42#include "scip/pub_message.h"
43#include "scip/pub_misc.h"
44#include "scip/pub_sol.h"
45#include "scip/pub_tree.h"
46#include "scip/pub_var.h"
47#include "scip/scip_branch.h"
48#include "scip/scip_cons.h"
49#include "scip/scip_general.h"
50#include "scip/scip_lp.h"
51#include "scip/scip_mem.h"
52#include "scip/scip_message.h"
53#include "scip/scip_nlp.h"
54#include "scip/scip_numerics.h"
55#include "scip/scip_param.h"
56#include "scip/scip_prob.h"
58#include "scip/scip_sol.h"
60#include "scip/scip_tree.h"
61#include "scip/scip_var.h"
62#include "scip/prop_symmetry.h"
63#include "scip/symmetry.h"
64#include <string.h>
65
66#define BRANCHRULE_NAME "relpscost"
67#define BRANCHRULE_DESC "reliability branching on pseudo cost values"
68#define BRANCHRULE_PRIORITY 10000
69#define BRANCHRULE_MAXDEPTH -1
70#define BRANCHRULE_MAXBOUNDDIST 1.0
71
72#define DEFAULT_CONFLICTWEIGHT 0.01 /**< weight in score calculations for conflict score */
73#define DEFAULT_CONFLENGTHWEIGHT 0.0 /**< weight in score calculations for conflict length score*/
74#define DEFAULT_INFERENCEWEIGHT 0.0001 /**< weight in score calculations for inference score */
75#define DEFAULT_CUTOFFWEIGHT 0.0001 /**< weight in score calculations for cutoff score */
76#define DEFAULT_GMIAVGEFFWEIGHT 0.0 /**< weight in score calculations of average GMI cut normed efficacies */
77#define DEFAULT_GMILASTEFFWEIGHT 0.00001 /**< weight in score calculations of last GMI cut normed efficacy */
78#define DEFAULT_PSCOSTWEIGHT 1.0 /**< weight in score calculations for pseudo cost score */
79#define DEFAULT_NLSCOREWEIGHT 0.1 /**< weight in score calculations for nlcount score */
80#define DEFAULT_MINRELIABLE 1.0 /**< minimal value for minimum pseudo cost size to regard pseudo cost value as reliable */
81#define DEFAULT_MAXRELIABLE 5.0 /**< maximal value for minimum pseudo cost size to regard pseudo cost value as reliable */
82#define DEFAULT_SBITERQUOT 0.5 /**< maximal fraction of strong branching LP iterations compared to normal iterations */
83#define DEFAULT_SBITEROFS 100000 /**< additional number of allowed strong branching LP iterations */
84#define DEFAULT_MAXLOOKAHEAD 9 /**< maximal number of further variables evaluated without better score */
85#define DEFAULT_INITCAND 100 /**< maximal number of candidates initialized with strong branching per node */
86#define DEFAULT_INITITER 0 /**< iteration limit for strong branching initialization of pseudo cost entries (0: auto) */
87#define DEFAULT_MAXBDCHGS 5 /**< maximal number of bound tightenings before the node is reevaluated (-1: unlimited) */
88#define DEFAULT_MAXPROPROUNDS -2 /**< maximum number of propagation rounds to be performed during strong branching
89 * before solving the LP (-1: no limit, -2: parameter settings) */
90#define DEFAULT_PROBINGBOUNDS TRUE /**< should valid bounds be identified in a probing-like fashion during strong
91 * branching (only with propagation)? */
92#define DEFAULT_USERELERRORFORRELIABILITY FALSE /**< should reliability be based on relative errors? */
93#define DEFAULT_LOWERRORTOL 0.05 /**< lowest tolerance beneath which relative errors are reliable */
94#define DEFAULT_HIGHERRORTOL 1.0 /**< highest tolerance beneath which relative errors are reliable */
95#define DEFAULT_USEHYPTESTFORRELIABILITY FALSE /**< should the strong branching decision be based on a hypothesis test? */
96#define DEFAULT_USEDYNAMICCONFIDENCE FALSE /**< should the confidence level be adjusted dynamically? */
97#define DEFAULT_STORESEMIINITCOSTS FALSE /**< should strong branching result be considered for pseudo costs if the other direction was infeasible? */
98#define DEFAULT_USESBLOCALINFO FALSE /**< should the scoring function use only local cutoff and inference information obtained for strong branching candidates? */
99#define DEFAULT_CONFIDENCELEVEL 2 /**< The confidence level for statistical methods, between 0 (Min) and 4 (Max). */
100#define DEFAULT_SKIPBADINITCANDS TRUE /**< should branching rule skip candidates that have a low probability to be
101 * better than the best strong-branching or pseudo-candidate? */
102#define DEFAULT_STARTRANDSEED 5 /**< start random seed for random number generation */
103#define DEFAULT_RANDINITORDER FALSE /**< should slight perturbation of scores be used to break ties in the prior scores? */
104#define DEFAULT_USESMALLWEIGHTSITLIM FALSE /**< should smaller weights be used for pseudo cost updates after hitting the LP iteration limit? */
105#define DEFAULT_DYNAMICWEIGHTS TRUE /**< should the weights of the branching rule be adjusted dynamically during solving based
106 * infeasible and objective leaf counters? */
107#define DEFAULT_DEGENERACYAWARE 1 /**< should degeneracy be taken into account to update weights and skip strong branching? (0: off, 1: after root, 2: always)*/
108
109/* symmetry handling */
110#define DEFAULT_FILTERCANDSSYM FALSE /**< Use symmetry to filter branching candidates? */
111#define DEFAULT_TRANSSYMPSCOST FALSE /**< Transfer pscost information to symmetric variables if filtering is performed? */
112
113/** branching rule data */
114struct SCIP_BranchruleData
115{
116 SCIP_Real conflictweight; /**< weight in score calculations for conflict score */
117 SCIP_Real conflengthweight; /**< weight in score calculations for conflict length score */
118 SCIP_Real inferenceweight; /**< weight in score calculations for inference score */
119 SCIP_Real cutoffweight; /**< weight in score calculations for cutoff score */
120 SCIP_Real gmiavgeffweight; /**< weight in score calculations of average GMI normed cut efficacies */
121 SCIP_Real gmilasteffweight; /**< weight in score calculations of last GMI cut normalized efficacy */
122 SCIP_Real pscostweight; /**< weight in score calculations for pseudo cost score */
123 SCIP_Real nlscoreweight; /**< weight in score calculations for nlcount score */
124 SCIP_Real minreliable; /**< minimal value for minimum pseudo cost size to regard pseudo cost value as reliable */
125 SCIP_Real maxreliable; /**< maximal value for minimum pseudo cost size to regard pseudo cost value as reliable */
126 SCIP_Real sbiterquot; /**< maximal fraction of strong branching LP iterations compared to normal iterations */
127 int sbiterofs; /**< additional number of allowed strong branching LP iterations */
128 int maxlookahead; /**< maximal number of further variables evaluated without better score */
129 int initcand; /**< maximal number of candidates initialized with strong branching per node */
130 int inititer; /**< iteration limit for strong branching initialization of pseudo cost entries (0: auto) */
131 int maxbdchgs; /**< maximal number of bound tightenings before the node is reevaluated (-1: unlimited) */
132 int maxproprounds; /**< maximum number of propagation rounds to be performed during strong branching
133 * before solving the LP (-1: no limit, -2: parameter settings) */
134 SCIP_Bool probingbounds; /**< should valid bounds be identified in a probing-like fashion during strong
135 * branching (only with propagation)? */
136 SCIP_Bool userelerrorforreliability; /**< should reliability be based on relative errors? */
137 SCIP_Real lowerrortol; /**< lowest tolerance beneath which relative errors are reliable */
138 SCIP_Real higherrortol; /**< highest tolerance beneath which relative errors are reliable */
139 SCIP_Bool usehyptestforreliability; /**< should the strong branching decision be based on a hypothesis test? */
140 SCIP_Bool usedynamicconfidence; /**< should the confidence level be adjusted dynamically? */
141 SCIP_Bool storesemiinitcosts; /**< should strong branching result be considered for pseudo costs if the
142 * other direction was infeasible? */
143 SCIP_Bool usesblocalinfo; /**< should the scoring function disregard cutoffs for variable if sb-lookahead was feasible ? */
144 SCIP_Bool skipbadinitcands; /**< should branching rule skip candidates that have a low probability to be
145 * better than the best strong-branching or pseudo-candidate? */
146 SCIP_Bool dynamicweights; /**< should the weights of the branching rule be adjusted dynamically during
147 * solving based on objective and infeasible leaf counters? */
148 int degeneracyaware; /**< should degeneracy be taken into account to update weights and skip strong branching? (0: off, 1: after root, 2: always) */
149 int confidencelevel; /**< The confidence level for statistical methods, between 0 (Min) and 4 (Max). */
150 int* nlcount; /**< array to store nonlinear count values */
151 int nlcountsize; /**< length of nlcount array */
152 int nlcountmax; /**< maximum entry in nlcount array or 1 if NULL */
153 SCIP_Bool randinitorder; /**< should slight perturbation of scores be used to break ties in the prior scores? */
154 SCIP_RANDNUMGEN* randnumgen; /**< random number generator */
155 int startrandseed; /**< start random seed for random number generation */
156 SCIP_Bool usesmallweightsitlim; /**< should smaller weights be used for pseudo cost updates after hitting the LP iteration limit? */
157 SCIP_TREEMODEL* treemodel; /**< Parameters for the Treemodel branching rules */
158
159 /* for symmetry */
160 SCIP_Bool filtercandssym; /**< Use symmetry to filter branching candidates? */
161 SCIP_Bool transsympscost; /**< Transfer pscost information to symmetric variables? */
162
163 SCIP_Bool nosymmetry; /**< No symmetry present? */
164 int* orbits; /**< array of non-trivial orbits */
165 int* orbitbegins; /**< array containing begin positions of new orbits in orbits array */
166 int norbits; /**< pointer to number of orbits currently stored in orbits */
167 int* varorbitmap; /**< array for storing indices of the containing orbit for each variable */
168 int* orbitrep; /**< representative variable of each orbit */
169 SCIP_VAR** permvars; /**< variables on which permutations act */
170 int npermvars; /**< number of variables for permutations */
171 SCIP_HASHMAP* permvarmap; /**< map of variables to indices in permvars array */
172};
173
174/*
175 * local methods
176 */
177
178/** initialize orbits */
179static
181 SCIP* scip, /**< SCIP data structure */
182 SCIP_BRANCHRULEDATA* branchruledata /**< branching rule data */
183 )
184{
185 int** permstrans = NULL;
186 int* components = NULL;
187 int* componentbegins = NULL;
188 int* vartocomponent = NULL;
189 int ncomponents = 0;
190 int nperms = -1;
191
192 assert( scip != NULL );
193 assert( branchruledata != NULL );
194 assert( branchruledata->filtercandssym );
195
196 /* exit if no symmetry or orbits already available */
197 if( branchruledata->nosymmetry || branchruledata->orbits != NULL )
198 return SCIP_OKAY;
199
200 assert( branchruledata->orbitbegins == NULL );
201 assert( branchruledata->varorbitmap == NULL );
202 assert( branchruledata->orbitrep == NULL );
203
204 /* obtain symmetry including permutations */
205 SCIP_CALL( SCIPgetSymmetry(scip, &branchruledata->npermvars, &branchruledata->permvars, &branchruledata->permvarmap,
206 &nperms, NULL, &permstrans, NULL, NULL, &components, &componentbegins, &vartocomponent, &ncomponents) );
207
208 /* turn off symmetry handling if there is no symmetry or the number of variables is not equal */
209 if( nperms <= 0 || branchruledata->npermvars != SCIPgetNVars(scip) )
210 {
211 branchruledata->nosymmetry = TRUE;
212 return SCIP_OKAY;
213 }
214 assert( branchruledata->permvars != NULL );
215 assert( branchruledata->permvarmap != NULL );
216 assert( branchruledata->npermvars > 0 );
217 assert( permstrans != NULL );
218 assert( components != NULL );
219 assert( componentbegins != NULL );
220 assert( vartocomponent != NULL );
221 assert( ncomponents > 0 );
222
223 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &branchruledata->orbits, branchruledata->npermvars) );
224 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &branchruledata->orbitbegins, branchruledata->npermvars) );
225 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &branchruledata->varorbitmap, branchruledata->npermvars) );
226 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &branchruledata->orbitrep, branchruledata->npermvars) );
227
228 /* Compute orbits on all variables, since this might help for branching and this computation is only done once. */
229 SCIP_CALL( SCIPcomputeOrbitsComponentsSym(scip, branchruledata->npermvars, permstrans, nperms,
230 components, componentbegins, vartocomponent, ncomponents,
231 branchruledata->orbits, branchruledata->orbitbegins, &branchruledata->norbits, branchruledata->varorbitmap) );
232 assert( branchruledata->norbits < branchruledata->npermvars );
233
234 return SCIP_OKAY;
235}
236
237/** filter out symmetric variables from branching variables */
238static
240 SCIP* scip, /**< SCIP data structure */
241 SCIP_BRANCHRULEDATA* branchruledata, /**< branching rule data */
242 SCIP_VAR** origbranchcands, /**< original branching candidates */
243 SCIP_Real* origbranchcandssol, /**< original solution value for the branching candidates */
244 SCIP_Real* origbranchcandsfrac,/**< original fractional part of the branching candidates */
245 int norigbranchcands, /**< original number of branching candidates */
246 SCIP_VAR** branchcands, /**< branching candidates */
247 SCIP_Real* branchcandssol, /**< solution value for the branching candidates */
248 SCIP_Real* branchcandsfrac, /**< fractional part of the branching candidates */
249 int* branchorbitidx, /**< array of indices of orbit of branching candidates */
250 int* nbranchcands /**< pointer to store number of branching candidates */
251 )
252{
253 int i;
254
255 assert( scip != NULL );
256 assert( branchruledata != NULL );
257 assert( origbranchcands != NULL );
258 assert( origbranchcandssol != NULL );
259 assert( origbranchcandsfrac != NULL );
260 assert( branchcands != NULL );
261 assert( branchcandssol != NULL );
262 assert( branchcandsfrac != NULL );
263 assert( nbranchcands != NULL );
264
265 assert( ! branchruledata->nosymmetry );
266 assert( branchruledata->orbitbegins != NULL );
267 assert( branchruledata->orbits != NULL );
268 assert( branchruledata->permvarmap != NULL );
269 assert( branchruledata->varorbitmap != NULL );
270 assert( branchruledata->orbitrep != NULL );
271 assert( branchruledata->norbits < branchruledata->npermvars );
272
273 /* init representatives (used to see whether variable is the first in an orbit) */
274 for( i = 0; i < branchruledata->norbits; ++i )
275 branchruledata->orbitrep[i] = -1;
276
277 /* loop through branching variables, determine orbit and whether they are the first ones */
278 *nbranchcands = 0;
279 for( i = 0; i < norigbranchcands; ++i )
280 {
281 int orbitidx = -1;
282 int varidx;
283
284 varidx = SCIPhashmapGetImageInt(branchruledata->permvarmap, (void*) origbranchcands[i]);
285 if( varidx != INT_MAX )
286 {
287 assert( 0 <= varidx && varidx < branchruledata->npermvars );
288 orbitidx = branchruledata->varorbitmap[varidx];
289 }
290 assert( -1 <= orbitidx && orbitidx < branchruledata->norbits );
291
292 /* Check whether the variable is not present (can happen if variable was added after computing symmetries or is in
293 * a singleton orbit). */
294 if( orbitidx == -1 )
295 {
296 branchcands[*nbranchcands] = origbranchcands[i];
297 branchcandssol[*nbranchcands] = origbranchcandssol[i];
298 branchcandsfrac[*nbranchcands] = origbranchcandsfrac[i];
299 branchorbitidx[*nbranchcands] = -1;
300 ++(*nbranchcands);
301 }
302 else if( branchruledata->orbitrep[orbitidx] == -1 )
303 {
304 /* if variable is the first in a nontrivial orbit */
305 assert( 0 <= varidx && varidx < branchruledata->npermvars );
306 branchruledata->orbitrep[orbitidx] = varidx;
307 branchcands[*nbranchcands] = origbranchcands[i];
308 branchcandssol[*nbranchcands] = origbranchcandssol[i];
309 branchcandsfrac[*nbranchcands] = origbranchcandsfrac[i];
310 branchorbitidx[*nbranchcands] = orbitidx;
311 ++(*nbranchcands);
312 }
313 }
314
315 SCIPdebugMsg(scip, "Filtered out %d variables by symmetry.\n", norigbranchcands - *nbranchcands);
316
317 return SCIP_OKAY;
318}
319
320/** updates the pseudo costs of the given variable and all its symmetric variables */
321static
323 SCIP* scip, /**< SCIP data structure */
324 SCIP_BRANCHRULEDATA* branchruledata, /**< branching rule data */
325 SCIP_VAR* branchvar, /**< branching variable candidate */
326 int* branchorbitidx, /**< array of orbit indices */
327 int branchvaridx, /**< index of variable in branchorbitidx */
328 SCIP_Real solvaldelta, /**< difference of variable's new LP value - old LP value */
329 SCIP_Real objdelta, /**< difference of new LP's objective value - old LP's objective value */
330 SCIP_Real weight /**< weight in (0,1] of this update in pseudo cost sum */
331 )
332{
333 int orbitidx;
334 int j;
335
336 assert( scip != NULL );
337 assert( branchruledata != NULL );
338
339 if( branchruledata->nosymmetry || ! branchruledata->transsympscost || branchorbitidx == NULL )
340 {
341 /* use original update function */
342 SCIP_CALL( SCIPupdateVarPseudocost(scip, branchvar, solvaldelta, objdelta, weight) );
343 return SCIP_OKAY;
344 }
345
346 assert( branchruledata->orbitbegins != NULL );
347 assert( branchruledata->orbits != NULL );
348 assert( 0 <= branchvaridx && branchvaridx < branchruledata->npermvars );
349
350 orbitidx = branchorbitidx[branchvaridx];
351 if( orbitidx < 0 )
352 {
353 /* only update given variable */
354 SCIP_CALL( SCIPupdateVarPseudocost(scip, branchvar, solvaldelta, objdelta, weight) );
355 return SCIP_OKAY;
356 }
357 assert( 0 <= orbitidx && orbitidx < branchruledata->norbits );
358
359 /* loop through orbit containing variable and update pseudo costs for all variables */
360 for( j = branchruledata->orbitbegins[orbitidx]; j < branchruledata->orbitbegins[orbitidx+1]; ++j )
361 {
362 SCIP_VAR* var;
363 int idx;
364
365 idx = branchruledata->orbits[j];
366 assert( 0 <= idx && idx < branchruledata->npermvars );
367
368 var = branchruledata->permvars[idx];
369 assert( var != NULL );
370
371 if( SCIPvarIsActive(var) )
372 {
373 SCIP_CALL( SCIPupdateVarPseudocost(scip, var, solvaldelta, objdelta, weight) );
374 }
375 }
376
377 return SCIP_OKAY;
378}
379
380/**! [SnippetCodeStyleStaticAsserts] */
381
382/** return probindex of variable or corresponding active variable (if negated or aggregated) or -1 (if multiaggregated) */
383static
385 SCIP* scip, /**< SCIP data structure */
386 SCIP_VAR* var, /**< binary variable */
387 int* probindex /**< buffer to store probindex */
388 )
389{
390 assert(scip != NULL);
391 assert(var != NULL);
392 assert(SCIPvarIsBinary(var));
393 assert(probindex != NULL);
394
395/**! [SnippetCodeStyleStaticAsserts] */
396
397 *probindex = SCIPvarGetProbindex(var);
398
399 /* if variable is not active, try to find active representative */
400 if( *probindex == -1 )
401 {
402 SCIP_VAR* repvar;
403 SCIP_Bool negated;
404
405 SCIP_CALL( SCIPgetBinvarRepresentative(scip, var, &repvar, &negated) );
406 assert(repvar != NULL);
407 assert(SCIPvarGetStatus(repvar) != SCIP_VARSTATUS_FIXED);
408
409 if( SCIPvarIsActive(repvar) )
410 *probindex = SCIPvarGetProbindex(repvar);
411 else if( SCIPvarIsNegated(repvar) )
412 *probindex = SCIPvarGetProbindex(SCIPvarGetNegationVar(repvar));
413 }
414
415 return SCIP_OKAY;
416}
417
418/**! [SnippetCodeStyleDeclaration] */
419
420/** counts number of nonlinear constraints in which each variable appears */
421static
423 SCIP* scip, /**< SCIP data structure */
424 int* nlcount, /**< pointer to array for storing count values */
425 int nlcountsize, /**< buffer for storing length of nlcount array */
426 int* nlcountmax /**< buffer for storing maximum value in nlcount array */
427 )
428{
429 SCIP_CONSHDLR* andconshdlr;
430 SCIP_VAR** vars;
431 int nvars;
432 int i;
433
434/**! [SnippetCodeStyleDeclaration] */
435
436 assert(scip != NULL);
437 assert(nlcount != NULL);
438 assert(nlcountmax != NULL);
439
440 SCIP_CALL( SCIPgetVarsData(scip, &vars, &nvars, NULL, NULL, NULL, NULL) );
441 assert(nlcountsize >= nvars);
442
443 /* get nonlinearity for constraints in NLP */
445 {
446 assert(SCIPgetNNLPVars(scip) == nvars);
448 }
449 else
450 {
451 BMSclearMemoryArray(nlcount, nvars);
452 }
453
454 /* increase counters for and constraints */
455 andconshdlr = SCIPfindConshdlr(scip, "and");
456 if( andconshdlr != NULL )
457 {
458 int c;
459
460 for( c = 0; c < SCIPconshdlrGetNActiveConss(andconshdlr); c++ )
461 {
462 SCIP_CONS* andcons;
463 SCIP_VAR** andvars;
464 SCIP_VAR* andres;
465 int probindex;
466 int nandvars;
467 int v;
468
469 /* get constraint and variables */
470 andcons = SCIPconshdlrGetConss(andconshdlr)[c];
471 nandvars = SCIPgetNVarsAnd(scip, andcons);
472 andvars = SCIPgetVarsAnd(scip, andcons);
473 andres = SCIPgetResultantAnd(scip, andcons);
474
475 probindex = -1;
476
477 /**! [SnippetCodeStyleIfFor] */
478
479 for( v = 0; v < nandvars; v++ )
480 {
481 /* don't rely on the and conshdlr removing fixed variables
482 * @todo fix the and conshdlr in that respect
483 */
484 if( SCIPvarGetStatus(andvars[v]) != SCIP_VARSTATUS_FIXED )
485 {
486 SCIP_CALL( binvarGetActiveProbindex(scip, andvars[v], &probindex) );
487 if( probindex >= 0 )
488 nlcount[probindex]++;
489 }
490 }
491
492 /**! [SnippetCodeStyleIfFor] */
493
494 SCIP_CALL( binvarGetActiveProbindex(scip, andres, &probindex) );
495 if( probindex >= 0 )
496 nlcount[probindex]++;
497 }
498 }
499
500 /* compute maximum count value */
501 *nlcountmax = 1;
502 for( i = 0; i < nvars; i++ )
503 {
504 if( *nlcountmax < nlcount[i] )
505 *nlcountmax = nlcount[i];
506 }
507
508 return SCIP_OKAY;
509}
510
511static
513 SCIP* scip, /**< SCIP data structure */
514 SCIP_BRANCHRULEDATA* branchruledata /**< branching rule data */
515 )
516{
517 int nvars;
518
519 assert(scip != NULL);
520 assert(branchruledata != NULL);
521
522 nvars = SCIPgetNVars(scip);
523
524 /**@todo test whether we want to apply this as if problem has only and constraints */
525 /**@todo update changes in and constraints */
526 if( branchruledata->nlscoreweight > 0.0 ) /* && SCIPisNLPConstructed(scip) */
527 {
528 if( branchruledata->nlcount == NULL )
529 {
530 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &branchruledata->nlcount, nvars) );
531 branchruledata->nlcountsize = nvars;
532
533 SCIP_CALL( countNonlinearities(scip, branchruledata->nlcount, branchruledata->nlcountsize, &branchruledata->nlcountmax) );
534 }
535 else if( branchruledata->nlcountsize < nvars )
536 {
537 SCIP_CALL( SCIPreallocBlockMemoryArray(scip, &branchruledata->nlcount, branchruledata->nlcountsize, nvars) );
538 /**@todo should we update nlcounts for new variables? */
539 BMSclearMemoryArray(&(branchruledata->nlcount[branchruledata->nlcountsize]), nvars - branchruledata->nlcountsize); /*lint !e866*/
540 branchruledata->nlcountsize = nvars;
541 }
542 assert(branchruledata->nlcount != NULL);
543 assert(branchruledata->nlcountsize == nvars);
544 assert(branchruledata->nlcountmax >= 1);
545 }
546 else
547 {
548 SCIPfreeBlockMemoryArrayNull(scip, &branchruledata->nlcount, branchruledata->nlcountsize);
549 branchruledata->nlcountsize = 0;
550 branchruledata->nlcountmax = 1;
551 }
552
553 return SCIP_OKAY;
554}
555
556
557/** calculates nlscore value between 0 and 1 */
558static
560 SCIP* scip, /**< SCIP data structure */
561 int* nlcount, /**< array to store count values */
562 int nlcountmax, /**< maximum value in nlcount array */
563 int probindex /**< index of branching candidate */
564 )
565{
566 if( nlcountmax >= 1 && nlcount != NULL )
567 {
568 SCIP_Real nlscore;
569
570 assert(scip != NULL);
571 assert(probindex >= 0);
572 assert(probindex < SCIPgetNVars(scip));
573
574 nlscore = nlcount[probindex] / (SCIP_Real)nlcountmax;
575
576 assert(nlscore >= 0.0);
577 assert(nlscore <= 1.0);
578 return nlscore;
579 }
580 else
581 return 0.0;
582}
583
584/** calculates an overall score value for the given individual score values */
585static
587 SCIP* scip, /**< SCIP data structure */
588 SCIP_BRANCHRULEDATA* branchruledata, /**< branching rule data */
589 SCIP_Real conflictscore, /**< conflict score of current variable */
590 SCIP_Real avgconflictscore, /**< average conflict score */
591 SCIP_Real conflengthscore, /**< conflict length score of current variable */
592 SCIP_Real avgconflengthscore, /**< average conflict length score */
593 SCIP_Real inferencescore, /**< inference score of current variable */
594 SCIP_Real avginferencescore, /**< average inference score */
595 SCIP_Real cutoffscore, /**< cutoff score of current variable */
596 SCIP_Real avgcutoffscore, /**< average cutoff score */
597 SCIP_Real gmieffscore, /**< normalized-eff of avg GMI cuts from row when var was frac and basic */
598 SCIP_Real lastgmieffscore, /**< last normalized gmieffscore when var was frac and basic */
599 SCIP_Real pscostscore, /**< pscost score of current variable */
600 SCIP_Real avgpscostscore, /**< average pscost score */
601 SCIP_Real nlscore, /**< nonlinear score of current variable between 0 and 1 */
602 SCIP_Real frac, /**< fractional value of variable in current solution */
603 SCIP_Real degeneracyfactor /**< factor to apply because of degeneracy */
604 )
605{
606 SCIP_Real score;
607 SCIP_Real dynamicfactor;
608
609 assert(branchruledata != NULL);
610 assert(0.0 < frac && frac < 1.0);
611
612 if( branchruledata->dynamicweights )
613 {
614 dynamicfactor = (SCIPgetNInfeasibleLeaves(scip) + 1.0) / (SCIPgetNObjlimLeaves(scip) + 1.0);
615 }
616 else
617 dynamicfactor = 1.0;
618
619 dynamicfactor *= degeneracyfactor;
620
621 score = dynamicfactor * (branchruledata->conflictweight * (1.0 - 1.0/(1.0+conflictscore/avgconflictscore))
622 + branchruledata->conflengthweight * (1.0 - 1.0/(1.0+conflengthscore/avgconflengthscore))
623 + branchruledata->inferenceweight * (1.0 - 1.0/(1.0+inferencescore/avginferencescore))
624 + branchruledata->cutoffweight * (1.0 - 1.0/(1.0+cutoffscore/avgcutoffscore))
625 + branchruledata->gmiavgeffweight * gmieffscore + branchruledata->gmilasteffweight * lastgmieffscore)
626 + branchruledata->pscostweight / dynamicfactor * (1.0 - 1.0/(1.0+pscostscore/avgpscostscore))
627 + branchruledata->nlscoreweight * nlscore;
628
629 /* avoid close to integral variables */
630 if( MIN(frac, 1.0 - frac) < 10.0 * SCIPfeastol(scip) )
631 score *= 1e-6;
632
633 return score;
634}
635
636/** adds given index and direction to bound change arrays */
637static
639 SCIP* scip, /**< SCIP data structure */
640 int** bdchginds, /**< pointer to bound change index array */
641 SCIP_BOUNDTYPE** bdchgtypes, /**< pointer to bound change types array */
642 SCIP_Real** bdchgbounds, /**< pointer to bound change new bounds array */
643 int* nbdchgs, /**< pointer to number of bound changes */
644 int ind, /**< index to store in bound change index array */
645 SCIP_BOUNDTYPE type, /**< type of the bound change to store in bound change type array */
646 SCIP_Real bound /**< new bound to store in bound change new bounds array */
647 )
648{
649 assert(bdchginds != NULL);
650 assert(bdchgtypes != NULL);
651 assert(bdchgbounds != NULL);
652 assert(nbdchgs != NULL);
653
654 SCIP_CALL( SCIPreallocBufferArray(scip, bdchginds, (*nbdchgs) + 1) );
655 SCIP_CALL( SCIPreallocBufferArray(scip, bdchgtypes, (*nbdchgs) + 1) );
656 SCIP_CALL( SCIPreallocBufferArray(scip, bdchgbounds, (*nbdchgs) + 1) );
657 assert(*bdchginds != NULL);
658 assert(*bdchgtypes != NULL);
659 assert(*bdchgbounds != NULL);
660
661 (*bdchginds)[*nbdchgs] = ind;
662 (*bdchgtypes)[*nbdchgs] = type;
663 (*bdchgbounds)[*nbdchgs] = bound;
664 (*nbdchgs)++;
665
666 return SCIP_OKAY;
667}
668
669/** frees bound change arrays */
670static
672 SCIP* scip, /**< SCIP data structure */
673 int** bdchginds, /**< pointer to bound change index array */
674 SCIP_BOUNDTYPE** bdchgtypes, /**< pointer to bound change types array */
675 SCIP_Real** bdchgbounds, /**< pointer to bound change new bounds array */
676 int* nbdchgs /**< pointer to number of bound changes */
677 )
678{
679 assert(bdchginds != NULL);
680 assert(bdchgtypes != NULL);
681 assert(bdchgbounds != NULL);
682 assert(nbdchgs != NULL);
683
684 SCIPfreeBufferArrayNull(scip, bdchgbounds);
685 SCIPfreeBufferArrayNull(scip, bdchgtypes);
686 SCIPfreeBufferArrayNull(scip, bdchginds);
687 *nbdchgs = 0;
688}
689
690/** applies bound changes stored in bound change arrays */
691static
693 SCIP* scip, /**< SCIP data structure */
694 SCIP_VAR** vars, /**< problem variables */
695 int* bdchginds, /**< bound change index array */
696 SCIP_BOUNDTYPE* bdchgtypes, /**< bound change types array */
697 SCIP_Real* bdchgbounds, /**< bound change new bound array */
698 int nbdchgs, /**< number of bound changes */
699 SCIP_RESULT* result /**< result pointer */
700 )
701{
702#ifndef NDEBUG
703 SCIP_BRANCHRULE* branchrule;
704 SCIP_BRANCHRULEDATA* branchruledata;
705#endif
706 SCIP_Bool infeasible;
707 SCIP_Bool tightened;
708 int i;
709
710 assert(vars != NULL);
711
712#ifndef NDEBUG
713 /* find branching rule */
715 assert(branchrule != NULL);
716
717 /* get branching rule data */
718 branchruledata = SCIPbranchruleGetData(branchrule);
719 assert(branchruledata != NULL);
720#endif
721
722 SCIPdebugMsg(scip, "applying %d bound changes\n", nbdchgs);
723
724 for( i = 0; i < nbdchgs; ++i )
725 {
726 int v;
727
728 v = bdchginds[i];
729
730 SCIPdebugMsg(scip, " -> <%s> [%g,%g]\n",
731 SCIPvarGetName(vars[v]), SCIPvarGetLbLocal(vars[v]), SCIPvarGetUbLocal(vars[v]));
732
733 if( bdchgtypes[i] == SCIP_BOUNDTYPE_LOWER )
734 {
735 /* change lower bound of variable to given bound */
736 SCIP_CALL( SCIPtightenVarLb(scip, vars[v], bdchgbounds[i], TRUE, &infeasible, &tightened) );
737 if( infeasible )
738 {
739 *result = SCIP_CUTOFF;
740 return SCIP_OKAY;
741 }
742
743 /* if we did propagation, the bound change might already have been added */
744 assert(tightened || (branchruledata->maxproprounds != 0));
745 }
746 else
747 {
748 assert(bdchgtypes[i] == SCIP_BOUNDTYPE_UPPER);
749
750 /* change upper bound of variable to given bound */
751 SCIP_CALL( SCIPtightenVarUb(scip, vars[v], bdchgbounds[i], TRUE, &infeasible, &tightened) );
752 if( infeasible )
753 {
754 *result = SCIP_CUTOFF;
755 return SCIP_OKAY;
756 }
757
758 /* if we did propagation, the bound change might already have been added */
759 assert(tightened || (branchruledata->maxproprounds != 0));
760 }
761 SCIPdebugMsg(scip, " -> [%g,%g]\n", SCIPvarGetLbLocal(vars[v]), SCIPvarGetUbLocal(vars[v]));
762 }
763
764 return SCIP_OKAY;
765}
766
767/** execute reliability pseudo cost branching */
768static
770 SCIP* scip, /**< SCIP data structure */
771 SCIP_BRANCHRULE* branchrule, /**< branching rule */
772 SCIP_VAR** branchcands, /**< branching candidates */
773 SCIP_Real* branchcandssol, /**< solution value for the branching candidates */
774 SCIP_Real* branchcandsfrac, /**< fractional part of the branching candidates */
775 int* branchorbitidx, /**< indices of orbit (or NULL) */
776 int nbranchcands, /**< number of branching candidates */
777 SCIP_Bool executebranch, /**< execute a branching step or run probing only */
778 SCIP_RESULT* result /**< pointer to the result of the execution */
779 )
780{ /*lint --e{715}*/
781 SCIP_BRANCHRULEDATA* branchruledata;
782 SCIP_Real lpobjval;
783 SCIP_Real bestsbdown;
784 SCIP_Real bestsbup;
785 SCIP_Real provedbound;
786 SCIP_Bool bestsbdownvalid;
787 SCIP_Bool bestsbupvalid;
788 SCIP_Bool bestsbdowncutoff;
789 SCIP_Bool bestsbupcutoff;
790 SCIP_Bool bestisstrongbranch;
791 SCIP_Bool allcolsinlp;
792 SCIP_Bool exactsolve;
793 int ninitcands;
794 int bestcand;
795
796 /* remember which variables strong branching is performed on, and the
797 * recorded lp bound changes that are observed */
798 SCIP_Bool* sbvars = NULL;
799 SCIP_Real* sbdown = NULL;
800 SCIP_Real* sbup = NULL;
801 SCIP_Bool* sbdownvalid = NULL;
802 SCIP_Bool* sbupvalid = NULL;
803
804 *result = SCIP_DIDNOTRUN;
805
807
808 /* get branching rule data */
809 branchruledata = SCIPbranchruleGetData(branchrule);
810 assert(branchruledata != NULL);
811
812 /* get current LP objective bound of the local sub problem and global cutoff bound */
813 lpobjval = SCIPgetLPObjval(scip);
814
815 /* check, if we want to solve the problem exactly, meaning that strong branching information is not useful
816 * for cutting off sub problems and improving lower bounds of children
817 */
818 exactsolve = SCIPisExactSolve(scip);
819
820 /* check, if all existing columns are in LP, and thus the strong branching results give lower bounds */
821 allcolsinlp = SCIPallColsInLP(scip);
822
823 bestcand = -1;
824 bestisstrongbranch = FALSE;
825 bestsbdown = SCIP_INVALID;
826 bestsbup = SCIP_INVALID;
827 bestsbdownvalid = FALSE;
828 bestsbupvalid = FALSE;
829 bestsbdowncutoff = FALSE;
830 bestsbupcutoff = FALSE;
831 provedbound = lpobjval;
832
833 /* Allocate memory to store the lp bounds of the up and down children
834 * for those of the variables that we performed sb on
835 */
836 SCIP_CALL( SCIPallocBufferArray(scip, &sbdown, nbranchcands) );
837 SCIP_CALL( SCIPallocBufferArray(scip, &sbup, nbranchcands) );
838 SCIP_CALL( SCIPallocBufferArray(scip, &sbdownvalid, nbranchcands) );
839 SCIP_CALL( SCIPallocBufferArray(scip, &sbupvalid, nbranchcands) );
840 SCIP_CALL( SCIPallocBufferArray(scip, &sbvars, nbranchcands) );
841
842 if( nbranchcands == 1 )
843 {
844 /* only one candidate: nothing has to be done */
845 bestcand = 0;
846 SCIPdebug(ninitcands = 0);
847 sbvars[0] = FALSE;
848 }
849 else
850 {
851 SCIP_VAR** vars;
852 int* initcands;
853 SCIP_Real* initcandscores;
854 SCIP_Real* newlbs = NULL;
855 SCIP_Real* newubs = NULL;
856 SCIP_Real* mingains = NULL;
857 SCIP_Real* maxgains = NULL;
858 /* scores computed from pseudocost branching */
859 SCIP_Real* scores = NULL;
860 SCIP_Real* scoresfrompc = NULL;
861 SCIP_Real* scoresfromothers = NULL;
862 int* bdchginds;
863 SCIP_BOUNDTYPE* bdchgtypes;
864 SCIP_Real* bdchgbounds;
865 int maxninitcands;
866 int nuninitcands;
867 int nbdchgs;
868 int nbdconflicts;
869 SCIP_Real avgconflictscore;
870 SCIP_Real avgconflengthscore;
871 SCIP_Real avginferencescore;
872 SCIP_Real avgcutoffscore;
873 SCIP_Real avgpscostscore;
874 SCIP_Real bestpsscore;
875 SCIP_Real bestpsfracscore;
876 SCIP_Real bestpsdomainscore;
877 SCIP_Real bestsbscore;
878 SCIP_Real bestuninitsbscore;
879 SCIP_Real bestsbfracscore;
880 SCIP_Real bestsbdomainscore;
881 SCIP_Real prio;
882 SCIP_Real reliable;
883 SCIP_Real maxlookahead;
884 SCIP_Real lookahead;
885 SCIP_Real relerrorthreshold;
886 SCIP_Bool initstrongbranching;
887 SCIP_Bool propagate;
888 SCIP_Bool probingbounds;
889 SCIP_Longint nodenum;
890 SCIP_Longint nlpiterationsquot;
891 SCIP_Longint nsblpiterations;
892 SCIP_Longint maxnsblpiterations;
893 int bestsolidx;
894 int maxbdchgs;
895 int bestpscand;
896 int bestsbcand;
897 int bestuninitsbcand;
898 int inititer;
899 int nvars;
900 int i;
901 int c;
903 SCIP_Real degeneracyfactor = 1.0;
904
905 /* get LP degeneracy information and compute a factor to change weighting of pseudo cost score vs. other scores */
906 if( branchruledata->degeneracyaware > 0 && (SCIPgetDepth(scip) > 0 || branchruledata->degeneracyaware > 1) )
907 {
908 SCIP_Real degeneracy;
909 SCIP_Real varconsratio;
910
911 /* get LP degeneracy information */
912 SCIP_CALL( SCIPgetLPDualDegeneracy(scip, &degeneracy, &varconsratio) );
913
914 assert(degeneracy >= 0.0);
915 assert(degeneracy <= 1.0);
916 assert(varconsratio >= 1.0);
917
918 /* increase factor for a degeneracy >= 80% */
919 if( degeneracy >= 0.8 )
920 {
921 degeneracy = 10.0 * (degeneracy - 0.7);
922 degeneracyfactor = degeneracyfactor * pow(10.0,degeneracy);
923 }
924 /* increase factor for a variable-constraint ratio >= 2.0 */
925 if( varconsratio >= 2.0 )
926 {
927 degeneracyfactor *= 10.0 * varconsratio;
928 }
929 }
930
931 vars = SCIPgetVars(scip);
932 nvars = SCIPgetNVars(scip);
933
934 bestsolidx = SCIPgetBestSol(scip) == NULL ? -1 : SCIPsolGetIndex(SCIPgetBestSol(scip));
935
936 /* get average conflict, inference, and pseudocost scores */
937 avgconflictscore = SCIPgetAvgConflictScore(scip);
938 avgconflictscore = MAX(avgconflictscore, 0.1);
939 avgconflengthscore = SCIPgetAvgConflictlengthScore(scip);
940 avgconflengthscore = MAX(avgconflengthscore, 0.1);
941 avginferencescore = SCIPgetAvgInferenceScore(scip);
942 avginferencescore = MAX(avginferencescore, 0.1);
943 avgcutoffscore = SCIPgetAvgCutoffScore(scip);
944 avgcutoffscore = MAX(avgcutoffscore, 0.1);
945 avgpscostscore = SCIPgetAvgPseudocostScore(scip);
946 avgpscostscore = MAX(avgpscostscore, 0.1);
947
948 /* get nonlinear counts according to parameters */
949 SCIP_CALL( branchruledataEnsureNlcount(scip, branchruledata) );
950
951 initstrongbranching = FALSE;
952
953 /* check whether propagation should be performed */
954 propagate = (branchruledata->maxproprounds != 0);
955
956 /* check whether valid bounds should be identified in probing-like fashion */
957 probingbounds = propagate && branchruledata->probingbounds;
958
959 /* get maximal number of candidates to initialize with strong branching; if the current solutions is not basic,
960 * we cannot warmstart the simplex algorithm and therefore don't initialize any candidates
961 */
962 maxninitcands = MIN(nbranchcands, branchruledata->initcand);
963 if( !SCIPisLPSolBasic(scip) )
964 maxninitcands = 0;
965
966 /* calculate maximal number of strong branching LP iterations; if we used too many, don't apply strong branching
967 * any more; also, if degeneracy is too high, don't run strong branching at this node
968 */
969 nlpiterationsquot = (SCIP_Longint)(branchruledata->sbiterquot * SCIPgetNNodeLPIterations(scip));
970 maxnsblpiterations = nlpiterationsquot + branchruledata->sbiterofs + SCIPgetNRootStrongbranchLPIterations(scip);
971 nsblpiterations = SCIPgetNStrongbranchLPIterations(scip);
972 if( nsblpiterations > maxnsblpiterations || degeneracyfactor >= 10.0 )
973 maxninitcands = 0;
974
975 /* get buffer for storing the unreliable candidates */
976 SCIP_CALL( SCIPallocBufferArray(scip, &initcands, maxninitcands+1) ); /* allocate one additional slot for convenience */
977 SCIP_CALL( SCIPallocBufferArray(scip, &initcandscores, maxninitcands+1) );
978 ninitcands = 0;
979
980 /* Allocate memory for the down and up gains, and the computed pseudocost scores */
981 SCIP_CALL( SCIPallocBufferArray(scip, &mingains, nbranchcands) );
982 SCIP_CALL( SCIPallocBufferArray(scip, &maxgains, nbranchcands) );
983 SCIP_CALL( SCIPallocBufferArray(scip, &scores, nbranchcands) );
984 SCIP_CALL( SCIPallocBufferArray(scip, &scoresfrompc, nbranchcands) );
985 SCIP_CALL( SCIPallocBufferArray(scip, &scoresfromothers, nbranchcands) );
986
987 /* get current node number */
988 nodenum = SCIPgetNNodes(scip);
989
990 /* initialize bound change arrays */
991 bdchginds = NULL;
992 bdchgtypes = NULL;
993 bdchgbounds = NULL;
994 nbdchgs = 0;
995 nbdconflicts = 0;
996 maxbdchgs = branchruledata->maxbdchgs;
997 if( maxbdchgs == -1 )
998 maxbdchgs = INT_MAX;
999
1000 /* calculate value used as reliability */
1001 prio = (maxnsblpiterations - nsblpiterations)/(nsblpiterations + 1.0);
1002 prio = MIN(prio, 1.0);
1003 prio = MAX(prio, (nlpiterationsquot - nsblpiterations)/(nsblpiterations + 1.0));
1004 reliable = (1.0-prio) * branchruledata->minreliable + prio * branchruledata->maxreliable;
1005
1006 /* calculate the threshold for the relative error in the same way; low tolerance is more strict than higher tolerance */
1007 relerrorthreshold = (1.0 - prio) * branchruledata->higherrortol + prio * branchruledata->lowerrortol;
1008
1009 clevel = (SCIP_CONFIDENCELEVEL)branchruledata->confidencelevel;
1010 /* determine the confidence level for hypothesis testing based on value of prio */
1011 if( branchruledata->usedynamicconfidence )
1012 {
1013 /* with decreasing priority, use a less strict confidence level */
1014 if( prio >= 0.9 )
1015 clevel = SCIP_CONFIDENCELEVEL_MAX;
1016 else if( prio >= 0.7 )
1018 else if( prio >= 0.5 )
1020 else if( prio >= 0.3 )
1021 clevel = SCIP_CONFIDENCELEVEL_LOW;
1022 else
1023 clevel = SCIP_CONFIDENCELEVEL_MIN;
1024 }
1025
1026 /* search for the best pseudo cost candidate, while remembering unreliable candidates in a sorted buffer */
1027 nuninitcands = 0;
1028 bestpscand = -1;
1029 bestpsscore = -SCIPinfinity(scip);
1030 bestpsfracscore = -SCIPinfinity(scip);
1031 bestpsdomainscore = -SCIPinfinity(scip);
1032
1033 /* search for the best candidate first */
1034 if( branchruledata->usehyptestforreliability )
1035 {
1036 for( c = 0; c < nbranchcands; ++c )
1037 {
1038 SCIP_Real conflictscore;
1039 SCIP_Real conflengthscore;
1040 SCIP_Real inferencescore;
1041 SCIP_Real cutoffscore;
1042 SCIP_Real gmieffscore;
1043 SCIP_Real lastgmieffscore;
1044 SCIP_Real pscostscore;
1045 SCIP_Real nlscore;
1046 SCIP_Real score;
1047
1048 conflictscore = SCIPgetVarConflictScore(scip, branchcands[c]);
1049 conflengthscore = SCIPgetVarConflictlengthScore(scip, branchcands[c]);
1050 inferencescore = SCIPgetVarAvgInferenceScore(scip, branchcands[c]);
1051 cutoffscore = SCIPgetVarAvgCutoffScore(scip, branchcands[c]);
1052 gmieffscore = SCIPgetVarAvgGMIScore(scip, branchcands[c]);
1053 lastgmieffscore = SCIPgetVarLastGMIScore(scip, branchcands[c]);
1054 nlscore = calcNlscore(scip, branchruledata->nlcount, branchruledata->nlcountmax, SCIPvarGetProbindex(branchcands[c]));
1055 pscostscore = SCIPgetVarPseudocostScore(scip, branchcands[c], branchcandssol[c]);
1056
1057 /* replace the pseudo cost score with the already calculated one;
1058 * @todo: use old data for strong branching with propagation?
1059 */
1060 if( SCIPgetVarStrongbranchNode(scip, branchcands[c]) == nodenum )
1061 {
1062 SCIP_Real down;
1063 SCIP_Real up;
1064 SCIP_Real lastlpobjval;
1065 SCIP_Real downgain;
1066 SCIP_Real upgain;
1067
1068 /* use the score of the strong branching call at the current node */
1069 SCIP_CALL( SCIPgetVarStrongbranchLast(scip, branchcands[c], &down, &up, NULL, NULL, NULL, &lastlpobjval) );
1070 downgain = MAX(down - lastlpobjval, 0.0);
1071 upgain = MAX(up - lastlpobjval, 0.0);
1072 pscostscore = SCIPgetBranchScore(scip, branchcands[c], downgain, upgain);
1073
1074 SCIPdebugMsg(scip, " -> strong branching on variable <%s> already performed (down=%g (%+g), up=%g (%+g), pscostscore=%g)\n",
1075 SCIPvarGetName(branchcands[c]), down, downgain, up, upgain, pscostscore);
1076 }
1077
1078 score = calcScore(scip, branchruledata, conflictscore, avgconflictscore, conflengthscore, avgconflengthscore,
1079 inferencescore, avginferencescore, cutoffscore, avgcutoffscore, gmieffscore, lastgmieffscore,
1080 pscostscore, avgpscostscore, nlscore, branchcandsfrac[c], degeneracyfactor);
1081
1082 /* check for better score of candidate */
1083 if( SCIPisSumGE(scip, score, bestpsscore) )
1084 {
1085 SCIP_Real fracscore;
1086 SCIP_Real domainscore;
1087
1088 fracscore = MIN(branchcandsfrac[c], 1.0 - branchcandsfrac[c]);
1089 domainscore = -(SCIPvarGetUbLocal(branchcands[c]) - SCIPvarGetLbLocal(branchcands[c]));
1090 if( SCIPisSumGT(scip, score, bestpsscore)
1091 || SCIPisSumGT(scip, fracscore, bestpsfracscore)
1092 || (SCIPisSumGE(scip, fracscore, bestpsfracscore) && domainscore > bestpsdomainscore) )
1093 {
1094 bestpscand = c;
1095 bestpsscore = score;
1096 bestpsfracscore = fracscore;
1097 bestpsdomainscore = domainscore;
1098 }
1099 }
1100 }
1101 }
1102
1103 for( c = 0; c < nbranchcands; ++c )
1104 {
1105 SCIP_Real conflictscore;
1106 SCIP_Real conflengthscore;
1107 SCIP_Real inferencescore;
1108 SCIP_Real cutoffscore;
1109 SCIP_Real gmieffscore;
1110 SCIP_Real lastgmieffscore;
1111 SCIP_Real pscostscore;
1112 SCIP_Real nlscore;
1113 SCIP_Real score;
1114 SCIP_Bool usesb;
1115 SCIP_Real downgain;
1116 SCIP_Real upgain;
1117 SCIP_Real fracpart;
1118
1119 assert(branchcands[c] != NULL);
1120 assert(!SCIPisFeasIntegral(scip, branchcandssol[c]));
1121 assert(!SCIPisFeasIntegral(scip, SCIPvarGetLPSol(branchcands[c])));
1122
1123 /* Record the variables current pseudocosts. These may be overwritten if
1124 * strong branching is performed.
1125 */
1126 sbvars[c] = FALSE;
1127 fracpart = SCIPfeasFrac(scip, SCIPvarGetLPSol(branchcands[c]));
1128 downgain = SCIPgetVarPseudocostVal(scip, branchcands[c], 0.0 - fracpart);
1129 upgain = SCIPgetVarPseudocostVal(scip, branchcands[c], 1.0 - fracpart);
1130 mingains[c] = MIN(downgain, upgain);
1131 maxgains[c] = MAX(downgain, upgain);
1132
1133 /* get conflict, inference, cutoff, nonlinear, and pseudo cost scores for candidate */
1134 conflictscore = SCIPgetVarConflictScore(scip, branchcands[c]);
1135 conflengthscore = SCIPgetVarConflictlengthScore(scip, branchcands[c]);
1136 inferencescore = SCIPgetVarAvgInferenceScore(scip, branchcands[c]);
1137 cutoffscore = SCIPgetVarAvgCutoffScore(scip, branchcands[c]);
1138 gmieffscore = SCIPgetVarAvgGMIScore(scip, branchcands[c]);
1139 lastgmieffscore = SCIPgetVarLastGMIScore(scip, branchcands[c]);
1140 nlscore = calcNlscore(scip, branchruledata->nlcount, branchruledata->nlcountmax, SCIPvarGetProbindex(branchcands[c]));
1141 pscostscore = SCIPgetVarPseudocostScore(scip, branchcands[c], branchcandssol[c]);
1142 usesb = FALSE;
1143
1144 /* don't use strong branching on variables that have already been initialized at the current node;
1145 * instead replace the pseudo cost score with the already calculated one;
1146 * @todo: use old data for strong branching with propagation?
1147 */
1148 if( SCIPgetVarStrongbranchNode(scip, branchcands[c]) == nodenum )
1149 {
1150 SCIP_Real down;
1151 SCIP_Real up;
1152 SCIP_Real lastlpobjval;
1153
1154 /* use the score of the strong branching call at the current node */
1155 SCIP_CALL( SCIPgetVarStrongbranchLast(scip, branchcands[c], &down, &up, NULL, NULL, NULL, &lastlpobjval) );
1156 downgain = MAX(down - lastlpobjval, 0.0);
1157 upgain = MAX(up - lastlpobjval, 0.0);
1158 pscostscore = SCIPgetBranchScore(scip, branchcands[c], downgain, upgain);
1159
1160 mingains[c] = MIN(downgain, upgain);
1161 maxgains[c] = MAX(downgain, upgain);
1162
1163 SCIPdebugMsg(scip, " -> strong branching on variable <%s> already performed (down=%g (%+g), up=%g (%+g), pscostscore=%g)\n",
1164 SCIPvarGetName(branchcands[c]), down, downgain, up, upgain, pscostscore);
1165 }
1166 else if( maxninitcands > 0 )
1167 {
1168 SCIP_Real downsize;
1169 SCIP_Real upsize;
1170 SCIP_Real size;
1171
1172 /* check, if the pseudo cost score of the variable is reliable */
1175 size = MIN(downsize, upsize);
1176
1177 /* determine if variable is considered reliable based on the current reliability setting */
1178 /* check fixed number threshold (aka original) reliability first */
1179 assert(!branchruledata->usehyptestforreliability || bestpscand >= 0);
1180 usesb = FALSE;
1181 if( size < reliable )
1182 usesb = TRUE;
1183 else if( branchruledata->userelerrorforreliability && branchruledata->usehyptestforreliability )
1184 {
1185 if( !SCIPisVarPscostRelerrorReliable(scip, branchcands[c], relerrorthreshold, clevel) &&
1186 !SCIPsignificantVarPscostDifference(scip, branchcands[bestpscand], branchcandsfrac[bestpscand],
1187 branchcands[c], branchcandsfrac[c], SCIP_BRANCHDIR_DOWNWARDS, clevel, TRUE) &&
1188 !SCIPsignificantVarPscostDifference(scip, branchcands[bestpscand], 1 - branchcandsfrac[bestpscand],
1189 branchcands[c], 1 - branchcandsfrac[c], SCIP_BRANCHDIR_UPWARDS, clevel, TRUE) )
1190 usesb = TRUE;
1191 }
1192 /* check if relative error is tolerable */
1193 else if( branchruledata->userelerrorforreliability &&
1194 !SCIPisVarPscostRelerrorReliable(scip, branchcands[c], relerrorthreshold, clevel))
1195 usesb = TRUE;
1196 /* check if best pseudo-candidate is significantly better in both directions, use strong-branching otherwise */
1197 else if( branchruledata->usehyptestforreliability &&
1198 !SCIPsignificantVarPscostDifference(scip, branchcands[bestpscand], branchcandsfrac[bestpscand],
1199 branchcands[c], branchcandsfrac[c], SCIP_BRANCHDIR_DOWNWARDS, clevel, TRUE) &&
1200 !SCIPsignificantVarPscostDifference(scip, branchcands[bestpscand], 1 - branchcandsfrac[bestpscand],
1201 branchcands[c], 1 - branchcandsfrac[c], SCIP_BRANCHDIR_UPWARDS, clevel, TRUE))
1202 usesb = TRUE;
1203
1204 /* count the number of variables that are completely uninitialized */
1205 if( size < 0.1 )
1206 nuninitcands++;
1207 }
1208
1209 /* combine the five score values */
1210 scoresfrompc[c] = calcScore(scip, branchruledata, 0.0, avgconflictscore, 0.0, avgconflengthscore,
1211 0.0, avginferencescore, 0.0, avgcutoffscore, 0.0, 0.0,
1212 pscostscore, avgpscostscore, 0.0, branchcandsfrac[c], degeneracyfactor);
1213 scoresfromothers[c] = calcScore(scip, branchruledata, conflictscore, avgconflictscore, conflengthscore, avgconflengthscore,
1214 inferencescore, avginferencescore, cutoffscore, avgcutoffscore, gmieffscore, lastgmieffscore,
1215 0.0, avgpscostscore, nlscore, branchcandsfrac[c], degeneracyfactor);
1216 score = scoresfrompc[c] + scoresfromothers[c];
1217 scores[c] = score;
1218 /*score = calcScore(scip, branchruledata, conflictscore, avgconflictscore, conflengthscore, avgconflengthscore,
1219 inferencescore, avginferencescore, cutoffscore, avgcutoffscore, gmieffscore, lastgmieffscore,
1220 pscostscore, avgpscostscore, nlscore, branchcandsfrac[c], degeneracyfactor);*/
1221 if( usesb )
1222 {
1223 int j;
1224
1225 mingains[c] = 0;
1226 maxgains[c] = 0;
1227 scoresfrompc[c] = 0;
1228 scoresfromothers[c] = 0;
1229
1230 /* assign a random score to this uninitialized candidate */
1231 if( branchruledata->randinitorder )
1232 score += SCIPrandomGetReal(branchruledata->randnumgen, 0.0, 1e-4);
1233
1234 /* pseudo cost of variable is not reliable: insert candidate in initcands buffer */
1235 for( j = ninitcands; j > 0 && score > initcandscores[j-1]; --j )
1236 {
1237 initcands[j] = initcands[j-1];
1238 initcandscores[j] = initcandscores[j-1];
1239 }
1240 initcands[j] = c;
1241 initcandscores[j] = score;
1242 ninitcands++;
1243 ninitcands = MIN(ninitcands, maxninitcands);
1244 }
1245 /* in the case of hypothesis reliability, the best pseudo candidate has been determined already */
1246 else if( !branchruledata->usehyptestforreliability )
1247 {
1248 /* variable will keep its pseudo cost value: check for better score of candidate */
1249 if( SCIPisSumGE(scip, score, bestpsscore) )
1250 {
1251 SCIP_Real fracscore;
1252 SCIP_Real domainscore;
1253
1254 fracscore = MIN(branchcandsfrac[c], 1.0 - branchcandsfrac[c]);
1255 domainscore = -(SCIPvarGetUbLocal(branchcands[c]) - SCIPvarGetLbLocal(branchcands[c]));
1256 if( SCIPisSumGT(scip, score, bestpsscore)
1257 || SCIPisSumGT(scip, fracscore, bestpsfracscore)
1258 || (SCIPisSumGE(scip, fracscore, bestpsfracscore) && domainscore > bestpsdomainscore) )
1259 {
1260 bestpscand = c;
1261 bestpsscore = score;
1262 bestpsfracscore = fracscore;
1263 bestpsdomainscore = domainscore;
1264 }
1265 }
1266 }
1267 }
1268
1269 /* in the special case that only the best pseudo candidate was selected for strong branching, skip the strong branching */
1270 if( branchruledata->usehyptestforreliability && ninitcands == 1 )
1271 {
1272 ninitcands = 0;
1273 SCIPdebugMsg(scip, "Only one single candidate for initialization-->Skipping strong branching\n");
1274 }
1275
1276 /* initialize unreliable candidates with strong branching until maxlookahead is reached,
1277 * search best strong branching candidate
1278 */
1279 maxlookahead = (SCIP_Real)branchruledata->maxlookahead * (1.0 + (SCIP_Real)nuninitcands/(SCIP_Real)nbranchcands);
1280 inititer = branchruledata->inititer;
1281 if( inititer == 0 )
1282 {
1283 SCIP_Longint nlpiterations;
1284 SCIP_Longint nlps;
1285
1286 /* iteration limit is set to twice the average number of iterations spent to resolve a dual feasible SCIP_LP;
1287 * at the first few nodes, this average is not very exact, so we better increase the iteration limit on
1288 * these very important nodes
1289 */
1290 nlpiterations = SCIPgetNDualResolveLPIterations(scip);
1292 if( nlps == 0 )
1293 {
1294 nlpiterations = SCIPgetNNodeInitLPIterations(scip);
1295 nlps = SCIPgetNNodeInitLPs(scip);
1296 if( nlps == 0 )
1297 {
1298 nlpiterations = 1000;
1299 nlps = 1;
1300 }
1301 }
1302 assert(nlps >= 1);
1303 inititer = (int)(2*nlpiterations / nlps);
1304 inititer = (int)((SCIP_Real)inititer * (1.0 + 20.0/nodenum));
1305 inititer = MAX(inititer, 10);
1306 inititer = MIN(inititer, 500);
1307 }
1308
1309 SCIPdebugMsg(scip, "strong branching (reliable=%g, %d/%d cands, %d uninit, maxcands=%d, maxlookahead=%g, maxbdchgs=%d, inititer=%d, iterations:%" SCIP_LONGINT_FORMAT "/%" SCIP_LONGINT_FORMAT ", basic:%u)\n",
1310 reliable, ninitcands, nbranchcands, nuninitcands, maxninitcands, maxlookahead, maxbdchgs, inititer,
1312
1313 bestsbcand = -1;
1314 bestsbscore = -SCIPinfinity(scip);
1315 bestsbfracscore = -SCIPinfinity(scip);
1316 bestsbdomainscore = -SCIPinfinity(scip);
1317 bestuninitsbscore = -SCIPinfinity(scip);
1318 bestuninitsbcand = -1;
1319 lookahead = 0.0;
1320 for( i = 0; i < ninitcands && lookahead < maxlookahead && nbdchgs + nbdconflicts < maxbdchgs
1321 && (i < (int) maxlookahead || SCIPgetNStrongbranchLPIterations(scip) < maxnsblpiterations); ++i )
1322 {
1323 SCIP_Real down;
1324 SCIP_Real up;
1325 SCIP_Real downgain;
1326 SCIP_Real upgain;
1327 SCIP_Bool downvalid;
1328 SCIP_Bool upvalid;
1329 SCIP_Longint ndomredsdown;
1330 SCIP_Longint ndomredsup;
1331 SCIP_Bool lperror;
1332 SCIP_Bool downinf;
1333 SCIP_Bool upinf;
1334 SCIP_Bool downconflict;
1335 SCIP_Bool upconflict;
1336
1337 /* get candidate number to initialize */
1338 c = initcands[i];
1339 assert(!SCIPisFeasIntegral(scip, branchcandssol[c]));
1340
1341 if( branchruledata->skipbadinitcands )
1342 {
1343 SCIP_Bool skipsb = FALSE;
1344 /* if the current best candidate is a candidate found by strong branching, determine if candidate pseudo-costs are
1345 * significantly smaller in at least one direction, in which case we safe the execution of strong-branching for now
1346 */
1347 if( bestsbscore > bestpsscore && bestsbscore > bestuninitsbscore && bestsbupvalid && bestsbdownvalid )
1348 {
1349 assert(bestsbcand != -1);
1350 assert(bestsbup != SCIP_INVALID && bestsbdown != SCIP_INVALID); /*lint !e777 lint doesn't like comparing floats */
1351
1352 /* test if the variable is unlikely to produce a better gain than the currently best one. Skip strong-branching
1353 * in such a case
1354 */
1355 if( SCIPpscostThresholdProbabilityTest(scip, branchcands[c], branchcandsfrac[c], bestsbdown,
1357 || SCIPpscostThresholdProbabilityTest(scip, branchcands[c], 1.0 - branchcandsfrac[c], bestsbup,
1358 SCIP_BRANCHDIR_UPWARDS, clevel) )
1359 skipsb = TRUE;
1360 }
1361 /* the currently best candidate is also a pseudo-candidate; apply significance test and skip candidate if it
1362 * is significantly worse in at least one direction
1363 */
1364 else if( bestpscand != -1 && bestpsscore > bestuninitsbscore )
1365 {
1366 if( SCIPsignificantVarPscostDifference(scip, branchcands[bestpscand], branchcandsfrac[bestpscand],
1367 branchcands[c], branchcandsfrac[c], SCIP_BRANCHDIR_DOWNWARDS, clevel, TRUE)
1368 || SCIPsignificantVarPscostDifference(scip, branchcands[bestpscand], 1.0 - branchcandsfrac[bestpscand],
1369 branchcands[c], 1.0 - branchcandsfrac[c], SCIP_BRANCHDIR_UPWARDS, clevel, TRUE) )
1370 skipsb = TRUE;
1371 }
1372 /* compare against the best init cand that has been skipped already */
1373 else if( bestuninitsbcand != -1 )
1374 {
1375 if( SCIPsignificantVarPscostDifference(scip, branchcands[bestuninitsbcand], branchcandsfrac[bestuninitsbcand],
1376 branchcands[c], branchcandsfrac[c], SCIP_BRANCHDIR_DOWNWARDS, clevel, TRUE)
1377 || SCIPsignificantVarPscostDifference(scip, branchcands[bestuninitsbcand], 1.0 - branchcandsfrac[bestuninitsbcand],
1378 branchcands[c], 1.0 - branchcandsfrac[c], SCIP_BRANCHDIR_UPWARDS, clevel, TRUE) )
1379 skipsb = TRUE;
1380 }
1381
1382 /* skip candidate, update the best score of an unitialized candidate */
1383 if( skipsb )
1384 {
1385 if( bestuninitsbcand == -1 )
1386 {
1387 bestuninitsbcand = c;
1388 bestuninitsbscore = initcandscores[i];
1389 }
1390 continue;
1391 }
1392 }
1393 SCIPdebugMsg(scip, "init pseudo cost (%g/%g) of <%s> at %g (score:%g) with strong branching (%d iterations) -- %" SCIP_LONGINT_FORMAT "/%" SCIP_LONGINT_FORMAT " iterations\n",
1396 SCIPvarGetName(branchcands[c]), branchcandssol[c], initcandscores[i],
1397 inititer, SCIPgetNStrongbranchLPIterations(scip), maxnsblpiterations);
1398
1399 /* use strong branching on candidate */
1400 if( !initstrongbranching )
1401 {
1402 initstrongbranching = TRUE;
1403
1404 SCIP_CALL( SCIPstartStrongbranch(scip, propagate) );
1405
1406 /* create arrays for probing-like bound tightening */
1407 if( probingbounds )
1408 {
1409 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &newlbs, nvars) );
1410 SCIP_CALL( SCIPallocBlockMemoryArray(scip, &newubs, nvars) );
1411 }
1412 }
1413
1414 if( propagate )
1415 {
1416 /* apply strong branching */
1417 SCIP_CALL( SCIPgetVarStrongbranchWithPropagation(scip, branchcands[c], branchcandssol[c], lpobjval, inititer,
1418 branchruledata->maxproprounds, &down, &up, &downvalid, &upvalid, &ndomredsdown, &ndomredsup, &downinf, &upinf,
1419 &downconflict, &upconflict, &lperror, newlbs, newubs) );
1420 }
1421 else
1422 {
1423 /* apply strong branching */
1424 SCIP_CALL( SCIPgetVarStrongbranchFrac(scip, branchcands[c], inititer, FALSE,
1425 &down, &up, &downvalid, &upvalid, &downinf, &upinf, &downconflict, &upconflict, &lperror) );
1426
1427 ndomredsdown = ndomredsup = 0;
1428 }
1429
1430 /* check for an error in strong branching */
1431 if( lperror )
1432 {
1433 if( !SCIPisStopped(scip) )
1434 {
1436 "(node %" SCIP_LONGINT_FORMAT ") error in strong branching call for variable <%s> with solution %g\n",
1437 SCIPgetNNodes(scip), SCIPvarGetName(branchcands[c]), branchcandssol[c]);
1438 }
1439 break;
1440 }
1441
1442 /* Strong branching might have found a new primal solution which updated the cutoff bound. In this case, the
1443 * provedbound computed before can be higher than the cutoffbound and the current node can be cut off.
1444 * Additionally, also if the value for the current best candidate is valid and exceeds the new cutoff bound,
1445 * we want to change the domain of this variable rather than branching on it.
1446 */
1447 if( SCIPgetBestSol(scip) != NULL && SCIPsolGetIndex(SCIPgetBestSol(scip)) != bestsolidx )
1448 {
1449 bestsolidx = SCIPsolGetIndex(SCIPgetBestSol(scip));
1450
1451 SCIPdebugMsg(scip, " -> strong branching on variable <%s> lead to a new incumbent\n",
1452 SCIPvarGetName(branchcands[c]));
1453
1454 /* proved bound for current node is larger than new cutoff bound -> cut off current node */
1455 if( SCIPisGE(scip, provedbound, SCIPgetCutoffbound(scip)) )
1456 {
1457 SCIPdebugMsg(scip, " -> node can be cut off (provedbound=%g, cutoff=%g)\n", provedbound, SCIPgetCutoffbound(scip));
1458
1459 *result = SCIP_CUTOFF;
1460 break; /* terminate initialization loop, because node is infeasible */
1461 }
1462 /* proved bound for down child of best candidate is larger than cutoff bound
1463 * -> increase lower bound of best candidate
1464 * we must only do this if the LP is complete, i.e., we are not doing column generation
1465 */
1466
1467 else if( bestsbcand != -1 && allcolsinlp )
1468 {
1469 if( !bestsbdowncutoff && bestsbdownvalid && SCIPisGE(scip, bestsbdown, SCIPgetCutoffbound(scip)) )
1470 {
1471 bestsbdowncutoff = TRUE;
1472
1473 SCIPdebugMsg(scip, " -> valid dual bound for down child of best candidate <%s> is higher than new cutoff bound (valid=%u, bestsbdown=%g, cutoff=%g)\n",
1474 SCIPvarGetName(branchcands[bestsbcand]), bestsbdownvalid, bestsbdown, SCIPgetCutoffbound(scip));
1475
1476 SCIPdebugMsg(scip, " -> increase lower bound of best candidate <%s> to %g\n",
1477 SCIPvarGetName(branchcands[bestsbcand]), SCIPfeasCeil(scip, branchcandssol[bestsbcand]));
1478
1479 SCIP_CALL( addBdchg(scip, &bdchginds, &bdchgtypes, &bdchgbounds, &nbdchgs, SCIPvarGetProbindex(branchcands[bestsbcand]),
1480 SCIP_BOUNDTYPE_LOWER, SCIPfeasCeil(scip, branchcandssol[bestsbcand])) );
1481 }
1482 /* proved bound for up child of best candidate is larger than cutoff bound -> decrease upper bound of best candidate */
1483 else if( !bestsbupcutoff && bestsbupvalid && SCIPisGE(scip, bestsbup, SCIPgetCutoffbound(scip)) )
1484 {
1485 bestsbupcutoff = TRUE;
1486
1487 SCIPdebugMsg(scip, " -> valid dual bound for up child of best candidate <%s> is higher than new cutoff bound (valid=%u, bestsbup=%g, cutoff=%g)\n",
1488 SCIPvarGetName(branchcands[bestsbcand]), bestsbupvalid, bestsbup, SCIPgetCutoffbound(scip));
1489
1490 SCIPdebugMsg(scip, " -> decrease upper bound of best candidate <%s> to %g\n",
1491 SCIPvarGetName(branchcands[bestsbcand]), SCIPfeasFloor(scip, branchcandssol[bestsbcand]));
1492
1493 SCIP_CALL( addBdchg(scip, &bdchginds, &bdchgtypes, &bdchgbounds, &nbdchgs, SCIPvarGetProbindex(branchcands[bestsbcand]),
1494 SCIP_BOUNDTYPE_UPPER, SCIPfeasFloor(scip, branchcandssol[bestsbcand])) );
1495 }
1496 }
1497 }
1498
1499 /* evaluate strong branching */
1500 down = MAX(down, lpobjval);
1501 up = MAX(up, lpobjval);
1502 downgain = down - lpobjval;
1503 upgain = up - lpobjval;
1504 assert(!allcolsinlp || exactsolve || !downvalid || downinf == SCIPisGE(scip, down, SCIPgetCutoffbound(scip)));
1505 assert(!allcolsinlp || exactsolve || !upvalid || upinf == SCIPisGE(scip, up, SCIPgetCutoffbound(scip)));
1506 assert(downinf || !downconflict);
1507 assert(upinf || !upconflict);
1508
1509 /* @todo: store pseudo cost only for valid bounds?
1510 * depending on the user parameter choice of storesemiinitcosts, pseudo costs are also updated in single directions,
1511 * if the node in the other direction was infeasible or cut off
1512 */
1513 if( !downinf
1514#ifdef WITH_LPSOLSTAT
1516#endif
1517 && (!upinf || branchruledata->storesemiinitcosts) )
1518 {
1519 SCIP_Real weight;
1520
1521 /* smaller weights are given if the strong branching hit the time limit in the corresponding direction */
1522 if( branchruledata->usesmallweightsitlim )
1524 else
1525 weight = 1.0;
1526
1527 /* update pseudo cost values */
1528 SCIP_CALL( SCIPupdateVarPseudocostSymmetric(scip, branchruledata, branchcands[c], branchorbitidx, c, 0.0 - branchcandsfrac[c], downgain, weight) );
1529 }
1530 if( !upinf
1531#ifdef WITH_LPSOLSTAT
1533#endif
1534 && (!downinf || branchruledata->storesemiinitcosts) )
1535 {
1536 SCIP_Real weight;
1537
1538 /* smaller weights are given if the strong branching hit the time limit in the corresponding direction */
1539 if( branchruledata->usesmallweightsitlim )
1541 else
1542 weight = 1.0;
1543
1544 /* update pseudo cost values */
1545 SCIP_CALL( SCIPupdateVarPseudocostSymmetric(scip, branchruledata, branchcands[c], branchorbitidx, c, 1.0 - branchcandsfrac[c], upgain, weight) );
1546 }
1547
1548 /* the minimal lower bound of both children is a proved lower bound of the current subtree */
1549 if( allcolsinlp && !exactsolve && downvalid && upvalid )
1550 {
1551 SCIP_Real minbound;
1552
1553 minbound = MIN(down, up);
1554 provedbound = MAX(provedbound, minbound);
1555 assert((downinf && upinf) || SCIPisLT(scip, provedbound, SCIPgetCutoffbound(scip)));
1556
1557 /* save probing-like bounds detected during strong branching */
1558 if( probingbounds )
1559 {
1560 int v;
1561
1562 assert(newlbs != NULL);
1563 assert(newubs != NULL);
1564
1565 for( v = 0; v < nvars; ++v )
1566 {
1567 if( SCIPisGT(scip, newlbs[v], SCIPvarGetLbLocal(vars[v])) )
1568 {
1569 SCIPdebugMsg(scip, "better lower bound for variable <%s>: %.9g -> %.9g (by strongbranching on <%s>)\n",
1570 SCIPvarGetName(vars[v]), SCIPvarGetLbLocal(vars[v]), newlbs[v], SCIPvarGetName(branchcands[c]));
1571
1572 SCIP_CALL( addBdchg(scip, &bdchginds, &bdchgtypes, &bdchgbounds, &nbdchgs, v,
1573 SCIP_BOUNDTYPE_LOWER, newlbs[v]) );
1574 }
1575 if( SCIPisLT(scip, newubs[v], SCIPvarGetUbLocal(vars[v])) )
1576 {
1577 SCIPdebugMsg(scip, "better upper bound for variable <%s>: %.9g -> %.9g (by strongbranching on <%s>)\n",
1578 SCIPvarGetName(vars[v]), SCIPvarGetUbLocal(vars[v]), newubs[v], SCIPvarGetName(branchcands[c]));
1579
1580 SCIP_CALL( addBdchg(scip, &bdchginds, &bdchgtypes, &bdchgbounds, &nbdchgs, v,
1581 SCIP_BOUNDTYPE_UPPER, newubs[v]) );
1582 }
1583 }
1584 }
1585 }
1586
1587 /* check if there are infeasible roundings */
1588 if( downinf || upinf )
1589 {
1590 assert(allcolsinlp || propagate);
1591 assert(!exactsolve);
1592
1593 if( downinf && upinf )
1594 {
1595 /* both roundings are infeasible -> node is infeasible */
1596 SCIPdebugMsg(scip, " -> variable <%s> is infeasible in both directions (conflict: %u/%u)\n",
1597 SCIPvarGetName(branchcands[c]), downconflict, upconflict);
1598 *result = SCIP_CUTOFF;
1599 break; /* terminate initialization loop, because node is infeasible */
1600 }
1601 else
1602 {
1603 /* rounding is infeasible in one direction -> round variable in other direction */
1604 SCIPdebugMsg(scip, " -> variable <%s> is infeasible in %s branch (conflict: %u/%u)\n",
1605 SCIPvarGetName(branchcands[c]), downinf ? "downward" : "upward", downconflict, upconflict);
1606 SCIP_CALL( addBdchg(scip, &bdchginds, &bdchgtypes, &bdchgbounds, &nbdchgs, SCIPvarGetProbindex(branchcands[c]),
1608 (downinf ? SCIPfeasCeil(scip, branchcandssol[c]) : SCIPfeasFloor(scip, branchcandssol[c]))) );
1609 if( nbdchgs + nbdconflicts >= maxbdchgs )
1610 break; /* terminate initialization loop, because enough roundings are performed */
1611 }
1612 }
1613 else
1614 {
1615 SCIP_Real conflictscore;
1616 SCIP_Real conflengthscore;
1617 SCIP_Real inferencescore;
1618 SCIP_Real cutoffscore;
1619 SCIP_Real gmieffscore;
1620 SCIP_Real lastgmieffscore;
1621 SCIP_Real pscostscore;
1622 SCIP_Real nlscore;
1623 SCIP_Real score;
1624
1625 mingains[c] = MIN(downgain, upgain);
1626 maxgains[c] = MAX(downgain, upgain);
1627
1628 sbdown[c] = down;
1629 sbup[c] = up;
1630 sbdownvalid[c] = downvalid;
1631 sbupvalid[c] = upvalid;
1632 sbvars[c] = TRUE;
1633
1634 /* check for a better score */
1635 conflictscore = SCIPgetVarConflictScore(scip, branchcands[c]);
1636 conflengthscore = SCIPgetVarConflictlengthScore(scip, branchcands[c]);
1637 nlscore = calcNlscore(scip, branchruledata->nlcount, branchruledata->nlcountmax, SCIPvarGetProbindex(branchcands[c]));
1638
1639 /* optionally, use only local information obtained via strong branching for this candidate, i.e., local
1640 * domain reductions and no cutoff score
1641 */
1642 inferencescore = branchruledata->usesblocalinfo ? SCIPgetBranchScore(scip, branchcands[c], (SCIP_Real)ndomredsdown, (SCIP_Real)ndomredsup)
1643 : SCIPgetVarAvgInferenceScore(scip, branchcands[c]);
1644 cutoffscore = branchruledata->usesblocalinfo ? 0.0 : SCIPgetVarAvgCutoffScore(scip, branchcands[c]);
1645 gmieffscore = branchruledata->usesblocalinfo ? 0.0 : SCIPgetVarAvgGMIScore(scip, branchcands[c]);
1646 lastgmieffscore = branchruledata->usesblocalinfo ? 0.0 : SCIPgetVarLastGMIScore(scip, branchcands[c]);
1647 pscostscore = SCIPgetBranchScore(scip, branchcands[c], downgain, upgain);
1648
1649 scoresfrompc[c] = calcScore(scip, branchruledata, 0.0, avgconflictscore, 0.0, avgconflengthscore,
1650 0.0, avginferencescore, 0.0, avgcutoffscore, 0.0, 0.0, pscostscore,
1651 avgpscostscore, 0.0, branchcandsfrac[c], degeneracyfactor);
1652 scoresfromothers[c] = calcScore(scip, branchruledata, conflictscore, avgconflictscore, conflengthscore, avgconflengthscore,
1653 inferencescore, avginferencescore, cutoffscore, avgcutoffscore, gmieffscore,
1654 lastgmieffscore, 0.0, avgpscostscore, nlscore, branchcandsfrac[c],
1655 degeneracyfactor);
1656 score = scoresfrompc[c] + scoresfromothers[c];
1657 scores[c] = score;
1658 /*score = calcScore(scip, branchruledata, conflictscore, avgconflictscore, conflengthscore, avgconflengthscore,
1659 inferencescore, avginferencescore, cutoffscore, avgcutoffscore, gmieffscore, lastgmieffscore,
1660 pscostscore, avgpscostscore, nlscore, branchcandsfrac[c], degeneracyfactor);*/
1661
1662 if( SCIPisSumGE(scip, score, bestsbscore) )
1663 {
1664 SCIP_Real fracscore;
1665 SCIP_Real domainscore;
1666
1667 fracscore = MIN(branchcandsfrac[c], 1.0 - branchcandsfrac[c]);
1668 domainscore = -(SCIPvarGetUbLocal(branchcands[c]) - SCIPvarGetLbLocal(branchcands[c]));
1669 if( SCIPisSumGT(scip, score, bestsbscore)
1670 || SCIPisSumGT(scip, fracscore, bestsbfracscore)
1671 || (SCIPisSumGE(scip, fracscore, bestsbfracscore) && domainscore > bestsbdomainscore) )
1672 {
1673 bestsbcand = c;
1674 bestsbscore = score;
1675 bestsbdown = down;
1676 bestsbup = up;
1677 bestsbdownvalid = downvalid;
1678 bestsbupvalid = upvalid;
1679 bestsbdowncutoff = FALSE;
1680 bestsbupcutoff = FALSE;
1681 bestsbfracscore = fracscore;
1682 bestsbdomainscore = domainscore;
1683 lookahead = 0.0;
1684 }
1685 else
1686 lookahead += 0.5;
1687 }
1688 else
1689 lookahead += 1.0;
1690
1691 SCIPdebugMsg(scip, " -> variable <%s> (solval=%g, down=%g (%+g,valid=%u), up=%g (%+g,valid=%u), score=%g/ %g/%g %g/%g %g -> %g)\n",
1692 SCIPvarGetName(branchcands[c]), branchcandssol[c], down, downgain, downvalid, up, upgain, upvalid,
1693 pscostscore, conflictscore, conflengthscore, inferencescore, cutoffscore, gmieffscore, score);
1694 }
1695 }
1696#ifdef SCIP_DEBUG
1697 if( bestsbcand >= 0 )
1698 {
1699 SCIPdebugMsg(scip, " -> best: <%s> (%g / %g / %g), lookahead=%g/%g\n",
1700 SCIPvarGetName(branchcands[bestsbcand]), bestsbscore, bestsbfracscore, bestsbdomainscore,
1701 lookahead, maxlookahead);
1702 }
1703#endif
1704
1705 if( initstrongbranching )
1706 {
1707 if( probingbounds )
1708 {
1709 assert(newlbs != NULL);
1710 assert(newubs != NULL);
1711
1712 SCIPfreeBlockMemoryArray(scip, &newubs, nvars);
1713 SCIPfreeBlockMemoryArray(scip, &newlbs, nvars);
1714 }
1715
1717
1719 {
1720 assert(SCIPhasCurrentNodeLP(scip));
1721 *result = SCIP_CUTOFF;
1722 }
1723 }
1724
1725 if( *result != SCIP_CUTOFF )
1726 {
1727 /* get the score of the best uninitialized strong branching candidate */
1728 if( i < ninitcands && bestuninitsbcand == -1 )
1729 bestuninitsbscore = initcandscores[i];
1730
1731 /* if the best pseudo cost candidate is better than the best uninitialized strong branching candidate,
1732 * compare it to the best initialized strong branching candidate
1733 */
1734 if( bestpsscore > bestuninitsbscore && SCIPisSumGT(scip, bestpsscore, bestsbscore) )
1735 {
1736 bestcand = bestpscand;
1737 bestisstrongbranch = FALSE;
1738 }
1739 else if( bestsbcand >= 0 )
1740 {
1741 bestcand = bestsbcand;
1742 bestisstrongbranch = TRUE;
1743 }
1744 else
1745 {
1746 /* no candidate was initialized, and the best score is the one of the first candidate in the initialization
1747 * queue
1748 */
1749 assert(ninitcands >= 1);
1750 bestcand = initcands[0];
1751 bestisstrongbranch = FALSE;
1752 }
1753
1754 /* update lower bound of current node */
1755 if( allcolsinlp && !exactsolve )
1756 {
1757 assert(SCIPisLT(scip, provedbound, SCIPgetCutoffbound(scip)));
1759 }
1760 }
1761
1762 /* apply domain reductions */
1763 if( nbdchgs > 0 )
1764 {
1765 if( *result != SCIP_CUTOFF )
1766 {
1767 SCIP_CALL( applyBdchgs(scip, vars, bdchginds, bdchgtypes, bdchgbounds, nbdchgs, result) );
1768 if( *result != SCIP_CUTOFF )
1769 *result = SCIP_REDUCEDDOM;
1770 }
1771 freeBdchgs(scip, &bdchginds, &bdchgtypes, &bdchgbounds, &nbdchgs);
1772 }
1773
1774 /* Apply the Treemodel branching rule to potentially select a better branching candidate than the current one. */
1775 if( *result != SCIP_CUTOFF && *result != SCIP_REDUCEDDOM && *result != SCIP_CONSADDED && SCIPtreemodelIsEnabled(scip, branchruledata->treemodel) )
1776 {
1777 SCIP_Real smallpscost;
1778 SCIP_Bool usetreemodel;
1779
1780 usetreemodel = TRUE;
1781
1782 /* If the pseudocosts are zero, use SCIPs best variable since the Treemodel is not applicable */
1783 if( SCIPisZero(scip, maxgains[bestcand]))
1784 {
1785 usetreemodel = FALSE;
1786 }
1787
1788 /* If SCIPs best candidate was selected due to hybrid branching scores
1789 * rather than because of pseudocosts, then we keep it.
1790 */
1791 SCIP_CALL( SCIPgetRealParam(scip, "branching/treemodel/smallpscost", &smallpscost) );
1792 if( usetreemodel == TRUE && avgpscostscore <= smallpscost )
1793 {
1794 int cand;
1795 for( cand = 0; cand < nbranchcands; ++cand )
1796 {
1797 if( scoresfrompc[cand] > scoresfrompc[bestcand] )
1798 {
1799 usetreemodel = FALSE;
1800 break;
1801 }
1802 }
1803 }
1804
1805 if( usetreemodel == TRUE )
1806 {
1808 scip, /* SCIP data structure */
1809 branchruledata->treemodel, /* branching rule */
1810 branchcands, /* branching candidate storage */
1811 mingains, /* minimum gain of rounding downwards or upwards */
1812 maxgains, /* maximum gain of rounding downwards or upwards */
1813 scoresfromothers, /* scores from other branching methods */
1814 nbranchcands, /* the number of branching candidates */
1815 &bestcand /* the best branching candidate found by SCIP */
1816 ) );
1817 }
1818 }
1819
1820 /* free buffer for the lp gains and pseudocost scores */
1821 SCIPfreeBufferArray(scip, &scoresfromothers);
1822 SCIPfreeBufferArray(scip, &scoresfrompc);
1823 SCIPfreeBufferArray(scip, &scores);
1824 SCIPfreeBufferArray(scip, &maxgains);
1825 SCIPfreeBufferArray(scip, &mingains);
1826
1827 /* free buffer for the unreliable candidates */
1828 SCIPfreeBufferArray(scip, &initcandscores);
1829 SCIPfreeBufferArray(scip, &initcands);
1830 }
1831
1832 /* if no domain could be reduced, create the branching */
1833 if( *result != SCIP_CUTOFF && *result != SCIP_REDUCEDDOM && *result != SCIP_CONSADDED && executebranch )
1834 {
1835 SCIP_NODE* downchild;
1836 SCIP_NODE* upchild;
1837 SCIP_VAR* var;
1838 SCIP_Real val;
1839
1840 assert(*result == SCIP_DIDNOTRUN);
1841 assert(0 <= bestcand && bestcand < nbranchcands);
1842 assert(!SCIPisFeasIntegral(scip, branchcandssol[bestcand]));
1843 assert(!allcolsinlp || SCIPisLT(scip, provedbound, SCIPgetCutoffbound(scip)));
1844 assert(!bestsbdowncutoff && !bestsbupcutoff);
1845
1846 var = branchcands[bestcand];
1847 val = branchcandssol[bestcand];
1848
1849 /* perform the branching */
1850 SCIPdebugMsg(scip, " -> %d (%d) cands, sel cand %d: var <%s> (sol=%g, down=%g (%+g), up=%g (%+g), sb=%u, psc=%g/%g [%g])\n",
1851 nbranchcands, ninitcands, bestcand, SCIPvarGetName(var), branchcandssol[bestcand],
1852 bestsbdown, bestsbdown - lpobjval, bestsbup, bestsbup - lpobjval, bestisstrongbranch,
1855 SCIPgetVarPseudocostScoreCurrentRun(scip, var, branchcandssol[bestcand]));
1856 SCIP_UNUSED(bestisstrongbranch);
1857 SCIP_CALL( SCIPbranchVarVal(scip, var, val, &downchild, NULL, &upchild) );
1858 assert(downchild != NULL);
1859 assert(upchild != NULL);
1860
1861 /* update the lower bounds in the children */
1862 if( sbvars[bestcand] && allcolsinlp && !exactsolve )
1863 {
1864 if( sbdownvalid[bestcand] )
1865 {
1866 assert(SCIPisLT(scip, sbdown[bestcand], SCIPgetCutoffbound(scip)));
1867 SCIP_CALL( SCIPupdateNodeLowerbound(scip, downchild, sbdown[bestcand]) );
1868 assert(SCIPisGE(scip, SCIPgetNodeLowerbound(scip, downchild), provedbound));
1869 }
1870 if( sbupvalid[bestcand] )
1871 {
1872 assert(SCIPisLT(scip, sbup[bestcand], SCIPgetCutoffbound(scip)));
1873 SCIP_CALL( SCIPupdateNodeLowerbound(scip, upchild, sbup[bestcand]) );
1874 assert(SCIPisGE(scip, SCIPgetNodeLowerbound(scip, upchild), provedbound));
1875 }
1876 }
1877
1878 SCIPdebugMsg(scip, " -> down child's lowerbound: %g\n", SCIPnodeGetLowerbound(downchild));
1879 SCIPdebugMsg(scip, " -> up child's lowerbound : %g\n", SCIPnodeGetLowerbound(upchild));
1880
1882
1883 *result = SCIP_BRANCHED;
1884 }
1885
1886 /* free buffer for the strong branching lp gains */
1887 SCIPfreeBufferArray(scip, &sbvars);
1888 SCIPfreeBufferArray(scip, &sbupvalid);
1889 SCIPfreeBufferArray(scip, &sbdownvalid);
1890 SCIPfreeBufferArray(scip, &sbup);
1891 SCIPfreeBufferArray(scip, &sbdown);
1892
1893 return SCIP_OKAY;
1894}
1895
1896
1897/*
1898 * Callback methods
1899 */
1900
1901/** copy method for branchrule plugins (called when SCIP copies plugins) */
1902static
1903SCIP_DECL_BRANCHCOPY(branchCopyRelpscost)
1904{ /*lint --e{715}*/
1905 assert(scip != NULL);
1906 assert(branchrule != NULL);
1907 assert(strcmp(SCIPbranchruleGetName(branchrule), BRANCHRULE_NAME) == 0);
1908
1909 /* call inclusion method of branchrule */
1911
1912 return SCIP_OKAY;
1913}
1914
1915/** destructor of branching rule to free user data (called when SCIP is exiting) */
1916static
1917SCIP_DECL_BRANCHFREE(branchFreeRelpscost)
1918{ /*lint --e{715}*/
1919 SCIP_BRANCHRULEDATA* branchruledata;
1920 branchruledata = SCIPbranchruleGetData(branchrule);
1921
1922 /* free Treemodel parameter data structure */
1923 SCIP_CALL( SCIPtreemodelFree(scip, &branchruledata->treemodel) );
1924
1925 /* free branching rule data */
1926 SCIPfreeBlockMemory(scip, &branchruledata);
1927 SCIPbranchruleSetData(branchrule, NULL);
1928
1929 return SCIP_OKAY;
1930}
1931
1932
1933/** solving process initialization method of branching rule (called when branch and bound process is about to begin) */
1934static
1935SCIP_DECL_BRANCHINITSOL(branchInitsolRelpscost)
1936{ /*lint --e{715}*/
1937 SCIP_BRANCHRULEDATA* branchruledata;
1938
1939 /* initialize branching rule data */
1940 branchruledata = SCIPbranchruleGetData(branchrule);
1941 branchruledata->nlcount = NULL;
1942 branchruledata->nlcountsize = 0;
1943 branchruledata->nlcountmax = 1;
1944 assert(branchruledata->startrandseed >= 0);
1945
1946 /* create a random number generator */
1947 SCIP_CALL( SCIPcreateRandom(scip, &branchruledata->randnumgen,
1948 (unsigned int)branchruledata->startrandseed, TRUE) );
1949
1950 return SCIP_OKAY;
1951}
1952
1953
1954/** solving process deinitialization method of branching rule (called before branch and bound process data is freed) */
1955static
1956SCIP_DECL_BRANCHEXITSOL(branchExitsolRelpscost)
1957{ /*lint --e{715}*/
1958 SCIP_BRANCHRULEDATA* branchruledata;
1959
1960 /* free memory in branching rule data */
1961 branchruledata = SCIPbranchruleGetData(branchrule);
1962 SCIPfreeBlockMemoryArrayNull(scip, &branchruledata->nlcount, branchruledata->nlcountsize);
1963
1964 /* free random number generator */
1965 SCIPfreeRandom(scip, &branchruledata->randnumgen);
1966
1967 SCIPfreeBlockMemoryArrayNull(scip, &branchruledata->orbitrep, branchruledata->npermvars);
1968 SCIPfreeBlockMemoryArrayNull(scip, &branchruledata->varorbitmap, branchruledata->npermvars);
1969 SCIPfreeBlockMemoryArrayNull(scip, &branchruledata->orbitbegins, branchruledata->npermvars);
1970 SCIPfreeBlockMemoryArrayNull(scip, &branchruledata->orbits, branchruledata->npermvars);
1971 branchruledata->nosymmetry = FALSE;
1972 branchruledata->norbits = 0;
1973 branchruledata->permvars = NULL;
1974 branchruledata->permvarmap = NULL;
1975 branchruledata->npermvars = 0;
1976
1977 return SCIP_OKAY;
1978}
1979
1980
1981/** branching execution method for fractional LP solutions */
1982static
1983SCIP_DECL_BRANCHEXECLP(branchExeclpRelpscost)
1984{ /*lint --e{715}*/
1985 SCIP_BRANCHRULEDATA* branchruledata;
1986 SCIP_VAR** lpcands;
1987 SCIP_Real* lpcandssol;
1988 SCIP_Real* lpcandsfrac;
1989 SCIP_VAR** filteredlpcands;
1990 SCIP_Real* filteredlpcandssol;
1991 SCIP_Real* filteredlpcandsfrac;
1992 SCIP_Bool runfiltering;
1993 int* filteredlpcandsorbitidx = NULL;
1994 int nfilteredlpcands;
1995 int nlpcands;
1996
1997 assert(branchrule != NULL);
1998 assert(strcmp(SCIPbranchruleGetName(branchrule), BRANCHRULE_NAME) == 0);
1999 assert(scip != NULL);
2000 assert(result != NULL);
2001
2002 SCIPdebugMsg(scip, "Execlp method of relpscost branching in node %" SCIP_LONGINT_FORMAT "\n", SCIPnodeGetNumber(SCIPgetCurrentNode(scip)));
2003
2005 {
2006 *result = SCIP_DIDNOTRUN;
2007 SCIPdebugMsg(scip, "Could not apply relpscost branching, as the current LP was not solved to optimality.\n");
2008
2009 return SCIP_OKAY;
2010 }
2011
2012 /* get branching candidates */
2013 SCIP_CALL( SCIPgetLPBranchCands(scip, &lpcands, &lpcandssol, &lpcandsfrac, NULL, &nlpcands, NULL) );
2014 assert(nlpcands > 0);
2015
2016 branchruledata = SCIPbranchruleGetData(branchrule);
2017 assert(branchruledata != NULL);
2018
2019 /* determine whether we should run filtering */
2020 runfiltering = ! branchruledata->nosymmetry && branchruledata->filtercandssym && SCIPgetSubscipDepth(scip) == 0 && ! SCIPinDive(scip) && ! SCIPinProbing(scip);
2021
2022 /* init orbits if necessary */
2023 if( runfiltering )
2024 {
2025 SCIP_CALL( initOrbits(scip, branchruledata) );
2026 }
2027
2028 /* determine fractional variables (possibly filter by using symmetries) */
2029 if( runfiltering && branchruledata->norbits != 0 )
2030 {
2031 SCIP_CALL( SCIPallocBufferArray(scip, &filteredlpcands, nlpcands) );
2032 SCIP_CALL( SCIPallocBufferArray(scip, &filteredlpcandssol, nlpcands) );
2033 SCIP_CALL( SCIPallocBufferArray(scip, &filteredlpcandsfrac, nlpcands) );
2034 SCIP_CALL( SCIPallocBufferArray(scip, &filteredlpcandsorbitidx, nlpcands) );
2035
2036 /* determine filtered fractional variables */
2037 SCIP_CALL( filterSymmetricVariables(scip, branchruledata, lpcands, lpcandssol, lpcandsfrac, nlpcands,
2038 filteredlpcands, filteredlpcandssol, filteredlpcandsfrac, filteredlpcandsorbitidx, &nfilteredlpcands) );
2039 }
2040 else
2041 {
2042 /* No orbits available. Copy all (unfiltered) branching candidates, because they will be updated w.r.t. the strong branching LP solution */
2043 SCIP_CALL( SCIPduplicateBufferArray(scip, &filteredlpcands, lpcands, nlpcands) );
2044 SCIP_CALL( SCIPduplicateBufferArray(scip, &filteredlpcandssol, lpcandssol, nlpcands) );
2045 SCIP_CALL( SCIPduplicateBufferArray(scip, &filteredlpcandsfrac, lpcandsfrac, nlpcands) );
2046 nfilteredlpcands = nlpcands;
2047 }
2048
2049 /* execute branching rule */
2050 SCIP_CALL( execRelpscost(scip, branchrule, filteredlpcands, filteredlpcandssol, filteredlpcandsfrac, filteredlpcandsorbitidx, nfilteredlpcands, TRUE, result) );
2051
2052 SCIPfreeBufferArrayNull(scip, &filteredlpcandsorbitidx);
2053 SCIPfreeBufferArray(scip, &filteredlpcandsfrac);
2054 SCIPfreeBufferArray(scip, &filteredlpcandssol);
2055 SCIPfreeBufferArray(scip, &filteredlpcands);
2056
2057 return SCIP_OKAY;
2058}
2059
2060
2061/*
2062 * branching specific interface methods
2063 */
2064
2065/** creates the reliable pseudo cost branching rule and includes it in SCIP */
2067 SCIP* scip /**< SCIP data structure */
2068 )
2069{
2070 SCIP_BRANCHRULEDATA* branchruledata;
2071 SCIP_BRANCHRULE* branchrule;
2072
2073 /* create relpscost branching rule data */
2074 SCIP_CALL( SCIPallocBlockMemory(scip, &branchruledata) );
2075
2076 branchruledata->nosymmetry = FALSE;
2077 branchruledata->orbits = NULL;
2078 branchruledata->orbitbegins = NULL;
2079 branchruledata->orbitrep = NULL;
2080 branchruledata->varorbitmap = NULL;
2081 branchruledata->norbits = 0;
2082 branchruledata->permvars = NULL;
2083 branchruledata->npermvars = 0;
2084 branchruledata->permvarmap = NULL;
2085
2086 /* include branching rule */
2088 BRANCHRULE_MAXDEPTH, BRANCHRULE_MAXBOUNDDIST, branchruledata) );
2089
2090 assert(branchrule != NULL);
2091
2092 /* set non-fundamental callbacks via specific setter functions*/
2093 SCIP_CALL( SCIPsetBranchruleCopy(scip, branchrule, branchCopyRelpscost) );
2094 SCIP_CALL( SCIPsetBranchruleFree(scip, branchrule, branchFreeRelpscost) );
2095 SCIP_CALL( SCIPsetBranchruleInitsol(scip, branchrule, branchInitsolRelpscost) );
2096 SCIP_CALL( SCIPsetBranchruleExitsol(scip, branchrule, branchExitsolRelpscost) );
2097 SCIP_CALL( SCIPsetBranchruleExecLp(scip, branchrule, branchExeclpRelpscost) );
2098
2099 /* relpscost branching rule parameters */
2101 "branching/relpscost/conflictweight",
2102 "weight in score calculations for conflict score",
2103 &branchruledata->conflictweight, TRUE, DEFAULT_CONFLICTWEIGHT, SCIP_REAL_MIN, SCIP_REAL_MAX, NULL, NULL) );
2105 "branching/relpscost/conflictlengthweight",
2106 "weight in score calculations for conflict length score",
2107 &branchruledata->conflengthweight, TRUE, DEFAULT_CONFLENGTHWEIGHT, SCIP_REAL_MIN, SCIP_REAL_MAX, NULL, NULL) );
2109 "branching/relpscost/inferenceweight",
2110 "weight in score calculations for inference score",
2111 &branchruledata->inferenceweight, TRUE, DEFAULT_INFERENCEWEIGHT, SCIP_REAL_MIN, SCIP_REAL_MAX, NULL, NULL) );
2113 "branching/relpscost/cutoffweight",
2114 "weight in score calculations for cutoff score",
2115 &branchruledata->cutoffweight, TRUE, DEFAULT_CUTOFFWEIGHT, SCIP_REAL_MIN, SCIP_REAL_MAX, NULL, NULL) );
2117 "branching/relpscost/gmiavgeffweight",
2118 "weight in score calculations for average GMI cuts normalized efficacy",
2119 &branchruledata->gmiavgeffweight, TRUE, DEFAULT_GMIAVGEFFWEIGHT, SCIP_REAL_MIN, SCIP_REAL_MAX, NULL, NULL) );
2121 "branching/relpscost/gmilasteffweight",
2122 "weight in score calculations for last GMI cuts normalized efficacy",
2123 &branchruledata->gmilasteffweight, TRUE, DEFAULT_GMILASTEFFWEIGHT, SCIP_REAL_MIN, SCIP_REAL_MAX, NULL, NULL) );
2125 "branching/relpscost/pscostweight",
2126 "weight in score calculations for pseudo cost score",
2127 &branchruledata->pscostweight, TRUE, DEFAULT_PSCOSTWEIGHT, SCIP_REAL_MIN, SCIP_REAL_MAX, NULL, NULL) );
2129 "branching/relpscost/nlscoreweight",
2130 "weight in score calculations for nlcount score",
2131 &branchruledata->nlscoreweight, TRUE, DEFAULT_NLSCOREWEIGHT, SCIP_REAL_MIN, SCIP_REAL_MAX, NULL, NULL) );
2133 "branching/relpscost/minreliable",
2134 "minimal value for minimum pseudo cost size to regard pseudo cost value as reliable",
2135 &branchruledata->minreliable, TRUE, DEFAULT_MINRELIABLE, 0.0, SCIP_REAL_MAX, NULL, NULL) );
2137 "branching/relpscost/maxreliable",
2138 "maximal value for minimum pseudo cost size to regard pseudo cost value as reliable",
2139 &branchruledata->maxreliable, TRUE, DEFAULT_MAXRELIABLE, 0.0, SCIP_REAL_MAX, NULL, NULL) );
2141 "branching/relpscost/sbiterquot",
2142 "maximal fraction of strong branching LP iterations compared to node relaxation LP iterations",
2143 &branchruledata->sbiterquot, FALSE, DEFAULT_SBITERQUOT, 0.0, SCIP_REAL_MAX, NULL, NULL) );
2145 "branching/relpscost/sbiterofs",
2146 "additional number of allowed strong branching LP iterations",
2147 &branchruledata->sbiterofs, FALSE, DEFAULT_SBITEROFS, 0, INT_MAX, NULL, NULL) );
2149 "branching/relpscost/maxlookahead",
2150 "maximal number of further variables evaluated without better score",
2151 &branchruledata->maxlookahead, TRUE, DEFAULT_MAXLOOKAHEAD, 1, INT_MAX, NULL, NULL) );
2153 "branching/relpscost/initcand",
2154 "maximal number of candidates initialized with strong branching per node",
2155 &branchruledata->initcand, FALSE, DEFAULT_INITCAND, 0, INT_MAX, NULL, NULL) );
2157 "branching/relpscost/inititer",
2158 "iteration limit for strong branching initializations of pseudo cost entries (0: auto)",
2159 &branchruledata->inititer, FALSE, DEFAULT_INITITER, 0, INT_MAX, NULL, NULL) );
2161 "branching/relpscost/maxbdchgs",
2162 "maximal number of bound tightenings before the node is reevaluated (-1: unlimited)",
2163 &branchruledata->maxbdchgs, TRUE, DEFAULT_MAXBDCHGS, -1, INT_MAX, NULL, NULL) );
2165 "branching/relpscost/maxproprounds",
2166 "maximum number of propagation rounds to be performed during strong branching before solving the LP (-1: no limit, -2: parameter settings)",
2167 &branchruledata->maxproprounds, TRUE, DEFAULT_MAXPROPROUNDS, -2, INT_MAX, NULL, NULL) );
2169 "branching/relpscost/probingbounds",
2170 "should valid bounds be identified in a probing-like fashion during strong branching (only with propagation)?",
2171 &branchruledata->probingbounds, TRUE, DEFAULT_PROBINGBOUNDS, NULL, NULL) );
2172
2173 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/userelerrorreliability",
2174 "should reliability be based on relative errors?", &branchruledata->userelerrorforreliability, TRUE, DEFAULT_USERELERRORFORRELIABILITY,
2175 NULL, NULL) );
2176
2177 SCIP_CALL( SCIPaddRealParam(scip, "branching/relpscost/lowerrortol", "low relative error tolerance for reliability",
2178 &branchruledata->lowerrortol, TRUE, DEFAULT_LOWERRORTOL, 0.0, SCIP_REAL_MAX, NULL, NULL) );
2179
2180 SCIP_CALL( SCIPaddRealParam(scip, "branching/relpscost/higherrortol", "high relative error tolerance for reliability",
2181 &branchruledata->higherrortol, TRUE, DEFAULT_HIGHERRORTOL, 0.0, SCIP_REAL_MAX, NULL, NULL) );
2182
2183/**! [SnippetCodeStyleParenIndent] */
2184
2185 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/storesemiinitcosts",
2186 "should strong branching result be considered for pseudo costs if the other direction was infeasible?",
2187 &branchruledata->storesemiinitcosts, TRUE, DEFAULT_STORESEMIINITCOSTS,
2188 NULL, NULL) );
2189
2190/**! [SnippetCodeStyleParenIndent] */
2191
2192 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/usesblocalinfo",
2193 "should the scoring function use only local cutoff and inference information obtained for strong branching candidates?",
2194 &branchruledata->usesblocalinfo, TRUE, DEFAULT_USESBLOCALINFO,
2195 NULL, NULL) );
2196
2197 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/usehyptestforreliability",
2198 "should the strong branching decision be based on a hypothesis test?",
2199 &branchruledata->usehyptestforreliability, TRUE, DEFAULT_USEHYPTESTFORRELIABILITY,
2200 NULL, NULL) );
2201
2202 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/usedynamicconfidence",
2203 "should the confidence level be adjusted dynamically?",
2204 &branchruledata->usedynamicconfidence, TRUE, DEFAULT_USEDYNAMICCONFIDENCE,
2205 NULL, NULL) );
2206 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/skipbadinitcands",
2207 "should branching rule skip candidates that have a low probability to "
2208 "be better than the best strong-branching or pseudo-candidate?",
2209 &branchruledata->skipbadinitcands, TRUE, DEFAULT_SKIPBADINITCANDS,
2210 NULL, NULL) );
2212 "branching/relpscost/confidencelevel",
2213 "the confidence level for statistical methods, between 0 (Min) and 4 (Max).",
2214 &branchruledata->confidencelevel, TRUE, DEFAULT_CONFIDENCELEVEL, 0, 4, NULL, NULL) );
2215
2216 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/randinitorder",
2217 "should candidates be initialized in randomized order?",
2218 &branchruledata->randinitorder, TRUE, DEFAULT_RANDINITORDER,
2219 NULL, NULL) );
2220
2221 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/usesmallweightsitlim",
2222 "should smaller weights be used for pseudo cost updates after hitting the LP iteration limit?",
2223 &branchruledata->usesmallweightsitlim, TRUE, DEFAULT_USESMALLWEIGHTSITLIM,
2224 NULL, NULL) );
2225
2226 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/dynamicweights",
2227 "should the weights of the branching rule be adjusted dynamically during solving based on objective and infeasible leaf counters?",
2228 &branchruledata->dynamicweights, TRUE, DEFAULT_DYNAMICWEIGHTS,
2229 NULL, NULL) );
2230 SCIP_CALL( SCIPaddIntParam(scip, "branching/relpscost/degeneracyaware",
2231 "should degeneracy be taken into account to update weights and skip strong branching? (0: off, 1: after root, 2: always)",
2232 &branchruledata->degeneracyaware, TRUE, DEFAULT_DEGENERACYAWARE, 0, 2,
2233 NULL, NULL) );
2234 SCIP_CALL( SCIPaddIntParam(scip, "branching/relpscost/startrandseed", "start seed for random number generation",
2235 &branchruledata->startrandseed, TRUE, DEFAULT_STARTRANDSEED, 0, INT_MAX, NULL, NULL) );
2236
2237 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/filtercandssym",
2238 "Use symmetry to filter branching candidates?",
2239 &branchruledata->filtercandssym, TRUE, DEFAULT_FILTERCANDSSYM, NULL, NULL) );
2240
2241 SCIP_CALL( SCIPaddBoolParam(scip, "branching/relpscost/transsympscost",
2242 "Transfer pscost information to symmetric variables?",
2243 &branchruledata->transsympscost, TRUE, DEFAULT_TRANSSYMPSCOST, NULL, NULL) );
2244
2245 /* initialise the Treemodel parameters */
2246 SCIP_CALL( SCIPtreemodelInit(scip, &branchruledata->treemodel) );
2247
2248 return SCIP_OKAY;
2249}
2250
2251/** execution reliability pseudo cost branching with the given branching candidates */
2253 SCIP* scip, /**< SCIP data structure */
2254 SCIP_VAR** branchcands, /**< branching candidates */
2255 SCIP_Real* branchcandssol, /**< solution value for the branching candidates */
2256 SCIP_Real* branchcandsfrac, /**< fractional part of the branching candidates */
2257 int nbranchcands, /**< number of branching candidates */
2258 SCIP_Bool executebranching, /**< perform a branching step after probing */
2259 SCIP_RESULT* result /**< pointer to the result of the execution */
2260 )
2261{
2262 SCIP_BRANCHRULE* branchrule;
2263
2264 assert(scip != NULL);
2265 assert(result != NULL);
2266
2267 /* find branching rule */
2269 assert(branchrule != NULL);
2270
2271 /* execute branching rule */
2272 SCIP_CALL( execRelpscost(scip, branchrule, branchcands, branchcandssol, branchcandsfrac, NULL, nbranchcands, executebranching, result) );
2273
2274 return SCIP_OKAY;
2275}
static long bound
#define BRANCHRULE_DESC
#define DEFAULT_DEGENERACYAWARE
#define DEFAULT_USESBLOCALINFO
#define DEFAULT_SKIPBADINITCANDS
#define DEFAULT_SBITERQUOT
static SCIP_Real calcNlscore(SCIP *scip, int *nlcount, int nlcountmax, int probindex)
#define DEFAULT_HIGHERRORTOL
#define BRANCHRULE_PRIORITY
#define DEFAULT_PROBINGBOUNDS
#define DEFAULT_INITITER
#define DEFAULT_USESMALLWEIGHTSITLIM
static SCIP_RETCODE applyBdchgs(SCIP *scip, SCIP_VAR **vars, int *bdchginds, SCIP_BOUNDTYPE *bdchgtypes, SCIP_Real *bdchgbounds, int nbdchgs, SCIP_RESULT *result)
#define DEFAULT_STARTRANDSEED
#define DEFAULT_FILTERCANDSSYM
#define DEFAULT_USEDYNAMICCONFIDENCE
#define DEFAULT_DYNAMICWEIGHTS
#define DEFAULT_MAXBDCHGS
#define DEFAULT_USEHYPTESTFORRELIABILITY
#define DEFAULT_GMIAVGEFFWEIGHT
static SCIP_RETCODE SCIPupdateVarPseudocostSymmetric(SCIP *scip, SCIP_BRANCHRULEDATA *branchruledata, SCIP_VAR *branchvar, int *branchorbitidx, int branchvaridx, SCIP_Real solvaldelta, SCIP_Real objdelta, SCIP_Real weight)
#define DEFAULT_INFERENCEWEIGHT
#define BRANCHRULE_NAME
#define DEFAULT_RANDINITORDER
static SCIP_RETCODE initOrbits(SCIP *scip, SCIP_BRANCHRULEDATA *branchruledata)
#define DEFAULT_MAXLOOKAHEAD
#define DEFAULT_SBITEROFS
static SCIP_DECL_BRANCHFREE(branchFreeRelpscost)
#define DEFAULT_NLSCOREWEIGHT
#define DEFAULT_CONFLICTWEIGHT
static SCIP_RETCODE filterSymmetricVariables(SCIP *scip, SCIP_BRANCHRULEDATA *branchruledata, SCIP_VAR **origbranchcands, SCIP_Real *origbranchcandssol, SCIP_Real *origbranchcandsfrac, int norigbranchcands, SCIP_VAR **branchcands, SCIP_Real *branchcandssol, SCIP_Real *branchcandsfrac, int *branchorbitidx, int *nbranchcands)
#define DEFAULT_GMILASTEFFWEIGHT
#define DEFAULT_USERELERRORFORRELIABILITY
#define DEFAULT_CUTOFFWEIGHT
static SCIP_DECL_BRANCHCOPY(branchCopyRelpscost)
static SCIP_RETCODE countNonlinearities(SCIP *scip, int *nlcount, int nlcountsize, int *nlcountmax)
static SCIP_DECL_BRANCHINITSOL(branchInitsolRelpscost)
#define DEFAULT_TRANSSYMPSCOST
static SCIP_RETCODE addBdchg(SCIP *scip, int **bdchginds, SCIP_BOUNDTYPE **bdchgtypes, SCIP_Real **bdchgbounds, int *nbdchgs, int ind, SCIP_BOUNDTYPE type, SCIP_Real bound)
#define DEFAULT_LOWERRORTOL
static SCIP_RETCODE execRelpscost(SCIP *scip, SCIP_BRANCHRULE *branchrule, SCIP_VAR **branchcands, SCIP_Real *branchcandssol, SCIP_Real *branchcandsfrac, int *branchorbitidx, int nbranchcands, SCIP_Bool executebranch, SCIP_RESULT *result)
static SCIP_RETCODE branchruledataEnsureNlcount(SCIP *scip, SCIP_BRANCHRULEDATA *branchruledata)
#define DEFAULT_CONFLENGTHWEIGHT
#define DEFAULT_CONFIDENCELEVEL
#define DEFAULT_MINRELIABLE
static SCIP_Real calcScore(SCIP *scip, SCIP_BRANCHRULEDATA *branchruledata, SCIP_Real conflictscore, SCIP_Real avgconflictscore, SCIP_Real conflengthscore, SCIP_Real avgconflengthscore, SCIP_Real inferencescore, SCIP_Real avginferencescore, SCIP_Real cutoffscore, SCIP_Real avgcutoffscore, SCIP_Real gmieffscore, SCIP_Real lastgmieffscore, SCIP_Real pscostscore, SCIP_Real avgpscostscore, SCIP_Real nlscore, SCIP_Real frac, SCIP_Real degeneracyfactor)
#define DEFAULT_STORESEMIINITCOSTS
static SCIP_RETCODE binvarGetActiveProbindex(SCIP *scip, SCIP_VAR *var, int *probindex)
#define DEFAULT_MAXPROPROUNDS
static SCIP_DECL_BRANCHEXITSOL(branchExitsolRelpscost)
#define DEFAULT_INITCAND
#define BRANCHRULE_MAXDEPTH
#define DEFAULT_MAXRELIABLE
#define DEFAULT_PSCOSTWEIGHT
#define BRANCHRULE_MAXBOUNDDIST
static void freeBdchgs(SCIP *scip, int **bdchginds, SCIP_BOUNDTYPE **bdchgtypes, SCIP_Real **bdchgbounds, int *nbdchgs)
static SCIP_DECL_BRANCHEXECLP(branchExeclpRelpscost)
reliable pseudo costs branching rule
Constraint handler for AND constraints, .
#define NULL
Definition: def.h:267
#define SCIP_Longint
Definition: def.h:158
#define SCIP_UNUSED(x)
Definition: def.h:428
#define SCIP_REAL_MAX
Definition: def.h:174
#define SCIP_INVALID
Definition: def.h:193
#define SCIP_Bool
Definition: def.h:91
#define MIN(x, y)
Definition: def.h:243
#define SCIP_Real
Definition: def.h:173
#define TRUE
Definition: def.h:93
#define FALSE
Definition: def.h:94
#define MAX(x, y)
Definition: def.h:239
#define SCIP_LONGINT_FORMAT
Definition: def.h:165
#define SCIP_REAL_MIN
Definition: def.h:175
#define SCIP_CALL(x)
Definition: def.h:374
SCIP_RETCODE SCIPexecRelpscostBranching(SCIP *scip, SCIP_VAR **branchcands, SCIP_Real *branchcandssol, SCIP_Real *branchcandsfrac, int nbranchcands, SCIP_Bool executebranching, SCIP_RESULT *result)
SCIP_RETCODE SCIPincludeBranchruleRelpscost(SCIP *scip)
SCIP_VAR * SCIPgetResultantAnd(SCIP *scip, SCIP_CONS *cons)
Definition: cons_and.c:5248
int SCIPgetNVarsAnd(SCIP *scip, SCIP_CONS *cons)
Definition: cons_and.c:5199
SCIP_VAR ** SCIPgetVarsAnd(SCIP *scip, SCIP_CONS *cons)
Definition: cons_and.c:5223
int SCIPgetSubscipDepth(SCIP *scip)
Definition: scip_copy.c:2605
SCIP_Bool SCIPisStopped(SCIP *scip)
Definition: scip_general.c:724
SCIP_Bool SCIPisExactSolve(SCIP *scip)
Definition: scip_general.c:611
SCIP_RETCODE SCIPgetVarsData(SCIP *scip, SCIP_VAR ***vars, int *nvars, int *nbinvars, int *nintvars, int *nimplvars, int *ncontvars)
Definition: scip_prob.c:1866
int SCIPgetNVars(SCIP *scip)
Definition: scip_prob.c:1992
SCIP_VAR ** SCIPgetVars(SCIP *scip)
Definition: scip_prob.c:1947
int SCIPhashmapGetImageInt(SCIP_HASHMAP *hashmap, void *origin)
Definition: misc.c:3281
SCIP_RETCODE SCIPupdateNodeLowerbound(SCIP *scip, SCIP_NODE *node, SCIP_Real newbound)
Definition: scip_prob.c:3757
SCIP_Real SCIPgetNodeLowerbound(SCIP *scip, SCIP_NODE *node)
Definition: scip_prob.c:3622
void SCIPverbMessage(SCIP *scip, SCIP_VERBLEVEL msgverblevel, FILE *file, const char *formatstr,...)
Definition: scip_message.c:225
#define SCIPdebugMsg
Definition: scip_message.h:78
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 SCIPgetRealParam(SCIP *scip, const char *name, SCIP_Real *value)
Definition: scip_param.c:307
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
SCIP_RETCODE SCIPsetBranchruleExecLp(SCIP *scip, SCIP_BRANCHRULE *branchrule, SCIP_DECL_BRANCHEXECLP((*branchexeclp)))
Definition: scip_branch.c:249
SCIP_BRANCHRULE * SCIPfindBranchrule(SCIP *scip, const char *name)
Definition: scip_branch.c:297
SCIP_RETCODE SCIPsetBranchruleCopy(SCIP *scip, SCIP_BRANCHRULE *branchrule, SCIP_DECL_BRANCHCOPY((*branchcopy)))
Definition: scip_branch.c:153
SCIP_RETCODE SCIPincludeBranchruleBasic(SCIP *scip, SCIP_BRANCHRULE **branchruleptr, const char *name, const char *desc, int priority, int maxdepth, SCIP_Real maxbounddist, SCIP_BRANCHRULEDATA *branchruledata)
Definition: scip_branch.c:116
const char * SCIPbranchruleGetName(SCIP_BRANCHRULE *branchrule)
Definition: branch.c:1971
SCIP_BRANCHRULEDATA * SCIPbranchruleGetData(SCIP_BRANCHRULE *branchrule)
Definition: branch.c:1849
SCIP_RETCODE SCIPsetBranchruleFree(SCIP *scip, SCIP_BRANCHRULE *branchrule, SCIP_DECL_BRANCHFREE((*branchfree)))
Definition: scip_branch.c:169
SCIP_RETCODE SCIPsetBranchruleExitsol(SCIP *scip, SCIP_BRANCHRULE *branchrule, SCIP_DECL_BRANCHEXITSOL((*branchexitsol)))
Definition: scip_branch.c:233
void SCIPbranchruleSetData(SCIP_BRANCHRULE *branchrule, SCIP_BRANCHRULEDATA *branchruledata)
Definition: branch.c:1859
SCIP_RETCODE SCIPsetBranchruleInitsol(SCIP *scip, SCIP_BRANCHRULE *branchrule, SCIP_DECL_BRANCHINITSOL((*branchinitsol)))
Definition: scip_branch.c:217
SCIP_RETCODE SCIPbranchVarVal(SCIP *scip, SCIP_VAR *var, SCIP_Real val, SCIP_NODE **downchild, SCIP_NODE **eqchild, SCIP_NODE **upchild)
Definition: scip_branch.c:1126
SCIP_RETCODE SCIPgetLPBranchCands(SCIP *scip, SCIP_VAR ***lpcands, SCIP_Real **lpcandssol, SCIP_Real **lpcandsfrac, int *nlpcands, int *npriolpcands, int *nfracimplvars)
Definition: scip_branch.c:395
SCIP_Real SCIPgetBranchScore(SCIP *scip, SCIP_VAR *var, SCIP_Real downgain, SCIP_Real upgain)
Definition: scip_branch.c:849
SCIP_CONSHDLR * SCIPfindConshdlr(SCIP *scip, const char *name)
Definition: scip_cons.c:941
int SCIPconshdlrGetNActiveConss(SCIP_CONSHDLR *conshdlr)
Definition: cons.c:4670
SCIP_CONS ** SCIPconshdlrGetConss(SCIP_CONSHDLR *conshdlr)
Definition: cons.c:4593
SCIP_RETCODE SCIPgetLPDualDegeneracy(SCIP *scip, SCIP_Real *degeneracy, SCIP_Real *varconsratio)
Definition: scip_lp.c:2792
SCIP_Bool SCIPinDive(SCIP *scip)
Definition: scip_lp.c:2775
SCIP_Bool SCIPhasCurrentNodeLP(SCIP *scip)
Definition: scip_lp.c:83
SCIP_LPSOLSTAT SCIPgetLPSolstat(SCIP *scip)
Definition: scip_lp.c:168
SCIP_Bool SCIPallColsInLP(SCIP *scip)
Definition: scip_lp.c:649
SCIP_Real SCIPgetLPObjval(SCIP *scip)
Definition: scip_lp.c:247
SCIP_Bool SCIPisLPSolBasic(SCIP *scip)
Definition: scip_lp.c:667
#define SCIPfreeBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:110
#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
#define SCIPduplicateBufferArray(scip, ptr, source, num)
Definition: scip_mem.h:132
#define SCIPallocBlockMemoryArray(scip, ptr, num)
Definition: scip_mem.h:93
#define SCIPreallocBlockMemoryArray(scip, ptr, oldnum, newnum)
Definition: scip_mem.h:99
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip_mem.h:108
#define SCIPfreeBlockMemoryArrayNull(scip, ptr, num)
Definition: scip_mem.h:111
#define SCIPfreeBufferArrayNull(scip, ptr)
Definition: scip_mem.h:137
#define SCIPallocBlockMemory(scip, ptr)
Definition: scip_mem.h:89
SCIP_Bool SCIPisNLPConstructed(SCIP *scip)
Definition: scip_nlp.c:110
int SCIPgetNNLPVars(SCIP *scip)
Definition: scip_nlp.c:201
SCIP_RETCODE SCIPgetNLPVarsNonlinearity(SCIP *scip, int *nlcount)
Definition: scip_nlp.c:223
SCIP_Real SCIPnodeGetLowerbound(SCIP_NODE *node)
Definition: tree.c:7513
SCIP_Longint SCIPnodeGetNumber(SCIP_NODE *node)
Definition: tree.c:7493
SCIP_Bool SCIPinProbing(SCIP *scip)
Definition: scip_probing.c:97
SCIP_SOL * SCIPgetBestSol(SCIP *scip)
Definition: scip_sol.c:2169
int SCIPsolGetIndex(SCIP_SOL *sol)
Definition: sol.c:2835
SCIP_Real SCIPgetAvgInferenceScore(SCIP *scip)
SCIP_Real SCIPgetAvgConflictlengthScore(SCIP *scip)
SCIP_Longint SCIPgetNInfeasibleLeaves(SCIP *scip)
SCIP_Longint SCIPgetNNodes(SCIP *scip)
SCIP_Longint SCIPgetNNodeInitLPs(SCIP *scip)
SCIP_Longint SCIPgetNStrongbranchLPIterations(SCIP *scip)
SCIP_Longint SCIPgetNNodeLPIterations(SCIP *scip)
SCIP_Real SCIPgetAvgConflictScore(SCIP *scip)
SCIP_Longint SCIPgetNDualResolveLPIterations(SCIP *scip)
SCIP_Real SCIPgetAvgPseudocostScore(SCIP *scip)
SCIP_Real SCIPgetCutoffbound(SCIP *scip)
SCIP_Longint SCIPgetNObjlimLeaves(SCIP *scip)
SCIP_Longint SCIPgetNNodeInitLPIterations(SCIP *scip)
SCIP_Longint SCIPgetNDualResolveLPs(SCIP *scip)
SCIP_Longint SCIPgetNRootStrongbranchLPIterations(SCIP *scip)
SCIP_Real SCIPgetAvgCutoffScore(SCIP *scip)
SCIP_RETCODE SCIPcomputeOrbitsComponentsSym(SCIP *scip, int npermvars, int **permstrans, int nperms, int *components, int *componentbegins, int *vartocomponent, int ncomponents, int *orbits, int *orbitbegins, int *norbits, int *varorbitmap)
Definition: symmetry.c:420
SCIP_Real SCIPinfinity(SCIP *scip)
SCIP_Bool SCIPisGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Real SCIPfeasFrac(SCIP *scip, SCIP_Real val)
SCIP_Real SCIPfeasCeil(SCIP *scip, SCIP_Real val)
SCIP_Real SCIPfeasFloor(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisFeasIntegral(SCIP *scip, SCIP_Real val)
SCIP_Real SCIPfeastol(SCIP *scip)
SCIP_Bool SCIPisSumGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPisGT(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)
SCIP_Bool SCIPisSumGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
int SCIPgetDepth(SCIP *scip)
Definition: scip_tree.c:670
SCIP_NODE * SCIPgetCurrentNode(SCIP *scip)
Definition: scip_tree.c:91
SCIP_RETCODE SCIPgetVarStrongbranchFrac(SCIP *scip, SCIP_VAR *var, int itlim, SCIP_Bool idempotent, SCIP_Real *down, SCIP_Real *up, SCIP_Bool *downvalid, SCIP_Bool *upvalid, SCIP_Bool *downinf, SCIP_Bool *upinf, SCIP_Bool *downconflict, SCIP_Bool *upconflict, SCIP_Bool *lperror)
Definition: scip_var.c:2919
SCIP_RETCODE SCIPtightenVarLb(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound, SCIP_Bool force, SCIP_Bool *infeasible, SCIP_Bool *tightened)
Definition: scip_var.c:5203
SCIP_Bool SCIPpscostThresholdProbabilityTest(SCIP *scip, SCIP_VAR *var, SCIP_Real frac, SCIP_Real threshold, SCIP_BRANCHDIR dir, SCIP_CONFIDENCELEVEL clevel)
Definition: scip_var.c:9059
SCIP_Bool SCIPvarIsActive(SCIP_VAR *var)
Definition: var.c:17748
SCIP_Bool SCIPvarIsBinary(SCIP_VAR *var)
Definition: var.c:17599
SCIP_Real SCIPgetVarPseudocostCountCurrentRun(SCIP *scip, SCIP_VAR *var, SCIP_BRANCHDIR dir)
Definition: scip_var.c:8950
SCIP_Real SCIPgetVarAvgInferenceScore(SCIP *scip, SCIP_VAR *var)
Definition: scip_var.c:9473
SCIP_RETCODE SCIPendStrongbranch(SCIP *scip)
Definition: scip_var.c:2744
SCIP_VARSTATUS SCIPvarGetStatus(SCIP_VAR *var)
Definition: var.c:17538
SCIP_Real SCIPgetVarPseudocostCurrentRun(SCIP *scip, SCIP_VAR *var, SCIP_BRANCHDIR dir)
Definition: scip_var.c:8896
SCIP_Real SCIPvarGetUbLocal(SCIP_VAR *var)
Definition: var.c:18144
SCIP_Bool SCIPisVarPscostRelerrorReliable(SCIP *scip, SCIP_VAR *var, SCIP_Real threshold, SCIP_CONFIDENCELEVEL clevel)
Definition: scip_var.c:9078
SCIP_RETCODE SCIPtightenVarUb(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound, SCIP_Bool force, SCIP_Bool *infeasible, SCIP_Bool *tightened)
Definition: scip_var.c:5320
int SCIPvarGetProbindex(SCIP_VAR *var)
Definition: var.c:17768
const char * SCIPvarGetName(SCIP_VAR *var)
Definition: var.c:17419
SCIP_Longint SCIPgetVarStrongbranchNode(SCIP *scip, SCIP_VAR *var)
Definition: scip_var.c:4160
SCIP_LPSOLSTAT SCIPgetLastStrongbranchLPSolStat(SCIP *scip, SCIP_BRANCHDIR branchdir)
Definition: scip_var.c:3988
SCIP_Real SCIPgetVarPseudocostVal(SCIP *scip, SCIP_VAR *var, SCIP_Real solvaldelta)
Definition: scip_var.c:8814
SCIP_Real SCIPgetVarConflictlengthScore(SCIP *scip, SCIP_VAR *var)
Definition: scip_var.c:9303
SCIP_Bool SCIPsignificantVarPscostDifference(SCIP *scip, SCIP_VAR *varx, SCIP_Real fracx, SCIP_VAR *vary, SCIP_Real fracy, SCIP_BRANCHDIR dir, SCIP_CONFIDENCELEVEL clevel, SCIP_Bool onesided)
Definition: scip_var.c:9029
SCIP_RETCODE SCIPgetVarStrongbranchWithPropagation(SCIP *scip, SCIP_VAR *var, SCIP_Real solval, SCIP_Real lpobjval, int itlim, int maxproprounds, SCIP_Real *down, SCIP_Real *up, SCIP_Bool *downvalid, SCIP_Bool *upvalid, SCIP_Longint *ndomredsdown, SCIP_Longint *ndomredsup, SCIP_Bool *downinf, SCIP_Bool *upinf, SCIP_Bool *downconflict, SCIP_Bool *upconflict, SCIP_Bool *lperror, SCIP_Real *newlbs, SCIP_Real *newubs)
Definition: scip_var.c:3352
SCIP_Real SCIPvarGetLPSol(SCIP_VAR *var)
Definition: var.c:18452
SCIP_Real SCIPgetVarPseudocostScoreCurrentRun(SCIP *scip, SCIP_VAR *var, SCIP_Real solval)
Definition: scip_var.c:9141
SCIP_Real SCIPvarGetLbLocal(SCIP_VAR *var)
Definition: var.c:18134
SCIP_Bool SCIPvarIsNegated(SCIP_VAR *var)
Definition: var.c:17574
SCIP_VAR * SCIPvarGetNegationVar(SCIP_VAR *var)
Definition: var.c:17904
SCIP_RETCODE SCIPupdateVarPseudocost(SCIP *scip, SCIP_VAR *var, SCIP_Real solvaldelta, SCIP_Real objdelta, SCIP_Real weight)
Definition: scip_var.c:8780
SCIP_Real SCIPgetVarAvgGMIScore(SCIP *scip, SCIP_VAR *var)
Definition: scip_var.c:9877
SCIP_Real SCIPgetVarLastGMIScore(SCIP *scip, SCIP_VAR *var)
Definition: scip_var.c:9932
SCIP_Real SCIPgetVarConflictScore(SCIP *scip, SCIP_VAR *var)
Definition: scip_var.c:9241
SCIP_Real SCIPgetVarAvgCutoffScore(SCIP *scip, SCIP_VAR *var)
Definition: scip_var.c:9727
SCIP_RETCODE SCIPgetBinvarRepresentative(SCIP *scip, SCIP_VAR *var, SCIP_VAR **repvar, SCIP_Bool *negated)
Definition: scip_var.c:1597
SCIP_RETCODE SCIPgetVarStrongbranchLast(SCIP *scip, SCIP_VAR *var, SCIP_Real *down, SCIP_Real *up, SCIP_Bool *downvalid, SCIP_Bool *upvalid, SCIP_Real *solval, SCIP_Real *lpobjval)
Definition: scip_var.c:4010
SCIP_Real SCIPgetVarPseudocostScore(SCIP *scip, SCIP_VAR *var, SCIP_Real solval)
Definition: scip_var.c:9103
SCIP_RETCODE SCIPstartStrongbranch(SCIP *scip, SCIP_Bool enablepropagation)
Definition: scip_var.c:2686
void SCIPfreeRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen)
SCIP_Real SCIPrandomGetReal(SCIP_RANDNUMGEN *randnumgen, SCIP_Real minrandval, SCIP_Real maxrandval)
Definition: misc.c:10130
SCIP_RETCODE SCIPcreateRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen, unsigned int initialseed, SCIP_Bool useglobalseed)
memory allocation routines
#define BMSclearMemoryArray(ptr, num)
Definition: memory.h:130
SCIP_RETCODE SCIPgetSymmetry(SCIP *scip, int *npermvars, SCIP_VAR ***permvars, SCIP_HASHMAP **permvarmap, int *nperms, int ***perms, int ***permstrans, SCIP_Real *log10groupsize, SCIP_Bool *binvaraffected, int **components, int **componentbegins, int **vartocomponent, int *ncomponents)
propagator for symmetry handling
public methods for branching rules
public methods for managing constraints
public methods for message output
#define SCIPdebug(x)
Definition: pub_message.h:93
public data structures and miscellaneous methods
public methods for primal CIP solutions
public methods for branch and bound tree
public methods for problem variables
public methods for branching rule plugins and branching
public methods for constraint handler plugins and constraints
general public methods
public methods for the LP relaxation, rows and columns
public methods for memory management
public methods for message handling
public methods for nonlinear relaxation
public methods for numerical tolerances
public methods for SCIP parameter handling
public methods for global and local (sub)problems
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
methods for handling symmetries
SCIP_RETCODE SCIPtreemodelSelectCandidate(SCIP *scip, SCIP_TREEMODEL *treemodel, SCIP_VAR **branchcands, SCIP_Real *mingains, SCIP_Real *maxgains, SCIP_Real *tiebreakerscore, int nbranchcands, int *bestcand)
Definition: treemodel.c:912
SCIP_RETCODE SCIPtreemodelInit(SCIP *scip, SCIP_TREEMODEL **treemodel)
Definition: treemodel.c:826
SCIP_RETCODE SCIPtreemodelFree(SCIP *scip, SCIP_TREEMODEL **treemodel)
Definition: treemodel.c:884
SCIP_Bool SCIPtreemodelIsEnabled(SCIP *scip, SCIP_TREEMODEL *treemodel)
Definition: treemodel.c:900
Branching rules based on the Single-Variable-Branching (SVB) model.
struct SCIP_BranchruleData SCIP_BRANCHRULEDATA
Definition: type_branch.h:57
@ SCIP_BRANCHDIR_DOWNWARDS
Definition: type_history.h:43
@ SCIP_BRANCHDIR_UPWARDS
Definition: type_history.h:44
@ SCIP_BOUNDTYPE_UPPER
Definition: type_lp.h:57
@ SCIP_BOUNDTYPE_LOWER
Definition: type_lp.h:56
enum SCIP_BoundType SCIP_BOUNDTYPE
Definition: type_lp.h:59
@ SCIP_LPSOLSTAT_OPTIMAL
Definition: type_lp.h:43
@ SCIP_LPSOLSTAT_INFEASIBLE
Definition: type_lp.h:44
@ SCIP_LPSOLSTAT_OBJLIMIT
Definition: type_lp.h:46
@ SCIP_LPSOLSTAT_ITERLIMIT
Definition: type_lp.h:47
@ SCIP_VERBLEVEL_HIGH
Definition: type_message.h:56
@ SCIP_CONFIDENCELEVEL_MAX
Definition: type_misc.h:51
@ SCIP_CONFIDENCELEVEL_MEDIUM
Definition: type_misc.h:49
@ SCIP_CONFIDENCELEVEL_HIGH
Definition: type_misc.h:50
@ SCIP_CONFIDENCELEVEL_MIN
Definition: type_misc.h:47
@ SCIP_CONFIDENCELEVEL_LOW
Definition: type_misc.h:48
enum SCIP_Confidencelevel SCIP_CONFIDENCELEVEL
Definition: type_misc.h:53
@ SCIP_DIDNOTRUN
Definition: type_result.h:42
@ SCIP_CUTOFF
Definition: type_result.h:48
@ SCIP_REDUCEDDOM
Definition: type_result.h:51
@ SCIP_CONSADDED
Definition: type_result.h:52
@ SCIP_BRANCHED
Definition: type_result.h:54
enum SCIP_Result SCIP_RESULT
Definition: type_result.h:61
@ SCIP_OKAY
Definition: type_retcode.h:42
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:63
@ SCIP_VARSTATUS_FIXED
Definition: type_var.h:52