presol_dualinfer.c
Go to the documentation of this file.
30 * This presolver does bound strengthening on continuous variables (columns) for getting bounds on dual variables y.
31 * The bounds of the dual variables are then used to fix primal variables or change the side of constraints.
42 * Chen et al. "Two-row and two-column mixed-integer presolve using hasing-based pairing methods".
45/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
72#define PRESOL_PRIORITY (-3000) /**< priority of the presolver (>= 0: before, < 0: after constraint handlers) */
73#define PRESOL_MAXROUNDS 0 /**< maximal number of presolving rounds the presolver participates in (-1: no limit) */
74#define PRESOL_TIMING SCIP_PRESOLTIMING_EXHAUSTIVE /* timing of the presolver (fast, medium, or exhaustive) */
76#define DEFAULT_TWOCOLUMN_COMBINE TRUE /**< should two column convex combination be used per default */
77#define DEFAULT_MAXLOOPS_DUALBNDSTR 12 /**< default maximal number of loops for dual bound strengthening */
78#define DEFAULT_MAXCONSIDEREDNONZEROS 100 /**< default maximal number of considered non-zeros within one row */
79#define DEFAULT_MAXRETRIEVEFAILS 1000 /**< default maximal number of consecutive useless hashtable retrieves */
80#define DEFAULT_MAXCOMBINEFAILS 1000 /**< default maximal number of consecutive useless row combines */
81#define DEFAULT_MAXHASHFAC 10 /**< default maximal number of hashlist entries as multiple of number of rows in the problem */
82#define DEFAULT_MAXPAIRFAC 1 /**< default maximal number of processed row pairs as multiple of the number of rows in the problem */
83#define DEFAULT_MAXROWSUPPORT 3 /**< default maximal number of non-zeros in one row for turning an inequality into an equality */
95 int maxpairfac; /**< maximal number of processed row pairs as multiple of the number of rows in the problem (-1: no limit) */
96 int maxhashfac; /**< maximal number of hashlist entries as multiple of number of rows in the problem (-1: no limit) */
99 int maxconsiderednonzeros; /**< maximal number of considered non-zeros within one row (-1: no limit) */
100 int maxrowsupport; /**< maximal number of non-zeros in one row for turning an inequality into an equality */
109};
118};
121/** Signum for convex-combined variable coefficients \f$(\lambda * A_{ri} + (1 - \lambda) * A_{si})\f$
223 * tries to derive upper and lower bounds for all variables via convex combinations of linear inequalities
315 /* We use 2.0 as default value for "no cancellation". For cancellations, this will be replaced by values in (0,1).
316 * A value larger than 1.0 is used because we sort the array and want to put non-cancellations to the end. */
409 /* The obvious preconditions for bound tightenings are met, so we try to calculate new bounds. */
478 SCIPdebugMsg(scip, "lambda = 0, l1 = %g, l2 = %g, ninfs = %d\n", i, breakpoints[i], l1, l2, ninfs);
580 SCIPdebugMsg(scip, "lambda_%d = %g, l1 = %g, l2 = %g, ninfs = %d\n", i, breakpoints[i], l1, l2, ninfs);
616 newlbs[j] = MAX(newlbs[j], (breakpoints[i] * l1 + (1 - breakpoints[i]) * l2 + coef * ubs[idx]) / coef);
623 newubs[j] = MIN(newubs[j], (breakpoints[i] * l1 + (1 - breakpoints[i]) * l2 + coef * lbs[idx]) / coef);
652 SCIPdebugMsg(scip, "lambda = 1, l1 = %g, l2 = %g, ninfs = %d\n", i, breakpoints[i], l1, l2, ninfs);
856 *rowub = (rhs - minresactivity) / val; // maybe one wants some kind of numerical guard of check that values is not too small for all these
998/** In the primal the residual activity of a constraint w.r.t. a variable is the activity of the constraint without the variable.
1196 SCIP_Bool ubinfchange, /**< flag indicating if the upper bound has changed from infinity to a finite value */
1197 SCIP_Bool lbinfchange /**< flag indicating if the lower bound has changed from -infinity to a finite value */
1246/** use LP calculations for determining the best dual variable bounds from a specific row index */
1297 if( SCIPvarGetType(var) == SCIP_VARTYPE_CONTINUOUS || SCIPvarGetType(var) == SCIP_VARTYPE_IMPLINT )
1360 if( SCIPvarGetType(var) == SCIP_VARTYPE_CONTINUOUS || SCIPvarGetType(var) == SCIP_VARTYPE_IMPLINT )
1556 FIXINGDIRECTION* varstofix, /**< array holding information for later upper/lower bound fixing */
1645 (SCIPvarGetType(var) != SCIP_VARTYPE_CONTINUOUS && SCIPvarGetType(var) != SCIP_VARTYPE_IMPLINT) )
1695 /* run convex combination on pairs of continuous variables (columns) using Belotti's algorithm */
1716 maxhashes = presoldata->maxhashfac == -1 ? SCIP_LONGINT_MAX : (((SCIP_Longint)ncols) * presoldata->maxhashfac);
1725 maxlen = MIN(presoldata->maxconsiderednonzeros, SCIPmatrixGetColNNonzs(matrix, implubvars[i])); /*lint !e666*/
1760 SCIPdebugMsg(scip, "hashlist sizes: pp %d, mm %d, pm %d, mp %d \n", pospp, posmm, pospm, posmp);
1786 maxcombines = presoldata->maxpairfac == -1 ? SCIP_LONGINT_MAX : (((SCIP_Longint)ncols) * presoldata->maxpairfac);
1797 // same as in the rworowbnd presolver - both while loops to basically the same with one using pp and mm and the other pm and mp
1838 SCIPdebugMsg(scip, "pm/mp: %d retrievefails before reset, %d combines\n", retrievefails, ncombines);
1888 maxcombines = presoldata->maxpairfac == -1 ? SCIP_LONGINT_MAX : (((SCIP_Longint)ncols) * presoldata->maxpairfac);
2209 /* the reductions made in this presolver apply to all optimal solutions because of complementary slackness */
2217 SCIPdebugMsg(scip, "DualInfer not executed because condition of existing dual solution is not fulfilled.\n");
2269 assert(varstofix[i] == NOFIX || (!SCIPmatrixUplockConflict(matrix, i) && !SCIPmatrixDownlockConflict(matrix, i)));
2311 if( SCIPvarGetType(var) == SCIP_VARTYPE_CONTINUOUS || SCIPvarGetType(var) == SCIP_VARTYPE_IMPLINT )
2428 SCIP_CALL( SCIPincludePresolBasic(scip, &presol, PRESOL_NAME, PRESOL_DESC, PRESOL_PRIORITY, PRESOL_MAXROUNDS,
2446 &presoldata->maxconsiderednonzeros, TRUE, DEFAULT_MAXCONSIDEREDNONZEROS, -1, INT_MAX, NULL, NULL) );
2460 "Maximum number of hashlist entries as multiple of number of columns in the problem (-1: no limit)",
2465 "Maximum number of processed column pairs as multiple of the number of columns in the problem (-1: no limit)",
Constraint handler for linear constraints in their most general form, .
static SCIP_RETCODE determineBestBounds(SCIP *scip, SCIP_VAR *var, SCIP_SOL *sol, SCIP_Real boundswitch, int usevbds, SCIP_Bool allowlocal, SCIP_Bool fixintegralrhs, SCIP_Bool ignoresol, int *boundsfortrans, SCIP_BOUNDTYPE *boundtypesfortrans, SCIP_Real *bestlb, SCIP_Real *bestub, int *bestlbtype, int *bestubtype, SCIP_BOUNDTYPE *selectedbound, SCIP_Bool *freevariable)
Definition: cuts.c:2718
SCIP_Real SCIPgetRhsLinear(SCIP *scip, SCIP_CONS *cons)
Definition: cons_linear.c:18490
SCIP_RETCODE SCIPchgRhsLinear(SCIP *scip, SCIP_CONS *cons, SCIP_Real rhs)
Definition: cons_linear.c:18535
SCIP_Real SCIPgetLhsLinear(SCIP *scip, SCIP_CONS *cons)
Definition: cons_linear.c:18466
SCIP_RETCODE SCIPcreateConsBasicLinear(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_Real *vals, SCIP_Real lhs, SCIP_Real rhs)
Definition: cons_linear.c:18066
SCIP_RETCODE SCIPchgLhsLinear(SCIP *scip, SCIP_CONS *cons, SCIP_Real lhs)
Definition: cons_linear.c:18514
SCIP_RETCODE SCIPwriteOrigProblem(SCIP *scip, const char *filename, const char *extension, SCIP_Bool genericnames)
Definition: scip_prob.c:601
SCIP_RETCODE SCIPsetObjsense(SCIP *scip, SCIP_OBJSENSE objsense)
Definition: scip_prob.c:1242
SCIP_RETCODE SCIPcreateProbBasic(SCIP *scip, const char *name)
Definition: scip_prob.c:180
void SCIPhashsetFree(SCIP_HASHSET **hashset, BMS_BLKMEM *blkmem)
Definition: misc.c:3793
SCIP_Bool SCIPhashsetExists(SCIP_HASHSET *hashset, void *element)
Definition: misc.c:3820
SCIP_RETCODE SCIPhashsetInsert(SCIP_HASHSET *hashset, BMS_BLKMEM *blkmem, void *element)
Definition: misc.c:3803
SCIP_RETCODE SCIPhashsetCreate(SCIP_HASHSET **hashset, BMS_BLKMEM *blkmem, int size)
Definition: misc.c:3762
void SCIPwarningMessage(SCIP *scip, const char *formatstr,...)
Definition: scip_message.c:120
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 SCIPsetIntParam(SCIP *scip, const char *name, int value)
Definition: scip_param.c:487
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 SCIPsetBoolParam(SCIP *scip, const char *name, SCIP_Bool value)
Definition: scip_param.c:429
SCIP_RETCODE SCIPincludePresolDualinfer(SCIP *scip)
Definition: presol_dualinfer.c:2417
SCIP_RETCODE SCIPreleaseCons(SCIP *scip, SCIP_CONS **cons)
Definition: scip_cons.c:1174
#define SCIPreallocBlockMemoryArray(scip, ptr, oldnum, newnum)
Definition: scip_mem.h:99
void SCIPpresolSetData(SCIP_PRESOL *presol, SCIP_PRESOLDATA *presoldata)
Definition: presol.c:522
SCIP_RETCODE SCIPsetPresolFree(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLFREE((*presolfree)))
Definition: scip_presol.c:156
SCIP_RETCODE SCIPsetPresolCopy(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLCOPY((*presolcopy)))
Definition: scip_presol.c:140
SCIP_RETCODE SCIPincludePresolBasic(SCIP *scip, SCIP_PRESOL **presolptr, const char *name, const char *desc, int priority, int maxrounds, SCIP_PRESOLTIMING timing, SCIP_DECL_PRESOLEXEC((*presolexec)), SCIP_PRESOLDATA *presoldata)
Definition: scip_presol.c:105
SCIP_RETCODE SCIPcheckSolOrig(SCIP *scip, SCIP_SOL *sol, SCIP_Bool *feasible, SCIP_Bool printreason, SCIP_Bool completely)
Definition: scip_sol.c:3305
SCIP_Bool SCIPisGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:497
SCIP_Bool SCIPisLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:471
SCIP_Bool SCIPisGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:484
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:445
SCIP_Bool SCIPisLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:458
SCIP_RETCODE SCIPfixVar(SCIP *scip, SCIP_VAR *var, SCIP_Real fixedval, SCIP_Bool *infeasible, SCIP_Bool *fixed)
Definition: scip_var.c:8399
SCIP_RETCODE SCIPcreateVarBasic(SCIP *scip, SCIP_VAR **var, const char *name, SCIP_Real lb, SCIP_Real ub, SCIP_Real obj, SCIP_VARTYPE vartype)
Definition: scip_var.c:194
SCIP_RETCODE SCIPchgVarObj(SCIP *scip, SCIP_VAR *var, SCIP_Real newobj)
Definition: scip_var.c:4636
void SCIPsortIntInt(int *intarray1, int *intarray2, int len)
void SCIPsortIntReal(int *intarray, SCIP_Real *realarray, int len)
void SCIPsortRealInt(SCIP_Real *realarray, int *intarray, int len)
static SCIP_RETCODE solveLP(SCIP *scip, SCIP_DIVESET *diveset, SCIP_Longint maxnlpiterations, SCIP_DIVECONTEXT divecontext, SCIP_Bool *lperror, SCIP_Bool *cutoff)
Definition: heuristics.c:48
SCIP_Bool SCIPmatrixUplockConflict(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1872
int * SCIPmatrixGetColIdxPtr(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1580
int SCIPmatrixGetRowNNonzs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1708
int SCIPmatrixGetColNNonzs(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1592
SCIP_Bool SCIPmatrixIsRowRhsInfinity(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1766
SCIP_Real SCIPmatrixGetRowLhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1742
SCIP_Real * SCIPmatrixGetRowValPtr(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1684
SCIP_Bool SCIPmatrixDownlockConflict(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1884
SCIP_Real SCIPmatrixGetRowRhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1754
SCIP_Real * SCIPmatrixGetColValPtr(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1568
SCIP_RETCODE SCIPmatrixCreate(SCIP *scip, SCIP_MATRIX **matrixptr, SCIP_Bool onlyifcomplete, SCIP_Bool *initialized, SCIP_Bool *complete, SCIP_Bool *infeasible, int *naddconss, int *ndelconss, int *nchgcoefs, int *nchgbds, int *nfixedvars)
Definition: matrix.c:454
SCIP_CONS * SCIPmatrixGetCons(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1860
int * SCIPmatrixGetRowIdxPtr(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1696
memory allocation routines
Definition: objbenders.h:44
static void calcMaxColActivity(SCIP *scip, SCIP_MATRIX *matrix, int col, SCIP_Real *lbdual, SCIP_Real *ubdual, SCIP_Real *maxcolact, int *maxcolactinf)
Definition: presol_dualinfer.c:1126
static void calcMinColActResidual(SCIP *scip, SCIP_MATRIX *matrix, int col, int row, SCIP_Real val, SCIP_Real *lbdual, SCIP_Real *ubdual, const SCIP_Real *mincolact, const int *mincolactinf, SCIP_Real *mincolresact)
Definition: presol_dualinfer.c:1003
static void getImpliedBounds(SCIP *scip, SCIP_MATRIX *matrix, int col, SCIP_Real *lbs, SCIP_Real *ubs, SCIP_Bool *ubimplied, SCIP_Bool *lbimplied)
Definition: presol_dualinfer.c:885
static void calcMinColActivity(SCIP *scip, SCIP_MATRIX *matrix, int col, SCIP_Real *lbdual, SCIP_Real *ubdual, SCIP_Real *mincolact, int *mincolactinf)
Definition: presol_dualinfer.c:1066
static SCIP_DECL_PRESOLFREE(presolFreeDualinfer)
Definition: presol_dualinfer.c:2167
static void getVarBoundsOfRow(SCIP *scip, SCIP_MATRIX *matrix, int col, int row, SCIP_Real val, SCIP_Real *lbs, SCIP_Real *ubs, SCIP_Real *rowub, SCIP_Bool *ubfound, SCIP_Real *rowlb, SCIP_Bool *lbfound)
Definition: presol_dualinfer.c:814
static SCIP_RETCODE combineCols(SCIP *scip, int *row1idxptr, int *row2idxptr, SCIP_Real *row1valptr, SCIP_Real *row2valptr, SCIP_Real b1, SCIP_Real b2, int row1len, int row2len, int ncols, SCIP_Bool swaprow1, SCIP_Bool swaprow2, SCIP_Real *lbs, SCIP_Real *ubs, SCIP_Bool *success)
Definition: presol_dualinfer.c:227
static SCIP_Real getMinColActWithoutRow(SCIP *scip, SCIP_MATRIX *matrix, int col, int withoutrow, SCIP_Real *lbdual, SCIP_Real *ubdual)
Definition: presol_dualinfer.c:947
static void updateDualBounds(SCIP *scip, SCIP_MATRIX *matrix, SCIP_Real objval, SCIP_Real val, int row, SCIP_Real mincolresact, SCIP_Real *lbdual, SCIP_Real *ubdual, int *boundchanges, SCIP_Bool *ubinfchange, SCIP_Bool *lbinfchange)
Definition: presol_dualinfer.c:1483
static void findNextBlock(const int *list, int len, int *start, int *end)
Definition: presol_dualinfer.c:205
static SCIP_RETCODE dualBoundStrengthening(SCIP *scip, SCIP_MATRIX *matrix, SCIP_PRESOLDATA *presoldata, FIXINGDIRECTION *varstofix, int *npossiblefixings, SIDECHANGE *sidestochange, int *npossiblesidechanges)
Definition: presol_dualinfer.c:1552
static SCIP_DECL_PRESOLEXEC(presolExecDualinfer)
Definition: presol_dualinfer.c:2183
static SCIP_DECL_PRESOLCOPY(presolCopyDualinfer)
Definition: presol_dualinfer.c:2153
static void getMinMaxActivityResiduals(SCIP *scip, SCIP_MATRIX *matrix, int withoutcol, int row, SCIP_Real *lbs, SCIP_Real *ubs, SCIP_Real *minresactivity, SCIP_Real *maxresactivity, SCIP_Bool *isminsettoinfinity, SCIP_Bool *ismaxsettoinfinity)
Definition: presol_dualinfer.c:719
static void infinityCountUpdate(SCIP *scip, SCIP_MATRIX *matrix, int row, SCIP_Real *lbdual, SCIP_Real *ubdual, const SCIP_Bool *isubimplied, SCIP_Real *mincolact, int *mincolactinf, SCIP_Bool ubinfchange, SCIP_Bool lbinfchange)
Definition: presol_dualinfer.c:1187
static SCIP_RETCODE addEntry(SCIP *scip, int *pos, int *listsize, int **hashlist, int **colidxlist, int hash, int colidx)
Definition: presol_dualinfer.c:173
dual inference presolver
public methods for managing constraints
public methods for matrix
public methods for message output
public methods for presolvers
public methods for problem variables
general public methods
public methods for memory management
public methods for message handling
public methods for numerical tolerances
public methods for presolving plugins
public methods for global and local (sub)problems
public methods for the probing mode
public methods for SCIP variables
SCIP_RETCODE SCIPincludeDefaultPlugins(SCIP *scip)
Definition: scipdefplugins.c:37
default SCIP plugins
Definition: presol_dualinfer.c:131
Definition: struct_cons.h:47
Definition: struct_cons.h:127
Definition: struct_misc.h:150
Definition: struct_matrix.h:48
Definition: struct_presol.h:47
Definition: struct_sol.h:74
Definition: struct_var.h:208
Definition: struct_scip.h:70