Lagromory separator.
This separator is based on the following article that discusses Lagromory separation using the relax-and-cut framework. Multiple enhancements have been implemented on top of the basic algorithm described in the article.
Fischetti M. and Salvagnin D. (2011).
A relax-and-cut framework for Gomory mixed-integer cuts.
Mathematical Programming Computation, 3, 79-102.
Consider the following linear relaxation at a node:
\[ \begin{array}{rrl} \min & c^T x &\\ & x & \in P, \end{array} \]
where \(P\) is the feasible region of the relaxation. Let the following be the cuts generated so far in the current separation round.
\[ {\alpha^i}^T x \leq \alpha^i_0, i = 1, 2, \ldots, M \]
Then, the following is the Lagrangian dual problem considered in the relax-and-cut framework used in the separator.
\[ z_D := \max\limits_{u \geq 0} \left\{L(u) := \min \left\{c^T x + \sum\limits_{i = 1}^{M} \left(u_i \left({\alpha^i}^T x - \alpha^i_0\right) \right) \mid x \in P\right\} \right\}, \]
where \(u\) are the Lagrangian multipliers (referred to as dualvector in this separator) used for penalizing the violation of the generated cuts, and \(z_D\) is the optimal objective value (which is approximated via ubparam in this separator). Then, the following are the steps of the relax-and-cut algorithm implemented in this separator.
\[ u^{j+1} = \left(u^j + \lambda_j s^k\right)_+, \]
where \(\lambda_j\) is the step length:\[ \lambda_j = \frac{\mu_j (UB - L(u^j))}{\|s^j\|^2_2}, \]
where \(mu_j\) is a factor (i.e., muparam) such that \(0 < \mu_j \leq 2\), UB isubparam, and \(s^j\) is the subgradient vector defined as: \[ s^j_k = \left({\alpha^k}^T x - \alpha^k_0\right), k = 1, 2, \ldots, M. \]
The factor \(mu_j\) is updated as below.\[ mu_j = \begin{cases} mubacktrackfactor * mu_j & \text{if } L(u^j) < bestLB - \delta\\ \begin{cases} muslab1factor * mu_j & \text{if } bestLB - avgLB < deltaslab1ub * delta\\ muslab2factor * mu_j & \text{if } deltaslab1ub * \delta \leq bestLB - avgLB < deltaslab2ub * delta\\ muslab3factor * mu_j & \text{otherwise} \end{cases} & \text{otherwise}, \end{cases} \]
where \(bestLB\) and \(avgLB\) are best and average Lagrangian values found so far, and \(\delta = UB - bestLB\).Definition in file sepa_lagromory.c.
Go to the source code of this file.
Functions | |
| static SCIP_RETCODE | createLPWithSoftCuts (SCIP *scip, SCIP_SEPADATA *sepadata) |
| static SCIP_RETCODE | deleteLPWithSoftCuts (SCIP *scip, SCIP_SEPADATA *sepadata) |
| static SCIP_RETCODE | createLPWithHardCuts (SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_ROW **cuts, int ncuts) |
| static SCIP_RETCODE | sepadataFree (SCIP *scip, SCIP_SEPADATA **sepadata) |
| static SCIP_RETCODE | updateMuSteplengthParam (SCIP *scip, SCIP_SEPADATA *sepadata, int subgradientiternum, SCIP_Real ubparam, SCIP_Real *lagrangianvals, SCIP_Real bestlagrangianval, SCIP_Real avglagrangianval, SCIP_Real *muparam, SCIP_Bool *backtrack) |
| static void | updateSubgradient (SCIP *scip, SCIP_SOL *sol, SCIP_ROW **cuts, int ncuts, SCIP_Real *subgradient, SCIP_Real *dualvector, SCIP_Bool *subgradientzero, int *ncutviols, SCIP_Real *maxcutviol, int *nnzsubgradientdualprod, SCIP_Real *maxnzsubgradientdualprod) |
| static void | updateLagrangianValue (SCIP *scip, SCIP_Real objval, SCIP_Real *dualvector, SCIP_ROW **cuts, int ncuts, SCIP_Real *lagrangianval) |
| static void | updateStepLength (SCIP *scip, SCIP_Real muparam, SCIP_Real ubparam, SCIP_Real lagrangianval, SCIP_Real *subgradient, int ncuts, SCIP_Real *steplength) |
| static void | updateBallRadius (SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_Real maxviolscore, SCIP_Real maxviolscoreold, SCIP_Real nviolscore, SCIP_Real nviolscoreold, int nlpiters, SCIP_Real *ballradius) |
| static SCIP_RETCODE | l1BallProjection (SCIP *scip, SCIP_Real *dualvector, int dualvectorlen, SCIP_Real radius) |
| static void | l2BallProjection (SCIP *scip, SCIP_Real *dualvector, int dualvectorlen, SCIP_Real radius) |
| static void | linfBallProjection (SCIP *scip, SCIP_Real *dualvector, int dualvectorlen, SCIP_Real radius) |
| static void | weightedDualVector (SCIP_SEPADATA *sepadata, SCIP_Real *dualvector, int dualvectorlen, SCIP_Real *stabilitycenter, int stabilitycenterlen, int nbestdualupdates, int totaliternum) |
| static SCIP_RETCODE | stabilizeDualVector (SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_Real *dualvector, int dualvectorlen, SCIP_Real *bestdualvector, int bestdualvectorlen, int nbestdualupdates, int subgradientiternum, int totaliternum, SCIP_Real maxviolscore, SCIP_Real maxviolscoreold, SCIP_Real nviolscore, SCIP_Real nviolscoreold, int nlpiters, SCIP_Real *ballradius) |
| static SCIP_RETCODE | updateDualVector (SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_Real *dualvector1, SCIP_Real *dualvector2, int dualvector2len, int ndualvector2updates, int subgradientiternum, int totaliternum, SCIP_Real steplength, SCIP_Real *subgradient, int ncuts, SCIP_Bool backtrack, SCIP_Real maxviolscore, SCIP_Real maxviolscoreold, SCIP_Real nviolscore, SCIP_Real nviolscoreold, int nlpiters, SCIP_Bool *dualvecsdiffer, SCIP_Real *ballradius) |
| static SCIP_RETCODE | checkLagrangianDualTermination (SCIP_SEPADATA *sepadata, int nnewaddedsoftcuts, int nyettoaddsoftcuts, SCIP_Bool objvecsdiffer, int ngeneratedcurrroundcuts, int nmaxgeneratedperroundcuts, int ncurrroundlpiters, int depth, SCIP_Bool *terminate) |
| static SCIP_RETCODE | solveLagromoryLP (SCIP *scip, SCIP_SEPADATA *sepadata, int depth, SCIP_Real origobjoffset, SCIP_Bool *solfound, SCIP_SOL *sol, SCIP_Real *solvals, SCIP_Real *objval, int *ncurrroundlpiters) |
| static SCIP_RETCODE | solveLPWithHardCuts (SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_Bool *solfound, SCIP_SOL *sol, SCIP_Real *solvals) |
| static SCIP_RETCODE | constructCutRow (SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, int mainiternum, int subgradientiternum, int cutnnz, int *cutinds, SCIP_Real *cutcoefs, SCIP_Real cutefficacy, SCIP_Real cutrhs, SCIP_Bool cutislocal, int cutrank, SCIP_ROW **generatedcuts, SCIP_Real *generatedcutefficacies, int ngeneratedcurrroundcuts, int *ngeneratednewcuts, SCIP_Bool *cutoff) |
| static SCIP_RETCODE | aggregateGeneratedCuts (SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, SCIP_ROW **generatedcuts, SCIP_Real *bestdualvector, int bestdualvectorlen, SCIP_ROW **aggrcuts, int *naggrcuts, SCIP_Bool *cutoff) |
| static SCIP_RETCODE | generateGMICuts (SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, int mainiternum, int subgradientiternum, SCIP_SOL *sol, SCIP_Real *solvals, int nmaxgeneratedperroundcuts, SCIP_Bool allowlocal, SCIP_ROW **generatedcurrroundcuts, SCIP_Real *generatedcutefficacies, int ngeneratedcurrroundcuts, int *ngeneratednewcuts, int depth, SCIP_Bool *cutoff) |
| static SCIP_RETCODE | updateObjectiveVector (SCIP *scip, SCIP_Real *dualvector, SCIP_ROW **cuts, int ncuts, SCIP_Real *origobjcoefs, SCIP_Bool *objvecsdiffer) |
| static SCIP_RETCODE | addGMICutsAsSoftConss (SCIP_Real *dualvector, int ngeneratedcuts, int *naddedcuts, int *nnewaddedsoftcuts) |
| static SCIP_RETCODE | solveLagrangianDual (SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, SCIP_SOL *sol, SCIP_Real *solvals, int mainiternum, SCIP_Real ubparam, int depth, SCIP_Bool allowlocal, int nmaxgeneratedperroundcuts, SCIP_Real *origobjcoefs, SCIP_Real origobjoffset, SCIP_Real *dualvector, int *nsoftcuts, SCIP_ROW **generatedcurrroundcuts, SCIP_Real *generatedcutefficacies, int *ngeneratedcutsperiter, int *ngeneratedcurrroundcuts, int *ncurrroundlpiters, SCIP_Bool *cutoff, SCIP_Real *bestlagrangianval, SCIP_Real *bestdualvector, int *bestdualvectorlen, int *nbestdualupdates, int *totaliternum) |
| static SCIP_RETCODE | generateInitCutPool (SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, int mainiternum, SCIP_SOL *sol, SCIP_Real *solvals, int nmaxgeneratedperroundcuts, SCIP_Bool allowlocal, SCIP_ROW **generatedcurrroundcuts, SCIP_Real *generatedcutefficacies, int *ngeneratedcutsperiter, int *ngeneratedcurrroundcuts, int depth, SCIP_Bool *cutoff) |
| static SCIP_RETCODE | addCuts (SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_ROW **cuts, int ncuts, SCIP_Longint maxdnom, SCIP_Real maxscale, int *naddedcuts, SCIP_Bool *cutoff) |
| static SCIP_RETCODE | checkMainLoopTermination (SCIP_SEPADATA *sepadata, SCIP_Bool cutoff, SCIP_Bool dualvecsdiffer, int ngeneratedcurrroundcuts, int nsoftcuts, int nmaxgeneratedperroundcuts, int ncurrroundlpiters, int depth, SCIP_Bool *terminate) |
| static SCIP_RETCODE | separateCuts (SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, SCIP_Real ubparam, int depth, SCIP_Bool allowlocal, SCIP_RESULT *result) |
| static SCIP_RETCODE | sepadataCreate (SCIP *scip, SCIP_SEPADATA **sepadata) |
| static | SCIP_DECL_SEPACOPY (sepaCopyLagromory) |
| static | SCIP_DECL_SEPAFREE (sepaFreeLagromory) |
| static | SCIP_DECL_SEPAINIT (sepaInitLagromory) |
| static | SCIP_DECL_SEPAEXIT (sepaExitLagromory) |
| static | SCIP_DECL_SEPAEXECLP (sepaExeclpLagromory) |
| SCIP_RETCODE | SCIPincludeSepaLagromory (SCIP *scip) |
| #define SEPA_NAME "lagromory" |
Definition at line 138 of file sepa_lagromory.c.
Definition at line 139 of file sepa_lagromory.c.
| #define SEPA_PRIORITY -8000 |
Definition at line 140 of file sepa_lagromory.c.
| #define SEPA_FREQ -1 |
Definition at line 141 of file sepa_lagromory.c.
| #define SEPA_MAXBOUNDDIST 1.0 |
Definition at line 142 of file sepa_lagromory.c.
| #define SEPA_USESSUBSCIP FALSE |
does the separator use a secondary SCIP instance?
Definition at line 143 of file sepa_lagromory.c.
| #define SEPA_DELAY FALSE |
should separation method be delayed, if other separators found cuts?
Definition at line 144 of file sepa_lagromory.c.
| #define DEFAULT_AWAY 0.01 |
minimal integrality violation of a basis variable to try separation
Definition at line 147 of file sepa_lagromory.c.
| #define DEFAULT_DELAYEDCUTS FALSE |
should cuts be added to the delayed cut pool?
Definition at line 148 of file sepa_lagromory.c.
| #define DEFAULT_SEPARATEROWS TRUE |
separate rows with integral slack?
Definition at line 149 of file sepa_lagromory.c.
| #define DEFAULT_SORTCUTOFFSOL TRUE |
sort fractional integer columns based on fractionality?
Definition at line 150 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_SIDETYPEBASIS TRUE |
choose side types of row (lhs/rhs) based on basis information?
Definition at line 151 of file sepa_lagromory.c.
| #define DEFAULT_DYNAMICCUTS TRUE |
should generated cuts be removed from the LP if they are no longer tight?
Definition at line 152 of file sepa_lagromory.c.
| #define DEFAULT_MAKEINTEGRAL FALSE |
try to scale all cuts to integral coefficients?
Definition at line 153 of file sepa_lagromory.c.
| #define DEFAULT_FORCECUTS FALSE |
force cuts to be added to the LP?
Definition at line 154 of file sepa_lagromory.c.
| #define DEFAULT_ALLOWLOCAL FALSE |
should locally valid cuts be generated?
Definition at line 155 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MAXROUNDSROOT 1 |
maximal number of separation rounds in the root node (-1: unlimited)
Definition at line 158 of file sepa_lagromory.c.
| #define DEFAULT_MAXROUNDS 1 |
maximal number of separation rounds per node (-1: unlimited)
Definition at line 159 of file sepa_lagromory.c.
| #define DEFAULT_DUALDEGENERACYRATETHRESHOLD 0.5 |
minimum dual degeneracy rate for separator execution
Definition at line 160 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_VARCONSRATIOTHRESHOLD 1.0 |
minimum variable-constraint ratio on optimal face for separator execution
Definition at line 161 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MINRESTART 1 |
minimum restart round for separator execution (0: from beginning of the instance solving, >= n with n >= 1: from restart round n)
Definition at line 162 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_PERLPMAXCUTSROOT 50 |
maximal number of cuts separated per Lagromory LP in the root node
Definition at line 164 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_PERLPMAXCUTS 10 |
maximal number of cuts separated per Lagromory LP in the non-root node
Definition at line 165 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_PERROUNDLPITERLIMITFACTOR -1.0 |
factor w.r.t. root node LP iterations for maximal separating LP iterations per separation round (negative for no limit)
Definition at line 166 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_ROOTLPITERLIMITFACTOR -1.0 |
factor w.r.t. root node LP iterations for maximal separating LP iterations in the root node (negative for no limit)
Definition at line 168 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_TOTALLPITERLIMITFACTOR -1.0 |
factor w.r.t. root node LP iterations for maximal separating LP iterations in the tree (negative for no limit)
Definition at line 170 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_PERROUNDMAXLPITERS 50000 |
maximal number of separating LP iterations per separation round (-1: unlimited)
Definition at line 172 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_PERROUNDCUTSFACTORROOT 1.0 |
factor w.r.t. number of integer columns for number of cuts separated per separation round in root node
Definition at line 173 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_PERROUNDCUTSFACTOR 0.5 |
factor w.r.t. number of integer columns for number of cuts separated per separation round at a non-root node
Definition at line 175 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_TOTALCUTSFACTOR 50.0 |
factor w.r.t. number of integer columns for total number of cuts separated
Definition at line 177 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MAXMAINITERS 4 |
maximal number of main loop iterations of the relax-and-cut algorithm
Definition at line 178 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MAXSUBGRADIENTITERS 6 |
maximal number of subgradient loop iterations of the relax-and-cut algorithm
Definition at line 179 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MUPARAMCONST TRUE |
is the mu parameter (factor for step length) constant?
Definition at line 182 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MUPARAMINIT 0.01 |
initial value of the mu parameter (factor for step length)
Definition at line 183 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MUPARAMLB 0.0 |
lower bound for the mu parameter (factor for step length)
Definition at line 184 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MUPARAMUB 2.0 |
upper bound for the mu parameter (factor for step length)
Definition at line 185 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MUBACKTRACKFACTOR 0.5 |
factor of mu while backtracking the mu parameter (factor for step length) - see updateMuSteplengthParam()
Definition at line 186 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MUSLAB1FACTOR 10.0 |
factor of mu parameter (factor for step length) for larger increment - see updateMuSteplengthParam()
Definition at line 188 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MUSLAB2FACTOR 2.0 |
factor of mu parameter (factor for step length) for smaller increment - see updateMuSteplengthParam()
Definition at line 190 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MUSLAB3FACTOR 0.5 |
factor of mu parameter (factor for step length) for reduction - see updateMuSteplengthParam()
Definition at line 192 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_DELTASLAB1UB 0.001 |
factor of delta deciding larger increment of mu parameter (factor for step length) - see updateMuSteplengthParam()
Definition at line 194 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_DELTASLAB2UB 0.01 |
factor of delta deciding smaller increment of mu parameter (factor for step length) - see updateMuSteplengthParam()
Definition at line 196 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_UBPARAMPOSFACTOR 2.0 |
factor for positive upper bound used as an estimate for the optimal Lagrangian dual value
Definition at line 198 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_UBPARAMNEGFACTOR 0.5 |
factor for negative upper bound used as an estimate for the optimal Lagrangian dual value
Definition at line 200 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MAXLAGRANGIANVALSFORAVG 2 |
maximal number of iterations for rolling average of Lagrangian value
Definition at line 202 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_MAXCONSECITERSFORMUUPDATE 10 |
consecutive number of iterations used to determine if mu needs to be backtracked
Definition at line 203 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_PERROOTLPITERFACTOR 0.2 |
factor w.r.t. root node LP iterations for iteration limit of each separating LP (negative for no limit)
Definition at line 204 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_PERLPITERFACTOR 0.1 |
factor w.r.t. non-root node LP iterations for iteration limit of each separating LP (negative for no limit)
Definition at line 206 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_CUTGENFREQ 1 |
frequency of subgradient iterations for generating cuts
Definition at line 208 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_CUTADDFREQ 1 |
frequency of subgradient iterations for adding cuts to objective function
Definition at line 209 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_CUTSFILTERFACTOR 1.0 |
fraction of generated cuts per explored basis to accept from separator
Definition at line 210 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_OPTIMALFACEPRIORITY 2 |
priority of the optimal face for separator execution (0: low priority, 1: medium priority, 2: high priority)
Definition at line 211 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_AGGREGATECUTS TRUE |
aggregate all generated cuts using the Lagrangian multipliers?
Definition at line 213 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_PROJECTIONTYPE 2 |
the ball into which the Lagrangian multipliers are projected for stabilization (0: no projection, 1: L1-norm ball projection, 2: L2-norm ball projection, 3: L_inf-norm ball projection)
Definition at line 216 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_STABILITYCENTERTYPE 1 |
type of stability center for taking weighted average of Lagrangian multipliers for stabilization (0: no weighted stabilization, 1: best Lagrangian multipliers)
Definition at line 219 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_RADIUSINIT 0.5 |
initial radius of the ball used in stabilization of Lagrangian multipliers
Definition at line 221 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_RADIUSMAX 20.0 |
maximum radius of the ball used in stabilization of Lagrangian multipliers
Definition at line 222 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_RADIUSMIN 1e-6 |
minimum radius of the ball used in stabilization of Lagrangian multipliers
Definition at line 223 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_CONST 2.0 |
a constant for stablity center based stabilization of Lagrangian multipliers
Definition at line 224 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define DEFAULT_RADIUSUPDATEWEIGHT 0.98 |
multiplier to evaluate cut violation score used for updating ball radius
Definition at line 225 of file sepa_lagromory.c.
Referenced by SCIPincludeSepaLagromory().
| #define RANDSEED 42 |
random seed
Definition at line 228 of file sepa_lagromory.c.
| #define MAKECONTINTEGRAL FALSE |
convert continuous variable to integral variables in SCIPmakeRowIntegral()?
Definition at line 229 of file sepa_lagromory.c.
| #define POSTPROCESS TRUE |
apply postprocessing after MIR calculation? - see SCIPcalcMIR()
Definition at line 230 of file sepa_lagromory.c.
| #define BOUNDSWITCH 0.9999 |
threshold for bound switching - see SCIPcalcMIR()
Definition at line 231 of file sepa_lagromory.c.
| #define VARTYPEUSEVBDS 2 |
We allow variable bound substitution for variables with continuous vartype only. See SCIPcalcMIR() for more information.
Definition at line 232 of file sepa_lagromory.c.
| #define FIXINTEGRALRHS FALSE |
try to generate an integral rhs? - see SCIPcalcMIR()
Definition at line 234 of file sepa_lagromory.c.
| #define MAXAGGRLEN | ( | ncols | ) |
maximal length of base inequality
Definition at line 235 of file sepa_lagromory.c.
|
static |
start the diving mode for solving LPs corresponding to the Lagrangian dual with fixed multipliers in the subgradient loop of the separator, and update some sepadata values
Definition at line 351 of file sepa_lagromory.c.
References assert(), i, MIN, NULL, SCIP_CALL, SCIP_Longint, SCIP_OKAY, SCIPcolIsIntegral(), SCIPgetLPColsData(), SCIPgetLPI(), SCIPgetLPRowsData(), SCIPgetNRootFirstLPIterations(), SCIPgetNRuns(), SCIPisInfinity(), SCIPstartDive(), and sepadata.
Referenced by separateCuts().
|
static |
end the diving mode that was used for solving LPs corresponding to the Lagrangian dual with fixed multipliers
Definition at line 432 of file sepa_lagromory.c.
References assert(), NULL, SCIP_CALL, SCIP_OKAY, SCIPendDive(), and sepadata.
Referenced by separateCuts().
|
static |
set up LP interface to solve LPs in the (outer) main loop of the relax-and-cut algorithm; these LPs are built by adding all the generated cuts to the node relaxation
| scip | SCIP data structure |
| sepadata | separator data structure |
| cuts | generated cuts to be added to the LP |
| ncuts | number of generated cuts to be added to the LP |
Definition at line 449 of file sepa_lagromory.c.
References assert(), i, NULL, SCIP_CALL, SCIP_OBJSEN_MINIMIZE, SCIP_OKAY, SCIP_Real, SCIPallocBufferArray, SCIPblkmem(), SCIPcolGetLb(), SCIPcolGetLPPos(), SCIPcolGetObj(), SCIPcolGetUb(), SCIPfreeBufferArray, SCIPgetLPColsData(), SCIPgetLPI(), SCIPgetLPRowsData(), SCIPgetMessagehdlr(), SCIPisInfinity(), SCIPlpiAddCols(), SCIPlpiAddRows(), SCIPlpiCreate(), SCIPlpiFree(), SCIPlpiFreeState(), SCIPlpiGetState(), SCIPlpiInfinity(), SCIPlpiSetState(), SCIProwGetCols(), SCIProwGetConstant(), SCIProwGetLhs(), SCIProwGetNLPNonz(), SCIProwGetNNonz(), SCIProwGetRhs(), SCIProwGetVals(), and sepadata.
Referenced by separateCuts().
|
static |
free separator data
Definition at line 648 of file sepa_lagromory.c.
References assert(), NULL, SCIP_CALL, SCIP_OKAY, SCIPfreeBlockMemory, SCIPlpiFree(), and sepadata.
Referenced by SCIP_DECL_SEPAFREE().
|
static |
update mu parameter which is used as a factor in the step length calculation; refer to the top of the file for a description of the formula.
| scip | SCIP data structure |
| sepadata | separator data structure |
| subgradientiternum | subgradient iteration number |
| ubparam | estimate of the optimal Lagrangian dual value |
| lagrangianvals | vector of Lagrangian values found so far |
| bestlagrangianval | best Lagrangian value found so far |
| avglagrangianval | rolling average of the Lagrangian values found so far |
| muparam | mu parameter to be updated |
| backtrack | whether mu parameter has been backtracked |
Definition at line 684 of file sepa_lagromory.c.
References assert(), FALSE, i, MAX, MIN, NULL, SCIP_Bool, SCIP_OKAY, SCIP_Real, SCIPisGE(), SCIPisPositive(), sepadata, and TRUE.
Referenced by solveLagrangianDual().
|
static |
update subgradient, i.e., residuals of generated cuts
| scip | SCIP data structure |
| sol | LP solution used in updating subgradient vector |
| cuts | cuts generated so far |
| ncuts | number of cuts generated so far |
| subgradient | vector of subgradients to be updated |
| dualvector | Lagrangian multipliers |
| subgradientzero | whether the subgradient vector is all zero |
| ncutviols | number of violations of generated cuts |
| maxcutviol | maximum violation of generated cuts |
| nnzsubgradientdualprod | number of nonzero products of subgradient vector and Lagrangian multipliers (i.e., number of complementarity slackness violations) |
| maxnzsubgradientdualprod | maximum value of nonzero products of subgradient vector and Lagrangian multipliers (i.e., maximum value of complementarity slackness violations) |
Definition at line 808 of file sepa_lagromory.c.
References assert(), FALSE, i, MAX, NULL, REALABS, SCIP_Bool, SCIP_Real, SCIPgetRowSolActivity(), SCIPisFeasPositive(), SCIPisFeasZero(), SCIPisInfinity(), SCIPisZero(), SCIProwGetConstant(), SCIProwGetLhs(), SCIProwGetRhs(), sol, and TRUE.
Referenced by solveLagrangianDual().
|
static |
update Lagrangian value, i.e., optimal value of the Lagrangian dual with fixed multipliers
| scip | SCIP data structure |
| objval | objective value of the Lagrangian dual with fixed multipliers |
| dualvector | Lagrangian multipliers |
| cuts | cuts generated so far |
| ncuts | number of cuts generated so far |
| lagrangianval | Lagrangian value to be updated |
Definition at line 879 of file sepa_lagromory.c.
References assert(), i, NULL, objval, SCIP_Real, SCIPisInfinity(), SCIProwGetConstant(), SCIProwGetLhs(), and SCIProwGetRhs().
Referenced by solveLagrangianDual().
|
static |
update step length based on various input arguments; refer to the top of the file for an expression
| scip | SCIP data structure |
| muparam | mu parameter used as a factor for step length |
| ubparam | estimate of the optimal Lagrangian dual value |
| lagrangianval | Lagrangian value |
| subgradient | subgradient vector |
| ncuts | number of cuts generated so far |
| steplength | step length to be updated |
Definition at line 904 of file sepa_lagromory.c.
References i, SCIP_Real, SCIPisFeasZero(), and SQR.
Referenced by solveLagrangianDual().
|
static |
update the ball radius (based on various violation metrics) that is used for stabilization of Lagrangian multipliers
| scip | SCIP data structure |
| sepadata | separator data structure |
| maxviolscore | weighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the current iteration |
| maxviolscoreold | weighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the previous iteration |
| nviolscore | weighted average of number of generated cut violations and number of complementarity slackness violations, in the current iteration |
| nviolscoreold | weighted average of number of generated cut violations and number of complementarity slackness violations, in the previous iteration |
| nlpiters | number of LP iterations taken for solving the Lagrangian dual with fixed multipliers in current iteration |
| ballradius | norm ball radius to be updated |
Definition at line 926 of file sepa_lagromory.c.
References assert(), MAX, MIN, NULL, SCIP_Bool, SCIP_Real, SCIPisNegative(), and sepadata.
Referenced by stabilizeDualVector().
|
static |
projection of Lagrangian multipliers onto L1-norm ball. This algorithm is based on the following article
Condat L. (2016).
Fast projection onto the simplex and the \(l_1\) ball.
Mathematical Programming, 158, 1-2, 575-585.
| scip | SCIP data structure |
| dualvector | Lagrangian multipliers to be projected onto L1-norm vall |
| dualvectorlen | length of the Lagrangian multipliers vector |
| radius | radius of the L1-norm ball |
Definition at line 995 of file sepa_lagromory.c.
References assert(), FALSE, i, MAX, REALABS, SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIP_Real, SCIPallocBufferArray, SCIPallocCleanBufferArray, SCIPfreeBufferArray, SCIPfreeCleanBufferArray, SCIPisGT(), SCIPisLE(), SCIPisNegative(), SCIPisPositive(), and TRUE.
Referenced by stabilizeDualVector().
|
static |
projection of Lagrangian multipliers onto L2-norm ball
| scip | SCIP data structure |
| dualvector | Lagrangian multipliers to be projected onto L2-norm vall |
| dualvectorlen | length of the Lagrangian multipliers vector |
| radius | radius of the L2-norm ball |
Definition at line 1135 of file sepa_lagromory.c.
References assert(), i, SCIP_Real, SCIPisGT(), SCIPisLT(), SCIPisNegative(), and SQR.
Referenced by stabilizeDualVector().
|
static |
projection of Lagrangian multipliers onto L_infinity-norm ball
| scip | SCIP data structure |
| dualvector | Lagrangian multipliers to be projected onto L_infinity-norm vall |
| dualvectorlen | length of the Lagrangian multipliers vector |
| radius | radius of the L_infinity-norm ball |
Definition at line 1166 of file sepa_lagromory.c.
References assert(), i, SCIP_Real, SCIPisGT(), SCIPisLT(), and SCIPisNegative().
Referenced by stabilizeDualVector().
|
static |
weighted Lagrangian multipliers based on a given vector as stability center
| sepadata | separator data structure |
| dualvector | Lagrangian multipliers |
| dualvectorlen | length of the Lagrangian multipliers vector |
| stabilitycenter | stability center (i.e., core vector of Lagrangian multipliers) |
| stabilitycenterlen | length of the stability center |
| nbestdualupdates | number of best Lagrangian values found so far |
| totaliternum | total number of iterations of the relax-and-cut algorithm performed so far |
Definition at line 1192 of file sepa_lagromory.c.
References alpha, assert(), i, MAX, MIN, SCIP_Real, and sepadata.
Referenced by stabilizeDualVector().
|
static |
stabilize Lagrangian multipliers
| scip | SCIP data structure |
| sepadata | separator data structure |
| dualvector | Lagrangian multipliers |
| dualvectorlen | length of the Lagrangian multipliers vector |
| bestdualvector | best Lagrangian multipliers found so far |
| bestdualvectorlen | length of the best Lagrangian multipliers vector |
| nbestdualupdates | number of best Lagrangian values found so far |
| subgradientiternum | iteration number of the subgradient algorithm |
| totaliternum | total number of iterations of the relax-and-cut algorithm performed so far |
| maxviolscore | weighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the current iteration |
| maxviolscoreold | weighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the previous iteration |
| nviolscore | weighted average of number of generated cut violations and number of complementarity slackness violations, in the current iteration |
| nviolscoreold | weighted average of number of generated cut violations and number of complementarity slackness violations, in the previous iteration |
| nlpiters | number of LP iterations taken for solving the Lagrangian dual with fixed multipliers in current iteration |
| ballradius | norm ball radius |
Definition at line 1225 of file sepa_lagromory.c.
References l1BallProjection(), l2BallProjection(), linfBallProjection(), SCIP_CALL, SCIP_OKAY, SCIP_Real, sepadata, updateBallRadius(), and weightedDualVector().
Referenced by updateDualVector().
|
static |
update Lagrangian multipliers
| scip | SCIP data structure |
| sepadata | separator data structure |
| dualvector1 | Lagrangian multipliers vector to be updated |
| dualvector2 | Lagrangian multipliers vector used for backtracking |
| dualvector2len | length of the Lagrangian multipliers vector used for backtracking |
| ndualvector2updates | number of best Lagrangian values found so far |
| subgradientiternum | iteration number of the subgradient algorithm |
| totaliternum | total number of iterations of the relax-and-cut algorithm performed so far |
| steplength | step length used for updating Lagrangian multipliers |
| subgradient | subgradient vector |
| ncuts | number of generated cuts so far |
| backtrack | whether the Lagrangian multipliers need to be backtracked |
| maxviolscore | weighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the current iteration |
| maxviolscoreold | weighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the previous iteration |
| nviolscore | weighted average of number of generated cut violations and number of complementarity slackness violations, in the current iteration |
| nviolscoreold | weighted average of number of generated cut violations and number of complementarity slackness violations, in the previous iteration |
| nlpiters | number of LP iterations taken for solving the Lagrangian dual with fixed multipliers in current iteration |
| dualvecsdiffer | whether the updated Lagrangian multipliers differ from the old one |
| ballradius | norm ball radius |
Definition at line 1285 of file sepa_lagromory.c.
References assert(), FALSE, i, MAX, NULL, SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIP_Real, SCIPallocCleanBufferArray, SCIPfreeCleanBufferArray, SCIPisEQ(), sepadata, stabilizeDualVector(), and TRUE.
Referenced by separateCuts(), and solveLagrangianDual().
|
static |
check different termination criteria
| sepadata | separator data structure |
| nnewaddedsoftcuts | number of cuts that were recently penalized and added to the Lagrangian dual's objective function |
| nyettoaddsoftcuts | number of cuts that are yet to be penalized and added to the Lagrangian dual's objective function |
| objvecsdiffer | whether the Lagrangian dual's objective function has changed |
| ngeneratedcurrroundcuts | number of cuts generated in the current separation round |
| nmaxgeneratedperroundcuts | maximal number of cuts allowed to generate per separation round |
| ncurrroundlpiters | number of separating LP iterations in the current separation round |
| depth | depth of the current node |
| terminate | whether to terminate the subgradient algorithm loop |
Definition at line 1394 of file sepa_lagromory.c.
References depth, FALSE, SCIP_Bool, SCIP_OKAY, sepadata, and TRUE.
Referenced by solveLagrangianDual().
|
static |
solve the LP corresponding to the Lagrangian dual with fixed Lagrangian multipliers
| scip | SCIP data structure |
| sepadata | separator data structure |
| depth | depth of the current node in the tree |
| origobjoffset | objective offset in the current node's relaxation |
| solfound | whether an LP optimal solution has been found |
| sol | data structure to store LP optimal solution, if found |
| solvals | values of the LP optimal solution, if found |
| objval | optimal objective value of the LP optimal solution, if found |
| ncurrroundlpiters | number of LP iterations taken for solving Lagrangian dual problems with fixed multipliers in the current separator round |
Definition at line 1441 of file sepa_lagromory.c.
References assert(), cutoff, depth, FALSE, i, lperror, MIN, NULL, objval, SCIP_Bool, SCIP_CALL, SCIP_Longint, SCIP_LPPAR_LPTILIM, SCIP_LPSOLSTAT_OPTIMAL, SCIP_OKAY, SCIP_Real, SCIPcolGetPrimsol(), SCIPcolGetVar(), SCIPdebugMsg, SCIPgetLPColsData(), SCIPgetLPObjval(), SCIPgetLPSolstat(), SCIPgetNLPIterations(), SCIPgetNNodeInitLPIterations(), SCIPgetNRootFirstLPIterations(), SCIPgetRealParam(), SCIPgetSolvingTime(), SCIPisInfinity(), SCIPisLPSolBasic(), SCIPlpiSetRealpar(), SCIPsetSolVal(), SCIPsolveDiveLP(), sepadata, sol, TRUE, and var.
Referenced by separateCuts(), and solveLagrangianDual().
|
static |
solve the LP corresponding to the node relaxation upon adding all the generated cuts
| scip | SCIP data structure |
| sepadata | separator data structure |
| solfound | whether an LP optimal solution has been found |
| sol | data structure to store LP optimal solution, if found |
| solvals | values of the LP optimal solution, if found |
Definition at line 1578 of file sepa_lagromory.c.
References assert(), FALSE, i, NULL, SCIP_Bool, SCIP_CALL, SCIP_LPPAR_LPTILIM, SCIP_OKAY, SCIP_Real, SCIPcolGetVar(), SCIPdebugMsg, SCIPgetLPColsData(), SCIPgetRealParam(), SCIPgetSolvingTime(), SCIPisInfinity(), SCIPlpiGetSol(), SCIPlpiIsPrimalFeasible(), SCIPlpiSetRealpar(), SCIPlpiSolvePrimal(), SCIPsetSolVal(), sepadata, sol, TRUE, and var.
Referenced by separateCuts().
|
static |
construct a cut based on the input cut coefficients, sides, etc
| scip | SCIP data structure |
| sepa | pointer to the separator |
| sepadata | separator data structure |
| mainiternum | iteration number of the outer loop of the relax-and-cut algorithm |
| subgradientiternum | iteration number of the subgradient algorithm |
| cutnnz | number of nonzeros in cut |
| cutinds | column indices in cut |
| cutcoefs | cut cofficients |
| cutefficacy | cut efficacy |
| cutrhs | RHS of cut |
| cutislocal | whether cut is local |
| cutrank | rank of cut |
| generatedcuts | array of generated cuts |
| generatedcutefficacies | array of generated cut efficacies w.r.t. the respective LP bases used for cut generations |
| ngeneratedcurrroundcuts | number of cuts generated until the previous basis in the current separation round |
| ngeneratednewcuts | number of new cuts generated using the current basis |
| cutoff | should the current node be cutoff? |
Definition at line 1642 of file sepa_lagromory.c.
References assert(), cutoff, FALSE, NULL, SCIP_Bool, SCIP_CALL, SCIP_LONGINT_FORMAT, SCIP_MAXSTRLEN, SCIP_OKAY, SCIP_Real, SCIPaddVarToRow(), SCIPcacheRowExtensions(), SCIPcolGetVar(), SCIPcreateEmptyRowSepa(), SCIPdebugMsg, SCIPflushRowExtensions(), SCIPgetLPCols(), SCIPgetRowMaxActivity(), SCIPgetRowMinActivity(), SCIPinfinity(), SCIPisEfficacious(), SCIPisFeasNegative(), SCIPisFeasPositive(), SCIPisInfinity(), SCIProwChgRank(), SCIProwGetLhs(), SCIProwGetName(), SCIProwGetNNonz(), SCIProwGetRhs(), SCIProwIsModifiable(), SCIPsepaGetName(), SCIPsnprintf(), sepadata, TRUE, and var.
Referenced by generateGMICuts().
|
static |
aggregated generated cuts based on the best Lagrangian multipliers
| scip | SCIP data structure |
| sepa | pointer to the separator |
| sepadata | separator data structure |
| generatedcuts | cuts generated in the current separation round |
| bestdualvector | best Lagrangian multipliers vector |
| bestdualvectorlen | length of the best Lagrangian multipliers vector |
| aggrcuts | aggregated cuts generated so far in the current separation round |
| naggrcuts | number of aggregated cuts generated so far in the current separation round |
| cutoff | should the current node be cutoff? |
Definition at line 1774 of file sepa_lagromory.c.
References assert(), cutoff, FALSE, i, MAX, NULL, QUAD, QUAD_ARRAY_LOAD, QUAD_ARRAY_SIZE, QUAD_ARRAY_STORE, QUAD_ASSIGN, QUAD_TO_DBL, SCIP_Bool, SCIP_CALL, SCIP_LONGINT_FORMAT, SCIP_MAXSTRLEN, SCIP_OKAY, SCIP_Real, SCIPaddVarToRow(), SCIPallocBufferArray, SCIPallocCleanBufferArray, SCIPcacheRowExtensions(), SCIPcolGetLPPos(), SCIPcolGetVar(), SCIPcreateEmptyRowSepa(), SCIPdebugMsg, SCIPflushRowExtensions(), SCIPfreeBufferArray, SCIPfreeCleanBufferArray, SCIPgetLPColsData(), SCIPgetRowMaxActivity(), SCIPgetRowMinActivity(), SCIPinfinity(), SCIPisFeasNegative(), SCIPisFeasPositive(), SCIPisGE(), SCIPisInfinity(), SCIPisZero(), SCIPquadprecProdDD, SCIPquadprecSumQQ, SCIProwChgRank(), SCIProwGetCols(), SCIProwGetConstant(), SCIProwGetLhs(), SCIProwGetName(), SCIProwGetNNonz(), SCIProwGetRank(), SCIProwGetRhs(), SCIProwGetVals(), SCIProwIsLocal(), SCIProwIsModifiable(), SCIPsepaGetName(), SCIPsnprintf(), sepadata, TRUE, and var.
Referenced by separateCuts().
|
static |
main method: LP solution separation method of separator
| scip | SCIP data structure |
| sepa | pointer to the separator |
| sepadata | separator data structure |
| mainiternum | iteration number of the outer loop of the relax-and-cut algorithm |
| subgradientiternum | iteration number of the subgradient algorithm |
| sol | LP solution to be used for cut generation |
| solvals | values of the LP solution to be used for cut generation |
| nmaxgeneratedperroundcuts | maximal number of cuts allowed to generate per separation round |
| allowlocal | should locally valid cuts be generated? |
| generatedcurrroundcuts | cuts generated in the current separation round |
| generatedcutefficacies | array of generated cut efficacies w.r.t. the respective LP bases used for cut generations |
| ngeneratedcurrroundcuts | number of cuts generated until the previous basis in the current separation round |
| ngeneratednewcuts | number of new cuts generated using the current basis |
| depth | depth of the current node in the tree |
| cutoff | should the current node be cutoff? |
Definition at line 1997 of file sepa_lagromory.c.
References assert(), BOUNDSWITCH, c, constructCutRow(), cutoff, depth, FIXINTEGRALRHS, frac, i, MAXAGGRLEN, MIN, NULL, POSTPROCESS, SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIP_Real, SCIPaggrRowCreate(), SCIPaggrRowFree(), SCIPaggrRowSumRows(), SCIPallocBufferArray, SCIPcalcMIR(), SCIPcolGetVar(), SCIPfeasFrac(), SCIPfreeBufferArray, SCIPgetLPBasisInd(), SCIPgetLPBInvRow(), SCIPgetLPColsData(), SCIPgetLPRowsData(), SCIPgetRowActivity(), SCIPisStopped(), SCIPrandomGetReal(), SCIProwIsIntegral(), SCIProwIsModifiable(), SCIPsortDownRealInt(), SCIPvarIsIntegral(), sepadata, sol, var, and VARTYPEUSEVBDS.
Referenced by generateInitCutPool(), and solveLagrangianDual().
|
static |
update objective vector w.r.t. the fixed Lagrangian multipliers
| scip | SCIP data structure |
| dualvector | Lagrangian multipliers vector |
| cuts | cuts generated so far in the current separation round |
| ncuts | number of cuts generated so far in the current separation round |
| origobjcoefs | original objective function coefficients of the node linear relaxation |
| objvecsdiffer | whether the updated objective function coefficients differ from the old ones |
Definition at line 2187 of file sepa_lagromory.c.
References assert(), FALSE, i, NULL, SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIP_Real, SCIPallocBufferArray, SCIPallocCleanBufferArray, SCIPchgVarObjDive(), SCIPcolGetLPPos(), SCIPcolGetVar(), SCIPfreeBufferArray, SCIPfreeCleanBufferArray, SCIPgetLPColsData(), SCIPgetVarObjDive(), SCIPisEQ(), SCIPisInfinity(), SCIPisZero(), SCIProwGetCols(), SCIProwGetLhs(), SCIProwGetNNonz(), SCIProwGetRhs(), SCIProwGetVals(), TRUE, and var.
Referenced by solveLagrangianDual().
|
static |
add GMI cuts to the objective function of the Lagrangian dual problem by introducing new Lagrangian multipliers
| dualvector | Lagrangian multipliers vector |
| ngeneratedcuts | number of cuts generated so far in the current separation round |
| naddedcuts | number of cuts added so far in the current separation round to the Lagrangian dual problem upon penalization |
| nnewaddedsoftcuts | number of cuts added newly to the Lagrangian dual problem upon penalization |
Definition at line 2272 of file sepa_lagromory.c.
References assert(), i, SCIP_OKAY, and SCIP_Real.
Referenced by solveLagrangianDual().
|
static |
solve the Lagrangian dual problem
| scip | SCIP data structure |
| sepa | pointer to the separator |
| sepadata | separator data structure |
| sol | data structure to store an LP solution upon solving a Lagrangian dual problem with fixed Lagrangian multipliers |
| solvals | values of the LP solution obtained upon solving a Lagrangian dual problem with fixed Lagrangian multipliers |
| mainiternum | iteration number of the outer loop of the relax-and-cut algorithm |
| ubparam | estimate of the optimal Lagrangian dual value |
| depth | depth of the current node in the tree |
| allowlocal | should locally valid cuts be generated? |
| nmaxgeneratedperroundcuts | maximal number of cuts allowed to generate per separation round |
| origobjcoefs | original objective function coefficients of the node linear relaxation |
| origobjoffset | original objective function offset of the node linear relaxation |
| dualvector | Lagrangian multipliers vector |
| nsoftcuts | number of generated cuts that were penalized and added to the Lagrangian dual problem |
| generatedcurrroundcuts | cuts generated in the current separation round |
| generatedcutefficacies | array of generated cut efficacies w.r.t. the respective LP bases used for cut generations |
| ngeneratedcutsperiter | number of cuts generated per subgradient iteration in the current separation round |
| ngeneratedcurrroundcuts | number of cuts generated so far in the current separation round |
| ncurrroundlpiters | number of LP iterations taken for solving Lagrangian dual problems with fixed multipliers in the current separator round |
| cutoff | should the current node be cutoff? |
| bestlagrangianval | best Lagrangian value found so far |
| bestdualvector | Lagrangian multipliers corresponding to the best Lagrangian value found so far |
| bestdualvectorlen | length of the Lagrangian multipliers corresponding to the best Lagrangian value found so far |
| nbestdualupdates | number of best Lagrangian values found so far |
| totaliternum | total number of iterations of the relax-and-cut algorithm performed so far |
Definition at line 2296 of file sepa_lagromory.c.
References addGMICutsAsSoftConss(), checkLagrangianDualTermination(), cutoff, depth, FALSE, generateGMICuts(), i, objval, SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIP_Real, SCIPallocBufferArray, SCIPallocCleanBufferArray, SCIPfreeBufferArray, SCIPfreeCleanBufferArray, SCIPisPositive(), SCIPisStopped(), sepadata, sol, solveLagromoryLP(), TRUE, updateDualVector(), updateLagrangianValue(), updateMuSteplengthParam(), updateObjectiveVector(), updateStepLength(), and updateSubgradient().
Referenced by separateCuts().
|
static |
generates initial cut pool before solving the Lagrangian dual
| scip | SCIP data structure |
| sepa | separator |
| sepadata | separator data structure |
| mainiternum | iteration number of the outer loop of the relax-and-cut algorithm |
| sol | LP solution to be used for cut generation |
| solvals | values of the LP solution to be used for cut generation |
| nmaxgeneratedperroundcuts | maximal number of cuts allowed to generate per separation round |
| allowlocal | should locally valid cuts be generated? |
| generatedcurrroundcuts | cuts generated in the current separation round |
| generatedcutefficacies | array of generated cut efficacies w.r.t. the respective LP bases used for cut generations |
| ngeneratedcutsperiter | number of cuts generated per subgradient iteration in the current separation round |
| ngeneratedcurrroundcuts | number of cuts generated so far in the current separation round |
| depth | depth of the current node in the tree |
| cutoff | should the current node be cutoff? |
Definition at line 2530 of file sepa_lagromory.c.
References cutoff, depth, generateGMICuts(), SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIP_Real, sepadata, and sol.
Referenced by separateCuts().
|
static |
add cuts to SCIP
| scip | SCIP data structure |
| sepadata | separator data structure |
| cuts | cuts generated so far in the current separation round |
| ncuts | number of cuts generated so far in the current separation round |
| maxdnom | maximum denominator in the rational representation of cuts |
| maxscale | maximal scale factor to scale the cuts to integral values |
| naddedcuts | number of cuts added to either global cutpool or sepastore |
| cutoff | should the current node be cutoff? |
Definition at line 2567 of file sepa_lagromory.c.
References assert(), cutoff, FALSE, i, MAKECONTINTEGRAL, NULL, SCIP_Bool, SCIP_CALL, SCIP_Longint, SCIP_OKAY, SCIP_Real, SCIPaddDelayedPoolCut(), SCIPaddPoolCut(), SCIPaddRow(), SCIPepsilon(), SCIPgetRowNumIntCols(), SCIPisCutNew(), SCIPisInfinity(), SCIPmakeRowIntegral(), SCIProwGetNNonz(), SCIProwGetRhs(), SCIProwIsLocal(), SCIPsumepsilon(), sepadata, and TRUE.
Referenced by separateCuts().
|
static |
check different termination criteria
| sepadata | separator data structure |
| cutoff | should the current node be cutoff? |
| dualvecsdiffer | whether the updated Lagrangian multipliers differ from the old one |
| ngeneratedcurrroundcuts | number of cuts generated in the current separation round |
| nsoftcuts | number of generated cuts that were penalized and added to the Lagrangian dual problem |
| nmaxgeneratedperroundcuts | maximal number of cuts allowed to generate per separation round |
| ncurrroundlpiters | number of LP iterations taken for solving Lagrangian dual problems with fixed multipliers in the current separator round |
| depth | depth of the current node in the tree |
| terminate | whether to terminate the relax-and-cut algorithm |
Definition at line 2644 of file sepa_lagromory.c.
References cutoff, depth, FALSE, SCIP_Bool, SCIP_OKAY, sepadata, and TRUE.
Referenced by separateCuts().
|
static |
searches and tries to add Lagromory cuts
| scip | SCIP data structure |
| sepa | separator |
| sepadata | separator data structure |
| ubparam | estimate of the optimal Lagrangian dual value |
| depth | depth of the current node in the tree |
| allowlocal | should locally valid cuts be generated? |
| result | final result of the separation round |
Definition at line 2692 of file sepa_lagromory.c.
References addCuts(), aggregateGeneratedCuts(), assert(), checkMainLoopTermination(), createLPWithHardCuts(), createLPWithSoftCuts(), cutoff, deleteLPWithSoftCuts(), depth, FALSE, generateInitCutPool(), i, maxdepth, NULL, objval, result, SCIP_Bool, SCIP_CALL, SCIP_CUTOFF, SCIP_DIDNOTFIND, SCIP_DIDNOTRUN, SCIP_Longint, SCIP_OKAY, SCIP_Real, SCIP_SEPARATED, SCIPallocBufferArray, SCIPallocCleanBufferArray, SCIPceil(), SCIPcolGetObj(), SCIPcreateSol(), SCIPdebugMsg, SCIPfreeBufferArray, SCIPfreeCleanBufferArray, SCIPfreeSol(), SCIPgetLPColsData(), SCIPgetMaxDepth(), SCIPgetNLPRows(), SCIPgetTransObjoffset(), SCIPheurPassSolTrySol(), SCIPinfinity(), SCIPisGE(), SCIPisPositive(), SCIPisStopped(), SCIPlpiGetSol(), SCIPlpiIsDualFeasible(), SCIPlpiWasSolved(), SCIPreleaseRow(), SCIPsepaGetNCalls(), SCIPsortDownRealInt(), sepadata, solveLagrangianDual(), solveLagromoryLP(), solveLPWithHardCuts(), TRUE, and updateDualVector().
Referenced by SCIP_DECL_SEPAEXECLP().
|
static |
creates separator data
Definition at line 3011 of file sepa_lagromory.c.
References assert(), BMSclearMemory, NULL, SCIP_CALL, SCIP_OKAY, SCIPallocBlockMemory, and sepadata.
Referenced by SCIPincludeSepaLagromory().
|
static |
copy method for separator plugins (called when SCIP copies plugins)
Definition at line 3032 of file sepa_lagromory.c.
References assert(), NULL, SCIP_CALL, SCIP_OKAY, SCIPincludeSepaLagromory(), SCIPsepaGetName(), and SEPA_NAME.
|
static |
destructor of separator to free user data (called when SCIP is exiting)
Definition at line 3047 of file sepa_lagromory.c.
References assert(), NULL, SCIP_CALL, SCIP_OKAY, SCIPsepaGetData(), SCIPsepaSetData(), sepadata, and sepadataFree().
|
static |
initialization method of separator (called after problem was transformed)
Definition at line 3063 of file sepa_lagromory.c.
References assert(), NULL, RANDSEED, SCIP_CALL, SCIP_OKAY, SCIPcreateRandom(), SCIPfindHeur(), SCIPsepaGetData(), sepadata, and TRUE.
|
static |
deinitialization method of separator (called before transformed problem is freed)
Definition at line 3083 of file sepa_lagromory.c.
References assert(), NULL, SCIP_OKAY, SCIPfreeRandom(), SCIPsepaGetData(), and sepadata.
|
static |
LP solution separation method of separator
Definition at line 3097 of file sepa_lagromory.c.
References assert(), depth, i, ncalls, NULL, result, SCIP_CALL, SCIP_DIDNOTRUN, SCIP_LPSOLSTAT_OPTIMAL, SCIP_OKAY, SCIP_Real, SCIPcolGetObj(), SCIPcolGetPrimsol(), SCIPcolGetVar(), SCIPgetBestSol(), SCIPgetLPColsData(), SCIPgetLPDualDegeneracy(), SCIPgetLPSolstat(), SCIPgetNLPRows(), SCIPgetNRuns(), SCIPgetSolVal(), SCIPgetTransObjoffset(), SCIPisLPSolBasic(), SCIPisPositive(), SCIPisStopped(), SCIPsepaGetData(), SCIPsepaGetName(), SCIPsepaGetNCallsAtNode(), SEPA_NAME, sepadata, separateCuts(), and var.