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-2025 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 * @author Alexander Hoen
29 *
30 * Calls the presolve library and communicates (multi-)aggregations, fixings, and bound
31 * changes to SCIP by utilizing the postsolve information. Constraint changes can currently
32 * only be communicated by deleting all constraints and adding new ones.
33 *
34 * @todo add infrastructure to SCIP for handling parallel columns
35 * @todo better communication of constraint changes by adding more information to the postsolve structure
36 * @todo allow to pass additional external locks to the presolve library that are considered when doing reductions
37 *
38 */
39
40/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
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"
71#include "scip/pub_matrix.h"
72#include "scip/pub_presol.h"
73#include "scip/pub_var.h"
74#include "scip/pub_cons.h"
75#include "scip/pub_message.h"
76#include "scip/scip_exact.h"
77#include "scip/scip_general.h"
78#include "scip/scip_presol.h"
79#include "scip/scip_var.h"
80#include "scip/scip_mem.h"
81#include "scip/scip_prob.h"
82#include "scip/scip_param.h"
83#include "scip/scip_cons.h"
84#include "scip/scip_numerics.h"
85#include "scip/scip_timing.h"
86#include "scip/scip_message.h"
88#if defined(SCIP_WITH_EXACTSOLVE)
90#endif
91#include "scip/rational.h"
92#include "papilo/core/Presolve.hpp"
93#include "papilo/core/ProblemBuilder.hpp"
94#include "papilo/Config.hpp"
95
96/* API since PaPILO 2.3.0 */
97#if !defined(PAPILO_API_VERSION)
98#define PAPILO_APIVERSION 0
99#elif !(PAPILO_API_VERSION + 0)
100#define PAPILO_APIVERSION 1
101#else
102#define PAPILO_APIVERSION PAPILO_API_VERSION
103#endif
104
105#if defined(SCIP_WITH_GMP) && defined(SCIP_WITH_EXACTSOLVE) && !defined(PAPILO_HAVE_GMP)
106#warning SCIP built with GMP and exact solving, but PaPILO without GMP disables exact presolving.
107#endif
108
109#if defined(SCIP_WITH_GMP) && defined(SCIP_WITH_EXACTSOLVE) && defined(PAPILO_HAVE_GMP)
110#define PAPILO_WITH_EXACTPRESOLVE
111#endif
112
113#define PRESOL_NAME "milp"
114#define PRESOL_DESC "MILP specific presolving methods"
115#define PRESOL_PRIORITY 9999999 /**< priority of the presolver (>= 0: before, < 0: after constraint handlers); combined with propagators */
116#define PRESOL_MAXROUNDS (-1) /**< maximal number of presolving rounds the presolver participates in (-1: no limit) */
117#define PRESOL_TIMING SCIP_PRESOLTIMING_MEDIUM /* timing of the presolver (fast, medium, or exhaustive) */
118
119/** general settings for PaPILO */
120#define DEFAULT_THREADS 1 /**< maximum number of threads presolving may use (0: automatic) */
121#define DEFAULT_ABORTFAC_EXHAUSTIVE 0.0008 /**< the abort factor for exhaustive presolving in PAPILO */
122#define DEFAULT_ABORTFAC_MEDIUM 0.0008 /**< the abort factor for medium presolving in PAPILO */
123#define DEFAULT_ABORTFAC_FAST 0.0008 /**< the abort factor for fast presolving in PAPILO */
124#define DEFAULT_DETECTLINDEP 0 /**< should linear dependent equations and free columns be removed? (0: never, 1: for LPs, 2: always) */
125#define DEFAULT_INTERNAL_MAXROUNDS (-1) /**< internal max rounds in PaPILO (-1: no limit, 0: model cleanup) */
126#define DEFAULT_MODIFYCONSFAC 0.8 /**< modify SCIP constraints when the number of nonzeros or rows is at most this
127 * factor times the number of nonzeros or rows before presolving */
128#define DEFAULT_RANDOMSEED 0 /**< the random seed used for randomization of tie breaking */
129
130/** numerics in PaPILO */
131#define DEFAULT_HUGEBOUND 1e8 /**< absolute bound value that is considered too huge for activitity based calculations */
132
133/** presolvers in PaPILO */
134#define DEFAULT_ENABLEDOMCOL TRUE /**< should the dominated column presolver be enabled within the presolve library? */
135#define DEFAULT_ENABLEDUALINFER TRUE /**< should the dualinfer presolver be enabled within the presolve library? */
136#define DEFAULT_ENABLEMULTIAGGR TRUE /**< should the multi-aggregation/substitution presolver be enabled within the presolve library? */
137#define DEFAULT_ENABLEPARALLELROWS TRUE /**< should the parallel rows presolver be enabled within the presolve library? */
138#define DEFAULT_ENABLEPROBING TRUE /**< should the probing presolver be enabled within the presolve library? */
139#define DEFAULT_ENABLESPARSIFY FALSE /**< should the sparsify presolver be enabled within the presolve library? */
140#define DEFAULT_ENABLECLIQUEMERGE FALSE /**< should the clique merging presolver be enabled within the presolve library? */
141
142/** parameters tied to a certain presolve technique in PaPILO */
143#define DEFAULT_MAXBADGESIZE_SEQ 15000 /**< the max badge size in Probing if PaPILO is executed in sequential mode */
144#define DEFAULT_MAXBADGESIZE_PAR (-1) /**< the max badge size in Probing if PaPILO is executed in parallel mode */
145#define DEFAULT_MARKOWITZTOLERANCE 0.01 /**< the markowitz tolerance used for substitutions */
146#define DEFAULT_MAXFILLINPERSUBST 3 /**< maximal possible fillin for substitutions to be considered */
147#define DEFAULT_MAXSHIFTPERROW 10 /**< maximal amount of nonzeros allowed to be shifted to make space for substitutions */
148#define DEFAULT_MAXEDGESPARALLEL 1000000 /**< maximal amount of edges in the parallel clique merging graph */
149#define DEFAULT_MAXEDGESSEQUENTIAL 100000 /**< maximal amount of edges in the sequential clique merging graph */
150#define DEFAULT_MAXCLIQUESIZE 100 /**< maximal size of clique considered for clique merging */
151#define DEFAULT_MAXGREEDYCALLS 10000 /**< maximal number of greedy max clique calls in a single thread */
152
153/** debug options for PaPILO */
154#define DEFAULT_FILENAME_PROBLEM "-" /**< default filename to store the instance before presolving */
155#define DEFAULT_VERBOSITY 0
156
157
158/*
159 * Data structures
160 */
161
162/** presolver data */
163struct SCIP_PresolMilpData
164{
165 int lastncols; /**< the number of columns from the last call */
166 int lastnrows; /**< the number of rows from the last call */
167 int threads; /**< maximum number of threads presolving may use (0: automatic) */
168 int maxfillinpersubstitution; /**< maximal possible fillin for substitutions to be considered */
169 int maxbadgesizeseq; /**< the max badge size in Probing if PaPILO is called in sequential mode */
170 int maxbadgesizepar; /**< the max badge size in Probing if PaPILO is called in parallel mode */
171 int internalmaxrounds; /**< internal max rounds in PaPILO (-1: no limit, 0: model cleanup) */
172 int maxshiftperrow; /**< maximal amount of nonzeros allowed to be shifted to make space for substitutions */
173 int detectlineardependency; /**< should linear dependent equations and free columns be removed? (0: never, 1: for LPs, 2: always) */
174#if PAPILO_APIVERSION >= 6
175 int maxedgesparallel; /**< maximal amount of edges in the parallel clique merging graph */
176 int maxedgessequential; /**< maximal amount of edges in the sequential clique merging graph */
177 int maxcliquesize; /**< maximal size of clique considered for clique merging */
178 int maxgreedycalls; /**< maximal number of greedy max clique calls in a single thread */
179#endif
180 int randomseed; /**< the random seed used for randomization of tie breaking */
181 int verbosity;
182
183 SCIP_Bool enablesparsify; /**< should the sparsify presolver be enabled within the presolve library? */
184 SCIP_Bool enabledomcol; /**< should the dominated column presolver be enabled within the presolve library? */
185 SCIP_Bool enableprobing; /**< should the probing presolver be enabled within the presolve library? */
186 SCIP_Bool enabledualinfer; /**< should the dualinfer presolver be enabled within the presolve library? */
187 SCIP_Bool enablemultiaggr; /**< should the multi-aggregation presolver be enabled within the presolve library? */
188 SCIP_Bool enableparallelrows; /**< should the parallel rows presolver be enabled within the presolve library? */
189#if PAPILO_APIVERSION >= 6
190 SCIP_Bool enablecliquemerging; /**< should the clique merging presolver be enabled within the presolve library? */
191#endif
192 SCIP_Real modifyconsfac; /**< modify SCIP constraints when the number of nonzeros or rows is at most this
193 * factor times the number of nonzeros or rows before presolving */
194 SCIP_Real markowitztolerance; /**< the markowitz tolerance used for substitutions */
195 SCIP_Real hugebound; /**< absolute bound value that is considered too huge for activitity based calculations */
196 SCIP_Real abortfacexhaustive; /**< abort factor for exhaustive presolving in PAPILO */
197 SCIP_Real abortfacmedium; /**< abort factor for medium presolving in PAPILO */
198 SCIP_Real abortfacfast; /**< abort factor for fast presolving in PAPILO */
199
200 char* filename = NULL; /**< filename to store the instance before presolving */
201};
202typedef struct SCIP_PresolMilpData SCIP_PRESOLMILPDATA;
203
204using namespace papilo;
205
206/*
207 * Local methods
208 */
209
210#if defined(PAPILO_WITH_EXACTPRESOLVE)
211/** casts rational value from PaPILO to SCIP */
212static
213void setRational(
214 SCIP* scip, /**< SCIP data structure */
215 SCIP_RATIONAL* res, /**< pointer to SCIP' rational to set */
216 papilo::Rational papiloval /**< rational value from PaPILO */
217 )
218{
219 assert(scip != NULL);
220 assert(res != NULL);
221
222 res->val = papilo::Rational(papiloval.backend().data());
225 {
227 }
228}
229
230/** builds PaPILO problem from SCIP matrix */
231static
232Problem<papilo::Rational> buildProblemRational(
233 SCIP* scip, /**< SCIP data structure */
234 SCIP_MATRIX* matrix /**< initialized SCIP_MATRIX data structure */
235 )
236{
237 ProblemBuilder<papilo::Rational> builder;
238
239 /* build problem from matrix */
240 int nnz = SCIPmatrixGetNNonzs(matrix);
241 int ncols = SCIPmatrixGetNColumns(matrix);
242 int nrows = SCIPmatrixGetNRows(matrix);
243 builder.reserve(nnz, nrows, ncols);
244
245 /* set up columns */
246 builder.setNumCols(ncols);
247 for( int i = 0; i != ncols; ++i )
248 {
249 SCIP_VAR* var = SCIPmatrixGetVar(matrix, i);
252 builder.setColLb(i, lb->val);
253 builder.setColUb(i, ub->val);
254 builder.setColLbInf(i, SCIPrationalIsNegInfinity(lb));
255 builder.setColUbInf(i, SCIPrationalIsInfinity(ub));
256
257 builder.setColIntegral(i, SCIPvarIsIntegral(var));
258 builder.setObj(i, SCIPvarGetObjExact(var)->val);
259 }
260
261 /* set up rows */
262 builder.setNumRows(nrows);
263 for( int i = 0; i != nrows; ++i )
264 {
265 int* rowcols = SCIPmatrixGetRowIdxPtr(matrix, i);
266 SCIP_RATIONAL** rowvalsscip = SCIPmatrixGetRowValPtrExact(matrix, i);
267 Vec<papilo::Rational> rowvals;
268 int rowlen = SCIPmatrixGetRowNNonzs(matrix, i);
269 for( int j = 0; j < rowlen; ++j )
270 rowvals.emplace_back(rowvalsscip[j]->val);
271 builder.addRowEntries(i, rowlen, rowcols, rowvals.data());
272
273 SCIP_RATIONAL* lhs = SCIPmatrixGetRowLhsExact(matrix, i);
274 SCIP_RATIONAL* rhs = SCIPmatrixGetRowRhsExact(matrix, i);
275 builder.setRowLhs(i, lhs->val);
276 builder.setRowRhs(i, rhs->val);
277 builder.setRowLhsInf(i, SCIPrationalIsNegInfinity(lhs));
278 builder.setRowRhsInf(i, SCIPrationalIsInfinity(rhs));
279 }
280
281 builder.setObjOffset(0);
282
283 return builder.build();
284}
285#endif
286
287/** builds PaPILO problem from SCIP matrix */
288static
289Problem<SCIP_Real> buildProblemReal(
290 SCIP* scip, /**< SCIP data structure */
291 SCIP_MATRIX* matrix /**< initialized SCIP_MATRIX data structure */
292 )
293{
294 ProblemBuilder<SCIP_Real> builder;
295
296 /* build problem from matrix */
297 int nnz = SCIPmatrixGetNNonzs(matrix);
298 int ncols = SCIPmatrixGetNColumns(matrix);
299 int nrows = SCIPmatrixGetNRows(matrix);
300 builder.reserve(nnz, nrows, ncols);
301
302 /* set up columns */
303 builder.setNumCols(ncols);
304 for( int i = 0; i != ncols; ++i )
305 {
306 SCIP_VAR* var = SCIPmatrixGetVar(matrix, i);
309 builder.setColLb(i, lb);
310 builder.setColUb(i, ub);
311 builder.setColLbInf(i, SCIPisInfinity(scip, -lb));
312 builder.setColUbInf(i, SCIPisInfinity(scip, ub));
313 builder.setColIntegral(i, SCIPvarGetType(var) != SCIP_VARTYPE_CONTINUOUS);
314#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 1)
315 builder.setColImplInt(i, SCIPvarIsImpliedIntegral(var));
316#endif
317 builder.setObj(i, SCIPvarGetObj(var));
318 }
319
320 /* set up rows */
321 builder.setNumRows(nrows);
322 for( int i = 0; i != nrows; ++i )
323 {
324 int* rowcols = SCIPmatrixGetRowIdxPtr(matrix, i);
325 SCIP_Real* rowvals = SCIPmatrixGetRowValPtr(matrix, i);
326 int rowlen = SCIPmatrixGetRowNNonzs(matrix, i);
327 builder.addRowEntries(i, rowlen, rowcols, rowvals);
328
329 SCIP_Real lhs = SCIPmatrixGetRowLhs(matrix, i);
330 SCIP_Real rhs = SCIPmatrixGetRowRhs(matrix, i);
331 builder.setRowLhs(i, lhs);
332 builder.setRowRhs(i, rhs);
333 builder.setRowLhsInf(i, SCIPisInfinity(scip, -lhs));
334 builder.setRowRhsInf(i, SCIPisInfinity(scip, rhs));
335 }
336
337 /* init objective offset - the value itself is irrelevant */
338 builder.setObjOffset(0);
339
340#ifdef SCIP_PRESOLLIB_ENABLE_OUTPUT
341 /* show problem name */
342 builder.setProblemName(SCIPgetProbName(scip));
343#endif
344
345 return builder.build();
346}
347
348/** sets up PaPILO's presolve object from the data */
349template <typename T>
350static
351SCIP_RETCODE setupPresolve(
352 SCIP* scip, /**< SCIP data structure */
353 Presolve<T>& presolve, /**< PaPILO's presolve object */
354 SCIP_PRESOLMILPDATA* data, /**< presolver data structure */
355 SCIP_Bool allowconsmodification /**< whether constraint modifications are allowed */
356 )
357{
358 SCIP_Real timelimit;
359
360 /* important so that SCIP does not throw an error, e.g. when an integer variable is substituted
361 * into a knapsack constraint */
362 presolve.getPresolveOptions().substitutebinarieswithints = false;
363
364 /* currently these changes cannot be communicated to SCIP correctly since a constraint needs
365 * to be modified in the cases where slackvariables are removed from constraints but for the
366 * presolve library those look like normal substitution on the postsolve stack */
367 presolve.getPresolveOptions().removeslackvars = false;
368
369 /* communicate the SCIP parameters to the presolve library */
370 presolve.getPresolveOptions().maxfillinpersubstitution = data->maxfillinpersubstitution;
371 presolve.getPresolveOptions().markowitz_tolerance = data->markowitztolerance;
372 presolve.getPresolveOptions().maxshiftperrow = data->maxshiftperrow;
373 presolve.getPresolveOptions().hugeval = data->hugebound;
374
375 /* removal of linear dependent equations has only an effect when constraint modifications are communicated */
376 presolve.getPresolveOptions().detectlindep = allowconsmodification ? data->detectlineardependency : 0;
377
378 /* communicate the random seed */
379 presolve.getPresolveOptions().randomseed = SCIPinitializeRandomSeed(scip, (unsigned int)data->randomseed);
380
381 /* set number of threads to be used for presolve */
382 presolve.getPresolveOptions().threads = data->threads;
383
384#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 3)
385 presolve.getPresolveOptions().maxrounds = data->internalmaxrounds;
386#endif
387
388 /* disable dual reductions that are not permitted */
390 presolve.getPresolveOptions().dualreds = 2;
391 else if( SCIPallowWeakDualReds(scip) )
392 presolve.getPresolveOptions().dualreds = 1;
393 else
394 presolve.getPresolveOptions().dualreds = 0;
395
396 /* set up the presolvers that shall participate */
397 using uptr = std::unique_ptr<PresolveMethod<T>>;
398
399 /* fast presolvers*/
400 presolve.addPresolveMethod( uptr( new SingletonCols<T>() ) );
401 presolve.addPresolveMethod( uptr( new CoefficientStrengthening<T>() ) );
402 presolve.addPresolveMethod( uptr( new ConstraintPropagation<T>() ) );
403
404 /* medium presolver */
405 presolve.addPresolveMethod( uptr( new SimpleProbing<T>() ) );
406 if( data->enableparallelrows )
407 presolve.addPresolveMethod( uptr( new ParallelRowDetection<T>() ) );
408 /* todo: parallel cols cannot be handled by SCIP currently
409 * addPresolveMethod( uptr( new ParallelColDetection<SCIP_Real>() ) ); */
410 presolve.addPresolveMethod( uptr( new SingletonStuffing<T>() ) );
411#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 1)
412 DualFix<T> *dualfix = new DualFix<T>();
413 dualfix->set_fix_to_infinity_allowed(false);
414 presolve.addPresolveMethod( uptr( dualfix ) );
415#else
416 presolve.addPresolveMethod( uptr( new DualFix<SCIP_Real>() ) );
417#endif
418 presolve.addPresolveMethod( uptr( new FixContinuous<T>() ) );
419 presolve.addPresolveMethod( uptr( new SimplifyInequalities<T>() ) );
420 presolve.addPresolveMethod( uptr( new SimpleSubstitution<T>() ) );
421#if PAPILO_APIVERSION >= 6
422 if( data->enablecliquemerging )
423 {
424 CliqueMerging<T>* cliquemerging = new CliqueMerging<T>();
425 cliquemerging->setParameters( data->maxedgesparallel, data->maxedgessequential,
426 data->maxcliquesize, data->maxgreedycalls );
427 presolve.addPresolveMethod( uptr( cliquemerging ) );
428 }
429#endif
430
431 /* exhaustive presolvers*/
432 presolve.addPresolveMethod( uptr( new ImplIntDetection<T>() ) );
433 if( data->enabledualinfer )
434 presolve.addPresolveMethod( uptr( new DualInfer<T>() ) );
435 if( data->enableprobing )
436 {
437#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 1)
438 Probing<T> *probing = new Probing<T>();
439 if( presolve.getPresolveOptions().runs_sequential() )
440 {
441 probing->set_max_badge_size( data->maxbadgesizeseq );
442 }
443 else
444 {
445 probing->set_max_badge_size( data->maxbadgesizepar );
446 }
447 presolve.addPresolveMethod( uptr( probing ) );
448#else
449 presolve.addPresolveMethod( uptr( new Probing<T>() ) );
450 if( data->maxbadgesizeseq != DEFAULT_MAXBADGESIZE_SEQ )
452 " The parameter 'presolving/milp/maxbadgesizeseq' can only be used with PaPILO 2.1.0 or later versions.\n");
453
454 if( data->maxbadgesizepar != DEFAULT_MAXBADGESIZE_PAR )
456 " The parameter 'presolving/milp/maxbadgesizepar' can only be used with PaPILO 2.1.0 or later versions.\n");
457#endif
458 }
459 if( data->enabledomcol )
460 presolve.addPresolveMethod( uptr( new DominatedCols<T>() ) );
461 if( data->enablemultiaggr )
462 presolve.addPresolveMethod( uptr( new Substitution<T>() ) );
463 if( data->enablesparsify )
464 presolve.addPresolveMethod( uptr( new Sparsify<T>() ) );
465
466 /* set numerical tolerances */
467#if PAPILO_APIVERSION >= 3
468 presolve.getPresolveOptions().useabsfeas = false;
469#endif
470 if( SCIPisExact(scip) )
471 {
472 presolve.getPresolveOptions().epsilon = 0.0;
473 presolve.getPresolveOptions().feastol = 0.0;
474 }
475 else
476 {
477 presolve.getPresolveOptions().epsilon = SCIPepsilon(scip);
478 presolve.getPresolveOptions().feastol = SCIPfeastol(scip);
479 }
480
481#ifndef SCIP_PRESOLLIB_ENABLE_OUTPUT
482 /* adjust output settings of presolve library */
483 presolve.setVerbosityLevel((VerbosityLevel) data->verbosity);
484#endif
485
486#if PAPILO_APIVERSION >= 2
487 presolve.getPresolveOptions().abortfac = data->abortfacexhaustive;
488 presolve.getPresolveOptions().abortfacmedium = data->abortfacmedium;
489 presolve.getPresolveOptions().abortfacfast = data->abortfacfast;
490#endif
491
492 /* communicate the time limit */
493 SCIPgetRealParam(scip, "limits/time", &timelimit);
494 if( !SCIPisInfinity(scip, timelimit) )
495 presolve.getPresolveOptions().tlim = timelimit - SCIPgetSolvingTime(scip);
496
497 return SCIP_OKAY;
498}
499
500
501#if defined(PAPILO_WITH_EXACTPRESOLVE)
502/** calls PaPILO presolving in rational arithmetic*/
503static
504SCIP_RETCODE performRationalPresolving(
505 SCIP* scip, /**< SCIP data structure */
506 SCIP_MATRIX* matrix, /**< initialized SCIP_MATRIX data structure */
507 SCIP_PRESOLMILPDATA* data, /**< plugin specific presol data */
508 SCIP_Bool initialized, /**< was the matrix initialized */
509 int* nfixedvars, /**< store number of fixed variables */
510 int* naggrvars, /**< store number of aggregated variables */
511 int* nchgvartypes, /**< store number of changed variable types */
512 int* nchgbds, /**< store number of changed bounds */
513 int* naddholes, /**< store number of added holes */
514 int* ndelconss, /**< store the number of deleted cons */
515 int* naddconss, /**< store the number of added cons */
516 int* nupgdconss, /**< store the number of upgraded cons */
517 int* nchgcoefs, /**< store the number of changed coefficients */
518 int* nchgsides, /**< store the number of changed sides */
519 SCIP_RESULT* result /**< result pointer */
520 )
521{
522 int nvars = SCIPgetNVars(scip);
523 int nconss = SCIPgetNConss(scip);
524
525 SCIP_CONSHDLR* linconshdlr = SCIPfindConshdlr(scip, "exactlinear");
526 assert(linconshdlr != NULL);
527 bool allowconsmodification = (SCIPconshdlrGetNCheckConss(linconshdlr) == SCIPmatrixGetNRows(matrix));
528
529 /* store current numbers of aggregations, fixings, and changed bounds for statistics */
530 int oldnaggrvars = *naggrvars;
531 int oldnfixedvars = *nfixedvars;
532 int oldnchgbds = *nchgbds;
533
534 Problem<papilo::Rational> problem = buildProblemRational(scip, matrix);
535 Presolve<papilo::Rational> presolve;
536 setupPresolve(scip, presolve, data, allowconsmodification);
537
538 /* call presolving (without storing information for dual postsolve) */
540 " (%.1fs) running MILP presolver%s\n", SCIPgetSolvingTime(scip),
541 presolve.getPresolveOptions().threads == 1 ? "" : " on multiple threads");
542
543 int oldnnz = problem.getConstraintMatrix().getNnz();
544
545#if (PAPILO_VERSION_MAJOR >= 2)
546 PresolveResult<papilo::Rational> res = presolve.apply(problem, false);
547#else
548 PresolveResult<papilo::Rational> res = presolve.apply(problem);
549#endif
550 data->lastncols = problem.getNCols();
551 data->lastnrows = problem.getNRows();
552
553 /* evaluate the result */
554 switch( res.status )
555 {
556 case PresolveStatus::kInfeasible:
557 *result = SCIP_CUTOFF;
559 " (%.1fs) MILP presolver detected infeasibility\n",
561 SCIPmatrixFree(scip, &matrix);
562 return SCIP_OKAY;
563 case PresolveStatus::kUnbndOrInfeas:
564 case PresolveStatus::kUnbounded:
565 *result = SCIP_UNBOUNDED;
567 " (%.1fs) MILP presolver detected unboundedness\n",
569 SCIPmatrixFree(scip, &matrix);
570 return SCIP_OKAY;
571 case PresolveStatus::kUnchanged:
572 *result = SCIP_DIDNOTFIND;
573 data->lastncols = nvars;
574 data->lastnrows = nconss;
576 " (%.1fs) MILP presolver found nothing\n",
578 SCIPmatrixFree(scip, &matrix);
579 return SCIP_OKAY;
580 case PresolveStatus::kReduced:
581 data->lastncols = problem.getNCols();
582 data->lastnrows = problem.getNRows();
583 *result = SCIP_SUCCESS;
584 }
585
586 /* result indicated success, now populate the changes into the SCIP structures */
587 Vec<SCIP_VAR*> tmpvars;
588 Vec<SCIP_Real> tmpvalsreal;
589
590 /* if the number of nonzeros decreased by a sufficient factor, rather create all constraints from scratch */
591 int newnnz = problem.getConstraintMatrix().getNnz();
592 bool constraintsReplaced = false;
593 if( newnnz == 0 || (allowconsmodification &&
594 (problem.getNRows() <= data->modifyconsfac * data->lastnrows ||
595 newnnz <= data->modifyconsfac * oldnnz)) )
596 {
597 int oldnrows = SCIPmatrixGetNRows(matrix);
598 int newnrows = problem.getNRows();
599
600 constraintsReplaced = true;
601
602 /* capture constraints that are still present in the problem after presolve */
603 for( int i = 0; i < newnrows; ++i )
604 {
605 SCIP_CONS* c = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
607 }
608
609 /* delete all constraints */
610 *ndelconss += oldnrows;
611 *naddconss += newnrows;
612
613 for( int i = 0; i < oldnrows; ++i )
614 {
616 }
617
618 /* now loop over rows of presolved problem and create them as new linear constraints,
619 * then release the old constraint after its name was passed to the new constraint
620 */
621 const Vec<RowFlags>& rflags = problem.getRowFlags();
622 const auto& consmatrix = problem.getConstraintMatrix();
623 for( int i = 0; i < newnrows; ++i )
624 {
625 auto rowvec = consmatrix.getRowCoefficients(i);
626 const int* rowcols = rowvec.getIndices();
627 /* SCIPcreateConsBasicLinear() requires a non const pointer */
628 papilo::Rational* rowvals = const_cast<papilo::Rational*>(rowvec.getValues());
629 int rowlen = rowvec.getLength();
630
631 /* retrieve SCIP compatible left and right hand sides */
632 papilo::Rational lhs = rflags[i].test(RowFlag::kLhsInf) ? - SCIPinfinity(scip) : consmatrix.getLeftHandSides()[i];
633 papilo::Rational rhs = rflags[i].test(RowFlag::kRhsInf) ? SCIPinfinity(scip) : consmatrix.getRightHandSides()[i];
634
635 /* create variable array matching the value array */
636 tmpvars.clear();
637 tmpvars.reserve(rowlen);
638 for( int j = 0; j < rowlen; ++j )
639 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[rowcols[j]]));
640
641 /* create and add new constraint with name of old constraint */
642 SCIP_CONS* oldcons = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
643 SCIP_CONS* cons;
644 SCIP_RATIONAL** tmpvals;
645 SCIP_RATIONAL* tmplhs;
646 SCIP_RATIONAL* tmprhs;
647
651
652 for( int j = 0; j < rowlen; j++ )
653 setRational(scip, tmpvals[j], rowvals[j]);
654
655 setRational(scip, tmprhs, rhs);
656 setRational(scip, tmplhs, lhs);
657 if( rflags[i].test(RowFlag::kLhsInf) )
659 if( rflags[i].test(RowFlag::kRhsInf) )
661
662 SCIP_CALL( SCIPcreateConsBasicExactLinear(scip, &cons, SCIPconsGetName(oldcons), rowlen, tmpvars.data(), tmpvals, tmplhs, tmprhs) );
663 SCIP_CALL( SCIPaddCons(scip, cons) );
664
665 /* release old and new constraint */
666 SCIP_CALL( SCIPreleaseCons(scip, &oldcons) );
667 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
668
671 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &tmpvals, rowlen);
672 }
673 }
674
675 /* loop over res.postsolve and add all fixed variables and aggregations to scip */
676 for( std::size_t i = 0; i != res.postsolve.types.size(); ++i )
677 {
678 ReductionType type =
679 res.postsolve.types[i];
680 int first = res.postsolve.start[i];
681 int last = res.postsolve.start[i + 1];
682
683 switch( type )
684 {
685 case ReductionType::kFixedCol:
686 {
687 SCIP_RATIONAL* tmpval;
688 SCIP_Bool infeas;
689 SCIP_Bool fixed;
690 int col = res.postsolve.indices[first];
691
692 SCIP_VAR* var = SCIPmatrixGetVar(matrix, col);
693
694 papilo::Rational value = res.postsolve.values[first];
695
697 setRational(scip, tmpval, value);
698
699 SCIPrationalDebugMessage("Papilo fix var %s to %q \n", SCIPvarGetName(var), tmpval);
700
701 /* SCIP has different rules for aggregation than PaPILO
702 * As a result, SCIP might have aggregated and replaced the variable that PaPILO now wants to fix
703 */
705 {
706 SCIP_RATIONAL* aggregatedScalar;
707 SCIP_RATIONAL* aggregatedConst;
708
709 aggregatedScalar = SCIPvarGetAggrScalarExact(var);
710 aggregatedConst = SCIPvarGetAggrConstantExact(var);
711
712 /* fix aggregation variable y in x = a*y + c, instead of fixing x directly */
714 assert( !SCIPrationalIsZero(aggregatedScalar));
715 if( SCIPrationalIsAbsInfinity(tmpval) )
716 SCIPrationalMultReal(tmpval, tmpval, SCIPrationalIsNegative(aggregatedScalar) ? -1 : 1);
717 else
718 {
719 SCIPrationalDiff(tmpval, tmpval, aggregatedConst);
720 SCIPrationalDiv(tmpval, tmpval, aggregatedScalar);
721 }
722 var = SCIPvarGetAggrVar(var);
723 }
724
725 /* SCIP might also have fixed the variable during aggregation */
727 break;
728
729 SCIP_CALL( SCIPfixVarExact(scip, var, tmpval, &infeas, &fixed) );
730
731 *nfixedvars += 1;
732
734
735 assert(!infeas);
736 assert(fixed || SCIPvarGetStatus(var) == SCIP_VARSTATUS_AGGREGATED ||
738 break;
739 }
740 /*
741 * Dual-postsolving in PaPILO required introducing a postsolve-type for substitution with additional information.
742 * Further, the different Substitution-postsolving types store the required postsolving data differently (in different order) in the postsolving stack.
743 * Therefore, we need to distinguish how to parse the required data (rowLength, col, side, startRowCoefficients, lastRowCoefficients) from the postsolving stack.
744 * If these values are accessed, the procedure is the same for both.
745 */
746#if (PAPILO_VERSION_MAJOR >= 2)
747 case ReductionType::kSubstitutedColWithDual:
748#endif
749 case ReductionType::kSubstitutedCol:
750 {
751 int col = 0;
752 papilo::Rational side = 0;
753
754 int rowlen = 0;
755 int startRowCoefficients = 0;
756 int lastRowCoefficients = 0;
757
758 if( type == ReductionType::kSubstitutedCol )
759 {
760 rowlen = last - first - 1;
761 col = res.postsolve.indices[first];
762 side = res.postsolve.values[first];
763
764 startRowCoefficients = first + 1;
765 lastRowCoefficients = last;
766 }
767#if (PAPILO_VERSION_MAJOR >= 2)
768 if( type == ReductionType::kSubstitutedColWithDual )
769 {
770 rowlen = (int) res.postsolve.values[first];
771 col = res.postsolve.indices[first + 3 + rowlen];
772 side = res.postsolve.values[first + 1];
773
774 startRowCoefficients = first + 3;
775 lastRowCoefficients = first + 3 + rowlen;
776
777 assert(side == res.postsolve.values[first + 2]);
778 assert(res.postsolve.indices[first + 1] == 0);
779 assert(res.postsolve.indices[first + 2] == 0);
780 }
781 assert( type == ReductionType::kSubstitutedCol || type == ReductionType::kSubstitutedColWithDual );
782#else
783 assert( type == ReductionType::kSubstitutedCol );
784#endif
785 SCIP_Bool infeas;
786 SCIP_Bool aggregated;
787 SCIP_Bool redundant = FALSE;
788 SCIP_RATIONAL* constant;
789 if( rowlen == 2 )
790 {
791 SCIP_VAR* varx = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients]);
792 SCIP_VAR* vary = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients + 1]);
793 papilo::Rational scalarx = res.postsolve.values[startRowCoefficients];
794 papilo::Rational scalary = res.postsolve.values[startRowCoefficients + 1];
795
796 SCIP_RATIONAL* tmpscalarx;
797 SCIP_RATIONAL* tmpscalary;
798 SCIP_RATIONAL* tmpside;
803
804 setRational(scip, tmpscalarx, scalarx);
805 setRational(scip, tmpscalary, scalary);
806
807 SCIP_CALL( SCIPgetProbvarSumExact(scip, &varx, tmpscalarx, constant) );
809
810 SCIP_CALL( SCIPgetProbvarSumExact(scip, &vary, tmpscalary, constant) );
812
813 setRational(scip, tmpside, side);
814 SCIPrationalDiff(tmpside, tmpside, constant);
815
816 SCIPrationalDebugMessage("Papilo aggregate vars %s, %s with scalars %q, %q and constant %q \n", SCIPvarGetName(varx), SCIPvarGetName(vary),
817 tmpscalarx, tmpscalary, constant);
818
819 SCIP_CALL( SCIPaggregateVarsExact(scip, varx, vary, tmpscalarx, tmpscalary, tmpside, &infeas, &redundant, &aggregated) );
820
825 }
826 else
827 {
828 SCIP_RATIONAL* colCoef;
829 SCIP_RATIONAL* updatedSide;
830 SCIP_RATIONAL** tmpvals;
831 int c = 0;
832
837
838 for( int j = startRowCoefficients; j < lastRowCoefficients; ++j )
839 {
840 if( res.postsolve.indices[j] == col )
841 {
842 setRational(scip, colCoef, res.postsolve.values[j]);
843 break;
844 }
845 }
846
847 tmpvars.clear();
848 tmpvars.reserve(rowlen);
849
850 assert(!SCIPrationalIsZero(colCoef));
851 SCIP_VAR* aggrvar = SCIPmatrixGetVar(matrix, col);
852
853 SCIP_CALL( SCIPgetProbvarSumExact(scip, &aggrvar, colCoef, constant) );
854 assert(SCIPvarGetStatus(aggrvar) != SCIP_VARSTATUS_MULTAGGR);
855
856 for( int j = startRowCoefficients; j < lastRowCoefficients; ++j )
857 {
858 if( res.postsolve.indices[j] == col )
859 continue;
860
861 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.indices[j]));
862 setRational(scip, tmpvals[c], -res.postsolve.values[j] / colCoef->val);
863 c++;
864 }
865 setRational(scip, updatedSide, side);
866 SCIPrationalDiff(updatedSide, updatedSide, constant);
867 SCIPrationalDiv(updatedSide, updatedSide, colCoef);
868
869 SCIPrationalDebugMessage("Papilo multiaggregate var %s, constant %q \n", SCIPvarGetName(aggrvar), updatedSide);
870
871 SCIP_CALL( SCIPmultiaggregateVarExact(scip, aggrvar, tmpvars.size(),
872 tmpvars.data(), tmpvals, updatedSide, &infeas, &aggregated) );
873
874 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &tmpvals, rowlen);
875 SCIPrationalFreeBuffer(SCIPbuffer(scip), &updatedSide);
878 }
879
880 if( aggregated )
881 *naggrvars += 1;
882 else if( constraintsReplaced && !redundant )
883 {
884 /* if the constraints where replaced, we need to add the failed substitution as an equality to SCIP */
885 SCIP_RATIONAL** tmpvals;
886 SCIP_RATIONAL* tmpside;
887
890
891 setRational(scip, tmpside, side);
892
893 tmpvars.clear();
894 for( int j = startRowCoefficients; j < lastRowCoefficients; ++j )
895 {
896 int idx = j - startRowCoefficients;
897 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.indices[j]));
898 setRational(scip, tmpvals[idx], res.postsolve.values[j]);
899 }
900
901 SCIP_CONS* cons;
902 String name = fmt::format("{}_failed_aggregation_equality", SCIPvarGetName(SCIPmatrixGetVar(matrix, col)));
903 SCIP_CALL( SCIPcreateConsBasicExactLinear(scip, &cons, name.c_str(),
904 tmpvars.size(), tmpvars.data(), tmpvals, tmpside, tmpside ) );
905
906 SCIPdebugMessage("Papilo adding failed aggregation equality: \n");
908 SCIP_CALL( SCIPaddCons(scip, cons) );
909 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
910 *naddconss += 1;
911
913 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &tmpvals, rowlen);
914 }
915
916 if( infeas )
917 {
918 *result = SCIP_CUTOFF;
919 break;
920 }
921
922 break;
923 }
924 case ReductionType::kParallelCol:
925 return SCIP_INVALIDRESULT;
926#if (PAPILO_VERSION_MAJOR <= 1 && PAPILO_VERSION_MINOR==0)
927#else
928 case ReductionType::kFixedInfCol: {
929 if(!constraintsReplaced)
930 continue;
931 SCIP_Bool infeas;
932 SCIP_Bool fixed;
933 SCIP_RATIONAL* value;
934
937
938 int column = res.postsolve.indices[first];
939 bool is_negative_infinity = res.postsolve.values[first] < 0;
940 SCIP_VAR* column_variable = SCIPmatrixGetVar(matrix, column);
941
942 if( is_negative_infinity )
943 {
945 }
946
947 SCIP_CALL( SCIPfixVarExact(scip, column_variable, value, &infeas, &fixed) );
948 *nfixedvars += 1;
949
950 assert(!infeas);
951 assert(fixed);
952
954
955 break;
956 }
957#endif
958#if (PAPILO_VERSION_MAJOR >= 2)
959 case ReductionType::kVarBoundChange :
960 case ReductionType::kRedundantRow :
961 case ReductionType::kRowBoundChange :
962 case ReductionType::kReasonForRowBoundChangeForcedByRow :
963 case ReductionType::kRowBoundChangeForcedByRow :
964 case ReductionType::kSaveRow :
965 case ReductionType::kReducedBoundsCost :
966 case ReductionType::kColumnDualValue :
967 case ReductionType::kRowDualValue :
968 case ReductionType::kCoefficientChange :
969 /* dual ReductionTypes should be only calculated for dual reductions and should not appear for MIP */
970 SCIPerrorMessage("PaPILO: PaPILO should not return dual postsolving reductions in SCIP!!\n");
971 SCIPABORT(); /*lint --e{527}*/
972 break;
973#endif
974 default:
975 SCIPdebugMsg(scip, "PaPILO returned unknown data type: \n" );
976 continue;
977 }
978 }
979
980 /* tighten bounds of variables that are still present after presolving */
981 if( *result != SCIP_CUTOFF )
982 {
983 VariableDomains<papilo::Rational>& varDomains = problem.getVariableDomains();
984 SCIP_RATIONAL* varbound;
987
988 for( int i = 0; i != problem.getNCols(); ++i )
989 {
990 SCIP_VAR* var = SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[i]);
991 if( !varDomains.flags[i].test(ColFlag::kLbInf) )
992 {
993 SCIP_Bool infeas;
994 SCIP_Bool tightened;
995
996 setRational(scip, varbound, varDomains.lower_bounds[i]);
997
998 SCIP_CALL( SCIPtightenVarLbExact(scip, var, varbound, &infeas, &tightened) );
999
1000 if( tightened )
1001 {
1002 *nchgbds += 1;
1003 SCIPrationalDebugMessage("Papilo tightened lb of variable %s \n", SCIPvarGetName(var));
1004 }
1005
1006 if( infeas )
1007 {
1008 *result = SCIP_CUTOFF;
1009 break;
1010 }
1011 }
1012
1013 if( !varDomains.flags[i].test(ColFlag::kUbInf) )
1014 {
1015 SCIP_Bool infeas;
1016 SCIP_Bool tightened;
1017 setRational(scip, varbound, varDomains.upper_bounds[i]);
1018
1019 SCIP_CALL( SCIPtightenVarUbExact(scip, var, varbound, &infeas, &tightened) );
1020
1021 if( tightened )
1022 {
1023 *nchgbds += 1;
1024 SCIPrationalDebugMessage("Papilo tightened ub of variable %s \n", SCIPvarGetName(var));
1025 }
1026
1027 if( infeas )
1028 {
1029 *result = SCIP_CUTOFF;
1030 break;
1031 }
1032 }
1033 }
1034
1036 }
1037
1038 /* finish with a final verb message and return */
1040 " (%.1fs) MILP presolver (%d rounds): %d aggregations, %d fixings, %d bound changes\n",
1041 SCIPgetSolvingTime(scip), presolve.getStatistics().nrounds, *naggrvars - oldnaggrvars,
1042 *nfixedvars - oldnfixedvars, *nchgbds - oldnchgbds);
1043
1044 /* free the matrix */
1045 assert(initialized);
1046 SCIPmatrixFree(scip, &matrix);
1047
1048 return SCIP_OKAY;
1049}
1050#endif
1051
1052/** calls PaPILO presolving in floating-point arithmetic */
1053static
1054SCIP_RETCODE performRealPresolving(
1055 SCIP* scip, /**< SCIP data structure */
1056 SCIP_MATRIX* matrix, /**< initialized SCIP_MATRIX data structure */
1057 SCIP_PRESOLMILPDATA* data, /**< plugin specific presol data */
1058 SCIP_Bool initialized, /**< was the matrix initialized */
1059 int* nfixedvars, /**< store number of fixed variables */
1060 int* naggrvars, /**< store number of aggregated variables */
1061 int* nchgvartypes, /**< store number of changed variable types */
1062 int* nchgbds, /**< store number of changed bounds */
1063 int* naddholes, /**< store number of added holes */
1064 int* ndelconss, /**< store the number of deleted cons */
1065 int* naddconss, /**< store the number of added cons */
1066 int* nupgdconss, /**< store the number of upgraded cons */
1067 int* nchgcoefs, /**< store the number of changed coefficients */
1068 int* nchgsides, /**< store the number of changed sides */
1069 SCIP_RESULT* result /**< result pointer */
1070 )
1071{
1072 int nvars = SCIPgetNVars(scip);
1073 int nconss = SCIPgetNConss(scip);
1074
1075 /* only allow communication of constraint modifications by deleting all constraints when some already have been upgraded */
1076 SCIP_CONSHDLR* linconshdlr = SCIPfindConshdlr(scip, "linear");
1077 assert(linconshdlr != NULL);
1078 SCIP_Bool allowconsmodification = (SCIPconshdlrGetNCheckConss(linconshdlr) == SCIPmatrixGetNRows(matrix));
1079
1080 /* store current numbers of aggregations, fixings, and changed bounds for statistics */
1081 int oldnaggrvars = *naggrvars;
1082 int oldnfixedvars = *nfixedvars;
1083 int oldnchgbds = *nchgbds;
1084
1085 /* create presolving objects */
1086 Problem<SCIP_Real> problem = buildProblemReal(scip, matrix);
1087 int oldnnz = problem.getConstraintMatrix().getNnz();
1088 Presolve<SCIP_Real> presolve;
1089 setupPresolve(scip, presolve, data, allowconsmodification);
1090
1091 /* call the presolving */
1093 " (%.1fs) running MILP presolver%s\n", SCIPgetSolvingTime(scip),
1094 presolve.getPresolveOptions().threads == 1 ? "" : " on multiple threads");
1095
1096 /* call presolving without storing information for dual postsolve */
1097#if (PAPILO_VERSION_MAJOR >= 2)
1098 PresolveResult<SCIP_Real> res = presolve.apply(problem, false);
1099#else
1100 PresolveResult<SCIP_Real> res = presolve.apply(problem);
1101#endif
1102 data->lastncols = problem.getNCols();
1103 data->lastnrows = problem.getNRows();
1104
1105 /* evaluate the result */
1106 switch( res.status )
1107 {
1108 case PresolveStatus::kInfeasible:
1109 *result = SCIP_CUTOFF;
1111 " (%.1fs) MILP presolver detected infeasibility\n",
1113 SCIPmatrixFree(scip, &matrix);
1114 return SCIP_OKAY;
1115 case PresolveStatus::kUnbndOrInfeas:
1116 case PresolveStatus::kUnbounded:
1117 *result = SCIP_UNBOUNDED;
1119 " (%.1fs) MILP presolver detected unboundedness\n",
1121 SCIPmatrixFree(scip, &matrix);
1122 return SCIP_OKAY;
1123 case PresolveStatus::kUnchanged:
1124 *result = SCIP_DIDNOTFIND;
1125 data->lastncols = nvars;
1126 data->lastnrows = nconss;
1128 " (%.1fs) MILP presolver found nothing\n",
1130 SCIPmatrixFree(scip, &matrix);
1131 return SCIP_OKAY;
1132 case PresolveStatus::kReduced:
1133 data->lastncols = problem.getNCols();
1134 data->lastnrows = problem.getNRows();
1135 *result = SCIP_SUCCESS;
1136 }
1137
1138 /* result indicated success, now populate the changes into the SCIP structures */
1139
1140 /* tighten bounds of variables that are still present after presolving */
1141 VariableDomains<SCIP_Real>& varDomains = problem.getVariableDomains();
1142 for( int i = 0; i != problem.getNCols(); ++i )
1143 {
1144 assert( ! varDomains.flags[i].test(ColFlag::kInactive) );
1145 SCIP_VAR* var = SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[i]);
1146 if( !varDomains.flags[i].test(ColFlag::kLbInf) )
1147 {
1148 SCIP_Bool infeas;
1149 SCIP_Bool tightened;
1150 SCIP_CALL( SCIPtightenVarLb(scip, var, varDomains.lower_bounds[i], TRUE, &infeas, &tightened) );
1151
1152 if( tightened )
1153 *nchgbds += 1;
1154
1155 if( infeas )
1156 {
1157 *result = SCIP_CUTOFF;
1158 break;
1159 }
1160 }
1161
1162 if( !varDomains.flags[i].test(ColFlag::kUbInf) )
1163 {
1164 SCIP_Bool infeas;
1165 SCIP_Bool tightened;
1166 SCIP_CALL( SCIPtightenVarUb(scip, var, varDomains.upper_bounds[i], TRUE, &infeas, &tightened) );
1167
1168 if( tightened )
1169 *nchgbds += 1;
1170
1171 if( infeas )
1172 {
1173 *result = SCIP_CUTOFF;
1174 break;
1175 }
1176 }
1177 }
1178
1179 if( *result == SCIP_CUTOFF )
1180 {
1182 " (%.1fs) MILP presolver detected infeasibility\n",
1184 SCIPmatrixFree(scip, &matrix);
1185 return SCIP_OKAY;
1186 }
1187
1188 /* transfer variable fixings and aggregations */
1189 Vec<SCIP_VAR*> tmpvars;
1190 Vec<SCIP_Real> tmpvals;
1191
1192 /* if the size of the problem decreased by a sufficient factor, create all constraints from scratch if allowed */
1193 int newnnz = problem.getConstraintMatrix().getNnz();
1194 bool constraintsReplaced = false;
1195 if( newnnz == 0 || (allowconsmodification &&
1196 (problem.getNCols() <= data->modifyconsfac * SCIPmatrixGetNColumns(matrix) ||
1197 problem.getNRows() <= data->modifyconsfac * SCIPmatrixGetNRows(matrix) ||
1198 newnnz <= data->modifyconsfac * oldnnz)) )
1199 {
1200 int oldnrows = SCIPmatrixGetNRows(matrix);
1201 int newnrows = problem.getNRows();
1202
1203 constraintsReplaced = true;
1204
1205 /* capture constraints that are still present in the problem after presolve */
1206 for( int i = 0; i < newnrows; ++i )
1207 {
1208 SCIP_CONS* c = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
1210 }
1211
1212 /* delete all constraints */
1213 *ndelconss += oldnrows;
1214 *naddconss += newnrows;
1215
1216 for( int i = 0; i < oldnrows; ++i )
1217 {
1219 }
1220
1221 /* now loop over rows of presolved problem and create them as new linear constraints,
1222 * then release the old constraint after its name was passed to the new constraint */
1223 const Vec<RowFlags>& rflags = problem.getRowFlags();
1224 const auto& consmatrix = problem.getConstraintMatrix();
1225 for( int i = 0; i < newnrows; ++i )
1226 {
1227 auto rowvec = consmatrix.getRowCoefficients(i);
1228 const int* rowcols = rowvec.getIndices();
1229 /* SCIPcreateConsBasicLinear() requires a non const pointer */
1230 SCIP_Real* rowvals = const_cast<SCIP_Real*>(rowvec.getValues());
1231 int rowlen = rowvec.getLength();
1232
1233 /* retrieve SCIP compatible left and right hand sides */
1234 SCIP_Real lhs = rflags[i].test(RowFlag::kLhsInf) ? - SCIPinfinity(scip) : consmatrix.getLeftHandSides()[i];
1235 SCIP_Real rhs = rflags[i].test(RowFlag::kRhsInf) ? SCIPinfinity(scip) : consmatrix.getRightHandSides()[i];
1236
1237 /* create variable array matching the value array */
1238 tmpvars.clear();
1239 tmpvars.reserve(rowlen);
1240 for( int j = 0; j < rowlen; ++j )
1241 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[rowcols[j]]));
1242
1243 /* create and add new constraint with name of old constraint */
1244 SCIP_CONS* oldcons = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
1245 SCIP_CONS* cons;
1246 SCIP_CALL( SCIPcreateConsBasicLinear(scip, &cons, SCIPconsGetName(oldcons), rowlen, tmpvars.data(), rowvals, lhs, rhs) );
1247 SCIP_CALL( SCIPaddCons(scip, cons) );
1248
1249 /* release old and new constraint */
1250 SCIP_CALL( SCIPreleaseCons(scip, &oldcons) );
1251 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
1252 }
1253 }
1254
1255 /* PaPILO's aggregations are valid regarding the constraints as they were presolved by PaPILO.
1256 * If coefficients were changed, but constraints in SCIP are not replaced by those from PaPILO,
1257 * then it can not be guaranteed that the bounds of multiaggregated variables will be enforced,
1258 * i.e., will be implied by the constraints in SCIP (see also #3704).
1259 * Only for variable aggregations, SCIP will ensure this by tightening the bounds on the aggregation
1260 * variable as part of SCIPaggregateVars(). For multiaggregations, we will only accept those
1261 * where we can be sure with a simple check that the bounds on the aggregated variable are implied.
1262 */
1263 bool checkmultaggr =
1264#if PAPILO_APIVERSION >= 1
1265 presolve.getStatistics().single_matrix_coefficient_changes > 0
1266#else
1267 presolve.getStatistics().ncoefchgs > 0
1268#endif
1269 && !constraintsReplaced;
1270
1271 /* loop over res.postsolve and add all fixed variables and aggregations to scip */
1272 for( std::size_t i = 0; i != res.postsolve.types.size(); ++i )
1273 {
1274 ReductionType type = res.postsolve.types[i];
1275 int first = res.postsolve.start[i];
1276 int last = res.postsolve.start[i + 1];
1277
1278 switch( type )
1279 {
1280 case ReductionType::kFixedCol:
1281 {
1282 SCIP_Bool infeas;
1283 SCIP_Bool fixed;
1284 int col = res.postsolve.indices[first];
1285
1286 SCIP_VAR* var = SCIPmatrixGetVar(matrix, col);
1287
1288 SCIP_Real value = res.postsolve.values[first];
1289
1290 SCIP_CALL( SCIPfixVar(scip, var, value, &infeas, &fixed) );
1291 *nfixedvars += 1;
1292
1293 assert(!infeas);
1294 /* SCIP has different rules for aggregating variables than PaPILO: therefore the variable PaPILO
1295 * tries to fix now may have been aggregated by SCIP before. Additionally, after aggregation SCIP
1296 * sometimes performs bound tightening resulting in possible fixings. These cases need to be excluded. */
1297 assert(fixed || SCIPvarGetStatus(var) == SCIP_VARSTATUS_AGGREGATED ||
1299 break;
1300 }
1301 /*
1302 * Dual-postsolving in PaPILO required introducing a postsolve-type for substitution with additional information.
1303 * Further, the different Substitution-postsolving types store the required postsolving data differently (in different order) in the postsolving stack.
1304 * Therefore, we need to distinguish how to parse the required data (rowLength, col, side, startRowCoefficients, lastRowCoefficients) from the postsolving stack.
1305 * If these values are accessed, the procedure is the same for both.
1306 */
1307#if (PAPILO_VERSION_MAJOR >= 2)
1308 case ReductionType::kSubstitutedColWithDual:
1309#endif
1310 case ReductionType::kSubstitutedCol:
1311 {
1312 int col = 0;
1313 SCIP_Real side = 0;
1314
1315 int rowlen = 0;
1316 int startRowCoefficients = 0;
1317 int lastRowCoefficients = 0;
1318
1319 if( type == ReductionType::kSubstitutedCol )
1320 {
1321 rowlen = last - first - 1;
1322 col = res.postsolve.indices[first];
1323 side = res.postsolve.values[first];
1324
1325 startRowCoefficients = first + 1;
1326 lastRowCoefficients = last;
1327 }
1328#if (PAPILO_VERSION_MAJOR >= 2)
1329 if( type == ReductionType::kSubstitutedColWithDual )
1330 {
1331 rowlen = (int) res.postsolve.values[first];
1332 col = res.postsolve.indices[first + 3 + rowlen];
1333 side = res.postsolve.values[first + 1];
1334
1335 startRowCoefficients = first + 3;
1336 lastRowCoefficients = first + 3 + rowlen;
1337
1338 assert(side == res.postsolve.values[first + 2]);
1339 assert(res.postsolve.indices[first + 1] == 0);
1340 assert(res.postsolve.indices[first + 2] == 0);
1341 }
1342 assert( type == ReductionType::kSubstitutedCol || type == ReductionType::kSubstitutedColWithDual );
1343#else
1344 assert( type == ReductionType::kSubstitutedCol );
1345#endif
1346 SCIP_Bool infeas;
1347 SCIP_Bool aggregated;
1348 SCIP_Bool redundant = FALSE;
1349 SCIP_Real constant = 0.0;
1350 if( rowlen == 2 )
1351 {
1352 SCIP_Real updatedSide;
1353 SCIP_VAR* varx = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients]);
1354 SCIP_VAR* vary = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients + 1]);
1355 SCIP_Real scalarx = res.postsolve.values[startRowCoefficients];
1356 SCIP_Real scalary = res.postsolve.values[startRowCoefficients + 1];
1357
1358 SCIP_CALL( SCIPgetProbvarSum(scip, &varx, &scalarx, &constant) );
1360
1361 SCIP_CALL( SCIPgetProbvarSum(scip, &vary, &scalary, &constant) );
1363
1364 /* If PaPILO tries to aggregate fixed variables then it missed some obvious fixings.
1365 * This might happen if another aggregation leads to fixings which are not applied immediately by PaPILO.
1366 * With the latest version of PaPILO, this should not occur.
1367 */
1369 {
1370 SCIPdebugMsg(scip, "Aggregation of <%s> and <%s> rejected because they are already fixed.\n",
1371 SCIPvarGetName(varx), SCIPvarGetName(vary));
1372
1373 break;
1374 }
1375
1376 updatedSide = side - constant;
1377
1378 SCIP_CALL( SCIPaggregateVars(scip, varx, vary, scalarx, scalary, updatedSide, &infeas, &redundant, &aggregated) );
1379 }
1380 else
1381 {
1382 SCIP_Real colCoef = 0.0;
1383 SCIP_Real updatedSide;
1384 SCIP_Bool checklbimplied;
1385 SCIP_Bool checkubimplied;
1386 SCIP_Real impliedlb;
1387 SCIP_Real impliedub;
1388 int j;
1389
1390 for( j = startRowCoefficients; j < lastRowCoefficients; ++j )
1391 {
1392 if( res.postsolve.indices[j] == col )
1393 {
1394 colCoef = res.postsolve.values[j];
1395 break;
1396 }
1397 }
1398
1399 tmpvars.clear();
1400 tmpvals.clear();
1401 tmpvars.reserve(rowlen);
1402 tmpvals.reserve(rowlen);
1403
1404 assert(colCoef != 0.0);
1405 SCIP_VAR* aggrvar = SCIPmatrixGetVar(matrix, col);
1406
1407 SCIP_CALL( SCIPgetProbvarSum(scip, &aggrvar, &colCoef, &constant) );
1408 assert(SCIPvarGetStatus(aggrvar) != SCIP_VARSTATUS_MULTAGGR);
1409
1410 /* If PaPILO tries to multi-aggregate a fixed variable, then it missed some obvious fixings.
1411 * This might happen if another aggregation leads to fixings which are not applied immediately by PaPILO.
1412 * With the latest version of PaPILO, this should not occur.
1413 */
1414 if( SCIPvarGetStatus(aggrvar) == SCIP_VARSTATUS_FIXED )
1415 {
1416 SCIPdebugMsg(scip, "Multi-aggregation of <%s> rejected because it is already fixed.\n",
1417 SCIPvarGetName(aggrvar));
1418
1419 break;
1420 }
1421
1422 updatedSide = side - constant;
1423
1424 /* we need to check whether lb/ub on aggrvar is implied by bounds of other variables in multiaggregation
1425 * if checkmultaggr is TRUE and the lb/ub is finite
1426 * it should be sufficient to ensure global bounds on aggrvar (and as we are in presolve, local=global anyway)
1427 */
1428 checklbimplied = checkmultaggr && !SCIPisInfinity(scip, -SCIPvarGetLbGlobal(aggrvar));
1429 checkubimplied = checkmultaggr && !SCIPisInfinity(scip, SCIPvarGetUbGlobal(aggrvar));
1430 impliedlb = impliedub = updatedSide / colCoef;
1431
1432 for( j = startRowCoefficients; j < lastRowCoefficients; ++j )
1433 {
1434 SCIP_Real coef;
1435 SCIP_VAR* var;
1436
1437 if( res.postsolve.indices[j] == col )
1438 continue;
1439
1440 coef = - res.postsolve.values[j] / colCoef;
1441 var = SCIPmatrixGetVar(matrix, res.postsolve.indices[j]);
1442
1443 if( checklbimplied )
1444 {
1445 if( coef > 0.0 )
1446 {
1447 /* if impliedlb will be -infinity, then we can give up: we cannot use this mutiaggregation */
1449 break;
1450 else
1451 impliedlb += coef * SCIPvarGetLbLocal(var);
1452 }
1453 else
1454 {
1456 break;
1457 else
1458 impliedlb += coef * SCIPvarGetUbLocal(var);
1459 }
1460 }
1461
1462 if( checkubimplied )
1463 {
1464 if( coef > 0.0 )
1465 {
1467 break;
1468 else
1469 impliedub += coef * SCIPvarGetUbLocal(var);
1470 }
1471 else
1472 {
1474 break;
1475 else
1476 impliedub += coef * SCIPvarGetLbLocal(var);
1477 }
1478 }
1479
1480 tmpvals.push_back(coef);
1481 tmpvars.push_back(var);
1482 }
1483
1484 /* if implied bounds are not sufficient to ensure bounds on aggrvar, then we cannot use the multiaggregation */
1485 if( j < lastRowCoefficients )
1486 break;
1487
1488 if( checklbimplied && SCIPisGT(scip, SCIPvarGetLbGlobal(aggrvar), impliedlb) )
1489 break;
1490
1491 if( checkubimplied && SCIPisLT(scip, SCIPvarGetUbGlobal(aggrvar), impliedub) )
1492 break;
1493
1494 SCIP_CALL( SCIPmultiaggregateVar(scip, aggrvar, tmpvars.size(),
1495 tmpvars.data(), tmpvals.data(), updatedSide / colCoef, &infeas, &aggregated) );
1496 }
1497
1498 if( aggregated )
1499 *naggrvars += 1;
1500 else if( constraintsReplaced && !redundant )
1501 {
1502 /* if the constraints where replaced, we need to add the failed substitution as an equality to SCIP */
1503 tmpvars.clear();
1504 tmpvals.clear();
1505 for( int j = startRowCoefficients; j < lastRowCoefficients; ++j )
1506 {
1507 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.indices[j]));
1508 tmpvals.push_back(res.postsolve.values[j]);
1509 }
1510
1511 SCIP_CONS* cons;
1512 String name = fmt::format("{}_failed_aggregation_equality", SCIPvarGetName(SCIPmatrixGetVar(matrix, col)));
1513 SCIP_CALL( SCIPcreateConsBasicLinear(scip, &cons, name.c_str(),
1514 tmpvars.size(), tmpvars.data(), tmpvals.data(), side, side ) );
1515 SCIP_CALL( SCIPaddCons(scip, cons) );
1516 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
1517 *naddconss += 1;
1518 }
1519
1520 if( infeas )
1521 {
1522 *result = SCIP_CUTOFF;
1523 break;
1524 }
1525
1526 break;
1527 }
1528 case ReductionType::kParallelCol:
1529 return SCIP_INVALIDRESULT;
1530#if PAPILO_VERSION_MAJOR > 1 || (PAPILO_VERSION_MAJOR == 1 && PAPILO_VERSION_MINOR >= 1)
1531 case ReductionType::kFixedInfCol: {
1532 /* todo: currently SCIP can not handle this kind of reduction (see issue #3391) */
1533 assert(false);
1534 if(!constraintsReplaced)
1535 continue;
1536 SCIP_Bool infeas;
1537 SCIP_Bool fixed;
1538 SCIP_Real value = SCIPinfinity(scip);
1539
1540 int column = res.postsolve.indices[first];
1541 bool is_negative_infinity = res.postsolve.values[first] < 0;
1542 SCIP_VAR* column_variable = SCIPmatrixGetVar(matrix, column);
1543
1544 if( is_negative_infinity )
1545 {
1546 value = -SCIPinfinity(scip);
1547 }
1548
1549 SCIP_CALL( SCIPfixVar(scip, column_variable, value, &infeas, &fixed) );
1550 *nfixedvars += 1;
1551
1552 assert(!infeas);
1553 assert(fixed);
1554 break;
1555 }
1556#endif
1557#if (PAPILO_VERSION_MAJOR >= 2)
1558 case ReductionType::kVarBoundChange :
1559 case ReductionType::kRedundantRow :
1560 case ReductionType::kRowBoundChange :
1561 case ReductionType::kReasonForRowBoundChangeForcedByRow :
1562 case ReductionType::kRowBoundChangeForcedByRow :
1563 case ReductionType::kSaveRow :
1564 case ReductionType::kReducedBoundsCost :
1565 case ReductionType::kColumnDualValue :
1566 case ReductionType::kRowDualValue :
1567 case ReductionType::kCoefficientChange :
1568 /* dual ReductionTypes should be only calculated for dual reductions and should not appear for MIP */
1569 SCIPerrorMessage("PaPILO: PaPILO should not return dual postsolving reductions in SCIP!!\n");
1570 SCIPABORT(); /*lint --e{527}*/
1571 break;
1572#endif
1573 default:
1574 SCIPdebugMsg(scip, "PaPILO returned unknown data type: \n" );
1575 continue;
1576 }
1577 }
1578
1579 /* finish with a final verb message and return */
1581 " (%.1fs) MILP presolver (%d rounds): %d aggregations, %d fixings, %d bound changes\n",
1582 SCIPgetSolvingTime(scip), presolve.getStatistics().nrounds, *naggrvars - oldnaggrvars,
1583 *nfixedvars - oldnfixedvars, *nchgbds - oldnchgbds);
1584
1585 /* free the matrix */
1586 assert(initialized);
1587 SCIPmatrixFree(scip, &matrix);
1588
1589 return SCIP_OKAY;
1590}
1591
1592
1593/*
1594 * Callback methods of presolver
1595 */
1596
1597/** copy method for constraint handler plugins (called when SCIP copies plugins) */
1598static
1599SCIP_DECL_PRESOLCOPY(presolCopyMILP)
1600{ /*lint --e{715}*/
1602
1603 return SCIP_OKAY;
1604}
1605
1606/** destructor of presolver to free user data (called when SCIP is exiting) */
1607static
1608SCIP_DECL_PRESOLFREE(presolFreeMILP)
1609{ /*lint --e{715}*/
1610 SCIP_PRESOLMILPDATA* data = (SCIP_PRESOLMILPDATA*)SCIPpresolGetData(presol);
1611 assert(data != NULL);
1612
1613 SCIPpresolSetData(presol, NULL);
1614 SCIPfreeBlockMemory(scip, &data);
1615 return SCIP_OKAY;
1616}
1617
1618/** initialization method of presolver (called after problem was transformed) */
1619static
1620SCIP_DECL_PRESOLINIT(presolInitMILP)
1621{ /*lint --e{715}*/
1622 SCIP_PRESOLMILPDATA* data = (SCIP_PRESOLMILPDATA*)SCIPpresolGetData(presol);
1623 assert(data != NULL);
1624
1625 data->lastncols = -1;
1626 data->lastnrows = -1;
1627
1628 return SCIP_OKAY;
1629}
1630
1631
1632/** execution method of presolver */
1633static
1634SCIP_DECL_PRESOLEXEC(presolExecMILP)
1635{ /*lint --e{715}*/
1636 SCIP_MATRIX* matrix;
1637 SCIP_PRESOLMILPDATA* data;
1638 SCIP_Bool initialized;
1639 SCIP_Bool complete;
1640 SCIP_Bool infeasible;
1641
1642 *result = SCIP_DIDNOTRUN;
1643
1644 data = (SCIP_PRESOLMILPDATA*)SCIPpresolGetData(presol);
1645
1646 int nvars = SCIPgetNVars(scip);
1647 int nconss = SCIPgetNConss(scip);
1648
1649 /* run only if the problem size reduced by some amount since the last call or if it is the first call */
1650 if( data->lastncols != -1 && data->lastnrows != -1 &&
1651 nvars > data->lastncols * 0.85 &&
1652 nconss > data->lastnrows * 0.85 )
1653 return SCIP_OKAY;
1654
1655 SCIP_CALL( SCIPmatrixCreate(scip, &matrix, TRUE, &initialized, &complete, &infeasible,
1656 naddconss, ndelconss, nchgcoefs, nchgbds, nfixedvars) );
1657
1658 /* if infeasibility was detected during matrix creation, return here */
1659 if( infeasible )
1660 {
1661 if( initialized )
1662 SCIPmatrixFree(scip, &matrix);
1663
1664 *result = SCIP_CUTOFF;
1665 return SCIP_OKAY;
1666 }
1667
1668 /* we only work on pure MIPs, also disable to try building the matrix again if it failed once */
1669 if( !initialized || !complete )
1670 {
1671 data->lastncols = 0;
1672 data->lastnrows = 0;
1673
1674 if( initialized )
1675 SCIPmatrixFree(scip, &matrix);
1676
1677 return SCIP_OKAY;
1678 }
1679
1680 if( 0 != strncmp(data->filename, DEFAULT_FILENAME_PROBLEM, strlen(DEFAULT_FILENAME_PROBLEM)) )
1681 {
1683 " writing transformed problem to %s (only enforced constraints)\n", data->filename);
1684 SCIP_CALL( SCIPwriteTransProblem(scip, data->filename, NULL, FALSE) );
1685 }
1686
1687 if( !SCIPisExact(scip) )
1688 return performRealPresolving(scip, matrix, data, initialized, nfixedvars, naggrvars, nchgvartypes, nchgbds,
1689 naddholes, ndelconss, naddconss, nupgdconss, nchgcoefs, nchgsides, result);
1690#if defined(PAPILO_WITH_EXACTPRESOLVE)
1691 else
1692 return performRationalPresolving(scip, matrix, data, initialized, nfixedvars, naggrvars, nchgvartypes, nchgbds,
1693 naddholes, ndelconss, naddconss, nupgdconss, nchgcoefs, nchgsides, result);
1694#endif
1695
1696 return SCIP_OKAY;
1697}
1698
1699
1700/*
1701 * presolver specific interface methods
1702 */
1703
1704/** creates the MILP presolver and includes it in SCIP */
1706 SCIP* scip /**< SCIP data structure */
1707 )
1708{
1709 SCIP_PRESOLMILPDATA* presoldata;
1710 SCIP_PRESOL* presol;
1711
1712#if defined(PAPILO_VERSION_TWEAK) && PAPILO_VERSION_TWEAK != 0
1713 String name = fmt::format("PaPILO {}.{}.{}.{}", PAPILO_VERSION_MAJOR, PAPILO_VERSION_MINOR, PAPILO_VERSION_PATCH, PAPILO_VERSION_TWEAK);
1714#else
1715 String name = fmt::format("PaPILO {}.{}.{}", PAPILO_VERSION_MAJOR, PAPILO_VERSION_MINOR, PAPILO_VERSION_PATCH);
1716#endif
1717
1718#if defined(PAPILO_GITHASH_AVAILABLE) && defined(PAPILO_TBB)
1719 String desc = fmt::format("parallel presolve for integer and linear optimization (github.com/scipopt/papilo) (built with TBB) [GitHash: {}]", PAPILO_GITHASH);
1720#elif !defined(PAPILO_GITHASH_AVAILABLE) && !defined(PAPILO_TBB)
1721 String desc("parallel presolve for integer and linear optimization (github.com/scipopt/papilo)");
1722#elif defined(PAPILO_GITHASH_AVAILABLE) && !defined(PAPILO_TBB)
1723 String desc = fmt::format("parallel presolve for integer and linear optimization (github.com/scipopt/papilo) [GitHash: {}]", PAPILO_GITHASH);
1724#elif !defined(PAPILO_GITHASH_AVAILABLE) && defined(PAPILO_TBB)
1725 String desc = fmt::format("parallel presolve for integer and linear optimization (github.com/scipopt/papilo) (built with TBB)");
1726#endif
1727
1728 /* add external code info for the presolve library */
1729 SCIP_CALL( SCIPincludeExternalCodeInformation(scip, name.c_str(), desc.c_str()) );
1730
1731 /* create MILP presolver data */
1732 presoldata = NULL;
1733 SCIP_CALL( SCIPallocBlockMemory(scip, &presoldata) );
1734 BMSclearMemory(presoldata);
1735
1736 presol = NULL;
1737
1738 /* include presolver */
1740 presolExecMILP,
1741 (SCIP_PRESOLDATA*)presoldata) );
1742
1743 assert(presol != NULL);
1744
1745 /* set non fundamental callbacks via setter functions */
1746 SCIP_CALL( SCIPsetPresolCopy(scip, presol, presolCopyMILP) );
1747 SCIP_CALL( SCIPsetPresolFree(scip, presol, presolFreeMILP) );
1748 SCIP_CALL( SCIPsetPresolInit(scip, presol, presolInitMILP) );
1749
1750#if defined(PAPILO_WITH_EXACTPRESOLVE)
1751 /* mark presolver as exact */
1752 SCIPpresolMarkExact(presol);
1753#endif
1754
1755 /* add MILP presolver parameters */
1756#ifdef PAPILO_TBB
1758 "presolving/" PRESOL_NAME "/threads",
1759 "maximum number of threads presolving may use (0: automatic)",
1760 &presoldata->threads, FALSE, DEFAULT_THREADS, 0, INT_MAX, NULL, NULL) );
1761#else
1762 presoldata->threads = 1;
1763#endif
1764
1766 "presolving/" PRESOL_NAME "/maxfillinpersubstitution",
1767 "maximal possible fillin for substitutions to be considered",
1768 &presoldata->maxfillinpersubstitution, FALSE, DEFAULT_MAXFILLINPERSUBST, INT_MIN, INT_MAX, NULL, NULL) );
1769
1771 "presolving/" PRESOL_NAME "/maxshiftperrow",
1772 "maximal amount of nonzeros allowed to be shifted to make space for substitutions",
1773 &presoldata->maxshiftperrow, TRUE, DEFAULT_MAXSHIFTPERROW, 0, INT_MAX, NULL, NULL) );
1774
1776 "presolving/" PRESOL_NAME "/randomseed",
1777 "the random seed used for randomization of tie breaking",
1778 &presoldata->randomseed, FALSE, DEFAULT_RANDOMSEED, INT_MIN, INT_MAX, NULL, NULL) );
1779
1780 if( DependentRows<double>::Enabled )
1781 {
1783 "presolving/" PRESOL_NAME "/detectlineardependency",
1784 "should linear dependent equations and free columns be removed? (0: never, 1: for LPs, 2: always)",
1785 &presoldata->detectlineardependency, TRUE, DEFAULT_DETECTLINDEP, 0, 2, NULL, NULL) );
1786 }
1787 else
1788 presoldata->detectlineardependency = DEFAULT_DETECTLINDEP;
1789
1791 "presolving/" PRESOL_NAME "/modifyconsfac",
1792 "modify SCIP constraints when the number of nonzeros or rows is at most this factor "
1793 "times the number of nonzeros or rows before presolving",
1794 &presoldata->modifyconsfac, FALSE, DEFAULT_MODIFYCONSFAC, 0.0, 1.0, NULL, NULL) );
1795
1797 "presolving/" PRESOL_NAME "/markowitztolerance",
1798 "the markowitz tolerance used for substitutions",
1799 &presoldata->markowitztolerance, FALSE, DEFAULT_MARKOWITZTOLERANCE, 0.0, 1.0, NULL, NULL) );
1800
1802 "presolving/" PRESOL_NAME "/hugebound",
1803 "absolute bound value that is considered too huge for activity based calculations",
1804 &presoldata->hugebound, FALSE, DEFAULT_HUGEBOUND, 0.0, SCIP_REAL_MAX, NULL, NULL) );
1805
1806#if PAPILO_APIVERSION >= 2
1807 SCIP_CALL( SCIPaddRealParam(scip, "presolving/" PRESOL_NAME "/abortfacexhaustive",
1808 "abort threshold for exhaustive presolving in PAPILO",
1809 &presoldata->abortfacexhaustive, TRUE, DEFAULT_ABORTFAC_EXHAUSTIVE, 0.0, 1.0, NULL, NULL) );
1810 SCIP_CALL( SCIPaddRealParam(scip, "presolving/" PRESOL_NAME "/abortfacmedium",
1811 "abort threshold for medium presolving in PAPILO",
1812 &presoldata->abortfacmedium, TRUE, DEFAULT_ABORTFAC_MEDIUM, 0.0, 1.0, NULL, NULL) );
1813 SCIP_CALL( SCIPaddRealParam(scip, "presolving/" PRESOL_NAME "/abortfacfast",
1814 "abort threshold for fast presolving in PAPILO",
1815 &presoldata->abortfacfast, TRUE, DEFAULT_ABORTFAC_FAST, 0.0, 1.0, NULL, NULL) );
1816#else
1817 presoldata->abortfacexhaustive = DEFAULT_ABORTFAC_EXHAUSTIVE;
1818 presoldata->abortfacmedium = DEFAULT_ABORTFAC_MEDIUM;
1819 presoldata->abortfacfast = DEFAULT_ABORTFAC_FAST;
1820#endif
1821
1822#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 1)
1823 SCIP_CALL( SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/maxbadgesizeseq",
1824 "maximal badge size in Probing in PaPILO if PaPILO is executed in sequential mode",
1825 &presoldata->maxbadgesizeseq, FALSE, DEFAULT_MAXBADGESIZE_SEQ, -1, INT_MAX, NULL, NULL) );
1826
1827 SCIP_CALL( SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/maxbadgesizepar",
1828 "maximal badge size in Probing in PaPILO if PaPILO is executed in parallel mode",
1829 &presoldata->maxbadgesizepar, FALSE, DEFAULT_MAXBADGESIZE_PAR, -1, INT_MAX, NULL, NULL) );
1830#else
1831 presoldata->maxbadgesizeseq = DEFAULT_MAXBADGESIZE_SEQ;
1832 presoldata->maxbadgesizepar = DEFAULT_MAXBADGESIZE_PAR;
1833#endif
1834
1835#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 3)
1836 SCIP_CALL( SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/internalmaxrounds",
1837 "internal maxrounds for each milp presolving (-1: no limit, 0: model cleanup)",
1838 &presoldata->internalmaxrounds, TRUE, DEFAULT_INTERNAL_MAXROUNDS, -1, INT_MAX, NULL, NULL) );
1839#else
1840 presoldata->internalmaxrounds = DEFAULT_INTERNAL_MAXROUNDS;
1841#endif
1842
1844 "presolving/" PRESOL_NAME "/enableparallelrows",
1845 "should the parallel rows presolver be enabled within the presolve library?",
1846 &presoldata->enableparallelrows, TRUE, DEFAULT_ENABLEPARALLELROWS, NULL, NULL) );
1847
1849 "presolving/" PRESOL_NAME "/enabledomcol",
1850 "should the dominated column presolver be enabled within the presolve library?",
1851 &presoldata->enabledomcol, TRUE, DEFAULT_ENABLEDOMCOL, NULL, NULL) );
1852
1854 "presolving/" PRESOL_NAME "/enabledualinfer",
1855 "should the dualinfer presolver be enabled within the presolve library?",
1856 &presoldata->enabledualinfer, TRUE, DEFAULT_ENABLEDUALINFER, NULL, NULL) );
1857
1859 "presolving/" PRESOL_NAME "/enablemultiaggr",
1860 "should the multi-aggregation presolver be enabled within the presolve library?",
1861 &presoldata->enablemultiaggr, TRUE, DEFAULT_ENABLEMULTIAGGR, NULL, NULL) );
1862
1864 "presolving/" PRESOL_NAME "/enableprobing",
1865 "should the probing presolver be enabled within the presolve library?",
1866 &presoldata->enableprobing, TRUE, DEFAULT_ENABLEPROBING, NULL, NULL) );
1867
1869 "presolving/" PRESOL_NAME "/enablesparsify",
1870 "should the sparsify presolver be enabled within the presolve library?",
1871 &presoldata->enablesparsify, TRUE, DEFAULT_ENABLESPARSIFY, NULL, NULL) );
1872
1873 SCIP_CALL( SCIPaddStringParam(scip, "presolving/" PRESOL_NAME "/probfilename",
1874 "filename to store the problem before MILP presolving starts (only enforced constraints)",
1875 &presoldata->filename, TRUE, DEFAULT_FILENAME_PROBLEM, NULL, NULL) );
1876
1877 SCIP_CALL( SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/verbosity",
1878 "verbosity level of PaPILO (0: quiet, 1: errors, 2: warnings, 3: normal, 4: detailed)",
1879 &presoldata->verbosity, FALSE, DEFAULT_VERBOSITY, 0, 4, NULL, NULL) );
1880#if PAPILO_APIVERSION >= 6
1882 "presolving/" PRESOL_NAME "/enablecliquemerging",
1883 "should the clique merging presolver be enabled within the presolve library?",
1884 &presoldata->enablecliquemerging, TRUE, DEFAULT_ENABLESPARSIFY, NULL, NULL) );
1886 "presolving/" PRESOL_NAME "/maxedgesparallel",
1887 "maximal amount of edges in the parallel clique merging graph",
1888 &presoldata->maxedgesparallel, FALSE, DEFAULT_MAXEDGESPARALLEL, -1, INT_MAX, NULL, NULL) );
1890 "presolving/" PRESOL_NAME "/maxedgessequential",
1891 "maximal amount of edges in the sequential clique merging graph",
1892 &presoldata->maxedgessequential, FALSE, DEFAULT_MAXEDGESSEQUENTIAL, -1, INT_MAX, NULL, NULL) );
1894 "presolving/" PRESOL_NAME "/maxcliquesize",
1895 "maximal size of clique considered for clique merging",
1896 &presoldata->maxcliquesize, FALSE, DEFAULT_MAXCLIQUESIZE, -1, INT_MAX, NULL, NULL) );
1898 "presolving/" PRESOL_NAME "/maxgreedycalls",
1899 "maximal number of greedy max clique calls in a single thread",
1900 &presoldata->maxgreedycalls, FALSE, DEFAULT_MAXGREEDYCALLS, -1, INT_MAX, NULL, NULL) );
1901#endif
1902
1903 return SCIP_OKAY;
1904}
1905
1906#endif
Constraint handler for linear constraints in their most general form, .
Constraint handler for linear constraints in their most general form, .
#define NULL
Definition: def.h:248
#define SCIP_REAL_MAX
Definition: def.h:158
#define SCIP_Bool
Definition: def.h:91
#define SCIP_Real
Definition: def.h:156
#define TRUE
Definition: def.h:93
#define FALSE
Definition: def.h:94
#define SCIPABORT()
Definition: def.h:327
#define REALABS(x)
Definition: def.h:182
#define SCIP_CALL(x)
Definition: def.h:355
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)
SCIP_RETCODE SCIPcreateConsBasicExactLinear(SCIP *scip, SCIP_CONS **cons, const char *name, int nvars, SCIP_VAR **vars, SCIP_RATIONAL **vals, SCIP_RATIONAL *lhs, SCIP_RATIONAL *rhs)
const char * SCIPgetProbName(SCIP *scip)
Definition: scip_prob.c:1242
int SCIPgetNVars(SCIP *scip)
Definition: scip_prob.c:2246
SCIP_RETCODE SCIPaddCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip_prob.c:3274
SCIP_RETCODE SCIPdelCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip_prob.c:3420
int SCIPgetNConss(SCIP *scip)
Definition: scip_prob.c:3620
SCIP_RETCODE SCIPwriteTransProblem(SCIP *scip, const char *filename, const char *extension, SCIP_Bool genericnames)
Definition: scip_prob.c:789
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:4798
SCIP_CONSHDLR * SCIPfindConshdlr(SCIP *scip, const char *name)
Definition: scip_cons.c:940
SCIP_RETCODE SCIPprintCons(SCIP *scip, SCIP_CONS *cons, FILE *file)
Definition: scip_cons.c:2536
const char * SCIPconsGetName(SCIP_CONS *cons)
Definition: cons.c:8389
SCIP_RETCODE SCIPreleaseCons(SCIP *scip, SCIP_CONS **cons)
Definition: scip_cons.c:1173
SCIP_RETCODE SCIPcaptureCons(SCIP *scip, SCIP_CONS *cons)
Definition: scip_cons.c:1138
SCIP_Bool SCIPisExact(SCIP *scip)
Definition: scip_exact.c:193
SCIP_RETCODE SCIPincludeExternalCodeInformation(SCIP *scip, const char *name, const char *description)
Definition: scip_general.c:769
BMS_BUFMEM * SCIPbuffer(SCIP *scip)
Definition: scip_mem.c:72
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip_mem.h:108
#define SCIPallocBlockMemory(scip, ptr)
Definition: scip_mem.h:89
void SCIPpresolMarkExact(SCIP_PRESOL *presol)
Definition: presol.c:615
void SCIPpresolSetData(SCIP_PRESOL *presol, SCIP_PRESOLDATA *presoldata)
Definition: presol.c:538
SCIP_PRESOLDATA * SCIPpresolGetData(SCIP_PRESOL *presol)
Definition: presol.c:528
SCIP_RETCODE SCIPsetPresolInit(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLINIT((*presolinit)))
Definition: scip_presol.c:180
SCIP_RETCODE SCIPsetPresolFree(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLFREE((*presolfree)))
Definition: scip_presol.c:164
SCIP_RETCODE SCIPsetPresolCopy(SCIP *scip, SCIP_PRESOL *presol, SCIP_DECL_PRESOLCOPY((*presolcopy)))
Definition: scip_presol.c:148
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:113
void SCIPrationalSetInfinity(SCIP_RATIONAL *res)
Definition: rational.cpp:618
SCIP_Real SCIPrationalGetReal(SCIP_RATIONAL *rational)
Definition: rational.cpp:2085
#define SCIPrationalDebugMessage
Definition: rational.h:641
void SCIPrationalDiv(SCIP_RATIONAL *res, SCIP_RATIONAL *op1, SCIP_RATIONAL *op2)
Definition: rational.cpp:1132
SCIP_Bool SCIPrationalIsAbsInfinity(SCIP_RATIONAL *rational)
Definition: rational.cpp:1680
void SCIPrationalFreeBuffer(BMS_BUFMEM *bufmem, SCIP_RATIONAL **rational)
Definition: rational.cpp:473
void SCIPrationalDiff(SCIP_RATIONAL *res, SCIP_RATIONAL *op1, SCIP_RATIONAL *op2)
Definition: rational.cpp:983
SCIP_RETCODE SCIPrationalCreateBuffer(BMS_BUFMEM *bufmem, SCIP_RATIONAL **rational)
Definition: rational.cpp:123
SCIP_Bool SCIPrationalIsZero(SCIP_RATIONAL *rational)
Definition: rational.cpp:1624
void SCIPrationalSetNegInfinity(SCIP_RATIONAL *res)
Definition: rational.cpp:630
SCIP_Bool SCIPrationalIsNegative(SCIP_RATIONAL *rational)
Definition: rational.cpp:1650
SCIP_Bool SCIPrationalIsInfinity(SCIP_RATIONAL *rational)
Definition: rational.cpp:1660
SCIP_RETCODE SCIPrationalCreateBufferArray(BMS_BUFMEM *mem, SCIP_RATIONAL ***rational, int size)
Definition: rational.cpp:214
SCIP_Bool SCIPrationalIsNegInfinity(SCIP_RATIONAL *rational)
Definition: rational.cpp:1670
void SCIPrationalMultReal(SCIP_RATIONAL *res, SCIP_RATIONAL *op1, SCIP_Real op2)
Definition: rational.cpp:1097
void SCIPrationalFreeBufferArray(BMS_BUFMEM *mem, SCIP_RATIONAL ***ratbufarray, int size)
Definition: rational.cpp:518
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:6401
SCIP_RETCODE SCIPtightenVarUbExact(SCIP *scip, SCIP_VAR *var, SCIP_RATIONAL *newbound, SCIP_Bool *infeasible, SCIP_Bool *tightened)
Definition: scip_var.c:6768
SCIP_RATIONAL * SCIPvarGetAggrScalarExact(SCIP_VAR *var)
Definition: var.c:23760
SCIP_VARSTATUS SCIPvarGetStatus(SCIP_VAR *var)
Definition: var.c:23386
SCIP_Bool SCIPvarIsImpliedIntegral(SCIP_VAR *var)
Definition: var.c:23498
SCIP_Real SCIPvarGetUbLocal(SCIP_VAR *var)
Definition: var.c:24268
SCIP_RETCODE SCIPaggregateVarsExact(SCIP *scip, SCIP_VAR *varx, SCIP_VAR *vary, SCIP_RATIONAL *scalarx, SCIP_RATIONAL *scalary, SCIP_RATIONAL *rhs, SCIP_Bool *infeasible, SCIP_Bool *redundant, SCIP_Bool *aggregated)
Definition: scip_var.c:10692
SCIP_RATIONAL * SCIPvarGetAggrConstantExact(SCIP_VAR *var)
Definition: var.c:23783
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:10550
SCIP_Real SCIPvarGetObj(SCIP_VAR *var)
Definition: var.c:23900
SCIP_RETCODE SCIPtightenVarUb(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound, SCIP_Bool force, SCIP_Bool *infeasible, SCIP_Bool *tightened)
Definition: scip_var.c:6651
SCIP_VARTYPE SCIPvarGetType(SCIP_VAR *var)
Definition: var.c:23453
SCIP_RETCODE SCIPgetProbvarSum(SCIP *scip, SCIP_VAR **var, SCIP_Real *scalar, SCIP_Real *constant)
Definition: scip_var.c:2499
SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:24142
SCIP_VARSTATUS SCIPvarGetStatusExact(SCIP_VAR *var)
Definition: var.c:23396
const char * SCIPvarGetName(SCIP_VAR *var)
Definition: var.c:23267
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:10834
SCIP_Bool SCIPvarIsIntegral(SCIP_VAR *var)
Definition: var.c:23490
SCIP_Real SCIPvarGetLbLocal(SCIP_VAR *var)
Definition: var.c:24234
SCIP_RATIONAL * SCIPvarGetLbGlobalExact(SCIP_VAR *var)
Definition: var.c:24130
SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:24120
SCIP_RETCODE SCIPfixVar(SCIP *scip, SCIP_VAR *var, SCIP_Real fixedval, SCIP_Bool *infeasible, SCIP_Bool *fixed)
Definition: scip_var.c:10318
SCIP_RETCODE SCIPgetProbvarSumExact(SCIP *scip, SCIP_VAR **var, SCIP_RATIONAL *scalar, SCIP_RATIONAL *constant)
Definition: scip_var.c:2538
SCIP_RATIONAL * SCIPvarGetObjExact(SCIP_VAR *var)
Definition: var.c:23910
SCIP_Bool SCIPallowWeakDualReds(SCIP *scip)
Definition: scip_var.c:10998
SCIP_RETCODE SCIPmultiaggregateVarExact(SCIP *scip, SCIP_VAR *var, int naggvars, SCIP_VAR **aggvars, SCIP_RATIONAL **scalars, SCIP_RATIONAL *constant, SCIP_Bool *infeasible, SCIP_Bool *aggregated)
Definition: scip_var.c:10879
SCIP_RETCODE SCIPtightenVarLbExact(SCIP *scip, SCIP_VAR *var, SCIP_RATIONAL *newbound, SCIP_Bool *infeasible, SCIP_Bool *tightened)
Definition: scip_var.c:6518
SCIP_Bool SCIPallowStrongDualReds(SCIP *scip)
Definition: scip_var.c:10984
SCIP_RATIONAL * SCIPvarGetUbGlobalExact(SCIP_VAR *var)
Definition: var.c:24152
SCIP_RETCODE SCIPfixVarExact(SCIP *scip, SCIP_VAR *var, SCIP_RATIONAL *fixedval, SCIP_Bool *infeasible, SCIP_Bool *fixed)
Definition: scip_var.c:10420
SCIP_VAR * SCIPvarGetAggrVar(SCIP_VAR *var)
Definition: var.c:23736
unsigned int SCIPinitializeRandomSeed(SCIP *scip, unsigned int initialseedvalue)
int SCIPmatrixGetNNonzs(SCIP_MATRIX *matrix)
Definition: matrix.c:2107
SCIP_RATIONAL * SCIPmatrixGetRowLhsExact(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2071
int SCIPmatrixGetRowNNonzs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2013
SCIP_Real SCIPmatrixGetRowLhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2047
SCIP_Real * SCIPmatrixGetRowValPtr(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1977
SCIP_RATIONAL ** SCIPmatrixGetRowValPtrExact(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1989
SCIP_Real SCIPmatrixGetRowRhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2059
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:703
int SCIPmatrixGetNColumns(SCIP_MATRIX *matrix)
Definition: matrix.c:1897
SCIP_CONS * SCIPmatrixGetCons(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2189
void SCIPmatrixFree(SCIP *scip, SCIP_MATRIX **matrix)
Definition: matrix.c:1348
SCIP_VAR * SCIPmatrixGetVar(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1953
SCIP_RATIONAL * SCIPmatrixGetRowRhsExact(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2083
int * SCIPmatrixGetRowIdxPtr(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2001
int SCIPmatrixGetNRows(SCIP_MATRIX *matrix)
Definition: matrix.c:2037
#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
#define SCIPdebug(x)
Definition: pub_message.h:93
#define SCIPdebugMessage
Definition: pub_message.h:96
public methods for presolvers
public methods for problem variables
wrapper for rational number arithmetic
#define DEFAULT_RANDOMSEED
Definition: reader_scflp.c:116
public methods for constraint handler plugins and constraints
public methods for exact solving
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:1138
public methods for timing
public methods for SCIP variables
scip::Rational val
unsigned int isfprepresentable
definition of wrapper class for rational numbers
@ SCIP_VERBLEVEL_HIGH
Definition: type_message.h:61
#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_ISFPREPRESENTABLE_UNKNOWN
Definition: type_rational.h:48
@ 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
enum SCIP_Result SCIP_RESULT
Definition: type_result.h:61
@ 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_CONTINUOUS
Definition: type_var.h:71
@ SCIP_VARSTATUS_FIXED
Definition: type_var.h:54
@ SCIP_VARSTATUS_MULTAGGR
Definition: type_var.h:56
@ SCIP_VARSTATUS_NEGATED
Definition: type_var.h:57
@ SCIP_VARSTATUS_AGGREGATED
Definition: type_var.h:55