Scippy

SCIP

Solving Constraint Integer Programs

Detailed Description

Lagromory separator.

Author
Suresh Bolusani
Mark Ruben Turner
Mathieu Besançon

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.

  • Generate an initial pool of cuts to build the initial Lagrangian dual problem.
  • Select initial values for Lagrangian multipliers \(u^0\) (e.g., all zeroes vector).
  • In the outer main loop \(i\) of the algorithm:
    1. Solve the Lagrangian dual problem until certain termination criterion is met. This results in an inner subgradient loop, whose iteration \(j\) is described below.
      1. Fix \(u^j\), and solve the LP corresponding to the Lagrangian dual with fixed multipliers. Gather its optimal simplex tableau and optimal objective value (i.e., the Lagrangian value) \(L(u^j)\).
      2. Update \(u^j\) to \(u^{j+1}\) as follows.

        \[ 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 is ubparam, 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\).
      3. Stabilize \(u^{j+1}\) by projecting onto a norm ball followed by taking a convex combination with a core vector of Lagrangian multipliers.
      4. Generate GMI cuts based on the optimal simplex tableau.
      5. Relax the newly generated cuts by penalizing and adding them to the objective function.
      6. Go to the next iteration \(j+1\).
    2. Gather all the generated cuts and build an LP by adding all these cuts to the node relaxation.
    3. Solve this LP to obtain its optimal primal and dual solutions.
    4. If this primal solution is MIP primal feasible, then add this solution to the solution pool, add all the generated cuts to the cutpool or sepastore as needed, and exit the separator.
    5. Otherwise, update the Lagrangian multipliers based on this optimal dual solution, and go to the next iteration \(i+1\).

Definition in file sepa_lagromory.c.

Go to the source code of this file.

Macros

#define SEPA_NAME   "lagromory"
 
#define SEPA_DESC   "separator for Lagromory cuts for MIP relaxations"
 
#define SEPA_PRIORITY   -8000
 
#define SEPA_FREQ   -1
 
#define SEPA_MAXBOUNDDIST   1.0
 
#define SEPA_USESSUBSCIP   FALSE
 
#define SEPA_DELAY   FALSE
 
#define DEFAULT_AWAY   0.01
 
#define DEFAULT_DELAYEDCUTS   FALSE
 
#define DEFAULT_SEPARATEROWS   TRUE
 
#define DEFAULT_SORTCUTOFFSOL   TRUE
 
#define DEFAULT_SIDETYPEBASIS   TRUE
 
#define DEFAULT_DYNAMICCUTS   TRUE
 
#define DEFAULT_MAKEINTEGRAL   FALSE
 
#define DEFAULT_FORCECUTS   FALSE
 
#define DEFAULT_ALLOWLOCAL   FALSE
 
#define DEFAULT_MAXROUNDSROOT   1
 
#define DEFAULT_MAXROUNDS   1
 
#define DEFAULT_DUALDEGENERACYRATETHRESHOLD   0.5
 
#define DEFAULT_VARCONSRATIOTHRESHOLD   1.0
 
#define DEFAULT_MINRESTART   1
 
#define DEFAULT_PERLPMAXCUTSROOT   50
 
#define DEFAULT_PERLPMAXCUTS   10
 
#define DEFAULT_PERROUNDLPITERLIMITFACTOR   -1.0
 
#define DEFAULT_ROOTLPITERLIMITFACTOR   -1.0
 
#define DEFAULT_TOTALLPITERLIMITFACTOR   -1.0
 
#define DEFAULT_PERROUNDMAXLPITERS   50000
 
#define DEFAULT_PERROUNDCUTSFACTORROOT   1.0
 
#define DEFAULT_PERROUNDCUTSFACTOR   0.5
 
#define DEFAULT_TOTALCUTSFACTOR   50.0
 
#define DEFAULT_MAXMAINITERS   4
 
#define DEFAULT_MAXSUBGRADIENTITERS   6
 
#define DEFAULT_MUPARAMCONST   TRUE
 
#define DEFAULT_MUPARAMINIT   0.01
 
#define DEFAULT_MUPARAMLB   0.0
 
#define DEFAULT_MUPARAMUB   2.0
 
#define DEFAULT_MUBACKTRACKFACTOR   0.5
 
#define DEFAULT_MUSLAB1FACTOR   10.0
 
#define DEFAULT_MUSLAB2FACTOR   2.0
 
#define DEFAULT_MUSLAB3FACTOR   0.5
 
#define DEFAULT_DELTASLAB1UB   0.001
 
#define DEFAULT_DELTASLAB2UB   0.01
 
#define DEFAULT_UBPARAMPOSFACTOR   2.0
 
#define DEFAULT_UBPARAMNEGFACTOR   0.5
 
#define DEFAULT_MAXLAGRANGIANVALSFORAVG   2
 
#define DEFAULT_MAXCONSECITERSFORMUUPDATE   10
 
#define DEFAULT_PERROOTLPITERFACTOR   0.2
 
#define DEFAULT_PERLPITERFACTOR   0.1
 
#define DEFAULT_CUTGENFREQ   1
 
#define DEFAULT_CUTADDFREQ   1
 
#define DEFAULT_CUTSFILTERFACTOR   1.0
 
#define DEFAULT_OPTIMALFACEPRIORITY   2
 
#define DEFAULT_AGGREGATECUTS   TRUE
 
#define DEFAULT_PROJECTIONTYPE   2
 
#define DEFAULT_STABILITYCENTERTYPE   1
 
#define DEFAULT_RADIUSINIT   0.5
 
#define DEFAULT_RADIUSMAX   20.0
 
#define DEFAULT_RADIUSMIN   1e-6
 
#define DEFAULT_CONST   2.0
 
#define DEFAULT_RADIUSUPDATEWEIGHT   0.98
 
#define RANDSEED   42
 
#define MAKECONTINTEGRAL   FALSE
 
#define POSTPROCESS   TRUE
 
#define BOUNDSWITCH   0.9999
 
#define USEVBDS   TRUE
 
#define FIXINTEGRALRHS   FALSE
 
#define MAXAGGRLEN(ncols)   (0.1*(ncols)+1000)
 

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)
 

Macro Definition Documentation

◆ SEPA_NAME

#define SEPA_NAME   "lagromory"

Definition at line 138 of file sepa_lagromory.c.

Referenced by sepadataCreate().

◆ SEPA_DESC

#define SEPA_DESC   "separator for Lagromory cuts for MIP relaxations"

Definition at line 139 of file sepa_lagromory.c.

◆ SEPA_PRIORITY

#define SEPA_PRIORITY   -8000

Definition at line 140 of file sepa_lagromory.c.

◆ SEPA_FREQ

#define SEPA_FREQ   -1

Definition at line 141 of file sepa_lagromory.c.

◆ SEPA_MAXBOUNDDIST

#define SEPA_MAXBOUNDDIST   1.0

Definition at line 142 of file sepa_lagromory.c.

◆ SEPA_USESSUBSCIP

#define SEPA_USESSUBSCIP   FALSE

does the separator use a secondary SCIP instance?

Definition at line 143 of file sepa_lagromory.c.

◆ SEPA_DELAY

#define SEPA_DELAY   FALSE

should separation method be delayed, if other separators found cuts?

Definition at line 144 of file sepa_lagromory.c.

◆ DEFAULT_AWAY

#define DEFAULT_AWAY   0.01

minimal integrality violation of a basis variable to try separation

Definition at line 147 of file sepa_lagromory.c.

◆ DEFAULT_DELAYEDCUTS

#define DEFAULT_DELAYEDCUTS   FALSE

should cuts be added to the delayed cut pool?

Definition at line 148 of file sepa_lagromory.c.

◆ DEFAULT_SEPARATEROWS

#define DEFAULT_SEPARATEROWS   TRUE

separate rows with integral slack?

Definition at line 149 of file sepa_lagromory.c.

◆ DEFAULT_SORTCUTOFFSOL

#define DEFAULT_SORTCUTOFFSOL   TRUE

sort fractional integer columns based on fractionality?

Definition at line 150 of file sepa_lagromory.c.

◆ DEFAULT_SIDETYPEBASIS

#define DEFAULT_SIDETYPEBASIS   TRUE

choose side types of row (lhs/rhs) based on basis information?

Definition at line 151 of file sepa_lagromory.c.

◆ DEFAULT_DYNAMICCUTS

#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.

◆ DEFAULT_MAKEINTEGRAL

#define DEFAULT_MAKEINTEGRAL   FALSE

try to scale all cuts to integral coefficients?

Definition at line 153 of file sepa_lagromory.c.

◆ DEFAULT_FORCECUTS

#define DEFAULT_FORCECUTS   FALSE

force cuts to be added to the LP?

Definition at line 154 of file sepa_lagromory.c.

◆ DEFAULT_ALLOWLOCAL

#define DEFAULT_ALLOWLOCAL   FALSE

should locally valid cuts be generated?

Definition at line 155 of file sepa_lagromory.c.

◆ DEFAULT_MAXROUNDSROOT

#define DEFAULT_MAXROUNDSROOT   1

maximal number of separation rounds in the root node (-1: unlimited)

Definition at line 158 of file sepa_lagromory.c.

◆ DEFAULT_MAXROUNDS

#define DEFAULT_MAXROUNDS   1

maximal number of separation rounds per node (-1: unlimited)

Definition at line 159 of file sepa_lagromory.c.

◆ DEFAULT_DUALDEGENERACYRATETHRESHOLD

#define DEFAULT_DUALDEGENERACYRATETHRESHOLD   0.5

minimum dual degeneracy rate for separator execution

Definition at line 160 of file sepa_lagromory.c.

◆ DEFAULT_VARCONSRATIOTHRESHOLD

#define DEFAULT_VARCONSRATIOTHRESHOLD   1.0

minimum variable-constraint ratio on optimal face for separator execution

Definition at line 161 of file sepa_lagromory.c.

◆ DEFAULT_MINRESTART

#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.

◆ DEFAULT_PERLPMAXCUTSROOT

#define DEFAULT_PERLPMAXCUTSROOT   50

maximal number of cuts separated per Lagromory LP in the root node

Definition at line 165 of file sepa_lagromory.c.

◆ DEFAULT_PERLPMAXCUTS

#define DEFAULT_PERLPMAXCUTS   10

maximal number of cuts separated per Lagromory LP in the non-root node

Definition at line 166 of file sepa_lagromory.c.

◆ DEFAULT_PERROUNDLPITERLIMITFACTOR

#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 167 of file sepa_lagromory.c.

◆ DEFAULT_ROOTLPITERLIMITFACTOR

#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 170 of file sepa_lagromory.c.

◆ DEFAULT_TOTALLPITERLIMITFACTOR

#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 173 of file sepa_lagromory.c.

◆ DEFAULT_PERROUNDMAXLPITERS

#define DEFAULT_PERROUNDMAXLPITERS   50000

maximal number of separating LP iterations per separation round (-1: unlimited)

Definition at line 176 of file sepa_lagromory.c.

◆ DEFAULT_PERROUNDCUTSFACTORROOT

#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 177 of file sepa_lagromory.c.

◆ DEFAULT_PERROUNDCUTSFACTOR

#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 180 of file sepa_lagromory.c.

◆ DEFAULT_TOTALCUTSFACTOR

#define DEFAULT_TOTALCUTSFACTOR   50.0

factor w.r.t. number of integer columns for total number of cuts separated

Definition at line 183 of file sepa_lagromory.c.

◆ DEFAULT_MAXMAINITERS

#define DEFAULT_MAXMAINITERS   4

maximal number of main loop iterations of the relax-and-cut algorithm

Definition at line 184 of file sepa_lagromory.c.

◆ DEFAULT_MAXSUBGRADIENTITERS

#define DEFAULT_MAXSUBGRADIENTITERS   6

maximal number of subgradient loop iterations of the relax-and-cut algorithm

Definition at line 185 of file sepa_lagromory.c.

◆ DEFAULT_MUPARAMCONST

#define DEFAULT_MUPARAMCONST   TRUE

is the mu parameter (factor for step length) constant?

Definition at line 188 of file sepa_lagromory.c.

◆ DEFAULT_MUPARAMINIT

#define DEFAULT_MUPARAMINIT   0.01

initial value of the mu parameter (factor for step length)

Definition at line 189 of file sepa_lagromory.c.

◆ DEFAULT_MUPARAMLB

#define DEFAULT_MUPARAMLB   0.0

lower bound for the mu parameter (factor for step length)

Definition at line 190 of file sepa_lagromory.c.

◆ DEFAULT_MUPARAMUB

#define DEFAULT_MUPARAMUB   2.0

upper bound for the mu parameter (factor for step length)

Definition at line 191 of file sepa_lagromory.c.

◆ DEFAULT_MUBACKTRACKFACTOR

#define DEFAULT_MUBACKTRACKFACTOR   0.5

factor of mu while backtracking the mu parameter (factor for step length) - see updateMuSteplengthParam()

Definition at line 192 of file sepa_lagromory.c.

◆ DEFAULT_MUSLAB1FACTOR

#define DEFAULT_MUSLAB1FACTOR   10.0

factor of mu parameter (factor for step length) for larger increment - see updateMuSteplengthParam()

Definition at line 195 of file sepa_lagromory.c.

◆ DEFAULT_MUSLAB2FACTOR

#define DEFAULT_MUSLAB2FACTOR   2.0

factor of mu parameter (factor for step length) for smaller increment - see updateMuSteplengthParam()

Definition at line 198 of file sepa_lagromory.c.

◆ DEFAULT_MUSLAB3FACTOR

#define DEFAULT_MUSLAB3FACTOR   0.5

factor of mu parameter (factor for step length) for reduction - see updateMuSteplengthParam()

Definition at line 201 of file sepa_lagromory.c.

◆ DEFAULT_DELTASLAB1UB

#define DEFAULT_DELTASLAB1UB   0.001

factor of delta deciding larger increment of mu parameter (factor for step length) - see updateMuSteplengthParam()

Definition at line 204 of file sepa_lagromory.c.

◆ DEFAULT_DELTASLAB2UB

#define DEFAULT_DELTASLAB2UB   0.01

factor of delta deciding smaller increment of mu parameter (factor for step length) - see updateMuSteplengthParam()

Definition at line 207 of file sepa_lagromory.c.

◆ DEFAULT_UBPARAMPOSFACTOR

#define DEFAULT_UBPARAMPOSFACTOR   2.0

factor for positive upper bound used as an estimate for the optimal Lagrangian dual value

Definition at line 210 of file sepa_lagromory.c.

◆ DEFAULT_UBPARAMNEGFACTOR

#define DEFAULT_UBPARAMNEGFACTOR   0.5

factor for negative upper bound used as an estimate for the optimal Lagrangian dual value

Definition at line 213 of file sepa_lagromory.c.

◆ DEFAULT_MAXLAGRANGIANVALSFORAVG

#define DEFAULT_MAXLAGRANGIANVALSFORAVG   2

maximal number of iterations for rolling average of Lagrangian value

Definition at line 216 of file sepa_lagromory.c.

◆ DEFAULT_MAXCONSECITERSFORMUUPDATE

#define DEFAULT_MAXCONSECITERSFORMUUPDATE   10

consecutive number of iterations used to determine if mu needs to be backtracked

Definition at line 217 of file sepa_lagromory.c.

◆ DEFAULT_PERROOTLPITERFACTOR

#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 218 of file sepa_lagromory.c.

◆ DEFAULT_PERLPITERFACTOR

#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 221 of file sepa_lagromory.c.

◆ DEFAULT_CUTGENFREQ

#define DEFAULT_CUTGENFREQ   1

frequency of subgradient iterations for generating cuts

Definition at line 224 of file sepa_lagromory.c.

◆ DEFAULT_CUTADDFREQ

#define DEFAULT_CUTADDFREQ   1

frequency of subgradient iterations for adding cuts to objective function

Definition at line 225 of file sepa_lagromory.c.

◆ DEFAULT_CUTSFILTERFACTOR

#define DEFAULT_CUTSFILTERFACTOR   1.0

fraction of generated cuts per explored basis to accept from separator

Definition at line 226 of file sepa_lagromory.c.

◆ DEFAULT_OPTIMALFACEPRIORITY

#define DEFAULT_OPTIMALFACEPRIORITY   2

priority of the optimal face for separator execution (0: low priority, 1: medium priority, 2: high priority)

Definition at line 227 of file sepa_lagromory.c.

◆ DEFAULT_AGGREGATECUTS

#define DEFAULT_AGGREGATECUTS   TRUE

aggregate all generated cuts using the Lagrangian multipliers?

Definition at line 230 of file sepa_lagromory.c.

◆ DEFAULT_PROJECTIONTYPE

#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 233 of file sepa_lagromory.c.

◆ DEFAULT_STABILITYCENTERTYPE

#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 238 of file sepa_lagromory.c.

◆ DEFAULT_RADIUSINIT

#define DEFAULT_RADIUSINIT   0.5

initial radius of the ball used in stabilization of Lagrangian multipliers

Definition at line 241 of file sepa_lagromory.c.

◆ DEFAULT_RADIUSMAX

#define DEFAULT_RADIUSMAX   20.0

maximum radius of the ball used in stabilization of Lagrangian multipliers

Definition at line 242 of file sepa_lagromory.c.

◆ DEFAULT_RADIUSMIN

#define DEFAULT_RADIUSMIN   1e-6

minimum radius of the ball used in stabilization of Lagrangian multipliers

Definition at line 243 of file sepa_lagromory.c.

◆ DEFAULT_CONST

#define DEFAULT_CONST   2.0

a constant for stablity center based stabilization of Lagrangian multipliers

Definition at line 244 of file sepa_lagromory.c.

◆ DEFAULT_RADIUSUPDATEWEIGHT

#define DEFAULT_RADIUSUPDATEWEIGHT   0.98

multiplier to evaluate cut violation score used for updating ball radius

Definition at line 245 of file sepa_lagromory.c.

◆ RANDSEED

#define RANDSEED   42

random seed

Definition at line 248 of file sepa_lagromory.c.

◆ MAKECONTINTEGRAL

#define MAKECONTINTEGRAL   FALSE

convert continuous variable to integral variables in SCIPmakeRowIntegral()?

Definition at line 249 of file sepa_lagromory.c.

Referenced by addCuts().

◆ POSTPROCESS

#define POSTPROCESS   TRUE

apply postprocessing after MIR calculation? - see SCIPcalcMIR()

Definition at line 250 of file sepa_lagromory.c.

◆ BOUNDSWITCH

#define BOUNDSWITCH   0.9999

threshold for bound switching - see SCIPcalcMIR()

Definition at line 251 of file sepa_lagromory.c.

◆ USEVBDS

#define USEVBDS   TRUE

use variable bounds? - see SCIPcalcMIR()

Definition at line 252 of file sepa_lagromory.c.

◆ FIXINTEGRALRHS

#define FIXINTEGRALRHS   FALSE

try to generate an integral rhs? - see SCIPcalcMIR()

Definition at line 253 of file sepa_lagromory.c.

◆ MAXAGGRLEN

#define MAXAGGRLEN (   ncols)    (0.1*(ncols)+1000)

maximal length of base inequality

Definition at line 254 of file sepa_lagromory.c.

Function Documentation

◆ createLPWithSoftCuts()

static SCIP_RETCODE createLPWithSoftCuts ( SCIP scip,
SCIP_SEPADATA sepadata 
)
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

Parameters
scipSCIP data structure
sepadataseparator data structure

Definition at line 370 of file sepa_lagromory.c.

◆ deleteLPWithSoftCuts()

static SCIP_RETCODE deleteLPWithSoftCuts ( SCIP scip,
SCIP_SEPADATA sepadata 
)
static

end the diving mode that was used for solving LPs corresponding to the Lagrangian dual with fixed multipliers

Parameters
scipSCIP data structure
sepadataseparator data structure

Definition at line 451 of file sepa_lagromory.c.

References SCIP_Real.

◆ createLPWithHardCuts()

static SCIP_RETCODE createLPWithHardCuts ( SCIP scip,
SCIP_SEPADATA sepadata,
SCIP_ROW **  cuts,
int  ncuts 
)
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

Parameters
scipSCIP data structure
sepadataseparator data structure
cutsgenerated cuts to be added to the LP
ncutsnumber of generated cuts to be added to the LP

Definition at line 468 of file sepa_lagromory.c.

◆ sepadataFree()

static SCIP_RETCODE sepadataFree ( SCIP scip,
SCIP_SEPADATA **  sepadata 
)
static

free separator data

Parameters
scipSCIP data structure
sepadataseparator data structure

Definition at line 667 of file sepa_lagromory.c.

◆ updateMuSteplengthParam()

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

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.

Parameters
scipSCIP data structure
sepadataseparator data structure
subgradientiternumsubgradient iteration number
ubparamestimate of the optimal Lagrangian dual value
lagrangianvalsvector of Lagrangian values found so far
bestlagrangianvalbest Lagrangian value found so far
avglagrangianvalrolling average of the Lagrangian values found so far
muparammu parameter to be updated
backtrackwhether mu parameter has been backtracked

Definition at line 703 of file sepa_lagromory.c.

Referenced by solveLagrangianDual().

◆ updateSubgradient()

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

update subgradient, i.e., residuals of generated cuts

Note
assuming that \(i^{th}\) cut is of the form \({\alpha^i}^T x \leq {\alpha^i_0}\).
Parameters
scipSCIP data structure
solLP solution used in updating subgradient vector
cutscuts generated so far
ncutsnumber of cuts generated so far
subgradientvector of subgradients to be updated
dualvectorLagrangian multipliers
subgradientzerowhether the subgradient vector is all zero
ncutviolsnumber of violations of generated cuts
maxcutviolmaximum violation of generated cuts
nnzsubgradientdualprodnumber of nonzero products of subgradient vector and Lagrangian multipliers (i.e., number of complementarity slackness violations)
maxnzsubgradientdualprodmaximum value of nonzero products of subgradient vector and Lagrangian multipliers (i.e., maximum value of complementarity slackness violations)

Definition at line 827 of file sepa_lagromory.c.

Referenced by solveLagrangianDual().

◆ updateLagrangianValue()

static void updateLagrangianValue ( SCIP scip,
SCIP_Real  objval,
SCIP_Real dualvector,
SCIP_ROW **  cuts,
int  ncuts,
SCIP_Real lagrangianval 
)
static

update Lagrangian value, i.e., optimal value of the Lagrangian dual with fixed multipliers

Parameters
scipSCIP data structure
objvalobjective value of the Lagrangian dual with fixed multipliers
dualvectorLagrangian multipliers
cutscuts generated so far
ncutsnumber of cuts generated so far
lagrangianvalLagrangian value to be updated

Definition at line 898 of file sepa_lagromory.c.

References SCIP_Real, SCIPisFeasZero(), and SQR.

Referenced by solveLagrangianDual().

◆ updateStepLength()

static void updateStepLength ( SCIP scip,
SCIP_Real  muparam,
SCIP_Real  ubparam,
SCIP_Real  lagrangianval,
SCIP_Real subgradient,
int  ncuts,
SCIP_Real steplength 
)
static

update step length based on various input arguments; refer to the top of the file for an expression

Parameters
scipSCIP data structure
muparammu parameter used as a factor for step length
ubparamestimate of the optimal Lagrangian dual value
lagrangianvalLagrangian value
subgradientsubgradient vector
ncutsnumber of cuts generated so far
steplengthstep length to be updated

Definition at line 923 of file sepa_lagromory.c.

References NULL, and SCIP_Bool.

Referenced by solveLagrangianDual().

◆ updateBallRadius()

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

update the ball radius (based on various violation metrics) that is used for stabilization of Lagrangian multipliers

Parameters
scipSCIP data structure
sepadataseparator data structure
maxviolscoreweighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the current iteration
maxviolscoreoldweighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the previous iteration
nviolscoreweighted average of number of generated cut violations and number of complementarity slackness violations, in the current iteration
nviolscoreoldweighted average of number of generated cut violations and number of complementarity slackness violations, in the previous iteration
nlpitersnumber of LP iterations taken for solving the Lagrangian dual with fixed multipliers in current iteration
ballradiusnorm ball radius to be updated

Definition at line 945 of file sepa_lagromory.c.

Referenced by stabilizeDualVector().

◆ l1BallProjection()

static SCIP_RETCODE l1BallProjection ( SCIP scip,
SCIP_Real dualvector,
int  dualvectorlen,
SCIP_Real  radius 
)
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.

Parameters
scipSCIP data structure
dualvectorLagrangian multipliers to be projected onto L1-norm vall
dualvectorlenlength of the Lagrangian multipliers vector
radiusradius of the L1-norm ball

Definition at line 1014 of file sepa_lagromory.c.

Referenced by stabilizeDualVector().

◆ l2BallProjection()

static void l2BallProjection ( SCIP scip,
SCIP_Real dualvector,
int  dualvectorlen,
SCIP_Real  radius 
)
static

projection of Lagrangian multipliers onto L2-norm ball

Parameters
scipSCIP data structure
dualvectorLagrangian multipliers to be projected onto L2-norm vall
dualvectorlenlength of the Lagrangian multipliers vector
radiusradius of the L2-norm ball

Definition at line 1154 of file sepa_lagromory.c.

References linfBallProjection(), SCIP_Real, SCIPisGT(), SCIPisLT(), and SCIPisNegative().

Referenced by stabilizeDualVector().

◆ linfBallProjection()

static void linfBallProjection ( SCIP scip,
SCIP_Real dualvector,
int  dualvectorlen,
SCIP_Real  radius 
)
static

projection of Lagrangian multipliers onto L_infinity-norm ball

Parameters
scipSCIP data structure
dualvectorLagrangian multipliers to be projected onto L_infinity-norm vall
dualvectorlenlength of the Lagrangian multipliers vector
radiusradius of the L_infinity-norm ball

Definition at line 1185 of file sepa_lagromory.c.

References MAX, MIN, and SCIP_Real.

Referenced by l2BallProjection(), and stabilizeDualVector().

◆ weightedDualVector()

static void weightedDualVector ( SCIP_SEPADATA sepadata,
SCIP_Real dualvector,
int  dualvectorlen,
SCIP_Real stabilitycenter,
int  stabilitycenterlen,
int  nbestdualupdates,
int  totaliternum 
)
static

weighted Lagrangian multipliers based on a given vector as stability center

Parameters
sepadataseparator data structure
dualvectorLagrangian multipliers
dualvectorlenlength of the Lagrangian multipliers vector
stabilitycenterstability center (i.e., core vector of Lagrangian multipliers)
stabilitycenterlenlength of the stability center
nbestdualupdatesnumber of best Lagrangian values found so far
totaliternumtotal number of iterations of the relax-and-cut algorithm performed so far

Definition at line 1211 of file sepa_lagromory.c.

Referenced by stabilizeDualVector().

◆ stabilizeDualVector()

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

stabilize Lagrangian multipliers

Parameters
scipSCIP data structure
sepadataseparator data structure
dualvectorLagrangian multipliers
dualvectorlenlength of the Lagrangian multipliers vector
bestdualvectorbest Lagrangian multipliers found so far
bestdualvectorlenlength of the best Lagrangian multipliers vector
nbestdualupdatesnumber of best Lagrangian values found so far
subgradientiternumiteration number of the subgradient algorithm
totaliternumtotal number of iterations of the relax-and-cut algorithm performed so far
maxviolscoreweighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the current iteration
maxviolscoreoldweighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the previous iteration
nviolscoreweighted average of number of generated cut violations and number of complementarity slackness violations, in the current iteration
nviolscoreoldweighted average of number of generated cut violations and number of complementarity slackness violations, in the previous iteration
nlpitersnumber of LP iterations taken for solving the Lagrangian dual with fixed multipliers in current iteration
ballradiusnorm ball radius

Definition at line 1244 of file sepa_lagromory.c.

References l1BallProjection(), l2BallProjection(), linfBallProjection(), SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIP_Real, updateBallRadius(), updateDualVector(), and weightedDualVector().

Referenced by updateDualVector().

◆ updateDualVector()

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

update Lagrangian multipliers

Parameters
scipSCIP data structure
sepadataseparator data structure
dualvector1Lagrangian multipliers vector to be updated
dualvector2Lagrangian multipliers vector used for backtracking
dualvector2lenlength of the Lagrangian multipliers vector used for backtracking
ndualvector2updatesnumber of best Lagrangian values found so far
subgradientiternumiteration number of the subgradient algorithm
totaliternumtotal number of iterations of the relax-and-cut algorithm performed so far
steplengthstep length used for updating Lagrangian multipliers
subgradientsubgradient vector
ncutsnumber of generated cuts so far
backtrackwhether the Lagrangian multipliers need to be backtracked
maxviolscoreweighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the current iteration
maxviolscoreoldweighted average of maximum value of generated cut violations and maximum value of complementarity slackness violations, in the previous iteration
nviolscoreweighted average of number of generated cut violations and number of complementarity slackness violations, in the current iteration
nviolscoreoldweighted average of number of generated cut violations and number of complementarity slackness violations, in the previous iteration
nlpitersnumber of LP iterations taken for solving the Lagrangian dual with fixed multipliers in current iteration
dualvecsdifferwhether the updated Lagrangian multipliers differ from the old one
ballradiusnorm ball radius

Definition at line 1304 of file sepa_lagromory.c.

References checkLagrangianDualTermination(), FALSE, MAX, NULL, SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIP_Real, SCIPallocCleanBufferArray, SCIPfreeCleanBufferArray, SCIPisEQ(), stabilizeDualVector(), and TRUE.

Referenced by solveLagrangianDual(), and stabilizeDualVector().

◆ checkLagrangianDualTermination()

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

check different termination criteria

Note
the criterion based on objvecsdiffer assumes deterministic solving process (i.e., we would get same LP solution for "Lagrangian dual with fixed Lagrangian multipliers" when the objective vector remains the same across iterations).
Parameters
sepadataseparator data structure
nnewaddedsoftcutsnumber of cuts that were recently penalized and added to the Lagrangian dual's objective function
nyettoaddsoftcutsnumber of cuts that are yet to be penalized and added to the Lagrangian dual's objective function
objvecsdifferwhether the Lagrangian dual's objective function has changed
ngeneratedcurrroundcutsnumber of cuts generated in the current separation round
nmaxgeneratedperroundcutsmaximal number of cuts allowed to generate per separation round
ncurrroundlpitersnumber of separating LP iterations in the current separation round
depthdepth of the current node
terminatewhether to terminate the subgradient algorithm loop

Definition at line 1413 of file sepa_lagromory.c.

References SCIP_Bool, and SCIP_Real.

Referenced by solveLagrangianDual(), and updateDualVector().

◆ solveLagromoryLP()

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

solve the LP corresponding to the Lagrangian dual with fixed Lagrangian multipliers

Parameters
scipSCIP data structure
sepadataseparator data structure
depthdepth of the current node in the tree
origobjoffsetobjective offset in the current node's relaxation
solfoundwhether an LP optimal solution has been found
soldata structure to store LP optimal solution, if found
solvalsvalues of the LP optimal solution, if found
objvaloptimal objective value of the LP optimal solution, if found
ncurrroundlpitersnumber of LP iterations taken for solving Lagrangian dual problems with fixed multipliers in the current separator round

Definition at line 1460 of file sepa_lagromory.c.

Referenced by solveLagrangianDual().

◆ solveLPWithHardCuts()

static SCIP_RETCODE solveLPWithHardCuts ( SCIP scip,
SCIP_SEPADATA sepadata,
SCIP_Bool solfound,
SCIP_SOL sol,
SCIP_Real solvals 
)
static

solve the LP corresponding to the node relaxation upon adding all the generated cuts

Parameters
scipSCIP data structure
sepadataseparator data structure
solfoundwhether an LP optimal solution has been found
soldata structure to store LP optimal solution, if found
solvalsvalues of the LP optimal solution, if found

Definition at line 1597 of file sepa_lagromory.c.

◆ constructCutRow()

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

construct a cut based on the input cut coefficients, sides, etc

Parameters
scipSCIP data structure
sepapointer to the separator
sepadataseparator data structure
mainiternumiteration number of the outer loop of the relax-and-cut algorithm
subgradientiternumiteration number of the subgradient algorithm
cutnnznumber of nonzeros in cut
cutindscolumn indices in cut
cutcoefscut cofficients
cutefficacycut efficacy
cutrhsRHS of cut
cutislocalwhether cut is local
cutrankrank of cut
generatedcutsarray of generated cuts
generatedcutefficaciesarray of generated cut efficacies w.r.t. the respective LP bases used for cut generations
ngeneratedcurrroundcutsnumber of cuts generated until the previous basis in the current separation round
ngeneratednewcutsnumber of new cuts generated using the current basis
cutoffshould the current node be cutoff?

Definition at line 1661 of file sepa_lagromory.c.

References aggregateGeneratedCuts(), 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(), and TRUE.

◆ aggregateGeneratedCuts()

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

aggregated generated cuts based on the best Lagrangian multipliers

Parameters
scipSCIP data structure
sepapointer to the separator
sepadataseparator data structure
generatedcutscuts generated in the current separation round
bestdualvectorbest Lagrangian multipliers vector
bestdualvectorlenlength of the best Lagrangian multipliers vector
aggrcutsaggregated cuts generated so far in the current separation round
naggrcutsnumber of aggregated cuts generated so far in the current separation round
cutoffshould the current node be cutoff?

Definition at line 1793 of file sepa_lagromory.c.

Referenced by constructCutRow().

◆ generateGMICuts()

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

main method: LP solution separation method of separator

Parameters
scipSCIP data structure
sepapointer to the separator
sepadataseparator data structure
mainiternumiteration number of the outer loop of the relax-and-cut algorithm
subgradientiternumiteration number of the subgradient algorithm
solLP solution to be used for cut generation
solvalsvalues of the LP solution to be used for cut generation
nmaxgeneratedperroundcutsmaximal number of cuts allowed to generate per separation round
allowlocalshould locally valid cuts be generated?
generatedcurrroundcutscuts generated in the current separation round
generatedcutefficaciesarray of generated cut efficacies w.r.t. the respective LP bases used for cut generations
ngeneratedcurrroundcutsnumber of cuts generated until the previous basis in the current separation round
ngeneratednewcutsnumber of new cuts generated using the current basis
depthdepth of the current node in the tree
cutoffshould the current node be cutoff?

Definition at line 2016 of file sepa_lagromory.c.

Referenced by solveLagrangianDual().

◆ updateObjectiveVector()

static SCIP_RETCODE updateObjectiveVector ( SCIP scip,
SCIP_Real dualvector,
SCIP_ROW **  cuts,
int  ncuts,
SCIP_Real origobjcoefs,
SCIP_Bool objvecsdiffer 
)
static

update objective vector w.r.t. the fixed Lagrangian multipliers

Parameters
scipSCIP data structure
dualvectorLagrangian multipliers vector
cutscuts generated so far in the current separation round
ncutsnumber of cuts generated so far in the current separation round
origobjcoefsoriginal objective function coefficients of the node linear relaxation
objvecsdifferwhether the updated objective function coefficients differ from the old ones

Definition at line 2206 of file sepa_lagromory.c.

Referenced by solveLagrangianDual().

◆ addGMICutsAsSoftConss()

static SCIP_RETCODE addGMICutsAsSoftConss ( SCIP_Real dualvector,
int  ngeneratedcuts,
int *  naddedcuts,
int *  nnewaddedsoftcuts 
)
static

add GMI cuts to the objective function of the Lagrangian dual problem by introducing new Lagrangian multipliers

Parameters
dualvectorLagrangian multipliers vector
ngeneratedcutsnumber of cuts generated so far in the current separation round
naddedcutsnumber of cuts added so far in the current separation round to the Lagrangian dual problem upon penalization
nnewaddedsoftcutsnumber of cuts added newly to the Lagrangian dual problem upon penalization

Definition at line 2291 of file sepa_lagromory.c.

Referenced by solveLagrangianDual().

◆ solveLagrangianDual()

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

solve the Lagrangian dual problem

Parameters
scipSCIP data structure
sepapointer to the separator
sepadataseparator data structure
soldata structure to store an LP solution upon solving a Lagrangian dual problem with fixed Lagrangian multipliers
solvalsvalues of the LP solution obtained upon solving a Lagrangian dual problem with fixed Lagrangian multipliers
mainiternumiteration number of the outer loop of the relax-and-cut algorithm
ubparamestimate of the optimal Lagrangian dual value
depthdepth of the current node in the tree
allowlocalshould locally valid cuts be generated?
nmaxgeneratedperroundcutsmaximal number of cuts allowed to generate per separation round
origobjcoefsoriginal objective function coefficients of the node linear relaxation
origobjoffsetoriginal objective function offset of the node linear relaxation
dualvectorLagrangian multipliers vector
nsoftcutsnumber of generated cuts that were penalized and added to the Lagrangian dual problem
generatedcurrroundcutscuts generated in the current separation round
generatedcutefficaciesarray of generated cut efficacies w.r.t. the respective LP bases used for cut generations
ngeneratedcutsperiternumber of cuts generated per subgradient iteration in the current separation round
ngeneratedcurrroundcutsnumber of cuts generated so far in the current separation round
ncurrroundlpitersnumber of LP iterations taken for solving Lagrangian dual problems with fixed multipliers in the current separator round
cutoffshould the current node be cutoff?
bestlagrangianvalbest Lagrangian value found so far
bestdualvectorLagrangian multipliers corresponding to the best Lagrangian value found so far
bestdualvectorlenlength of the Lagrangian multipliers corresponding to the best Lagrangian value found so far
nbestdualupdatesnumber of best Lagrangian values found so far
totaliternumtotal number of iterations of the relax-and-cut algorithm performed so far

Definition at line 2315 of file sepa_lagromory.c.

References addGMICutsAsSoftConss(), checkLagrangianDualTermination(), FALSE, generateGMICuts(), generateInitCutPool(), SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIP_Real, SCIPallocBufferArray, SCIPallocCleanBufferArray, SCIPfreeBufferArray, SCIPfreeCleanBufferArray, SCIPisPositive(), SCIPisStopped(), solveLagromoryLP(), TRUE, updateDualVector(), updateLagrangianValue(), updateMuSteplengthParam(), updateObjectiveVector(), updateStepLength(), and updateSubgradient().

◆ generateInitCutPool()

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

generates initial cut pool before solving the Lagrangian dual

Parameters
scipSCIP data structure
sepaseparator
sepadataseparator data structure
mainiternumiteration number of the outer loop of the relax-and-cut algorithm
solLP solution to be used for cut generation
solvalsvalues of the LP solution to be used for cut generation
nmaxgeneratedperroundcutsmaximal number of cuts allowed to generate per separation round
allowlocalshould locally valid cuts be generated?
generatedcurrroundcutscuts generated in the current separation round
generatedcutefficaciesarray of generated cut efficacies w.r.t. the respective LP bases used for cut generations
ngeneratedcutsperiternumber of cuts generated per subgradient iteration in the current separation round
ngeneratedcurrroundcutsnumber of cuts generated so far in the current separation round
depthdepth of the current node in the tree
cutoffshould the current node be cutoff?

Definition at line 2553 of file sepa_lagromory.c.

Referenced by solveLagrangianDual().

◆ addCuts()

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

add cuts to SCIP

Parameters
scipSCIP data structure
sepadataseparator data structure
cutscuts generated so far in the current separation round
ncutsnumber of cuts generated so far in the current separation round
maxdnommaximum denominator in the rational representation of cuts
maxscalemaximal scale factor to scale the cuts to integral values
naddedcutsnumber of cuts added to either global cutpool or sepastore
cutoffshould the current node be cutoff?

Definition at line 2590 of file sepa_lagromory.c.

References checkMainLoopTermination(), FALSE, MAKECONTINTEGRAL, SCIP_Bool, SCIP_CALL, SCIP_OKAY, SCIPaddDelayedPoolCut(), SCIPaddPoolCut(), SCIPaddRow(), SCIPepsilon(), SCIPgetRowNumIntCols(), SCIPisCutNew(), SCIPisInfinity(), SCIPmakeRowIntegral(), SCIProwGetNNonz(), SCIProwGetRhs(), SCIProwIsLocal(), SCIPsumepsilon(), and TRUE.

◆ checkMainLoopTermination()

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

check different termination criteria

Parameters
sepadataseparator data structure
cutoffshould the current node be cutoff?
dualvecsdifferwhether the updated Lagrangian multipliers differ from the old one
ngeneratedcurrroundcutsnumber of cuts generated in the current separation round
nsoftcutsnumber of generated cuts that were penalized and added to the Lagrangian dual problem
nmaxgeneratedperroundcutsmaximal number of cuts allowed to generate per separation round
ncurrroundlpitersnumber of LP iterations taken for solving Lagrangian dual problems with fixed multipliers in the current separator round
depthdepth of the current node in the tree
terminatewhether to terminate the relax-and-cut algorithm

Definition at line 2667 of file sepa_lagromory.c.

References SCIP_Bool, and SCIP_Real.

Referenced by addCuts().

◆ separateCuts()

static SCIP_RETCODE separateCuts ( SCIP scip,
SCIP_SEPA sepa,
SCIP_SEPADATA sepadata,
SCIP_Real  ubparam,
int  depth,
SCIP_Bool  allowlocal,
SCIP_RESULT result 
)
static

searches and tries to add Lagromory cuts

Parameters
scipSCIP data structure
sepaseparator
sepadataseparator data structure
ubparamestimate of the optimal Lagrangian dual value
depthdepth of the current node in the tree
allowlocalshould locally valid cuts be generated?
resultfinal result of the separation round

Definition at line 2715 of file sepa_lagromory.c.

◆ sepadataCreate()

static SCIP_RETCODE sepadataCreate ( SCIP scip,
SCIP_SEPADATA **  sepadata 
)
static

creates separator data

Parameters
scipSCIP data structure
sepadataseparator data structure

Definition at line 3035 of file sepa_lagromory.c.

References NULL, SCIP_CALL, SCIP_DECL_SEPAFREE(), SCIP_OKAY, SCIPincludeSepaLagromory(), SCIPsepaGetData(), SCIPsepaGetName(), and SEPA_NAME.

◆ SCIP_DECL_SEPACOPY()

static SCIP_DECL_SEPACOPY ( sepaCopyLagromory  )
static

copy method for separator plugins (called when SCIP copies plugins)

Definition at line 3056 of file sepa_lagromory.c.

◆ SCIP_DECL_SEPAFREE()

static SCIP_DECL_SEPAFREE ( sepaFreeLagromory  )
static

destructor of separator to free user data (called when SCIP is exiting)

Definition at line 3071 of file sepa_lagromory.c.

Referenced by sepadataCreate().

◆ SCIP_DECL_SEPAINIT()

static SCIP_DECL_SEPAINIT ( sepaInitLagromory  )
static

initialization method of separator (called after problem was transformed)

Definition at line 3087 of file sepa_lagromory.c.

References NULL, SCIP_DECL_SEPAEXECLP(), SCIP_OKAY, SCIPfreeRandom(), and SCIPsepaGetData().

◆ SCIP_DECL_SEPAEXIT()

static SCIP_DECL_SEPAEXIT ( sepaExitLagromory  )
static

deinitialization method of separator (called before transformed problem is freed)

Definition at line 3107 of file sepa_lagromory.c.

◆ SCIP_DECL_SEPAEXECLP()

static SCIP_DECL_SEPAEXECLP ( sepaExeclpLagromory  )
static

LP solution separation method of separator

Definition at line 3121 of file sepa_lagromory.c.

Referenced by SCIP_DECL_SEPAINIT().