SCIP Doxygen Documentation
Loading...
Searching...
No Matches
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-2026 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
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#if PAPILO_APIVERSION >= 12
448 // TODO: enable this after performance test. On MIPLIB this brought 3% on instances with cliques.
449 probing->set_numcliquefails(0);
450#endif
451 presolve.addPresolveMethod( uptr( probing ) );
452#else
453 presolve.addPresolveMethod( uptr( new Probing<T>() ) );
454 if( data->maxbadgesizeseq != DEFAULT_MAXBADGESIZE_SEQ )
456 " The parameter 'presolving/milp/maxbadgesizeseq' can only be used with PaPILO 2.1.0 or later versions.\n");
457
458 if( data->maxbadgesizepar != DEFAULT_MAXBADGESIZE_PAR )
460 " The parameter 'presolving/milp/maxbadgesizepar' can only be used with PaPILO 2.1.0 or later versions.\n");
461#endif
462 }
463 if( data->enabledomcol )
464 presolve.addPresolveMethod( uptr( new DominatedCols<T>() ) );
465 if( data->enablemultiaggr )
466 presolve.addPresolveMethod( uptr( new Substitution<T>() ) );
467 if( data->enablesparsify )
468 presolve.addPresolveMethod( uptr( new Sparsify<T>() ) );
469
470 /* set numerical tolerances */
471#if PAPILO_APIVERSION >= 3
472 presolve.getPresolveOptions().useabsfeas = false;
473#endif
474 if( SCIPisExact(scip) )
475 {
476 presolve.getPresolveOptions().epsilon = 0.0;
477 presolve.getPresolveOptions().feastol = 0.0;
478 }
479 else
480 {
481 presolve.getPresolveOptions().epsilon = SCIPepsilon(scip);
482 presolve.getPresolveOptions().feastol = SCIPfeastol(scip);
483 }
484
485#ifndef SCIP_PRESOLLIB_ENABLE_OUTPUT
486 /* adjust output settings of presolve library */
487 presolve.setVerbosityLevel((VerbosityLevel) data->verbosity);
488#endif
489
490#if PAPILO_APIVERSION >= 2
491 presolve.getPresolveOptions().abortfac = data->abortfacexhaustive;
492 presolve.getPresolveOptions().abortfacmedium = data->abortfacmedium;
493 presolve.getPresolveOptions().abortfacfast = data->abortfacfast;
494#endif
495
496 /* communicate the time limit */
497 SCIPgetRealParam(scip, "limits/time", &timelimit);
498 if( !SCIPisInfinity(scip, timelimit) )
499 presolve.getPresolveOptions().tlim = timelimit - SCIPgetSolvingTime(scip);
500
501 return SCIP_OKAY;
502}
503
504
505#if defined(PAPILO_WITH_EXACTPRESOLVE)
506/** calls PaPILO presolving in rational arithmetic*/
507static
508SCIP_RETCODE performRationalPresolving(
509 SCIP* scip, /**< SCIP data structure */
510 SCIP_MATRIX* matrix, /**< initialized SCIP_MATRIX data structure */
511 SCIP_PRESOLMILPDATA* data, /**< plugin specific presol data */
512 SCIP_Bool initialized, /**< was the matrix initialized */
513 int* nfixedvars, /**< store number of fixed variables */
514 int* naggrvars, /**< store number of aggregated variables */
515 int* nchgvartypes, /**< store number of changed variable types */
516 int* nchgbds, /**< store number of changed bounds */
517 int* naddholes, /**< store number of added holes */
518 int* ndelconss, /**< store the number of deleted cons */
519 int* naddconss, /**< store the number of added cons */
520 int* nupgdconss, /**< store the number of upgraded cons */
521 int* nchgcoefs, /**< store the number of changed coefficients */
522 int* nchgsides, /**< store the number of changed sides */
523 SCIP_RESULT* result /**< result pointer */
524 )
525{
526 int nvars = SCIPgetNVars(scip);
527 int nconss = SCIPgetNConss(scip);
528
529 SCIP_CONSHDLR* linconshdlr = SCIPfindConshdlr(scip, "exactlinear");
530 assert(linconshdlr != NULL);
531 bool allowconsmodification = (SCIPconshdlrGetNCheckConss(linconshdlr) == SCIPmatrixGetNRows(matrix));
532
533 /* store current numbers of aggregations, fixings, and changed bounds for statistics */
534 int oldnaggrvars = *naggrvars;
535 int oldnfixedvars = *nfixedvars;
536 int oldnchgbds = *nchgbds;
537
538 Problem<papilo::Rational> problem = buildProblemRational(scip, matrix);
539 Presolve<papilo::Rational> presolve;
540 setupPresolve(scip, presolve, data, allowconsmodification);
541
542 /* call presolving (without storing information for dual postsolve) */
544 " (%.1fs) running MILP presolver%s\n", SCIPgetSolvingTime(scip),
545 presolve.getPresolveOptions().threads == 1 ? "" : " on multiple threads");
546
547 int oldnnz = problem.getConstraintMatrix().getNnz();
548
549#if (PAPILO_VERSION_MAJOR >= 2)
550 PresolveResult<papilo::Rational> res = presolve.apply(problem, false);
551#else
552 PresolveResult<papilo::Rational> res = presolve.apply(problem);
553#endif
554 data->lastncols = problem.getNCols();
555 data->lastnrows = problem.getNRows();
556
557 /* evaluate the result */
558 switch( res.status )
559 {
560 case PresolveStatus::kInfeasible:
563 " (%.1fs) MILP presolver detected infeasibility\n",
565 SCIPmatrixFree(scip, &matrix);
566 return SCIP_OKAY;
567 case PresolveStatus::kUnbndOrInfeas:
568 case PresolveStatus::kUnbounded:
571 " (%.1fs) MILP presolver detected unboundedness\n",
573 SCIPmatrixFree(scip, &matrix);
574 return SCIP_OKAY;
575 case PresolveStatus::kUnchanged:
577 data->lastncols = nvars;
578 data->lastnrows = nconss;
580 " (%.1fs) MILP presolver found nothing\n",
582 SCIPmatrixFree(scip, &matrix);
583 return SCIP_OKAY;
584 case PresolveStatus::kReduced:
585 data->lastncols = problem.getNCols();
586 data->lastnrows = problem.getNRows();
588 }
589
590 /* result indicated success, now populate the changes into the SCIP structures */
591 Vec<SCIP_VAR*> tmpvars;
592 Vec<SCIP_Real> tmpvalsreal;
593
594 /* if the number of nonzeros decreased by a sufficient factor, rather create all constraints from scratch */
595 int newnnz = problem.getConstraintMatrix().getNnz();
596 bool constraintsReplaced = false;
597 if( newnnz == 0 || (allowconsmodification &&
598 (problem.getNRows() <= data->modifyconsfac * data->lastnrows ||
599 newnnz <= data->modifyconsfac * oldnnz)) )
600 {
601 int oldnrows = SCIPmatrixGetNRows(matrix);
602 int newnrows = problem.getNRows();
603
604 constraintsReplaced = true;
605
606 /* capture constraints that are still present in the problem after presolve */
607 for( int i = 0; i < newnrows; ++i )
608 {
609 SCIP_CONS* c = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
611 }
612
613 /* delete all constraints */
614 *ndelconss += oldnrows;
615 *naddconss += newnrows;
616
617 for( int i = 0; i < oldnrows; ++i )
618 {
620 }
621
622 /* now loop over rows of presolved problem and create them as new linear constraints,
623 * then release the old constraint after its name was passed to the new constraint
624 */
625 const Vec<RowFlags>& rflags = problem.getRowFlags();
626 const auto& consmatrix = problem.getConstraintMatrix();
627 for( int i = 0; i < newnrows; ++i )
628 {
629 auto rowvec = consmatrix.getRowCoefficients(i);
630 const int* rowcols = rowvec.getIndices();
631 /* SCIPcreateConsBasicLinear() requires a non const pointer */
632 papilo::Rational* rowvals = const_cast<papilo::Rational*>(rowvec.getValues());
633 int rowlen = rowvec.getLength();
634
635 /* retrieve SCIP compatible left and right hand sides */
636 papilo::Rational lhs = rflags[i].test(RowFlag::kLhsInf) ? - SCIPinfinity(scip) : consmatrix.getLeftHandSides()[i];
637 papilo::Rational rhs = rflags[i].test(RowFlag::kRhsInf) ? SCIPinfinity(scip) : consmatrix.getRightHandSides()[i];
638
639 /* create variable array matching the value array */
640 tmpvars.clear();
641 tmpvars.reserve(rowlen);
642 for( int j = 0; j < rowlen; ++j )
643 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[rowcols[j]]));
644
645 /* create and add new constraint with name of old constraint */
646 SCIP_CONS* oldcons = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
647 SCIP_CONS* cons;
648 SCIP_RATIONAL** tmpvals;
649 SCIP_RATIONAL* tmplhs;
650 SCIP_RATIONAL* tmprhs;
651
655
656 for( int j = 0; j < rowlen; j++ )
657 setRational(scip, tmpvals[j], rowvals[j]);
658
659 setRational(scip, tmprhs, rhs);
660 setRational(scip, tmplhs, lhs);
661 if( rflags[i].test(RowFlag::kLhsInf) )
663 if( rflags[i].test(RowFlag::kRhsInf) )
665
666 SCIP_CALL( SCIPcreateConsBasicExactLinear(scip, &cons, SCIPconsGetName(oldcons), rowlen, tmpvars.data(), tmpvals, tmplhs, tmprhs) );
667 SCIP_CALL( SCIPaddCons(scip, cons) );
668
669 /* release old and new constraint */
670 SCIP_CALL( SCIPreleaseCons(scip, &oldcons) );
671 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
672
675 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &tmpvals, rowlen);
676 }
677 }
678
679 /* loop over res.postsolve and add all fixed variables and aggregations to scip */
680 for( std::size_t i = 0; i != res.postsolve.types.size(); ++i )
681 {
682 ReductionType type =
683 res.postsolve.types[i];
684 int first = res.postsolve.start[i];
685 int last = res.postsolve.start[i + 1];
686
687 switch( type )
688 {
689 case ReductionType::kFixedCol:
690 {
691 SCIP_RATIONAL* tmpval;
692 SCIP_Bool infeas;
693 SCIP_Bool fixed;
694 int col = res.postsolve.indices[first];
695
696 SCIP_VAR* var = SCIPmatrixGetVar(matrix, col);
697
698 papilo::Rational value = res.postsolve.values[first];
699
701 setRational(scip, tmpval, value);
702
703 SCIPrationalDebugMessage("Papilo fix var %s to %q \n", SCIPvarGetName(var), tmpval);
704
705 /* SCIP has different rules for aggregation than PaPILO
706 * As a result, SCIP might have aggregated and replaced the variable that PaPILO now wants to fix
707 */
709 {
710 SCIP_RATIONAL* aggregatedScalar;
711 SCIP_RATIONAL* aggregatedConst;
712
713 aggregatedScalar = SCIPvarGetAggrScalarExact(var);
714 aggregatedConst = SCIPvarGetAggrConstantExact(var);
715
716 /* fix aggregation variable y in x = a*y + c, instead of fixing x directly */
718 assert( !SCIPrationalIsZero(aggregatedScalar));
719 if( SCIPrationalIsAbsInfinity(tmpval) )
720 SCIPrationalMultReal(tmpval, tmpval, SCIPrationalIsNegative(aggregatedScalar) ? -1 : 1);
721 else
722 {
723 SCIPrationalDiff(tmpval, tmpval, aggregatedConst);
724 SCIPrationalDiv(tmpval, tmpval, aggregatedScalar);
725 }
727 }
728
729 /* SCIP might also have fixed the variable during aggregation */
731 break;
732
733 SCIP_CALL( SCIPfixVarExact(scip, var, tmpval, &infeas, &fixed) );
734
735 *nfixedvars += 1;
736
738
739 assert(!infeas);
742 break;
743 }
744 /*
745 * Dual-postsolving in PaPILO required introducing a postsolve-type for substitution with additional information.
746 * Further, the different Substitution-postsolving types store the required postsolving data differently (in different order) in the postsolving stack.
747 * Therefore, we need to distinguish how to parse the required data (rowLength, col, side, startRowCoefficients, lastRowCoefficients) from the postsolving stack.
748 * If these values are accessed, the procedure is the same for both.
749 */
750#if (PAPILO_VERSION_MAJOR >= 2)
751 case ReductionType::kSubstitutedColWithDual:
752#endif
753 case ReductionType::kSubstitutedCol:
754 {
755 int col = 0;
756 papilo::Rational side = 0;
757
758 int rowlen = 0;
759 int startRowCoefficients = 0;
760 int lastRowCoefficients = 0;
761
762 if( type == ReductionType::kSubstitutedCol )
763 {
764 rowlen = last - first - 1;
765 col = res.postsolve.indices[first];
766 side = res.postsolve.values[first];
767
768 startRowCoefficients = first + 1;
769 lastRowCoefficients = last;
770 }
771#if (PAPILO_VERSION_MAJOR >= 2)
772 if( type == ReductionType::kSubstitutedColWithDual )
773 {
774 rowlen = (int) res.postsolve.values[first];
775 col = res.postsolve.indices[first + 3 + rowlen];
776 side = res.postsolve.values[first + 1];
777
778 startRowCoefficients = first + 3;
779 lastRowCoefficients = first + 3 + rowlen;
780
781 assert(side == res.postsolve.values[first + 2]);
782 assert(res.postsolve.indices[first + 1] == 0);
783 assert(res.postsolve.indices[first + 2] == 0);
784 }
785 assert( type == ReductionType::kSubstitutedCol || type == ReductionType::kSubstitutedColWithDual );
786#else
787 assert( type == ReductionType::kSubstitutedCol );
788#endif
789 SCIP_Bool infeas;
790 SCIP_Bool aggregated;
791 SCIP_Bool redundant = FALSE;
792 SCIP_RATIONAL* constant;
793 if( rowlen == 2 )
794 {
795 SCIP_VAR* varx = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients]);
796 SCIP_VAR* vary = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients + 1]);
797 papilo::Rational scalarx = res.postsolve.values[startRowCoefficients];
798 papilo::Rational scalary = res.postsolve.values[startRowCoefficients + 1];
799
800 SCIP_RATIONAL* tmpscalarx;
801 SCIP_RATIONAL* tmpscalary;
802 SCIP_RATIONAL* tmpside;
807
808 setRational(scip, tmpscalarx, scalarx);
809 setRational(scip, tmpscalary, scalary);
810
811 SCIP_CALL( SCIPgetProbvarSumExact(scip, &varx, tmpscalarx, constant) );
813
814 SCIP_CALL( SCIPgetProbvarSumExact(scip, &vary, tmpscalary, constant) );
816
817 setRational(scip, tmpside, side);
818 SCIPrationalDiff(tmpside, tmpside, constant);
819
820 SCIPrationalDebugMessage("Papilo aggregate vars %s, %s with scalars %q, %q and constant %q \n", SCIPvarGetName(varx), SCIPvarGetName(vary),
821 tmpscalarx, tmpscalary, constant);
822
823 SCIP_CALL( SCIPaggregateVarsExact(scip, varx, vary, tmpscalarx, tmpscalary, tmpside, &infeas, &redundant, &aggregated) );
824
829 }
830 else
831 {
832 SCIP_RATIONAL* colCoef;
833 SCIP_RATIONAL* updatedSide;
834 SCIP_RATIONAL** tmpvals;
835 int c = 0;
836
841
842 for( int j = startRowCoefficients; j < lastRowCoefficients; ++j )
843 {
844 if( res.postsolve.indices[j] == col )
845 {
846 setRational(scip, colCoef, res.postsolve.values[j]);
847 break;
848 }
849 }
850
851 tmpvars.clear();
852 tmpvars.reserve(rowlen);
853
854 assert(!SCIPrationalIsZero(colCoef));
855 SCIP_VAR* aggrvar = SCIPmatrixGetVar(matrix, col);
856
857 SCIP_CALL( SCIPgetProbvarSumExact(scip, &aggrvar, colCoef, constant) );
859
860 for( int j = startRowCoefficients; j < lastRowCoefficients; ++j )
861 {
862 if( res.postsolve.indices[j] == col )
863 continue;
864
865 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.indices[j]));
866 setRational(scip, tmpvals[c], -res.postsolve.values[j] / colCoef->val);
867 c++;
868 }
869 setRational(scip, updatedSide, side);
870 SCIPrationalDiff(updatedSide, updatedSide, constant);
871 SCIPrationalDiv(updatedSide, updatedSide, colCoef);
872
873 SCIPrationalDebugMessage("Papilo multiaggregate var %s, constant %q \n", SCIPvarGetName(aggrvar), updatedSide);
874
875 SCIP_CALL( SCIPmultiaggregateVarExact(scip, aggrvar, tmpvars.size(),
876 tmpvars.data(), tmpvals, updatedSide, &infeas, &aggregated) );
877
878 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &tmpvals, rowlen);
879 SCIPrationalFreeBuffer(SCIPbuffer(scip), &updatedSide);
882 }
883
884 if( aggregated )
885 *naggrvars += 1;
886 else if( constraintsReplaced && !redundant )
887 {
888 /* if the constraints where replaced, we need to add the failed substitution as an equality to SCIP */
889 SCIP_RATIONAL** tmpvals;
890 SCIP_RATIONAL* tmpside;
891
894
895 setRational(scip, tmpside, side);
896
897 tmpvars.clear();
898 for( int j = startRowCoefficients; j < lastRowCoefficients; ++j )
899 {
900 int idx = j - startRowCoefficients;
901 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.indices[j]));
902 setRational(scip, tmpvals[idx], res.postsolve.values[j]);
903 }
904
905 SCIP_CONS* cons;
906 String name = fmt::format("{}_failed_aggregation_equality", SCIPvarGetName(SCIPmatrixGetVar(matrix, col)));
907 SCIP_CALL( SCIPcreateConsBasicExactLinear(scip, &cons, name.c_str(),
908 tmpvars.size(), tmpvars.data(), tmpvals, tmpside, tmpside ) );
909
910 SCIPdebugMessage("Papilo adding failed aggregation equality: \n");
912 SCIP_CALL( SCIPaddCons(scip, cons) );
913 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
914 *naddconss += 1;
915
917 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &tmpvals, rowlen);
918 }
919
920 if( infeas )
921 {
923 break;
924 }
925
926 break;
927 }
928 case ReductionType::kParallelCol:
929 return SCIP_INVALIDRESULT;
930#if (PAPILO_VERSION_MAJOR <= 1 && PAPILO_VERSION_MINOR==0)
931#else
932 case ReductionType::kFixedInfCol: {
933 if(!constraintsReplaced)
934 continue;
935 SCIP_Bool infeas;
936 SCIP_Bool fixed;
937 SCIP_RATIONAL* value;
938
941
942 int column = res.postsolve.indices[first];
943 bool is_negative_infinity = res.postsolve.values[first] < 0;
944 SCIP_VAR* column_variable = SCIPmatrixGetVar(matrix, column);
945
946 if( is_negative_infinity )
947 {
949 }
950
951 SCIP_CALL( SCIPfixVarExact(scip, column_variable, value, &infeas, &fixed) );
952 *nfixedvars += 1;
953
954 assert(!infeas);
955 assert(fixed);
956
958
959 break;
960 }
961#endif
962#if (PAPILO_VERSION_MAJOR >= 2)
963 case ReductionType::kVarBoundChange :
964 case ReductionType::kRedundantRow :
965 case ReductionType::kRowBoundChange :
966 case ReductionType::kReasonForRowBoundChangeForcedByRow :
967 case ReductionType::kRowBoundChangeForcedByRow :
968 case ReductionType::kSaveRow :
969 case ReductionType::kReducedBoundsCost :
970 case ReductionType::kColumnDualValue :
971 case ReductionType::kRowDualValue :
972 case ReductionType::kCoefficientChange :
973 /* dual ReductionTypes should be only calculated for dual reductions and should not appear for MIP */
974 SCIPerrorMessage("PaPILO: PaPILO should not return dual postsolving reductions in SCIP!!\n");
975 SCIPABORT(); /*lint --e{527}*/
976 break;
977#endif
978 default:
979 SCIPdebugMsg(scip, "PaPILO returned unknown data type: \n" );
980 continue;
981 }
982 }
983
984 /* tighten bounds of variables that are still present after presolving */
985 if( *result != SCIP_CUTOFF )
986 {
987 VariableDomains<papilo::Rational>& varDomains = problem.getVariableDomains();
988 SCIP_RATIONAL* varbound;
991
992 for( int i = 0; i != problem.getNCols(); ++i )
993 {
994 SCIP_VAR* var = SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[i]);
995 if( !varDomains.flags[i].test(ColFlag::kLbInf) )
996 {
997 SCIP_Bool infeas;
998 SCIP_Bool tightened;
999
1000 setRational(scip, varbound, varDomains.lower_bounds[i]);
1001
1002 SCIP_CALL( SCIPtightenVarLbExact(scip, var, varbound, &infeas, &tightened) );
1003
1004 if( tightened )
1005 {
1006 *nchgbds += 1;
1007 SCIPrationalDebugMessage("Papilo tightened lb of variable %s \n", SCIPvarGetName(var));
1008 }
1009
1010 if( infeas )
1011 {
1013 break;
1014 }
1015 }
1016
1017 if( !varDomains.flags[i].test(ColFlag::kUbInf) )
1018 {
1019 SCIP_Bool infeas;
1020 SCIP_Bool tightened;
1021 setRational(scip, varbound, varDomains.upper_bounds[i]);
1022
1023 SCIP_CALL( SCIPtightenVarUbExact(scip, var, varbound, &infeas, &tightened) );
1024
1025 if( tightened )
1026 {
1027 *nchgbds += 1;
1028 SCIPrationalDebugMessage("Papilo tightened ub of variable %s \n", SCIPvarGetName(var));
1029 }
1030
1031 if( infeas )
1032 {
1034 break;
1035 }
1036 }
1037 }
1038
1040 }
1041
1042 /* finish with a final verb message and return */
1044 " (%.1fs) MILP presolver (%d rounds): %d aggregations, %d fixings, %d bound changes\n",
1045 SCIPgetSolvingTime(scip), presolve.getStatistics().nrounds, *naggrvars - oldnaggrvars,
1046 *nfixedvars - oldnfixedvars, *nchgbds - oldnchgbds);
1047
1048 /* free the matrix */
1049 assert(initialized);
1050 SCIPmatrixFree(scip, &matrix);
1051
1052 return SCIP_OKAY;
1053}
1054#endif
1055
1056/** calls PaPILO presolving in floating-point arithmetic */
1057static
1058SCIP_RETCODE performRealPresolving(
1059 SCIP* scip, /**< SCIP data structure */
1060 SCIP_MATRIX* matrix, /**< initialized SCIP_MATRIX data structure */
1061 SCIP_PRESOLMILPDATA* data, /**< plugin specific presol data */
1062 SCIP_Bool initialized, /**< was the matrix initialized */
1063 int* nfixedvars, /**< store number of fixed variables */
1064 int* naggrvars, /**< store number of aggregated variables */
1065 int* nchgvartypes, /**< store number of changed variable types */
1066 int* nchgbds, /**< store number of changed bounds */
1067 int* naddholes, /**< store number of added holes */
1068 int* ndelconss, /**< store the number of deleted cons */
1069 int* naddconss, /**< store the number of added cons */
1070 int* nupgdconss, /**< store the number of upgraded cons */
1071 int* nchgcoefs, /**< store the number of changed coefficients */
1072 int* nchgsides, /**< store the number of changed sides */
1073 SCIP_RESULT* result /**< result pointer */
1074 )
1075{
1076 int nvars = SCIPgetNVars(scip);
1077 int nconss = SCIPgetNConss(scip);
1078
1079 /* only allow communication of constraint modifications by deleting all constraints when some already have been upgraded */
1080 SCIP_CONSHDLR* linconshdlr = SCIPfindConshdlr(scip, "linear");
1081 assert(linconshdlr != NULL);
1082 SCIP_Bool allowconsmodification = (SCIPconshdlrGetNCheckConss(linconshdlr) == SCIPmatrixGetNRows(matrix));
1083
1084 /* store current numbers of aggregations, fixings, and changed bounds for statistics */
1085 int oldnaggrvars = *naggrvars;
1086 int oldnfixedvars = *nfixedvars;
1087 int oldnchgbds = *nchgbds;
1088
1089 /* create presolving objects */
1090 Problem<SCIP_Real> problem = buildProblemReal(scip, matrix);
1091 int oldnnz = problem.getConstraintMatrix().getNnz();
1092 Presolve<SCIP_Real> presolve;
1093 setupPresolve(scip, presolve, data, allowconsmodification);
1094
1095 /* call the presolving */
1097 " (%.1fs) running MILP presolver%s\n", SCIPgetSolvingTime(scip),
1098 presolve.getPresolveOptions().threads == 1 ? "" : " on multiple threads");
1099
1100 /* call presolving without storing information for dual postsolve */
1101#if (PAPILO_VERSION_MAJOR >= 2)
1102 PresolveResult<SCIP_Real> res = presolve.apply(problem, false);
1103#else
1104 PresolveResult<SCIP_Real> res = presolve.apply(problem);
1105#endif
1106 data->lastncols = problem.getNCols();
1107 data->lastnrows = problem.getNRows();
1108
1109 /* evaluate the result */
1110 switch( res.status )
1111 {
1112 case PresolveStatus::kInfeasible:
1115 " (%.1fs) MILP presolver detected infeasibility\n",
1117 SCIPmatrixFree(scip, &matrix);
1118 return SCIP_OKAY;
1119 case PresolveStatus::kUnbndOrInfeas:
1120 case PresolveStatus::kUnbounded:
1123 " (%.1fs) MILP presolver detected unboundedness\n",
1125 SCIPmatrixFree(scip, &matrix);
1126 return SCIP_OKAY;
1127 case PresolveStatus::kUnchanged:
1129 data->lastncols = nvars;
1130 data->lastnrows = nconss;
1132 " (%.1fs) MILP presolver found nothing\n",
1134 SCIPmatrixFree(scip, &matrix);
1135 return SCIP_OKAY;
1136 case PresolveStatus::kReduced:
1137 data->lastncols = problem.getNCols();
1138 data->lastnrows = problem.getNRows();
1140 }
1141
1142 /* result indicated success, now populate the changes into the SCIP structures */
1143
1144 /* tighten bounds of variables that are still present after presolving */
1145 VariableDomains<SCIP_Real>& varDomains = problem.getVariableDomains();
1146 for( int i = 0; i != problem.getNCols(); ++i )
1147 {
1148 assert( ! varDomains.flags[i].test(ColFlag::kInactive) );
1149 SCIP_VAR* var = SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[i]);
1150 if( !varDomains.flags[i].test(ColFlag::kLbInf) )
1151 {
1152 SCIP_Bool infeas;
1153 SCIP_Bool tightened;
1154 SCIP_CALL( SCIPtightenVarLb(scip, var, varDomains.lower_bounds[i], TRUE, &infeas, &tightened) );
1155
1156 if( tightened )
1157 *nchgbds += 1;
1158
1159 if( infeas )
1160 {
1162 break;
1163 }
1164 }
1165
1166 if( !varDomains.flags[i].test(ColFlag::kUbInf) )
1167 {
1168 SCIP_Bool infeas;
1169 SCIP_Bool tightened;
1170 SCIP_CALL( SCIPtightenVarUb(scip, var, varDomains.upper_bounds[i], TRUE, &infeas, &tightened) );
1171
1172 if( tightened )
1173 *nchgbds += 1;
1174
1175 if( infeas )
1176 {
1178 break;
1179 }
1180 }
1181 }
1182
1183 if( *result == SCIP_CUTOFF )
1184 {
1186 " (%.1fs) MILP presolver detected infeasibility\n",
1188 SCIPmatrixFree(scip, &matrix);
1189 return SCIP_OKAY;
1190 }
1191
1192 /* transfer variable fixings and aggregations */
1193 Vec<SCIP_VAR*> tmpvars;
1194 Vec<SCIP_Real> tmpvals;
1195
1196 /* if the size of the problem decreased by a sufficient factor, create all constraints from scratch if allowed */
1197 int newnnz = problem.getConstraintMatrix().getNnz();
1198 bool constraintsReplaced = false;
1199 if( newnnz == 0 || (allowconsmodification &&
1200 (problem.getNCols() <= data->modifyconsfac * SCIPmatrixGetNColumns(matrix) ||
1201 problem.getNRows() <= data->modifyconsfac * SCIPmatrixGetNRows(matrix) ||
1202 newnnz <= data->modifyconsfac * oldnnz)) )
1203 {
1204 int oldnrows = SCIPmatrixGetNRows(matrix);
1205 int newnrows = problem.getNRows();
1206
1207 constraintsReplaced = true;
1208
1209 /* capture constraints that are still present in the problem after presolve */
1210 for( int i = 0; i < newnrows; ++i )
1211 {
1212 SCIP_CONS* c = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
1214 }
1215
1216 /* delete all constraints */
1217 *ndelconss += oldnrows;
1218 *naddconss += newnrows;
1219
1220 for( int i = 0; i < oldnrows; ++i )
1221 {
1223 }
1224
1225 /* now loop over rows of presolved problem and create them as new linear constraints,
1226 * then release the old constraint after its name was passed to the new constraint */
1227 const Vec<RowFlags>& rflags = problem.getRowFlags();
1228 const auto& consmatrix = problem.getConstraintMatrix();
1229 for( int i = 0; i < newnrows; ++i )
1230 {
1231 auto rowvec = consmatrix.getRowCoefficients(i);
1232 const int* rowcols = rowvec.getIndices();
1233 /* SCIPcreateConsBasicLinear() requires a non const pointer */
1234 SCIP_Real* rowvals = const_cast<SCIP_Real*>(rowvec.getValues());
1235 int rowlen = rowvec.getLength();
1236
1237 /* retrieve SCIP compatible left and right hand sides */
1238 SCIP_Real lhs = rflags[i].test(RowFlag::kLhsInf) ? - SCIPinfinity(scip) : consmatrix.getLeftHandSides()[i];
1239 SCIP_Real rhs = rflags[i].test(RowFlag::kRhsInf) ? SCIPinfinity(scip) : consmatrix.getRightHandSides()[i];
1240
1241 /* create variable array matching the value array */
1242 tmpvars.clear();
1243 tmpvars.reserve(rowlen);
1244 for( int j = 0; j < rowlen; ++j )
1245 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.origcol_mapping[rowcols[j]]));
1246
1247 /* create and add new constraint with name of old constraint */
1248 SCIP_CONS* oldcons = SCIPmatrixGetCons(matrix, res.postsolve.origrow_mapping[i]);
1249 SCIP_CONS* cons;
1250 SCIP_CALL( SCIPcreateConsBasicLinear(scip, &cons, SCIPconsGetName(oldcons), rowlen, tmpvars.data(), rowvals, lhs, rhs) );
1251 SCIP_CALL( SCIPaddCons(scip, cons) );
1252
1253 /* release old and new constraint */
1254 SCIP_CALL( SCIPreleaseCons(scip, &oldcons) );
1255 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
1256 }
1257 }
1258
1259 /* PaPILO's aggregations are valid regarding the constraints as they were presolved by PaPILO.
1260 * If coefficients were changed, but constraints in SCIP are not replaced by those from PaPILO,
1261 * then it can not be guaranteed that the bounds of multiaggregated variables will be enforced,
1262 * i.e., will be implied by the constraints in SCIP (see also #3704).
1263 * Only for variable aggregations, SCIP will ensure this by tightening the bounds on the aggregation
1264 * variable as part of SCIPaggregateVars(). For multiaggregations, we will only accept those
1265 * where we can be sure with a simple check that the bounds on the aggregated variable are implied.
1266 */
1267 bool checkmultaggr =
1268#if PAPILO_APIVERSION >= 1
1269 presolve.getStatistics().single_matrix_coefficient_changes > 0
1270#else
1271 presolve.getStatistics().ncoefchgs > 0
1272#endif
1273 && !constraintsReplaced;
1274
1275 /* loop over res.postsolve and add all fixed variables and aggregations to scip */
1276 for( std::size_t i = 0; i != res.postsolve.types.size(); ++i )
1277 {
1278 ReductionType type = res.postsolve.types[i];
1279 int first = res.postsolve.start[i];
1280 int last = res.postsolve.start[i + 1];
1281
1282 switch( type )
1283 {
1284 case ReductionType::kFixedCol:
1285 {
1286 SCIP_Bool infeas;
1287 SCIP_Bool fixed;
1288 int col = res.postsolve.indices[first];
1289
1290 SCIP_VAR* var = SCIPmatrixGetVar(matrix, col);
1291
1292 SCIP_Real value = res.postsolve.values[first];
1293
1294 SCIP_CALL( SCIPfixVar(scip, var, value, &infeas, &fixed) );
1295 *nfixedvars += 1;
1296
1297 assert(!infeas);
1298 /* SCIP has different rules for aggregating variables than PaPILO: therefore the variable PaPILO
1299 * tries to fix now may have been aggregated by SCIP before. Additionally, after aggregation SCIP
1300 * sometimes performs bound tightening resulting in possible fixings. These cases need to be excluded. */
1303 break;
1304 }
1305 /*
1306 * Dual-postsolving in PaPILO required introducing a postsolve-type for substitution with additional information.
1307 * Further, the different Substitution-postsolving types store the required postsolving data differently (in different order) in the postsolving stack.
1308 * Therefore, we need to distinguish how to parse the required data (rowLength, col, side, startRowCoefficients, lastRowCoefficients) from the postsolving stack.
1309 * If these values are accessed, the procedure is the same for both.
1310 */
1311#if (PAPILO_VERSION_MAJOR >= 2)
1312 case ReductionType::kSubstitutedColWithDual:
1313#endif
1314 case ReductionType::kSubstitutedCol:
1315 {
1316 int col = 0;
1317 SCIP_Real side = 0;
1318
1319 int rowlen = 0;
1320 int startRowCoefficients = 0;
1321 int lastRowCoefficients = 0;
1322
1323 if( type == ReductionType::kSubstitutedCol )
1324 {
1325 rowlen = last - first - 1;
1326 col = res.postsolve.indices[first];
1327 side = res.postsolve.values[first];
1328
1329 startRowCoefficients = first + 1;
1330 lastRowCoefficients = last;
1331 }
1332#if (PAPILO_VERSION_MAJOR >= 2)
1333 if( type == ReductionType::kSubstitutedColWithDual )
1334 {
1335 rowlen = (int) res.postsolve.values[first];
1336 col = res.postsolve.indices[first + 3 + rowlen];
1337 side = res.postsolve.values[first + 1];
1338
1339 startRowCoefficients = first + 3;
1340 lastRowCoefficients = first + 3 + rowlen;
1341
1342 assert(side == res.postsolve.values[first + 2]);
1343 assert(res.postsolve.indices[first + 1] == 0);
1344 assert(res.postsolve.indices[first + 2] == 0);
1345 }
1346 assert( type == ReductionType::kSubstitutedCol || type == ReductionType::kSubstitutedColWithDual );
1347#else
1348 assert( type == ReductionType::kSubstitutedCol );
1349#endif
1350 SCIP_Bool infeas;
1351 SCIP_Bool aggregated;
1352 SCIP_Bool redundant = FALSE;
1353 SCIP_Real constant = 0.0;
1354 if( rowlen == 2 )
1355 {
1356 SCIP_Real updatedSide;
1357 SCIP_VAR* varx = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients]);
1358 SCIP_VAR* vary = SCIPmatrixGetVar(matrix, res.postsolve.indices[startRowCoefficients + 1]);
1359 SCIP_Real scalarx = res.postsolve.values[startRowCoefficients];
1360 SCIP_Real scalary = res.postsolve.values[startRowCoefficients + 1];
1361
1362 SCIP_CALL( SCIPgetProbvarSum(scip, &varx, &scalarx, &constant) );
1364
1365 SCIP_CALL( SCIPgetProbvarSum(scip, &vary, &scalary, &constant) );
1367
1368 /* If PaPILO tries to aggregate fixed variables then it missed some obvious fixings.
1369 * This might happen if another aggregation leads to fixings which are not applied immediately by PaPILO.
1370 * With the latest version of PaPILO, this should not occur.
1371 */
1373 {
1374 SCIPdebugMsg(scip, "Aggregation of <%s> and <%s> rejected because they are already fixed.\n",
1375 SCIPvarGetName(varx), SCIPvarGetName(vary));
1376
1377 break;
1378 }
1379
1380 updatedSide = side - constant;
1381
1382 SCIP_CALL( SCIPaggregateVars(scip, varx, vary, scalarx, scalary, updatedSide, &infeas, &redundant, &aggregated) );
1383 }
1384 else
1385 {
1386 SCIP_Real colCoef = 0.0;
1387 SCIP_Real updatedSide;
1388 SCIP_Bool checklbimplied;
1389 SCIP_Bool checkubimplied;
1390 SCIP_Real impliedlb;
1391 SCIP_Real impliedub;
1392 int j;
1393
1394 for( j = startRowCoefficients; j < lastRowCoefficients; ++j )
1395 {
1396 if( res.postsolve.indices[j] == col )
1397 {
1398 colCoef = res.postsolve.values[j];
1399 break;
1400 }
1401 }
1402
1403 tmpvars.clear();
1404 tmpvals.clear();
1405 tmpvars.reserve(rowlen);
1406 tmpvals.reserve(rowlen);
1407
1408 assert(colCoef != 0.0);
1409 SCIP_VAR* aggrvar = SCIPmatrixGetVar(matrix, col);
1410
1411 SCIP_CALL( SCIPgetProbvarSum(scip, &aggrvar, &colCoef, &constant) );
1413
1414 /* If PaPILO tries to multi-aggregate a fixed variable, then it missed some obvious fixings.
1415 * This might happen if another aggregation leads to fixings which are not applied immediately by PaPILO.
1416 * With the latest version of PaPILO, this should not occur.
1417 */
1418 if( SCIPvarGetStatus(aggrvar) == SCIP_VARSTATUS_FIXED )
1419 {
1420 SCIPdebugMsg(scip, "Multi-aggregation of <%s> rejected because it is already fixed.\n",
1421 SCIPvarGetName(aggrvar));
1422
1423 break;
1424 }
1425
1426 updatedSide = side - constant;
1427
1428 /* we need to check whether lb/ub on aggrvar is implied by bounds of other variables in multiaggregation
1429 * if checkmultaggr is TRUE and the lb/ub is finite
1430 * it should be sufficient to ensure global bounds on aggrvar (and as we are in presolve, local=global anyway)
1431 */
1432 checklbimplied = checkmultaggr && !SCIPisInfinity(scip, -SCIPvarGetLbGlobal(aggrvar));
1433 checkubimplied = checkmultaggr && !SCIPisInfinity(scip, SCIPvarGetUbGlobal(aggrvar));
1434 impliedlb = impliedub = updatedSide / colCoef;
1435
1436 for( j = startRowCoefficients; j < lastRowCoefficients; ++j )
1437 {
1438 SCIP_Real coef;
1439 SCIP_VAR* var;
1440
1441 if( res.postsolve.indices[j] == col )
1442 continue;
1443
1444 coef = - res.postsolve.values[j] / colCoef;
1445 var = SCIPmatrixGetVar(matrix, res.postsolve.indices[j]);
1446
1447 if( checklbimplied )
1448 {
1449 if( coef > 0.0 )
1450 {
1451 /* if impliedlb will be -infinity, then we can give up: we cannot use this mutiaggregation */
1453 break;
1454 else
1455 impliedlb += coef * SCIPvarGetLbLocal(var);
1456 }
1457 else
1458 {
1460 break;
1461 else
1462 impliedlb += coef * SCIPvarGetUbLocal(var);
1463 }
1464 }
1465
1466 if( checkubimplied )
1467 {
1468 if( coef > 0.0 )
1469 {
1471 break;
1472 else
1473 impliedub += coef * SCIPvarGetUbLocal(var);
1474 }
1475 else
1476 {
1478 break;
1479 else
1480 impliedub += coef * SCIPvarGetLbLocal(var);
1481 }
1482 }
1483
1484 tmpvals.push_back(coef);
1485 tmpvars.push_back(var);
1486 }
1487
1488 /* if implied bounds are not sufficient to ensure bounds on aggrvar, then we cannot use the multiaggregation */
1489 if( j < lastRowCoefficients )
1490 break;
1491
1492 if( checklbimplied && SCIPisGT(scip, SCIPvarGetLbGlobal(aggrvar), impliedlb) )
1493 break;
1494
1495 if( checkubimplied && SCIPisLT(scip, SCIPvarGetUbGlobal(aggrvar), impliedub) )
1496 break;
1497
1498 SCIP_CALL( SCIPmultiaggregateVar(scip, aggrvar, tmpvars.size(),
1499 tmpvars.data(), tmpvals.data(), updatedSide / colCoef, &infeas, &aggregated) );
1500 }
1501
1502 if( aggregated )
1503 *naggrvars += 1;
1504 else if( constraintsReplaced && !redundant )
1505 {
1506 /* if the constraints where replaced, we need to add the failed substitution as an equality to SCIP */
1507 tmpvars.clear();
1508 tmpvals.clear();
1509 for( int j = startRowCoefficients; j < lastRowCoefficients; ++j )
1510 {
1511 tmpvars.push_back(SCIPmatrixGetVar(matrix, res.postsolve.indices[j]));
1512 tmpvals.push_back(res.postsolve.values[j]);
1513 }
1514
1515 SCIP_CONS* cons;
1516 String name = fmt::format("{}_failed_aggregation_equality", SCIPvarGetName(SCIPmatrixGetVar(matrix, col)));
1517 SCIP_CALL( SCIPcreateConsBasicLinear(scip, &cons, name.c_str(),
1518 tmpvars.size(), tmpvars.data(), tmpvals.data(), side, side ) );
1519 SCIP_CALL( SCIPaddCons(scip, cons) );
1520 SCIP_CALL( SCIPreleaseCons(scip, &cons) );
1521 *naddconss += 1;
1522 }
1523
1524 if( infeas )
1525 {
1527 break;
1528 }
1529
1530 break;
1531 }
1532 case ReductionType::kParallelCol:
1533 return SCIP_INVALIDRESULT;
1534#if PAPILO_VERSION_MAJOR > 1 || (PAPILO_VERSION_MAJOR == 1 && PAPILO_VERSION_MINOR >= 1)
1535 case ReductionType::kFixedInfCol: {
1536 /* todo: currently SCIP can not handle this kind of reduction (see issue #3391) */
1537 assert(false);
1538 if(!constraintsReplaced)
1539 continue;
1540 SCIP_Bool infeas;
1541 SCIP_Bool fixed;
1542 SCIP_Real value = SCIPinfinity(scip);
1543
1544 int column = res.postsolve.indices[first];
1545 bool is_negative_infinity = res.postsolve.values[first] < 0;
1546 SCIP_VAR* column_variable = SCIPmatrixGetVar(matrix, column);
1547
1548 if( is_negative_infinity )
1549 {
1550 value = -SCIPinfinity(scip);
1551 }
1552
1553 SCIP_CALL( SCIPfixVar(scip, column_variable, value, &infeas, &fixed) );
1554 *nfixedvars += 1;
1555
1556 assert(!infeas);
1557 assert(fixed);
1558 break;
1559 }
1560#endif
1561#if (PAPILO_VERSION_MAJOR >= 2)
1562 case ReductionType::kVarBoundChange :
1563 case ReductionType::kRedundantRow :
1564 case ReductionType::kRowBoundChange :
1565 case ReductionType::kReasonForRowBoundChangeForcedByRow :
1566 case ReductionType::kRowBoundChangeForcedByRow :
1567 case ReductionType::kSaveRow :
1568 case ReductionType::kReducedBoundsCost :
1569 case ReductionType::kColumnDualValue :
1570 case ReductionType::kRowDualValue :
1571 case ReductionType::kCoefficientChange :
1572 /* dual ReductionTypes should be only calculated for dual reductions and should not appear for MIP */
1573 SCIPerrorMessage("PaPILO: PaPILO should not return dual postsolving reductions in SCIP!!\n");
1574 SCIPABORT(); /*lint --e{527}*/
1575 break;
1576#endif
1577 default:
1578 SCIPdebugMsg(scip, "PaPILO returned unknown data type: \n" );
1579 continue;
1580 }
1581 }
1582
1583 /* finish with a final verb message and return */
1585 " (%.1fs) MILP presolver (%d rounds): %d aggregations, %d fixings, %d bound changes\n",
1586 SCIPgetSolvingTime(scip), presolve.getStatistics().nrounds, *naggrvars - oldnaggrvars,
1587 *nfixedvars - oldnfixedvars, *nchgbds - oldnchgbds);
1588
1589 /* free the matrix */
1590 assert(initialized);
1591 SCIPmatrixFree(scip, &matrix);
1592
1593 return SCIP_OKAY;
1594}
1595
1596
1597/*
1598 * Callback methods of presolver
1599 */
1600
1601/** copy method for constraint handler plugins (called when SCIP copies plugins) */
1602static
1603SCIP_DECL_PRESOLCOPY(presolCopyMILP)
1604{ /*lint --e{715}*/
1606
1607 return SCIP_OKAY;
1608}
1609
1610/** destructor of presolver to free user data (called when SCIP is exiting) */
1611static
1612SCIP_DECL_PRESOLFREE(presolFreeMILP)
1613{ /*lint --e{715}*/
1614 SCIP_PRESOLMILPDATA* data = (SCIP_PRESOLMILPDATA*)SCIPpresolGetData(presol);
1615 assert(data != NULL);
1616
1617 SCIPpresolSetData(presol, NULL);
1618 SCIPfreeBlockMemory(scip, &data);
1619 return SCIP_OKAY;
1620}
1621
1622/** initialization method of presolver (called after problem was transformed) */
1623static
1624SCIP_DECL_PRESOLINIT(presolInitMILP)
1625{ /*lint --e{715}*/
1626 SCIP_PRESOLMILPDATA* data = (SCIP_PRESOLMILPDATA*)SCIPpresolGetData(presol);
1627 assert(data != NULL);
1628
1629 data->lastncols = -1;
1630 data->lastnrows = -1;
1631
1632 return SCIP_OKAY;
1633}
1634
1635
1636/** execution method of presolver */
1637static
1638SCIP_DECL_PRESOLEXEC(presolExecMILP)
1639{ /*lint --e{715}*/
1640 SCIP_MATRIX* matrix;
1641 SCIP_PRESOLMILPDATA* data;
1642 SCIP_Bool initialized;
1643 SCIP_Bool complete;
1644 SCIP_Bool infeasible;
1645
1647
1648 data = (SCIP_PRESOLMILPDATA*)SCIPpresolGetData(presol);
1649
1650 int nvars = SCIPgetNVars(scip);
1651 int nconss = SCIPgetNConss(scip);
1652
1653 /* run only if the problem size reduced by some amount since the last call or if it is the first call */
1654 if( data->lastncols != -1 && data->lastnrows != -1 &&
1655 nvars > data->lastncols * 0.85 &&
1656 nconss > data->lastnrows * 0.85 )
1657 return SCIP_OKAY;
1658
1659 SCIP_CALL( SCIPmatrixCreate(scip, &matrix, TRUE, &initialized, &complete, &infeasible,
1660 naddconss, ndelconss, nchgcoefs, nchgbds, nfixedvars) );
1661
1662 /* if infeasibility was detected during matrix creation, return here */
1663 if( infeasible )
1664 {
1665 if( initialized )
1666 SCIPmatrixFree(scip, &matrix);
1667
1669 return SCIP_OKAY;
1670 }
1671
1672 /* we only work on pure MIPs, also disable to try building the matrix again if it failed once */
1673 if( !initialized || !complete )
1674 {
1675 data->lastncols = 0;
1676 data->lastnrows = 0;
1677
1678 if( initialized )
1679 SCIPmatrixFree(scip, &matrix);
1680
1681 return SCIP_OKAY;
1682 }
1683
1684 if( 0 != strncmp(data->filename, DEFAULT_FILENAME_PROBLEM, strlen(DEFAULT_FILENAME_PROBLEM)) )
1685 {
1687 " writing transformed problem to %s (only enforced constraints)\n", data->filename);
1688 SCIP_CALL( SCIPwriteTransProblem(scip, data->filename, NULL, FALSE) );
1689 }
1690
1691 if( !SCIPisExact(scip) )
1692 return performRealPresolving(scip, matrix, data, initialized, nfixedvars, naggrvars, nchgvartypes, nchgbds,
1693 naddholes, ndelconss, naddconss, nupgdconss, nchgcoefs, nchgsides, result);
1694#if defined(PAPILO_WITH_EXACTPRESOLVE)
1695 else
1696 return performRationalPresolving(scip, matrix, data, initialized, nfixedvars, naggrvars, nchgvartypes, nchgbds,
1697 naddholes, ndelconss, naddconss, nupgdconss, nchgcoefs, nchgsides, result);
1698#endif
1699
1700 return SCIP_OKAY;
1701}
1702
1703
1704/*
1705 * presolver specific interface methods
1706 */
1707
1708/** creates the MILP presolver and includes it in SCIP */
1710 SCIP* scip /**< SCIP data structure */
1711 )
1712{
1713 SCIP_PRESOLMILPDATA* presoldata;
1714 SCIP_PRESOL* presol;
1715
1716#if defined(PAPILO_VERSION_TWEAK) && PAPILO_VERSION_TWEAK != 0
1717 String name = fmt::format("PaPILO {}.{}.{}.{}", PAPILO_VERSION_MAJOR, PAPILO_VERSION_MINOR, PAPILO_VERSION_PATCH, PAPILO_VERSION_TWEAK);
1718#else
1719 String name = fmt::format("PaPILO {}.{}.{}", PAPILO_VERSION_MAJOR, PAPILO_VERSION_MINOR, PAPILO_VERSION_PATCH);
1720#endif
1721
1722#if defined(PAPILO_GITHASH_AVAILABLE) && defined(PAPILO_TBB)
1723 String desc = fmt::format("parallel presolve for integer and linear optimization (github.com/scipopt/papilo) (built with TBB) [GitHash: {}]", PAPILO_GITHASH);
1724#elif !defined(PAPILO_GITHASH_AVAILABLE) && !defined(PAPILO_TBB)
1725 String desc("parallel presolve for integer and linear optimization (github.com/scipopt/papilo)");
1726#elif defined(PAPILO_GITHASH_AVAILABLE) && !defined(PAPILO_TBB)
1727 String desc = fmt::format("parallel presolve for integer and linear optimization (github.com/scipopt/papilo) [GitHash: {}]", PAPILO_GITHASH);
1728#elif !defined(PAPILO_GITHASH_AVAILABLE) && defined(PAPILO_TBB)
1729 String desc = fmt::format("parallel presolve for integer and linear optimization (github.com/scipopt/papilo) (built with TBB)");
1730#endif
1731
1732 /* add external code info for the presolve library */
1733 SCIP_CALL( SCIPincludeExternalCodeInformation(scip, name.c_str(), desc.c_str()) );
1734
1735 /* create MILP presolver data */
1736 presoldata = NULL;
1737 SCIP_CALL( SCIPallocBlockMemory(scip, &presoldata) );
1738 BMSclearMemory(presoldata);
1739
1740 presol = NULL;
1741
1742 /* include presolver */
1744 presolExecMILP,
1745 (SCIP_PRESOLDATA*)presoldata) );
1746
1747 assert(presol != NULL);
1748
1749 /* set non fundamental callbacks via setter functions */
1750 SCIP_CALL( SCIPsetPresolCopy(scip, presol, presolCopyMILP) );
1751 SCIP_CALL( SCIPsetPresolFree(scip, presol, presolFreeMILP) );
1752 SCIP_CALL( SCIPsetPresolInit(scip, presol, presolInitMILP) );
1753
1754#if defined(PAPILO_WITH_EXACTPRESOLVE)
1755 /* mark presolver as exact */
1756 SCIPpresolMarkExact(presol);
1757#endif
1758
1759 /* add MILP presolver parameters */
1760#ifdef PAPILO_TBB
1762 "presolving/" PRESOL_NAME "/threads",
1763 "maximum number of threads presolving may use (0: automatic)",
1764 &presoldata->threads, FALSE, DEFAULT_THREADS, 0, INT_MAX, NULL, NULL) );
1765#else
1766 presoldata->threads = 1;
1767#endif
1768
1770 "presolving/" PRESOL_NAME "/maxfillinpersubstitution",
1771 "maximal possible fillin for substitutions to be considered",
1772 &presoldata->maxfillinpersubstitution, FALSE, DEFAULT_MAXFILLINPERSUBST, INT_MIN, INT_MAX, NULL, NULL) );
1773
1775 "presolving/" PRESOL_NAME "/maxshiftperrow",
1776 "maximal amount of nonzeros allowed to be shifted to make space for substitutions",
1777 &presoldata->maxshiftperrow, TRUE, DEFAULT_MAXSHIFTPERROW, 0, INT_MAX, NULL, NULL) );
1778
1780 "presolving/" PRESOL_NAME "/randomseed",
1781 "the random seed used for randomization of tie breaking",
1782 &presoldata->randomseed, FALSE, DEFAULT_RANDOMSEED, INT_MIN, INT_MAX, NULL, NULL) );
1783
1784 if( DependentRows<double>::Enabled )
1785 {
1787 "presolving/" PRESOL_NAME "/detectlineardependency",
1788 "should linear dependent equations and free columns be removed? (0: never, 1: for LPs, 2: always)",
1789 &presoldata->detectlineardependency, TRUE, DEFAULT_DETECTLINDEP, 0, 2, NULL, NULL) );
1790 }
1791 else
1792 presoldata->detectlineardependency = DEFAULT_DETECTLINDEP;
1793
1795 "presolving/" PRESOL_NAME "/modifyconsfac",
1796 "modify SCIP constraints when the number of nonzeros or rows is at most this factor "
1797 "times the number of nonzeros or rows before presolving",
1798 &presoldata->modifyconsfac, FALSE, DEFAULT_MODIFYCONSFAC, 0.0, 1.0, NULL, NULL) );
1799
1801 "presolving/" PRESOL_NAME "/markowitztolerance",
1802 "the markowitz tolerance used for substitutions",
1803 &presoldata->markowitztolerance, FALSE, DEFAULT_MARKOWITZTOLERANCE, 0.0, 1.0, NULL, NULL) );
1804
1806 "presolving/" PRESOL_NAME "/hugebound",
1807 "absolute bound value that is considered too huge for activity based calculations",
1808 &presoldata->hugebound, FALSE, DEFAULT_HUGEBOUND, 0.0, SCIP_REAL_MAX, NULL, NULL) );
1809
1810#if PAPILO_APIVERSION >= 2
1811 SCIP_CALL( SCIPaddRealParam(scip, "presolving/" PRESOL_NAME "/abortfacexhaustive",
1812 "abort threshold for exhaustive presolving in PAPILO",
1813 &presoldata->abortfacexhaustive, TRUE, DEFAULT_ABORTFAC_EXHAUSTIVE, 0.0, 1.0, NULL, NULL) );
1814 SCIP_CALL( SCIPaddRealParam(scip, "presolving/" PRESOL_NAME "/abortfacmedium",
1815 "abort threshold for medium presolving in PAPILO",
1816 &presoldata->abortfacmedium, TRUE, DEFAULT_ABORTFAC_MEDIUM, 0.0, 1.0, NULL, NULL) );
1817 SCIP_CALL( SCIPaddRealParam(scip, "presolving/" PRESOL_NAME "/abortfacfast",
1818 "abort threshold for fast presolving in PAPILO",
1819 &presoldata->abortfacfast, TRUE, DEFAULT_ABORTFAC_FAST, 0.0, 1.0, NULL, NULL) );
1820#else
1821 presoldata->abortfacexhaustive = DEFAULT_ABORTFAC_EXHAUSTIVE;
1822 presoldata->abortfacmedium = DEFAULT_ABORTFAC_MEDIUM;
1823 presoldata->abortfacfast = DEFAULT_ABORTFAC_FAST;
1824#endif
1825
1826#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 1)
1827 SCIP_CALL( SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/maxbadgesizeseq",
1828 "maximal badge size in Probing in PaPILO if PaPILO is executed in sequential mode",
1829 &presoldata->maxbadgesizeseq, FALSE, DEFAULT_MAXBADGESIZE_SEQ, -1, INT_MAX, NULL, NULL) );
1830
1831 SCIP_CALL( SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/maxbadgesizepar",
1832 "maximal badge size in Probing in PaPILO if PaPILO is executed in parallel mode",
1833 &presoldata->maxbadgesizepar, FALSE, DEFAULT_MAXBADGESIZE_PAR, -1, INT_MAX, NULL, NULL) );
1834#else
1835 presoldata->maxbadgesizeseq = DEFAULT_MAXBADGESIZE_SEQ;
1836 presoldata->maxbadgesizepar = DEFAULT_MAXBADGESIZE_PAR;
1837#endif
1838
1839#if PAPILO_VERSION_MAJOR > 2 || (PAPILO_VERSION_MAJOR == 2 && PAPILO_VERSION_MINOR >= 3)
1840 SCIP_CALL( SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/internalmaxrounds",
1841 "internal maxrounds for each milp presolving (-1: no limit, 0: model cleanup)",
1842 &presoldata->internalmaxrounds, TRUE, DEFAULT_INTERNAL_MAXROUNDS, -1, INT_MAX, NULL, NULL) );
1843#else
1844 presoldata->internalmaxrounds = DEFAULT_INTERNAL_MAXROUNDS;
1845#endif
1846
1848 "presolving/" PRESOL_NAME "/enableparallelrows",
1849 "should the parallel rows presolver be enabled within the presolve library?",
1850 &presoldata->enableparallelrows, TRUE, DEFAULT_ENABLEPARALLELROWS, NULL, NULL) );
1851
1853 "presolving/" PRESOL_NAME "/enabledomcol",
1854 "should the dominated column presolver be enabled within the presolve library?",
1855 &presoldata->enabledomcol, TRUE, DEFAULT_ENABLEDOMCOL, NULL, NULL) );
1856
1858 "presolving/" PRESOL_NAME "/enabledualinfer",
1859 "should the dualinfer presolver be enabled within the presolve library?",
1860 &presoldata->enabledualinfer, TRUE, DEFAULT_ENABLEDUALINFER, NULL, NULL) );
1861
1863 "presolving/" PRESOL_NAME "/enablemultiaggr",
1864 "should the multi-aggregation presolver be enabled within the presolve library?",
1865 &presoldata->enablemultiaggr, TRUE, DEFAULT_ENABLEMULTIAGGR, NULL, NULL) );
1866
1868 "presolving/" PRESOL_NAME "/enableprobing",
1869 "should the probing presolver be enabled within the presolve library?",
1870 &presoldata->enableprobing, TRUE, DEFAULT_ENABLEPROBING, NULL, NULL) );
1871
1873 "presolving/" PRESOL_NAME "/enablesparsify",
1874 "should the sparsify presolver be enabled within the presolve library?",
1875 &presoldata->enablesparsify, TRUE, DEFAULT_ENABLESPARSIFY, NULL, NULL) );
1876
1877 SCIP_CALL( SCIPaddStringParam(scip, "presolving/" PRESOL_NAME "/probfilename",
1878 "filename to store the problem before MILP presolving starts (only enforced constraints)",
1879 &presoldata->filename, TRUE, DEFAULT_FILENAME_PROBLEM, NULL, NULL) );
1880
1881 SCIP_CALL( SCIPaddIntParam(scip, "presolving/" PRESOL_NAME "/verbosity",
1882 "verbosity level of PaPILO (0: quiet, 1: errors, 2: warnings, 3: normal, 4: detailed)",
1883 &presoldata->verbosity, FALSE, DEFAULT_VERBOSITY, 0, 4, NULL, NULL) );
1884#if PAPILO_APIVERSION >= 6
1886 "presolving/" PRESOL_NAME "/enablecliquemerging",
1887 "should the clique merging presolver be enabled within the presolve library?",
1888 &presoldata->enablecliquemerging, TRUE, DEFAULT_ENABLECLIQUEMERGE, NULL, NULL) );
1890 "presolving/" PRESOL_NAME "/maxedgesparallel",
1891 "maximal amount of edges in the parallel clique merging graph",
1892 &presoldata->maxedgesparallel, FALSE, DEFAULT_MAXEDGESPARALLEL, -1, INT_MAX, NULL, NULL) );
1894 "presolving/" PRESOL_NAME "/maxedgessequential",
1895 "maximal amount of edges in the sequential clique merging graph",
1896 &presoldata->maxedgessequential, FALSE, DEFAULT_MAXEDGESSEQUENTIAL, -1, INT_MAX, NULL, NULL) );
1898 "presolving/" PRESOL_NAME "/maxcliquesize",
1899 "maximal size of clique considered for clique merging",
1900 &presoldata->maxcliquesize, FALSE, DEFAULT_MAXCLIQUESIZE, -1, INT_MAX, NULL, NULL) );
1902 "presolving/" PRESOL_NAME "/maxgreedycalls",
1903 "maximal number of greedy max clique calls in a single thread",
1904 &presoldata->maxgreedycalls, FALSE, DEFAULT_MAXGREEDYCALLS, -1, INT_MAX, NULL, NULL) );
1905#endif
1906
1907 return SCIP_OKAY;
1908}
1909
1910#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:255
#define SCIP_REAL_MAX
Definition def.h:165
#define SCIP_Bool
Definition def.h:98
#define SCIP_Real
Definition def.h:163
#define TRUE
Definition def.h:100
#define FALSE
Definition def.h:101
#define SCIPABORT()
Definition def.h:334
#define REALABS(x)
Definition def.h:189
#define SCIP_CALL(x)
Definition def.h:362
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,...)
#define SCIPdebugMsg
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)
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)
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
SCIP_RETCODE SCIPsetPresolFree(SCIP *scip, SCIP_PRESOL *presol,)
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 SCIPsetPresolCopy(SCIP *scip, SCIP_PRESOL *presol,)
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)
SCIP_RETCODE SCIPsetPresolInit(SCIP *scip, SCIP_PRESOL *presol,)
void SCIPrationalSetInfinity(SCIP_RATIONAL *res)
Definition rational.cpp:619
SCIP_Real SCIPrationalGetReal(SCIP_RATIONAL *rational)
#define SCIPrationalDebugMessage
Definition rational.h:641
void SCIPrationalDiv(SCIP_RATIONAL *res, SCIP_RATIONAL *op1, SCIP_RATIONAL *op2)
SCIP_Bool SCIPrationalIsAbsInfinity(SCIP_RATIONAL *rational)
void SCIPrationalFreeBuffer(BMS_BUFMEM *bufmem, SCIP_RATIONAL **rational)
Definition rational.cpp:474
void SCIPrationalDiff(SCIP_RATIONAL *res, SCIP_RATIONAL *op1, SCIP_RATIONAL *op2)
Definition rational.cpp:984
SCIP_RETCODE SCIPrationalCreateBuffer(BMS_BUFMEM *bufmem, SCIP_RATIONAL **rational)
Definition rational.cpp:124
SCIP_Bool SCIPrationalIsZero(SCIP_RATIONAL *rational)
void SCIPrationalSetNegInfinity(SCIP_RATIONAL *res)
Definition rational.cpp:631
SCIP_Bool SCIPrationalIsNegative(SCIP_RATIONAL *rational)
SCIP_Bool SCIPrationalIsInfinity(SCIP_RATIONAL *rational)
SCIP_RETCODE SCIPrationalCreateBufferArray(BMS_BUFMEM *mem, SCIP_RATIONAL ***rational, int size)
Definition rational.cpp:215
SCIP_Bool SCIPrationalIsNegInfinity(SCIP_RATIONAL *rational)
void SCIPrationalMultReal(SCIP_RATIONAL *res, SCIP_RATIONAL *op1, SCIP_Real op2)
void SCIPrationalFreeBufferArray(BMS_BUFMEM *mem, SCIP_RATIONAL ***ratbufarray, int size)
Definition rational.cpp:519
SCIP_Real SCIPgetSolvingTime(SCIP *scip)
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)
return SCIP_OKAY
int c
assert(minobj< SCIPgetCutoffbound(scip))
int nvars
SCIP_VAR * var
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
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)
public methods for timing
public methods for SCIP variables
scip::Rational val
unsigned int isfprepresentable
definition of wrapper class for rational numbers
struct SCIP_Cons SCIP_CONS
Definition type_cons.h:63
struct SCIP_Conshdlr SCIP_CONSHDLR
Definition type_cons.h:62
struct SCIP_Matrix SCIP_MATRIX
Definition type_matrix.h:42
@ SCIP_VERBLEVEL_HIGH
#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
struct SCIP_Presol SCIP_PRESOL
Definition type_presol.h:50
#define SCIP_DECL_PRESOLINIT(x)
Definition type_presol.h:76
#define SCIP_DECL_PRESOLEXEC(x)
struct SCIP_Rational SCIP_RATIONAL
@ SCIP_ISFPREPRESENTABLE_UNKNOWN
@ 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
enum SCIP_Retcode SCIP_RETCODE
struct Scip SCIP
Definition type_scip.h:39
struct SCIP_Var SCIP_VAR
Definition type_var.h:166
@ 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