sepa_lagromory.c
Go to the documentation of this file.
32 * This separator is based on the following article that discusses Lagromory separation using the relax-and-cut
33 * framework. Multiple enhancements have been implemented on top of the basic algorithm described in the article.
48 * where \f$P\f$ is the feasible region of the relaxation. Let the following be the cuts generated so far in the current
55 * Then, the following is the Lagrangian dual problem considered in the relax-and-cut framework used in the separator.
58 * z_D := \max\limits_{u \geq 0} \left\{L(u) := \min \left\{c^T x + \sum\limits_{i = 1}^{M} \left(u_i
62 * where \f$u\f$ are the Lagrangian multipliers (referred to as \a dualvector in this separator) used for penalizing the
63 * violation of the generated cuts, and \f$z_D\f$ is the optimal objective value (which is approximated via \a ubparam in this separator).
64 * Then, the following are the steps of the relax-and-cut algorithm implemented in this separator.
71 * \item Solve the Lagrangian dual problem until certain termination criterion is met. This results in an inner
74 * \item Fix \f$u^j\f$, and solve the LP corresponding to the Lagrangian dual with fixed multipliers.
85 * where \f$mu_j\f$ is a factor (i.e., \a muparam) such that \f$0 < \mu_j \leq 2\f$, UB is \p ubparam,
96 * \p muslab2factor * mu_j & \text{if } \p deltaslab1ub * \delta \leq bestLB - avgLB < \p deltaslab2ub * delta\\
103 * \item Stabilize \f$u^{j+1}\f$ by projecting onto a norm ball followed by taking a convex combination
109 * \item Gather all the generated cuts and build an LP by adding all these cuts to the node relaxation.
111 * \item If this primal solution is MIP primal feasible, then add this solution to the solution pool, add all
113 * \item Otherwise, update the Lagrangian multipliers based on this optimal dual solution, and go to the next
128 * @todo for warm starting, if there are additional rows/cols, set their basis status to non-basic and then set WS info.
139 /*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
150 #define SEPA_DELAY FALSE /**< should separation method be delayed, if other separators found cuts? */
153 #define DEFAULT_AWAY 0.01 /**< minimal integrality violation of a basis variable to try separation */
156 #define DEFAULT_SORTCUTOFFSOL TRUE /**< sort fractional integer columns based on fractionality? */
157 #define DEFAULT_SIDETYPEBASIS TRUE /**< choose side types of row (lhs/rhs) based on basis information? */
158 #define DEFAULT_DYNAMICCUTS TRUE /**< should generated cuts be removed from the LP if they are no longer tight? */
164 #define DEFAULT_MAXROUNDSROOT 1 /**< maximal number of separation rounds in the root node (-1: unlimited) */
165 #define DEFAULT_MAXROUNDS 1 /**< maximal number of separation rounds per node (-1: unlimited) */
166 #define DEFAULT_DUALDEGENERACYRATETHRESHOLD 0.5 /**< minimum dual degeneracy rate for separator execution */
167 #define DEFAULT_VARCONSRATIOTHRESHOLD 1.0 /**< minimum variable-constraint ratio on optimal face for separator execution */
168 #define DEFAULT_MINRESTART 1 /**< minimum restart round for separator execution (0: from beginning of the
170 #define DEFAULT_PERLPMAXCUTSROOT 50 /**< maximal number of cuts separated per Lagromory LP in the root node */
171 #define DEFAULT_PERLPMAXCUTS 10 /**< maximal number of cuts separated per Lagromory LP in the non-root node */
172 #define DEFAULT_PERROUNDLPITERLIMITFACTOR -1.0 /**< factor w.r.t. root node LP iterations for maximal separating LP iterations
174 #define DEFAULT_ROOTLPITERLIMITFACTOR -1.0 /**< factor w.r.t. root node LP iterations for maximal separating LP iterations
176 #define DEFAULT_TOTALLPITERLIMITFACTOR -1.0 /**< factor w.r.t. root node LP iterations for maximal separating LP iterations
178 #define DEFAULT_PERROUNDMAXLPITERS 50000 /**< maximal number of separating LP iterations per separation round (-1: unlimited) */
179 #define DEFAULT_PERROUNDCUTSFACTORROOT 1.0 /**< factor w.r.t. number of integer columns for number of cuts separated per
181 #define DEFAULT_PERROUNDCUTSFACTOR 0.5 /**< factor w.r.t. number of integer columns for number of cuts separated per
183 #define DEFAULT_TOTALCUTSFACTOR 50.0 /**< factor w.r.t. number of integer columns for total number of cuts separated */
184 #define DEFAULT_MAXMAINITERS 4 /**< maximal number of main loop iterations of the relax-and-cut algorithm */
185 #define DEFAULT_MAXSUBGRADIENTITERS 6 /**< maximal number of subgradient loop iterations of the relax-and-cut algorithm */
188 #define DEFAULT_MUPARAMCONST TRUE /**< is the mu parameter (factor for step length) constant? */
189 #define DEFAULT_MUPARAMINIT 0.01 /**< initial value of the mu parameter (factor for step length) */
192 #define DEFAULT_MUBACKTRACKFACTOR 0.5 /**< factor of mu while backtracking the mu parameter (factor for step length) -
194 #define DEFAULT_MUSLAB1FACTOR 10.0 /**< factor of mu parameter (factor for step length) for larger increment - see
196 #define DEFAULT_MUSLAB2FACTOR 2.0 /**< factor of mu parameter (factor for step length) for smaller increment - see
198 #define DEFAULT_MUSLAB3FACTOR 0.5 /**< factor of mu parameter (factor for step length) for reduction - see
200 #define DEFAULT_DELTASLAB1UB 0.001 /**< factor of delta deciding larger increment of mu parameter (factor for step
202 #define DEFAULT_DELTASLAB2UB 0.01 /**< factor of delta deciding smaller increment of mu parameter (factor for step
204 #define DEFAULT_UBPARAMPOSFACTOR 2.0 /**< factor for positive upper bound used as an estimate for the optimal
206 #define DEFAULT_UBPARAMNEGFACTOR 0.5 /**< factor for negative upper bound used as an estimate for the optimal
208 #define DEFAULT_MAXLAGRANGIANVALSFORAVG 2 /**< maximal number of iterations for rolling average of Lagrangian value */
209 #define DEFAULT_MAXCONSECITERSFORMUUPDATE 10 /**< consecutive number of iterations used to determine if mu needs to be backtracked */
210 #define DEFAULT_PERROOTLPITERFACTOR 0.2 /**< factor w.r.t. root node LP iterations for iteration limit of each separating
212 #define DEFAULT_PERLPITERFACTOR 0.1 /**< factor w.r.t. non-root node LP iterations for iteration limit of each
215 #define DEFAULT_CUTADDFREQ 1 /**< frequency of subgradient iterations for adding cuts to objective function */
216 #define DEFAULT_CUTSFILTERFACTOR 1.0 /**< fraction of generated cuts per explored basis to accept from separator */
217 #define DEFAULT_OPTIMALFACEPRIORITY 2 /**< priority of the optimal face for separator execution (0: low priority, 1:
219 #define DEFAULT_AGGREGATECUTS TRUE /**< aggregate all generated cuts using the Lagrangian multipliers? */
222 #define DEFAULT_PROJECTIONTYPE 2 /**< the ball into which the Lagrangian multipliers are projected for
225 #define DEFAULT_STABILITYCENTERTYPE 1 /**< type of stability center for taking weighted average of Lagrangian multipliers for
227 #define DEFAULT_RADIUSINIT 0.5 /**< initial radius of the ball used in stabilization of Lagrangian multipliers */
228 #define DEFAULT_RADIUSMAX 20.0 /**< maximum radius of the ball used in stabilization of Lagrangian multipliers */
229 #define DEFAULT_RADIUSMIN 1e-6 /**< minimum radius of the ball used in stabilization of Lagrangian multipliers */
230 #define DEFAULT_CONST 2.0 /**< a constant for stablity center based stabilization of Lagrangian multipliers */
231 #define DEFAULT_RADIUSUPDATEWEIGHT 0.98 /**< multiplier to evaluate cut violation score used for updating ball radius */
235 #define MAKECONTINTEGRAL FALSE /**< convert continuous variable to integral variables in SCIPmakeRowIntegral()? */
236 #define POSTPROCESS TRUE /**< apply postprocessing after MIR calculation? - see SCIPcalcMIR() */
248 {
255 SCIP_Bool dynamiccuts; /**< should generated cuts be removed from the LP if they are no longer tight? */
263 SCIP_LPI* lpiwithsoftcuts; /**< pointer to the lpi interface of Lagrangian dual with fixed multipliers */
272 SCIP_Real dualdegeneracyratethreshold; /**< minimum dual degeneracy rate for separator execution */
273 SCIP_Real varconsratiothreshold; /**< minimum variable-constraint ratio on optimal face for separator execution */
274 int minrestart; /**< minimum restart round for separator execution (0: from beginning of the instance
276 int nmaxcutsperlproot; /**< maximal number of cuts separated per Lagromory LP in the root node */
277 int nmaxcutsperlp; /**< maximal number of cuts separated per Lagromory LP in the non-root node */
278 SCIP_Real perroundlpiterlimitfactor; /**< factor w.r.t. root node LP iterations for maximal separating LP iterations per
280 SCIP_Real rootlpiterlimitfactor; /**< factor w.r.t. root node LP iterations for maximal separating LP iterations in the
282 SCIP_Real totallpiterlimitfactor; /**< factor w.r.t. root node LP iterations for maximal separating LP iterations in the
284 int perroundnmaxlpiters; /**< maximal number of separating LP iterations per separation round (-1: unlimited) */
285 SCIP_Real perroundcutsfactorroot; /**< factor w.r.t. number of integer columns for number of cuts separated per
287 SCIP_Real perroundcutsfactor; /**< factor w.r.t. number of integer columns for number of cuts separated per
289 SCIP_Real totalcutsfactor; /**< factor w.r.t. number of integer columns for total number of cuts separated */
290 int nmaxmainiters; /**< maximal number of main loop iterations of the relax-and-cut algorithm */
291 int nmaxsubgradientiters; /**< maximal number of subgradient loop iterations of the relax-and-cut algorithm */
292 int nmaxperroundlpiters; /**< maximal number of separating LP iterations per separation round */
297 int nmaxperroundcutsroot; /**< maximal number of cuts separated per separation round in root node */
307 SCIP_Real mubacktrackfactor; /**< factor of mu while backtracking the mu parameter (factor for step length) - see
309 SCIP_Real muslab1factor; /**< factor of mu parameter (factor for step length) for larger increment - see
311 SCIP_Real muslab2factor; /**< factor of mu parameter (factor for step length) for smaller increment - see
313 SCIP_Real muslab3factor; /**< factor of mu parameter (factor for step length) for reduction - see updateMuSteplengthParam() */
314 SCIP_Real deltaslab1ub; /**< factor of delta deciding larger increment of mu parameter (factor for step
316 SCIP_Real deltaslab2ub; /**< factor of delta deciding smaller increment of mu parameter (factor for step
318 SCIP_Real ubparamposfactor; /**< factor for positive upper bound used as an estimate for the optimal Lagrangian
320 SCIP_Real ubparamnegfactor; /**< factor for negative upper bound used as an estimate for the optimal Lagrangian
322 int nmaxlagrangianvalsforavg; /**< maximal number of iterations for rolling average of Lagrangian value */
323 int nmaxconsecitersformuupdate; /**< consecutive number of iterations used to determine if mu needs to be backtracked */
324 SCIP_Real perrootlpiterfactor; /**< factor w.r.t. root node LP iterations for iteration limit of each separating LP
326 SCIP_Real perlpiterfactor; /**< factor w.r.t. non-root node LP iterations for iteration limit of each separating
329 int cutaddfreq; /**< frequency of subgradient iterations for adding cuts to objective function */
330 SCIP_Real cutsfilterfactor; /**< fraction of generated cuts per explored basis to accept from separator */
331 int optimalfacepriority; /**< priority of the optimal face for separator execution (0: low priority, 1: medium
336 int projectiontype; /**< the ball into which the Lagrangian multipliers are projected for stabilization
339 int stabilitycentertype; /**< type of stability center for taking weighted average of Lagrangian multipliers for
341 SCIP_Real radiusinit; /**< initial radius of the ball used in stabilization of Lagrangian multipliers */
342 SCIP_Real radiusmax; /**< maximum radius of the ball used in stabilization of Lagrangian multipliers */
343 SCIP_Real radiusmin; /**< minimum radius of the ball used in stabilization of Lagrangian multipliers */
344 SCIP_Real constant; /**< a constant for stablity center based stabilization of Lagrangian multipliers */
345 SCIP_Real radiusupdateweight; /**< multiplier to evaluate cut violation score used for updating ball radius */
353 /** start the diving mode for solving LPs corresponding to the Lagrangian dual with fixed multipliers in the subgradient
376 runnum = SCIPgetNRuns(scip); /* current run number of SCIP (starts from 1) indicating restarts */
377 nruns = sepadata->nrunsinsoftcutslp; /* previous run number of SCIP in which the diving LP was created */
386 /* if called for the first time in a restart (including the initial run), set certain sepadata values */
392 /* calculate maximum number of LP iterations allowed for all separation calls in the root node */
393 if( (sepadata->rootlpiterlimitfactor >= 0.0) && !SCIPisInfinity(scip, sepadata->rootlpiterlimitfactor) )
398 /* calculate maximum number of LP iterations allowed for all separation calls in the entire tree */
399 if( (sepadata->totallpiterlimitfactor >= 0.0) && !SCIPisInfinity(scip, sepadata->totallpiterlimitfactor) )
405 if( (sepadata->perroundlpiterlimitfactor >= 0.0) && !SCIPisInfinity(scip, sepadata->perroundlpiterlimitfactor) )
412 sepadata->nmaxperroundlpiters = ((sepadata->nmaxperroundlpiters >= 0) ? MIN(sepadata->nmaxperroundlpiters,
415 /* set maximum number of cuts allowed to generate per round in root and non-root nodes as well as the total tree */
428 /* update the run number of solving to represent the restart number in which the above limits were set */
435 /** end the diving mode that was used for solving LPs corresponding to the Lagrangian dual with fixed multipliers */
450 /** set up LP interface to solve LPs in the (outer) main loop of the relax-and-cut algorithm; these LPs are built by
495 /* if this function is called for the first time in this separation round, then create an LPI and add cols & rows */
506 SCIP_CALL( SCIPlpiCreate(lpi, SCIPgetMessagehdlr(scip), "node LP with generated cuts", SCIP_OBJSEN_MINIMIZE) );
531 assert(!(SCIPisInfinity(scip, -SCIProwGetLhs(rows[i])) && SCIPisInfinity(scip, SCIProwGetRhs(rows[i]))));
553 rowlhs[i] = SCIPisInfinity(scip, -SCIProwGetLhs(rows[i])) ? ninf : SCIProwGetLhs(rows[i]) - rowconst;
554 rowrhs[i] = SCIPisInfinity(scip, SCIProwGetRhs(rows[i])) ? pinf : SCIProwGetRhs(rows[i]) - rowconst;
568 SCIP_CALL( SCIPlpiAddRows(*lpi, nrows, rowlhs, rowrhs, NULL, rowbegs[nrows], rowbegs, rowcolinds, rowvals) );
590 assert(!(SCIPisInfinity(scip, -SCIProwGetLhs(cuts[i])) && SCIPisInfinity(scip, SCIProwGetRhs(cuts[i]))));
597 SCIP_CALL( SCIPallocBufferArray(scip, &rowbegs, (ncuts - sepadata->nrowsinhardcutslp + nrows) + 1) );
598 SCIP_CALL( SCIPallocBufferArray(scip, &rowlhs, (ncuts - sepadata->nrowsinhardcutslp + nrows)) );
599 SCIP_CALL( SCIPallocBufferArray(scip, &rowrhs, (ncuts - sepadata->nrowsinhardcutslp + nrows)) );
610 rowbegs[i - sepadata->nrowsinhardcutslp + nrows + 1] = rowbegs[i - sepadata->nrowsinhardcutslp + nrows] +
613 rowlhs[i - sepadata->nrowsinhardcutslp + nrows] = SCIPisInfinity(scip, -SCIProwGetLhs(cuts[i])) ? ninf :
615 rowrhs[i - sepadata->nrowsinhardcutslp + nrows] = SCIPisInfinity(scip, SCIProwGetRhs(cuts[i])) ? pinf :
630 SCIP_CALL( SCIPlpiAddRows(*lpi, (ncuts - sepadata->nrowsinhardcutslp + nrows), rowlhs, rowrhs, NULL,
673 (*sepadata)->nmaxtotallpiters = 0;
685 /** update mu parameter which is used as a factor in the step length calculation; refer to the top of the file for a
774 /* if certain number of iterations are done, then check for a possibility of backtracking and apply accordingly */
823 int* nnzsubgradientdualprod, /**< number of nonzero products of subgradient vector and Lagrangian multipliers (i.e.,
825 SCIP_Real* maxnzsubgradientdualprod /**< maximum value of nonzero products of subgradient vector and Lagrangian multipliers
851 subgradient[i] = SCIPgetRowSolActivity(scip, cuts[i], sol) + SCIProwGetConstant(cuts[i]) - SCIProwGetRhs(cuts[i]);
882 /** update Lagrangian value, i.e., optimal value of the Lagrangian dual with fixed multipliers */
904 }
907 /** update step length based on various input arguments; refer to the top of the file for an expression */
929 /** update the ball radius (based on various violation metrics) that is used for stabilization of Lagrangian multipliers */
934 SCIP_Real maxviolscore, /**< weighted average of maximum value of generated cut violations and maximum value of
936 SCIP_Real maxviolscoreold, /**< weighted average of maximum value of generated cut violations and maximum value of
938 SCIP_Real nviolscore, /**< weighted average of number of generated cut violations and number of complementarity
940 SCIP_Real nviolscoreold, /**< weighted average of number of generated cut violations and number of complementarity
942 int nlpiters, /**< number of LP iterations taken for solving the Lagrangian dual with fixed multipliers in
957 /* both the maximum violation and number of violations scores have become better, so, increase the radius */
971 /* both the maximum violation and number of violations scores have become worse, so, decrease the radius */
977 /* only one among the maximum violation and number of violations scores has become better, and the LP basis did
992 /** projection of Lagrangian multipliers onto L1-norm ball. This algorithm is based on the following article
1028 /* if the vector of Lagrangian multipliers lies outside the L1-norm ball, then do the projection */
1093 /* @note the third condition (temp1len - ntemp1removed > 0) is true as long as the first condition
1161 /* if the vector of Lagrangian multipliers is outside the L2-norm ball, then do the projection */
1173 SCIP_Real* dualvector, /**< Lagrangian multipliers to be projected onto L_infinity-norm vall */
1182 /* if the vector of Lagrangian multipliers is outside the L_infinity-norm ball, then do the projection */
1193 /* @todo calculate weight outside this function and pass it (so that this function becomes generic and independent of
1201 SCIP_Real* stabilitycenter, /**< stability center (i.e., core vector of Lagrangian multipliers) */
1204 int totaliternum /**< total number of iterations of the relax-and-cut algorithm performed so far */
1239 int totaliternum, /**< total number of iterations of the relax-and-cut algorithm performed so far */
1240 SCIP_Real maxviolscore, /**< weighted average of maximum value of generated cut violations and maximum value of
1242 SCIP_Real maxviolscoreold, /**< weighted average of maximum value of generated cut violations and maximum value of
1244 SCIP_Real nviolscore, /**< weighted average of number of generated cut violations and number of complementarity
1246 SCIP_Real nviolscoreold, /**< weighted average of number of generated cut violations and number of complementarity
1248 int nlpiters, /**< number of LP iterations taken for solving the Lagrangian dual with fixed multipliers in
1258 updateBallRadius(scip, sepadata, maxviolscore, maxviolscoreold, nviolscore, nviolscoreold, nlpiters,
1281 /* weighted Lagrangian multipliers based on best Langrangian multipliers as stability center */
1282 weightedDualVector(sepadata, dualvector, dualvectorlen, bestdualvector, bestdualvectorlen, nbestdualupdates, totaliternum);
1298 int totaliternum, /**< total number of iterations of the relax-and-cut algorithm performed so far */
1303 SCIP_Real maxviolscore, /**< weighted average of maximum value of generated cut violations and maximum value of
1305 SCIP_Real maxviolscoreold, /**< weighted average of maximum value of generated cut violations and maximum value of
1307 SCIP_Real nviolscore, /**< weighted average of number of generated cut violations and number of complementarity
1309 SCIP_Real nviolscoreold, /**< weighted average of number of generated cut violations and number of complementarity
1311 int nlpiters, /**< number of LP iterations taken for solving the Lagrangian dual with fixed multipliers in
1313 SCIP_Bool* dualvecsdiffer, /**< whether the updated Lagrangian multipliers differ from the old one */
1325 /* @todo do allocation and free operations outside only once instead of every time this function is called? */
1326 /* copy of the Lagrangian multipliers to be used to check if the updated vector is different than this */
1331 /* if backtracking was not identified at the time of the mu parameter update, then update the Lagrangian multipliers
1354 SCIP_CALL( stabilizeDualVector(scip, sepadata, dualvector1, ncuts, dualvector2, dualvector2len,
1355 ndualvector2updates, subgradientiternum, totaliternum, maxviolscore, maxviolscoreold, nviolscore,
1358 /* projection onto non-negative orthant again in case stabilization changed some components negative*/
1362 /* if backtracking was identified at the time of the mu parameter update, then backtrack the Lagrangian multipliers
1394 * @note the criterion based on objvecsdiffer assumes deterministic solving process (i.e., we would get same LP solution
1395 * for "Lagrangian dual with fixed Lagrangian multipliers" when the objective vector remains the same across iterations).
1401 int nnewaddedsoftcuts, /**< number of cuts that were recently penalized and added to the Lagrangian dual's
1403 int nyettoaddsoftcuts, /**< number of cuts that are yet to be penalized and added to the Lagrangian dual's
1407 int nmaxgeneratedperroundcuts, /**< maximal number of cuts allowed to generate per separation round */
1408 int ncurrroundlpiters, /**< number of separating LP iterations in the current separation round */
1415 /* check if no new cuts were added to the Lagrangian dual, no cuts are remaining to be added, and the objective
1416 * function of the Lagrangian dual with fixed multipliers had not changed from the previous iteration
1429 /* check if allowed number of simplex iterations in this separation round have already been used up */
1430 if( (sepadata->nmaxperroundlpiters >= 0) && (ncurrroundlpiters >= sepadata->nmaxperroundlpiters) )
1434 if( (depth == 0) && (sepadata->nmaxrootlpiters >= 0) && (sepadata->nrootlpiters >= sepadata->nmaxrootlpiters) )
1438 if( (sepadata->nmaxtotallpiters >= 0) && (sepadata->ntotallpiters >= sepadata->nmaxtotallpiters) )
1455 int* ncurrroundlpiters /**< number of LP iterations taken for solving Lagrangian dual problems with fixed multipliers
1497 /* @note the following direct LPI call is being used because of the lack of an equivalent function call in
1505 (sepadata->perrootlpiterfactor >= 0.0 && !SCIPisInfinity(scip, sepadata->perrootlpiterfactor)) )
1526 /* @todo impose a finite iteration limit only when the dualvector changes from zero to non-zero for the first time because
1527 * many simplex pivots are performed in this case even with warm starting (compared to the case when the
1541 /* @todo is there any way to accept terminations due to iterlimit and timelimit as well? It is not possible
1544 /* @note ideally, only primal feasibility is sufficient. But, there is no such option with SCIPgetLPSolstat. */
1661 SCIP_Real* generatedcutefficacies, /**< array of generated cut efficacies w.r.t. the respective LP bases used for cut
1663 int ngeneratedcurrroundcuts, /**< number of cuts generated until the previous basis in the current separation round */
1667 {
1683 SCIPdebugMsg(scip, " -> Lagromory cut detected node infeasibility with cut 0 <= %g.\n", cutrhs);
1699 (void) SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "%s_%" SCIP_LONGINT_FORMAT "_%d" "_%d" "_%d", SCIPsepaGetName(sepa),
1700 sepadata->ncalls, mainiternum, subgradientiternum, ngeneratedcurrroundcuts + *ngeneratednewcuts);
1704 (void) SCIPsnprintf(cutname, SCIP_MAXSTRLEN, "%s_%" SCIP_LONGINT_FORMAT "_%d" "_%d", SCIPsepaGetName(sepa),
1709 SCIP_CALL( SCIPcreateEmptyRowSepa(scip, &cut, sepa, cutname, -SCIPinfinity(scip), cutrhs, cutislocal, FALSE,
1733 SCIPdebugMsg(scip, " -> Lagromory cut detected node infeasibility with cut 0 <= %g.\n", cutrhs);
1739 /* gathering lhs and rhs again in case the separator is extended later to add other cuts/constraints that may
1746 SCIPdebugMsg(scip, " -> %s cut <%s>: rhs=%f, eff=%f\n", "lagromory", cutname, cutrhs, cutefficacy);
1787 int* naggrcuts, /**< number of aggregated cuts generated so far in the current separation round */
1911 SCIPdebugMsg(scip, " -> Lagromory cut detected node infeasibility with cut 0 <= %g.\n", aggrcutrhs);
1917 (void) SCIPsnprintf(aggrcutname, SCIP_MAXSTRLEN, "%s_%" SCIP_LONGINT_FORMAT "_aggregated", SCIPsepaGetName(sepa),
1921 SCIP_CALL( SCIPcreateEmptyRowSepa(scip, &aggrcut, sepa, aggrcutname, aggrcutlhs - aggrcutconst, aggrcutrhs -
1943 SCIPdebugMsg(scip, " -> Lagromory cut detected node infeasibility with cut 0 <= %g.\n", aggrcutrhs);
1949 /* gathering lhs and rhs again in case the separator is extended later to add other cuts/constraints that may
2010 int nmaxgeneratedperroundcuts, /**< maximal number of cuts allowed to generate per separation round */
2013 SCIP_Real* generatedcutefficacies, /**< array of generated cut efficacies w.r.t. the respective LP bases used for cut
2015 int ngeneratedcurrroundcuts, /**< number of cuts generated until the previous basis in the current separation round */
2071 /* check if the rows of the simplex tableau are suitable for cut generation and build an array of fractions */
2148 /* get the row of B^-1 for this basic integer variable with fractional solution value and call aggregate function */
2151 SCIP_CALL( SCIPaggrRowSumRows(scip, aggrrow, binvrow, inds, ninds, sepadata->sidetypebasis, allowlocal, 2,
2157 /* @todo Currently we are using the SCIPcalcMIR() function to compute the coefficients of the Gomory
2158 * cut. Alternatively, we could use the direct version (see thesis of Achterberg formula (8.4)) which
2163 SCIP_CALL( SCIPcalcMIR(scip, sol, POSTPROCESS, BOUNDSWITCH, USEVBDS, allowlocal, FIXINTEGRALRHS, NULL,
2164 NULL, minfrac, maxfrac, 1.0, aggrrow, cutcoefs, &cutrhs, cutinds, &cutnnz, &cutefficacy, &cutrank,
2171 SCIP_CALL( constructCutRow(scip, sepa, sepadata, mainiternum, subgradientiternum, cutnnz, cutinds, cutcoefs,
2197 SCIP_Real* origobjcoefs, /**< original objective function coefficients of the node linear relaxation */
2198 SCIP_Bool* objvecsdiffer /**< whether the updated objective function coefficients differ from the old ones */
2230 assert(!(SCIPisInfinity(scip, -SCIProwGetLhs(cuts[i])) && SCIPisInfinity(scip, SCIProwGetRhs(cuts[i]))));
2275 /** add GMI cuts to the objective function of the Lagrangian dual problem by introducing new Lagrangian multipliers */
2280 int* naddedcuts, /**< number of cuts added so far in the current separation round to the Lagrangian dual problem
2282 int* nnewaddedsoftcuts /**< number of cuts added newly to the Lagrangian dual problem upon penalization */
2297 }
2305 SCIP_SOL* sol, /**< data structure to store an LP solution upon solving a Lagrangian dual problem with
2307 SCIP_Real* solvals, /**< values of the LP solution obtained upon solving a Lagrangian dual problem with fixed
2313 int nmaxgeneratedperroundcuts,/**< maximal number of cuts allowed to generate per separation round */
2314 SCIP_Real* origobjcoefs, /**< original objective function coefficients of the node linear relaxation */
2315 SCIP_Real origobjoffset, /**< original objective function offset of the node linear relaxation */
2317 int* nsoftcuts, /**< number of generated cuts that were penalized and added to the Lagrangian dual problem */
2319 SCIP_Real* generatedcutefficacies, /**< array of generated cut efficacies w.r.t. the respective LP bases used for cut
2321 int* ngeneratedcutsperiter, /**< number of cuts generated per subgradient iteration in the current separation round */
2322 int* ngeneratedcurrroundcuts, /**< number of cuts generated so far in the current separation round */
2323 int* ncurrroundlpiters, /**< number of LP iterations taken for solving Lagrangian dual problems with fixed
2327 SCIP_Real* bestdualvector, /**< Lagrangian multipliers corresponding to the best Lagrangian value found so far */
2328 int* bestdualvectorlen, /**< length of the Lagrangian multipliers corresponding to the best Lagrangian value
2331 int* totaliternum /**< total number of iterations of the relax-and-cut algorithm performed so far */
2390 SCIP_CALL( updateObjectiveVector(scip, dualvector, generatedcurrroundcuts, *nsoftcuts, origobjcoefs, &objvecsdiffer) );
2411 SCIP_CALL( solveLagromoryLP(scip, sepadata, depth, origobjoffset, &solfound, sol, solvals, &objval,
2428 ngeneratedcutsperiter[mainiternum * sepadata->nmaxsubgradientiters + i + 1] = ngeneratednewcuts;
2431 /* update subgradient, i.e., find the residuals of the penalized cuts, and determine various violations */
2432 updateSubgradient(scip, sol, generatedcurrroundcuts, *nsoftcuts, subgradient, dualvector, &subgradientzero,
2435 /* calculate Lagrangian value for the fixed Lagrangian multipliers, and update best and avg values */
2436 updateLagrangianValue(scip, objval, dualvector, generatedcurrroundcuts, *nsoftcuts, &lagrangianval);
2457 /* if the subgradient vector is non-zero, then update the mu parameter and the Lagrangian multipliers */
2461 SCIP_CALL( updateMuSteplengthParam(scip, sepadata, i, ubparam, lagrangianvals, *bestlagrangianval, avglagrangianval,
2475 *nbestdualupdates, i, *totaliternum, steplength, subgradient, *nsoftcuts, backtrack, maxviolscore,
2481 SCIP_CALL( updateObjectiveVector(scip, dualvector, generatedcurrroundcuts, *nsoftcuts, origobjcoefs, &objvecsdiffer) );
2484 /* if the subgradient vector if zero, then simply mark that the Lagrangian multipliers and the objective
2492 /* add generated GMI cuts to the objective function of the Lagrangian dual problem by introducing new
2497 SCIP_CALL( addGMICutsAsSoftConss(dualvector, *ngeneratedcurrroundcuts, nsoftcuts, &nnewaddedsoftcuts) );
2502 /* add any remaining generated GMI cuts to the objective function of the Lagrangian dual problem by introducing
2506 SCIP_CALL( addGMICutsAsSoftConss(dualvector, *ngeneratedcurrroundcuts, nsoftcuts, &nnewaddedsoftcuts) );
2513 SCIP_CALL( checkLagrangianDualTermination(sepadata, nnewaddedsoftcuts, *ngeneratedcurrroundcuts - *nsoftcuts,
2520 /* add any remaining generated GMI cuts to the objective function of the Lagrangian dual problem by introducing new
2524 SCIP_CALL( addGMICutsAsSoftConss(dualvector, *ngeneratedcurrroundcuts, nsoftcuts, &nnewaddedsoftcuts) );
2546 int nmaxgeneratedperroundcuts,/**< maximal number of cuts allowed to generate per separation round */
2549 SCIP_Real* generatedcutefficacies, /**< array of generated cut efficacies w.r.t. the respective LP bases used for cut
2551 int* ngeneratedcutsperiter, /**< number of cuts generated per subgradient iteration in the current separation round */
2552 int* ngeneratedcurrroundcuts, /**< number of cuts generated so far in the current separation round */
2559 ngeneratednewcuts = 0;
2562 SCIP_CALL( generateGMICuts(scip, sepa, sepadata, mainiternum, -1, sol, solvals, nmaxgeneratedperroundcuts,
2563 allowlocal, generatedcurrroundcuts, generatedcutefficacies, *ngeneratedcurrroundcuts, &ngeneratednewcuts,
2605 /* Add the bound change as cut to avoid that the LP gets modified. This would mean that the LP is not flushed
2606 * and the method SCIPgetLPBInvRow() fails; SCIP internally will apply this bound change automatically. */
2656 SCIP_Bool dualvecsdiffer, /**< whether the updated Lagrangian multipliers differ from the old one */
2658 int nsoftcuts, /**< number of generated cuts that were penalized and added to the Lagrangian dual problem */
2659 int nmaxgeneratedperroundcuts, /**< maximal number of cuts allowed to generate per separation round */
2660 int ncurrroundlpiters, /**< number of LP iterations taken for solving Lagrangian dual problems with fixed
2672 /* check if the Lagrangian multipliers do not differ from the previous iteration and no new cuts exist for penalizing */
2684 /* check if allowed number of simplex iterations in this separation round have already been used up */
2685 if( (sepadata->nmaxperroundlpiters >= 0) && (ncurrroundlpiters >= sepadata->nmaxperroundlpiters) )
2689 if( (depth == 0) && (sepadata->nmaxrootlpiters >= 0) && (sepadata->nrootlpiters >= sepadata->nmaxrootlpiters) )
2693 if( (sepadata->nmaxtotallpiters >= 0) && (sepadata->ntotallpiters >= sepadata->nmaxtotallpiters) )
2776 /* initialize the LP that will be used to solve the Lagrangian dual with fixed Lagrangian multipliers */
2778 /* create the LP that represents the node relaxation including all the generated cuts in this separator */
2785 nmaxgeneratedperroundcuts = ((depth == 0) ? sepadata->nmaxperroundcutsroot : sepadata->nmaxperroundcuts);
2810 SCIP_CALL( solveLagromoryLP(scip, sepadata, depth, origobjoffset, &solfound, softcutslpsol, softcutslpsolvals, &objval,
2814 SCIP_CALL( generateInitCutPool(scip, sepa, sepadata, 0, softcutslpsol, softcutslpsolvals, nmaxgeneratedperroundcuts, allowlocal,
2815 generatedcurrroundcuts, generatedcutefficacies, ngeneratedcutsperiter, &ngeneratedcurrroundcuts, depth,
2819 SCIP_CALL( checkMainLoopTermination(sepadata, cutoff, TRUE, ngeneratedcurrroundcuts, nsoftcuts,
2830 SCIP_CALL( solveLagrangianDual(scip, sepa, sepadata, softcutslpsol, softcutslpsolvals, i, ubparam, depth, allowlocal,
2831 nmaxgeneratedperroundcuts, origobjcoefs, origobjoffset, dualvector, &nsoftcuts, generatedcurrroundcuts,
2832 generatedcutefficacies, ngeneratedcutsperiter, &ngeneratedcurrroundcuts, &ncurrroundlpiters, &cutoff,
2835 /* @todo filter cuts before adding them to the new LP that was created based on the node relaxation? */
2840 SCIP_CALL( createLPWithHardCuts(scip, sepadata, generatedcurrroundcuts, ngeneratedcurrroundcuts) );
2846 /* @note if trysol heuristic is was not present, then the solution found above, which can potentially be a MIP
2854 /* if dual feasible, then fetch dual solution and reset Lagrangian multipliers based on it. otherwise, retain the
2860 ngeneratedcurrroundcuts, 0, -1, -1, 0.0, NULL, ngeneratedcurrroundcuts, TRUE, 0.0, 0.0, 0.0, 0.0, -1,
2865 SCIP_CALL( updateDualVector(scip, sepadata, dualvector, dualvector, nsoftcuts, 0, -1, -1, 0.0, NULL,
2871 SCIP_CALL( checkMainLoopTermination(sepadata, cutoff, dualvecsdiffer, ngeneratedcurrroundcuts, nsoftcuts,
2875 /* set the maximal denominator in rational representation of gomory cut and the maximal scale factor to
2878 /* @todo find better but still stable gomory cut settings: look at dcmulti, gesa3, khb0525, misc06, p2756 */
2879 /* @note above todo was copied from sepa_gomory.c. So, if gomory code is changed, same changes can be done here. */
2910 SCIP_CALL( addCuts(scip, sepadata, generatedcurrroundcuts, nselectedcurrroundcuts, maxdnom, maxscale,
2925 SCIPsortDownRealInt(&generatedcutefficacies[nprocessedcuts], cutindsperm, ngeneratedcutsperiter[i]);
2929 SCIP_CALL( addCuts(scip, sepadata, &generatedcurrroundcuts[nprocessedcuts], nselectedcuts, maxdnom,
2941 SCIP_CALL( addCuts(scip, sepadata, generatedcurrroundcuts, nselectedcurrroundcuts, maxdnom, maxscale,
2949 SCIP_CALL( aggregateGeneratedCuts(scip, sepa, sepadata, generatedcurrroundcuts, bestdualvector, bestdualvectorlen,
2954 SCIP_CALL( addCuts(scip, sepadata, aggregatedcurrroundcuts, naggregatedcurrroundcuts, maxdnom, maxscale,
3183 /* if a MIP primal solution exists, use it to estimate the optimal value of the Lagrangian dual problem */
3198 /* if a MIP primal solution does not exist, then use the node relaxation's LP solutin to estimate the optimal value
3212 ubparam *= SCIPisPositive(scip, ubparam) ? sepadata->ubparamposfactor : sepadata->ubparamnegfactor;
3216 SCIP_CALL( separateCuts(scip, sepa, sepadata, ubparam, depth, (allowlocal && sepadata->allowlocal), result) );
3240 SCIP_CALL( SCIPincludeSepaBasic(scip, &sepa, SEPA_NAME, SEPA_DESC, SEPA_PRIORITY, SEPA_FREQ, SEPA_MAXBOUNDDIST,
3252 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/away", "minimal integrality violation of a basis "
3255 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/rootlpiterlimitfactor", "factor w.r.t. root node LP "
3257 &sepadata->rootlpiterlimitfactor, TRUE, DEFAULT_ROOTLPITERLIMITFACTOR, -1.0, SCIP_REAL_MAX, NULL, NULL) );
3259 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/totallpiterlimitfactor", "factor w.r.t. root node LP "
3261 &sepadata->totallpiterlimitfactor, TRUE, DEFAULT_TOTALLPITERLIMITFACTOR, -1.0, SCIP_REAL_MAX, NULL, NULL) );
3263 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/perroundlpiterlimitfactor", "factor w.r.t. root node LP "
3264 "iterations for maximal separating LP iterations per separation round (negative for no limit)",
3265 &sepadata->perroundlpiterlimitfactor, TRUE, DEFAULT_PERROUNDLPITERLIMITFACTOR, -1.0, SCIP_REAL_MAX, NULL,
3268 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/perroundcutsfactorroot", "factor w.r.t. number of integer "
3269 "columns for number of cuts separated per separation round in root node", &sepadata->perroundcutsfactorroot,
3272 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/perroundcutsfactor", "factor w.r.t. number of integer "
3274 &sepadata->perroundcutsfactor, TRUE, DEFAULT_PERROUNDCUTSFACTOR, 0.0, SCIP_REAL_MAX, NULL, NULL) );
3276 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/totalcutsfactor", "factor w.r.t. number of integer "
3277 "columns for total number of cuts separated", &sepadata->totalcutsfactor, TRUE, DEFAULT_TOTALCUTSFACTOR, 0.0,
3280 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/muparaminit", "initial value of the mu parameter (factor "
3281 "for step length)", &sepadata->muparaminit, TRUE, DEFAULT_MUPARAMINIT, 0.0, 100.0, NULL, NULL) );
3283 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/muparamlb", "lower bound of the mu parameter (factor for "
3286 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/muparamub", "upper bound of the mu parameter (factor for "
3289 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/mubacktrackfactor", "factor of mu while backtracking the "
3290 "mu parameter (factor for step length)", &sepadata->mubacktrackfactor, TRUE, DEFAULT_MUBACKTRACKFACTOR,
3293 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/muslab1factor", "factor of mu parameter (factor for step "
3294 "length) for larger increment" , &sepadata->muslab1factor, TRUE, DEFAULT_MUSLAB1FACTOR, 0.0, SCIP_REAL_MAX, NULL,
3297 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/muslab2factor", "factor of mu parameter (factor for step "
3298 "length) for smaller increment", &sepadata->muslab2factor, TRUE, DEFAULT_MUSLAB2FACTOR, 0.0, SCIP_REAL_MAX, NULL,
3301 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/muslab3factor", "factor of mu parameter (factor for step "
3302 "length) for reduction", &sepadata->muslab3factor, TRUE, DEFAULT_MUSLAB3FACTOR, 0.0, SCIP_REAL_MAX, NULL, NULL) );
3304 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/deltaslab1ub", "factor of delta deciding larger increment "
3305 "of mu parameter (factor for step length)", &sepadata->deltaslab1ub, TRUE, DEFAULT_DELTASLAB1UB, 0.0, 1.0,
3308 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/deltaslab2ub", "factor of delta deciding smaller "
3309 "increment of mu parameter (factor for step length)", &sepadata->deltaslab2ub, TRUE, DEFAULT_DELTASLAB2UB,
3312 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/ubparamposfactor", "factor for positive upper bound used "
3313 "as an estimate for the optimal Lagrangian dual value", &sepadata->ubparamposfactor, TRUE, DEFAULT_UBPARAMPOSFACTOR,
3316 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/ubparamnegfactor", "factor for negative upper bound used "
3317 "as an estimate for the optimal Lagrangian dual value", &sepadata->ubparamnegfactor, TRUE, DEFAULT_UBPARAMNEGFACTOR,
3320 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/perrootlpiterfactor", "factor w.r.t. root node LP "
3322 &sepadata->perrootlpiterfactor, TRUE, DEFAULT_PERROOTLPITERFACTOR, -1.0, SCIP_REAL_MAX, NULL, NULL) );
3324 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/perlpiterfactor", "factor w.r.t. non-root node LP "
3325 "iterations for iteration limit of each separating LP (negative for no limit)", &sepadata->perlpiterfactor, TRUE,
3328 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/cutsfilterfactor", "fraction of generated cuts per "
3329 "explored basis to accept from separator", &sepadata->cutsfilterfactor, TRUE, DEFAULT_CUTSFILTERFACTOR, 0.0,
3332 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/radiusinit", "initial radius of the ball used in "
3333 "stabilization of Lagrangian multipliers", &sepadata->radiusinit, TRUE, DEFAULT_RADIUSINIT, 0.0, 1.0, NULL,
3336 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/radiusmax", "maximum radius of the ball used in "
3337 "stabilization of Lagrangian multipliers", &sepadata->radiusmax, TRUE, DEFAULT_RADIUSMAX, 0.0, SCIP_REAL_MAX,
3340 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/radiusmin", "minimum radius of the ball used in "
3341 "stabilization of Lagrangian multipliers", &sepadata->radiusmin, TRUE, DEFAULT_RADIUSMIN, 0.0, SCIP_REAL_MAX,
3344 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/constant", "a constant for stablity center based "
3345 "stabilization of Lagrangian multipliers", &sepadata->constant, TRUE, DEFAULT_CONST, 2.0, SCIP_REAL_MAX,
3348 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/radiusupdateweight", "multiplier to evaluate cut "
3352 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/dualdegeneracyratethreshold", "minimum dual degeneracy "
3356 SCIP_CALL( SCIPaddRealParam(scip, "separating/" SEPA_NAME "/varconsratiothreshold", "minimum variable-constraint "
3360 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/muparamconst", "is the mu parameter (factor for step "
3363 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/separaterows", "separate rows with integral slack?",
3366 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/sortcutoffsol", "sort fractional integer columns"
3367 "based on fractionality?" , &sepadata->sortcutoffsol, TRUE, DEFAULT_SORTCUTOFFSOL, NULL, NULL) );
3369 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/sidetypebasis", "choose side types of row (lhs/rhs) "
3370 "based on basis information?", &sepadata->sidetypebasis, TRUE, DEFAULT_SIDETYPEBASIS, NULL, NULL) );
3372 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/dynamiccuts", "should generated cuts be removed from "
3373 "LP if they are no longer tight?", &sepadata->dynamiccuts, FALSE, DEFAULT_DYNAMICCUTS, NULL, NULL) );
3375 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/makeintegral", "try to scale all cuts to integral "
3378 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/forcecuts", "force cuts to be added to the LP?",
3381 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/delayedcuts", "should cuts be added to the delayed cut "
3384 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/allowlocal", "should locally valid cuts be generated?",
3387 SCIP_CALL( SCIPaddBoolParam(scip, "separating/" SEPA_NAME "/aggregatecuts", "aggregate all generated cuts using the "
3388 "Lagrangian multipliers?", &sepadata->aggregatecuts, TRUE, DEFAULT_AGGREGATECUTS, NULL, NULL) );
3390 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/maxrounds", "maximal number of separation rounds per node "
3393 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/maxroundsroot", "maximal number of separation rounds in "
3394 "the root node (-1: unlimited)", &sepadata->maxroundsroot, FALSE, DEFAULT_MAXROUNDSROOT, -1, INT_MAX, NULL,
3397 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/perroundnmaxlpiters", "maximal number of separating LP "
3401 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/nmaxcutsperlp", "maximal number of cuts separated per "
3402 "Lagromory LP in the non-root node", &sepadata->nmaxcutsperlp, FALSE, DEFAULT_PERLPMAXCUTS, 0, INT_MAX, NULL,
3405 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/nmaxcutsperlproot", "maximal number of cuts separated per "
3406 "Lagromory LP in the root node", &sepadata->nmaxcutsperlproot, FALSE, DEFAULT_PERLPMAXCUTSROOT, 0, INT_MAX,
3409 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/nmaxmainiters", "maximal number of main loop iterations of "
3410 "the relax-and-cut algorithm", &sepadata->nmaxmainiters, TRUE, DEFAULT_MAXMAINITERS, 0, INT_MAX, NULL, NULL) );
3412 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/nmaxsubgradientiters", "maximal number of subgradient loop "
3416 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/cutgenfreq", "frequency of subgradient iterations for "
3419 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/cutaddfreq", "frequency of subgradient iterations for "
3420 "adding cuts to objective function", &sepadata->cutaddfreq, TRUE, DEFAULT_CUTADDFREQ, 0, INT_MAX, NULL, NULL) );
3422 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/nmaxlagrangianvalsforavg", "maximal number of iterations "
3426 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/nmaxconsecitersformuupdate", "consecutive number of "
3427 "iterations used to determine if mu needs to be backtracked", &sepadata->nmaxconsecitersformuupdate, TRUE,
3430 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/projectiontype", "the ball into which the Lagrangian multipliers "
3431 "are projected for stabilization (0: no projection, 1: L1-norm ball projection, 2: L2-norm ball projection, 3: "
3432 "L_inf-norm ball projection)", &sepadata->projectiontype, TRUE, DEFAULT_PROJECTIONTYPE, 0, 3, NULL, NULL));
3434 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/stabilitycentertype", "type of stability center for "
3435 "taking weighted average of Lagrangian multipliers for stabilization (0: no weighted stabilization, 1: best "
3436 "Lagrangian multipliers)", &sepadata->stabilitycentertype, TRUE, DEFAULT_STABILITYCENTERTYPE, 0, 1, NULL,
3439 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/optimalfacepriority", "priority of the optimal face for "
3443 SCIP_CALL( SCIPaddIntParam(scip, "separating/" SEPA_NAME "/minrestart", "minimum restart round for separator "
3444 "execution (0: from beginning of the instance solving, >= n with n >= 1: from restart round n)",
void SCIPfreeRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen)
Definition: scip_randnumgen.c:79
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)
Definition: sepa_lagromory.c:1419
Definition: type_result.h:42
SCIP_RETCODE SCIPgetLPBInvRow(SCIP *scip, int r, SCIP_Real *coefs, int *inds, int *ninds)
Definition: scip_lp.c:714
primal heuristic that tries a given solution
SCIP_RETCODE SCIPcacheRowExtensions(SCIP *scip, SCIP_ROW *row)
Definition: scip_lp.c:1635
SCIP_RETCODE SCIPlpiSetState(SCIP_LPI *lpi, BMS_BLKMEM *blkmem, const SCIP_LPISTATE *lpistate)
Definition: lpi_clp.cpp:3429
static SCIP_RETCODE l1BallProjection(SCIP *scip, SCIP_Real *dualvector, int dualvectorlen, SCIP_Real radius)
Definition: sepa_lagromory.c:1020
SCIP_Longint SCIPgetNLPIterations(SCIP *scip)
Definition: scip_solvingstats.c:448
Definition: struct_scip.h:69
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)
Definition: sepa_lagromory.c:833
#define DEFAULT_PERROUNDCUTSFACTORROOT
Definition: sepa_lagromory.c:183
static SCIP_RETCODE addGMICutsAsSoftConss(SCIP_Real *dualvector, int ngeneratedcuts, int *naddedcuts, int *nnewaddedsoftcuts)
Definition: sepa_lagromory.c:2297
SCIP_RETCODE SCIPflushRowExtensions(SCIP *scip, SCIP_ROW *row)
Definition: scip_lp.c:1658
SCIP_RETCODE SCIPlpiGetSol(SCIP_LPI *lpi, SCIP_Real *objval, SCIP_Real *primsol, SCIP_Real *dualsol, SCIP_Real *activity, SCIP_Real *redcost)
Definition: lpi_clp.cpp:2788
SCIP_RETCODE SCIPgetRealParam(SCIP *scip, const char *name, SCIP_Real *value)
Definition: scip_param.c:307
Definition: struct_cuts.h:40
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)
Definition: sepa_lagromory.c:2022
SCIP_RETCODE SCIPaddVarToRow(SCIP *scip, SCIP_ROW *row, SCIP_VAR *var, SCIP_Real val)
Definition: scip_lp.c:1701
SCIP_RETCODE SCIPgetLPDualDegeneracy(SCIP *scip, SCIP_Real *degeneracy, SCIP_Real *varconsratio)
Definition: scip_lp.c:2792
static SCIP_RETCODE createLPWithSoftCuts(SCIP *scip, SCIP_SEPADATA *sepadata)
Definition: sepa_lagromory.c:376
SCIP_Bool SCIPisGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:497
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)
Definition: sepa_lagromory.c:2321
Definition: struct_var.h:207
SCIP_Bool SCIPisFeasNegative(SCIP *scip, SCIP_Real val)
Definition: scip_numerics.c:869
SCIP_RETCODE SCIPaddDelayedPoolCut(SCIP *scip, SCIP_ROW *row)
Definition: scip_cut.c:641
static SCIP_RETCODE sepadataCreate(SCIP *scip, SCIP_SEPADATA **sepadata)
Definition: sepa_lagromory.c:3041
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)
Definition: sepa_lagromory.c:2559
Definition: struct_misc.h:268
SCIP_RETCODE SCIPlpiSetRealpar(SCIP_LPI *lpi, SCIP_LPPARAM type, SCIP_Real dval)
Definition: lpi_clp.cpp:3833
SCIP_Longint SCIPgetNRootFirstLPIterations(SCIP *scip)
Definition: scip_solvingstats.c:511
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)
Definition: sepa_lagromory.c:1667
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)
Definition: sepa_lagromory.c:951
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)
Definition: sepa_lagromory.c:1466
Definition: type_result.h:49
static SCIP_RETCODE updateObjectiveVector(SCIP *scip, SCIP_Real *dualvector, SCIP_ROW **cuts, int ncuts, SCIP_Real *origobjcoefs, SCIP_Bool *objvecsdiffer)
Definition: sepa_lagromory.c:2212
void SCIPsortDownRealInt(SCIP_Real *realarray, int *intarray, int len)
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:445
Definition: struct_sepa.h:46
SCIP_RETCODE SCIPgetLPColsData(SCIP *scip, SCIP_COL ***cols, int *ncols)
Definition: scip_lp.c:471
SCIP_RETCODE SCIPsetSepaCopy(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPACOPY((*sepacopy)))
Definition: scip_sepa.c:151
SCIP_RETCODE SCIPaddIntParam(SCIP *scip, const char *name, const char *desc, int *valueptr, SCIP_Bool isadvanced, int defaultvalue, int minvalue, int maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:83
SCIP_RETCODE SCIPlpiCreate(SCIP_LPI **lpi, SCIP_MESSAGEHDLR *messagehdlr, const char *name, SCIP_OBJSEN objsen)
Definition: lpi_clp.cpp:531
SCIP_RETCODE SCIPlpiAddCols(SCIP_LPI *lpi, int ncols, const SCIP_Real *obj, const SCIP_Real *lb, const SCIP_Real *ub, char **colnames, int nnonz, const int *beg, const int *ind, const SCIP_Real *val)
Definition: lpi_clp.cpp:758
SCIP_Real SCIPgetRowMaxActivity(SCIP *scip, SCIP_ROW *row)
Definition: scip_lp.c:1956
Definition: struct_lp.h:135
Definition: struct_sol.h:73
Definition: type_lpi.h:61
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)
Definition: sepa_lagromory.c:1250
SCIP_RETCODE SCIPheurPassSolTrySol(SCIP *scip, SCIP_HEUR *heur, SCIP_SOL *sol)
Definition: heur_trysol.c:252
#define SCIPallocCleanBufferArray(scip, ptr, num)
Definition: scip_mem.h:142
SCIP_RETCODE SCIPlpiAddRows(SCIP_LPI *lpi, int nrows, const SCIP_Real *lhs, const SCIP_Real *rhs, char **rownames, int nnonz, const int *beg, const int *ind, const SCIP_Real *val)
Definition: lpi_clp.cpp:914
#define DEFAULT_MAXCONSECITERSFORMUUPDATE
Definition: sepa_lagromory.c:223
Definition: type_result.h:44
static void updateLagrangianValue(SCIP *scip, SCIP_Real objval, SCIP_Real *dualvector, SCIP_ROW **cuts, int ncuts, SCIP_Real *lagrangianval)
Definition: sepa_lagromory.c:904
SCIP_RETCODE SCIPsolveDiveLP(SCIP *scip, int itlim, SCIP_Bool *lperror, SCIP_Bool *cutoff)
Definition: scip_lp.c:2678
SCIP_Bool SCIPisLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:458
SCIP_RETCODE SCIPsetSepaExit(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPAEXIT((*sepaexit)))
Definition: scip_sepa.c:199
SCIP_RETCODE SCIPincludeSepaLagromory(SCIP *scip)
Definition: sepa_lagromory.c:3247
#define DEFAULT_MAXLAGRANGIANVALSFORAVG
Definition: sepa_lagromory.c:222
Definition: type_lpi.h:43
SCIP_RETCODE SCIPlpiFreeState(SCIP_LPI *lpi, BMS_BLKMEM *blkmem, SCIP_LPISTATE **lpistate)
Definition: lpi_clp.cpp:3503
SCIP_Bool SCIPisEfficacious(SCIP *scip, SCIP_Real efficacy)
Definition: scip_cut.c:135
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)
Definition: sepa_lagromory.c:1799
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)
Definition: sepa_lagromory.c:2673
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)
Definition: sepa_lagromory.c:1310
Definition: type_retcode.h:42
static SCIP_RETCODE solveLPWithHardCuts(SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_Bool *solfound, SCIP_SOL *sol, SCIP_Real *solvals)
Definition: sepa_lagromory.c:1603
SCIP_RETCODE SCIPaddRow(SCIP *scip, SCIP_ROW *row, SCIP_Bool forcecut, SCIP_Bool *infeasible)
Definition: scip_cut.c:250
Definition: struct_heur.h:97
Definition: type_lp.h:43
#define DEFAULT_TOTALLPITERLIMITFACTOR
Definition: sepa_lagromory.c:179
SCIP_RETCODE SCIPincludeSepaBasic(SCIP *scip, SCIP_SEPA **sepa, const char *name, const char *desc, int priority, int freq, SCIP_Real maxbounddist, SCIP_Bool usessubscip, SCIP_Bool delay, SCIP_DECL_SEPAEXECLP((*sepaexeclp)), SCIP_DECL_SEPAEXECSOL((*sepaexecsol)), SCIP_SEPADATA *sepadata)
Definition: scip_sepa.c:109
static void weightedDualVector(SCIP_SEPADATA *sepadata, SCIP_Real *dualvector, int dualvectorlen, SCIP_Real *stabilitycenter, int stabilitycenterlen, int nbestdualupdates, int totaliternum)
Definition: sepa_lagromory.c:1217
SCIP_RETCODE SCIPcreateRandom(SCIP *scip, SCIP_RANDNUMGEN **randnumgen, unsigned int initialseed, SCIP_Bool useglobalseed)
Definition: scip_randnumgen.c:56
SCIP_RETCODE SCIPcalcMIR(SCIP *scip, SCIP_SOL *sol, SCIP_Bool postprocess, SCIP_Real boundswitch, SCIP_Bool usevbds, SCIP_Bool allowlocal, SCIP_Bool fixintegralrhs, int *boundsfortrans, SCIP_BOUNDTYPE *boundtypesfortrans, SCIP_Real minfrac, SCIP_Real maxfrac, SCIP_Real scale, SCIP_AGGRROW *aggrrow, SCIP_Real *cutcoefs, SCIP_Real *cutrhs, int *cutinds, int *cutnnz, SCIP_Real *cutefficacy, int *cutrank, SCIP_Bool *cutislocal, SCIP_Bool *success)
Definition: cuts.c:3879
SCIP_RETCODE SCIPsetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var, SCIP_Real val)
Definition: scip_sol.c:1077
#define DEFAULT_DUALDEGENERACYRATETHRESHOLD
Definition: sepa_lagromory.c:166
static SCIP_RETCODE separateCuts(SCIP *scip, SCIP_SEPA *sepa, SCIP_SEPADATA *sepadata, SCIP_Real ubparam, int depth, SCIP_Bool allowlocal, SCIP_RESULT *result)
Definition: sepa_lagromory.c:2721
lagromory separator
SCIP_RETCODE SCIPlpiGetState(SCIP_LPI *lpi, BMS_BLKMEM *blkmem, SCIP_LPISTATE **lpistate)
Definition: lpi_clp.cpp:3389
Definition: struct_lp.h:201
SCIP_RETCODE SCIPcreateEmptyRowSepa(SCIP *scip, SCIP_ROW **row, SCIP_SEPA *sepa, const char *name, SCIP_Real lhs, SCIP_Real rhs, SCIP_Bool local, SCIP_Bool modifiable, SCIP_Bool removable)
Definition: scip_lp.c:1453
static void l2BallProjection(SCIP *scip, SCIP_Real *dualvector, int dualvectorlen, SCIP_Real radius)
Definition: sepa_lagromory.c:1160
static SCIP_RETCODE createLPWithHardCuts(SCIP *scip, SCIP_SEPADATA *sepadata, SCIP_ROW **cuts, int ncuts)
Definition: sepa_lagromory.c:474
SCIP_Real SCIPgetRowSolActivity(SCIP *scip, SCIP_ROW *row, SCIP_SOL *sol)
Definition: scip_lp.c:2144
static void linfBallProjection(SCIP *scip, SCIP_Real *dualvector, int dualvectorlen, SCIP_Real radius)
Definition: sepa_lagromory.c:1191
SCIP_Real SCIPrandomGetReal(SCIP_RANDNUMGEN *randnumgen, SCIP_Real minrandval, SCIP_Real maxrandval)
Definition: misc.c:10130
SCIP_Longint SCIPgetNNodeInitLPIterations(SCIP *scip)
Definition: scip_solvingstats.c:823
SCIP_RETCODE SCIPsetSepaFree(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPAFREE((*sepafree)))
Definition: scip_sepa.c:167
SCIP_RETCODE SCIPsetSepaInit(SCIP *scip, SCIP_SEPA *sepa, SCIP_DECL_SEPAINIT((*sepainit)))
Definition: scip_sepa.c:183
SCIP_Bool SCIPisGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:484
SCIP_RETCODE SCIPaggrRowSumRows(SCIP *scip, SCIP_AGGRROW *aggrrow, SCIP_Real *weights, int *rowinds, int nrowinds, SCIP_Bool sidetypebasis, SCIP_Bool allowlocal, int negslack, int maxaggrlen, SCIP_Bool *valid)
Definition: cuts.c:2287
SCIP_RETCODE SCIPchgVarObjDive(SCIP *scip, SCIP_VAR *var, SCIP_Real newobj)
Definition: scip_lp.c:2378
Definition: lpi_clp.cpp:132
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)
Definition: sepa_lagromory.c:2596
static SCIP_RETCODE sepadataFree(SCIP *scip, SCIP_SEPADATA **sepadata)
Definition: sepa_lagromory.c:673
SCIP_Bool SCIPisFeasPositive(SCIP *scip, SCIP_Real val)
Definition: scip_numerics.c:857
static void updateStepLength(SCIP *scip, SCIP_Real muparam, SCIP_Real ubparam, SCIP_Real lagrangianval, SCIP_Real *subgradient, int ncuts, SCIP_Real *steplength)
Definition: sepa_lagromory.c:929
Definition: lpi_clp.cpp:104
SCIP_Real SCIPgetRowMinActivity(SCIP *scip, SCIP_ROW *row)
Definition: scip_lp.c:1939
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)
Definition: sepa_lagromory.c:709
SCIP_RETCODE SCIPaggrRowCreate(SCIP *scip, SCIP_AGGRROW **aggrrow)
Definition: cuts.c:1731
#define DEFAULT_PERROUNDLPITERLIMITFACTOR
Definition: sepa_lagromory.c:173
SCIP_Bool SCIPisLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
Definition: scip_numerics.c:471
SCIP_RETCODE SCIPgetLPRowsData(SCIP *scip, SCIP_ROW ***rows, int *nrows)
Definition: scip_lp.c:570
static SCIP_RETCODE deleteLPWithSoftCuts(SCIP *scip, SCIP_SEPADATA *sepadata)
Definition: sepa_lagromory.c:457
Definition: objbenders.h:43
SCIP_RETCODE SCIPmakeRowIntegral(SCIP *scip, SCIP_ROW *row, SCIP_Real mindelta, SCIP_Real maxdelta, SCIP_Longint maxdnom, SCIP_Real maxscale, SCIP_Bool usecontvars, SCIP_Bool *success)
Definition: scip_lp.c:1844
SCIP_Real SCIPgetSolVal(SCIP *scip, SCIP_SOL *sol, SCIP_VAR *var)
Definition: scip_sol.c:1217
static SCIP_DECL_SEPAEXECLP(sepaExeclpLagromory)
Definition: sepa_lagromory.c:3127
SCIP_RETCODE SCIPaddRealParam(SCIP *scip, const char *name, const char *desc, SCIP_Real *valueptr, SCIP_Bool isadvanced, SCIP_Real defaultvalue, SCIP_Real minvalue, SCIP_Real maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:139
Definition: type_result.h:48
SCIP_RETCODE SCIPaddBoolParam(SCIP *scip, const char *name, const char *desc, SCIP_Bool *valueptr, SCIP_Bool isadvanced, SCIP_Bool defaultvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:57
SCIP_RETCODE SCIPcreateSol(SCIP *scip, SCIP_SOL **sol, SCIP_HEUR *heur)
Definition: scip_sol.c:184
Definition: type_var.h:71