Scippy

SCIP

Solving Constraint Integer Programs

presol_domcol.c
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_domcol.c
26 * @ingroup DEFPLUGINS_PRESOL
27 * @brief dominated column presolver
28 * @author Dieter Weninger
29 * @author Gerald Gamrath
30 *
31 * This presolver looks for dominance relations between variable pairs.
32 * From a dominance relation and certain bound/clique-constellations
33 * variable fixings mostly at the lower bound of the dominated variable can be derived.
34 * Additionally it is possible to improve bounds by predictive bound strengthening.
35 *
36 * @todo Also run on general CIPs, if the number of locks of the investigated variables comes only from (upgraded)
37 * linear constraints.
38 *
39 * @todo Instead of choosing variables for comparison out of one row, we should try to use 'hashvalues' for columns that
40 * indicate in which constraint type (<=, >=, or ranged row / ==) they are existing. Then sort the variables (and
41 * corresponding data) after the ranged row/equation hashvalue and only try to derive dominance on variables with
42 * the same hashvalue on ranged row/equation and fitting hashvalues for the other constraint types.
43 * @todo run on incomplete matrices; in order to do so, check at the time when dominance is detected that the locks are
44 * consistent; probably, it is sufficient to check one lock direction for each of the two variables
45 *
46 */
47
48/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
49
51#include "scip/presol_domcol.h"
52#include "scip/pub_matrix.h"
53#include "scip/pub_message.h"
54#include "scip/pub_misc_sort.h"
55#include "scip/pub_presol.h"
56#include "scip/pub_var.h"
57#include "scip/scip_general.h"
58#include "scip/scip_mem.h"
59#include "scip/scip_message.h"
60#include "scip/scip_nlp.h"
61#include "scip/scip_numerics.h"
62#include "scip/scip_param.h"
63#include "scip/scip_presol.h"
64#include "scip/scip_pricer.h"
65#include "scip/scip_prob.h"
66#include "scip/scip_probing.h"
67#include "scip/scip_var.h"
68#include <string.h>
69
70#define PRESOL_NAME "domcol"
71#define PRESOL_DESC "dominated column presolver"
72#define PRESOL_PRIORITY -1000 /**< priority of the presolver (>= 0: before, < 0: after constraint handlers) */
73#define PRESOL_MAXROUNDS -1 /**< maximal number of presolving rounds the presolver participates in (-1: no limit) */
74#define PRESOL_TIMING SCIP_PRESOLTIMING_EXHAUSTIVE /* timing of the presolver (fast, medium, or exhaustive) */
75
76#define DEFAULT_NUMMINPAIRS 1024 /**< minimal number of pair comparisons */
77#define DEFAULT_NUMMAXPAIRS 1048576 /**< maximal number of pair comparisons */
78
79#define DEFAULT_PREDBNDSTR FALSE /**< should predictive bound strengthening be applied? */
80#define DEFAULT_CONTINUOUS_RED TRUE /**< should reductions for continuous variables be carried out? */
81
82
83
84/*
85 * Data structures
86 */
87
88/** control parameters */
89struct SCIP_PresolData
90{
91 int numminpairs; /**< minimal number of pair comparisons */
92 int nummaxpairs; /**< maximal number of pair comparisons */
93 int numcurrentpairs; /**< current number of pair comparisons */
94 SCIP_Bool predbndstr; /**< flag indicating if predictive bound strengthening should be applied */
95 SCIP_Bool continuousred; /**< flag indicating if reductions for continuous variables should be performed */
96};
97
98/** type of fixing direction */
100{
101 FIXATLB = -1, /**< fix variable at lower bound */
102 NOFIX = 0, /**< do not fix variable */
103 FIXATUB = 1 /**< fix variable at upper bound */
106
107
108/*
109 * Local methods
110 */
111#ifdef SCIP_DEBUG
112/** print a row from the constraint matrix */
113static
114void printRow(
115 SCIP* scip, /**< SCIP main data structure */
116 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
117 int row /**< row index for printing */
118 )
119{
120 int* rowpnt;
121 int* rowend;
122 int c;
123 SCIP_Real val;
124 SCIP_Real* valpnt;
125 char relation;
126 SCIP_VAR* var;
127
128 if( !SCIPmatrixIsRowRhsInfinity(matrix, row) &&
129 SCIPisEQ(scip, SCIPmatrixGetRowLhs(matrix, row), SCIPmatrixGetRowRhs(matrix, row)) )
130 {
131 relation='e';
132 }
133 else if( !SCIPmatrixIsRowRhsInfinity(matrix, row) &&
134 !SCIPisEQ(scip, SCIPmatrixGetRowLhs(matrix, row), SCIPmatrixGetRowRhs(matrix, row)) )
135 {
136 relation='r';
137 }
138 else
139 {
140 relation='g';
141 }
142
143 rowpnt = SCIPmatrixGetRowIdxPtr(matrix, row);
144 rowend = rowpnt + SCIPmatrixGetRowNNonzs(matrix, row);
145 valpnt = SCIPmatrixGetRowValPtr(matrix, row);
146
147 SCIPdebugMsgPrint(scip, "\n(L:%g) [%c] %g <=",
148 (SCIPmatrixGetRowNMinActPosInf(matrix, row) + SCIPmatrixGetRowNMinActNegInf(matrix,row) > 0) ?
150 SCIPmatrixGetRowMinActivity(matrix, row), relation, SCIPmatrixGetRowLhs(matrix, row));
151 for(; (rowpnt < rowend); rowpnt++, valpnt++)
152 {
153 c = *rowpnt;
154 val = *valpnt;
155 var = SCIPmatrixGetVar(matrix, c);
156 SCIPdebugMsgPrint(scip, " %g{%s[idx:%d][bnd:%g,%g]}",
157 val, SCIPvarGetName(var), c, SCIPvarGetLbGlobal(var), SCIPvarGetUbGlobal(var));
158 }
159 SCIPdebugMsgPrint(scip, " <= %g (U:%g)", (SCIPmatrixGetRowNMaxActPosInf(matrix, row) + SCIPmatrixGetRowNMaxActNegInf(matrix, row) > 0) ?
161 SCIPmatrixGetRowRhs(matrix, row), SCIPmatrixGetRowMaxActivity(matrix , row));
162}
163
164/** print all rows from the constraint matrix containing a variable */
165static
166SCIP_RETCODE printRowsOfCol(
167 SCIP* scip, /**< SCIP main data structure */
168 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
169 int col /**< column index for printing */
170 )
171{
172 int numrows;
173 int* rows;
174 int i;
175 int* colpnt;
176 int* colend;
177
178 numrows = SCIPmatrixGetColNNonzs(matrix, col);
179
180 SCIP_CALL( SCIPallocBufferArray(scip, &rows, numrows) );
181
182 colpnt = SCIPmatrixGetColIdxPtr(matrix, col);
183 colend = colpnt + SCIPmatrixGetColNNonzs(matrix, col);
184 for( i = 0; (colpnt < colend); colpnt++, i++ )
185 {
186 rows[i] = *colpnt;
187 }
188
189 SCIPdebugMsgPrint(scip, "\n-------");
190 SCIPdebugMsgPrint(scip, "\ncol %d number rows: %d",col,numrows);
191 for( i = 0; i < numrows; i++ )
192 {
193 printRow(scip, matrix, rows[i]);
194 }
195 SCIPdebugMsgPrint(scip, "\n-------");
196
198
199 return SCIP_OKAY;
200}
201
202/** print information about a dominance relation */
203static
204SCIP_RETCODE printDomRelInfo(
205 SCIP* scip, /**< SCIP main data structure */
206 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
207 SCIP_VAR* dominatingvar, /**< dominating variable */
208 int dominatingidx, /**< index of dominating variable */
209 SCIP_VAR* dominatedvar, /**< dominated variable */
210 int dominatedidx, /**< index of dominated variable */
211 SCIP_Real dominatingub, /**< predicted upper bound of dominating variable */
212 SCIP_Real dominatingwclb /**< worst case lower bound of dominating variable */
213 )
214{
215 char type;
216
217 assert(SCIPvarGetType(dominatingvar)==SCIPvarGetType(dominatedvar));
218
219 switch(SCIPvarGetType(dominatingvar))
220 {
222 type='C';
223 break;
225 type='B';
226 break;
229 type='I';
230 break;
231 default:
232 SCIPerrorMessage("unknown variable type\n");
233 SCIPABORT();
234 return SCIP_INVALIDDATA; /*lint !e527*/
235 }
236
237 SCIPdebugMsgPrint(scip, "\n\n### [%c], obj:%g->%g,\t%s[idx:%d](nrows:%d)->%s[idx:%d](nrows:%d)\twclb=%g, ub'=%g, ub=%g",
238 type, SCIPvarGetObj(dominatingvar), SCIPvarGetObj(dominatedvar),
239 SCIPvarGetName(dominatingvar), dominatingidx, SCIPmatrixGetColNNonzs(matrix, dominatingidx),
240 SCIPvarGetName(dominatedvar), dominatedidx, SCIPmatrixGetColNNonzs(matrix, dominatedidx),
241 dominatingwclb, dominatingub, SCIPvarGetUbGlobal(dominatingvar));
242
243 SCIP_CALL( printRowsOfCol(scip, matrix, dominatingidx) );
244
245 return SCIP_OKAY;
246}
247#endif
248
249
250/** get minimum/maximum residual activity for the specified variable and with another variable set to its upper bound */
251static
253 SCIP* scip,
254 SCIP_MATRIX* matrix,
255 int row,
256 int col,
257 SCIP_Real coef,
258 int upperboundcol,
259 SCIP_Real upperboundcoef,
260 SCIP_Real* minresactivity,
261 SCIP_Real* maxresactivity,
262 SCIP_Bool* success
263 )
264{
265 SCIP_VAR* var;
266 SCIP_VAR* upperboundvar;
267 SCIP_Real lb;
268 SCIP_Real ub;
269 SCIP_Real lbupperboundvar;
270 SCIP_Real ubupperboundvar;
271 SCIP_Real tmpmaxact;
272 SCIP_Real tmpminact;
273 int maxactinf;
274 int minactinf;
275
276 assert(scip != NULL);
277 assert(matrix != NULL);
278 assert(row < SCIPmatrixGetNRows(matrix));
279 assert(minresactivity != NULL);
280 assert(maxresactivity != NULL);
281
282 var = SCIPmatrixGetVar(matrix, col);
283 upperboundvar = SCIPmatrixGetVar(matrix, upperboundcol);
284 assert(var != NULL);
285 assert(upperboundvar != NULL);
286
287 lb = SCIPvarGetLbGlobal(var);
288 ub = SCIPvarGetUbGlobal(var);
289
290 maxactinf = SCIPmatrixGetRowNMaxActPosInf(matrix, row) + SCIPmatrixGetRowNMaxActNegInf(matrix, row);
291 minactinf = SCIPmatrixGetRowNMinActPosInf(matrix, row) + SCIPmatrixGetRowNMinActNegInf(matrix, row);
292
293 assert(!SCIPisInfinity(scip, lb));
294 assert(!SCIPisInfinity(scip, -ub));
295
296 lbupperboundvar = SCIPvarGetLbGlobal(upperboundvar);
297 ubupperboundvar = SCIPvarGetUbGlobal(upperboundvar);
298
299 assert(!SCIPisInfinity(scip, lbupperboundvar));
300 assert(!SCIPisInfinity(scip, -ubupperboundvar));
301
302 tmpminact = SCIPmatrixGetRowMinActivity(matrix, row);
303 tmpmaxact = SCIPmatrixGetRowMaxActivity(matrix, row);
304
305 /* first, adjust min and max activity such that upperboundvar is always set to its upper bound */
306
307 /* abort if the upper bound is infty */
308 if( SCIPisInfinity(scip, ubupperboundvar) )
309 {
310 *success = FALSE;
311 return;
312 }
313 else
314 {
315 /* coef > 0 --> the lower bound is used for the minactivity --> update the minactivity */
316 if( upperboundcoef > 0.0 )
317 {
318 if( SCIPisInfinity(scip, -lbupperboundvar) )
319 {
320 assert(minactinf > 0);
321 minactinf--;
322 tmpminact += (upperboundcoef * ubupperboundvar);
323 }
324 else
325 {
326 tmpminact = tmpminact - (upperboundcoef * lbupperboundvar) + (upperboundcoef * ubupperboundvar);
327 }
328 }
329 /* coef < 0 --> the lower bound is used for the maxactivity --> update the maxactivity */
330 else
331 {
332 if( SCIPisInfinity(scip, -lbupperboundvar) )
333 {
334 assert(maxactinf > 0);
335 maxactinf--;
336 tmpmaxact += (upperboundcoef * ubupperboundvar);
337 }
338 else
339 {
340 tmpmaxact = tmpmaxact - (upperboundcoef * lbupperboundvar) + (upperboundcoef * ubupperboundvar);
341 }
342 }
343 }
344
345 *success = TRUE;
346
347 /* the coefficient is positive, so the upper bound contributed to the maxactivity and the lower bound to the minactivity */
348 if( coef >= 0.0 )
349 {
350 if( SCIPisInfinity(scip, ub) )
351 {
352 assert(maxactinf >= 1);
353 if( maxactinf == 1 )
354 {
355 *maxresactivity = tmpmaxact;
356 }
357 else
358 *maxresactivity = SCIPinfinity(scip);
359 }
360 else
361 {
362 if( maxactinf > 0 )
363 *maxresactivity = SCIPinfinity(scip);
364 else
365 *maxresactivity = tmpmaxact - coef * ub;
366 }
367
368 if( SCIPisInfinity(scip, -lb) )
369 {
370 assert(minactinf >= 1);
371 if( minactinf == 1 )
372 {
373 *minresactivity = tmpminact;
374 }
375 else
376 *minresactivity = -SCIPinfinity(scip);
377 }
378 else
379 {
380 if( minactinf > 0 )
381 *minresactivity = -SCIPinfinity(scip);
382 else
383 *minresactivity = tmpminact - coef * lb;
384 }
385 }
386 /* the coefficient is negative, so the lower bound contributed to the maxactivity and the upper bound to the minactivity */
387 else
388 {
389 if( SCIPisInfinity(scip, -lb) )
390 {
391 assert(maxactinf >= 1);
392 if( maxactinf == 1 )
393 {
394 *maxresactivity = tmpmaxact;
395 }
396 else
397 *maxresactivity = SCIPinfinity(scip);
398 }
399 else
400 {
401 if( maxactinf > 0 )
402 *maxresactivity = SCIPinfinity(scip);
403 else
404 *maxresactivity = tmpmaxact - coef * lb;
405 }
406
407 if( SCIPisInfinity(scip, ub) )
408 {
409 assert(minactinf >= 1);
410 if( minactinf == 1 )
411 {
412 *minresactivity = tmpminact;
413 }
414 else
415 *minresactivity = -SCIPinfinity(scip);
416 }
417 else
418 {
419 if( minactinf > 0 )
420 *minresactivity = -SCIPinfinity(scip);
421 else
422 *minresactivity = tmpminact - coef * ub;
423 }
424 }
425}
426
427/** get minimum/maximum residual activity for the specified variable and with another variable set to its lower bound */
428static
430 SCIP* scip, /**< SCIP main data structure */
431 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
432 int row, /**< row index */
433 int col, /**< column index */
434 SCIP_Real coef, /**< coefficient of the column in this row */
435 int lowerboundcol, /**< column index of variable to set to its lower bound */
436 SCIP_Real lowerboundcoef, /**< coefficient in this row of the column to be set to its lower bound */
437 SCIP_Real* minresactivity, /**< minimum residual activity of this row */
438 SCIP_Real* maxresactivity, /**< maximum residual activity of this row */
439 SCIP_Bool* success /**< pointer to store whether the computation was successful */
440 )
441{
442 SCIP_VAR* var;
443 SCIP_VAR* lowerboundvar;
444 SCIP_Real lb;
445 SCIP_Real ub;
446 SCIP_Real lblowerboundvar;
447 SCIP_Real ublowerboundvar;
448 SCIP_Real tmpmaxact;
449 SCIP_Real tmpminact;
450 int maxactinf;
451 int minactinf;
452
453 assert(scip != NULL);
454 assert(matrix != NULL);
455 assert(row < SCIPmatrixGetNRows(matrix));
456 assert(minresactivity != NULL);
457 assert(maxresactivity != NULL);
458
459 var = SCIPmatrixGetVar(matrix, col);;
460 lowerboundvar = SCIPmatrixGetVar(matrix, lowerboundcol);
461 assert(var != NULL);
462 assert(lowerboundvar != NULL);
463
464 lb = SCIPvarGetLbGlobal(var);
465 ub = SCIPvarGetUbGlobal(var);
466
467 maxactinf = SCIPmatrixGetRowNMaxActPosInf(matrix, row) + SCIPmatrixGetRowNMaxActNegInf(matrix, row);
468 minactinf = SCIPmatrixGetRowNMinActPosInf(matrix, row) + SCIPmatrixGetRowNMinActNegInf(matrix, row);
469
470 assert(!SCIPisInfinity(scip, lb));
471 assert(!SCIPisInfinity(scip, -ub));
472
473 lblowerboundvar = SCIPvarGetLbGlobal(lowerboundvar);
474 ublowerboundvar = SCIPvarGetUbGlobal(lowerboundvar);
475
476 assert(!SCIPisInfinity(scip, lblowerboundvar));
477 assert(!SCIPisInfinity(scip, -ublowerboundvar));
478
479 tmpminact = SCIPmatrixGetRowMinActivity(matrix, row);
480 tmpmaxact = SCIPmatrixGetRowMaxActivity(matrix, row);
481
482 /* first, adjust min and max activity such that lowerboundvar is always set to its lower bound */
483
484 /* abort if the lower bound is -infty */
485 if( SCIPisInfinity(scip, -lblowerboundvar) )
486 {
487 *success = FALSE;
488 return;
489 }
490 else
491 {
492 /* coef > 0 --> the upper bound is used for the maxactivity --> update the maxactivity */
493 if( lowerboundcoef > 0.0 )
494 {
495 if( SCIPisInfinity(scip, ublowerboundvar) )
496 {
497 assert(maxactinf > 0);
498 maxactinf--;
499 tmpmaxact += (lowerboundcoef * lblowerboundvar);
500 }
501 else
502 {
503 tmpmaxact = tmpmaxact - (lowerboundcoef * ublowerboundvar) + (lowerboundcoef * lblowerboundvar);
504 }
505 }
506 /* coef < 0 --> the upper bound is used for the minactivity --> update the minactivity */
507 else
508 {
509 if( SCIPisInfinity(scip, ublowerboundvar) )
510 {
511 assert(minactinf > 0);
512 minactinf--;
513 tmpminact += (lowerboundcoef * lblowerboundvar);
514 }
515 else
516 {
517 tmpminact = tmpminact - (lowerboundcoef * ublowerboundvar) + (lowerboundcoef * lblowerboundvar);
518 }
519 }
520 }
521
522 *success = TRUE;
523
524 /* the coefficient is positive, so the upper bound contributed to the maxactivity and the lower bound to the minactivity */
525 if( coef >= 0.0 )
526 {
527 if( SCIPisInfinity(scip, ub) )
528 {
529 assert(maxactinf >= 1);
530 if( maxactinf == 1 )
531 {
532 *maxresactivity = tmpmaxact;
533 }
534 else
535 *maxresactivity = SCIPinfinity(scip);
536 }
537 else
538 {
539 if( maxactinf > 0 )
540 *maxresactivity = SCIPinfinity(scip);
541 else
542 *maxresactivity = tmpmaxact - coef * ub;
543 }
544
545 if( SCIPisInfinity(scip, -lb) )
546 {
547 assert(minactinf >= 1);
548 if( minactinf == 1 )
549 {
550 *minresactivity = tmpminact;
551 }
552 else
553 *minresactivity = -SCIPinfinity(scip);
554 }
555 else
556 {
557 if( minactinf > 0 )
558 *minresactivity = -SCIPinfinity(scip);
559 else
560 *minresactivity = tmpminact - coef * lb;
561 }
562 }
563 /* the coefficient is negative, so the lower bound contributed to the maxactivity and the upper bound to the minactivity */
564 else
565 {
566 if( SCIPisInfinity(scip, -lb) )
567 {
568 assert(maxactinf >= 1);
569 if( maxactinf == 1 )
570 {
571 *maxresactivity = tmpmaxact;
572 }
573 else
574 *maxresactivity = SCIPinfinity(scip);
575 }
576 else
577 {
578 if( maxactinf > 0 )
579 *maxresactivity = SCIPinfinity(scip);
580 else
581 *maxresactivity = tmpmaxact - coef * lb;
582 }
583
584 if( SCIPisInfinity(scip, ub) )
585 {
586 assert(minactinf >= 1);
587 if( minactinf == 1 )
588 {
589 *minresactivity = tmpminact;
590 }
591 else
592 *minresactivity = -SCIPinfinity(scip);
593 }
594 else
595 {
596 if( minactinf > 0 )
597 *minresactivity = -SCIPinfinity(scip);
598 else
599 *minresactivity = tmpminact - coef * ub;
600 }
601 }
602}
603
604/** Calculate bounds of the dominated variable by rowbound analysis.
605 * We use a special kind of predictive rowbound analysis by first setting the dominating variable to its upper bound.
606 */
607static
609 SCIP* scip, /**< SCIP main data structure */
610 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
611 int row, /**< current row index */
612 int coldominating, /**< column index of dominating variable */
613 SCIP_Real valdominating, /**< row coefficient of dominating variable */
614 int coldominated, /**< column index of dominated variable */
615 SCIP_Real valdominated, /**< row coefficient of dominated variable */
616 SCIP_Bool* ubcalculated, /**< was a upper bound calculated? */
617 SCIP_Real* calculatedub, /**< predicted upper bound */
618 SCIP_Bool* wclbcalculated, /**< was a lower worst case bound calculated? */
619 SCIP_Real* calculatedwclb, /**< predicted worst case lower bound */
620 SCIP_Bool* lbcalculated, /**< was a lower bound calculated? */
621 SCIP_Real* calculatedlb, /**< predicted lower bound */
622 SCIP_Bool* wcubcalculated, /**< was a worst case upper bound calculated? */
623 SCIP_Real* calculatedwcub /**< calculated worst case upper bound */
624 )
625{
626 SCIP_Real minresactivity;
627 SCIP_Real maxresactivity;
628 SCIP_Real lhs;
629 SCIP_Real rhs;
630 SCIP_Bool success;
631
632 assert(scip != NULL);
633 assert(matrix != NULL);
634 assert(0 <= row && row < SCIPmatrixGetNRows(matrix));
635 assert(0 <= coldominating && coldominating < SCIPmatrixGetNColumns(matrix));
636 assert(0 <= coldominated && coldominated < SCIPmatrixGetNColumns(matrix));
637
638 assert(ubcalculated != NULL);
639 assert(calculatedub != NULL);
640 assert(wclbcalculated != NULL);
641 assert(calculatedwclb != NULL);
642 assert(lbcalculated != NULL);
643 assert(calculatedlb != NULL);
644 assert(wcubcalculated != NULL);
645 assert(calculatedwcub != NULL);
646
647 assert(!SCIPisZero(scip, valdominated));
648 assert(SCIPmatrixGetVar(matrix, coldominated) != NULL);
649
650 *ubcalculated = FALSE;
651 *wclbcalculated = FALSE;
652 *lbcalculated = FALSE;
653 *wcubcalculated = FALSE;
654
655 /* no rowbound analysis for multiaggregated variables, which should not exist, because the matrix only consists of
656 * active variables
657 */
658 assert(SCIPvarGetStatus(SCIPmatrixGetVar(matrix, coldominating)) != SCIP_VARSTATUS_MULTAGGR);
659 assert(SCIPvarGetStatus(SCIPmatrixGetVar(matrix, coldominated)) != SCIP_VARSTATUS_MULTAGGR);
660
661 lhs = SCIPmatrixGetRowLhs(matrix, row);
662 rhs = SCIPmatrixGetRowRhs(matrix, row);
663 assert(!SCIPisInfinity(scip, lhs));
664 assert(!SCIPisInfinity(scip, -rhs));
665
666 /* get minimum/maximum activity of this row without the dominated variable */
667 getActivityResidualsUpperBound(scip, matrix, row, coldominated, valdominated, coldominating, valdominating,
668 &minresactivity, &maxresactivity, &success);
669
670 if( !success )
671 return SCIP_OKAY;
672
673 assert(!SCIPisInfinity(scip, minresactivity));
674 assert(!SCIPisInfinity(scip, -maxresactivity));
675
676 *calculatedub = SCIPinfinity(scip);
677 *calculatedwclb = -SCIPinfinity(scip);
678 *calculatedlb = -SCIPinfinity(scip);
679 *calculatedwcub = SCIPinfinity(scip);
680
681 /* predictive rowbound analysis */
682
683 if( valdominated > 0.0 )
684 {
685 /* lower bound calculation */
686 if( !SCIPisInfinity(scip, maxresactivity) )
687 {
688 *calculatedlb = (lhs - maxresactivity)/valdominated;
689 *lbcalculated = TRUE;
690 }
691
692 /* worst case calculation of lower bound */
693 if( !SCIPisInfinity(scip, -minresactivity) )
694 {
695 *calculatedwclb = (lhs - minresactivity)/valdominated;
696 *wclbcalculated = TRUE;
697 }
698 else
699 {
700 /* worst case lower bound is infinity */
701 *calculatedwclb = SCIPinfinity(scip);
702 *wclbcalculated = TRUE;
703 }
704
705 /* we can only calculate upper bounds, of the right hand side is finite */
706 if( !SCIPmatrixIsRowRhsInfinity(matrix, row) )
707 {
708 /* upper bound calculation */
709 if( !SCIPisInfinity(scip, -minresactivity) )
710 {
711 *calculatedub = (rhs - minresactivity)/valdominated;
712 *ubcalculated = TRUE;
713 }
714
715 /* worst case calculation of upper bound */
716 if( !SCIPisInfinity(scip, maxresactivity) )
717 {
718 *calculatedwcub = (rhs - maxresactivity)/valdominated;
719 *wcubcalculated = TRUE;
720 }
721 else
722 {
723 /* worst case upper bound is -infinity */
724 *calculatedwcub = -SCIPinfinity(scip);
725 *wcubcalculated = TRUE;
726 }
727 }
728 }
729 else
730 {
731 /* upper bound calculation */
732 if( !SCIPisInfinity(scip, maxresactivity) )
733 {
734 *calculatedub = (lhs - maxresactivity)/valdominated;
735 *ubcalculated = TRUE;
736 }
737
738 /* worst case calculation of upper bound */
739 if( !SCIPisInfinity(scip, -minresactivity) )
740 {
741 *calculatedwcub = (lhs - minresactivity)/valdominated;
742 *wcubcalculated = TRUE;
743 }
744 else
745 {
746 /* worst case upper bound is -infinity */
747 *calculatedwcub = -SCIPinfinity(scip);
748 *wcubcalculated = TRUE;
749 }
750
751 /* we can only calculate lower bounds, of the right hand side is finite */
752 if( !SCIPmatrixIsRowRhsInfinity(matrix, row) )
753 {
754 /* lower bound calculation */
755 if( !SCIPisInfinity(scip, -minresactivity) )
756 {
757 *calculatedlb = (rhs - minresactivity)/valdominated;
758 *lbcalculated = TRUE;
759 }
760
761 /* worst case calculation of lower bound */
762 if( !SCIPisInfinity(scip, maxresactivity) )
763 {
764 *calculatedwclb = (rhs - maxresactivity)/valdominated;
765 *wclbcalculated = TRUE;
766 }
767 else
768 {
769 /* worst case lower bound is infinity */
770 *calculatedwclb = SCIPinfinity(scip);
771 *wclbcalculated = TRUE;
772 }
773 }
774 }
775
776 return SCIP_OKAY;
777}
778
779/** Calculate bounds of the dominating variable by rowbound analysis.
780 * We use a special kind of predictive rowbound analysis by first setting the dominated variable to its lower bound.
781 */
782static
784 SCIP* scip, /**< SCIP main data structure */
785 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
786 int row, /**< current row index */
787 int coldominating, /**< column index of dominating variable */
788 SCIP_Real valdominating, /**< row coefficient of dominating variable */
789 int coldominated, /**< column index of dominated variable */
790 SCIP_Real valdominated, /**< row coefficient of dominated variable */
791 SCIP_Bool* ubcalculated, /**< was a upper bound calculated? */
792 SCIP_Real* calculatedub, /**< predicted upper bound */
793 SCIP_Bool* wclbcalculated, /**< was a lower worst case bound calculated? */
794 SCIP_Real* calculatedwclb, /**< predicted worst case lower bound */
795 SCIP_Bool* lbcalculated, /**< was a lower bound calculated? */
796 SCIP_Real* calculatedlb, /**< predicted lower bound */
797 SCIP_Bool* wcubcalculated, /**< was a worst case upper bound calculated? */
798 SCIP_Real* calculatedwcub /**< calculated worst case upper bound */
799 )
800{
801 SCIP_Real minresactivity;
802 SCIP_Real maxresactivity;
803 SCIP_Real lhs;
804 SCIP_Real rhs;
805 SCIP_Bool success;
806
807 assert(scip != NULL);
808 assert(matrix != NULL);
809 assert(0 <= row && row < SCIPmatrixGetNRows(matrix) );
810 assert(0 <= coldominating && coldominating < SCIPmatrixGetNColumns(matrix));
811 assert(0 <= coldominated && coldominated < SCIPmatrixGetNColumns(matrix));
812
813 assert(ubcalculated != NULL);
814 assert(calculatedub != NULL);
815 assert(wclbcalculated != NULL);
816 assert(calculatedwclb != NULL);
817 assert(lbcalculated != NULL);
818 assert(calculatedlb != NULL);
819 assert(wcubcalculated != NULL);
820 assert(calculatedwcub != NULL);
821
822 assert(!SCIPisZero(scip, valdominating));
823 assert(SCIPmatrixGetVar(matrix, coldominating) != NULL);
824
825 *ubcalculated = FALSE;
826 *wclbcalculated = FALSE;
827 *lbcalculated = FALSE;
828 *wcubcalculated = FALSE;
829
830 /* no rowbound analysis for multiaggregated variables, which should not exist, because the matrix only consists of
831 * active variables
832 */
833 assert(SCIPvarGetStatus(SCIPmatrixGetVar(matrix, coldominating)) != SCIP_VARSTATUS_MULTAGGR);
834 assert(SCIPvarGetStatus(SCIPmatrixGetVar(matrix, coldominated)) != SCIP_VARSTATUS_MULTAGGR);
835
836 lhs = SCIPmatrixGetRowLhs(matrix, row);
837 rhs = SCIPmatrixGetRowRhs(matrix, row);
838 assert(!SCIPisInfinity(scip, lhs));
839 assert(!SCIPisInfinity(scip, -rhs));
840
841 /* get minimum/maximum activity of this row without the dominating variable */
842 getActivityResidualsLowerBound(scip, matrix, row, coldominating, valdominating, coldominated, valdominated,
843 &minresactivity, &maxresactivity, &success);
844
845 if( !success )
846 return SCIP_OKAY;
847
848 assert(!SCIPisInfinity(scip, minresactivity));
849 assert(!SCIPisInfinity(scip, -maxresactivity));
850
851 *calculatedub = SCIPinfinity(scip);
852 *calculatedwclb = -SCIPinfinity(scip);
853 *calculatedlb = -SCIPinfinity(scip);
854 *calculatedwcub = SCIPinfinity(scip);
855
856 /* predictive rowbound analysis */
857
858 if( valdominating > 0.0 )
859 {
860 /* lower bound calculation */
861 if( !SCIPisInfinity(scip, maxresactivity) )
862 {
863 *calculatedlb = (lhs - maxresactivity)/valdominating;
864 *lbcalculated = TRUE;
865 }
866
867 /* worst case calculation of lower bound */
868 if( !SCIPisInfinity(scip, -minresactivity) )
869 {
870 *calculatedwclb = (lhs - minresactivity)/valdominating;
871 *wclbcalculated = TRUE;
872 }
873 else
874 {
875 /* worst case lower bound is infinity */
876 *calculatedwclb = SCIPinfinity(scip);
877 *wclbcalculated = TRUE;
878 }
879
880 /* we can only calculate upper bounds, of the right hand side is finite */
881 if( !SCIPmatrixIsRowRhsInfinity(matrix, row) )
882 {
883 /* upper bound calculation */
884 if( !SCIPisInfinity(scip, -minresactivity) )
885 {
886 *calculatedub = (rhs - minresactivity)/valdominating;
887 *ubcalculated = TRUE;
888 }
889
890 /* worst case calculation of upper bound */
891 if( !SCIPisInfinity(scip, maxresactivity) )
892 {
893 *calculatedwcub = (rhs - maxresactivity)/valdominating;
894 *wcubcalculated = TRUE;
895 }
896 else
897 {
898 /* worst case upper bound is -infinity */
899 *calculatedwcub = -SCIPinfinity(scip);
900 *wcubcalculated = TRUE;
901 }
902 }
903 }
904 else
905 {
906 /* upper bound calculation */
907 if( !SCIPisInfinity(scip, maxresactivity) )
908 {
909 *calculatedub = (lhs - maxresactivity)/valdominating;
910 *ubcalculated = TRUE;
911 }
912
913 /* worst case calculation of upper bound */
914 if( !SCIPisInfinity(scip, -minresactivity) )
915 {
916 *calculatedwcub = (lhs - minresactivity)/valdominating;
917 *wcubcalculated = TRUE;
918 }
919 else
920 {
921 /* worst case upper bound is -infinity */
922 *calculatedwcub = -SCIPinfinity(scip);
923 *wcubcalculated = TRUE;
924 }
925
926 /* we can only calculate lower bounds, of the right hand side is finite */
927 if( !SCIPmatrixIsRowRhsInfinity(matrix, row) )
928 {
929 /* lower bound calculation */
930 if( !SCIPisInfinity(scip, -minresactivity) )
931 {
932 *calculatedlb = (rhs - minresactivity)/valdominating;
933 *lbcalculated = TRUE;
934 }
935
936 /* worst case calculation of lower bound */
937 if( !SCIPisInfinity(scip, maxresactivity) )
938 {
939 *calculatedwclb = (rhs - maxresactivity)/valdominating;
940 *wclbcalculated = TRUE;
941 }
942 else
943 {
944 /* worst case lower bound is infinity */
945 *calculatedwclb = SCIPinfinity(scip);
946 *wclbcalculated = TRUE;
947 }
948 }
949 }
950
951 return SCIP_OKAY;
952}
953
954/** try to find new variable bounds and update them when they are better then the old bounds */
955static
957 SCIP* scip, /**< SCIP main data structure */
958 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
959 int row, /**< row index */
960 int col1, /**< dominating variable index */
961 SCIP_Real val1, /**< dominating variable coefficient */
962 int col2, /**< dominated variable index */
963 SCIP_Real val2, /**< dominated variable coefficient */
964 SCIP_Bool predictdominating, /**< flag indicating if bounds of dominating or dominated variable are predicted */
965 SCIP_Real* upperbound, /**< predicted upper bound */
966 SCIP_Real* wclowerbound, /**< predicted worst case lower bound */
967 SCIP_Real* lowerbound, /**< predicted lower bound */
968 SCIP_Real* wcupperbound /**< predicted worst case upper bound */
969 )
970{
971 SCIP_Bool ubcalculated;
972 SCIP_Bool wclbcalculated;
973 SCIP_Bool lbcalculated;
974 SCIP_Bool wcubcalculated;
975 SCIP_Real newub;
976 SCIP_Real newwclb;
977 SCIP_Real newlb;
978 SCIP_Real newwcub;
979
980 assert(scip != NULL);
981 assert(matrix != NULL);
982 assert(upperbound != NULL);
983 assert(wclowerbound != NULL);
984 assert(lowerbound != NULL);
985 assert(wcupperbound != NULL);
986
987 newub = SCIPinfinity(scip);
988 newlb = -SCIPinfinity(scip);
989 newwcub = newub;
990 newwclb = newlb;
991
992 if( predictdominating )
993 {
994 /* do predictive rowbound analysis for the dominating variable */
995 SCIP_CALL( calcVarBoundsDominating(scip, matrix, row, col1, val1, col2, val2,
996 &ubcalculated, &newub, &wclbcalculated, &newwclb,
997 &lbcalculated, &newlb, &wcubcalculated, &newwcub) );
998 }
999 else
1000 {
1001 /* do predictive rowbound analysis for the dominated variable */
1002 SCIP_CALL( calcVarBoundsDominated(scip, matrix, row, col1, val1, col2, val2,
1003 &ubcalculated, &newub, &wclbcalculated, &newwclb,
1004 &lbcalculated, &newlb, &wcubcalculated, &newwcub) );
1005 }
1006
1007 /* update bounds in case if they are better */
1008 if( ubcalculated )
1009 {
1010 if( newub < *upperbound )
1011 *upperbound = newub;
1012 }
1013 if( wclbcalculated )
1014 {
1015 if( newwclb > *wclowerbound )
1016 *wclowerbound = newwclb;
1017 }
1018 if( lbcalculated )
1019 {
1020 if( newlb > *lowerbound )
1021 *lowerbound = newlb;
1022 }
1023 if( wcubcalculated )
1024 {
1025 if( newwcub < *wcupperbound )
1026 *wcupperbound = newwcub;
1027 }
1028
1029 return SCIP_OKAY;
1030}
1031
1032/** detect parallel columns by using the algorithm of Bixby and Wagner
1033 * see paper: "A note on Detecting Simple Redundancies in Linear Systems", June 1986
1034 */
1035static
1037 SCIP* scip, /**< SCIP main data structure */
1038 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
1039 int* pclass, /**< parallel column classes */
1040 SCIP_Bool* varineq /**< indicating if variable is within an equation */
1041 )
1042{
1043 SCIP_Real* valpnt;
1044 SCIP_Real* values;
1045 SCIP_Real* scale;
1046 int* classsizes;
1047 int* pcset;
1048 int* rowpnt;
1049 int* rowend;
1050 int* colindices;
1051 int* pcs;
1052 SCIP_Real startval;
1053 SCIP_Real aij;
1054 int startpc;
1055 int startk;
1056 int startt;
1057 int pcsetfill;
1058 int colidx;
1059 int k;
1060 int t;
1061 int m;
1062 int i;
1063 int r;
1064 int newpclass;
1065 int pc;
1066 int nrows;
1067 int ncols;
1068
1069 assert(scip != NULL);
1070 assert(matrix != NULL);
1071 assert(pclass != NULL);
1072 assert(varineq != NULL);
1073
1074 nrows = SCIPmatrixGetNRows(matrix);
1075 ncols = SCIPmatrixGetNColumns(matrix);
1076
1077 SCIP_CALL( SCIPallocBufferArray(scip, &classsizes, ncols) );
1078 SCIP_CALL( SCIPallocBufferArray(scip, &scale, ncols) );
1079 SCIP_CALL( SCIPallocBufferArray(scip, &pcset, ncols) );
1080 SCIP_CALL( SCIPallocBufferArray(scip, &values, ncols) );
1081 SCIP_CALL( SCIPallocBufferArray(scip, &colindices, ncols) );
1082 SCIP_CALL( SCIPallocBufferArray(scip, &pcs, ncols) );
1083
1084 BMSclearMemoryArray(scale, ncols);
1085 BMSclearMemoryArray(pclass, ncols);
1086 BMSclearMemoryArray(classsizes, ncols);
1087
1088 classsizes[0] = ncols;
1089 pcsetfill = 0;
1090 for( t = 1; t < ncols; ++t )
1091 pcset[pcsetfill++] = t;
1092
1093 /* loop over all rows */
1094 for( r = 0; r < nrows; ++r )
1095 {
1096 /* we consider only non-empty equations or ranged rows */
1097 if( !SCIPmatrixIsRowRhsInfinity(matrix, r) && SCIPmatrixGetRowNNonzs(matrix, r) > 0 )
1098 {
1099 rowpnt = SCIPmatrixGetRowIdxPtr(matrix, r);
1100 rowend = rowpnt + SCIPmatrixGetRowNNonzs(matrix, r);
1101 valpnt = SCIPmatrixGetRowValPtr(matrix, r);
1102
1103 i = 0;
1104 for( ; (rowpnt < rowend); rowpnt++, valpnt++ )
1105 {
1106 aij = *valpnt;
1107 colidx = *rowpnt;
1108
1109 /* remember variable was part of an equation or ranged row */
1110 varineq[colidx] = TRUE;
1111
1112 if( scale[colidx] == 0.0 )
1113 scale[colidx] = aij;
1114 assert(scale[colidx] != 0.0);
1115
1116 colindices[i] = colidx;
1117 values[i] = aij / scale[colidx];
1118 pc = pclass[colidx];
1119 assert(pc < ncols);
1120
1121 /* update class sizes and pclass set */
1122 assert(classsizes[pc] > 0);
1123 classsizes[pc]--;
1124 if( classsizes[pc] == 0 )
1125 {
1126 assert(pcsetfill < ncols);
1127 pcset[pcsetfill++] = pc;
1128 }
1129 pcs[i] = pc;
1130
1131 i++;
1132 }
1133
1134 assert(i > 0);
1135
1136 /* sort on the pclass values */
1137 if( i > 1 )
1138 {
1139 SCIPsortIntIntReal(pcs, colindices, values, i);
1140 }
1141
1142 k = 0;
1143 while( TRUE ) /*lint !e716*/
1144 {
1145 assert(k < i);
1146 startpc = pcs[k];
1147 startk = k;
1148
1149 /* find pclass-sets */
1150 while( k < i && pcs[k] == startpc )
1151 k++;
1152
1153 /* sort on the A values which have equal pclass values */
1154 if( k - startk > 1 )
1155 SCIPsortRealInt(&(values[startk]), &(colindices[startk]), k - startk);
1156
1157 t = 0;
1158 while( TRUE ) /*lint !e716*/
1159 {
1160 assert(startk + t < i);
1161 startval = values[startk + t];
1162 startt = t;
1163
1164 /* find A-sets */
1165 while( t < k - startk && SCIPisEQ(scip, startval, values[startk + t]) )
1166 t++;
1167
1168 /* get new pclass */
1169 newpclass = pcset[0];
1170 assert(pcsetfill > 0);
1171 pcset[0] = pcset[--pcsetfill];
1172
1173 /* renumbering */
1174 for( m = startk + startt; m < startk + t; m++ )
1175 {
1176 assert(m < i);
1177 assert(colindices[m] < ncols);
1178 assert(newpclass < ncols);
1179
1180 pclass[colindices[m]] = newpclass;
1181 classsizes[newpclass]++;
1182 }
1183
1184 if( t == k - startk )
1185 break;
1186 }
1187
1188 if( k == SCIPmatrixGetRowNNonzs(matrix, r) )
1189 break;
1190 }
1191 }
1192 }
1193
1195 SCIPfreeBufferArray(scip, &colindices);
1196 SCIPfreeBufferArray(scip, &values);
1197 SCIPfreeBufferArray(scip, &pcset);
1198 SCIPfreeBufferArray(scip, &scale);
1199 SCIPfreeBufferArray(scip, &classsizes);
1200
1201 return SCIP_OKAY;
1202}
1203
1204
1205/** try to improve variable bounds by predictive bound strengthening */
1206static
1208 SCIP* scip, /**< SCIP main data structure */
1209 SCIP_VAR* dominatingvar, /**< dominating variable */
1210 int dominatingidx, /**< column index of the dominating variable */
1211 SCIP_Real dominatingub, /**< predicted upper bound of the dominating variable */
1212 SCIP_Real dominatinglb, /**< predicted lower bound of the dominating variable */
1213 SCIP_Real dominatingwcub, /**< predicted worst case upper bound of the dominating variable */
1214 SCIP_VAR* dominatedvar, /**< dominated variable */
1215 int dominatedidx, /**< column index of the dominated variable */
1216 SCIP_Real dominatedub, /**< predicted upper bound of the dominated variable */
1217 SCIP_Real dominatedwclb, /**< predicted worst case lower bound of the dominated variable */
1218 SCIP_Real dominatedlb, /**< predicted lower bound of the dominated variable */
1219 FIXINGDIRECTION* varstofix, /**< array holding fixing information */
1220 int* nchgbds /**< count number of bound changes */
1221 )
1222{
1223 /* we compare only variables from the same type */
1224 if( !(SCIPvarGetType(dominatingvar) == SCIPvarGetType(dominatedvar) ||
1225 SCIPvarIsBinary(dominatingvar) == SCIPvarIsBinary(dominatedvar) ||
1226 (SCIPvarGetType(dominatingvar) == SCIP_VARTYPE_INTEGER && SCIPvarGetType(dominatedvar) == SCIP_VARTYPE_IMPLINT) ||
1227 (SCIPvarGetType(dominatedvar) == SCIP_VARTYPE_INTEGER && SCIPvarGetType(dominatingvar) == SCIP_VARTYPE_IMPLINT)) )
1228 {
1229 return SCIP_OKAY;
1230 }
1231
1232 if( varstofix[dominatingidx] == NOFIX )
1233 {
1234 /* assume x dominates y (x->y). we get this bound from a positive coefficient
1235 * of the dominating variable. because x->y the dominated variable y has
1236 * a positive coefficient too. thus y contributes to the minactivity with its
1237 * lower bound. but this case is considered within predictive bound analysis.
1238 * thus the dominating upper bound is a common upper bound.
1239 */
1240 if( !SCIPisInfinity(scip, dominatingub) )
1241 {
1242 SCIP_Real newub;
1243 SCIP_Real oldub;
1244 SCIP_Real lb;
1245
1246 newub = dominatingub;
1247 oldub = SCIPvarGetUbGlobal(dominatingvar);
1248 lb = SCIPvarGetLbGlobal(dominatingvar);
1249
1250 /* if( SCIPvarGetType(dominatingvar) != SCIP_VARTYPE_CONTINUOUS )
1251 newub = SCIPceil(scip, newub); */
1252
1253 if( SCIPisLE(scip, lb, newub) && SCIPisLT(scip, newub, oldub) )
1254 {
1255 SCIPdebugMsg(scip, "[ub]\tupper bound for dominating variable <%s> changed: [%.17f,%.17f] -> [%.17f,%.17f]\n",
1256 SCIPvarGetName(dominatingvar), lb, oldub, lb, newub);
1257 SCIP_CALL( SCIPchgVarUb(scip, dominatingvar, newub) );
1258 (*nchgbds)++;
1259 }
1260 }
1261
1262 /* assume x dominates y (x->y). we get this lower bound of the dominating variable
1263 * from a negative coefficient within a <= relation. if y has a positive coefficient
1264 * we get a common lower bound of x from predictive bound analysis. if y has a
1265 * negative coefficient we get an improved lower bound of x because the minactivity
1266 * is greater. for discrete variables we have to round down the lower bound.
1267 */
1268 if( !SCIPisInfinity(scip, -dominatinglb) )
1269 {
1270 SCIP_Real newlb;
1271 SCIP_Real oldlb;
1272 SCIP_Real ub;
1273
1274 newlb = dominatinglb;
1275 oldlb = SCIPvarGetLbGlobal(dominatingvar);
1276 ub = SCIPvarGetUbGlobal(dominatingvar);
1277
1278 if( SCIPvarGetType(dominatingvar) != SCIP_VARTYPE_CONTINUOUS )
1279 newlb = SCIPfloor(scip, newlb);
1280
1281 if( SCIPisLT(scip, oldlb, newlb) && SCIPisLE(scip, newlb, ub) )
1282 {
1283 SCIPdebugMsg(scip, "[lb]\tlower bound of dominating variable <%s> changed: [%.17f,%.17f] -> [%.17f,%.17f]\n",
1284 SCIPvarGetName(dominatingvar), oldlb, ub, newlb, ub);
1285 SCIP_CALL( SCIPchgVarLb(scip, dominatingvar, newlb) );
1286 (*nchgbds)++;
1287 }
1288 }
1289
1290 /* assume x dominates y (x->y). we get this bound from a positive coefficient
1291 * of x within a <= relation. from x->y it follows, that y has a positive
1292 * coefficient in this row too. the worst case upper bound of x is an estimation
1293 * of how great x can be in every case. if the objective coefficient of x is
1294 * negative we get thus a lower bound of x. for discrete variables we have
1295 * to round down the lower bound before.
1296 */
1297 if( !SCIPisInfinity(scip, dominatingwcub) && SCIPisNegative(scip, SCIPvarGetObj(dominatingvar)))
1298 {
1299 SCIP_Real newlb;
1300 SCIP_Real oldlb;
1301 SCIP_Real ub;
1302
1303 newlb = dominatingwcub;
1304 oldlb = SCIPvarGetLbGlobal(dominatingvar);
1305 ub = SCIPvarGetUbGlobal(dominatingvar);
1306
1307 if( SCIPvarGetType(dominatingvar) != SCIP_VARTYPE_CONTINUOUS )
1308 newlb = SCIPfloor(scip, newlb);
1309
1310 if( SCIPisLT(scip, oldlb, newlb) && SCIPisLE(scip, newlb, ub) )
1311 {
1312 SCIPdebugMsg(scip, "[wcub]\tlower bound of dominating variable <%s> changed: [%.17f,%.17f] -> [%.17f,%.17f]\n",
1313 SCIPvarGetName(dominatingvar), oldlb, ub, newlb, ub);
1314 SCIP_CALL( SCIPchgVarLb(scip, dominatingvar, newlb) );
1315 (*nchgbds)++;
1316 }
1317 }
1318 }
1319
1320 if( varstofix[dominatedidx] == NOFIX )
1321 {
1322 /* assume x dominates y (x->y). we get this bound for a positive coefficient of y
1323 * within a <= relation. if x has a negative coefficient we get a common upper
1324 * bound of y. if x has a positive coefficient we get an improved upper bound
1325 * of y because the minactivity is greater.
1326 */
1327 if( !SCIPisInfinity(scip, dominatedub) )
1328 {
1329 SCIP_Real newub;
1330 SCIP_Real oldub;
1331 SCIP_Real lb;
1332
1333 newub = dominatedub;
1334 oldub = SCIPvarGetUbGlobal(dominatedvar);
1335 lb = SCIPvarGetLbGlobal(dominatedvar);
1336
1337 if( SCIPisLE(scip, lb, newub) && SCIPisLT(scip, newub, oldub) )
1338 {
1339 SCIPdebugMsg(scip, "[ub]\tupper bound of dominated variable <%s> changed: [%.17f,%.17f] -> [%.17f,%.17f]\n",
1340 SCIPvarGetName(dominatedvar), lb, oldub, lb, newub);
1341 SCIP_CALL( SCIPchgVarUb(scip, dominatedvar, newub) );
1342 (*nchgbds)++;
1343 }
1344 }
1345
1346 /* assume x dominates y (x->y). we get this bound only from a negative
1347 * coefficient of y within a <= relation. because of x->y then x has a negative
1348 * coefficient too. the worst case lower bound is an estimation what property
1349 * the dominated variable must have if the dominating variable is at its upper bound.
1350 * to get an upper bound of the dominated variable we need to consider a positive
1351 * objective coefficient. for discrete variables we have to round up the upper bound.
1352 */
1353 if( !SCIPisInfinity(scip, -dominatedwclb) && SCIPisPositive(scip, SCIPvarGetObj(dominatedvar)) )
1354 {
1355 SCIP_Real newub;
1356 SCIP_Real oldub;
1357 SCIP_Real lb;
1358
1359 newub = dominatedwclb;
1360 oldub = SCIPvarGetUbGlobal(dominatedvar);
1361 lb = SCIPvarGetLbGlobal(dominatedvar);
1362
1363 if( SCIPvarGetType(dominatedvar) != SCIP_VARTYPE_CONTINUOUS )
1364 newub = SCIPceil(scip, newub);
1365
1366 if( SCIPisLE(scip, lb, newub) && SCIPisLT(scip, newub, oldub) )
1367 {
1368 SCIPdebugMsg(scip, "[wclb]\tupper bound of dominated variable <%s> changed: [%.17f,%.17f] -> [%.17f,%.17f]\n",
1369 SCIPvarGetName(dominatedvar), lb, oldub, lb, newub);
1370 SCIP_CALL( SCIPchgVarUb(scip, dominatedvar, newub) );
1371 (*nchgbds)++;
1372 }
1373 }
1374
1375 /* assume x dominates y (x->y). we get a lower bound of y from a negative
1376 * coefficient within a <= relation. but if x->y then x has a negative
1377 * coefficient too and contributes with its upper bound to the minactivity.
1378 * thus in all we have a common lower bound calculation and no rounding is necessary.
1379 */
1380 if( !SCIPisInfinity(scip, -dominatedlb) )
1381 {
1382 SCIP_Real newlb;
1383 SCIP_Real oldlb;
1384 SCIP_Real ub;
1385
1386 newlb = dominatedlb;
1387 oldlb = SCIPvarGetLbGlobal(dominatedvar);
1388 ub = SCIPvarGetUbGlobal(dominatedvar);
1389
1390 if( SCIPisLT(scip, oldlb, newlb) && SCIPisLE(scip, newlb, ub) )
1391 {
1392 SCIPdebugMsg(scip, "[lb]\tlower bound of dominated variable <%s> changed: [%.17f,%.17f] -> [%.17f,%.17f]\n",
1393 SCIPvarGetName(dominatedvar), oldlb, ub, newlb, ub);
1394 SCIP_CALL( SCIPchgVarLb(scip, dominatedvar, newlb) );
1395 (*nchgbds)++;
1396 }
1397 }
1398 }
1399
1400 return SCIP_OKAY;
1401}
1402
1403/** try to find variable fixings */
1404static
1406 SCIP* scip, /**< SCIP main data structure */
1407 SCIP_MATRIX* matrix, /**< constraint matrix structure */
1408 SCIP_VAR* dominatingvar, /**< dominating variable */
1409 int dominatingidx, /**< column index of the dominating variable */
1410 SCIP_Real dominatingub, /**< predicted upper bound of the dominating variable */
1411 SCIP_Real dominatingwclb, /**< predicted worst case lower bound of the dominating variable */
1412 SCIP_Real dominatinglb, /**< predicted lower bound of the dominating variable */
1413 SCIP_Real dominatingwcub, /**< predicted worst case upper bound of the dominating variable */
1414 SCIP_VAR* dominatedvar, /**< dominated variable */
1415 int dominatedidx, /**< column index of the dominated variable */
1416 FIXINGDIRECTION* varstofix, /**< array holding fixing information */
1417 SCIP_Bool onlybinvars, /**< flag indicating only binary variables are present */
1418 SCIP_Bool onlyoneone, /**< when onlybinvars is TRUE, flag indicates if both binary variables are in clique */
1419 int* nfixings /**< counter for possible fixings */
1420 )
1421{
1422 /* we compare only variables from the same type */
1423 if( !(SCIPvarGetType(dominatingvar) == SCIPvarGetType(dominatedvar) ||
1424 SCIPvarIsBinary(dominatingvar) == SCIPvarIsBinary(dominatedvar) ||
1425 (SCIPvarGetType(dominatingvar) == SCIP_VARTYPE_INTEGER && SCIPvarGetType(dominatedvar) == SCIP_VARTYPE_IMPLINT) ||
1426 (SCIPvarGetType(dominatedvar) == SCIP_VARTYPE_INTEGER && SCIPvarGetType(dominatingvar) == SCIP_VARTYPE_IMPLINT)) )
1427 {
1428 return SCIP_OKAY;
1429 }
1430
1431 if( varstofix[dominatedidx] == NOFIX && SCIPmatrixGetColNNonzs(matrix, dominatingidx) == 1
1432 && SCIPmatrixGetColNNonzs(matrix, dominatedidx) == 1 )
1433 {
1434 /* We have a x->y dominance relation and only one equality constraint
1435 * where the dominating variable x with an infinity upper bound and the
1436 * dominated variable y are present. Then the dominated variable y
1437 * can be fixed at its lower bound.
1438 */
1439 int row;
1440 row = *(SCIPmatrixGetColIdxPtr(matrix, dominatedidx));
1441
1442 if( SCIPisEQ(scip, SCIPmatrixGetRowLhs(matrix, row), SCIPmatrixGetRowRhs(matrix, row)) &&
1443 SCIPisInfinity(scip, SCIPvarGetUbGlobal(dominatingvar)) )
1444 {
1445 varstofix[dominatedidx] = FIXATLB;
1446 (*nfixings)++;
1447
1448 return SCIP_OKAY;
1449 }
1450 }
1451
1452 if( varstofix[dominatedidx] == NOFIX && !SCIPisNegative(scip, SCIPvarGetObj(dominatedvar)) )
1453 {
1454 if( !SCIPisInfinity(scip, -dominatingwclb) &&
1455 SCIPisLE(scip, dominatingwclb, SCIPvarGetUbGlobal(dominatingvar)) )
1456 {
1457 /* we have a x->y dominance relation with a positive obj coefficient
1458 * of the dominated variable y. we need to secure feasibility
1459 * by testing if the predicted lower worst case bound is less equal the
1460 * current upper bound. it is possible, that the lower worst case bound
1461 * is infinity and the upper bound of the dominating variable x is
1462 * infinity too.
1463 */
1464 varstofix[dominatedidx] = FIXATLB;
1465 (*nfixings)++;
1466 }
1467 }
1468
1469 if( varstofix[dominatedidx] == NOFIX && !SCIPisInfinity(scip, dominatingub) &&
1470 SCIPisLE(scip, dominatingub, SCIPvarGetUbGlobal(dominatingvar)) )
1471 {
1472 /* we have a x->y dominance relation with an arbitrary obj coefficient
1473 * of the dominating variable x. in all cases we have to look
1474 * if the predicted upper bound of the dominating variable is great enough.
1475 * by testing, that the predicted upper bound is not infinity we avoid problems
1476 * with x->y e.g.
1477 * min -x -y
1478 * s.t. -x -y <= -1
1479 * 0<=x<=1, 0<=y<=1
1480 * where y is not at their lower bound.
1481 */
1482 varstofix[dominatedidx] = FIXATLB;
1483 (*nfixings)++;
1484 }
1485
1486 if( varstofix[dominatingidx] == NOFIX && !SCIPisPositive(scip, SCIPvarGetObj(dominatingvar)) )
1487 {
1488 /* we have a x->y dominance relation with a negative obj coefficient
1489 * of the dominating variable x. if the worst case upper bound is
1490 * greater equal than upper bound, we fix x at the upper bound
1491 */
1492 if( !SCIPisInfinity(scip, dominatingwcub) &&
1493 SCIPisGE(scip, dominatingwcub, SCIPvarGetUbGlobal(dominatingvar)) )
1494 {
1495 varstofix[dominatingidx] = FIXATUB;
1496 (*nfixings)++;
1497 }
1498 }
1499
1500 if( varstofix[dominatingidx] == NOFIX && !SCIPisInfinity(scip, -dominatinglb) &&
1501 SCIPisGE(scip, dominatinglb, SCIPvarGetUbGlobal(dominatingvar)) )
1502 {
1503 /* we have a x->y dominance relation with an arbitrary obj coefficient
1504 * of the dominating variable x. if the predicted lower bound is greater
1505 * equal than upper bound, we fix x at the upper bound.
1506 */
1507 varstofix[dominatingidx] = FIXATUB;
1508 (*nfixings)++;
1509 }
1510
1511 if( onlybinvars )
1512 {
1513 if( varstofix[dominatedidx] == NOFIX && (onlyoneone || SCIPvarsHaveCommonClique(dominatingvar, TRUE, dominatedvar, TRUE, TRUE)) )
1514 {
1515 /* We have a (1->1)-clique with dominance relation (x->y) (x dominates y).
1516 * From this dominance relation, we know (1->0) is possible and not worse than (0->1)
1517 * concerning the objective function. It follows that only (1->0) or (0->0) are possible,
1518 * but in both cases y has the value 0 => y=0.
1519 */
1520 varstofix[dominatedidx] = FIXATLB;
1521 (*nfixings)++;
1522 }
1523
1524 if( varstofix[dominatingidx] == NOFIX && SCIPvarsHaveCommonClique(dominatingvar, FALSE, dominatedvar, FALSE, TRUE) )
1525 {
1526 /* We have a (0->0)-clique with dominance relation x->y (x dominates y).
1527 * From this dominance relation, we know (1->0) is possible and not worse than (0->1)
1528 * concerning the objective function. It follows that only (1->0) or (1->1) are possible,
1529 * but in both cases x has the value 1 => x=1
1530 */
1531 varstofix[dominatingidx] = FIXATUB;
1532 (*nfixings)++;
1533 }
1534 }
1535 else
1536 assert(!onlyoneone);
1537
1538 return SCIP_OKAY;
1539}
1540
1541/** find dominance relation between variable pairs */
1542static
1544 SCIP* scip, /**< SCIP main data structure */
1545 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
1546 SCIP_PRESOLDATA* presoldata, /**< presolver data */
1547 int* searchcols, /**< indexes of variables for pair comparisons */
1548 int searchsize, /**< number of variables for pair comparisons */
1549 SCIP_Bool onlybinvars, /**< flag indicating searchcols contains only binary variable indexes */
1550 FIXINGDIRECTION* varstofix, /**< array holding information for later upper/lower bound fixing */
1551 int* nfixings, /**< found number of possible fixings */
1552 SCIP_Longint* ndomrelations, /**< found number of dominance relations */
1553 int* nchgbds /**< number of changed bounds */
1554 )
1555{
1556 SCIP_Real* vals1;
1557 SCIP_Real* vals2;
1558 SCIP_Real tmpupperbounddominatingcol1;
1559 SCIP_Real tmpupperbounddominatingcol2;
1560 SCIP_Real tmpwclowerbounddominatingcol1;
1561 SCIP_Real tmpwclowerbounddominatingcol2;
1562 SCIP_Real tmplowerbounddominatingcol1;
1563 SCIP_Real tmplowerbounddominatingcol2;
1564 SCIP_Real tmpwcupperbounddominatingcol1;
1565 SCIP_Real tmpwcupperbounddominatingcol2;
1566 int* rows1;
1567 int* rows2;
1568 int nrows1;
1569 int nrows2;
1570 SCIP_Real tmpupperbounddominatedcol1;
1571 SCIP_Real tmpupperbounddominatedcol2;
1572 SCIP_Real tmpwclowerbounddominatedcol1;
1573 SCIP_Real tmpwclowerbounddominatedcol2;
1574 SCIP_Real tmplowerbounddominatedcol1;
1575 SCIP_Real tmplowerbounddominatedcol2;
1576 SCIP_Real tmpwcupperbounddominatedcol1;
1577 SCIP_Real tmpwcupperbounddominatedcol2;
1578 SCIP_Real obj1;
1579 SCIP_Bool col1domcol2;
1580 SCIP_Bool col2domcol1;
1581 SCIP_Bool onlyoneone;
1582 int cnt1;
1583 int cnt2;
1584 int col1;
1585 int col2;
1586 int r1;
1587 int r2;
1588 int paircnt;
1589 int oldnfixings;
1590
1591 assert(scip != NULL);
1592 assert(matrix != NULL);
1593 assert(presoldata != NULL);
1594 assert(searchcols != NULL);
1595 assert(varstofix != NULL);
1596 assert(nfixings != NULL);
1597 assert(ndomrelations != NULL);
1598 assert(nchgbds != NULL);
1599
1600 paircnt = 0;
1601 oldnfixings = *nfixings;
1602
1603 /* pair comparisons */
1604 for( cnt1 = 0; cnt1 < searchsize; cnt1++ )
1605 {
1606 SCIP_VAR* varcol1;
1607 SCIP_VAR* varcol2;
1608
1609 /* get index of the first variable */
1610 col1 = searchcols[cnt1];
1611
1612 if( varstofix[col1] == FIXATLB )
1613 continue;
1614
1615 varcol1 = SCIPmatrixGetVar(matrix, col1);
1616 obj1 = SCIPvarGetObj(varcol1);
1617
1618 for( cnt2 = cnt1+1; cnt2 < searchsize; cnt2++ )
1619 {
1620 /* get index of the second variable */
1621 col2 = searchcols[cnt2];
1622 varcol2 = SCIPmatrixGetVar(matrix, col2);
1623 onlyoneone = FALSE;
1624
1625 /* we always have minimize as obj sense */
1626
1627 /* column 1 dominating column 2 */
1628 col1domcol2 = (obj1 <= SCIPvarGetObj(varcol2));
1629
1630 /* column 2 dominating column 1 */
1631 col2domcol1 = (SCIPvarGetObj(varcol2) <= obj1);
1632
1633 /* search only if nothing was found yet */
1634 col1domcol2 = col1domcol2 && (varstofix[col2] == NOFIX);
1635 col2domcol1 = col2domcol1 && (varstofix[col1] == NOFIX);
1636
1637 /* we only search for a dominance relation if the lower bounds are not negative */
1638 if( !onlybinvars )
1639 {
1640 if( SCIPisLT(scip, SCIPvarGetLbGlobal(varcol1), 0.0) ||
1641 SCIPisLT(scip, SCIPvarGetLbGlobal(varcol2), 0.0) )
1642 {
1643 col1domcol2 = FALSE;
1644 col2domcol1 = FALSE;
1645 }
1646 }
1647
1648 /* pair comparison control */
1649 if( paircnt == presoldata->numcurrentpairs )
1650 {
1651 assert(*nfixings >= oldnfixings);
1652 if( *nfixings == oldnfixings )
1653 {
1654 /* not enough fixings found, decrement number of comparisons */
1655 presoldata->numcurrentpairs >>= 1; /*lint !e702*/
1656 if( presoldata->numcurrentpairs < presoldata->numminpairs )
1657 presoldata->numcurrentpairs = presoldata->numminpairs;
1658
1659 /* stop searching in this row */
1660 return SCIP_OKAY;
1661 }
1662 oldnfixings = *nfixings;
1663 paircnt = 0;
1664
1665 /* increment number of comparisons */
1666 presoldata->numcurrentpairs <<= 1; /*lint !e701*/
1667 if( presoldata->numcurrentpairs > presoldata->nummaxpairs )
1668 presoldata->numcurrentpairs = presoldata->nummaxpairs;
1669 }
1670 paircnt++;
1671
1672 if( !col1domcol2 && !col2domcol1 )
1673 continue;
1674
1675 /* get the data for both columns */
1676 vals1 = SCIPmatrixGetColValPtr(matrix, col1);
1677 rows1 = SCIPmatrixGetColIdxPtr(matrix, col1);
1678 nrows1 = SCIPmatrixGetColNNonzs(matrix, col1);
1679 r1 = 0;
1680 vals2 = SCIPmatrixGetColValPtr(matrix, col2);
1681 rows2 = SCIPmatrixGetColIdxPtr(matrix, col2);
1682 nrows2 = SCIPmatrixGetColNNonzs(matrix, col2);
1683 r2 = 0;
1684
1685 /* do we have a obj constant? */
1686 if( nrows1 == 0 || nrows2 == 0 )
1687 continue;
1688
1689 /* initialize temporary bounds of dominating variable */
1690 tmpupperbounddominatingcol1 = SCIPinfinity(scip);
1691 tmpupperbounddominatingcol2 = tmpupperbounddominatingcol1;
1692 tmpwclowerbounddominatingcol1 = -SCIPinfinity(scip);
1693 tmpwclowerbounddominatingcol2 = tmpwclowerbounddominatingcol1;
1694 tmplowerbounddominatingcol1 = -SCIPinfinity(scip);
1695 tmplowerbounddominatingcol2 = tmplowerbounddominatingcol1;
1696 tmpwcupperbounddominatingcol1 = SCIPinfinity(scip);
1697 tmpwcupperbounddominatingcol2 = tmpwcupperbounddominatingcol1;
1698
1699 /* initialize temporary bounds of dominated variable */
1700 tmpupperbounddominatedcol1 = SCIPinfinity(scip);
1701 tmpupperbounddominatedcol2 = tmpupperbounddominatedcol1;
1702 tmpwclowerbounddominatedcol1 = -SCIPinfinity(scip);
1703 tmpwclowerbounddominatedcol2 = tmpwclowerbounddominatedcol1;
1704 tmplowerbounddominatedcol1 = -SCIPinfinity(scip);
1705 tmplowerbounddominatedcol2 = tmplowerbounddominatedcol1;
1706 tmpwcupperbounddominatedcol1 = SCIPinfinity(scip);
1707 tmpwcupperbounddominatedcol2 = tmpwcupperbounddominatedcol1;
1708
1709 /* compare rows of this column pair */
1710 while( (col1domcol2 || col2domcol1) && (r1 < nrows1 || r2 < nrows2) )
1711 {
1712 assert((r1 >= nrows1-1) || (rows1[r1] < rows1[r1+1]));
1713 assert((r2 >= nrows2-1) || (rows2[r2] < rows2[r2+1]));
1714
1715 /* there is a nonredundant row containing column 1 but not column 2 */
1716 if( r1 < nrows1 && (r2 == nrows2 || rows1[r1] < rows2[r2]) )
1717 {
1718 /* dominance depends on the type of relation */
1719 if( !SCIPmatrixIsRowRhsInfinity(matrix, rows1[r1]) )
1720 {
1721 /* no dominance relation for equations or ranged rows */
1722 col2domcol1 = FALSE;
1723 col1domcol2 = FALSE;
1724 }
1725 else
1726 {
1727 /* >= relation, larger coefficients dominate smaller ones */
1728 if( vals1[r1] > 0.0 )
1729 col2domcol1 = FALSE;
1730 else
1731 col1domcol2 = FALSE;
1732 }
1733
1734 r1++;
1735 }
1736 /* there is a nonredundant row containing column 2, but not column 1 */
1737 else if( r2 < nrows2 && (r1 == nrows1 || rows1[r1] > rows2[r2]) )
1738 {
1739 /* dominance depends on the type of relation */
1740 if( !SCIPmatrixIsRowRhsInfinity(matrix, rows2[r2]) )
1741 {
1742 /* no dominance relation for equations or ranged rows */
1743 col2domcol1 = FALSE;
1744 col1domcol2 = FALSE;
1745 }
1746 else
1747 {
1748 /* >= relation, larger coefficients dominate smaller ones */
1749 if( vals2[r2] < 0.0 )
1750 col2domcol1 = FALSE;
1751 else
1752 col1domcol2 = FALSE;
1753 }
1754
1755 r2++;
1756 }
1757 /* if both columns appear in a common row, compare the coefficients */
1758 else
1759 {
1760 assert(r1 < nrows1 && r2 < nrows2);
1761 assert(rows1[r1] == rows2[r2]);
1762
1763 /* if both columns are binary variables we check if they have a common clique
1764 * and do not calculate any bounds
1765 */
1766 if( onlybinvars && !onlyoneone )
1767 {
1768 if( vals1[r1] < 0 && vals2[r2] < 0 )
1769 {
1770 if( (SCIPmatrixGetRowNMaxActPosInf(matrix, rows1[r1]) + SCIPmatrixGetRowNMaxActNegInf(matrix, rows1[r1]) == 0)
1771 && SCIPisFeasLE(scip, SCIPmatrixGetRowMaxActivity(matrix, rows1[r1]) + MAX(vals1[r1], vals2[r2]),
1772 SCIPmatrixGetRowLhs(matrix, rows1[r1])) )
1773 {
1774 onlyoneone = TRUE;
1775 }
1776 }
1777
1778 if( !onlyoneone && !SCIPmatrixIsRowRhsInfinity(matrix, rows1[r1]) )
1779 {
1780 if ( vals1[r1] > 0 && vals2[r2] > 0 )
1781 {
1782 if( (SCIPmatrixGetRowNMinActPosInf(matrix, rows1[r1]) + SCIPmatrixGetRowNMinActNegInf(matrix, rows1[r1]) == 0)
1783 && SCIPisFeasGE(scip, SCIPmatrixGetRowMinActivity(matrix, rows1[r1]) + MIN(vals1[r1], vals2[r2]),
1784 SCIPmatrixGetRowRhs(matrix, rows1[r1])) )
1785 {
1786 onlyoneone = TRUE;
1787 }
1788 }
1789 }
1790
1791 if( onlyoneone )
1792 {
1793 /* reset bounds */
1794 tmpupperbounddominatingcol1 = SCIPinfinity(scip);
1795 tmpupperbounddominatingcol2 = tmpupperbounddominatingcol1;
1796 tmpwclowerbounddominatingcol1 = -SCIPinfinity(scip);
1797 tmpwclowerbounddominatingcol2 = tmpwclowerbounddominatingcol1;
1798 tmplowerbounddominatingcol1 = -SCIPinfinity(scip);
1799 tmplowerbounddominatingcol2 = tmplowerbounddominatingcol1;
1800 tmpwcupperbounddominatingcol1 = SCIPinfinity(scip);
1801 tmpwcupperbounddominatingcol2 = tmpwcupperbounddominatingcol1;
1802
1803 tmpupperbounddominatedcol1 = SCIPinfinity(scip);
1804 tmpupperbounddominatedcol2 = tmpupperbounddominatedcol1;
1805 tmpwclowerbounddominatedcol1 = -SCIPinfinity(scip);
1806 tmpwclowerbounddominatedcol2 = tmpwclowerbounddominatedcol1;
1807 tmplowerbounddominatedcol1 = -SCIPinfinity(scip);
1808 tmplowerbounddominatedcol2 = tmplowerbounddominatedcol1;
1809 tmpwcupperbounddominatedcol1 = SCIPinfinity(scip);
1810 tmpwcupperbounddominatedcol2 = tmpwcupperbounddominatedcol1;
1811 }
1812 }
1813
1814 /* dominance depends on the type of inequality */
1815 if( !SCIPmatrixIsRowRhsInfinity(matrix, rows1[r1]) )
1816 {
1817 /* no dominance relation if coefficients differ for equations or ranged rows */
1818 if( !SCIPisEQ(scip, vals1[r1], vals2[r2]) )
1819 {
1820 col2domcol1 = FALSE;
1821 col1domcol2 = FALSE;
1822 }
1823 }
1824 else
1825 {
1826 /* >= relation, larger coefficients dominate smaller ones */
1827 if( vals1[r1] > vals2[r2] )
1828 col2domcol1 = FALSE;
1829 else if( vals1[r1] < vals2[r2] )
1830 col1domcol2 = FALSE;
1831 }
1832
1833 /* we do not use bound calulations if two binary variable are in one common clique.
1834 * for the other cases we claim the same sign for the coefficients to
1835 * achieve monotonically decreasing predictive bound functions.
1836 */
1837 if( !onlyoneone &&
1838 ((vals1[r1] < 0 && vals2[r2] < 0) || (vals1[r1] > 0 && vals2[r2] > 0)) )
1839 {
1840 if( col1domcol2 )
1841 {
1842 /* update bounds of dominating variable for column 1 */
1843 SCIP_CALL( updateBounds(scip, matrix, rows1[r1],
1844 col1, vals1[r1], col2, vals2[r2], TRUE,
1845 &tmpupperbounddominatingcol1, &tmpwclowerbounddominatingcol1,
1846 &tmplowerbounddominatingcol1, &tmpwcupperbounddominatingcol1) );
1847
1848 /* update bounds of dominated variable for column 1 */
1849 SCIP_CALL( updateBounds(scip, matrix, rows1[r1],
1850 col1, vals1[r1], col2, vals2[r2], FALSE,
1851 &tmpupperbounddominatedcol1, &tmpwclowerbounddominatedcol1,
1852 &tmplowerbounddominatedcol1, &tmpwcupperbounddominatedcol1) );
1853 }
1854
1855 if( col2domcol1 )
1856 {
1857 /* update bounds of dominating variable for column 2 */
1858 SCIP_CALL( updateBounds(scip, matrix, rows2[r2],
1859 col2, vals2[r2], col1, vals1[r1], TRUE,
1860 &tmpupperbounddominatingcol2, &tmpwclowerbounddominatingcol2,
1861 &tmplowerbounddominatingcol2, &tmpwcupperbounddominatingcol2) );
1862
1863 /* update bounds of dominated variable for column 2 */
1864 SCIP_CALL( updateBounds(scip, matrix, rows2[r2],
1865 col2, vals2[r2], col1, vals1[r1], FALSE,
1866 &tmpupperbounddominatedcol2, &tmpwclowerbounddominatedcol2,
1867 &tmplowerbounddominatedcol2, &tmpwcupperbounddominatedcol2) );
1868 }
1869 }
1870
1871 r1++;
1872 r2++;
1873 }
1874 }
1875
1876 /* a column is only dominated, if there are no more rows in which it is contained */
1877 col1domcol2 = col1domcol2 && r2 == nrows2;
1878 col2domcol1 = col2domcol1 && r1 == nrows1;
1879
1880 if( !col1domcol2 && !col2domcol1 )
1881 continue;
1882
1883 /* no dominance relation for left equations or ranged rows */
1884 while( r1 < nrows1 )
1885 {
1886 if( !SCIPmatrixIsRowRhsInfinity(matrix, rows1[r1]) )
1887 {
1888 col2domcol1 = FALSE;
1889 col1domcol2 = FALSE;
1890 break;
1891 }
1892 r1++;
1893 }
1894 if( !col1domcol2 && !col2domcol1 )
1895 continue;
1896 while( r2 < nrows2 )
1897 {
1898 if( !SCIPmatrixIsRowRhsInfinity(matrix, rows2[r2]) )
1899 {
1900 col2domcol1 = FALSE;
1901 col1domcol2 = FALSE;
1902 break;
1903 }
1904 r2++;
1905 }
1906
1907 if( col1domcol2 || col2domcol1 )
1908 (*ndomrelations)++;
1909
1910 if( col1domcol2 && col2domcol1 )
1911 {
1912 /* prefer the variable as dominating variable with the greater upper bound */
1913 if( SCIPisGE(scip, SCIPvarGetUbGlobal(varcol1), SCIPvarGetUbGlobal(varcol2)) )
1914 col2domcol1 = FALSE;
1915 else
1916 col1domcol2 = FALSE;
1917 }
1918
1919 /* use dominance relation and clique/bound-information
1920 * to find variable fixings. additionally try to strengthen
1921 * variable bounds by predictive bound strengthening.
1922 */
1923 if( col1domcol2 )
1924 {
1925 SCIP_CALL( findFixings(scip, matrix, varcol1, col1,
1926 tmpupperbounddominatingcol1, tmpwclowerbounddominatingcol1,
1927 tmplowerbounddominatingcol1, tmpwcupperbounddominatingcol1,
1928 varcol2, col2,
1929 varstofix, onlybinvars, onlyoneone, nfixings) );
1930
1931 if( presoldata->predbndstr )
1932 {
1933 SCIP_CALL( predBndStr(scip, varcol1, col1,
1934 tmpupperbounddominatingcol1,
1935 tmplowerbounddominatingcol1, tmpwcupperbounddominatingcol1,
1936 varcol2, col2,
1937 tmpupperbounddominatedcol1, tmpwclowerbounddominatedcol1,
1938 tmplowerbounddominatedcol1,
1939 varstofix, nchgbds) );
1940 }
1941 }
1942 else if( col2domcol1 )
1943 {
1944 SCIP_CALL( findFixings(scip, matrix, varcol2, col2,
1945 tmpupperbounddominatingcol2, tmpwclowerbounddominatingcol2,
1946 tmplowerbounddominatingcol2, tmpwcupperbounddominatingcol2,
1947 varcol1, col1,
1948 varstofix, onlybinvars, onlyoneone, nfixings) );
1949
1950 if( presoldata->predbndstr )
1951 {
1952 SCIP_CALL( predBndStr(scip, varcol2, col2,
1953 tmpupperbounddominatingcol2,
1954 tmplowerbounddominatingcol2, tmpwcupperbounddominatingcol2,
1955 varcol1, col1,
1956 tmpupperbounddominatedcol2, tmpwclowerbounddominatedcol2,
1957 tmplowerbounddominatedcol2,
1958 varstofix, nchgbds) );
1959 }
1960 }
1961 if( varstofix[col1] == FIXATLB )
1962 break;
1963 }
1964 }
1965
1966 return SCIP_OKAY;
1967}
1968
1969
1970/*
1971 * Callback methods of presolver
1972 */
1973
1974
1975/** copy method for constraint handler plugins (called when SCIP copies plugins) */
1976static
1977SCIP_DECL_PRESOLCOPY(presolCopyDomcol)
1978{ /*lint --e{715}*/
1979 assert(scip != NULL);
1980 assert(presol != NULL);
1981 assert(strcmp(SCIPpresolGetName(presol), PRESOL_NAME) == 0);
1982
1983 /* call inclusion method of presolver */
1985
1986 return SCIP_OKAY;
1987}
1988
1989/** destructor of presolver to free user data (called when SCIP is exiting) */
1990static
1991SCIP_DECL_PRESOLFREE(presolFreeDomcol)
1992{ /*lint --e{715}*/
1993 SCIP_PRESOLDATA* presoldata;
1994
1995 /* free presolver data */
1996 presoldata = SCIPpresolGetData(presol);
1997 assert(presoldata != NULL);
1998
1999 SCIPfreeBlockMemory(scip, &presoldata);
2000 SCIPpresolSetData(presol, NULL);
2001
2002 return SCIP_OKAY;
2003}
2004
2005/** execution method of presolver */
2006static
2007SCIP_DECL_PRESOLEXEC(presolExecDomcol)
2008{ /*lint --e{715}*/
2009 SCIP_PRESOLDATA* presoldata;
2010 SCIP_MATRIX* matrix;
2011 SCIP_Bool initialized;
2012 SCIP_Bool complete;
2013 SCIP_Bool infeasible;
2014 int nfixings;
2015 SCIP_Longint ndomrelations;
2016 int v;
2017 int r;
2018 FIXINGDIRECTION* varstofix;
2019 SCIP_Bool* varsprocessed;
2020 int nrows;
2021 int ncols;
2022 int* rowidxsorted;
2023 int* rowsparsity;
2024 int varcount;
2025 int* consearchcols;
2026 int* intsearchcols;
2027 int* binsearchcols;
2028 int nconfill;
2029 int nintfill;
2030 int nbinfill;
2031#ifdef SCIP_DEBUG
2032 int nconvarsfixed = 0;
2033 int nintvarsfixed = 0;
2034 int nbinvarsfixed = 0;
2035#endif
2036 int* pclass;
2037 int* colidx;
2038 int pclassstart;
2039 int pc;
2040 SCIP_Bool* varineq;
2041
2042 assert(result != NULL);
2043 *result = SCIP_DIDNOTRUN;
2044
2046 return SCIP_OKAY;
2047
2048 presoldata = SCIPpresolGetData(presol);
2049 assert(presoldata != NULL);
2050
2051 /* don't run for pure LPs */
2052 if( !presoldata->continuousred && (SCIPgetNBinVars(scip) + SCIPgetNIntVars(scip) == 0) )
2053 return SCIP_OKAY;
2054
2055 *result = SCIP_DIDNOTFIND;
2056
2057 matrix = NULL;
2058 SCIP_CALL( SCIPmatrixCreate(scip, &matrix, TRUE, &initialized, &complete, &infeasible,
2059 naddconss, ndelconss, nchgcoefs, nchgbds, nfixedvars) );
2060
2061 /* if infeasibility was detected during matrix creation, return here */
2062 if( infeasible )
2063 {
2064 if( initialized )
2065 SCIPmatrixFree(scip, &matrix);
2066
2067 *result = SCIP_CUTOFF;
2068 return SCIP_OKAY;
2069 }
2070
2071 if( !initialized )
2072 return SCIP_OKAY;
2073
2074 nfixings = 0;
2075 ndomrelations = 0;
2076 nrows = SCIPmatrixGetNRows(matrix);
2077 ncols = SCIPmatrixGetNColumns(matrix);
2078 assert(SCIPgetNVars(scip) == ncols);
2079
2080 SCIP_CALL( SCIPallocBufferArray(scip, &varstofix, ncols) );
2081 BMSclearMemoryArray(varstofix, ncols);
2082
2083 SCIP_CALL( SCIPallocBufferArray(scip, &varsprocessed, ncols) );
2084
2085 /* do not process columns that do not have all locks represented in the matrix */
2086 for( v = 0; v < ncols; ++v )
2087 varsprocessed[v] = SCIPmatrixUplockConflict(matrix, v) || SCIPmatrixDownlockConflict(matrix, v);
2088
2089 SCIP_CALL( SCIPallocBufferArray(scip, &consearchcols, ncols) );
2090 SCIP_CALL( SCIPallocBufferArray(scip, &intsearchcols, ncols) );
2091 SCIP_CALL( SCIPallocBufferArray(scip, &binsearchcols, ncols) );
2092
2093 SCIP_CALL( SCIPallocBufferArray(scip, &rowidxsorted, nrows) );
2094 SCIP_CALL( SCIPallocBufferArray(scip, &rowsparsity, nrows) );
2095 for( r = 0; r < nrows; ++r )
2096 {
2097 rowidxsorted[r] = r;
2098 rowsparsity[r] = SCIPmatrixGetRowNNonzs(matrix, r);
2099 }
2100
2101 SCIP_CALL( SCIPallocBufferArray(scip, &pclass, ncols) );
2102 SCIP_CALL( SCIPallocBufferArray(scip, &colidx, ncols) );
2103 SCIP_CALL( SCIPallocBufferArray(scip, &varineq, ncols) );
2104 for( v = 0; v < ncols; v++ )
2105 {
2106 colidx[v] = v;
2107 varineq[v] = FALSE;
2108 }
2109
2110 /* init pair comparision control */
2111 presoldata->numcurrentpairs = presoldata->nummaxpairs;
2112
2113 varcount = 0;
2114
2115 if( (presoltiming & SCIP_PRESOLTIMING_EXHAUSTIVE) != 0 )
2116 {
2117 /* 1.stage: search dominance relations of parallel columns
2118 * within equalities and ranged rows
2119 */
2120 SCIP_CALL( detectParallelCols(scip, matrix, pclass, varineq) );
2121 SCIPsortIntInt(pclass, colidx, ncols);
2122
2123 pc = 0;
2124 while( pc < ncols )
2125 {
2126 int varidx;
2127
2128 varidx = 0;
2129 nconfill = 0;
2130 nintfill = 0;
2131 nbinfill = 0;
2132
2133 pclassstart = pclass[pc];
2134 while( pc < ncols && pclassstart == pclass[pc] )
2135 {
2136 SCIP_VAR* var;
2137
2138 varidx = colidx[pc];
2139 var = SCIPmatrixGetVar(matrix, varidx);
2140
2141 /* we only regard variables which were not processed yet and
2142 * are present within equalities or ranged rows
2143 */
2144 if( !varsprocessed[varidx] && varineq[varidx] )
2145 {
2146 /* we search only for dominance relations between the same variable type */
2148 {
2149 consearchcols[nconfill++] = varidx;
2150 }
2151 else if( SCIPvarIsBinary(var) )
2152 {
2153 binsearchcols[nbinfill++] = varidx;
2154 }
2155 else
2156 {
2158 intsearchcols[nintfill++] = varidx;
2159 }
2160 }
2161 ++pc;
2162 }
2163
2164 /* continuous variables */
2165 if( nconfill > 1 && presoldata->continuousred )
2166 {
2167 SCIP_CALL( findDominancePairs(scip, matrix, presoldata, consearchcols, nconfill, FALSE,
2168 varstofix, &nfixings, &ndomrelations, nchgbds) );
2169
2170 for( v = 0; v < nconfill; ++v )
2171 varsprocessed[consearchcols[v]] = TRUE;
2172
2173 varcount += nconfill;
2174 }
2175 else if( nconfill == 1 )
2176 {
2177 if( varineq[varidx] )
2178 varsprocessed[consearchcols[0]] = TRUE;
2179 }
2180
2181 /* integer and impl-integer variables */
2182 if( nintfill > 1 )
2183 {
2184 SCIP_CALL( findDominancePairs(scip, matrix, presoldata, intsearchcols, nintfill, FALSE,
2185 varstofix, &nfixings, &ndomrelations, nchgbds) );
2186
2187 for( v = 0; v < nintfill; ++v )
2188 varsprocessed[intsearchcols[v]] = TRUE;
2189
2190 varcount += nintfill;
2191 }
2192 else if( nintfill == 1 )
2193 {
2194 if( varineq[varidx] )
2195 varsprocessed[intsearchcols[0]] = TRUE;
2196 }
2197
2198 /* binary variables */
2199 if( nbinfill > 1 )
2200 {
2201 SCIP_CALL( findDominancePairs(scip, matrix, presoldata, binsearchcols, nbinfill, TRUE,
2202 varstofix, &nfixings, &ndomrelations, nchgbds) );
2203
2204 for( v = 0; v < nbinfill; ++v )
2205 varsprocessed[binsearchcols[v]] = TRUE;
2206
2207 varcount += nbinfill;
2208 }
2209 else if( nbinfill == 1 )
2210 {
2211 if( varineq[varidx] )
2212 varsprocessed[binsearchcols[0]] = TRUE;
2213 }
2214
2215 if( varcount >= ncols )
2216 break;
2217 }
2218
2219 /* 2.stage: search dominance relations for the remaining columns
2220 * by increasing row-sparsity
2221 */
2222 SCIPsortIntInt(rowsparsity, rowidxsorted, nrows);
2223
2224 for( r = 0; r < nrows; ++r )
2225 {
2226 int rowidx;
2227 int* rowpnt;
2228 int* rowend;
2229
2230 /* break if the time limit was reached; since the check is expensive,
2231 * we only check all 1000 constraints
2232 */
2233 if( (r % 1000 == 0) && SCIPisStopped(scip) )
2234 break;
2235
2236 rowidx = rowidxsorted[r];
2237 rowpnt = SCIPmatrixGetRowIdxPtr(matrix, rowidx);
2238 rowend = rowpnt + SCIPmatrixGetRowNNonzs(matrix, rowidx);
2239
2240 if( SCIPmatrixGetRowNNonzs(matrix, rowidx) == 1 )
2241 continue;
2242
2243 nconfill = 0;
2244 nintfill = 0;
2245 nbinfill = 0;
2246
2247 for( ; rowpnt < rowend; rowpnt++ )
2248 {
2249 if( !(varsprocessed[*rowpnt]) )
2250 {
2251 int varidx;
2252 SCIP_VAR* var;
2253
2254 varidx = *rowpnt;
2255 var = SCIPmatrixGetVar(matrix, varidx);
2256
2257 /* we search only for dominance relations between the same variable type */
2259 {
2260 consearchcols[nconfill++] = varidx;
2261 }
2262 else if( SCIPvarIsBinary(var) )
2263 {
2264 binsearchcols[nbinfill++] = varidx;
2265 }
2266 else
2267 {
2269 intsearchcols[nintfill++] = varidx;
2270 }
2271 }
2272 }
2273
2274 /* continuous variables */
2275 if( nconfill > 1 && presoldata->continuousred )
2276 {
2277 SCIP_CALL( findDominancePairs(scip, matrix, presoldata, consearchcols, nconfill, FALSE,
2278 varstofix, &nfixings, &ndomrelations, nchgbds) );
2279
2280 for( v = 0; v < nconfill; ++v )
2281 varsprocessed[consearchcols[v]] = TRUE;
2282
2283 varcount += nconfill;
2284 }
2285
2286 /* integer and impl-integer variables */
2287 if( nintfill > 1 )
2288 {
2289 SCIP_CALL( findDominancePairs(scip, matrix, presoldata, intsearchcols, nintfill, FALSE,
2290 varstofix, &nfixings, &ndomrelations, nchgbds) );
2291
2292 for( v = 0; v < nintfill; ++v )
2293 varsprocessed[intsearchcols[v]] = TRUE;
2294
2295 varcount += nintfill;
2296 }
2297
2298 /* binary variables */
2299 if( nbinfill > 1 )
2300 {
2301 SCIP_CALL( findDominancePairs(scip, matrix, presoldata, binsearchcols, nbinfill, TRUE,
2302 varstofix, &nfixings, &ndomrelations, nchgbds) );
2303
2304 for( v = 0; v < nbinfill; ++v )
2305 varsprocessed[binsearchcols[v]] = TRUE;
2306
2307 varcount += nbinfill;
2308 }
2309
2310 if( varcount >= ncols )
2311 break;
2312 }
2313 }
2314
2315 if( nfixings > 0 )
2316 {
2317 int oldnfixedvars;
2318
2319 oldnfixedvars = *nfixedvars;
2320
2321 for( v = ncols - 1; v >= 0; --v )
2322 {
2323 SCIP_Bool fixed;
2324 SCIP_VAR* var;
2325
2326 var = SCIPmatrixGetVar(matrix,v);
2327
2328 /* there should be no fixings for variables with lock conflicts,
2329 * since they where marked as processed and therefore should be skipped
2330 */
2331 assert(varstofix[v] == NOFIX || (!SCIPmatrixUplockConflict(matrix, v) && !SCIPmatrixDownlockConflict(matrix, v)) );
2332
2333 if( varstofix[v] == FIXATLB )
2334 {
2335 SCIP_Real lb;
2336
2337 lb = SCIPvarGetLbGlobal(var);
2339
2340 /* fix at lower bound */
2341 SCIP_CALL( SCIPfixVar(scip, var, lb, &infeasible, &fixed) );
2342 if( infeasible )
2343 {
2344 SCIPdebugMsg(scip, " -> infeasible fixing\n");
2345 *result = SCIP_CUTOFF;
2346
2347 break;
2348 }
2349 assert(fixed);
2350 (*nfixedvars)++;
2351
2352#ifdef SCIP_DEBUG
2354 nconvarsfixed++;
2355 else if( SCIPvarIsBinary(var) )
2356 nbinvarsfixed++;
2357 else
2358 nintvarsfixed++;
2359#endif
2360 }
2361 else if( varstofix[v] == FIXATUB )
2362 {
2363 SCIP_Real ub;
2364
2365 ub = SCIPvarGetUbGlobal(var);
2367
2368 /* fix at upper bound */
2369 SCIP_CALL( SCIPfixVar(scip, var, ub, &infeasible, &fixed) );
2370 if( infeasible )
2371 {
2372 SCIPdebugMsg(scip, " -> infeasible fixing\n");
2373 *result = SCIP_CUTOFF;
2374
2375 break;
2376 }
2377 assert(fixed);
2378 (*nfixedvars)++;
2379
2380#ifdef SCIP_DEBUG
2382 nconvarsfixed++;
2383 else if( SCIPvarIsBinary(var) )
2384 nbinvarsfixed++;
2385 else
2386 nintvarsfixed++;
2387#endif
2388 }
2389 }
2390
2391 if( *result != SCIP_CUTOFF && *nfixedvars > oldnfixedvars )
2392 *result = SCIP_SUCCESS;
2393 }
2394
2395 SCIPfreeBufferArray(scip, &varineq);
2396 SCIPfreeBufferArray(scip, &colidx);
2397 SCIPfreeBufferArray(scip, &pclass);
2398 SCIPfreeBufferArray(scip, &rowsparsity);
2399 SCIPfreeBufferArray(scip, &rowidxsorted);
2400 SCIPfreeBufferArray(scip, &binsearchcols);
2401 SCIPfreeBufferArray(scip, &intsearchcols);
2402 SCIPfreeBufferArray(scip, &consearchcols);
2403 SCIPfreeBufferArray(scip, &varsprocessed);
2404 SCIPfreeBufferArray(scip, &varstofix);
2405
2406#ifdef SCIP_DEBUG
2407 if( (nconvarsfixed + nintvarsfixed + nbinvarsfixed) > 0 )
2408 {
2409 SCIPdebugMsg(scip, "### %d vars [%" SCIP_LONGINT_FORMAT " dom] => fixed [cont: %d, int: %d, bin: %d], %scutoff detected\n",
2410 ncols, ndomrelations, nconvarsfixed, nintvarsfixed, nbinvarsfixed, (*result != SCIP_CUTOFF) ? "no " : "");
2411 }
2412#endif
2413
2414 SCIPmatrixFree(scip, &matrix);
2415
2416 return SCIP_OKAY;
2417}
2418
2419/*
2420 * presolver specific interface methods
2421 */
2422
2423/** creates the domcol presolver and includes it in SCIP */
2425 SCIP* scip /**< SCIP data structure */
2426 )
2427{
2428 SCIP_PRESOLDATA* presoldata;
2429 SCIP_PRESOL* presol;
2430
2431 /* create domcol presolver data */
2432 SCIP_CALL( SCIPallocBlockMemory(scip, &presoldata) );
2433
2434 /* include presolver */
2436 PRESOL_TIMING, presolExecDomcol, presoldata) );
2437 SCIP_CALL( SCIPsetPresolCopy(scip, presol, presolCopyDomcol) );
2438 SCIP_CALL( SCIPsetPresolFree(scip, presol, presolFreeDomcol) );
2439
2441 "presolving/domcol/numminpairs",
2442 "minimal number of pair comparisons",
2443 &presoldata->numminpairs, FALSE, DEFAULT_NUMMINPAIRS, 100, DEFAULT_NUMMAXPAIRS, NULL, NULL) );
2444
2446 "presolving/domcol/nummaxpairs",
2447 "maximal number of pair comparisons",
2448 &presoldata->nummaxpairs, FALSE, DEFAULT_NUMMAXPAIRS, DEFAULT_NUMMINPAIRS, 1000000000, NULL, NULL) );
2449
2451 "presolving/domcol/predbndstr",
2452 "should predictive bound strengthening be applied?",
2453 &presoldata->predbndstr, FALSE, DEFAULT_PREDBNDSTR, NULL, NULL) );
2454
2456 "presolving/domcol/continuousred",
2457 "should reductions for continuous variables be performed?",
2458 &presoldata->continuousred, FALSE, DEFAULT_CONTINUOUS_RED, NULL, NULL) );
2459
2460 return SCIP_OKAY;
2461}
SCIP_Real * r
Definition: circlepacking.c:59
#define NULL
Definition: def.h:262
#define SCIP_Longint
Definition: def.h:157
#define SCIP_Bool
Definition: def.h:91
#define MIN(x, y)
Definition: def.h:238
#define SCIP_Real
Definition: def.h:172
#define TRUE
Definition: def.h:93
#define FALSE
Definition: def.h:94
#define MAX(x, y)
Definition: def.h:234
#define SCIP_LONGINT_FORMAT
Definition: def.h:164
#define SCIPABORT()
Definition: def.h:341
#define SCIP_CALL(x)
Definition: def.h:369
SCIP_Bool SCIPisStopped(SCIP *scip)
Definition: scip_general.c:734
int SCIPgetNIntVars(SCIP *scip)
Definition: scip_prob.c:2082
int SCIPgetNVars(SCIP *scip)
Definition: scip_prob.c:1992
int SCIPgetNBinVars(SCIP *scip)
Definition: scip_prob.c:2037
#define SCIPdebugMsgPrint
Definition: scip_message.h:79
#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 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 SCIPincludePresolDomcol(SCIP *scip)
#define SCIPallocBufferArray(scip, ptr, num)
Definition: scip_mem.h:124
#define SCIPfreeBufferArray(scip, ptr)
Definition: scip_mem.h:136
#define SCIPfreeBlockMemory(scip, ptr)
Definition: scip_mem.h:108
#define SCIPallocBlockMemory(scip, ptr)
Definition: scip_mem.h:89
void SCIPpresolSetData(SCIP_PRESOL *presol, SCIP_PRESOLDATA *presoldata)
Definition: presol.c:533
SCIP_PRESOLDATA * SCIPpresolGetData(SCIP_PRESOL *presol)
Definition: presol.c:523
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
const char * SCIPpresolGetName(SCIP_PRESOL *presol)
Definition: presol.c:610
SCIP_Bool SCIPisFeasGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Real SCIPinfinity(SCIP *scip)
SCIP_Bool SCIPisGE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPisPositive(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Real SCIPfloor(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisInfinity(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisFeasLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPisFeasIntegral(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisNegative(SCIP *scip, SCIP_Real val)
SCIP_Real SCIPceil(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPisZero(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisLT(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPvarIsBinary(SCIP_VAR *var)
Definition: var.c:17617
SCIP_RETCODE SCIPchgVarLb(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound)
Definition: scip_var.c:4799
SCIP_VARSTATUS SCIPvarGetStatus(SCIP_VAR *var)
Definition: var.c:17556
SCIP_RETCODE SCIPchgVarUb(SCIP *scip, SCIP_VAR *var, SCIP_Real newbound)
Definition: scip_var.c:4889
SCIP_Real SCIPvarGetObj(SCIP_VAR *var)
Definition: var.c:17944
SCIP_VARTYPE SCIPvarGetType(SCIP_VAR *var)
Definition: var.c:17602
SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:18106
const char * SCIPvarGetName(SCIP_VAR *var)
Definition: var.c:17437
SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:18096
SCIP_RETCODE SCIPfixVar(SCIP *scip, SCIP_VAR *var, SCIP_Real fixedval, SCIP_Bool *infeasible, SCIP_Bool *fixed)
Definition: scip_var.c:8381
SCIP_Bool SCIPvarsHaveCommonClique(SCIP_VAR *var1, SCIP_Bool value1, SCIP_VAR *var2, SCIP_Bool value2, SCIP_Bool regardimplics)
Definition: var.c:11493
SCIP_Bool SCIPallowStrongDualReds(SCIP *scip)
Definition: scip_var.c:8722
void SCIPsortIntInt(int *intarray1, int *intarray2, int len)
void SCIPsortRealInt(SCIP_Real *realarray, int *intarray, int len)
void SCIPsortIntIntReal(int *intarray1, int *intarray2, SCIP_Real *realarray, int len)
SCIP_Bool SCIPmatrixUplockConflict(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1872
int SCIPmatrixGetRowNMinActNegInf(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1812
int * SCIPmatrixGetColIdxPtr(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1580
int SCIPmatrixGetRowNNonzs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1708
int SCIPmatrixGetColNNonzs(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1592
SCIP_Bool SCIPmatrixIsRowRhsInfinity(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1766
SCIP_Real SCIPmatrixGetRowMaxActivity(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1800
SCIP_Real SCIPmatrixGetRowLhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1742
SCIP_Real * SCIPmatrixGetRowValPtr(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1684
SCIP_Bool SCIPmatrixDownlockConflict(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1884
SCIP_Real SCIPmatrixGetRowRhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1754
SCIP_Real * SCIPmatrixGetColValPtr(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1568
int SCIPmatrixGetRowNMinActPosInf(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1824
SCIP_RETCODE SCIPmatrixCreate(SCIP *scip, SCIP_MATRIX **matrixptr, SCIP_Bool onlyifcomplete, SCIP_Bool *initialized, SCIP_Bool *complete, SCIP_Bool *infeasible, int *naddconss, int *ndelconss, int *nchgcoefs, int *nchgbds, int *nfixedvars)
Definition: matrix.c:454
int SCIPmatrixGetNColumns(SCIP_MATRIX *matrix)
Definition: matrix.c:1604
SCIP_Real SCIPmatrixGetRowMinActivity(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1788
void SCIPmatrixFree(SCIP *scip, SCIP_MATRIX **matrix)
Definition: matrix.c:1072
int SCIPmatrixGetRowNMaxActPosInf(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1848
int SCIPmatrixGetRowNMaxActNegInf(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1836
SCIP_VAR * SCIPmatrixGetVar(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1660
int * SCIPmatrixGetRowIdxPtr(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1696
int SCIPmatrixGetNRows(SCIP_MATRIX *matrix)
Definition: matrix.c:1732
memory allocation routines
#define BMSclearMemoryArray(ptr, num)
Definition: memory.h:130
static void getActivityResidualsUpperBound(SCIP *scip, SCIP_MATRIX *matrix, int row, int col, SCIP_Real coef, int upperboundcol, SCIP_Real upperboundcoef, SCIP_Real *minresactivity, SCIP_Real *maxresactivity, SCIP_Bool *success)
static SCIP_RETCODE findFixings(SCIP *scip, SCIP_MATRIX *matrix, SCIP_VAR *dominatingvar, int dominatingidx, SCIP_Real dominatingub, SCIP_Real dominatingwclb, SCIP_Real dominatinglb, SCIP_Real dominatingwcub, SCIP_VAR *dominatedvar, int dominatedidx, FIXINGDIRECTION *varstofix, SCIP_Bool onlybinvars, SCIP_Bool onlyoneone, int *nfixings)
static SCIP_RETCODE predBndStr(SCIP *scip, SCIP_VAR *dominatingvar, int dominatingidx, SCIP_Real dominatingub, SCIP_Real dominatinglb, SCIP_Real dominatingwcub, SCIP_VAR *dominatedvar, int dominatedidx, SCIP_Real dominatedub, SCIP_Real dominatedwclb, SCIP_Real dominatedlb, FIXINGDIRECTION *varstofix, int *nchgbds)
#define PRESOL_NAME
Definition: presol_domcol.c:70
Fixingdirection
@ FIXATUB
@ FIXATLB
@ NOFIX
static SCIP_DECL_PRESOLEXEC(presolExecDomcol)
enum Fixingdirection FIXINGDIRECTION
static SCIP_RETCODE updateBounds(SCIP *scip, SCIP_MATRIX *matrix, int row, int col1, SCIP_Real val1, int col2, SCIP_Real val2, SCIP_Bool predictdominating, SCIP_Real *upperbound, SCIP_Real *wclowerbound, SCIP_Real *lowerbound, SCIP_Real *wcupperbound)
#define DEFAULT_CONTINUOUS_RED
Definition: presol_domcol.c:80
#define DEFAULT_NUMMINPAIRS
Definition: presol_domcol.c:76
#define DEFAULT_NUMMAXPAIRS
Definition: presol_domcol.c:77
#define PRESOL_PRIORITY
Definition: presol_domcol.c:72
static SCIP_DECL_PRESOLFREE(presolFreeDomcol)
static SCIP_RETCODE detectParallelCols(SCIP *scip, SCIP_MATRIX *matrix, int *pclass, SCIP_Bool *varineq)
static SCIP_RETCODE calcVarBoundsDominated(SCIP *scip, SCIP_MATRIX *matrix, int row, int coldominating, SCIP_Real valdominating, int coldominated, SCIP_Real valdominated, SCIP_Bool *ubcalculated, SCIP_Real *calculatedub, SCIP_Bool *wclbcalculated, SCIP_Real *calculatedwclb, SCIP_Bool *lbcalculated, SCIP_Real *calculatedlb, SCIP_Bool *wcubcalculated, SCIP_Real *calculatedwcub)
static SCIP_RETCODE calcVarBoundsDominating(SCIP *scip, SCIP_MATRIX *matrix, int row, int coldominating, SCIP_Real valdominating, int coldominated, SCIP_Real valdominated, SCIP_Bool *ubcalculated, SCIP_Real *calculatedub, SCIP_Bool *wclbcalculated, SCIP_Real *calculatedwclb, SCIP_Bool *lbcalculated, SCIP_Real *calculatedlb, SCIP_Bool *wcubcalculated, SCIP_Real *calculatedwcub)
static SCIP_RETCODE findDominancePairs(SCIP *scip, SCIP_MATRIX *matrix, SCIP_PRESOLDATA *presoldata, int *searchcols, int searchsize, SCIP_Bool onlybinvars, FIXINGDIRECTION *varstofix, int *nfixings, SCIP_Longint *ndomrelations, int *nchgbds)
#define PRESOL_MAXROUNDS
Definition: presol_domcol.c:73
static SCIP_DECL_PRESOLCOPY(presolCopyDomcol)
#define PRESOL_TIMING
Definition: presol_domcol.c:74
static void getActivityResidualsLowerBound(SCIP *scip, SCIP_MATRIX *matrix, int row, int col, SCIP_Real coef, int lowerboundcol, SCIP_Real lowerboundcoef, SCIP_Real *minresactivity, SCIP_Real *maxresactivity, SCIP_Bool *success)
#define DEFAULT_PREDBNDSTR
Definition: presol_domcol.c:79
#define PRESOL_DESC
Definition: presol_domcol.c:71
dominated column presolver
public methods for matrix
public methods for message output
#define SCIPerrorMessage
Definition: pub_message.h:64
methods for sorting joint arrays of various types
public methods for presolvers
public methods for problem variables
static SCIP_RETCODE printRow(SCIP *scip, FZNOUTPUT *fznoutput, const char *type, SCIP_VAR **vars, SCIP_Real *vals, int nvars, SCIP_Real rhs, SCIP_Bool hasfloats)
Definition: reader_fzn.c:3996
general public methods
public methods for memory management
public methods for message handling
public methods for nonlinear relaxation
public methods for numerical tolerances
public methods for SCIP parameter handling
public methods for presolving plugins
public methods for variable pricer plugins
public methods for global and local (sub)problems
public methods for the probing mode
public methods for SCIP variables
struct SCIP_PresolData SCIP_PRESOLDATA
Definition: type_presol.h:51
@ SCIP_DIDNOTRUN
Definition: type_result.h:42
@ SCIP_CUTOFF
Definition: type_result.h:48
@ SCIP_DIDNOTFIND
Definition: type_result.h:44
@ SCIP_SUCCESS
Definition: type_result.h:58
@ SCIP_INVALIDDATA
Definition: type_retcode.h:52
@ SCIP_OKAY
Definition: type_retcode.h:42
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:63
#define SCIP_PRESOLTIMING_EXHAUSTIVE
Definition: type_timing.h:54
@ SCIP_VARTYPE_INTEGER
Definition: type_var.h:63
@ SCIP_VARTYPE_CONTINUOUS
Definition: type_var.h:71
@ SCIP_VARTYPE_IMPLINT
Definition: type_var.h:64
@ SCIP_VARTYPE_BINARY
Definition: type_var.h:62
@ SCIP_VARSTATUS_MULTAGGR
Definition: type_var.h:54