Scippy

SCIP

Solving Constraint Integer Programs

matrix.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 matrix.c
26 * @ingroup OTHER_CFILES
27 * @brief methods for MIP matrix data structure
28 * @author Dieter Weninger
29 * @author Gerald Gamrath
30 *
31 * The MIP matrix is organized as sparse data structure in row and
32 * and column major format.
33 *
34 * @todo disregard relaxation-only variables in lock check and don't copy them to the matrix
35 */
36
37/*---+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
38
40#include "scip/cons_knapsack.h"
41#include "scip/cons_linear.h"
43#include "scip/cons_logicor.h"
44#include "scip/cons_setppc.h"
45#include "scip/cons_varbound.h"
46#include "scip/pub_matrix.h"
47#include "scip/pub_cons.h"
48#include "scip/pub_message.h"
49#include "scip/pub_misc_sort.h"
50#include "scip/pub_var.h"
51#include "scip/scip_cons.h"
52#include "scip/scip_exact.h"
53#include "scip/scip_general.h"
54#include "scip/scip_mem.h"
55#include "scip/scip_message.h"
56#include "scip/scip_numerics.h"
57#include "scip/scip_pricer.h"
58#include "scip/scip_prob.h"
59#include "scip/scip_var.h"
60#include "scip/struct_matrix.h"
61#include "scip/rational.h"
62#include <string.h>
63
64/*
65 * private functions
66 */
67
68/** transforms given variables, scalars and constant to the corresponding active variables, scalars and constant */
69static
71 SCIP* scip, /**< SCIP instance */
72 SCIP_VAR*** vars, /**< vars array to get active variables for */
73 SCIP_Real** scalars, /**< scalars a_1, ..., a_n in linear sum a_1*x_1 + ... + a_n*x_n + c */
74 int* nvars, /**< pointer to number of variables and values in vars and vals array */
75 SCIP_Real* constant /**< pointer to constant c in linear sum a_1*x_1 + ... + a_n*x_n + c */
76 )
77{
78 int requiredsize;
79
80 assert(scip != NULL);
81 assert(vars != NULL);
82 assert(scalars != NULL);
83 assert(*vars != NULL);
84 assert(*scalars != NULL);
85 assert(nvars != NULL);
86 assert(constant != NULL);
87
88 SCIP_CALL( SCIPgetProbvarLinearSum(scip, *vars, *scalars, nvars, *nvars, constant, &requiredsize) );
89
90 if( requiredsize > *nvars )
91 {
92 SCIP_CALL( SCIPreallocBufferArray(scip, vars, requiredsize) );
94
95 /* call function a second time with enough memory */
96 SCIP_CALL( SCIPgetProbvarLinearSum(scip, *vars, *scalars, nvars, requiredsize, constant, &requiredsize) );
97 }
98 assert(requiredsize == *nvars);
99
100 return SCIP_OKAY;
101}
102
103/** transforms given variables, scalars and constant to the corresponding active variables, scalars and constant */
104static
106 SCIP* scip, /**< SCIP instance */
107 SCIP_VAR*** vars, /**< vars array to get active variables for */
108 SCIP_RATIONAL** scalars, /**< scalars a_1, ..., a_n in linear sum a_1*x_1 + ... + a_n*x_n + c */
109 int* nvars, /**< pointer to number of variables and values in vars and vals array */
110 SCIP_RATIONAL* constant /**< pointer to constant c in linear sum a_1*x_1 + ... + a_n*x_n + c */
111 )
112{
113 int requiredsize;
114
115 assert(scip != NULL);
116 assert(vars != NULL);
117 assert(scalars != NULL);
118 assert(*vars != NULL);
119 assert(*scalars != NULL);
120 assert(nvars != NULL);
121 assert(constant != NULL);
122
123 SCIP_CALL( SCIPgetProbvarLinearSumExact(scip, *vars, scalars, nvars, *nvars, constant, &requiredsize, TRUE) );
124
125 if( requiredsize > *nvars )
126 {
127 SCIP_CALL( SCIPreallocBufferArray(scip, vars, requiredsize) );
129
130 /* call function a second time with enough memory */
131 SCIP_CALL( SCIPgetProbvarLinearSumExact(scip, *vars, scalars, nvars, requiredsize, constant, &requiredsize, TRUE) );
132 assert(requiredsize <= *nvars);
133 }
134
135 assert(requiredsize <= *nvars);
136
137 return SCIP_OKAY;
138}
139
140/** add one row to the constraint matrix */
141static
143 SCIP* scip, /**< SCIP data structure */
144 SCIP_MATRIX* matrix, /**< constraint matrix */
145 SCIP_VAR** vars, /**< variables of this row */
146 SCIP_Real* vals, /**< coefficients of this row */
147 int nvars, /**< number of variables of this row */
148 SCIP_Real lhs, /**< left hand side */
149 SCIP_Real rhs, /**< right hand side */
150 int maxnnonzsmem, /**< maximal number of fillable elements */
151 SCIP_Bool* rowadded /**< flag indicating if constraint was added to matrix */
152 )
153{
154 int j;
155 int probindex;
156 int rowidx;
157 SCIP_Real factor;
158 SCIP_Bool rangedorequality;
159
160 assert(vars != NULL);
161 assert(vals != NULL);
162
163 rowidx = matrix->nrows;
164 rangedorequality = FALSE;
165
166 if( SCIPisInfinity(scip, -lhs) )
167 {
168 factor = -1.0;
169 matrix->lhs[rowidx] = -rhs;
170 matrix->rhs[rowidx] = SCIPinfinity(scip);
171 matrix->isrhsinfinite[rowidx] = TRUE;
172 }
173 else
174 {
175 factor = 1.0;
176 matrix->lhs[rowidx] = lhs;
177 matrix->rhs[rowidx] = rhs;
178 matrix->isrhsinfinite[rowidx] = SCIPisInfinity(scip, matrix->rhs[rowidx]);
179
180 if( !SCIPisInfinity(scip, rhs) )
181 rangedorequality = TRUE;
182 }
183
184 if(SCIPisInfinity(scip, -matrix->lhs[rowidx]))
185 {
186 /* ignore redundant constraint */
187 *rowadded = FALSE;
188 return SCIP_OKAY;
189 }
190
191 matrix->rowmatbeg[rowidx] = matrix->nnonzs;
192
193 /* = or ranged */
194 if( rangedorequality )
195 {
196 assert(factor > 0);
197
198 for( j = 0; j < nvars; j++ )
199 {
200 assert(maxnnonzsmem > matrix->nnonzs);
201
202 /* ignore variables with very small coefficients */
203 if( SCIPisZero(scip, vals[j]) )
204 continue;
205
206 matrix->rowmatval[matrix->nnonzs] = factor * vals[j];
207 probindex = SCIPvarGetProbindex(vars[j]);
208 assert(matrix->vars[probindex] == vars[j]);
209
210 matrix->nuplocks[probindex]++;
211 matrix->ndownlocks[probindex]++;
212
213 assert(0 <= probindex && probindex < matrix->ncols);
214 matrix->rowmatind[matrix->nnonzs] = probindex;
215
216 (matrix->nnonzs)++;
217 }
218 }
219 /* >= or <= */
220 else
221 {
222 for( j = 0; j < nvars; j++ )
223 {
224 assert(maxnnonzsmem > matrix->nnonzs);
225
226 /* ignore variables with very small coefficients */
227 if( SCIPisZero(scip, vals[j]) )
228 continue;
229
230 /* due to the factor, <= constraints will be transfered to >= */
231 matrix->rowmatval[matrix->nnonzs] = factor * vals[j];
232 probindex = SCIPvarGetProbindex(vars[j]);
233 assert(matrix->vars[probindex] == vars[j]);
234
235 if( matrix->rowmatval[matrix->nnonzs] > 0 )
236 matrix->ndownlocks[probindex]++;
237 else
238 {
239 assert(matrix->rowmatval[matrix->nnonzs] < 0);
240 matrix->nuplocks[probindex]++;
241 }
242
243 assert(0 <= probindex && probindex < matrix->ncols);
244 matrix->rowmatind[matrix->nnonzs] = probindex;
245
246 (matrix->nnonzs)++;
247 }
248 }
249
250 matrix->rowmatcnt[rowidx] = matrix->nnonzs - matrix->rowmatbeg[rowidx];
251
252 ++(matrix->nrows);
253 *rowadded = TRUE;
254
255 return SCIP_OKAY;
256}
257
258/** add one row to the constraint matrix */
259static
261 SCIP_MATRIX* matrix, /**< constraint matrix */
262 SCIP_VAR** vars, /**< variables of this row */
263 SCIP_RATIONAL** vals, /**< coefficients of this row */
264 int nvars, /**< number of variables of this row */
265 SCIP_RATIONAL* lhs, /**< left hand side */
266 SCIP_RATIONAL* rhs, /**< right hand side */
267 int maxnnonzsmem, /**< maximal number of fillable elements */
268 SCIP_Bool* rowadded /**< flag indicating if constraint was added to matrix */
269 )
270{
271 int j;
272 int probindex;
273 int rowidx;
274 SCIP_Real factor;
275 SCIP_Bool rangedorequality;
276
277 assert(vars != NULL);
278 assert(vals != NULL);
279
280 rowidx = matrix->nrows;
281 rangedorequality = FALSE;
282
284 {
285 factor = -1.0;
286 SCIPrationalNegate(matrix->matrixvalsexact->lhsexact[rowidx], rhs);
288 matrix->isrhsinfinite[rowidx] = TRUE;
289 }
290 else
291 {
292 factor = 1.0;
293 SCIPrationalSetRational(matrix->matrixvalsexact->lhsexact[rowidx], lhs);
294 SCIPrationalSetRational(matrix->matrixvalsexact->rhsexact[rowidx], rhs);
295 matrix->isrhsinfinite[rowidx] = SCIPrationalIsInfinity(matrix->matrixvalsexact->rhsexact[rowidx]);
296
297 if( !SCIPrationalIsInfinity(rhs) )
298 rangedorequality = TRUE;
299 }
300
302 {
303 /* ignore redundant constraint */
304 *rowadded = FALSE;
305 return SCIP_OKAY;
306 }
307
308 matrix->rowmatbeg[rowidx] = matrix->nnonzs;
309
310 /* = or ranged */
311 if( rangedorequality )
312 {
313 assert(factor > 0);
314
315 for( j = 0; j < nvars; j++ )
316 {
317 assert(maxnnonzsmem > matrix->nnonzs);
318
319 /* ignore variables 0 - coefficients */
320 if( SCIPrationalIsZero(vals[j]) )
321 continue;
322
323 SCIPrationalMultReal(matrix->matrixvalsexact->rowmatvalexact[matrix->nnonzs], vals[j], factor);
324 matrix->rowmatval[matrix->nnonzs] = SCIPrationalGetReal(matrix->matrixvalsexact->rowmatvalexact[matrix->nnonzs]);
325 probindex = SCIPvarGetProbindex(vars[j]);
326 assert(matrix->vars[probindex] == vars[j]);
327
328 matrix->nuplocks[probindex]++;
329 matrix->ndownlocks[probindex]++;
330
331 assert(0 <= probindex && probindex < matrix->ncols);
332 matrix->rowmatind[matrix->nnonzs] = probindex;
333
334 (matrix->nnonzs)++;
335 }
336 }
337 /* >= or <= */
338 else
339 {
340 for( j = 0; j < nvars; j++ )
341 {
342 assert(maxnnonzsmem > matrix->nnonzs);
343
344 /* ignore variables with very small coefficients */
345 if( SCIPrationalIsZero(vals[j]) )
346 continue;
347
348 /* due to the factor, <= constraints will be transfered to >= */
349 SCIPrationalMultReal(matrix->matrixvalsexact->rowmatvalexact[matrix->nnonzs], vals[j], factor);
350 matrix->rowmatval[matrix->nnonzs] = SCIPrationalGetReal(matrix->matrixvalsexact->rowmatvalexact[matrix->nnonzs]);
351 probindex = SCIPvarGetProbindex(vars[j]);
352 assert(matrix->vars[probindex] == vars[j]);
353
355 matrix->ndownlocks[probindex]++;
356 else
357 {
359 matrix->nuplocks[probindex]++;
360 }
361
362 assert(0 <= probindex && probindex < matrix->ncols);
363 matrix->rowmatind[matrix->nnonzs] = probindex;
364
365 (matrix->nnonzs)++;
366 }
367 }
368
369 matrix->rowmatcnt[rowidx] = matrix->nnonzs - matrix->rowmatbeg[rowidx];
370
371 ++(matrix->nrows);
372 *rowadded = TRUE;
373
374 return SCIP_OKAY;
375}
376
377/** add one constraint to matrix */
378static
380 SCIP* scip, /**< current scip instance */
381 SCIP_MATRIX* matrix, /**< constraint matrix */
382 SCIP_VAR** vars, /**< variables of this constraint */
383 SCIP_Real* vals, /**< variable coefficients of this constraint */
384 int nvars, /**< number of variables */
385 SCIP_Real lhs, /**< left hand side */
386 SCIP_Real rhs, /**< right hand side */
387 int maxnnonzsmem, /**< maximal number of fillable elements */
388 SCIP_Bool* rowadded /**< flag indicating of row was added to matrix */
389 )
390{
391 SCIP_VAR** activevars;
392 SCIP_Real* activevals;
393 SCIP_Real activeconstant;
394 int nactivevars;
395 int v;
396
397 assert(scip != NULL);
398 assert(matrix != NULL);
399 assert(vars != NULL || nvars == 0);
400 assert(SCIPisLE(scip, lhs, rhs));
401 assert(rowadded != NULL);
402
403 *rowadded = FALSE;
404
405 /* constraint is redundant */
406 if( SCIPisInfinity(scip, -lhs) && SCIPisInfinity(scip, rhs) )
407 return SCIP_OKAY;
408
409 /* we do not add empty constraints to the matrix */
410 if( nvars == 0 )
411 return SCIP_OKAY;
412
413 activevars = NULL;
414 activevals = NULL;
415 nactivevars = nvars;
416 activeconstant = 0.0;
417
418 /* duplicate variable and value array */
419 SCIP_CALL( SCIPduplicateBufferArray(scip, &activevars, vars, nactivevars ) );
420 if( vals != NULL )
421 {
422 SCIP_CALL( SCIPduplicateBufferArray(scip, &activevals, vals, nactivevars ) );
423 }
424 else
425 {
426 SCIP_CALL( SCIPallocBufferArray(scip, &activevals, nactivevars) );
427
428 for( v = 0; v < nactivevars; v++ )
429 activevals[v] = 1.0;
430 }
431
432 /* retransform given variables to active variables */
433 SCIP_CALL( getActiveVariables(scip, &activevars, &activevals, &nactivevars, &activeconstant) );
434
435 /* adapt left and right hand side */
436 if( !SCIPisInfinity(scip, -lhs) )
437 lhs -= activeconstant;
438 if( !SCIPisInfinity(scip, rhs) )
439 rhs -= activeconstant;
440
441 /* add single row to matrix */
442 if( nactivevars > 0 )
443 {
444 SCIP_CALL( addRow(scip, matrix, activevars, activevals, nactivevars, lhs, rhs, maxnnonzsmem, rowadded) );
445 }
446
447 /* free buffer arrays */
448 SCIPfreeBufferArray(scip, &activevals);
449 SCIPfreeBufferArray(scip, &activevars);
450
451 return SCIP_OKAY;
452}
453
454/** add one constraint to matrix */
455static
457 SCIP* scip, /**< current scip instance */
458 SCIP_MATRIX* matrix, /**< constraint matrix */
459 SCIP_VAR** vars, /**< variables of this constraint */
460 SCIP_RATIONAL** vals, /**< variable coefficients of this constraint */
461 int nvars, /**< number of variables */
462 SCIP_RATIONAL* lhs, /**< left hand side */
463 SCIP_RATIONAL* rhs, /**< right hand side */
464 int maxnnonzsmem, /**< maximal number of fillable elements */
465 SCIP_Bool* rowadded /**< flag indicating of row was added to matrix */
466 )
467{
468 SCIP_VAR** activevars;
469 SCIP_RATIONAL** activevals;
470 SCIP_RATIONAL* activeconstant;
471 SCIP_RATIONAL* tmplhs; /* need these due to the constant */
472 SCIP_RATIONAL* tmprhs;
473 int nactivevars;
474 int v;
475
476 assert(scip != NULL);
477 assert(matrix != NULL);
478 assert(vars != NULL || nvars == 0);
479 assert(SCIPrationalIsLE(lhs, rhs));
480 assert(rowadded != NULL);
481
482 *rowadded = FALSE;
483
484 /* constraint is redundant */
486 return SCIP_OKAY;
487
488 /* we do not add empty constraints to the matrix */
489 if( nvars == 0 )
490 return SCIP_OKAY;
491
492 activevars = NULL;
493 activevals = NULL;
494 nactivevars = nvars;
498
499 /* duplicate variable and value array */
500 SCIP_CALL( SCIPduplicateBufferArray(scip, &activevars, vars, nactivevars ) );
501 if( vals != NULL )
502 {
503 SCIP_CALL( SCIPrationalCopyBufferArray(SCIPbuffer(scip), &activevals, vals, nactivevars ) );
504 }
505 else
506 {
507 SCIP_CALL( SCIPrationalCreateBufferArray(SCIPbuffer(scip), &activevals, nactivevars) );
508
509 for( v = 0; v < nactivevars; v++ )
510 SCIPrationalSetFraction(activevals[v], 1LL, 1LL);
511 }
512
513 /* retransform given variables to active variables */
514 SCIP_CALL( getActiveVariablesExact(scip, &activevars, activevals, &nactivevars, activeconstant) );
515
516 /* adapt left and right hand side */
517 if( !SCIPrationalIsNegInfinity(lhs) )
518 SCIPrationalDiff(tmplhs, lhs, activeconstant);
519 if( !SCIPrationalIsInfinity(rhs) )
520 SCIPrationalDiff(tmprhs, rhs, activeconstant);
521
522 /* add single row to matrix */
523 if( nactivevars > 0 )
524 {
525 SCIP_CALL( addRowExact(matrix, activevars, activevals, nactivevars, tmplhs, tmprhs, maxnnonzsmem, rowadded) );
526 }
527
528 /* free buffer arrays */
529 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &activevals, nvars);
530 SCIPfreeBufferArray(scip, &activevars);
531
534 SCIPrationalFreeBuffer(SCIPbuffer(scip), &activeconstant);
535
536 return SCIP_OKAY;
537}
538
539/** transform row major format into column major format */
540static
542 SCIP* scip, /**< current scip instance */
543 SCIP_MATRIX* matrix /**< constraint matrix */
544 )
545{
546 int colidx;
547 int i;
548 int* rowpnt;
549 int* rowend;
550 SCIP_Real* valpnt;
551 SCIP_RATIONAL* valpntrational;
552 int* fillidx;
553
554 assert(scip != NULL);
555 assert(matrix != NULL);
556 assert(matrix->colmatval != NULL);
557 assert(matrix->colmatind != NULL);
558 assert(matrix->colmatbeg != NULL);
559 assert(matrix->colmatcnt != NULL);
560 assert(matrix->rowmatval != NULL);
561 assert(matrix->rowmatind != NULL);
562 assert(matrix->rowmatbeg != NULL);
563 assert(matrix->rowmatcnt != NULL);
564
565 SCIP_CALL( SCIPallocBufferArray(scip, &fillidx, matrix->ncols) );
566 BMSclearMemoryArray(fillidx, matrix->ncols);
567 BMSclearMemoryArray(matrix->colmatcnt, matrix->ncols);
568
569 for( i = 0; i < matrix->nrows; i++ )
570 {
571 rowpnt = matrix->rowmatind + matrix->rowmatbeg[i];
572 rowend = rowpnt + matrix->rowmatcnt[i];
573 for( ; rowpnt < rowend; rowpnt++ )
574 {
575 colidx = *rowpnt;
576 (matrix->colmatcnt[colidx])++;
577 }
578 }
579
580 matrix->colmatbeg[0] = 0;
581 for( i = 0; i < matrix->ncols-1; i++ )
582 {
583 matrix->colmatbeg[i+1] = matrix->colmatbeg[i] + matrix->colmatcnt[i];
584 }
585
586 for( i = 0; i < matrix->nrows; i++ )
587 {
588 rowpnt = matrix->rowmatind + matrix->rowmatbeg[i];
589 rowend = rowpnt + matrix->rowmatcnt[i];
590 valpnt = matrix->rowmatval + matrix->rowmatbeg[i];
591 if( SCIPisExact(scip) )
592 valpntrational = matrix->matrixvalsexact->rowmatvalexact[matrix->rowmatbeg[i]];
593
594 for( ; rowpnt < rowend; rowpnt++, valpnt++ )
595 {
596 assert(*rowpnt < matrix->ncols);
597 colidx = *rowpnt;
598 matrix->colmatval[matrix->colmatbeg[colidx] + fillidx[colidx]] = *valpnt;
599 if( SCIPisExact(scip) )
600 SCIPrationalSetRational(matrix->matrixvalsexact->colmatvalexact[matrix->colmatbeg[colidx] + fillidx[colidx]], valpntrational); /*lint !e644*/
601 matrix->colmatind[matrix->colmatbeg[colidx] + fillidx[colidx]] = i;
602 fillidx[colidx]++;
603 }
604 }
605
606 SCIPfreeBufferArray(scip, &fillidx);
607
608 return SCIP_OKAY;
609}
610
611/** calculate min/max activity per row */
612static
614 SCIP* scip, /**< current scip instance */
615 SCIP_MATRIX* matrix /**< constraint matrix */
616 )
617{
618 SCIP_Real val;
619 int* rowpnt;
620 int* rowend;
621 SCIP_Real* valpnt;
622 int col;
623 int row;
624
625 assert(scip != NULL);
626 assert(matrix != NULL);
627
628 for( row = 0; row < matrix->nrows; row++ )
629 {
630 matrix->minactivity[row] = 0;
631 matrix->maxactivity[row] = 0;
632 matrix->minactivityneginf[row] = 0;
633 matrix->minactivityposinf[row] = 0;
634 matrix->maxactivityneginf[row] = 0;
635 matrix->maxactivityposinf[row] = 0;
636
637 rowpnt = matrix->rowmatind + matrix->rowmatbeg[row];
638 rowend = rowpnt + matrix->rowmatcnt[row];
639 valpnt = matrix->rowmatval + matrix->rowmatbeg[row];
640
641 for( ; rowpnt < rowend; rowpnt++, valpnt++ )
642 {
643 /* get column index */
644 col = *rowpnt;
645
646 /* get variable coefficient */
647 val = *valpnt;
648 assert(!SCIPisZero(scip, val) || SCIPisExact(scip));
649
650 assert(matrix->ncols > col);
651
652 assert(!SCIPisInfinity(scip, matrix->lb[col]));
653 assert(!SCIPisInfinity(scip, -matrix->ub[col]));
654
655 /* positive coefficient */
656 if( val > 0.0 )
657 {
658 if( SCIPisInfinity(scip, matrix->ub[col]) )
659 matrix->maxactivityposinf[row]++;
660 else
661 matrix->maxactivity[row] += val * matrix->ub[col];
662
663 if( SCIPisInfinity(scip, -matrix->lb[col]) )
664 matrix->minactivityneginf[row]++;
665 else
666 matrix->minactivity[row] += val * matrix->lb[col];
667 }
668 /* negative coefficient */
669 else
670 {
671 if( SCIPisInfinity(scip, -matrix->lb[col]) )
672 matrix->maxactivityneginf[row]++;
673 else
674 matrix->maxactivity[row] += val * matrix->lb[col];
675
676 if( SCIPisInfinity(scip, matrix->ub[col]) )
677 matrix->minactivityposinf[row]++;
678 else
679 matrix->minactivity[row] += val * matrix->ub[col];
680 }
681 }
682
683 /* consider infinite bound contributions for the activities */
684 if( matrix->maxactivityneginf[row] + matrix->maxactivityposinf[row] > 0 )
685 matrix->maxactivity[row] = SCIPinfinity(scip);
686
687 if( matrix->minactivityneginf[row] + matrix->minactivityposinf[row] > 0 )
688 matrix->minactivity[row] = -SCIPinfinity(scip);
689 }
690
691 return SCIP_OKAY;
692}
693
694/*
695 * public functions
696 */
697
698/** initialize matrix by copying all check constraints
699 *
700 * @note Completeness is checked by testing whether all check constraints are from a list of linear constraint handlers
701 * that can be represented.
702 */
704 SCIP* scip, /**< current scip instance */
705 SCIP_MATRIX** matrixptr, /**< pointer to constraint matrix object to be initialized */
706 SCIP_Bool onlyifcomplete, /**< should matrix creation be skipped if matrix will not be complete? */
707 SCIP_Bool* initialized, /**< was the initialization successful? */
708 SCIP_Bool* complete, /**< are all constraint represented within the matrix? */
709 SCIP_Bool* infeasible, /**< pointer to return whether problem was detected to be infeasible during matrix creation */
710 int* naddconss, /**< pointer to count number of added (linear) constraints during matrix creation */
711 int* ndelconss, /**< pointer to count number of deleted specialized linear constraints during matrix creation */
712 int* nchgcoefs, /**< pointer to count number of changed coefficients during matrix creation */
713 int* nchgbds, /**< pointer to count number of changed bounds during matrix creation */
714 int* nfixedvars /**< pointer to count number of fixed variables during matrix creation */
715 )
716{
717 SCIP_MATRIX* matrix;
718 SCIP_CONSHDLR** conshdlrs;
719 const char* conshdlrname;
720 SCIP_Bool stopped;
721 SCIP_VAR** vars;
722 SCIP_VAR* var;
723 SCIP_CONS** conshdlrconss;
724 SCIP_CONS* cons;
725 int nconshdlrconss;
726 int nconshdlrs;
727 int nconss;
728 int nconssall;
729 int nnonzstmp;
730 int nvars;
731 int c;
732 int i;
733 int v;
734 int cnt;
735
736 nnonzstmp = 0;
737
738 assert(scip != NULL);
739 assert(matrixptr != NULL);
740 assert(initialized != NULL);
741 assert(complete != NULL);
742
743 *initialized = FALSE;
744 *complete = FALSE;
745 *infeasible = FALSE;
746
747 /* return if no variables or constraints are present */
748 if( SCIPgetNVars(scip) == 0 || SCIPgetNConss(scip) == 0 )
749 return SCIP_OKAY;
750
751 /* return if pricers are present and the matrix should only be built when complete */
752 if( onlyifcomplete && SCIPgetNActivePricers(scip) != 0 )
753 return SCIP_OKAY;
754
755 /* loop over all constraint handlers and collect the number of checked constraints */
756 nconshdlrs = SCIPgetNConshdlrs(scip);
757 conshdlrs = SCIPgetConshdlrs(scip);
758 nconss = 0;
759 nconssall = 0;
760
761 for( i = 0; i < nconshdlrs; ++i )
762 {
763 nconshdlrconss = SCIPconshdlrGetNCheckConss(conshdlrs[i]);
764
765 if( nconshdlrconss > 0 )
766 {
767 conshdlrname = SCIPconshdlrGetName(conshdlrs[i]);
768
769 if( strcmp(conshdlrname, "linear") == 0 || strcmp(conshdlrname, "knapsack") == 0
770 || strcmp(conshdlrname, "setppc") == 0 || strcmp(conshdlrname, "logicor") == 0
771 || strcmp(conshdlrname, "varbound") == 0 || (strcmp(conshdlrname, "exactlinear") == 0) )
772 {
773 /* increment number of supported constraints */
774 nconss += nconshdlrconss;
775 }
776/* disabled because some of the presolvers can currently only handle 1-1 row-cons relationships */
777#ifdef SCIP_DISABLED_CODE
778 else if( strcmp(conshdlrname, "linking") == 0 )
779 {
780 /* the linear representation of linking constraints involves two linear constraints */
781 nconss += 2* nconshdlrconss;
782 }
783#endif
784 /* increment number of supported and unsupported constraints */
785 nconssall += nconshdlrconss;
786 }
787 }
788
789 /* print warning if we have unsupported constraint types; we only abort the matrix creation process if requested,
790 * because it makes sometimes sense to work on an incomplete matrix as long as the number of interesting variable
791 * uplocks or downlocks of the matrix and scip are the same
792 */
793 if( nconss < nconssall )
794 {
795 SCIPdebugMsg(scip, "Warning: milp matrix not complete!\n");
796 if( onlyifcomplete )
797 return SCIP_OKAY;
798 }
799 else
800 {
801 /* all constraints represented within the matrix */
802 *complete = TRUE;
803 }
804
805 /* do nothing if we have no checked constraints */
806 if( nconss == 0 )
807 return SCIP_OKAY;
808
809 stopped = FALSE;
810
811 /* first, clean up aggregations and fixings in varbound costraints, since this can lead
812 * to boundchanges and the varbound constraint can get downgraded to a linear constraint
813 */
814 SCIP_CALL( SCIPcleanupConssVarbound(scip, TRUE, infeasible, naddconss, ndelconss, nchgbds ) );
815 if( *infeasible )
816 return SCIP_OKAY;
817
818 /* next, clean up aggregations and fixings in setppc costraints, since this can lead
819 * to fixings and the setppc constraint can get downgraded to a linear constraint
820 */
821 SCIP_CALL( SCIPcleanupConssSetppc(scip, TRUE, infeasible, naddconss, ndelconss, nchgcoefs, nfixedvars ) );
822 if( *infeasible )
823 return SCIP_OKAY;
824
825 /* next, clean up aggregations and fixings in logicor costraints, since this cannot lead
826 * to further fixings but the logicor constraint can also get downgraded to a linear constraint
827 */
828 SCIP_CALL( SCIPcleanupConssLogicor(scip, TRUE, naddconss, ndelconss, nchgcoefs) );
829
830 /* finally, clean up aggregations and fixings in knapsack and linear constraints since now no new linaer constraints
831 * can come up due to downgrading and the remaining cleanup methods cannot fix any more variables
832 */
833
834 SCIP_CALL( SCIPcleanupConssKnapsack(scip, TRUE, infeasible, ndelconss) );
835 if( *infeasible )
836 return SCIP_OKAY;
837
838 SCIP_CALL( SCIPcleanupConssLinear(scip, TRUE, infeasible, ndelconss) );
839 if( *infeasible )
840 return SCIP_OKAY;
841
842 vars = SCIPgetVars(scip);
843 nvars = SCIPgetNVars(scip);
844
845 /* approximate number of nonzeros by taking for each variable the number of up- and downlocks;
846 * this counts nonzeros in equalities twice, but can be at most two times as high as the exact number
847 */
848 for( i = nvars - 1; i >= 0; --i )
849 {
850 nnonzstmp += SCIPvarGetNLocksDownType(vars[i], SCIP_LOCKTYPE_MODEL);
851 nnonzstmp += SCIPvarGetNLocksUpType(vars[i], SCIP_LOCKTYPE_MODEL);
852 }
853
854 /* do nothing if we have no entries */
855 if( nnonzstmp == 0 )
856 return SCIP_OKAY;
857
858 /* build the matrix structure */
859 SCIP_CALL( SCIPallocBuffer(scip, matrixptr) );
860 matrix = *matrixptr;
861
862 /* copy vars array and set number of variables */
863 SCIP_CALL( SCIPduplicateBufferArray(scip, &matrix->vars, vars, nvars) );
864 matrix->ncols = nvars;
865
866 matrix->nrows = 0;
867 matrix->nnonzs = 0;
868
869 /* allocate memory */
870 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->colmatval, nnonzstmp) );
871 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->colmatind, nnonzstmp) );
872 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->colmatbeg, matrix->ncols) );
873 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->colmatcnt, matrix->ncols) );
874 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->lb, matrix->ncols) );
875 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->ub, matrix->ncols) );
876 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->nuplocks, matrix->ncols) );
877 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->ndownlocks, matrix->ncols) );
878 if( SCIPisExact(scip) )
879 {
887 matrix->matrixvalsexact->buffersize = nnonzstmp;
888 matrix->matrixvalsexact->buffersizenconss = nconss;
889 }
890
891 BMSclearMemoryArray(matrix->nuplocks, matrix->ncols);
892 BMSclearMemoryArray(matrix->ndownlocks, matrix->ncols);
893
894 /* init bounds */
895 for( v = 0; v < matrix->ncols; v++ )
896 {
897 var = matrix->vars[v];
898 assert(var != NULL);
899
900 matrix->lb[v] = SCIPvarGetLbGlobal(var);
901 matrix->ub[v] = SCIPvarGetUbGlobal(var);
902
903 if( SCIPisExact(scip) )
904 {
907 }
908 }
909
910 /* allocate memory */
911 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->rowmatval, nnonzstmp) );
912 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->rowmatind, nnonzstmp) );
913 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->rowmatbeg, nconss) );
914 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->rowmatcnt, nconss) );
915 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->lhs, nconss) );
916 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->rhs, nconss) );
917 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->cons, nconss) );
919 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->minactivity, nconss) );
920 SCIP_CALL( SCIPallocBufferArray(scip, &matrix->maxactivity, nconss) );
925
926 cnt = 0;
927
928 /* loop a second time over constraints handlers and add supported constraints to the matrix */
929 for( i = 0; i < nconshdlrs; ++i )
930 {
931 SCIP_Bool rowadded;
932
933 if( SCIPisStopped(scip) || (onlyifcomplete && !(*complete)) )
934 {
935 stopped = TRUE;
936 break;
937 }
938
939 conshdlrname = SCIPconshdlrGetName(conshdlrs[i]);
940 conshdlrconss = SCIPconshdlrGetCheckConss(conshdlrs[i]);
941 nconshdlrconss = SCIPconshdlrGetNCheckConss(conshdlrs[i]);
942
943 if( strcmp(conshdlrname, "linear") == 0 )
944 {
945 for( c = 0; c < nconshdlrconss && (c % 1000 != 0 || !SCIPisStopped(scip)); ++c )
946 {
947 cons = conshdlrconss[c];
948 assert(SCIPconsIsTransformed(cons));
949
950 /* do not include constraints that can be altered due to column generation */
951 if( SCIPconsIsModifiable(cons) )
952 {
953 *complete = FALSE;
954
955 if( onlyifcomplete )
956 break;
957
958 continue;
959 }
960
963 SCIPgetLhsLinear(scip, cons), SCIPgetRhsLinear(scip, cons), nnonzstmp, &rowadded) );
964
965 if(rowadded)
966 {
967 assert(cnt < nconss);
968 matrix->cons[cnt] = cons;
969 cnt++;
970 }
971 }
972 }
973 else if( strcmp(conshdlrname, "exactlinear") == 0 )
974 {
975 for( c = 0; c < nconshdlrconss && (c % 1000 != 0 || !SCIPisStopped(scip)); ++c )
976 {
977 cons = conshdlrconss[c];
978 assert(SCIPconsIsTransformed(cons));
979
980 /* do not include constraints that can be altered due to column generation */
981 if( SCIPconsIsModifiable(cons) )
982 {
983 *complete = FALSE;
984
985 if( onlyifcomplete )
986 break;
987
988 continue;
989 }
990
993 SCIPgetLhsExactLinear(scip, cons), SCIPgetRhsExactLinear(scip, cons), nnonzstmp, &rowadded) );
994
995 if(rowadded)
996 {
997 assert(cnt < nconss);
998 matrix->cons[cnt] = cons;
999 cnt++;
1000 }
1001 }
1002 }
1003 else if( strcmp(conshdlrname, "knapsack") == 0 )
1004 {
1005 if( nconshdlrconss > 0 )
1006 {
1007 SCIP_Real* consvals;
1008 int valssize;
1009
1010 valssize = 100;
1011 SCIP_CALL( SCIPallocBufferArray(scip, &consvals, valssize) );
1012
1013 for( c = 0; c < nconshdlrconss && (c % 1000 != 0 || !SCIPisStopped(scip)); ++c )
1014 {
1015 SCIP_Longint* weights;
1016
1017 cons = conshdlrconss[c];
1018 assert(SCIPconsIsTransformed(cons));
1019
1020 /* do not include constraints that can be altered due to column generation */
1021 if( SCIPconsIsModifiable(cons) )
1022 {
1023 *complete = FALSE;
1024
1025 if( onlyifcomplete )
1026 break;
1027
1028 continue;
1029 }
1030
1031 weights = SCIPgetWeightsKnapsack(scip, cons);
1032 nvars = SCIPgetNVarsKnapsack(scip, cons);
1033
1034 if( nvars > valssize )
1035 {
1036 valssize = (int) (1.5 * nvars);
1037 SCIP_CALL( SCIPreallocBufferArray(scip, &consvals, valssize) );
1038 }
1039
1040 for( v = 0; v < nvars; v++ )
1041 consvals[v] = (SCIP_Real)weights[v];
1042
1043 SCIP_CALL( addConstraint(scip, matrix, SCIPgetVarsKnapsack(scip, cons), consvals,
1045 (SCIP_Real)SCIPgetCapacityKnapsack(scip, cons), nnonzstmp, &rowadded) );
1046
1047 if(rowadded)
1048 {
1049 assert(cnt < nconss);
1050 matrix->cons[cnt] = cons;
1051 cnt++;
1052 }
1053 }
1054
1055 SCIPfreeBufferArray(scip, &consvals);
1056 }
1057 }
1058 else if( strcmp(conshdlrname, "setppc") == 0 )
1059 {
1060 for( c = 0; c < nconshdlrconss && (c % 1000 != 0 || !SCIPisStopped(scip)); ++c )
1061 {
1062 SCIP_Real lhs;
1063 SCIP_Real rhs;
1064
1065 cons = conshdlrconss[c];
1066 assert(SCIPconsIsTransformed(cons));
1067
1068 /* do not include constraints that can be altered due to column generation */
1069 if( SCIPconsIsModifiable(cons) )
1070 {
1071 *complete = FALSE;
1072
1073 if( onlyifcomplete )
1074 break;
1075
1076 continue;
1077 }
1078
1079 switch( SCIPgetTypeSetppc(scip, cons) )
1080 {
1082 lhs = 1.0;
1083 rhs = 1.0;
1084 break;
1086 lhs = -SCIPinfinity(scip);
1087 rhs = 1.0;
1088 break;
1090 lhs = 1.0;
1091 rhs = SCIPinfinity(scip);
1092 break;
1093 default:
1094 return SCIP_ERROR;
1095 }
1096
1098 SCIPgetNVarsSetppc(scip, cons), lhs, rhs, nnonzstmp, &rowadded) );
1099
1100 if(rowadded)
1101 {
1102 assert(cnt < nconss);
1103 matrix->cons[cnt] = cons;
1104 cnt++;
1105 }
1106 }
1107 }
1108 else if( strcmp(conshdlrname, "logicor") == 0 )
1109 {
1110 for( c = 0; c < nconshdlrconss && (c % 1000 != 0 || !SCIPisStopped(scip)); ++c )
1111 {
1112 cons = conshdlrconss[c];
1113 assert(SCIPconsIsTransformed(cons));
1114
1115 /* do not include constraints that can be altered due to column generation */
1116 if( SCIPconsIsModifiable(cons) )
1117 {
1118 *complete = FALSE;
1119
1120 if( onlyifcomplete )
1121 break;
1122
1123 continue;
1124 }
1125
1127 NULL, SCIPgetNVarsLogicor(scip, cons), 1.0, SCIPinfinity(scip), nnonzstmp, &rowadded) );
1128
1129 if(rowadded)
1130 {
1131 assert(cnt < nconss);
1132 matrix->cons[cnt] = cons;
1133 cnt++;
1134 }
1135 }
1136 }
1137 else if( strcmp(conshdlrname, "varbound") == 0 )
1138 {
1139 if( nconshdlrconss > 0 )
1140 {
1141 SCIP_VAR** consvars;
1142 SCIP_Real* consvals;
1143
1144 SCIP_CALL( SCIPallocBufferArray(scip, &consvars, 2) );
1145 SCIP_CALL( SCIPallocBufferArray(scip, &consvals, 2) );
1146 consvals[0] = 1.0;
1147
1148 for( c = 0; c < nconshdlrconss && (c % 1000 != 0 || !SCIPisStopped(scip)); ++c )
1149 {
1150 cons = conshdlrconss[c];
1151 assert(SCIPconsIsTransformed(cons));
1152
1153 /* do not include constraints that can be altered due to column generation */
1154 if( SCIPconsIsModifiable(cons) )
1155 {
1156 *complete = FALSE;
1157
1158 if( onlyifcomplete )
1159 break;
1160
1161 continue;
1162 }
1163
1164 consvars[0] = SCIPgetVarVarbound(scip, cons);
1165 consvars[1] = SCIPgetVbdvarVarbound(scip, cons);
1166
1167 consvals[1] = SCIPgetVbdcoefVarbound(scip, cons);
1168
1169 SCIP_CALL( addConstraint(scip, matrix, consvars, consvals, 2, SCIPgetLhsVarbound(scip, cons),
1170 SCIPgetRhsVarbound(scip, cons), nnonzstmp, &rowadded) );
1171
1172 if(rowadded)
1173 {
1174 assert(cnt < nconss);
1175 matrix->cons[cnt] = cons;
1176 cnt++;
1177 }
1178 }
1179
1180 SCIPfreeBufferArray(scip, &consvals);
1181 SCIPfreeBufferArray(scip, &consvars);
1182 }
1183 }
1184/* the code below is correct. However, it needs to be disabled
1185 * because some of the presolvers can currently only handle 1-1 row-cons relationships,
1186 * while the linking constraint handler requires a representation as 2 linear constraints.
1187 */
1188#ifdef SCIP_DISABLED_CODE
1189 else if( strcmp(conshdlrname, "linking") == 0 )
1190 {
1191 if( nconshdlrconss > 0 )
1192 {
1193 SCIP_VAR** consvars;
1194 SCIP_VAR** curconsvars;
1195 SCIP_Real* consvals;
1196 int* curconsvals;
1197 int valssize;
1198 int nconsvars;
1199 int j;
1200
1201 valssize = 100;
1202 SCIP_CALL( SCIPallocBufferArray(scip, &consvars, valssize) );
1203 SCIP_CALL( SCIPallocBufferArray(scip, &consvals, valssize) );
1204
1205 for( c = 0; c < nconshdlrconss && (c % 1000 != 0 || !SCIPisStopped(scip)); ++c )
1206 {
1207 cons = conshdlrconss[c];
1208 assert(SCIPconsIsTransformed(cons));
1209
1210 /* do not include constraints that can be altered due to column generation */
1211 if( SCIPconsIsModifiable(cons) )
1212 {
1213 *complete = FALSE;
1214
1215 if( onlyifcomplete )
1216 break;
1217
1218 continue;
1219 }
1220
1221 /* get constraint variables and their amount */
1222 SCIP_CALL( SCIPgetBinvarsLinking(scip, cons, &curconsvars, &nconsvars) );
1223 curconsvals = SCIPgetValsLinking(scip, cons);
1224
1225 /* SCIPgetBinVarsLinking returns the number of binary variables, but we also need the integer variable */
1226 nconsvars++;
1227
1228 if( nconsvars > valssize )
1229 {
1230 valssize = (int) (1.5 * nconsvars);
1231 SCIP_CALL( SCIPreallocBufferArray(scip, &consvars, valssize) );
1232 SCIP_CALL( SCIPreallocBufferArray(scip, &consvals, valssize) );
1233 }
1234
1235 /* copy vars and vals for binary variables */
1236 for( j = 0; j < nconsvars - 1; j++ )
1237 {
1238 consvars[j] = curconsvars[j];
1239 consvals[j] = (SCIP_Real) curconsvals[j];
1240 }
1241
1242 /* set final entry of vars and vals to the linking variable and its coefficient, respectively */
1243 consvars[nconsvars - 1] = SCIPgetIntvarLinking(scip, cons);
1244 consvals[nconsvars - 1] = -1;
1245
1246 SCIP_CALL( addConstraint(scip, matrix, consvars, consvals, nconsvars, 0.0, 0.0, nnonzstmp, &rowadded) );
1247 SCIP_CALL( addConstraint(scip, matrix, consvars, NULL, nconsvars - 1, 1.0, 1.0, nnonzstmp, &rowadded) );
1248
1249 if(rowadded)
1250 {
1251 assert(cnt < nconss);
1252 matrix->cons[cnt] = cons;
1253 matrix->cons[cnt + 1] = cons;
1254 cnt += 2;
1255 }
1256 }
1257
1258 SCIPfreeBufferArray(scip, &consvals);
1259 SCIPfreeBufferArray(scip, &consvars);
1260 }
1261 }
1262#endif
1263 }
1264 assert(matrix->nrows == cnt);
1265 assert(matrix->nrows <= nconss);
1266 assert(matrix->nnonzs <= nnonzstmp);
1267
1268 if( *complete )
1269 {
1270 SCIP_Bool lockmismatch = FALSE;
1271
1272 for( i = 0; i < matrix->ncols; ++i )
1273 {
1274 if( SCIPmatrixUplockConflict(matrix, i) || SCIPmatrixDownlockConflict(matrix, i) )
1275 {
1276 lockmismatch = TRUE;
1277 break;
1278 }
1279 }
1280
1281 if( lockmismatch )
1282 {
1283 *complete = FALSE;
1284 if( onlyifcomplete )
1285 stopped = TRUE;
1286 }
1287 }
1288
1289 if( !stopped )
1290 {
1291 /* calculate row activity bounds */
1292 SCIP_CALL( calcActivityBounds(scip, matrix) );
1293
1294 /* transform row major format into column major format */
1296
1297 *initialized = TRUE;
1298 }
1299 else
1300 {
1307
1309 SCIPfreeBufferArray(scip, &matrix->cons);
1310
1311 SCIPfreeBufferArray(scip, &matrix->rhs);
1312 SCIPfreeBufferArray(scip, &matrix->lhs);
1317
1320
1321 if( SCIPisExact(scip) )
1322 {
1330 }
1331
1332 SCIPfreeBufferArray(scip, &matrix->ub);
1333 SCIPfreeBufferArray(scip, &matrix->lb);
1339
1340 SCIPfreeBuffer(scip, matrixptr);
1341 }
1342
1343 return SCIP_OKAY;
1344}
1345
1346
1347/** frees the constraint matrix */
1349 SCIP* scip, /**< current SCIP instance */
1350 SCIP_MATRIX** matrix /**< constraint matrix object */
1351 )
1352{
1353 assert(scip != NULL);
1354 assert(matrix != NULL);
1355
1356 if( (*matrix) != NULL )
1357 {
1358 assert((*matrix)->colmatval != NULL);
1359 assert((*matrix)->colmatind != NULL);
1360 assert((*matrix)->colmatbeg != NULL);
1361 assert((*matrix)->colmatcnt != NULL);
1362 assert((*matrix)->lb != NULL);
1363 assert((*matrix)->ub != NULL);
1364 assert((*matrix)->nuplocks != NULL);
1365 assert((*matrix)->ndownlocks != NULL);
1366
1367 assert((*matrix)->rowmatval != NULL);
1368 assert((*matrix)->rowmatind != NULL);
1369 assert((*matrix)->rowmatbeg != NULL);
1370 assert((*matrix)->rowmatcnt != NULL);
1371 assert((*matrix)->lhs != NULL);
1372 assert((*matrix)->rhs != NULL);
1373
1374 SCIPfreeBufferArray(scip, &((*matrix)->maxactivityposinf));
1375 SCIPfreeBufferArray(scip, &((*matrix)->maxactivityneginf));
1376 SCIPfreeBufferArray(scip, &((*matrix)->minactivityposinf));
1377 SCIPfreeBufferArray(scip, &((*matrix)->minactivityneginf));
1378 SCIPfreeBufferArray(scip, &((*matrix)->maxactivity));
1379 SCIPfreeBufferArray(scip, &((*matrix)->minactivity));
1380
1381 SCIPfreeMemoryArray(scip, &((*matrix)->isrhsinfinite));
1382 SCIPfreeBufferArray(scip, &((*matrix)->cons));
1383
1384 SCIPfreeBufferArray(scip, &((*matrix)->rhs));
1385 SCIPfreeBufferArray(scip, &((*matrix)->lhs));
1386 SCIPfreeBufferArray(scip, &((*matrix)->rowmatcnt));
1387 SCIPfreeBufferArray(scip, &((*matrix)->rowmatbeg));
1388 SCIPfreeBufferArray(scip, &((*matrix)->rowmatind));
1389 SCIPfreeBufferArray(scip, &((*matrix)->rowmatval));
1390
1391 SCIPfreeBufferArray(scip, &((*matrix)->ndownlocks));
1392 SCIPfreeBufferArray(scip, &((*matrix)->nuplocks));
1393
1394 if( SCIPisExact(scip) )
1395 {
1396 assert((*matrix)->matrixvalsexact != NULL);
1397 SCIPfreeBufferArray(scip, &(*matrix)->matrixvalsexact->ubexact);
1398 SCIPfreeBufferArray(scip, &(*matrix)->matrixvalsexact->lbexact);
1399 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &(*matrix)->matrixvalsexact->rhsexact, (*matrix)->matrixvalsexact->buffersizenconss);
1400 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &(*matrix)->matrixvalsexact->lhsexact, (*matrix)->matrixvalsexact->buffersizenconss);
1401
1402 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &(*matrix)->matrixvalsexact->rowmatvalexact, (*matrix)->matrixvalsexact->buffersize);
1403 SCIPrationalFreeBufferArray(SCIPbuffer(scip), &(*matrix)->matrixvalsexact->colmatvalexact, (*matrix)->matrixvalsexact->buffersize);
1404 SCIPfreeBuffer(scip, &(*matrix)->matrixvalsexact);
1405 }
1406
1407 SCIPfreeBufferArray(scip, &((*matrix)->ub));
1408 SCIPfreeBufferArray(scip, &((*matrix)->lb));
1409 SCIPfreeBufferArray(scip, &((*matrix)->colmatcnt));
1410 SCIPfreeBufferArray(scip, &((*matrix)->colmatbeg));
1411 SCIPfreeBufferArray(scip, &((*matrix)->colmatind));
1412 SCIPfreeBufferArray(scip, &((*matrix)->colmatval));
1413
1414 (*matrix)->nrows = 0;
1415 (*matrix)->ncols = 0;
1416 (*matrix)->nnonzs = 0;
1417
1418 SCIPfreeBufferArrayNull(scip, &((*matrix)->vars));
1419
1420 SCIPfreeBuffer(scip, matrix);
1421 }
1422}
1423
1424/** print one row of the matrix */
1426 SCIP* scip, /**< current SCIP instance */
1427 SCIP_MATRIX* matrix, /**< constraint matrix object */
1428 int row /**< row index */
1429 )
1430{
1431 int* rowpnt;
1432 int* rowend;
1433 int col;
1434 SCIP_Real val;
1435 SCIP_Real* valpnt;
1436
1438
1439 rowpnt = matrix->rowmatind + matrix->rowmatbeg[row];
1440 rowend = rowpnt + matrix->rowmatcnt[row];
1441 valpnt = matrix->rowmatval + matrix->rowmatbeg[row];
1442
1443 printf("### %s: %.15g <=", SCIPconsGetName(matrix->cons[row]), matrix->lhs[row]);
1444 for(; (rowpnt < rowend); rowpnt++, valpnt++)
1445 {
1446 col = *rowpnt;
1447 val = *valpnt;
1448 if( val < 0 )
1449 printf(" %.15g %s [%.15g,%.15g]", val, SCIPvarGetName(matrix->vars[col]),
1450 SCIPvarGetLbGlobal(matrix->vars[col]), SCIPvarGetUbGlobal(matrix->vars[col]));
1451 else
1452 printf(" +%.15g %s [%.15g,%.15g]", val, SCIPvarGetName(matrix->vars[col]),
1453 SCIPvarGetLbGlobal(matrix->vars[col]), SCIPvarGetUbGlobal(matrix->vars[col]));
1454 }
1455 printf(" <= %.15g ###\n", matrix->rhs[row]);
1456}
1457
1458/** removes the bounds of a column and updates the activities accordingly */
1460 SCIP* scip, /**< current scip instance */
1461 SCIP_MATRIX* matrix, /**< constraint matrix */
1462 int col /**< column variable to remove bounds from */
1463 )
1464{
1465 int colmatend = matrix->colmatbeg[col] + matrix->colmatcnt[col];
1466 int i;
1467
1468 for( i = matrix->colmatbeg[col]; i != colmatend; ++i )
1469 {
1470 int row = matrix->colmatind[i];
1471 SCIP_Real val = matrix->colmatval[i];
1472
1473 /* set lower bound to -infinity if necessary */
1474 if( !SCIPisInfinity(scip, -matrix->lb[col]) )
1475 {
1476 if( val > 0.0 )
1477 matrix->minactivityneginf[row]++;
1478 else
1479 matrix->maxactivityneginf[row]++;
1480 }
1481
1482 /* set upper bound to infinity if necessary */
1483 if( !SCIPisInfinity(scip, matrix->ub[col]) )
1484 {
1485 if( val > 0.0 )
1486 matrix->maxactivityposinf[row]++;
1487 else
1488 matrix->minactivityposinf[row]++;
1489 }
1490
1491 assert(matrix->maxactivityneginf[row] + matrix->maxactivityposinf[row] > 0);
1492 assert(matrix->minactivityneginf[row] + matrix->minactivityposinf[row] > 0);
1493
1494 /* mark the activities of the rows to be infinite */
1495 matrix->maxactivity[row] = SCIPinfinity(scip);
1496 matrix->minactivity[row] = -SCIPinfinity(scip);
1497 }
1498
1499 matrix->lb[col] = -SCIPinfinity(scip);
1500 matrix->ub[col] = SCIPinfinity(scip);
1501}
1502
1503/** detect parallel rows of matrix. rhs/lhs are ignored. */
1505 SCIP* scip, /**< SCIP instance */
1506 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
1507 SCIP_Real* scale, /**< scale factors of rows */
1508 int* pclass /**< parallel row classes */
1509 )
1510{
1511 SCIP_Real* valpnt;
1512 SCIP_Real* values;
1513 int* classsizes;
1514 int* pcset;
1515 int* colpnt;
1516 int* colend;
1517 int* rowindices;
1518 int* pcs;
1519 SCIP_Real startval;
1520 SCIP_Real aij;
1521 int startpc;
1522 int startk;
1523 int startt;
1524 int pcsetfill;
1525 int rowidx;
1526 int k;
1527 int t;
1528 int m;
1529 int i;
1530 int c;
1531 int newpclass;
1532 int pc;
1533
1534 assert(scip != NULL);
1535 assert(matrix != NULL);
1536 assert(pclass != NULL);
1537
1538 SCIP_CALL( SCIPallocBufferArray(scip, &classsizes, matrix->nrows) );
1539 SCIP_CALL( SCIPallocBufferArray(scip, &pcset, matrix->nrows) );
1540 SCIP_CALL( SCIPallocBufferArray(scip, &values, matrix->nrows) );
1541 SCIP_CALL( SCIPallocBufferArray(scip, &rowindices, matrix->nrows) );
1542 SCIP_CALL( SCIPallocBufferArray(scip, &pcs, matrix->nrows) );
1543
1544 /* init */
1545 BMSclearMemoryArray(scale, matrix->nrows);
1546 BMSclearMemoryArray(pclass, matrix->nrows);
1547 BMSclearMemoryArray(classsizes, matrix->nrows);
1548 classsizes[0] = matrix->nrows;
1549 pcsetfill = 0;
1550 for( t = 1; t < matrix->nrows; ++t )
1551 pcset[pcsetfill++] = t;
1552
1553 /* loop over all columns */
1554 for( c = 0; c < matrix->ncols; ++c )
1555 {
1556 if( matrix->colmatcnt[c] == 0 )
1557 continue;
1558
1559 colpnt = matrix->colmatind + matrix->colmatbeg[c];
1560 colend = colpnt + matrix->colmatcnt[c];
1561 valpnt = matrix->colmatval + matrix->colmatbeg[c];
1562
1563 i = 0;
1564 for( ; (colpnt < colend); colpnt++, valpnt++ )
1565 {
1566 aij = *valpnt;
1567 rowidx = *colpnt;
1568
1569 if( scale[rowidx] == 0.0 )
1570 scale[rowidx] = aij;
1571 assert(scale[rowidx] != 0.0);
1572
1573 rowindices[i] = rowidx;
1574 values[i] = aij / scale[rowidx];
1575 pc = pclass[rowidx];
1576 assert(pc < matrix->nrows);
1577
1578 /* update class sizes and pclass set */
1579 assert(classsizes[pc] > 0);
1580 classsizes[pc]--;
1581 if( classsizes[pc] == 0 )
1582 {
1583 assert(pcsetfill < matrix->nrows);
1584 pcset[pcsetfill++] = pc;
1585 }
1586 pcs[i] = pc;
1587
1588 i++;
1589 }
1590
1591 /* sort on the pclass values */
1592 if( i > 1 )
1593 {
1594 SCIPsortIntIntReal(pcs, rowindices, values, i);
1595 }
1596
1597 k = 0;
1598 while( TRUE ) /*lint !e716*/
1599 {
1600 assert(k < i);
1601 startpc = pcs[k];
1602 startk = k;
1603
1604 /* find pclass-sets */
1605 while( k < i && pcs[k] == startpc )
1606 k++;
1607
1608 /* sort on the A values which have equal pclass values */
1609 if( k - startk > 1 )
1610 SCIPsortRealInt(&(values[startk]), &(rowindices[startk]), k - startk);
1611
1612 t = 0;
1613 while( TRUE ) /*lint !e716*/
1614 {
1615 assert(startk + t < i);
1616 startval = values[startk + t];
1617 startt = t;
1618
1619 /* find A-sets */
1620 while( t < k - startk && SCIPisEQ(scip, startval, values[startk + t]) )
1621 t++;
1622
1623 /* get new pclass */
1624 newpclass = pcset[0];
1625 assert(pcsetfill > 0);
1626 pcset[0] = pcset[--pcsetfill];
1627
1628 /* renumbering */
1629 for( m = startk + startt; m < startk + t; m++ )
1630 {
1631 assert(m < i);
1632 assert(rowindices[m] < matrix->nrows);
1633 assert(newpclass < matrix->nrows);
1634
1635 pclass[rowindices[m]] = newpclass;
1636 classsizes[newpclass]++;
1637 }
1638
1639 if( t == k - startk )
1640 break;
1641 }
1642
1643 if( k == matrix->colmatcnt[c] )
1644 break;
1645 }
1646 }
1647
1649 SCIPfreeBufferArray(scip, &rowindices);
1650 SCIPfreeBufferArray(scip, &values);
1651 SCIPfreeBufferArray(scip, &pcset);
1652 SCIPfreeBufferArray(scip, &classsizes);
1653
1654 return SCIP_OKAY;
1655}
1656
1657/** detect parallel rows of matrix.
1658 * obj coefficients are ignored.
1659 */
1661 SCIP* scip, /**< SCIP instance */
1662 SCIP_MATRIX* matrix, /**< matrix containing the constraints */
1663 SCIP_Real* scale, /**< scale factors of cols */
1664 int* pclass, /**< parallel column classes */
1665 SCIP_Bool* varineq /**< indicating if variable is within an equation */
1666 )
1667{
1668 SCIP_Real* valpnt;
1669 SCIP_Real* values;
1670 int* classsizes;
1671 int* pcset;
1672 int* rowpnt;
1673 int* rowend;
1674 int* colindices;
1675 int* pcs;
1676 SCIP_Real startval;
1677 SCIP_Real aij;
1678 int startpc;
1679 int startk;
1680 int startt;
1681 int pcsetfill;
1682 int colidx;
1683 int k;
1684 int t;
1685 int m;
1686 int i;
1687 int r;
1688 int newpclass;
1689 int pc;
1690
1691 assert(scip != NULL);
1692 assert(matrix != NULL);
1693 assert(pclass != NULL);
1694 assert(varineq != NULL);
1695
1696 SCIP_CALL( SCIPallocBufferArray(scip, &classsizes, matrix->ncols) );
1697 SCIP_CALL( SCIPallocBufferArray(scip, &pcset, matrix->ncols) );
1698 SCIP_CALL( SCIPallocBufferArray(scip, &values, matrix->ncols) );
1699 SCIP_CALL( SCIPallocBufferArray(scip, &colindices, matrix->ncols) );
1700 SCIP_CALL( SCIPallocBufferArray(scip, &pcs, matrix->ncols) );
1701
1702 /* init */
1703 BMSclearMemoryArray(scale, matrix->ncols);
1704 BMSclearMemoryArray(pclass, matrix->ncols);
1705 BMSclearMemoryArray(classsizes, matrix->ncols);
1706 classsizes[0] = matrix->ncols;
1707 pcsetfill = 0;
1708 for( t = 1; t < matrix->ncols; ++t )
1709 pcset[pcsetfill++] = t;
1710
1711 /* loop over all rows */
1712 for( r = 0; r < matrix->nrows; ++r )
1713 {
1714 /* we consider only equations or ranged rows */
1715 if( !matrix->isrhsinfinite[r] )
1716 {
1717 rowpnt = matrix->rowmatind + matrix->rowmatbeg[r];
1718 rowend = rowpnt + matrix->rowmatcnt[r];
1719 valpnt = matrix->rowmatval + matrix->rowmatbeg[r];
1720
1721 i = 0;
1722 for( ; (rowpnt < rowend); rowpnt++, valpnt++ )
1723 {
1724 aij = *valpnt;
1725 colidx = *rowpnt;
1726
1727 /* remember variable was part of an equation or ranged row */
1728 varineq[colidx] = TRUE;
1729
1730 if( scale[colidx] == 0.0 )
1731 scale[colidx] = aij;
1732 assert(scale[colidx] != 0.0);
1733
1734 colindices[i] = colidx;
1735 values[i] = aij / scale[colidx];
1736 pc = pclass[colidx];
1737 assert(pc < matrix->ncols);
1738
1739 /* update class sizes and pclass set */
1740 assert(classsizes[pc] > 0);
1741 classsizes[pc]--;
1742 if( classsizes[pc] == 0 )
1743 {
1744 assert(pcsetfill < matrix->ncols);
1745 pcset[pcsetfill++] = pc;
1746 }
1747 pcs[i] = pc;
1748
1749 i++;
1750 }
1751
1752 /* sort on the pclass values */
1753 if( i > 1 )
1754 {
1755 SCIPsortIntIntReal(pcs, colindices, values, i);
1756 }
1757
1758 k = 0;
1759 while( TRUE ) /*lint !e716*/
1760 {
1761 assert(k < i);
1762 startpc = pcs[k];
1763 startk = k;
1764
1765 /* find pclass-sets */
1766 while( k < i && pcs[k] == startpc )
1767 k++;
1768
1769 /* sort on the A values which have equal pclass values */
1770 if( k - startk > 1 )
1771 SCIPsortRealInt(&(values[startk]), &(colindices[startk]), k - startk);
1772
1773 t = 0;
1774 while( TRUE ) /*lint !e716*/
1775 {
1776 assert(startk + t < i);
1777 startval = values[startk + t];
1778 startt = t;
1779
1780 /* find A-sets */
1781 while( t < k - startk && SCIPisEQ(scip, startval, values[startk + t]) )
1782 t++;
1783
1784 /* get new pclass */
1785 newpclass = pcset[0];
1786 assert(pcsetfill > 0);
1787 pcset[0] = pcset[--pcsetfill];
1788
1789 /* renumbering */
1790 for( m = startk + startt; m < startk + t; m++ )
1791 {
1792 assert(m < i);
1793 assert(colindices[m] < matrix->ncols);
1794 assert(newpclass < matrix->ncols);
1795
1796 pclass[colindices[m]] = newpclass;
1797 classsizes[newpclass]++;
1798 }
1799
1800 if( t == k - startk )
1801 break;
1802 }
1803
1804 if( k == matrix->rowmatcnt[r] )
1805 break;
1806 }
1807 }
1808 }
1809
1811 SCIPfreeBufferArray(scip, &colindices);
1812 SCIPfreeBufferArray(scip, &values);
1813 SCIPfreeBufferArray(scip, &pcset);
1814 SCIPfreeBufferArray(scip, &classsizes);
1815
1816 return SCIP_OKAY;
1817}
1818
1819
1820/*
1821 * access functions implemented as defines
1822 */
1823
1824/* In debug mode, the following methods are implemented as function calls to ensure
1825 * type validity.
1826 * In optimized mode, the methods are implemented as defines to improve performance.
1827 * However, we want to have them in the library anyways, so we have to undef the defines.
1828 */
1829
1830#undef SCIPmatrixGetColValPtr
1831#undef SCIPmatrixGetColIdxPtr
1832#undef SCIPmatrixGetColNNonzs
1833#undef SCIPmatrixGetNColumns
1834#undef SCIPmatrixGetColUb
1835#undef SCIPmatrixGetColLb
1836#undef SCIPmatrixGetColNUplocks
1837#undef SCIPmatrixGetColNDownlocks
1838#undef SCIPmatrixGetVar
1839#undef SCIPmatrixGetColName
1840#undef SCIPmatrixGetRowValPtr
1841#undef SCIPmatrixGetRowValPtrExact
1842#undef SCIPmatrixGetRowIdxPtr
1843#undef SCIPmatrixGetRowNNonzs
1844#undef SCIPmatrixGetRowName
1845#undef SCIPmatrixGetNRows
1846#undef SCIPmatrixGetRowLhs
1847#undef SCIPmatrixGetRowLhsExact
1848#undef SCIPmatrixGetRowRhs
1849#undef SCIPmatrixGetRowRhsExact
1850#undef SCIPmatrixIsRowRhsInfinity
1851#undef SCIPmatrixGetNNonzs
1852#undef SCIPmatrixGetRowMinActivity
1853#undef SCIPmatrixGetRowMaxActivity
1854#undef SCIPmatrixGetRowNMinActNegInf
1855#undef SCIPmatrixGetRowNMinActPosInf
1856#undef SCIPmatrixGetRowNMaxActNegInf
1857#undef SCIPmatrixGetRowNMaxActPosInf
1858#undef SCIPmatrixGetCons
1859
1860/** get column based start pointer of values */
1862 SCIP_MATRIX* matrix, /**< matrix instance */
1863 int col /**< column index */
1864 )
1865{
1866 assert(matrix != NULL);
1867 assert(0 <= col && col < matrix->ncols);
1868
1869 return matrix->colmatval + matrix->colmatbeg[col];
1870}
1871
1872/** get column based start pointer of row indices */
1874 SCIP_MATRIX* matrix, /**< matrix instance */
1875 int col /**< column index */
1876 )
1877{
1878 assert(matrix != NULL);
1879 assert(0 <= col && col < matrix->ncols);
1880
1881 return matrix->colmatind + matrix->colmatbeg[col];
1882}
1883
1884/** get the number of non-zero entries of this column */
1886 SCIP_MATRIX* matrix, /**< matrix instance */
1887 int col /**< column index */
1888 )
1889{
1890 assert(matrix != NULL);
1891 assert(0 <= col && col < matrix->ncols);
1892
1893 return matrix->colmatcnt[col];
1894}
1895
1896/** get number of columns of the matrix */
1898 SCIP_MATRIX* matrix /**< matrix instance */
1899 )
1900{
1901 assert(matrix != NULL);
1902
1903 return matrix->ncols;
1904}
1905
1906/** get upper bound of column */
1908 SCIP_MATRIX* matrix, /**< matrix instance */
1909 int col /**< column index */
1910 )
1911{
1912 assert(matrix != NULL);
1913
1914 return matrix->ub[col];
1915}
1916
1917/** get lower bound of column */
1919 SCIP_MATRIX* matrix, /**< matrix instance */
1920 int col /**< column index */
1921 )
1922{
1923 assert(matrix != NULL);
1924
1925 return matrix->lb[col];
1926}
1927
1928/** get number of uplocks of column */
1930 SCIP_MATRIX* matrix, /**< matrix instance */
1931 int col /**< column index */
1932 )
1933{
1934 assert(matrix != NULL);
1935 assert(0 <= col && col < matrix->ncols);
1936
1937 return matrix->nuplocks[col];
1938}
1939
1940/** get number of downlocks of column */
1942 SCIP_MATRIX* matrix, /**< matrix instance */
1943 int col /**< column index */
1944 )
1945{
1946 assert(matrix != NULL);
1947 assert(0 <= col && col < matrix->ncols);
1948
1949 return matrix->ndownlocks[col];
1950}
1951
1952/** get variable pointer of column */
1954 SCIP_MATRIX* matrix, /**< matrix instance */
1955 int col /**< column index */
1956 )
1957{
1958 assert(matrix != NULL);
1959 assert(0 <= col && col < matrix->ncols);
1960
1961 return matrix->vars[col];
1962}
1963
1964/** get name of column/variable */
1966 SCIP_MATRIX* matrix, /**< matrix instance */
1967 int col /**< column index */
1968 )
1969{
1970 assert(matrix != NULL);
1971 assert(0 <= col && col < matrix->ncols);
1972
1973 return SCIPvarGetName(matrix->vars[col]);
1974}
1975
1976/** get row based start pointer of values */
1978 SCIP_MATRIX* matrix, /**< matrix instance */
1979 int row /**< row index */
1980 )
1981{
1982 assert(matrix != NULL);
1983 assert(0 <= row && row < matrix->nrows);
1984
1985 return matrix->rowmatval + matrix->rowmatbeg[row];
1986}
1987
1988/** get row based start pointer of values */
1990 SCIP_MATRIX* matrix, /**< matrix instance */
1991 int row /**< row index */
1992 )
1993{
1994 assert(matrix != NULL);
1995 assert(0 <= row && row < matrix->nrows);
1996
1997 return matrix->matrixvalsexact->rowmatvalexact + matrix->rowmatbeg[row];
1998}
1999
2000/** get row based start pointer of column indices */
2002 SCIP_MATRIX* matrix, /**< matrix instance */
2003 int row /**< row index */
2004 )
2005{
2006 assert(matrix != NULL);
2007 assert(0 <= row && row < matrix->nrows);
2008
2009 return matrix->rowmatind + matrix->rowmatbeg[row];
2010}
2011
2012/** get number of non-zeros of this row */
2014 SCIP_MATRIX* matrix, /**< matrix instance */
2015 int row /**< row index */
2016 )
2017{
2018 assert(matrix != NULL);
2019 assert(0 <= row && row < matrix->nrows);
2020
2021 return matrix->rowmatcnt[row];
2022}
2023
2024/** get name of row */
2026 SCIP_MATRIX* matrix, /**< matrix instance */
2027 int row /**< row index */
2028 )
2029{
2030 assert(matrix != NULL);
2031 assert(0 <= row && row < matrix->nrows);
2032
2033 return SCIPconsGetName(matrix->cons[row]);
2034}
2035
2036/** get number of rows of the matrix */
2038 SCIP_MATRIX* matrix /**< matrix instance */
2039 )
2040{
2041 assert(matrix != NULL);
2042
2043 return matrix->nrows;
2044}
2045
2046/** get left-hand-side of row */
2048 SCIP_MATRIX* matrix, /**< matrix instance */
2049 int row /**< row index */
2050 )
2051{
2052 assert(matrix != NULL);
2053 assert(0 <= row && row < matrix->nrows);
2054
2055 return matrix->lhs[row];
2056}
2057
2058/** get right-hand-side of row */
2060 SCIP_MATRIX* matrix, /**< matrix instance */
2061 int row /**< row index */
2062 )
2063{
2064 assert(matrix != NULL);
2065 assert(0 <= row && row < matrix->nrows);
2066
2067 return matrix->rhs[row];
2068}
2069
2070/** get left-hand-side of row */
2072 SCIP_MATRIX* matrix, /**< matrix instance */
2073 int row /**< row index */
2074 )
2075{
2076 assert(matrix != NULL);
2077 assert(0 <= row && row < matrix->nrows);
2078
2079 return matrix->matrixvalsexact->lhsexact[row];
2080}
2081
2082/** get right-hand-side of row */
2084 SCIP_MATRIX* matrix, /**< matrix instance */
2085 int row /**< row index */
2086 )
2087{
2088 assert(matrix != NULL);
2089 assert(0 <= row && row < matrix->nrows);
2090
2091 return matrix->matrixvalsexact->rhsexact[row];
2092}
2093
2094/** flag indicating if right-hand-side of row is infinity */
2096 SCIP_MATRIX* matrix, /**< matrix instance */
2097 int row /**< row index */
2098 )
2099{
2100 assert(matrix != NULL);
2101 assert(0 <= row && row < matrix->nrows);
2102
2103 return matrix->isrhsinfinite[row];
2104}
2105
2106/** get number of non-zeros of matrix */
2108 SCIP_MATRIX* matrix /**< matrix instance */
2109 )
2110{
2111 assert(matrix != NULL);
2112
2113 return matrix->nnonzs;
2114}
2115
2116/** get minimal activity of row */
2118 SCIP_MATRIX* matrix, /**< matrix instance */
2119 int row /**< row index */
2120 )
2121{
2122 assert(matrix != NULL);
2123 assert(0 <= row && row < matrix->nrows);
2124
2125 return matrix->minactivity[row];
2126}
2127
2128/** get maximal activity of row */
2130 SCIP_MATRIX* matrix, /**< matrix instance */
2131 int row /**< row index */
2132 )
2133{
2134 assert(matrix != NULL);
2135 assert(0 <= row && row < matrix->nrows);
2136
2137 return matrix->maxactivity[row];
2138}
2139
2140/** get number of negative infinities present within minimal activity */
2142 SCIP_MATRIX* matrix, /**< matrix instance */
2143 int row /**< row index */
2144 )
2145{
2146 assert(matrix != NULL);
2147 assert(0 <= row && row < matrix->nrows);
2148
2149 return matrix->minactivityneginf[row];
2150}
2151
2152/** get number of positive infinities present within minimal activity */
2154 SCIP_MATRIX* matrix, /**< matrix instance */
2155 int row /**< row index */
2156 )
2157{
2158 assert(matrix != NULL);
2159 assert(0 <= row && row < matrix->nrows);
2160
2161 return matrix->minactivityposinf[row];
2162}
2163
2164/** get number of negative infinities present within maximal activity */
2166 SCIP_MATRIX* matrix, /**< matrix instance */
2167 int row /**< row index */
2168 )
2169{
2170 assert(matrix != NULL);
2171 assert(0 <= row && row < matrix->nrows);
2172
2173 return matrix->maxactivityneginf[row];
2174}
2175
2176/** get number of positive infinities present within maximal activity */
2178 SCIP_MATRIX* matrix, /**< matrix instance */
2179 int row /**< row index */
2180 )
2181{
2182 assert(matrix != NULL);
2183 assert(0 <= row && row < matrix->nrows);
2184
2185 return matrix->maxactivityposinf[row];
2186}
2187
2188/** get constraint pointer for constraint representing row */
2190 SCIP_MATRIX* matrix, /**< matrix instance */
2191 int row /**< row index */
2192 )
2193{
2194 assert(matrix != NULL);
2195 assert(0 <= row && row < matrix->nrows);
2196
2197 return matrix->cons[row];
2198}
2199
2200/** get if conflicting uplocks of a specific variable present */
2202 SCIP_MATRIX* matrix, /**< matrix instance */
2203 int col /**< column index */
2204 )
2205{
2206 assert(matrix != NULL);
2207 assert(0 <= col && col < matrix->ncols);
2208
2209 return (SCIPvarGetNLocksUpType(matrix->vars[col], SCIP_LOCKTYPE_MODEL) != matrix->nuplocks[col]);
2210}
2211
2212/** get if conflicting downlocks of a specific variable present */
2214 SCIP_MATRIX* matrix, /**< matrix instance */
2215 int col /**< column index */
2216 )
2217{
2218 assert(matrix != NULL);
2219 assert(0 <= col && col < matrix->ncols);
2220
2221 return (SCIPvarGetNLocksDownType(matrix->vars[col], SCIP_LOCKTYPE_MODEL) != matrix->ndownlocks[col]);
2222}
SCIP_Real * r
Definition: circlepacking.c:59
Constraint handler for linear constraints in their most general form, .
Constraint handler for knapsack constraints of the form , x binary and .
Constraint handler for linear constraints in their most general form, .
Constraint handler for logicor constraints (equivalent to set covering, but algorithms are suited fo...
Constraint handler for the set partitioning / packing / covering constraints .
Constraint handler for variable bound constraints .
#define NULL
Definition: def.h:248
#define SCIP_Longint
Definition: def.h:141
#define SCIP_UNUSED(x)
Definition: def.h:409
#define SCIP_Bool
Definition: def.h:91
#define SCIP_Real
Definition: def.h:156
#define TRUE
Definition: def.h:93
#define FALSE
Definition: def.h:94
#define SCIP_CALL(x)
Definition: def.h:355
SCIP_RETCODE SCIPgetBinvarsLinking(SCIP *scip, SCIP_CONS *cons, SCIP_VAR ***binvars, int *nbinvars)
int SCIPgetNVarsKnapsack(SCIP *scip, SCIP_CONS *cons)
SCIP_Real SCIPgetVbdcoefVarbound(SCIP *scip, SCIP_CONS *cons)
SCIP_RATIONAL * SCIPgetLhsExactLinear(SCIP *scip, SCIP_CONS *cons)
int SCIPgetNVarsLogicor(SCIP *scip, SCIP_CONS *cons)
SCIP_RETCODE SCIPcleanupConssSetppc(SCIP *scip, SCIP_Bool onlychecked, SCIP_Bool *infeasible, int *naddconss, int *ndelconss, int *nchgcoefs, int *nfixedvars)
Definition: cons_setppc.c:9839
SCIP_Real SCIPgetRhsLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_VAR ** SCIPgetVarsLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_RATIONAL * SCIPgetRhsExactLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_RETCODE SCIPcleanupConssLinear(SCIP *scip, SCIP_Bool onlychecked, SCIP_Bool *infeasible, int *ndelconss)
SCIP_VAR ** SCIPgetVarsExactLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_Real SCIPgetLhsLinear(SCIP *scip, SCIP_CONS *cons)
int SCIPgetNVarsLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_RETCODE SCIPcleanupConssKnapsack(SCIP *scip, SCIP_Bool onlychecked, SCIP_Bool *infeasible, int *ndelconss)
int SCIPgetNVarsExactLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_Real * SCIPgetValsLinear(SCIP *scip, SCIP_CONS *cons)
SCIP_VAR * SCIPgetVbdvarVarbound(SCIP *scip, SCIP_CONS *cons)
int SCIPgetNVarsSetppc(SCIP *scip, SCIP_CONS *cons)
Definition: cons_setppc.c:9596
SCIP_RETCODE SCIPcleanupConssVarbound(SCIP *scip, SCIP_Bool onlychecked, SCIP_Bool *infeasible, int *naddconss, int *ndelconss, int *nchgbds)
SCIP_RETCODE SCIPcleanupConssLogicor(SCIP *scip, SCIP_Bool onlychecked, int *naddconss, int *ndelconss, int *nchgcoefs)
SCIP_VAR ** SCIPgetVarsSetppc(SCIP *scip, SCIP_CONS *cons)
Definition: cons_setppc.c:9619
SCIP_VAR * SCIPgetVarVarbound(SCIP *scip, SCIP_CONS *cons)
SCIP_Longint * SCIPgetWeightsKnapsack(SCIP *scip, SCIP_CONS *cons)
SCIP_Longint SCIPgetCapacityKnapsack(SCIP *scip, SCIP_CONS *cons)
SCIP_Real SCIPgetLhsVarbound(SCIP *scip, SCIP_CONS *cons)
SCIP_SETPPCTYPE SCIPgetTypeSetppc(SCIP *scip, SCIP_CONS *cons)
Definition: cons_setppc.c:9642
SCIP_VAR ** SCIPgetVarsLogicor(SCIP *scip, SCIP_CONS *cons)
SCIP_Real SCIPgetRhsVarbound(SCIP *scip, SCIP_CONS *cons)
SCIP_VAR ** SCIPgetVarsKnapsack(SCIP *scip, SCIP_CONS *cons)
SCIP_Real * SCIPgetValsLinking(SCIP *scip, SCIP_CONS *cons)
SCIP_RATIONAL ** SCIPgetValsExactLinear(SCIP *scip, SCIP_CONS *cons)
@ SCIP_SETPPCTYPE_PARTITIONING
Definition: cons_setppc.h:87
@ SCIP_SETPPCTYPE_COVERING
Definition: cons_setppc.h:89
@ SCIP_SETPPCTYPE_PACKING
Definition: cons_setppc.h:88
SCIP_Bool SCIPisStopped(SCIP *scip)
Definition: scip_general.c:759
int SCIPgetNVars(SCIP *scip)
Definition: scip_prob.c:2246
int SCIPgetNConss(SCIP *scip)
Definition: scip_prob.c:3620
SCIP_VAR ** SCIPgetVars(SCIP *scip)
Definition: scip_prob.c:2201
#define SCIPdebugMsg
Definition: scip_message.h:78
int SCIPconshdlrGetNCheckConss(SCIP_CONSHDLR *conshdlr)
Definition: cons.c:4798
SCIP_CONS ** SCIPconshdlrGetCheckConss(SCIP_CONSHDLR *conshdlr)
Definition: cons.c:4755
int SCIPgetNConshdlrs(SCIP *scip)
Definition: scip_cons.c:964
const char * SCIPconshdlrGetName(SCIP_CONSHDLR *conshdlr)
Definition: cons.c:4316
SCIP_CONSHDLR ** SCIPgetConshdlrs(SCIP *scip)
Definition: scip_cons.c:953
SCIP_Bool SCIPconsIsTransformed(SCIP_CONS *cons)
Definition: cons.c:8698
const char * SCIPconsGetName(SCIP_CONS *cons)
Definition: cons.c:8389
SCIP_Bool SCIPconsIsModifiable(SCIP_CONS *cons)
Definition: cons.c:8638
SCIP_Bool SCIPisExact(SCIP *scip)
Definition: scip_exact.c:193
#define SCIPfreeBuffer(scip, ptr)
Definition: scip_mem.h:134
BMS_BUFMEM * SCIPbuffer(SCIP *scip)
Definition: scip_mem.c:72
#define SCIPallocBufferArray(scip, ptr, num)
Definition: scip_mem.h:124
#define SCIPreallocBufferArray(scip, ptr, num)
Definition: scip_mem.h:128
#define SCIPallocClearMemoryArray(scip, ptr, num)
Definition: scip_mem.h:66
#define SCIPfreeBufferArray(scip, ptr)
Definition: scip_mem.h:136
#define SCIPfreeMemoryArray(scip, ptr)
Definition: scip_mem.h:80
#define SCIPduplicateBufferArray(scip, ptr, source, num)
Definition: scip_mem.h:132
#define SCIPallocBuffer(scip, ptr)
Definition: scip_mem.h:122
#define SCIPfreeBufferArrayNull(scip, ptr)
Definition: scip_mem.h:137
int SCIPgetNActivePricers(SCIP *scip)
Definition: scip_pricer.c:348
void SCIPrationalSetInfinity(SCIP_RATIONAL *res)
Definition: rational.cpp:618
SCIP_Real SCIPrationalGetReal(SCIP_RATIONAL *rational)
Definition: rational.cpp:2085
void SCIPrationalFreeBuffer(BMS_BUFMEM *bufmem, SCIP_RATIONAL **rational)
Definition: rational.cpp:473
void SCIPrationalDiff(SCIP_RATIONAL *res, SCIP_RATIONAL *op1, SCIP_RATIONAL *op2)
Definition: rational.cpp:983
SCIP_RETCODE SCIPrationalCopyBufferArray(BMS_BUFMEM *mem, SCIP_RATIONAL ***result, SCIP_RATIONAL **src, int len)
Definition: rational.cpp:267
SCIP_Bool SCIPrationalIsPositive(SCIP_RATIONAL *rational)
Definition: rational.cpp:1640
SCIP_RETCODE SCIPrationalCreateBuffer(BMS_BUFMEM *bufmem, SCIP_RATIONAL **rational)
Definition: rational.cpp:123
SCIP_Bool SCIPrationalIsZero(SCIP_RATIONAL *rational)
Definition: rational.cpp:1624
void SCIPrationalSetRational(SCIP_RATIONAL *res, SCIP_RATIONAL *src)
Definition: rational.cpp:569
void SCIPrationalSetFraction(SCIP_RATIONAL *res, SCIP_Longint nom, SCIP_Longint denom)
Definition: rational.cpp:582
void SCIPrationalNegate(SCIP_RATIONAL *res, SCIP_RATIONAL *op)
Definition: rational.cpp:1297
SCIP_Bool SCIPrationalIsNegative(SCIP_RATIONAL *rational)
Definition: rational.cpp:1650
SCIP_Bool SCIPrationalIsInfinity(SCIP_RATIONAL *rational)
Definition: rational.cpp:1660
SCIP_RETCODE SCIPrationalCreateBufferArray(BMS_BUFMEM *mem, SCIP_RATIONAL ***rational, int size)
Definition: rational.cpp:214
SCIP_Bool SCIPrationalIsNegInfinity(SCIP_RATIONAL *rational)
Definition: rational.cpp:1670
SCIP_RETCODE SCIPrationalReallocBufferArray(BMS_BUFMEM *mem, SCIP_RATIONAL ***result, int oldlen, int newlen)
Definition: rational.cpp:314
void SCIPrationalMultReal(SCIP_RATIONAL *res, SCIP_RATIONAL *op1, SCIP_Real op2)
Definition: rational.cpp:1097
void SCIPrationalFreeBufferArray(BMS_BUFMEM *mem, SCIP_RATIONAL ***ratbufarray, int size)
Definition: rational.cpp:518
SCIP_RETCODE SCIPrationalCopyBuffer(BMS_BUFMEM *bufmem, SCIP_RATIONAL **result, SCIP_RATIONAL *src)
Definition: rational.cpp:165
SCIP_Bool SCIPrationalIsLE(SCIP_RATIONAL *rat1, SCIP_RATIONAL *rat2)
Definition: rational.cpp:1521
SCIP_Real SCIPinfinity(SCIP *scip)
SCIP_Bool SCIPisLE(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPisInfinity(SCIP *scip, SCIP_Real val)
SCIP_Bool SCIPisEQ(SCIP *scip, SCIP_Real val1, SCIP_Real val2)
SCIP_Bool SCIPisZero(SCIP *scip, SCIP_Real val)
int SCIPvarGetNLocksUpType(SCIP_VAR *var, SCIP_LOCKTYPE locktype)
Definition: var.c:4386
SCIP_Real SCIPvarGetUbGlobal(SCIP_VAR *var)
Definition: var.c:24142
int SCIPvarGetProbindex(SCIP_VAR *var)
Definition: var.c:23662
const char * SCIPvarGetName(SCIP_VAR *var)
Definition: var.c:23267
SCIP_RETCODE SCIPgetProbvarLinearSumExact(SCIP *scip, SCIP_VAR **vars, SCIP_RATIONAL **scalars, int *nvars, int varssize, SCIP_RATIONAL *constant, int *requiredsize, SCIP_Bool mergemultiples)
Definition: scip_var.c:2443
SCIP_RETCODE SCIPgetProbvarLinearSum(SCIP *scip, SCIP_VAR **vars, SCIP_Real *scalars, int *nvars, int varssize, SCIP_Real *constant, int *requiredsize)
Definition: scip_var.c:2378
SCIP_RATIONAL * SCIPvarGetLbGlobalExact(SCIP_VAR *var)
Definition: var.c:24130
SCIP_Real SCIPvarGetLbGlobal(SCIP_VAR *var)
Definition: var.c:24120
int SCIPvarGetNLocksDownType(SCIP_VAR *var, SCIP_LOCKTYPE locktype)
Definition: var.c:4328
SCIP_RATIONAL * SCIPvarGetUbGlobalExact(SCIP_VAR *var)
Definition: var.c:24152
void SCIPsortRealInt(SCIP_Real *realarray, int *intarray, int len)
void SCIPsortIntIntReal(int *intarray1, int *intarray2, SCIP_Real *realarray, int len)
static const SCIP_Real scalars[]
Definition: lp.c:5959
SCIP_Bool SCIPmatrixUplockConflict(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:2201
int SCIPmatrixGetRowNMinActNegInf(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2141
int * SCIPmatrixGetColIdxPtr(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1873
int SCIPmatrixGetNNonzs(SCIP_MATRIX *matrix)
Definition: matrix.c:2107
SCIP_RATIONAL * SCIPmatrixGetRowLhsExact(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2071
static SCIP_RETCODE setColumnMajorFormat(SCIP *scip, SCIP_MATRIX *matrix)
Definition: matrix.c:541
int SCIPmatrixGetRowNNonzs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2013
const char * SCIPmatrixGetRowName(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2025
int SCIPmatrixGetColNDownlocks(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1941
int SCIPmatrixGetColNNonzs(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1885
SCIP_Bool SCIPmatrixIsRowRhsInfinity(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2095
int SCIPmatrixGetColNUplocks(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1929
SCIP_Real SCIPmatrixGetRowMaxActivity(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2129
static SCIP_RETCODE calcActivityBounds(SCIP *scip, SCIP_MATRIX *matrix)
Definition: matrix.c:613
SCIP_Real SCIPmatrixGetColLb(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1918
static SCIP_RETCODE addRowExact(SCIP_MATRIX *matrix, SCIP_VAR **vars, SCIP_RATIONAL **vals, int nvars, SCIP_RATIONAL *lhs, SCIP_RATIONAL *rhs, int maxnnonzsmem, SCIP_Bool *rowadded)
Definition: matrix.c:260
SCIP_Real SCIPmatrixGetRowLhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2047
const char * SCIPmatrixGetColName(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1965
SCIP_Real * SCIPmatrixGetRowValPtr(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1977
SCIP_Bool SCIPmatrixDownlockConflict(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:2213
SCIP_RATIONAL ** SCIPmatrixGetRowValPtrExact(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1989
static SCIP_RETCODE addRow(SCIP *scip, SCIP_MATRIX *matrix, SCIP_VAR **vars, SCIP_Real *vals, int nvars, SCIP_Real lhs, SCIP_Real rhs, int maxnnonzsmem, SCIP_Bool *rowadded)
Definition: matrix.c:142
SCIP_Real SCIPmatrixGetRowRhs(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2059
SCIP_Real * SCIPmatrixGetColValPtr(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1861
int SCIPmatrixGetRowNMinActPosInf(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2153
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_Real SCIPmatrixGetRowMinActivity(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2117
SCIP_CONS * SCIPmatrixGetCons(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2189
SCIP_RETCODE SCIPmatrixGetParallelRows(SCIP *scip, SCIP_MATRIX *matrix, SCIP_Real *scale, int *pclass)
Definition: matrix.c:1504
static SCIP_RETCODE getActiveVariables(SCIP *scip, SCIP_VAR ***vars, SCIP_Real **scalars, int *nvars, SCIP_Real *constant)
Definition: matrix.c:70
static SCIP_RETCODE getActiveVariablesExact(SCIP *scip, SCIP_VAR ***vars, SCIP_RATIONAL **scalars, int *nvars, SCIP_RATIONAL *constant)
Definition: matrix.c:105
void SCIPmatrixFree(SCIP *scip, SCIP_MATRIX **matrix)
Definition: matrix.c:1348
int SCIPmatrixGetRowNMaxActPosInf(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2177
int SCIPmatrixGetRowNMaxActNegInf(SCIP_MATRIX *matrix, int row)
Definition: matrix.c:2165
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
static SCIP_RETCODE addConstraintExact(SCIP *scip, SCIP_MATRIX *matrix, SCIP_VAR **vars, SCIP_RATIONAL **vals, int nvars, SCIP_RATIONAL *lhs, SCIP_RATIONAL *rhs, int maxnnonzsmem, SCIP_Bool *rowadded)
Definition: matrix.c:456
static SCIP_RETCODE addConstraint(SCIP *scip, SCIP_MATRIX *matrix, SCIP_VAR **vars, SCIP_Real *vals, int nvars, SCIP_Real lhs, SCIP_Real rhs, int maxnnonzsmem, SCIP_Bool *rowadded)
Definition: matrix.c:379
void SCIPmatrixPrintRow(SCIP *scip, SCIP_MATRIX *matrix, int row)
Definition: matrix.c:1425
SCIP_RETCODE SCIPmatrixGetParallelCols(SCIP *scip, SCIP_MATRIX *matrix, SCIP_Real *scale, int *pclass, SCIP_Bool *varineq)
Definition: matrix.c:1660
int SCIPmatrixGetNRows(SCIP_MATRIX *matrix)
Definition: matrix.c:2037
void SCIPmatrixRemoveColumnBounds(SCIP *scip, SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1459
SCIP_Real SCIPmatrixGetColUb(SCIP_MATRIX *matrix, int col)
Definition: matrix.c:1907
memory allocation routines
#define BMSclearMemoryArray(ptr, num)
Definition: memory.h:130
public methods for managing constraints
public methods for matrix
public methods for message output
methods for sorting joint arrays of various types
public methods for problem variables
wrapper for rational number arithmetic
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 variable pricer plugins
public methods for global and local (sub)problems
public methods for SCIP variables
SCIP_RATIONAL ** rhsexact
Definition: struct_matrix.h:51
SCIP_RATIONAL ** lbexact
Definition: struct_matrix.h:48
SCIP_RATIONAL ** ubexact
Definition: struct_matrix.h:49
SCIP_RATIONAL ** lhsexact
Definition: struct_matrix.h:50
SCIP_RATIONAL ** colmatvalexact
Definition: struct_matrix.h:52
SCIP_RATIONAL ** rowmatvalexact
Definition: struct_matrix.h:53
SCIP_Real * lb
Definition: struct_matrix.h:68
SCIP_Bool * isrhsinfinite
Definition: struct_matrix.h:86
int * maxactivityneginf
Definition: struct_matrix.h:92
SCIP_VAR ** vars
Definition: struct_matrix.h:73
int * rowmatcnt
Definition: struct_matrix.h:78
int * minactivityposinf
Definition: struct_matrix.h:91
int * colmatbeg
Definition: struct_matrix.h:65
SCIP_MATRIXVALSEXACT * matrixvalsexact
Definition: struct_matrix.h:94
int * rowmatbeg
Definition: struct_matrix.h:77
SCIP_Real * lhs
Definition: struct_matrix.h:81
SCIP_Real * minactivity
Definition: struct_matrix.h:88
int * colmatind
Definition: struct_matrix.h:64
int * rowmatind
Definition: struct_matrix.h:76
SCIP_Real * colmatval
Definition: struct_matrix.h:63
int * nuplocks
Definition: struct_matrix.h:70
SCIP_CONS ** cons
Definition: struct_matrix.h:84
SCIP_Real * maxactivity
Definition: struct_matrix.h:89
int * maxactivityposinf
Definition: struct_matrix.h:93
SCIP_Real * rowmatval
Definition: struct_matrix.h:75
int * colmatcnt
Definition: struct_matrix.h:66
SCIP_Real * rhs
Definition: struct_matrix.h:82
int * minactivityneginf
Definition: struct_matrix.h:90
SCIP_Real * ub
Definition: struct_matrix.h:69
int * ndownlocks
Definition: struct_matrix.h:71
data structure for MIP matrix
@ SCIP_OKAY
Definition: type_retcode.h:42
@ SCIP_ERROR
Definition: type_retcode.h:43
enum SCIP_Retcode SCIP_RETCODE
Definition: type_retcode.h:63
@ SCIP_LOCKTYPE_MODEL
Definition: type_var.h:141