Scippy

SCIP

Solving Constraint Integer Programs

presol_milp.cpp
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the program and library */
4/* SCIP --- Solving Constraint Integer Programs */
5/* */
6/* Copyright (c) 2002-2024 Zuse Institute Berlin (ZIB) */
7/* */
8/* Licensed under the Apache License, Version 2.0 (the "License"); */
9/* you may not use this file except in compliance with the License. */
10/* You may obtain a copy of the License at */
11/* */
12/* http://www.apache.org/licenses/LICENSE-2.0 */
13/* */
14/* Unless required by applicable law or agreed to in writing, software */
15/* distributed under the License is distributed on an "AS IS" BASIS, */
16/* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. */
17/* See the License for the specific language governing permissions and */
18/* limitations under the License. */
19/* */
20/* You should have received a copy of the Apache-2.0 license */
21/* along with SCIP; see the file LICENSE. If not visit scipopt.org. */
22/* */
23/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
24
25/**@file presol_milp.cpp
26 * @brief MILP presolver
27 * @author Leona Gottwald
28 *
29 * Calls the presolve library and communicates (multi-)aggregations, fixings, and bound
30 * changes to SCIP by utilizing the postsolve information. Constraint changes can currently
31 * only be communicated by deleting all constraints and adding new ones.
32 *
33 * @todo add infrastructure to SCIP for handling parallel columns
34 * @todo better communication of constraint changes by adding more information to the postsolve structure
35 * @todo allow to pass additional external locks to the presolve library that are considered when doing reductions
36 *
37 */
38
39/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
40
41#include "scip/presol_milp.h"
42
43#ifndef SCIP_WITH_PAPILO
44
45/** creates the MILP presolver and includes it in SCIP */
47 SCIP* scip /**< SCIP data structure */
48 )
49{
50 assert(scip != NULL);
51 return SCIP_OKAY;
52}
53
54#else
55
56/* disable some warnings that come up in header files of PAPILOs dependencies */
57#ifdef __GNUC__
58#pragma GCC diagnostic ignored "-Wshadow"
59#pragma GCC diagnostic ignored "-Wctor-dtor-privacy"
60#pragma GCC diagnostic ignored "-Wredundant-decls"
61
62/* disable false warning, !3076, https://gcc.gnu.org/bugzilla/show_bug.cgi?id=106199 */
63#if __GNUC__ == 12 && __GNUC__MINOR__ <= 2
64#pragma GCC diagnostic ignored "-Wstringop-overflow"
65#endif
66#endif
67
68#include <assert.h>
69#include "scip/cons_linear.h"
70#include "scip/pub_matrix.h"
71#include "scip/pub_presol.h"
72#include "scip/pub_var.h"
73#include "scip/pub_cons.h"
74#include "scip/pub_message.h"
75#include "scip/scip_general.h"
76#include "scip/scip_presol.h"
77#include "scip/scip_var.h"
78#include "scip/scip_mem.h"
79#include "scip/scip_prob.h"
80#include "scip/scip_param.h"
81#include "scip/scip_cons.h"
82#include "scip/scip_numerics.h"
83#include "scip/scip_timing.h"
84#include "scip/scip_message.h"
86#include "papilo/core/Presolve.hpp"
87#include "papilo/core/ProblemBuilder.hpp"
88#include "papilo/Config.hpp"
89
90#define PRESOL_NAME "milp"
91#define PRESOL_DESC "MILP specific presolving methods"
92#define PRESOL_PRIORITY 9999999 /**< priority of the presolver (>= 0: before, < 0: after constraint handlers); combined with propagators */
93#define PRESOL_MAXROUNDS -1 /**< maximal number of presolving rounds the presolver participates in (-1: no limit) */
94#define PRESOL_TIMING SCIP_PRESOLTIMING_MEDIUM /* timing of the presolver (fast, medium, or exhaustive) */
95
96/* default parameter values */
97#define DEFAULT_THREADS 1 /**< maximum number of threads presolving may use (0: automatic) */
98#define DEFAULT_MAXFILLINPERSUBST 3 /**< maximal possible fillin for substitutions to be considered */
99#define DEFAULT_MAXSHIFTPERROW 10 /**< maximal amount of nonzeros allowed to be shifted to make space for substitutions */
100#define DEFAULT_DETECTLINDEP 0 /**< should linear dependent equations and free columns be removed? (0: never, 1: for LPs, 2: always) */
101#define DEFAULT_MAXBADGESIZE_SEQ 15000 /**< the max badge size in Probing if PaPILO is executed in sequential mode */
102#define DEFAULT_MAXBADGESIZE_PAR -1 /**< the max badge size in Probing if PaPILO is executed in parallel mode */
103#define DEFAULT_INTERNAL_MAXROUNDS -1 /**< internal max rounds in PaPILO (-1: no limit, 0: model cleanup) */
104#define DEFAULT_RANDOMSEED 0 /**< the random seed used for randomization of tie breaking */
105#define DEFAULT_MODIFYCONSFAC 0.8 /**< modify SCIP constraints when the number of nonzeros or rows is at most this
106 * factor times the number of nonzeros or rows before presolving */
107#define DEFAULT_MARKOWITZTOLERANCE 0.01 /**< the markowitz tolerance used for substitutions */
108#define DEFAULT_VERBOSITY 0
109#define DEFAULT_HUGEBOUND 1e8 /**< absolute bound value that is considered too huge for activitity based calculations */
110#define DEFAULT_ENABLEPARALLELROWS TRUE /**< should the parallel rows presolver be enabled within the presolve library? */
111#define DEFAULT_ENABLEDOMCOL TRUE /**< should the dominated column presolver be enabled within the presolve library? */
112#define DEFAULT_ENABLEDUALINFER TRUE /**< should the dualinfer presolver be enabled within the presolve library? */
113#define DEFAULT_ENABLEMULTIAGGR TRUE /**< should the multi-aggregation presolver be enabled within the presolve library? */
114#define DEFAULT_ENABLEPROBING TRUE /**< should the probing presolver be enabled within the presolve library? */
115#define DEFAULT_ENABLESPARSIFY FALSE /**< should the sparsify presolver be enabled within the presolve library? */
116#define DEFAULT_FILENAME_PROBLEM "-" /**< default filename to store the instance before presolving */
117
118/*
119 * Data structures
120 */
121
122/** presolver data */
123struct SCIP_PresolData
124{
125 int lastncols; /**< the number of columns from the last call */
126 int lastnrows; /**< the number of rows from the last call */
127 int threads; /**< maximum number of threads presolving may use (0: automatic) */
128 int maxfillinpersubstitution; /**< maximal possible fillin for substitutions to be considered */
129 int maxbadgesizeseq; /**< the max badge size in Probing if PaPILO is called in sequential mode */
130 int maxbadgesizepar; /**< the max badge size in Probing if PaPILO is called in parallel mode */
131 int internalmaxrounds; /**< internal max rounds in PaPILO (-1: no limit, 0: model cleanup) */
132 int maxshiftperrow; /**< maximal amount of nonzeros allowed to be shifted to make space for substitutions */
133 int detectlineardependency; /**< should linear dependent equations and free columns be removed? (0: never, 1: for LPs, 2: always) */
134 int randomseed; /**< the random seed used for randomization of tie breaking */
135 int verbosity;
136
137 SCIP_Bool enablesparsify; /**< should the sparsify presolver be enabled within the presolve library? */
138 SCIP_Bool enabledomcol; /**< should the dominated column presolver be enabled within the presolve library? */
139 SCIP_Bool enableprobing; /**< should the probing presolver be enabled within the presolve library? */
140 SCIP_Bool enabledualinfer; /**< should the dualinfer presolver be enabled within the presolve library? */
141 SCIP_Bool enablemultiaggr; /**< should the multi-aggregation presolver be enabled within the presolve library? */
142 SCIP_Bool enableparallelrows; /**< should the parallel rows presolver be enabled within the presolve library? */
143 SCIP_Real modifyconsfac; /**< modify SCIP constraints when the number of nonzeros or rows is at most this
144 * factor times the number of nonzeros or rows before presolving */
145 SCIP_Real markowitztolerance; /**< the markowitz tolerance used for substitutions */
146 SCIP_Real hugebound; /**< absolute bound value that is considered too huge for activitity based calculations */
147
148 char* filename = NULL; /**< filename to store the instance before presolving */
149};
150
151using namespace papilo;
152
153/*
154 * Local methods
155 */
156
157/** builds the presolvelib problem datastructure from the matrix */
158static
159Problem<SCIP_Real> buildProblem(
160 SCIP* scip, /**< SCIP data structure */
161 SCIP_MATRIX* matrix /**< initialized SCIP_MATRIX data structure */
162 )
163{
164 ProblemBuilder<SCIP_Real> builder;
165
166 /* build problem from matrix */
167 int nnz = SCIPmatrixGetNNonzs(matrix);
168 int ncols = SCIPmatrixGetNColumns(matrix);
169 int nrows = SCIPmatrixGetNRows(matrix);
170 builder.reserve(nnz, nrows, ncols);
171
172 /* set up columns */
173 builder.setNumCols(ncols);
174 for(int i = 0; i != ncols; ++i)
175 {
176 SCIP_VAR* var = SCIPmatrixGetVar(matrix, i);
179 builder.setColLb(i, lb);
180 builder.setColUb(i, ub);
181 builder.setColLbInf(i, SCIPisInfinity(scip, -lb));
182 builder.setColUbInf(i, SCIPisInfinity(scip, ub));
183#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 1)
185 builder.setColImplInt(i, TRUE);
186 else
187 builder.setColIntegral(i, SCIPvarIsIntegral(var));
188#else
189 builder.setColIntegral(i, SCIPvarIsIntegral(var));
190#endif
191 builder.setObj(i, SCIPvarGetObj(var));
192 }
193
194 /* set up rows */
195 builder.setNumRows(nrows);
196 for(int i = 0; i != nrows; ++i)
197 {
198 int* rowcols = SCIPmatrixGetRowIdxPtr(matrix, i);
199 SCIP_Real* rowvals = SCIPmatrixGetRowValPtr(matrix, i);
200 int rowlen = SCIPmatrixGetRowNNonzs(matrix, i);
201 builder.addRowEntries(i, rowlen, rowcols, rowvals);
202
203 SCIP_Real lhs = SCIPmatrixGetRowLhs(matrix, i);
204 SCIP_Real rhs = SCIPmatrixGetRowRhs(matrix, i);
205 builder.setRowLhs(i, lhs);
206 builder.setRowRhs(i, rhs);
207 builder.setRowLhsInf(i, SCIPisInfinity(scip, -lhs));
208 builder.setRowRhsInf(i, SCIPisInfinity(scip, rhs));
209 }
210
211 /* init objective offset - the value itself is irrelevant */
212 builder.setObjOffset(0);
213
214#ifdef SCIP_PRESOLLIB_ENABLE_OUTPUT
215 /* show problem name */
216 builder.setProblemName(SCIPgetProbName(scip));
217#endif
218
219 return builder.build();
220}
221
222/** sets up the presolvelib presolve datastructure from the data */
223static
224Presolve<SCIP_Real> setupPresolve(
225 SCIP* scip, /**< SCIP data structure */
226 SCIP_PRESOLDATA* data, /**< presolver data structure */
227 SCIP_Bool allowconsmodification /**< whether constraint modifications are allowed */
228 )
229{
230 Presolve<SCIP_Real> presolve;
231 SCIP_Real timelimit;
232
233 /* important so that SCIP does not throw an error, e.g. when an integer variable is substituted
234 * into a knapsack constraint */
235 presolve.getPresolveOptions().substitutebinarieswithints = false;
236
237 /* currently these changes cannot be communicated to SCIP correctly since a constraint needs
238 * to be modified in the cases where slackvariables are removed from constraints but for the
239 * presolve library those look like normal substitution on the postsolve stack */
240 presolve.getPresolveOptions().removeslackvars = false;
241
242 /* communicate the SCIP parameters to the presolve library */
243 presolve.getPresolveOptions().maxfillinpersubstitution = data->maxfillinpersubstitution;
244 presolve.getPresolveOptions().markowitz_tolerance = data->markowitztolerance;
245 presolve.getPresolveOptions().maxshiftperrow = data->maxshiftperrow;
246 presolve.getPresolveOptions().hugeval = data->hugebound;
247
248 /* removal of linear dependent equations has only an effect when constraint modifications are communicated */
249 presolve.getPresolveOptions().detectlindep = allowconsmodification ? data->detectlineardependency : 0;
250
251 /* communicate the random seed */
252 presolve.getPresolveOptions().randomseed = SCIPinitializeRandomSeed(scip, (unsigned int)data->randomseed);
253
254#ifdef PAPILO_TBB
255 /* set number of threads to be used for presolve */
256 presolve.getPresolveOptions().threads = data->threads;
257#else
258 if (data->threads != DEFAULT_THREADS)
260 "PaPILO can utilize only multiple threads if it is build with TBB.\n");
261 presolve.getPresolveOptions().threads = 1;
262#endif
263
264#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 3)
265 presolve.getPresolveOptions().maxrounds = data->internalmaxrounds;
266#endif
267
268 /* disable dual reductions that are not permitted */
270 presolve.getPresolveOptions().dualreds = 2;
271 else if( SCIPallowWeakDualReds(scip) )
272 presolve.getPresolveOptions().dualreds = 1;
273 else
274 presolve.getPresolveOptions().dualreds = 0;
275
276 /* set up the presolvers that shall participate */
277 using uptr = std::unique_ptr<PresolveMethod<SCIP_Real>>;
278
279 /* fast presolvers*/
280 presolve.addPresolveMethod( uptr( new SingletonCols<SCIP_Real>() ) );
281 presolve.addPresolveMethod( uptr( new CoefficientStrengthening<SCIP_Real>() ) );
282 presolve.addPresolveMethod( uptr( new ConstraintPropagation<SCIP_Real>() ) );
283
284 /* medium presolver */
285 presolve.addPresolveMethod( uptr( new SimpleProbing<SCIP_Real>() ) );
286 if( data->enableparallelrows )
287 presolve.addPresolveMethod( uptr( new ParallelRowDetection<SCIP_Real>() ) );
288 /* todo: parallel cols cannot be handled by SCIP currently
289 * addPresolveMethod( uptr( new ParallelColDetection<SCIP_Real>() ) ); */
290 presolve.addPresolveMethod( uptr( new SingletonStuffing<SCIP_Real>() ) );
291#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 1)
292 DualFix<SCIP_Real> *dualfix = new DualFix<SCIP_Real>();
293 dualfix->set_fix_to_infinity_allowed(false);
294 presolve.addPresolveMethod( uptr( dualfix ) );
295#else
296 presolve.addPresolveMethod( uptr( new DualFix<SCIP_Real>() ) );
297#endif
298 presolve.addPresolveMethod( uptr( new FixContinuous<SCIP_Real>() ) );
299 presolve.addPresolveMethod( uptr( new SimplifyInequalities<SCIP_Real>() ) );
300 presolve.addPresolveMethod( uptr( new SimpleSubstitution<SCIP_Real>() ) );
301
302 /* exhaustive presolvers*/
303 presolve.addPresolveMethod( uptr( new ImplIntDetection<SCIP_Real>() ) );
304 if( data->enabledualinfer )
305 presolve.addPresolveMethod( uptr( new DualInfer<SCIP_Real>() ) );
306 if( data->enableprobing )
307 {
308#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 1)
309 Probing<SCIP_Real> *probing = new Probing<SCIP_Real>();
310 if( presolve.getPresolveOptions().runs_sequential() )
311 {
312 probing->set_max_badge_size( data->maxbadgesizeseq );
313 }
314 else
315 {
316 probing->set_max_badge_size( data->maxbadgesizepar );
317 }
318 presolve.addPresolveMethod( uptr( probing ) );
319
320#else
321 presolve.addPresolveMethod( uptr( new Probing<SCIP_Real>() ) );
322 if( data->maxbadgesizeseq != DEFAULT_MAXBADGESIZE_SEQ )
324 " The parameter 'presolving/milp/maxbadgesizeseq' can only be used with PaPILO 2.1.0 or later versions.\n");
325
326 if( data->maxbadgesizepar != DEFAULT_MAXBADGESIZE_PAR )
328 " The parameter 'presolving/milp/maxbadgesizepar' can only be used with PaPILO 2.1.0 or later versions.\n");
329
330#endif
331 }
332 if( data->enabledomcol )
333 presolve.addPresolveMethod( uptr( new DominatedCols<SCIP_Real>() ) );
334 if( data->enablemultiaggr )
335 presolve.addPresolveMethod( uptr( new Substitution<SCIP_Real>() ) );
336 if( data->enablesparsify )
337 presolve.addPresolveMethod( uptr( new Sparsify<SCIP_Real>() ) );
338
339 /* set tolerances */
340 presolve.getPresolveOptions().feastol = SCIPfeastol(scip);
341 presolve.getPresolveOptions().epsilon = SCIPepsilon(scip);
342
343#ifndef SCIP_PRESOLLIB_ENABLE_OUTPUT
344 /* adjust output settings of presolve library */
345 presolve.setVerbosityLevel((VerbosityLevel) data->verbosity);
346#endif
347
348 /* communicate the time limit */
349 SCIPgetRealParam(scip, "limits/time", &timelimit);
350 if( !SCIPisInfinity(scip, timelimit) )
351 presolve.getPresolveOptions().tlim = timelimit - SCIPgetSolvingTime(scip);
352
353 return presolve;
354}
355
356/*
357 * Callback methods of presolver
358 */
359
360/** copy method for constraint handler plugins (called when SCIP copies plugins) */
361static
362SCIP_DECL_PRESOLCOPY(presolCopyMILP)
363{ /*lint --e{715}*/
365
366 return SCIP_OKAY;
367}
368
369/** destructor of presolver to free user data (called when SCIP is exiting) */
370static
371SCIP_DECL_PRESOLFREE(presolFreeMILP)
372{ /*lint --e{715}*/
373 SCIP_PRESOLDATA* data = SCIPpresolGetData(presol);
374 assert(data != NULL);
375
376 SCIPpresolSetData(presol, NULL);
378 return SCIP_OKAY;
379}
380
381/** initialization method of presolver (called after problem was transformed) */
382static
383SCIP_DECL_PRESOLINIT(presolInitMILP)
384{ /*lint --e{715}*/
385 SCIP_PRESOLDATA* data = SCIPpresolGetData(presol);
386 assert(data != NULL);
387
388 data->lastncols = -1;
389 data->lastnrows = -1;
390
391 return SCIP_OKAY;
392}
393
394/** execution method of presolver */
395static
396SCIP_DECL_PRESOLEXEC(presolExecMILP)
397{ /*lint --e{715}*/
398 SCIP_MATRIX* matrix;
399 SCIP_PRESOLDATA* data;
400 SCIP_Bool initialized;
401 SCIP_Bool complete;
402 SCIP_Bool infeasible;
403
404 *result = SCIP_DIDNOTRUN;
405
406 data = SCIPpresolGetData(presol);
407
408 int nvars = SCIPgetNVars(scip);
409 int nconss = SCIPgetNConss(scip);
410
411 /* run only if the problem size reduced by some amount since the last call or if it is the first call */
412 if( data->lastncols != -1 && data->lastnrows != -1 &&
413 nvars > data->lastncols * 0.85 &&
414 nconss > data->lastnrows * 0.85 )
415 return SCIP_OKAY;
416
417 SCIP_CALL( SCIPmatrixCreate(scip, &matrix, TRUE, &initialized, &complete, &infeasible,
418 naddconss, ndelconss, nchgcoefs, nchgbds, nfixedvars) );
419
420 /* if infeasibility was detected during matrix creation, return here */
421 if( infeasible )
422 {
423 if( initialized )
424 SCIPmatrixFree(scip, &matrix);
425
426 *result = SCIP_CUTOFF;
427 return SCIP_OKAY;
428 }
429
430 /* we only work on pure MIPs, also disable to try building the matrix again if it failed once */
431 if( !initialized || !complete )
432 {
433 data->lastncols = 0;
434 data->lastnrows = 0;
435
436 if( initialized )
437 SCIPmatrixFree(scip, &matrix);
438
439 return SCIP_OKAY;
440 }
441
442 if( 0 != strncmp(data->filename, DEFAULT_FILENAME_PROBLEM, strlen(DEFAULT_FILENAME_PROBLEM)) )
443 {
445 " writing transformed problem to %s (only enforced constraints)\n", data->filename);
446 SCIP_CALL(SCIPwriteTransProblem(scip, data->filename, NULL, FALSE));
447 }
448
449 /* only allow communication of constraint modifications by deleting all constraints when some already have been upgraded */
450 SCIP_CONSHDLR* linconshdlr = SCIPfindConshdlr(scip, "linear");
451 assert(linconshdlr != NULL);
452 SCIP_Bool allowconsmodification = (SCIPconshdlrGetNCheckConss(linconshdlr) == SCIPmatrixGetNRows(matrix));
453
454 /* store current numbers of aggregations, fixings, and changed bounds for statistics */
455 int oldnaggrvars = *naggrvars;
456 int oldnfixedvars = *nfixedvars;
457 int oldnchgbds = *nchgbds;
458
459 /* create presolving objects */
460 Problem<SCIP_Real> problem = buildProblem(scip, matrix);
461 int oldnnz = problem.getConstraintMatrix().getNnz();
462 Presolve<SCIP_Real> presolve = setupPresolve(scip, data, allowconsmodification);
463
464 /* call the presolving */
466 " (%.1fs) running MILP presolver\n", SCIPgetSolvingTime(scip));
467
468 /*call presolving without storing information for dual postsolve*/
469#if (PAPILO_VERSION_MAJOR >= 2)
470 PresolveResult<SCIP_Real> res = presolve.apply(problem, false);
471#else
472 PresolveResult<SCIP_Real> res = presolve.apply(problem);
473#endif
474 data->lastncols = problem.getNCols();
475 data->lastnrows = problem.getNRows();
476
477 /* evaluate the result */
478 switch(res.status)
479 {
480 case PresolveStatus::kInfeasible:
481 *result = SCIP_CUTOFF;
483 " (%.1fs) MILP presolver detected infeasibility\n",
485 SCIPmatrixFree(scip, &matrix);
486 return SCIP_OKAY;
487 case PresolveStatus::kUnbndOrInfeas:
488 case PresolveStatus::kUnbounded:
489 *result = SCIP_UNBOUNDED;
491 " (%.1fs) MILP presolver detected unboundedness\n",
493 SCIPmatrixFree(scip, &matrix);
494 return SCIP_OKAY;
495 case PresolveStatus::kUnchanged:
496 *result = SCIP_DIDNOTFIND;
497 data->lastncols = nvars;
498 data->lastnrows = nconss;
500 " (%.1fs) MILP presolver found nothing\n",
502 SCIPmatrixFree(scip, &matrix);
503 return SCIP_OKAY;
504 case PresolveStatus::kReduced:
505 data->lastncols = problem.getNCols();
506 data->lastnrows = problem.getNRows();
507 *result = SCIP_SUCCESS;
508 }
509
510 /* result indicated success, now populate the changes into the SCIP structures */
511
512 /* tighten bounds of variables that are still present after presolving */
513 VariableDomains<SCIP_Real>& varDomains = problem.getVariableDomains();
514 for( int i = 0; i != problem.getNCols(); ++i )
515 {
516 assert( ! varDomains.flags[i].test(ColFlag::kInactive) );
517 SCIP_VAR* var = SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[i]);
518 if( !varDomains.flags[i].test(ColFlag::kLbInf) )
519 {
520 SCIP_Bool infeas;
521 SCIP_Bool tightened;
522 SCIP_CALL( SCIPtightenVarLb(scip, var, varDomains.lower_bounds[i], TRUE, &infeas, &tightened) );
523
524 if( tightened )
525 *nchgbds += 1;
526
527 if( infeas )
528 {
529 *result = SCIP_CUTOFF;
530 break;
531 }
532 }
533
534 if( !varDomains.flags[i].test(ColFlag::kUbInf) )
535 {
536 SCIP_Bool infeas;
537 SCIP_Bool tightened;
538 SCIP_CALL( SCIPtightenVarUb(scip, var, varDomains.upper_bounds[i], TRUE, &infeas, &tightened) );
539
540 if( tightened )
541 *nchgbds += 1;
542
543 if( infeas )
544 {
545 *result = SCIP_CUTOFF;
546 break;
547 }
548 }
549 }
550
551 if( *result == SCIP_CUTOFF )
552 {
554 " (%.1fs) MILP presolver detected infeasibility\n",
556 SCIPmatrixFree(scip, &matrix);
557 return SCIP_OKAY;
558 }
559
560 /* transfer variable fixings and aggregations */
561 std::vector<SCIP_VAR*> tmpvars;
562 std::vector<SCIP_Real> tmpvals;
563
564 /* if the size of the problem decreased by a sufficient factor, create all constraints from scratch if allowed */
565 int newnnz = problem.getConstraintMatrix().getNnz();
566 bool constraintsReplaced = false;
567 if( newnnz == 0 || (allowconsmodification &&
568 (problem.getNCols() <= data->modifyconsfac * SCIPmatrixGetNColumns(matrix) ||
569 problem.getNRows() <= data->modifyconsfac * SCIPmatrixGetNRows(matrix) ||
570 newnnz <= data->modifyconsfac * oldnnz)) )
571 {
572 int oldnrows = SCIPmatrixGetNRows(matrix);
573 int newnrows = problem.getNRows();
574
575 constraintsReplaced = true;
576
577 /* capture constraints that are still present in the problem after presolve */
578 for( int i = 0; i < newnrows; ++i )
579 {
580 SCIP_CONS* c = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
582 }
583
584 /* delete all constraints */
585 *ndelconss += oldnrows;
586 *naddconss += newnrows;
587
588 for( int i = 0; i < oldnrows; ++i )
589 {
591 }
592
593 /* now loop over rows of presolved problem and create them as new linear constraints,
594 * then release the old constraint after its name was passed to the new constraint */
595 const Vec<RowFlags>& rflags = problem.getRowFlags();
596 const auto& consmatrix = problem.getConstraintMatrix();
597 for( int i = 0; i < newnrows; ++i )
598 {
599 auto rowvec = consmatrix.getRowCoefficients(i);
600 const int* rowcols = rowvec.getIndices();
601 /* SCIPcreateConsBasicLinear() requires a non const pointer */
602 SCIP_Real* rowvals = const_cast<SCIP_Real*>(rowvec.getValues());
603 int rowlen = rowvec.getLength();
604
605 /* retrieve SCIP compatible left and right hand sides */
606 SCIP_Real lhs = rflags[i].test(RowFlag::kLhsInf) ? - SCIPinfinity(scip) : consmatrix.getLeftHandSides()[i];
607 SCIP_Real rhs = rflags[i].test(RowFlag::kRhsInf) ? SCIPinfinity(scip) : consmatrix.getRightHandSides()[i];
608
609 /* create variable array matching the value array */
610 tmpvars.clear();
611 tmpvars.reserve(rowlen);
612 for( int j = 0; j < rowlen; ++j )
613 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[rowcols[j]]));
614
615 /* create and add new constraint with name of old constraint */
616 SCIP_CONS* oldcons = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
617 SCIP_CONS* cons;
618 SCIP_CALL( SCIPcreateConsBasicLinear(scip, &cons, SCIPconsGetName(oldcons), rowlen, tmpvars.data(), rowvals, lhs, rhs) );
619 SCIP_CALL( SCIPaddCons(scip, cons) );
620
621 /* release old and new constraint */
622 SCIP_CALL( SCIPreleaseCons(scip, &oldcons) );
623 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
624 }
625 }
626
627 /* PaPILO's aggregations are valid regarding the constraints as they were presolved by PaPILO.
628 * If coefficients were changed, but constraints in SCIP are not replaced by those from PaPILO,
629 * then it can not be guaranteed that the bounds of multiaggregated variables will be enforced,
630 * i.e., will be implied by the constraints in SCIP (see also #3704).
631 * Only for variable aggregations, SCIP will ensure this by tightening the bounds on the aggregation
632 * variable as part of SCIPaggregateVars(). For multiaggregations, we will only accept those
633 * where we can be sure with a simple check that the bounds on the aggregated variable are implied.
634 */
635 bool checkmultaggr =
636#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 3)
637 presolve.getStatistics().single_matrix_coefficient_changes > 0
638#else
639 presolve.getStatistics().ncoefchgs > 0
640#endif
641 && !constraintsReplaced;
642
643 /* loop over res.postsolve and add all fixed variables and aggregations to scip */
644 for( std::size_t i = 0; i != res.postsolve.types.size(); ++i )
645 {
646 ReductionType type = res.postsolve.types[i];
647 int first = res.postsolve.start[i];
648 int last = res.postsolve.start[i + 1];
649
650 switch( type )
651 {
652 case ReductionType::kFixedCol:
653 {
654 SCIP_Bool infeas;
655 SCIP_Bool fixed;
656 int col = res.postsolve.indices[first];
657
658 SCIP_VAR* var = SCIPmatrixGetVar(matrix, col);
659
660 SCIP_Real value = res.postsolve.values[first];
661
662 SCIP_CALL( SCIPfixVar(scip, var, value, &infeas, &fixed) );
663 *nfixedvars += 1;
664
665 assert(!infeas);
666 /* SCIP has different rules for aggregating variables than PaPILO: therefore the variable PaPILO
667 * tries to fix now may have been aggregated by SCIP before. Additionally, after aggregation SCIP
668 * sometimes performs bound tightening resulting in possible fixings. These cases need to be excluded. */
669 assert(fixed || SCIPvarGetStatus(var) == SCIP_VARSTATUS_AGGREGATED ||
671 break;
672 }
673/*
674 * Dual-postsolving in PaPILO required introducing a postsolve-type for substitution with additional information.
675 * Further, the different Substitution-postsolving types store the required postsolving data differently (in different order) in the postsolving stack.
676 * Therefore, we need to distinguish how to parse the required data (rowLength, col, side, startRowCoefficients, lastRowCoefficients) from the postsolving stack.
677 * If these values are accessed, the procedure is the same for both.
678 */
679#if (PAPILO_VERSION_MAJOR >= 2)
680 case ReductionType::kSubstitutedColWithDual:
681#endif
682 case ReductionType::kSubstitutedCol:
683 {
684 int col = 0;
685 SCIP_Real side = 0;
686
687 int rowlen = 0;
688 int startRowCoefficients = 0;
689 int lastRowCoefficients = 0;
690
691 if( type == ReductionType::kSubstitutedCol )
692 {
693 rowlen = last - first - 1;
694 col = res.postsolve.indices[first];
695 side = res.postsolve.values[first];
696
697 startRowCoefficients = first + 1;
698 lastRowCoefficients = last;
699 }
700#if (PAPILO_VERSION_MAJOR >= 2)
701 if( type == ReductionType::kSubstitutedColWithDual )
702 {
703 rowlen = (int) res.postsolve.values[first];
704 col = res.postsolve.indices[first + 3 + rowlen];
705 side = res.postsolve.values[first + 1];
706
707 startRowCoefficients = first + 3;
708 lastRowCoefficients = first + 3 + rowlen;
709
710 assert(side == res.postsolve.values[first + 2]);
711 assert(res.postsolve.indices[first + 1] == 0);
712 assert(res.postsolve.indices[first + 2] == 0);
713 }
714 assert( type == ReductionType::kSubstitutedCol || type == ReductionType::kSubstitutedColWithDual );
715#else
716 assert( type == ReductionType::kSubstitutedCol );
717#endif
718 SCIP_Bool infeas;
719 SCIP_Bool aggregated;
720 SCIP_Bool redundant = FALSE;
721 SCIP_Real constant = 0.0;
722 if( rowlen == 2 )
723 {
724 SCIP_Real updatedSide;
725 SCIP_VAR* varx = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients]);
726 SCIP_VAR* vary = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients + 1]);
727 SCIP_Real scalarx = res.postsolve.values[startRowCoefficients];
728 SCIP_Real scalary = res.postsolve.values[startRowCoefficients + 1];
729
730 SCIP_CALL( SCIPgetProbvarSum(scip, &varx, &scalarx, &constant) );
732
733 SCIP_CALL( SCIPgetProbvarSum(scip, &vary, &scalary, &constant) );
735
736 updatedSide = side - constant;
737
738 SCIP_CALL( SCIPaggregateVars(scip, varx, vary, scalarx, scalary, updatedSide, &infeas, &redundant, &aggregated) );
739 }
740 else
741 {
742 SCIP_Real colCoef = 0.0;
743 SCIP_Real updatedSide;
744 SCIP_Bool checklbimplied;
745 SCIP_Bool checkubimplied;
746 SCIP_Real impliedlb;
747 SCIP_Real impliedub;
748 int j;
749
750 for( j = startRowCoefficients; j < lastRowCoefficients; ++j )
751 {
752 if( res.postsolve.indices[j] == col )
753 {
754 colCoef = res.postsolve.values[j];
755 break;
756 }
757 }
758
759 tmpvars.clear();
760 tmpvals.clear();
761 tmpvars.reserve(rowlen);
762 tmpvals.reserve(rowlen);
763
764 assert(colCoef != 0.0);
765 SCIP_VAR* aggrvar = SCIPmatrixGetVar(matrix, col);
766
767 SCIP_CALL( SCIPgetProbvarSum(scip, &aggrvar, &colCoef, &constant) );
768 assert(SCIPvarGetStatus(aggrvar) != SCIP_VARSTATUS_MULTAGGR);
769
770 updatedSide = side - constant;
771
772 /* we need to check whether lb/ub on aggrvar is implied by bounds of other variables in multiaggregation
773 * if checkmultaggr is TRUE and the lb/ub is finite
774 * it should be sufficient to ensure global bounds on aggrvar (and as we are in presolve, local=global anyway)
775 */
776 checklbimplied = checkmultaggr && !SCIPisInfinity(scip, -SCIPvarGetLbGlobal(aggrvar));
777 checkubimplied = checkmultaggr && !SCIPisInfinity(scip, SCIPvarGetUbGlobal(aggrvar));
778 impliedlb = impliedub = updatedSide / colCoef;
779
780 for( j = startRowCoefficients; j < lastRowCoefficients; ++j )
781 {
782 SCIP_Real coef;
783 SCIP_VAR* var;
784
785 if( res.postsolve.indices[j] == col )
786 continue;
787
788 coef = - res.postsolve.values[j] / colCoef;
789 var = SCIPmatrixGetVar(matrix, res.postsolve.indices[j]);
790
791 if( checklbimplied )
792 {
793 if( coef > 0.0 )
794 {
795 /* if impliedlb will be -infinity, then we can give up: we cannot use this mutiaggregation */
797 break;
798 else
799 impliedlb += coef * SCIPvarGetLbLocal(var);
800 }
801 else
802 {
804 break;
805 else
806 impliedlb += coef * SCIPvarGetUbLocal(var);
807 }
808 }
809
810 if( checkubimplied )
811 {
812 if( coef > 0.0 )
813 {
815 break;
816 else
817 impliedub += coef * SCIPvarGetUbLocal(var);
818 }
819 else
820 {
822 break;
823 else
824 impliedub += coef * SCIPvarGetLbLocal(var);
825 }
826 }
827
828 tmpvals.push_back(coef);
829 tmpvars.push_back(var);
830 }
831
832 /* if implied bounds are not sufficient to ensure bounds on aggrvar, then we cannot use the multiaggregation */
833 if( j < lastRowCoefficients )
834 break;
835
836 if( checklbimplied && SCIPisGT(scip, SCIPvarGetLbGlobal(aggrvar), impliedlb) )
837 break;
838
839 if( checkubimplied && SCIPisLT(scip, SCIPvarGetUbGlobal(aggrvar), impliedub) )
840 break;
841
842 SCIP_CALL( SCIPmultiaggregateVar(scip, aggrvar, tmpvars.size(),
843 tmpvars.data(), tmpvals.data(), updatedSide / colCoef, &infeas, &aggregated) );
844 }
845
846 if( aggregated )
847 *naggrvars += 1;
848 else if( constraintsReplaced && !redundant )
849 {
850 /* if the constraints where replaced, we need to add the failed substitution as an equality to SCIP */
851 tmpvars.clear();
852 tmpvals.clear();
853 for( int j = startRowCoefficients; j < lastRowCoefficients; ++j )
854 {
855 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.indices[j]));
856 tmpvals.push_back(res.postsolve.values[j]);
857 }
858
859 SCIP_CONS* cons;
860 String name = fmt::format("{}_failed_aggregation_equality", SCIPvarGetName(SCIPmatrixGetVar(matrix, col)));
861 SCIP_CALL( SCIPcreateConsBasicLinear(scip, &cons, name.c_str(),
862 tmpvars.size(), tmpvars.data(), tmpvals.data(), side, side ) );
863 SCIP_CALL( SCIPaddCons(scip, cons) );
864 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
865 *naddconss += 1;
866 }
867
868 if( infeas )
869 {
870 *result = SCIP_CUTOFF;
871 break;
872 }
873
874 break;
875 }
876 case ReductionType::kParallelCol:
877 return SCIP_INVALIDRESULT;
878#if (PAPILO_VERSION_MAJOR <= 1 && PAPILO_VERSION_MINOR==0)
879#else
880 case ReductionType::kFixedInfCol: {
881 /* todo: currently SCIP can not handle this kind of reduction (see issue #3391) */
882 assert(false);
883 if(!constraintsReplaced)
884 continue;
885 SCIP_Bool infeas;
886 SCIP_Bool fixed;
887 SCIP_Real value = SCIPinfinity(scip);
888
889 int column = res.postsolve.indices[first];
890 bool is_negative_infinity = res.postsolve.values[first] < 0;
891 SCIP_VAR* column_variable = SCIPmatrixGetVar(matrix, column);
892
893 if( is_negative_infinity )
894 {
895 value = -SCIPinfinity(scip);
896 }
897
898 SCIP_CALL( SCIPfixVar(scip, column_variable, value, &infeas, &fixed) );
899 *nfixedvars += 1;
900
901 assert(!infeas);
902 assert(fixed);
903 break;
904 }
905#endif
906#if (PAPILO_VERSION_MAJOR >= 2)
907 case ReductionType::kVarBoundChange :
908 case ReductionType::kRedundantRow :
909 case ReductionType::kRowBoundChange :
910 case ReductionType::kReasonForRowBoundChangeForcedByRow :
911 case ReductionType::kRowBoundChangeForcedByRow :
912 case ReductionType::kSaveRow :
913 case ReductionType::kReducedBoundsCost :
914 case ReductionType::kColumnDualValue :
915 case ReductionType::kRowDualValue :
916 case ReductionType::kCoefficientChange :
917 // dual ReductionTypes should be only calculated for dual reductions and should not appear for MIP
918 SCIPerrorMessage("PaPILO: PaPILO should not return dual postsolving reductions in SCIP!!\n");
919 SCIPABORT(); /*lint --e{527}*/
920 break;
921#endif
922 default:
923 SCIPdebugMsg(scip, "PaPILO returned unknown data type: \n" );
924 continue;
925 }
926 }
927
928 /* finish with a final verb message and return */
930 " (%.1fs) MILP presolver (%d rounds): %d aggregations, %d fixings, %d bound changes\n",
931 SCIPgetSolvingTime(scip), presolve.getStatistics().nrounds, *naggrvars - oldnaggrvars,
932 *nfixedvars - oldnfixedvars, *nchgbds - oldnchgbds);
933
934 /* free the matrix */
935 assert(initialized);
936 SCIPmatrixFree(scip, &matrix);
937
938 return SCIP_OKAY;
939}
940
941
942/*
943 * presolver specific interface methods
944 */
945
946/** creates the xyz presolver and includes it in SCIP */
948 SCIP* scip /**< SCIP data structure */
949 )
950{
951 SCIP_PRESOLDATA* presoldata;
952 SCIP_PRESOL* presol;
953
954#if defined(PAPILO_VERSION_TWEAK) && PAPILO_VERSION_TWEAK != 0
955 String name = fmt::format("PaPILO {}.{}.{}.{}", PAPILO_VERSION_MAJOR, PAPILO_VERSION_MINOR, PAPILO_VERSION_PATCH, PAPILO_VERSION_TWEAK);
956#else
957 String name = fmt::format("PaPILO {}.{}.{}", PAPILO_VERSION_MAJOR, PAPILO_VERSION_MINOR, PAPILO_VERSION_PATCH);
958#endif
959
960#if defined(PAPILO_GITHASH_AVAILABLE) && defined(PAPILO_TBB)
961 String desc = fmt::format("parallel presolve for integer and linear optimization (github.com/scipopt/papilo) (built with TBB) [GitHash: {}]", PAPILO_GITHASH);
962#elif !defined(PAPILO_GITHASH_AVAILABLE) && !defined(PAPILO_TBB)
963 String desc("parallel presolve for integer and linear optimization (github.com/scipopt/papilo)");
964#elif defined(PAPILO_GITHASH_AVAILABLE) && !defined(PAPILO_TBB)
965 String desc = fmt::format("parallel presolve for integer and linear optimization (github.com/scipopt/papilo) [GitHash: {}]", PAPILO_GITHASH);
966#elif !defined(PAPILO_GITHASH_AVAILABLE) && defined(PAPILO_TBB)
967 String desc = fmt::format("parallel presolve for integer and linear optimization (github.com/scipopt/papilo) (built with TBB)");
968#endif
969
970 /* add external code info for the presolve library */
971 SCIP_CALL( SCIPincludeExternalCodeInformation(scip, name.c_str(), desc.c_str()) );
972
973 /* create MILP presolver data */
974 presoldata = NULL;
975 SCIP_CALL( SCIPallocBlockMemory(scip, &presoldata) );
976 BMSclearMemory(presoldata);
977
978 presol = NULL;
979
980 /* include presolver */
982 presolExecMILP,
983 presoldata) );
984
985 assert(presol != NULL);
986
987 /* set non fundamental callbacks via setter functions */
988 SCIP_CALL( SCIPsetPresolCopy(scip, presol, presolCopyMILP) );
989 SCIP_CALL( SCIPsetPresolFree(scip, presol, presolFreeMILP) );
990 SCIP_CALL( SCIPsetPresolInit(scip, presol, presolInitMILP) );
991
992 /* add MILP presolver parameters */
994 "presolving/" PRESOL_NAME "/threads",
995 "maximum number of threads presolving may use (0: automatic)",
996 &presoldata->threads, FALSE, DEFAULT_THREADS, 0, INT_MAX, NULL, NULL) );
997
999 "presolving/" PRESOL_NAME "/maxfillinpersubstitution",
1000 "maximal possible fillin for substitutions to be considered",
1001 &presoldata->maxfillinpersubstitution, FALSE, DEFAULT_MAXFILLINPERSUBST, INT_MIN, INT_MAX, NULL, NULL) );
1002
1004 "presolving/" PRESOL_NAME "/maxshiftperrow",
1005 "maximal amount of nonzeros allowed to be shifted to make space for substitutions",
1006 &presoldata->maxshiftperrow, TRUE, DEFAULT_MAXSHIFTPERROW, 0, INT_MAX, NULL, NULL) );
1007
1009 "presolving/" PRESOL_NAME "/randomseed",
1010 "the random seed used for randomization of tie breaking",
1011 &presoldata->randomseed, FALSE, DEFAULT_RANDOMSEED, INT_MIN, INT_MAX, NULL, NULL) );
1012
1013 if( DependentRows<double>::Enabled )
1014 {
1016 "presolving/" PRESOL_NAME "/detectlineardependency",
1017 "should linear dependent equations and free columns be removed? (0: never, 1: for LPs, 2: always)",
1018 &presoldata->detectlineardependency, TRUE, DEFAULT_DETECTLINDEP, 0, 2, NULL, NULL) );
1019 }
1020 else
1021 presoldata->detectlineardependency = 0;
1022
1024 "presolving/" PRESOL_NAME "/modifyconsfac",
1025 "modify SCIP constraints when the number of nonzeros or rows is at most this factor "
1026 "times the number of nonzeros or rows before presolving",
1027 &presoldata->modifyconsfac, FALSE, DEFAULT_MODIFYCONSFAC, 0.0, 1.0, NULL, NULL) );
1028
1030 "presolving/" PRESOL_NAME "/markowitztolerance",
1031 "the markowitz tolerance used for substitutions",
1032 &presoldata->markowitztolerance, FALSE, DEFAULT_MARKOWITZTOLERANCE, 0.0, 1.0, NULL, NULL) );
1033
1035 "presolving/" PRESOL_NAME "/hugebound",
1036 "absolute bound value that is considered too huge for activitity based calculations",
1037 &presoldata->hugebound, FALSE, DEFAULT_HUGEBOUND, 0.0, SCIP_REAL_MAX, NULL, NULL) );
1038
1039#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 1)
1040 SCIP_CALL(SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/maxbadgesizeseq",
1041 "maximal badge size in Probing in PaPILO if PaPILO is executed in sequential mode",
1042 &presoldata->maxbadgesizeseq, FALSE, DEFAULT_MAXBADGESIZE_SEQ, -1, INT_MAX, NULL, NULL));
1043
1044 SCIP_CALL(SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/maxbadgesizepar",
1045 "maximal badge size in Probing in PaPILO if PaPILO is executed in parallel mode",
1046 &presoldata->maxbadgesizepar, FALSE, DEFAULT_MAXBADGESIZE_PAR, -1, INT_MAX, NULL, NULL));
1047#else
1048 presoldata->maxbadgesizeseq = DEFAULT_MAXBADGESIZE_SEQ;
1049 presoldata->maxbadgesizepar = DEFAULT_MAXBADGESIZE_PAR;
1050#endif
1051
1052#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 3)
1053 SCIP_CALL(SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/internalmaxrounds",
1054 "internal maxrounds for each milp presolving (-1: no limit, 0: model cleanup)",
1055 &presoldata->internalmaxrounds, TRUE, DEFAULT_INTERNAL_MAXROUNDS, -1, INT_MAX, NULL, NULL));
1056#else
1057 presoldata->internalmaxrounds = DEFAULT_INTERNAL_MAXROUNDS;
1058#endif
1059
1061 "presolving/" PRESOL_NAME "/enableparallelrows",
1062 "should the parallel rows presolver be enabled within the presolve library?",
1063 &presoldata->enableparallelrows, TRUE, DEFAULT_ENABLEPARALLELROWS, NULL, NULL) );
1064
1066 "presolving/" PRESOL_NAME "/enabledomcol",
1067 "should the dominated column presolver be enabled within the presolve library?",
1068 &presoldata->enabledomcol, TRUE, DEFAULT_ENABLEDOMCOL, NULL, NULL) );
1069
1071 "presolving/" PRESOL_NAME "/enabledualinfer",
1072 "should the dualinfer presolver be enabled within the presolve library?",
1073 &presoldata->enabledualinfer, TRUE, DEFAULT_ENABLEDUALINFER, NULL, NULL) );
1074
1076 "presolving/" PRESOL_NAME "/enablemultiaggr",
1077 "should the multi-aggregation presolver be enabled within the presolve library?",
1078 &presoldata->enablemultiaggr, TRUE, DEFAULT_ENABLEMULTIAGGR, NULL, NULL) );
1079
1081 "presolving/" PRESOL_NAME "/enableprobing",
1082 "should the probing presolver be enabled within the presolve library?",
1083 &presoldata->enableprobing, TRUE, DEFAULT_ENABLEPROBING, NULL, NULL) );
1084
1086 "presolving/" PRESOL_NAME "/enablesparsify",
1087 "should the sparsify presolver be enabled within the presolve library?",
1088 &presoldata->enablesparsify, TRUE, DEFAULT_ENABLESPARSIFY, NULL, NULL) );
1089
1090 SCIP_CALL( SCIPaddStringParam(scip, "presolving/" PRESOL_NAME "/probfilename",
1091 "filename to store the problem before MILP presolving starts (only enforced constraints)",
1092 &presoldata->filename, TRUE, DEFAULT_FILENAME_PROBLEM, NULL, NULL) );
1093
1094 SCIP_CALL(SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/verbosity",
1095 "verbosity level of PaPILO (0: quiet, 1: errors, 2: warnings, 3: normal, 4: detailed)",
1096 &presoldata->verbosity, FALSE, DEFAULT_VERBOSITY, 0, 4, NULL, NULL));
1097
1098 return SCIP_OKAY;
1099}
1100
1101#endif
Constraint handler for linear constraints in their most general form, .
#define NULL
Definition: def.h:267
#define SCIP_REAL_MAX
Definition: def.h:174
#define SCIP_Bool
Definition: def.h:91
#define SCIP_Real
Definition: def.h:173
#define TRUE
Definition: def.h:93
#define FALSE
Definition: def.h:94
#define SCIPABORT()
Definition: def.h:346
#define SCIP_CALL(x)
Definition: def.h:374
SCIP_RETCODE SCIPcreateConsBasicLinear(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_Real *vals, SCIP_Real lhs, SCIP_Real rhs)
const char * SCIPgetProbName(SCIP *scip)
Definition: scip_prob.c:1067
int SCIPgetNVars(SCIP *scip)
Definition: scip_prob.c:1992
SCIP_RETCODE SCIPaddCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip_prob.c:2770
SCIP_RETCODE SCIPdelCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip_prob.c:2843
int SCIPgetNConss(SCIP *scip)
Definition: scip_prob.c:3042
SCIP_RETCODE SCIPwriteTransProblem(SCIP *scip, const char *filename, const char *extension, SCIP_Bool genericnames)
Definition: scip_prob.c:648
void SCIPverbMessage(SCIP *scip, SCIP_VERBLEVEL msgverblevel, FILE *file, const char *formatstr,...)
Definition: scip_message.c:225
#define SCIPdebugMsg
Definition: scip_message.h:78
SCIP_RETCODE SCIPaddIntParam(SCIP *scip, const char *name, const char *desc, int *valueptr, SCIP_Bool isadvanced, int defaultvalue, int minvalue, int maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:83
SCIP_RETCODE SCIPaddStringParam(SCIP *scip, const char *name, const char *desc, char **valueptr, SCIP_Bool isadvanced, const char *defaultvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:194
SCIP_RETCODE SCIPaddRealParam(SCIP *scip, const char *name, const char *desc, SCIP_Real *valueptr, SCIP_Bool isadvanced, SCIP_Real defaultvalue, SCIP_Real minvalue, SCIP_Real maxvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:139
SCIP_RETCODE SCIPgetRealParam(SCIP *scip, const char *name, SCIP_Real *value)
Definition: scip_param.c:307
SCIP_RETCODE SCIPaddBoolParam(SCIP *scip, const char *name, const char *desc, SCIP_Bool *valueptr, SCIP_Bool isadvanced, SCIP_Bool defaultvalue, SCIP_DECL_PARAMCHGD((*paramchgd)), SCIP_PARAMDATA *paramdata)
Definition: scip_param.c:57
SCIP_RETCODE SCIPincludePresolMILP(SCIP *scip)
Definition: presol_milp.cpp:46
int SCIPconshdlrGetNCheckConss(SCIP_CONSHDLR *conshdlr)
Definition: cons.c:4656
SCIP_CONSHDLR * SCIPfindConshdlr(SCIP *scip, const char *name)
Definition: scip_cons.c:941
const char * SCIPconsGetName(SCIP_CONS *cons)
Definition: cons.c:8214
SCIP_RETCODE SCIPreleaseCons(SCIP *scip, SCIP_CONS **cons)
Definition: scip_cons.c:1174
SCIP_RETCODE SCIPcaptureCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip_cons.c:1139
SCIP_RETCODE SCIPincludeExternalCodeInformation(SCIP *scip, const char *name, const char *description)
Definition: scip_general.c:734
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip_mem.h:108
#define SCIPallocBlockMemory(scip, ptr)
Definition: scip_mem.h:89
void SCIPpresolSetData(SCIP_PRESOL *presol, SCIP_PRESOLDATA *presoldata)
Definition: presol.c:522
SCIP_PRESOLDATA * SCIPpresolGetData(SCIP_PRESOL *presol)
Definition: presol.c:512
SCIP_RETCODE SCIPsetPresolInit(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLINIT((*presolinit)))
Definition: scip_presol.c:172
SCIP_RETCODE SCIPsetPresolFree(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLFREE((*presolfree)))
Definition: scip_presol.c:156
SCIP_RETCODE SCIPsetPresolCopy(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLCOPY((*presolcopy)))
Definition: scip_presol.c:140
SCIP_RETCODE SCIPincludePresolBasic(SCIP *scip, SCIP_PRESOL **presolptr, const char *name, const char *desc, int priority, int maxrounds, SCIP_PRESOLTIMING timing, SCIP_DECL_PRESOLEXEC((*presolexec)), SCIP_PRESOLDATA *presoldata)
Definition: scip_presol.c:105
SCIP_Real SCIPgetSolvingTime(SCIP *scip)
Definition: scip_timing.c:378
SCIP_Real SCIPinfinity(SCIP *scip)
SCIP_Bool SCIPisInfinity(SCIP *scip, SCIP_Real val)
SCIP_Real SCIPfeastol(SCIP *scip)
SCIP_Bool SCIPisGT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Real SCIPepsilon(SCIP *scip)
SCIP_Bool SCIPisLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_RETCODE SCIPtightenVarLb(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound, SCIP_Bool force, SCIP_Bool *infeasible, SCIP_Bool *tightened)
Definition: scip_var.c:5203
SCIP_VARSTATUS SCIPvarGetStatus(SCIP_VAR *var)
Definition: var.c:17538
SCIP_Real SCIPvarGetUbLocal(SCIP_VAR *var)
Definition: var.c:18144
SCIP_RETCODE SCIPaggregateVars(SCIP *scip, SCIP_VAR *varx, SCIP_VAR *vary, SCIP_Real scalarx, SCIP_Real scalary, SCIP_Real rhs, SCIP_Bool *infeasible, SCIP_Bool *redundant, SCIP_Bool *aggregated)
Definition: scip_var.c:8401
SCIP_Real SCIPvarGetObj(SCIP_VAR *var)
Definition: var.c:17926
SCIP_RETCODE SCIPtightenVarUb(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound, SCIP_Bool force, SCIP_Bool *infeasible, SCIP_Bool *tightened)
Definition: scip_var.c:5320
SCIP_VARTYPE SCIPvarGetType(SCIP_VAR *var)
Definition: var.c:17584
SCIP_RETCODE SCIPgetProbvarSum(SCIP *scip, SCIP_VAR **var, SCIP_Real *scalar, SCIP_Real *constant)
Definition: scip_var.c:1794
SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:18088
const char * SCIPvarGetName(SCIP_VAR *var)
Definition: var.c:17419
SCIP_RETCODE SCIPmultiaggregateVar(SCIP *scip, SCIP_VAR *var, int naggvars, SCIP_VAR **aggvars, SCIP_Real *scalars, SCIP_Real constant, SCIP_Bool *infeasible, SCIP_Bool *aggregated)
Definition: scip_var.c:8535
SCIP_Bool SCIPvarIsIntegral(SCIP_VAR *var)
Definition: var.c:17610
SCIP_Real SCIPvarGetLbLocal(SCIP_VAR *var)
Definition: var.c:18134
SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:18078
SCIP_RETCODE SCIPfixVar(SCIP *scip, SCIP_VAR *var, SCIP_Real fixedval, SCIP_Bool *infeasible, SCIP_Bool *fixed)
Definition: scip_var.c:8276
SCIP_Bool SCIPallowWeakDualReds(SCIP *scip)
Definition: scip_var.c:8656
SCIP_Bool SCIPallowStrongDualReds(SCIP *scip)
Definition: scip_var.c:8629
unsigned int SCIPinitializeRandomSeed(SCIP *scip, unsigned int initialseedvalue)
int SCIPmatrixGetNNonzs(SCIP_MATRIX *matrix)
Definition: matrix.c:1778
int SCIPmatrixGetRowNNonzs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1708
SCIP_Real SCIPmatrixGetRowLhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1742
SCIP_Real * SCIPmatrixGetRowValPtr(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1684
SCIP_Real SCIPmatrixGetRowRhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1754
SCIP_RETCODE SCIPmatrixCreate(SCIP *scip, SCIP_MATRIX **matrixptr, SCIP_Bool onlyifcomplete, SCIP_Bool *initialized, SCIP_Bool *complete, SCIP_Bool *infeasible, int *naddconss, int *ndelconss, int *nchgcoefs, int *nchgbds, int *nfixedvars)
Definition: matrix.c:454
int SCIPmatrixGetNColumns(SCIP_MATRIX *matrix)
Definition: matrix.c:1604
SCIP_CONS * SCIPmatrixGetCons(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1860
void SCIPmatrixFree(SCIP *scip, SCIP_MATRIX **matrix)
Definition: matrix.c:1072
SCIP_VAR * SCIPmatrixGetVar(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1660
int * SCIPmatrixGetRowIdxPtr(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1696
int SCIPmatrixGetNRows(SCIP_MATRIX *matrix)
Definition: matrix.c:1732
#define BMSclearMemory(ptr)
Definition: memory.h:129
#define PRESOL_NAME
#define PRESOL_PRIORITY
#define PRESOL_MAXROUNDS
#define PRESOL_TIMING
#define PRESOL_DESC
MILP presolver that calls the presolve library on the constraint matrix.
public methods for managing constraints
public methods for matrix
public methods for message output
#define SCIPerrorMessage
Definition: pub_message.h:64
public methods for presolvers
public methods for problem variables
#define DEFAULT_RANDOMSEED
Definition: reader_scflp.c:116
public methods for constraint handler plugins and constraints
general public methods
public methods for memory management
public methods for message handling
public methods for numerical tolerances
public methods for SCIP parameter handling
public methods for presolving plugins
public methods for global and local (sub)problems
public methods for random numbers
static SCIP_RETCODE presolve(SCIP *scip, SCIP_Bool *unbounded, SCIP_Bool *infeasible, SCIP_Bool *vanished)
Definition: scip_solve.c:1110
public methods for timing
public methods for SCIP variables
@ SCIP_VERBLEVEL_HIGH
Definition: type_message.h:56
#define SCIP_DECL_PRESOLCOPY(x)
Definition: type_presol.h:60
struct SCIP_PresolData SCIP_PRESOLDATA
Definition: type_presol.h:51
#define SCIP_DECL_PRESOLFREE(x)
Definition: type_presol.h:68
#define SCIP_DECL_PRESOLINIT(x)
Definition: type_presol.h:76
#define SCIP_DECL_PRESOLEXEC(x)
Definition: type_presol.h:167
@ SCIP_DIDNOTRUN
Definition: type_result.h:42
@ SCIP_CUTOFF
Definition: type_result.h:48
@ SCIP_DIDNOTFIND
Definition: type_result.h:44
@ SCIP_UNBOUNDED
Definition: type_result.h:47
@ SCIP_SUCCESS
Definition: type_result.h:58
@ SCIP_INVALIDRESULT
Definition: type_retcode.h:53
@ SCIP_OKAY
Definition: type_retcode.h:42
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:63
@ SCIP_VARTYPE_IMPLINT
Definition: type_var.h:64
@ SCIP_VARSTATUS_FIXED
Definition: type_var.h:52
@ SCIP_VARSTATUS_MULTAGGR
Definition: type_var.h:54
@ SCIP_VARSTATUS_NEGATED
Definition: type_var.h:55
@ SCIP_VARSTATUS_AGGREGATED
Definition: type_var.h:53