Scippy

    SCIP

    Solving Constraint Integer Programs

    lpi_glop.cpp
    Go to the documentation of this file.
    1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
    2/* */
    3/* This file is part of the program and library */
    4/* SCIP --- Solving Constraint Integer Programs */
    5/* */
    6/* Copyright (c) 2002-2025 Zuse Institute Berlin (ZIB) */
    7/* */
    8/* Licensed under the Apache License, Version 2.0 (the "License"); */
    9/* you may not use this file except in compliance with the License. */
    10/* You may obtain a copy of the License at */
    11/* */
    12/* http://www.apache.org/licenses/LICENSE-2.0 */
    13/* */
    14/* Unless required by applicable law or agreed to in writing, software */
    15/* distributed under the License is distributed on an "AS IS" BASIS, */
    16/* WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. */
    17/* See the License for the specific language governing permissions and */
    18/* limitations under the License. */
    19/* */
    20/* You should have received a copy of the Apache-2.0 license */
    21/* along with SCIP; see the file LICENSE. If not visit scipopt.org. */
    22/* */
    23/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
    24
    25/**@file lpi_glop.cpp
    26 * @ingroup LPIS
    27 * @brief LP interface for Glop
    28 * @author Frederic Didier
    29 * @author Marc Pfetsch
    30 */
    31
    32/*--+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
    33
    34/* turn off some warnings from includes */
    35#pragma GCC diagnostic ignored "-Wsign-compare"
    36#pragma GCC diagnostic ignored "-Wpedantic"
    37#pragma GCC diagnostic ignored "-Wignored-qualifiers"
    38#pragma GCC diagnostic ignored "-Wshadow"
    39#pragma GCC diagnostic ignored "-Wnon-virtual-dtor"
    40#pragma GCC diagnostic ignored "-Wctor-dtor-privacy"
    41#pragma GCC diagnostic ignored "-Woverflow"
    42
    43#include "ortools/base/version.h"
    44#include "ortools/glop/lp_solver.h"
    45#include "ortools/glop/revised_simplex.h"
    46#include "ortools/lp_data/lp_print_utils.h"
    47#include "ortools/lp_data/lp_data_utils.h"
    48#include "ortools/lp_data/proto_utils.h"
    49#include "ortools/util/file_util.h"
    50#include "ortools/util/stats.h"
    51#include "ortools/util/time_limit.h"
    52
    53#include "ortools/base/logging.h"
    54#include "ortools/base/vlog_is_on.h"
    55
    56#include "lpi/lpi.h"
    57#include "scip/pub_message.h"
    58
    59#include <assert.h>
    60
    61/* turn warnings on again */
    62#pragma GCC diagnostic warning "-Wsign-compare"
    63#pragma GCC diagnostic warning "-Wpedantic"
    64#pragma GCC diagnostic warning "-Wignored-qualifiers"
    65#pragma GCC diagnostic warning "-Wshadow"
    66#pragma GCC diagnostic warning "-Wnon-virtual-dtor"
    67#pragma GCC diagnostic warning "-Wctor-dtor-privacy"
    68#pragma GCC diagnostic warning "-Woverflow"
    69
    70
    71using operations_research::TimeLimit;
    72using operations_research::glop::BasisState;
    73using operations_research::glop::ColIndex;
    74using operations_research::glop::ColIndexVector;
    75using operations_research::glop::ConstraintStatus;
    76using operations_research::glop::ConstraintStatusColumn;
    77using operations_research::glop::DenseBooleanColumn;
    78using operations_research::glop::DenseBooleanRow;
    79using operations_research::glop::DenseColumn;
    80using operations_research::glop::DenseRow;
    81using operations_research::glop::SparseColumn;
    82using operations_research::glop::ScatteredColumn;
    83using operations_research::glop::ScatteredColumnIterator;
    84using operations_research::glop::SparseMatrix;
    85using operations_research::glop::Fractional;
    86using operations_research::glop::GetProblemStatusString;
    87using operations_research::glop::ProblemStatus;
    88using operations_research::glop::RowIndex;
    89using operations_research::glop::ScatteredRow;
    90using operations_research::glop::ScatteredRowIterator;
    91using operations_research::glop::VariableStatus;
    92using operations_research::glop::VariableStatusRow;
    93using operations_research::MPModelProto;
    94
    95/** LP interface */
    96struct SCIP_LPi
    97{
    98 operations_research::glop::LinearProgram* linear_program; /**< the linear program */
    99 operations_research::glop::LinearProgram* scaled_lp; /**< scaled linear program */
    100 operations_research::glop::RevisedSimplex* solver; /**< direct reference to the revised simplex, not passing through lp_solver */
    101 operations_research::glop::GlopParameters* parameters; /**< parameters */
    102 operations_research::glop::LpScalingHelper* scaler; /**< scaler auxiliary class */
    103
    104 /* the following is used by SCIPlpiWasSolved() */
    107
    108 /* store the values of some parameters in order to be able to return them */
    109 bool lp_info; /**< whether additional output is turned on */
    110 SCIP_PRICING pricing; /**< SCIP pricing setting */
    111 bool from_scratch; /**< store whether basis is ignored for next solving call */
    112 int numthreads; /**< number of threads used to solve the LP (0 = automatic) */
    113 SCIP_Real conditionlimit; /**< maximum condition number of LP basis counted as stable (-1.0: no limit) */
    114 bool checkcondition; /**< Should condition number of LP basis be checked for stability? */
    115 int timing; /**< type of timer (1 - cpu, 2 - wallclock, 0 - off) */
    116
    117 /* other data */
    118 SCIP_Longint niterations; /**< number of iterations used */
    119
    120 /* Temporary vectors allocated here for speed. This gain is non-negligible
    121 * because in many situations, only a few entries of these vectors are
    122 * inspected (hypersparsity) and allocating them is in O(num_rows) or
    123 * O(num_cols) instead of O(num_non_zeros) to read/clear them. */
    124 ScatteredRow* tmp_row; /**< temporary vector */
    125 ScatteredColumn* tmp_column; /**< temporary vector */
    126};
    127
    128/** uncomment to turn off scaling */
    129/* #define NOSCALING */
    130
    131/** define feasibility check to possibly reoptimize: 0: no check, 1: completely new check, 2: check unscaled variable and activity values */
    132#ifdef NOSCALING
    133#define UNSCALEDFEAS_CHECK 0
    134#else
    135#define UNSCALEDFEAS_CHECK 2
    136#endif
    137
    138
    139/*
    140 * LP Interface Methods
    141 */
    142
    143char* initGlopName( );
    144
    145static char* glopname = initGlopName( );
    146
    148{
    149 glopname = new char[100];
    150 (void) snprintf(glopname, 100, "Glop %d.%d", operations_research::OrToolsMajorVersion(), operations_research::OrToolsMinorVersion());
    151 return glopname;
    152}
    153
    154/** gets name and version of LP solver */
    156 void
    157 )
    158{
    159 return glopname;
    160}
    161
    162/** gets description of LP solver (developer, webpage, ...) */
    164 void
    165 )
    166{
    167 return "Glop Linear Solver, developed by Google, part of OR-Tools (developers.google.com/optimization)";
    168}
    169
    170/** gets pointer for LP solver - use only with great care */
    172 SCIP_LPI* lpi /**< pointer to an LP interface structure */
    173 )
    174{
    175 assert( lpi != NULL );
    176
    177 SCIPerrorMessage("SCIPlpiGetSolverPointer() has not been implemented yet.\n");
    178
    179 return NULL;
    180}
    181
    182/** pass integrality information to LP solver */
    184 SCIP_LPI* lpi, /**< pointer to an LP interface structure */
    185 int ncols, /**< length of integrality array */
    186 int* intInfo /**< integrality array (0: continuous, 1: integer). May be NULL iff ncols is 0. */
    187 )
    188{
    189 assert( lpi != NULL );
    190 assert( lpi->linear_program != NULL );
    191 assert( ncols == 0 || ncols == lpi->linear_program->num_variables().value() );
    192
    193 /* Pass on integrality information (currently not used by Glop). */
    194 for (ColIndex col(0); col < ColIndex(ncols); ++col)
    195 {
    196 assert( intInfo != NULL );
    197 int info = intInfo[col.value()];
    198 assert( info == 0 || info == 1 );
    199 if ( info == 0 )
    200 lpi->linear_program->SetVariableType(col, operations_research::glop::LinearProgram::VariableType::CONTINUOUS);
    201 else
    202 lpi->linear_program->SetVariableType(col, operations_research::glop::LinearProgram::VariableType::INTEGER);
    203 }
    204
    205 return SCIP_OKAY;
    206}
    207
    208/** informs about availability of a primal simplex solving method */
    210 void
    211 )
    212{
    213 return TRUE;
    214}
    215
    216/** informs about availability of a dual simplex solving method */
    218 void
    219 )
    220{
    221 return TRUE;
    222}
    223
    224/** informs about availability of a barrier solving method */
    226 void
    227 )
    228{
    229 return FALSE;
    230}
    231
    232
    233/*
    234 * LPI Creation and Destruction Methods
    235 */
    236
    237/**@name LPI Creation and Destruction Methods */
    238/**@{ */
    239
    240/** creates an LP problem object */
    242 SCIP_LPI** lpi, /**< pointer to an LP interface structure */
    243 SCIP_MESSAGEHDLR* messagehdlr, /**< message handler to use for printing messages, or NULL */
    244 const char* name, /**< problem name */
    245 SCIP_OBJSEN objsen /**< objective sense */
    246 )
    247{
    248 assert( lpi != NULL );
    249 assert( name != NULL );
    250
    251 /* Initilialize memory. */
    253 (*lpi)->linear_program = new operations_research::glop::LinearProgram();
    254 (*lpi)->scaled_lp = new operations_research::glop::LinearProgram();
    255 (*lpi)->solver = new operations_research::glop::RevisedSimplex();
    256 (*lpi)->parameters = new operations_research::glop::GlopParameters();
    257 (*lpi)->scaler = new operations_research::glop::LpScalingHelper();
    258
    259 /* Set problem name and objective direction. */
    260 (*lpi)->linear_program->SetName(std::string(name));
    261 SCIP_CALL( SCIPlpiChgObjsen(*lpi, objsen) );
    262
    263 (*lpi)->from_scratch = false;
    264 (*lpi)->lp_info = false;
    265 (*lpi)->pricing = SCIP_PRICING_LPIDEFAULT;
    266 (*lpi)->lp_modified_since_last_solve = true;
    267 (*lpi)->lp_time_limit_was_reached = false;
    268 (*lpi)->conditionlimit = -1.0;
    269 (*lpi)->checkcondition = false;
    270 (*lpi)->niterations = 0LL;
    271
    272 (*lpi)->tmp_row = new ScatteredRow();
    273 (*lpi)->tmp_column = new ScatteredColumn();
    274
    275#ifdef NOSCALING
    276 (*lpi)->parameters->set_use_scaling(false);
    277#endif
    278
    279 return SCIP_OKAY;
    280}
    281
    282/** deletes an LP problem object */
    284 SCIP_LPI** lpi /**< pointer to an LP interface structure */
    285 )
    286{
    287 SCIPdebugMessage("SCIPlpiFree\n");
    288
    289 delete (*lpi)->scaler;
    290 delete (*lpi)->parameters;
    291 delete (*lpi)->solver;
    292 delete (*lpi)->scaled_lp;
    293 delete (*lpi)->linear_program;
    294
    295 delete (*lpi)->tmp_row;
    296 delete (*lpi)->tmp_column;
    297
    298 BMSfreeMemory(lpi);
    299
    300 return SCIP_OKAY;
    301}
    302
    303/**@} */
    304
    305
    306
    307
    308/*
    309 * Modification Methods
    310 */
    311
    312/**@name Modification Methods */
    313/**@{ */
    314
    315/** copies LP data with column matrix into LP solver */
    317 SCIP_LPI* lpi, /**< LP interface structure */
    318 SCIP_OBJSEN objsen, /**< objective sense */
    319 int ncols, /**< number of columns */
    320 const SCIP_Real* obj, /**< objective function values of columns */
    321 const SCIP_Real* lb, /**< lower bounds of columns */
    322 const SCIP_Real* ub, /**< upper bounds of columns */
    323 char** colnames, /**< column names, or NULL */
    324 int nrows, /**< number of rows */
    325 const SCIP_Real* lhs, /**< left hand sides of rows */
    326 const SCIP_Real* rhs, /**< right hand sides of rows */
    327 char** rownames, /**< row names, or NULL */
    328 int nnonz, /**< number of nonzero elements in the constraint matrix */
    329 const int* beg, /**< start index of each column in ind- and val-array */
    330 const int* ind, /**< row indices of constraint matrix entries */
    331 const SCIP_Real* val /**< values of constraint matrix entries */
    332 )
    333{
    334 assert( lpi != NULL );
    335 assert( lpi->linear_program != NULL );
    336 assert( obj != NULL );
    337 assert( lb != NULL );
    338 assert( ub != NULL );
    339 assert( beg != NULL );
    340 assert( ind != NULL );
    341 assert( val != NULL );
    342
    343 lpi->linear_program->Clear();
    344 SCIP_CALL( SCIPlpiAddRows(lpi, nrows, lhs, rhs, rownames, 0, NULL, NULL, NULL) );
    345 SCIP_CALL( SCIPlpiAddCols(lpi, ncols, obj, lb, ub, colnames, nnonz, beg, ind, val) );
    346 SCIP_CALL( SCIPlpiChgObjsen(lpi, objsen) );
    347
    348 return SCIP_OKAY;
    349}
    350
    351/** adds columns to the LP */
    353 SCIP_LPI* lpi, /**< LP interface structure */
    354 int ncols, /**< number of columns to be added */
    355 const SCIP_Real* obj, /**< objective function values of new columns */
    356 const SCIP_Real* lb, /**< lower bounds of new columns */
    357 const SCIP_Real* ub, /**< upper bounds of new columns */
    358 char** colnames, /**< column names, or NULL */
    359 int nnonz, /**< number of nonzero elements to be added to the constraint matrix */
    360 const int* beg, /**< start index of each column in ind- and val-array, or NULL if nnonz == 0 */
    361 const int* ind, /**< row indices of constraint matrix entries, or NULL if nnonz == 0 */
    362 const SCIP_Real* val /**< values of constraint matrix entries, or NULL if nnonz == 0 */
    363 )
    364{
    365 assert( lpi != NULL );
    366 assert( lpi->linear_program != NULL );
    367 assert( obj != NULL );
    368 assert( lb != NULL );
    369 assert( ub != NULL );
    370 assert( nnonz >= 0) ;
    371 assert( ncols >= 0) ;
    372
    373 SCIPdebugMessage("adding %d columns with %d nonzeros.\n", ncols, nnonz);
    374
    375 /* @todo add names */
    376 if ( nnonz > 0 )
    377 {
    378 assert( beg != NULL );
    379 assert( ind != NULL );
    380 assert( val != NULL );
    381 assert( ncols > 0 );
    382
    383#ifndef NDEBUG
    384 /* perform check that no new rows are added */
    385 RowIndex num_rows = lpi->linear_program->num_constraints();
    386 for (int j = 0; j < nnonz; ++j)
    387 {
    388 assert( 0 <= ind[j] && ind[j] < num_rows.value() );
    389 assert( val[j] != 0.0 );
    390 }
    391#endif
    392
    393 int nz = 0;
    394 for (int i = 0; i < ncols; ++i)
    395 {
    396 const ColIndex col = lpi->linear_program->CreateNewVariable();
    397 lpi->linear_program->SetVariableBounds(col, lb[i], ub[i]);
    398 lpi->linear_program->SetObjectiveCoefficient(col, obj[i]);
    399 const int end = (nnonz == 0 || i == ncols - 1) ? nnonz : beg[i + 1];
    400 while ( nz < end )
    401 {
    402 lpi->linear_program->SetCoefficient(RowIndex(ind[nz]), col, val[nz]);
    403 ++nz;
    404 }
    405 }
    406 assert( nz == nnonz );
    407 }
    408 else
    409 {
    410 for (int i = 0; i < ncols; ++i)
    411 {
    412 const ColIndex col = lpi->linear_program->CreateNewVariable();
    413 lpi->linear_program->SetVariableBounds(col, lb[i], ub[i]);
    414 lpi->linear_program->SetObjectiveCoefficient(col, obj[i]);
    415 }
    416 }
    417
    419
    420 return SCIP_OKAY;
    421}
    422
    423/** deletes all columns in the given range from LP */
    425 SCIP_LPI* lpi, /**< LP interface structure */
    426 int firstcol, /**< first column to be deleted */
    427 int lastcol /**< last column to be deleted */
    428 )
    429{
    430 assert(lpi != NULL);
    431 assert(lpi->linear_program != NULL);
    432 assert(firstcol >= 0);
    433 assert(lastcol < lpi->linear_program->num_variables());
    434 assert(firstcol <= lastcol + 1);
    435
    436 SCIPdebugMessage("deleting columns %d to %d.\n", firstcol, lastcol);
    437
    438 // handle empty range
    439 if( firstcol > lastcol )
    440 return SCIP_OKAY;
    441
    442 const ColIndex num_cols = lpi->linear_program->num_variables();
    443 DenseBooleanRow columns_to_delete(num_cols, false);
    444 for( int i = firstcol; i <= lastcol; ++i )
    445 columns_to_delete[ColIndex(i)] = true;
    446
    447 lpi->linear_program->DeleteColumns(columns_to_delete);
    449
    450 return SCIP_OKAY;
    451}
    452
    453/** deletes columns from SCIP_LP; the new position of a column must not be greater that its old position */
    455 SCIP_LPI* lpi, /**< LP interface structure */
    456 int* dstat /**< deletion status of columns
    457 * input: 1 if column should be deleted, 0 if not
    458 * output: new position of column, -1 if column was deleted */
    459 )
    460{
    461 assert( lpi != NULL );
    462 assert( lpi->linear_program != NULL );
    463 assert( dstat != NULL );
    464
    465 const ColIndex num_cols = lpi->linear_program->num_variables();
    466 DenseBooleanRow columns_to_delete(num_cols, false);
    467 int new_index = 0;
    468 int num_deleted_columns = 0;
    469 for (ColIndex col(0); col < num_cols; ++col)
    470 {
    471 int i = col.value();
    472 if ( dstat[i] == 1 )
    473 {
    474 columns_to_delete[col] = true;
    475 dstat[i] = -1;
    476 ++num_deleted_columns;
    477 }
    478 else
    479 dstat[i] = new_index++;
    480 }
    481 SCIPdebugMessage("SCIPlpiDelColset: deleting %d columns.\n", num_deleted_columns);
    482 lpi->linear_program->DeleteColumns(columns_to_delete);
    484
    485 return SCIP_OKAY;
    486}
    487
    488/** adds rows to the LP */
    490 SCIP_LPI* lpi, /**< LP interface structure */
    491 int nrows, /**< number of rows to be added */
    492 const SCIP_Real* lhs, /**< left hand sides of new rows */
    493 const SCIP_Real* rhs, /**< right hand sides of new rows */
    494 char** rownames, /**< row names, or NULL */
    495 int nnonz, /**< number of nonzero elements to be added to the constraint matrix */
    496 const int* beg, /**< start index of each row in ind- and val-array, or NULL if nnonz == 0 */
    497 const int* ind, /**< column indices of constraint matrix entries, or NULL if nnonz == 0 */
    498 const SCIP_Real* val /**< values of constraint matrix entries, or NULL if nnonz == 0 */
    499 )
    500{
    501 assert( lpi != NULL );
    502 assert( lpi->linear_program != NULL );
    503 assert( lhs != NULL );
    504 assert( rhs != NULL );
    505 assert( nnonz >= 0) ;
    506 assert( nrows >= 0) ;
    507
    508 SCIPdebugMessage("adding %d rows with %d nonzeros.\n", nrows, nnonz);
    509
    510 /* @todo add names */
    511 if ( nnonz > 0 )
    512 {
    513 assert( beg != NULL );
    514 assert( ind != NULL );
    515 assert( val != NULL );
    516 assert( nrows > 0 );
    517
    518#ifndef NDEBUG
    519 /* perform check that no new columns are added - this is likely to be a mistake */
    520 const ColIndex num_cols = lpi->linear_program->num_variables();
    521 for (int j = 0; j < nnonz; ++j)
    522 {
    523 assert( val[j] != 0.0 );
    524 assert( 0 <= ind[j] && ind[j] < num_cols.value() );
    525 }
    526#endif
    527
    528 int nz = 0;
    529 for (int i = 0; i < nrows; ++i)
    530 {
    531 const RowIndex row = lpi->linear_program->CreateNewConstraint();
    532 lpi->linear_program->SetConstraintBounds(row, lhs[i], rhs[i]);
    533 const int end = (nnonz == 0 || i == nrows - 1) ? nnonz : beg[i + 1];
    534 while ( nz < end )
    535 {
    536 lpi->linear_program->SetCoefficient(row, ColIndex(ind[nz]), val[nz]);
    537 ++nz;
    538 }
    539 }
    540 assert( nz == nnonz );
    541 }
    542 else
    543 {
    544 for (int i = 0; i < nrows; ++i)
    545 {
    546 const RowIndex row = lpi->linear_program->CreateNewConstraint();
    547 lpi->linear_program->SetConstraintBounds(row, lhs[i], rhs[i]);
    548 }
    549 }
    550
    552
    553 return SCIP_OKAY;
    554}
    555
    556/** delete rows from LP and update the current basis */
    557static
    559 SCIP_LPI* lpi, /**< LP interface structure */
    560 const DenseBooleanColumn& rows_to_delete /**< array to mark rows that should be deleted */
    561 )
    562{
    563 const RowIndex num_rows = lpi->linear_program->num_constraints();
    564 const ColIndex num_cols = lpi->linear_program->num_variables();
    565
    566 /* try to repair basis status if problem size has not changed before */
    567 BasisState state = lpi->solver->GetState();
    568 if ( state.statuses.size() == num_cols.value() + num_rows.value() )
    569 {
    570 /* Shift the status of the non-deleted rows. Note that if the deleted rows where part of the basis (i.e., constraint
    571 * not tight), then we should be left with a correct basis afterward. This should be the most common use case in SCIP. */
    572 ColIndex new_size = num_cols;
    573 for (RowIndex row(0); row < num_rows; ++row)
    574 {
    575 if ( rows_to_delete[row] )
    576 continue;
    577 state.statuses[new_size++] = state.statuses[num_cols + RowToColIndex(row)];
    578 }
    579 state.statuses.resize(new_size);
    580 lpi->solver->LoadStateForNextSolve(state);
    581 }
    582
    583 lpi->linear_program->DeleteRows(rows_to_delete);
    585}
    586
    587/** deletes all rows in the given range from LP */
    589 SCIP_LPI* lpi, /**< LP interface structure */
    590 int firstrow, /**< first row to be deleted */
    591 int lastrow /**< last row to be deleted */
    592 )
    593{
    594 assert(lpi != NULL);
    595 assert(lpi->linear_program != NULL);
    596 assert(firstrow >= 0);
    597 assert(lastrow < lpi->linear_program->num_constraints());
    598 assert(firstrow <= lastrow + 1);
    599
    600 SCIPdebugMessage("deleting rows %d to %d.\n", firstrow, lastrow);
    601
    602 // handle empty range
    603 if( firstrow > lastrow )
    604 return SCIP_OKAY;
    605
    606 const RowIndex num_rows = lpi->linear_program->num_constraints();
    607 DenseBooleanColumn rows_to_delete(num_rows, false);
    608 for( int i = firstrow; i <= lastrow; ++i )
    609 rows_to_delete[RowIndex(i)] = true;
    610
    611 deleteRowsAndUpdateCurrentBasis(lpi, rows_to_delete);
    612
    613 return SCIP_OKAY;
    614}
    615
    616/** deletes rows from SCIP_LP; the new position of a row must not be greater that its old position */
    618 SCIP_LPI* lpi, /**< LP interface structure */
    619 int* dstat /**< deletion status of rows
    620 * input: 1 if row should be deleted, 0 if not
    621 * output: new position of row, -1 if row was deleted */
    622 )
    623{
    624 assert( lpi != NULL );
    625 assert( lpi->linear_program != NULL );
    626
    627 const RowIndex num_rows = lpi->linear_program->num_constraints();
    628 DenseBooleanColumn rows_to_delete(num_rows, false);
    629 int new_index = 0;
    630 int num_deleted_rows = 0;
    631 for (RowIndex row(0); row < num_rows; ++row)
    632 {
    633 int i = row.value();
    634 if (dstat[i] == 1)
    635 {
    636 rows_to_delete[row] = true;
    637 dstat[i] = -1;
    638 ++num_deleted_rows;
    639 }
    640 else
    641 dstat[i] = new_index++;
    642 }
    643
    644 SCIPdebugMessage("SCIPlpiDelRowset: deleting %d rows.\n", num_deleted_rows);
    645 deleteRowsAndUpdateCurrentBasis(lpi, rows_to_delete);
    646
    647 return SCIP_OKAY;
    648}
    649
    650/** clears the whole LP */
    652 SCIP_LPI* lpi /**< LP interface structure */
    653 )
    654{
    655 assert( lpi != NULL );
    656 assert( lpi->linear_program != NULL );
    657
    658 SCIPdebugMessage("SCIPlpiClear\n");
    659
    660 lpi->linear_program->Clear();
    662
    663 return SCIP_OKAY;
    664}
    665
    666/** changes lower and upper bounds of columns */
    668 SCIP_LPI* lpi, /**< LP interface structure */
    669 int ncols, /**< number of columns to change bounds for */
    670 const int* ind, /**< column indices or NULL if ncols is zero */
    671 const SCIP_Real* lb, /**< values for the new lower bounds or NULL if ncols is zero */
    672 const SCIP_Real* ub /**< values for the new upper bounds or NULL if ncols is zero */
    673 )
    674{
    675 assert( lpi != NULL );
    676 assert( lpi->linear_program != NULL );
    677 assert( ncols == 0 || (ind != NULL && lb != NULL && ub != NULL) );
    678
    679 SCIPdebugMessage("changing %d bounds.\n", ncols);
    680 if ( ncols <= 0 )
    681 return SCIP_OKAY;
    682
    683 for (int i = 0; i < ncols; ++i)
    684 {
    685 SCIPdebugMessage(" col %d: [%g,%g]\n", ind[i], lb[i], ub[i]);
    686
    687 if ( SCIPlpiIsInfinity(lpi, lb[i]) )
    688 {
    689 SCIPerrorMessage("LP Error: fixing lower bound for variable %d to infinity.\n", ind[i]);
    690 return SCIP_LPERROR;
    691 }
    692 if ( SCIPlpiIsInfinity(lpi, -ub[i]) )
    693 {
    694 SCIPerrorMessage("LP Error: fixing upper bound for variable %d to -infinity.\n", ind[i]);
    695 return SCIP_LPERROR;
    696 }
    697
    698 lpi->linear_program->SetVariableBounds(ColIndex(ind[i]), lb[i], ub[i]);
    699 }
    701
    702 return SCIP_OKAY;
    703}
    704
    705/** changes left and right hand sides of rows */
    707 SCIP_LPI* lpi, /**< LP interface structure */
    708 int nrows, /**< number of rows to change sides for */
    709 const int* ind, /**< row indices */
    710 const SCIP_Real* lhs, /**< new values for left hand sides */
    711 const SCIP_Real* rhs /**< new values for right hand sides */
    712 )
    713{
    714 assert( lpi != NULL );
    715 assert( lpi->linear_program != NULL );
    716
    717 if ( nrows <= 0 )
    718 return SCIP_OKAY;
    719
    720 SCIPdebugMessage("changing %d sides\n", nrows);
    721
    722 for (int i = 0; i < nrows; ++i)
    723 lpi->linear_program->SetConstraintBounds(RowIndex(ind[i]), lhs[i], rhs[i]);
    724
    726
    727 return SCIP_OKAY;
    728}
    729
    730/** changes a single coefficient */
    732 SCIP_LPI* lpi, /**< LP interface structure */
    733 int row, /**< row number of coefficient to change */
    734 int col, /**< column number of coefficient to change */
    735 SCIP_Real newval /**< new value of coefficient */
    736 )
    737{
    738 assert( lpi != NULL );
    739 assert( lpi->linear_program != NULL );
    740 assert( 0 <= row && row < lpi->linear_program->num_constraints().value() );
    741 assert( 0 <= col && col < lpi->linear_program->num_variables().value() );
    742
    743 SCIPdebugMessage("Set coefficient (%d,%d) to %f.\n", row, col, newval);
    744 lpi->linear_program->CleanUp();
    745 assert( lpi->linear_program->IsCleanedUp() );
    746 lpi->linear_program->SetCoefficient(RowIndex(row), ColIndex(col), newval);
    747
    749
    750 return SCIP_OKAY;
    751}
    752
    753/** changes the objective sense */
    755 SCIP_LPI* lpi, /**< LP interface structure */
    756 SCIP_OBJSEN objsen /**< new objective sense */
    757 )
    758{
    759 assert( lpi != NULL);
    760 assert( lpi->linear_program != NULL);
    761
    762 SCIPdebugMessage("changing objective sense to %d\n", objsen);
    763
    764 switch (objsen)
    765 {
    767 lpi->linear_program->SetMaximizationProblem(true);
    768 break;
    770 lpi->linear_program->SetMaximizationProblem(false);
    771 break;
    772 }
    774
    775 return SCIP_OKAY;
    776}
    777
    778/** changes objective values of columns in the LP */
    780 SCIP_LPI* lpi, /**< LP interface structure */
    781 int ncols, /**< number of columns to change objective value for */
    782 const int* ind, /**< column indices to change objective value for */
    783 const SCIP_Real* obj /**< new objective values for columns */
    784 )
    785{
    786 assert( lpi != NULL );
    787 assert( lpi->linear_program != NULL );
    788 assert( ind != NULL );
    789 assert( obj != NULL );
    790
    791 SCIPdebugMessage("changing %d objective values\n", ncols);
    792
    793 for (int i = 0; i < ncols; ++i)
    794 lpi->linear_program->SetObjectiveCoefficient(ColIndex(ind[i]), obj[i]);
    795
    797
    798 return SCIP_OKAY;
    799}
    800
    801/** multiplies a row with a non-zero scalar; for negative scalars, the row's sense is switched accordingly */
    803 SCIP_LPI* lpi, /**< LP interface structure */
    804 int row, /**< row number to scale */
    805 SCIP_Real scaleval /**< scaling multiplier */
    806 )
    807{
    808 SCIP_Real* vals;
    809 SCIP_Real lhs;
    810 SCIP_Real rhs;
    811 int nnonz;
    812 int* inds;
    813 int beg;
    814
    815 assert( lpi != NULL );
    816 assert( lpi->linear_program != NULL );
    817
    818 SCIPdebugMessage("Scale row %d by %f.\n", row, scaleval);
    819
    820 /* alloc memory */
    821 ColIndex num_cols = lpi->linear_program->num_variables();
    822
    823 SCIP_ALLOC( BMSallocMemoryArray(&inds, num_cols.value()) );
    824 SCIP_ALLOC( BMSallocMemoryArray(&vals, num_cols.value()) );
    825
    826 /* get the row */
    827 SCIP_CALL( SCIPlpiGetRows(lpi, row, row, &lhs, &rhs, &nnonz, &beg, inds, vals) );
    828
    829 /* scale row coefficients */
    830 for (int j = 0; j < nnonz; ++j)
    831 {
    832 SCIP_CALL( SCIPlpiChgCoef(lpi, row, inds[j], vals[j] * scaleval) );
    833 }
    834 BMSfreeMemoryArray(&vals);
    835 BMSfreeMemoryArray(&inds);
    836
    837 /* scale row sides */
    838 if ( ! SCIPlpiIsInfinity(lpi, -lhs) )
    839 lhs *= scaleval;
    840 else if ( scaleval < 0.0 )
    841 lhs = SCIPlpiInfinity(lpi);
    842
    843 if ( ! SCIPlpiIsInfinity(lpi, rhs) )
    844 rhs *= scaleval;
    845 else if ( scaleval < 0.0 )
    846 rhs = -SCIPlpiInfinity(lpi);
    847
    848 if ( scaleval > 0.0 )
    849 {
    850 SCIP_CALL( SCIPlpiChgSides(lpi, 1, &row, &lhs, &rhs) );
    851 }
    852 else
    853 {
    854 SCIP_CALL( SCIPlpiChgSides(lpi, 1, &row, &rhs, &lhs) );
    855 }
    856
    858
    859 return SCIP_OKAY;
    860}
    861
    862/** multiplies a column with a non-zero scalar; the objective value is multiplied with the scalar, and the bounds
    863 * are divided by the scalar; for negative scalars, the column's bounds are switched
    864 */
    866 SCIP_LPI* lpi, /**< LP interface structure */
    867 int col, /**< column number to scale */
    868 SCIP_Real scaleval /**< scaling multiplier */
    869 )
    870{
    871 SCIP_Real* vals;
    872 SCIP_Real lb;
    873 SCIP_Real ub;
    874 SCIP_Real obj;
    875 int nnonz;
    876 int* inds;
    877 int beg;
    878
    879 assert( lpi != NULL );
    880 assert( lpi->linear_program != NULL );
    881
    882 SCIPdebugMessage("Scale column %d by %f.\n", col, scaleval);
    883
    884 /* alloc memory */
    885 RowIndex num_rows = lpi->linear_program->num_constraints();
    886
    887 SCIP_ALLOC( BMSallocMemoryArray(&inds, num_rows.value()) );
    888 SCIP_ALLOC( BMSallocMemoryArray(&vals, num_rows.value()) );
    889
    890 /* get the column */
    891 SCIP_CALL( SCIPlpiGetCols(lpi, col, col, &lb, &ub, &nnonz, &beg, inds, vals) );
    892
    893 /* scale column coefficients */
    894 for (int j = 0; j < nnonz; ++j)
    895 {
    896 SCIP_CALL( SCIPlpiChgCoef(lpi, col, inds[j], vals[j] * scaleval) );
    897 }
    898 BMSfreeMemoryArray(&vals);
    899 BMSfreeMemoryArray(&inds);
    900
    901 /* scale objective value */
    902 SCIP_CALL( SCIPlpiGetObj(lpi, col, col, &obj) );
    903 obj *= scaleval;
    904 SCIP_CALL( SCIPlpiChgObj(lpi, 1, &col, &obj) );
    905
    906 /* scale bound */
    907 if ( ! SCIPlpiIsInfinity(lpi, -lb) )
    908 lb *= scaleval;
    909 else if ( scaleval < 0.0 )
    910 lb = SCIPlpiInfinity(lpi);
    911
    912 if ( ! SCIPlpiIsInfinity(lpi, ub) )
    913 ub *= scaleval;
    914 else if ( scaleval < 0.0 )
    915 ub = -SCIPlpiInfinity(lpi);
    916
    917 if ( scaleval > 0.0 )
    918 {
    919 SCIP_CALL( SCIPlpiChgBounds(lpi, 1, &col, &lb, &ub) );
    920 }
    921 else
    922 {
    923 SCIP_CALL( SCIPlpiChgBounds(lpi, 1, &col, &ub, &lb) );
    924 }
    925
    926 return SCIP_OKAY;
    927}
    928
    929
    930/**@} */
    931
    932
    933
    934
    935/*
    936 * Data Accessing Methods
    937 */
    938
    939/**@name Data Accessing Methods */
    940/**@{ */
    941
    942/** gets the number of rows in the LP */
    944 SCIP_LPI* lpi, /**< LP interface structure */
    945 int* nrows /**< pointer to store the number of rows */
    946 )
    947{
    948 assert( lpi != NULL );
    949 assert( lpi->linear_program != NULL );
    950 assert( nrows != NULL );
    951
    952 SCIPdebugMessage("getting number of rows.\n");
    953
    954 *nrows = lpi->linear_program->num_constraints().value();
    955
    956 return SCIP_OKAY;
    957}
    958
    959/** gets the number of columns in the LP */
    961 SCIP_LPI* lpi, /**< LP interface structure */
    962 int* ncols /**< pointer to store the number of cols */
    963 )
    964{
    965 assert( lpi != NULL );
    966 assert( lpi->linear_program != NULL );
    967 assert( ncols != NULL );
    968
    969 SCIPdebugMessage("getting number of columns.\n");
    970
    971 *ncols = lpi->linear_program->num_variables().value();
    972
    973 return SCIP_OKAY;
    974}
    975
    976/** gets objective sense of the LP */
    978 SCIP_LPI* lpi, /**< LP interface structure */
    979 SCIP_OBJSEN* objsen /**< pointer to store objective sense */
    980 )
    981{
    982 assert( lpi != NULL );
    983 assert( lpi->linear_program != NULL );
    984 assert( objsen != NULL );
    985
    986 SCIPdebugMessage("getting objective sense.\n");
    987
    988 *objsen = lpi->linear_program->IsMaximizationProblem() ? SCIP_OBJSEN_MAXIMIZE : SCIP_OBJSEN_MINIMIZE;
    989
    990 return SCIP_OKAY;
    991}
    992
    993/** gets the number of nonzero elements in the LP constraint matrix */
    995 SCIP_LPI* lpi, /**< LP interface structure */
    996 int* nnonz /**< pointer to store the number of nonzeros */
    997 )
    998{
    999 assert( lpi != NULL );
    1000 assert( lpi->linear_program != NULL );
    1001 assert( nnonz != NULL );
    1002
    1003 SCIPdebugMessage("getting number of non-zeros.\n");
    1004
    1005 *nnonz = (int) lpi->linear_program->num_entries().value();
    1006
    1007 return SCIP_OKAY;
    1008}
    1009
    1010/** gets columns from LP problem object; the arrays have to be large enough to store all values
    1011 * Either both, lb and ub, have to be NULL, or both have to be non-NULL,
    1012 * either nnonz, beg, ind, and val have to be NULL, or all of them have to be non-NULL.
    1013 */
    1015 SCIP_LPI* lpi, /**< LP interface structure */
    1016 int firstcol, /**< first column to get from LP */
    1017 int lastcol, /**< last column to get from LP */
    1018 SCIP_Real* lb, /**< buffer to store the lower bound vector, or NULL */
    1019 SCIP_Real* ub, /**< buffer to store the upper bound vector, or NULL */
    1020 int* nnonz, /**< pointer to store the number of nonzero elements returned, or NULL */
    1021 int* beg, /**< buffer to store start index of each column in ind- and val-array, or NULL */
    1022 int* ind, /**< buffer to store row indices of constraint matrix entries, or NULL */
    1023 SCIP_Real* val /**< buffer to store values of constraint matrix entries, or NULL */
    1024 )
    1025{
    1026 assert(lpi != NULL);
    1027 assert(lpi->linear_program != NULL);
    1028 assert((lb != NULL && ub != NULL) || (lb == NULL && ub == NULL));
    1029 assert((nnonz != NULL && beg != NULL && ind != NULL && val != NULL) || (nnonz == NULL && beg == NULL && ind == NULL && val == NULL));
    1030 assert(firstcol >= 0);
    1031 assert(lastcol < lpi->linear_program->num_variables());
    1032 assert(firstcol <= lastcol + 1);
    1033
    1034 const DenseRow& tmplb = lpi->linear_program->variable_lower_bounds();
    1035 const DenseRow& tmpub = lpi->linear_program->variable_upper_bounds();
    1036
    1037 if( nnonz != NULL )
    1038 {
    1039 assert(beg != NULL);
    1040 assert(ind != NULL);
    1041 assert(val != NULL);
    1042
    1043 *nnonz = 0;
    1044 int index = 0;
    1045 for( ColIndex col(firstcol); col <= ColIndex(lastcol); ++col, ++index )
    1046 {
    1047 if( lb != NULL )
    1048 lb[index] = tmplb[col];
    1049 if( ub != NULL )
    1050 ub[index] = tmpub[col];
    1051
    1052 beg[index] = *nnonz;
    1053 const SparseColumn& column = lpi->linear_program->GetSparseColumn(col);
    1054 for( const SparseColumn::Entry& entry : column )
    1055 {
    1056 const RowIndex row = entry.row();
    1057 ind[*nnonz] = row.value();
    1058 val[*nnonz] = entry.coefficient();
    1059 ++(*nnonz);
    1060 }
    1061 }
    1062 }
    1063 else
    1064 {
    1065 int index = 0;
    1066 for( ColIndex col(firstcol); col <= ColIndex(lastcol); ++col, ++index )
    1067 {
    1068 if( lb != NULL )
    1069 lb[index] = tmplb[col];
    1070 if( ub != NULL )
    1071 ub[index] = tmpub[col];
    1072 }
    1073 }
    1074
    1075 return SCIP_OKAY;
    1076}
    1077
    1078/** gets rows from LP problem object; the arrays have to be large enough to store all values.
    1079 * Either both, lhs and rhs, have to be NULL, or both have to be non-NULL,
    1080 * either nnonz, beg, ind, and val have to be NULL, or all of them have to be non-NULL.
    1081 */
    1083 SCIP_LPI* lpi, /**< LP interface structure */
    1084 int firstrow, /**< first row to get from LP */
    1085 int lastrow, /**< last row to get from LP */
    1086 SCIP_Real* lhs, /**< buffer to store left hand side vector, or NULL */
    1087 SCIP_Real* rhs, /**< buffer to store right hand side vector, or NULL */
    1088 int* nnonz, /**< pointer to store the number of nonzero elements returned, or NULL */
    1089 int* beg, /**< buffer to store start index of each row in ind- and val-array, or NULL */
    1090 int* ind, /**< buffer to store column indices of constraint matrix entries, or NULL */
    1091 SCIP_Real* val /**< buffer to store values of constraint matrix entries, or NULL */
    1092 )
    1093{
    1094 assert(lpi != NULL);
    1095 assert(lpi->linear_program != NULL);
    1096 assert((lhs == NULL && rhs == NULL) || (rhs != NULL && lhs != NULL));
    1097 assert((nnonz != NULL && beg != NULL && ind != NULL && val != NULL) || (nnonz == NULL && beg == NULL && ind == NULL && val == NULL));
    1098 assert(firstrow >= 0);
    1099 assert(lastrow < lpi->linear_program->num_constraints());
    1100 assert(firstrow <= lastrow + 1);
    1101
    1102 const DenseColumn& tmplhs = lpi->linear_program->constraint_lower_bounds();
    1103 const DenseColumn& tmprhs = lpi->linear_program->constraint_upper_bounds();
    1104
    1105 if( nnonz != NULL )
    1106 {
    1107 assert(beg != NULL);
    1108 assert(ind != NULL);
    1109 assert(val != NULL);
    1110
    1111 const SparseMatrix& matrixtrans = lpi->linear_program->GetTransposeSparseMatrix();
    1112
    1113 *nnonz = 0;
    1114 int index = 0;
    1115 for( RowIndex row(firstrow); row <= RowIndex(lastrow); ++row, ++index )
    1116 {
    1117 if( lhs != NULL )
    1118 lhs[index] = tmplhs[row];
    1119 if( rhs != NULL )
    1120 rhs[index] = tmprhs[row];
    1121
    1122 beg[index] = *nnonz;
    1123 const SparseColumn& column = matrixtrans.column(ColIndex(row.value()));
    1124 for( const SparseColumn::Entry& entry : column )
    1125 {
    1126 const RowIndex rowidx = entry.row();
    1127 ind[*nnonz] = rowidx.value();
    1128 val[*nnonz] = entry.coefficient();
    1129 ++(*nnonz);
    1130 }
    1131 }
    1132 }
    1133 else
    1134 {
    1135 int index = 0;
    1136 for( RowIndex row(firstrow); row <= RowIndex(lastrow); ++row, ++index )
    1137 {
    1138 if( lhs != NULL )
    1139 lhs[index] = tmplhs[row];
    1140 if( rhs != NULL )
    1141 rhs[index] = tmprhs[row];
    1142 }
    1143 }
    1144
    1145 return SCIP_OKAY;
    1146}
    1147
    1148/** gets column names */
    1150 SCIP_LPI* lpi, /**< LP interface structure */
    1151 int firstcol, /**< first column to get name from LP */
    1152 int lastcol, /**< last column to get name from LP */
    1153 char** colnames, /**< pointers to column names (of size at least lastcol-firstcol+1) or NULL if namestoragesize is zero */
    1154 char* namestorage, /**< storage for col names or NULL if namestoragesize is zero */
    1155 int namestoragesize, /**< size of namestorage (if 0, storageleft returns the storage needed) */
    1156 int* storageleft /**< amount of storage left (if < 0 the namestorage was not big enough) or NULL if namestoragesize is zero */
    1157 )
    1158{
    1159 assert(lpi != NULL);
    1160 assert(lpi->linear_program != NULL);
    1161 assert(colnames != NULL || namestoragesize == 0);
    1162 assert(namestorage != NULL || namestoragesize == 0);
    1163 assert(namestoragesize >= 0);
    1164 assert(storageleft != NULL);
    1165 assert(firstcol >= 0);
    1166 assert(lastcol < lpi->linear_program->num_variables());
    1167 assert(firstcol <= lastcol + 1);
    1168
    1169 SCIPerrorMessage("SCIPlpiGetColNames() has not been implemented yet.\n");
    1170
    1171 return SCIP_NOTIMPLEMENTED;
    1172}
    1173
    1174/** gets row names */
    1176 SCIP_LPI* lpi, /**< LP interface structure */
    1177 int firstrow, /**< first row to get name from LP */
    1178 int lastrow, /**< last row to get name from LP */
    1179 char** rownames, /**< pointers to row names (of size at least lastrow-firstrow+1) or NULL if namestoragesize is zero */
    1180 char* namestorage, /**< storage for row names or NULL if namestoragesize is zero */
    1181 int namestoragesize, /**< size of namestorage (if 0, -storageleft returns the storage needed) */
    1182 int* storageleft /**< amount of storage left (if < 0 the namestorage was not big enough) or NULL if namestoragesize is zero */
    1183 )
    1184{
    1185 assert(lpi != NULL);
    1186 assert(lpi->linear_program != NULL);
    1187 assert(rownames != NULL || namestoragesize == 0);
    1188 assert(namestorage != NULL || namestoragesize == 0);
    1189 assert(namestoragesize >= 0);
    1190 assert(storageleft != NULL);
    1191 assert(firstrow >= 0);
    1192 assert(lastrow < lpi->linear_program->num_constraints());
    1193 assert(firstrow <= lastrow + 1);
    1194
    1195 SCIPerrorMessage("SCIPlpiGetRowNames() has not been implemented yet.\n");
    1196
    1197 return SCIP_NOTIMPLEMENTED;
    1198}
    1199
    1200/** gets objective coefficients from LP problem object */
    1202 SCIP_LPI* lpi, /**< LP interface structure */
    1203 int firstcol, /**< first column to get objective coefficient for */
    1204 int lastcol, /**< last column to get objective coefficient for */
    1205 SCIP_Real* vals /**< array to store objective coefficients */
    1206 )
    1207{
    1208 assert(lpi != NULL);
    1209 assert(lpi->linear_program != NULL);
    1210 assert(vals != NULL);
    1211 assert(firstcol >= 0);
    1212 assert(lastcol < lpi->linear_program->num_variables());
    1213 assert(firstcol <= lastcol + 1);
    1214
    1215 SCIPdebugMessage("getting objective values %d to %d\n", firstcol, lastcol);
    1216
    1217 int index = 0;
    1218 for( ColIndex col(firstcol); col <= ColIndex(lastcol); ++col )
    1219 {
    1220 vals[index] = lpi->linear_program->objective_coefficients()[col];
    1221 ++index;
    1222 }
    1223
    1224 return SCIP_OKAY;
    1225}
    1226
    1227/** gets current bounds from LP problem object */
    1229 SCIP_LPI* lpi, /**< LP interface structure */
    1230 int firstcol, /**< first column to get bounds for */
    1231 int lastcol, /**< last column to get bounds for */
    1232 SCIP_Real* lbs, /**< array to store lower bound values, or NULL */
    1233 SCIP_Real* ubs /**< array to store upper bound values, or NULL */
    1234 )
    1235{
    1236 assert(lpi != NULL);
    1237 assert(lpi->linear_program != NULL);
    1238 assert(firstcol >= 0);
    1239 assert(lastcol < lpi->linear_program->num_variables());
    1240 assert(firstcol <= lastcol + 1);
    1241
    1242 SCIPdebugMessage("getting bounds %d to %d\n", firstcol, lastcol);
    1243
    1244 int index = 0;
    1245 for( ColIndex col(firstcol); col <= ColIndex(lastcol); ++col )
    1246 {
    1247 if( lbs != NULL )
    1248 lbs[index] = lpi->linear_program->variable_lower_bounds()[col];
    1249
    1250 if( ubs != NULL )
    1251 ubs[index] = lpi->linear_program->variable_upper_bounds()[col];
    1252
    1253 ++index;
    1254 }
    1255
    1256 return SCIP_OKAY;
    1257}
    1258
    1259/** gets current row sides from LP problem object */
    1261 SCIP_LPI* lpi, /**< LP interface structure */
    1262 int firstrow, /**< first row to get sides for */
    1263 int lastrow, /**< last row to get sides for */
    1264 SCIP_Real* lhss, /**< array to store left hand side values, or NULL */
    1265 SCIP_Real* rhss /**< array to store right hand side values, or NULL */
    1266 )
    1267{
    1268 assert(lpi != NULL);
    1269 assert(lpi->linear_program != NULL);
    1270 assert(firstrow >= 0);
    1271 assert(lastrow < lpi->linear_program->num_constraints());
    1272 assert(firstrow <= lastrow + 1);
    1273
    1274 SCIPdebugMessage("getting row sides %d to %d\n", firstrow, lastrow);
    1275
    1276 int index = 0;
    1277 for( RowIndex row(firstrow); row <= RowIndex(lastrow); ++row )
    1278 {
    1279 if( lhss != NULL )
    1280 lhss[index] = lpi->linear_program->constraint_lower_bounds()[row];
    1281
    1282 if( rhss != NULL )
    1283 rhss[index] = lpi->linear_program->constraint_upper_bounds()[row];
    1284
    1285 ++index;
    1286 }
    1287
    1288 return SCIP_OKAY;
    1289}
    1290
    1291/** gets a single coefficient */
    1293 SCIP_LPI* lpi, /**< LP interface structure */
    1294 int row, /**< row number of coefficient */
    1295 int col, /**< column number of coefficient */
    1296 SCIP_Real* val /**< pointer to store the value of the coefficient */
    1297 )
    1298{
    1299 assert( lpi != NULL );
    1300 assert( lpi->linear_program != NULL );
    1301 assert( val != NULL );
    1302
    1303 /* quite slow method: possibly needs linear time if matrix is not sorted */
    1304 const SparseMatrix& matrix = lpi->linear_program->GetSparseMatrix();
    1305 *val = matrix.LookUpValue(RowIndex(row), ColIndex(col));
    1306
    1307 return SCIP_OKAY;
    1308}
    1309
    1310/**@} */
    1311
    1312
    1313
    1314
    1315/*
    1316 * Solving Methods
    1317 */
    1318
    1319/**@name Solving Methods */
    1320/**@{ */
    1321
    1322/** update scaled linear program */
    1323static
    1325 SCIP_LPI* lpi /**< LP interface structure */
    1326 )
    1327{
    1328 if ( ! lpi->lp_modified_since_last_solve )
    1329 return;
    1330
    1331 lpi->scaled_lp->PopulateFromLinearProgram(*lpi->linear_program);
    1332 lpi->scaled_lp->AddSlackVariablesWhereNecessary(false);
    1333
    1334 /* @todo: Avoid doing a copy if there is no scaling. */
    1335 /* @todo: Avoid rescaling if not much changed. */
    1336 if ( lpi->parameters->use_scaling() )
    1337 lpi->scaler->Scale(lpi->scaled_lp);
    1338 else
    1339 lpi->scaler->Clear();
    1340}
    1341
    1342/** check primal feasibility */
    1343static
    1345 SCIP_LPI* lpi /**< LP interface structure */
    1346 )
    1347{
    1348 assert( lpi != NULL );
    1349 assert( lpi->solver != NULL );
    1350 assert( lpi->linear_program != NULL );
    1351
    1352#if UNSCALEDFEAS_CHECK == 1
    1353 /* get unscaled solution */
    1354 const ColIndex num_cols = lpi->linear_program->num_variables();
    1355 DenseRow unscaledsol(num_cols);
    1356 for (ColIndex col = ColIndex(0); col < num_cols; ++col)
    1357 unscaledsol[col] = lpi->scaler->UnscaleVariableValue(col, lpi->solver->GetVariableValue(col));
    1358
    1359 /* if the solution is not feasible w.r.t. absolute tolerances, try to fix it in the unscaled problem */
    1360 const double feastol = lpi->parameters->primal_feasibility_tolerance();
    1361 return lpi->linear_program->SolutionIsLPFeasible(unscaledsol, feastol);
    1362
    1363#elif UNSCALEDFEAS_CHECK == 2
    1364 const double feastol = lpi->parameters->primal_feasibility_tolerance();
    1365
    1366 /* check bounds of unscaled solution */
    1367 const ColIndex num_cols = lpi->linear_program->num_variables();
    1368 for (ColIndex col = ColIndex(0); col < num_cols; ++col)
    1369 {
    1370 const Fractional val = lpi->scaler->UnscaleVariableValue(col, lpi->solver->GetVariableValue(col));
    1371 const Fractional lb = lpi->linear_program->variable_lower_bounds()[col];
    1372 if ( val < lb - feastol )
    1373 return false;
    1374 const Fractional ub = lpi->linear_program->variable_upper_bounds()[col];
    1375 if ( val > ub + feastol )
    1376 return false;
    1377 }
    1378
    1379 /* check activities of unscaled solution */
    1380 const RowIndex num_rows = lpi->linear_program->num_constraints();
    1381 for (RowIndex row(0); row < num_rows; ++row)
    1382 {
    1383 const Fractional val = lpi->scaler->UnscaleConstraintActivity(row, lpi->solver->GetConstraintActivity(row));
    1384 const Fractional lhs = lpi->linear_program->constraint_lower_bounds()[row];
    1385 if ( val < lhs - feastol )
    1386 return false;
    1387 const Fractional rhs = lpi->linear_program->constraint_upper_bounds()[row];
    1388 if ( val > rhs + feastol )
    1389 return false;
    1390 }
    1391#endif
    1392
    1393 return true;
    1394}
    1395
    1396/** common function between the two LPI Solve() functions */
    1397static
    1399 SCIP_LPI* lpi, /**< LP interface structure */
    1400 bool recursive, /**< Is this a recursive call? */
    1401 std::unique_ptr<TimeLimit>& time_limit /**< time limit */
    1402 )
    1403{
    1404 assert( lpi != NULL );
    1405 assert( lpi->solver != NULL );
    1406 assert( lpi->parameters != NULL );
    1407
    1408 updateScaledLP(lpi);
    1409
    1410 lpi->solver->SetParameters(*(lpi->parameters));
    1411 lpi->lp_time_limit_was_reached = false;
    1412
    1413 /* possibly ignore warm start information for next solve */
    1414 if ( lpi->from_scratch )
    1415 lpi->solver->ClearStateForNextSolve();
    1416
    1417 if ( ! lpi->solver->Solve(*(lpi->scaled_lp), time_limit.get()).ok() )
    1418 {
    1419 return SCIP_LPERROR;
    1420 }
    1421 lpi->lp_time_limit_was_reached = time_limit->LimitReached();
    1422 if ( recursive )
    1423 lpi->niterations += (SCIP_Longint) lpi->solver->GetNumberOfIterations();
    1424 else
    1425 lpi->niterations = (SCIP_Longint) lpi->solver->GetNumberOfIterations();
    1426
    1427 SCIPdebugMessage("status=%s obj=%f iter=%ld.\n", GetProblemStatusString(lpi->solver->GetProblemStatus()).c_str(),
    1428 lpi->solver->GetObjectiveValue(), lpi->solver->GetNumberOfIterations());
    1429
    1430 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1431 if ( (status == ProblemStatus::PRIMAL_FEASIBLE || status == ProblemStatus::OPTIMAL) && lpi->parameters->use_scaling() )
    1432 {
    1433 if ( ! checkUnscaledPrimalFeasibility(lpi) )
    1434 {
    1435 SCIPdebugMessage("Solution not feasible w.r.t. absolute tolerance %g -> reoptimize.\n", lpi->parameters->primal_feasibility_tolerance());
    1436
    1437 /* Re-solve without scaling to try to fix the infeasibility. */
    1438 lpi->parameters->set_use_scaling(false);
    1439 lpi->lp_modified_since_last_solve = true;
    1440 SolveInternal(lpi, true, time_limit); /* inherit time limit, so used time is not reset; do not change iteration limit for resolve */
    1441 lpi->parameters->set_use_scaling(true);
    1442
    1443#ifdef SCIP_DEBUG
    1444 if ( ! checkUnscaledPrimalFeasibility(lpi) )
    1445 SCIPdebugMessage("Solution still not feasible after turning off scaling.\n");
    1446#endif
    1447 }
    1448 }
    1449
    1450 lpi->lp_modified_since_last_solve = false;
    1451
    1452 return SCIP_OKAY;
    1453}
    1454
    1455/** calls primal simplex to solve the LP */
    1457 SCIP_LPI* lpi /**< LP interface structure */
    1458 )
    1459{
    1460 assert( lpi != NULL );
    1461 assert( lpi->solver != NULL );
    1462 assert( lpi->linear_program != NULL );
    1463 assert( lpi->parameters != NULL );
    1464
    1465 SCIPdebugMessage("SCIPlpiSolvePrimal: %d rows, %d cols.\n", lpi->linear_program->num_constraints().value(), lpi->linear_program->num_variables().value());
    1466 std::unique_ptr<TimeLimit> time_limit = TimeLimit::FromParameters(*lpi->parameters);
    1467 lpi->niterations = 0;
    1468
    1469 lpi->parameters->set_use_dual_simplex(false);
    1470 return SolveInternal(lpi, false, time_limit);
    1471}
    1472
    1473/** calls dual simplex to solve the LP */
    1475 SCIP_LPI* lpi /**< LP interface structure */
    1476 )
    1477{
    1478 assert( lpi != NULL );
    1479 assert( lpi->solver != NULL );
    1480 assert( lpi->linear_program != NULL );
    1481 assert( lpi->parameters != NULL );
    1482
    1483 SCIPdebugMessage("SCIPlpiSolveDual: %d rows, %d cols.\n", lpi->linear_program->num_constraints().value(), lpi->linear_program->num_variables().value());
    1484 std::unique_ptr<TimeLimit> time_limit = TimeLimit::FromParameters(*lpi->parameters);
    1485 lpi->niterations = 0;
    1486
    1487 lpi->parameters->set_use_dual_simplex(true);
    1488 return SolveInternal(lpi, false, time_limit);
    1489}
    1490
    1491/** calls barrier or interior point algorithm to solve the LP with crossover to simplex basis */
    1493 SCIP_LPI* lpi, /**< LP interface structure */
    1494 SCIP_Bool crossover /**< perform crossover */
    1495 )
    1496{
    1497 assert( lpi != NULL );
    1498 assert( lpi->solver != NULL );
    1499 assert( lpi->linear_program != NULL );
    1500 assert( lpi->parameters != NULL );
    1501
    1502 SCIPerrorMessage("SCIPlpiSolveBarrier - Not supported.\n");
    1503
    1504 return SCIP_NOTIMPLEMENTED;
    1505}
    1506
    1507/** start strong branching */
    1509 SCIP_LPI* lpi /**< LP interface structure */
    1510 )
    1511{ /*lint --e{715}*/
    1512 assert( lpi != NULL );
    1513 assert( lpi->linear_program != NULL );
    1514 assert( lpi->solver != NULL );
    1515
    1516 updateScaledLP(lpi);
    1517
    1518 /* @todo Save state and do all the branching from there. */
    1519 return SCIP_OKAY;
    1520}
    1521
    1522/** end strong branching */
    1524 SCIP_LPI* lpi /**< LP interface structure */
    1525 )
    1526{ /*lint --e{715}*/
    1527 assert( lpi != NULL );
    1528 assert( lpi->linear_program != NULL );
    1529 assert( lpi->solver != NULL );
    1530
    1531 /* @todo Restore the saved state. */
    1532 return SCIP_OKAY;
    1533}
    1534
    1535/** determine whether the dual bound is valid */
    1536static
    1538 ProblemStatus status /**< status to be checked */
    1539 )
    1540{
    1541 return status == ProblemStatus::OPTIMAL || status == ProblemStatus::DUAL_FEASIBLE || status == ProblemStatus::DUAL_UNBOUNDED;
    1542}
    1543
    1544/** performs strong branching iterations */
    1545static
    1547 SCIP_LPI* lpi, /**< LP interface structure */
    1548 int col_index, /**< column to apply strong branching on */
    1549 SCIP_Real psol, /**< fractional current primal solution value of column */
    1550 int itlim, /**< iteration limit for strong branchings */
    1551 SCIP_Real* down, /**< stores dual bound after branching column down */
    1552 SCIP_Real* up, /**< stores dual bound after branching column up */
    1553 SCIP_Bool* downvalid, /**< stores whether the returned down value is a valid dual bound;
    1554 * otherwise, it can only be used as an estimate value */
    1555 SCIP_Bool* upvalid, /**< stores whether the returned up value is a valid dual bound;
    1556 * otherwise, it can only be used as an estimate value */
    1557 int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
    1558 )
    1559{
    1560 assert( lpi != NULL );
    1561 assert( lpi->scaled_lp != NULL );
    1562 assert( down != NULL );
    1563 assert( up != NULL );
    1564 assert( downvalid != NULL );
    1565 assert( upvalid != NULL );
    1566
    1567 SCIPdebugMessage("calling strongbranching on variable %d (%d iterations)\n", col_index, itlim);
    1568
    1569 /* We work on the scaled problem. */
    1570 const ColIndex col(col_index);
    1571 const Fractional lb = lpi->scaled_lp->variable_lower_bounds()[col];
    1572 const Fractional ub = lpi->scaled_lp->variable_upper_bounds()[col];
    1573 const double value = psol * lpi->scaler->VariableScalingFactor(col);
    1574
    1575 /* Configure solver. */
    1576
    1577 /* @todo Use the iteration limit once glop supports incrementality. */
    1578 int num_iterations = 0;
    1579 lpi->parameters->set_use_dual_simplex(true);
    1580 lpi->solver->SetParameters(*(lpi->parameters));
    1581 const Fractional eps = lpi->parameters->primal_feasibility_tolerance();
    1582
    1583 std::unique_ptr<TimeLimit> time_limit = TimeLimit::FromParameters(*lpi->parameters);
    1584
    1585 /* Down branch. */
    1586 const Fractional newub = EPSCEIL(value - 1.0, eps);
    1587 if ( newub >= lb - 0.5 )
    1588 {
    1589 lpi->scaled_lp->SetVariableBounds(col, lb, newub);
    1590
    1591 if ( lpi->solver->Solve(*(lpi->scaled_lp), time_limit.get()).ok() )
    1592 {
    1593 num_iterations += (int) lpi->solver->GetNumberOfIterations();
    1594 *down = lpi->solver->GetObjectiveValue();
    1595 *downvalid = IsDualBoundValid(lpi->solver->GetProblemStatus()) ? TRUE : FALSE;
    1596
    1597 SCIPdebugMessage("down: itlim=%d col=%d [%f,%f] obj=%f status=%d iter=%ld.\n", itlim, col_index, lb, EPSCEIL(value - 1.0, eps),
    1598 lpi->solver->GetObjectiveValue(), (int) lpi->solver->GetProblemStatus(), lpi->solver->GetNumberOfIterations());
    1599 }
    1600 else
    1601 {
    1602 SCIPerrorMessage("error during solve");
    1603 *down = 0.0;
    1604 *downvalid = FALSE;
    1605 }
    1606 }
    1607 else
    1608 {
    1609 if ( lpi->linear_program->IsMaximizationProblem() )
    1610 *down = lpi->parameters->objective_lower_limit();
    1611 else
    1612 *down = lpi->parameters->objective_upper_limit();
    1613 *downvalid = TRUE;
    1614 }
    1615
    1616 /* Up branch. */
    1617 const Fractional newlb = EPSFLOOR(value + 1.0, eps);
    1618 if ( newlb <= ub + 0.5 )
    1619 {
    1620 lpi->scaled_lp->SetVariableBounds(col, newlb, ub);
    1621
    1622 if ( lpi->solver->Solve(*(lpi->scaled_lp), time_limit.get()).ok() )
    1623 {
    1624 num_iterations += (int) lpi->solver->GetNumberOfIterations();
    1625 *up = lpi->solver->GetObjectiveValue();
    1626 *upvalid = IsDualBoundValid(lpi->solver->GetProblemStatus()) ? TRUE : FALSE;
    1627
    1628 SCIPdebugMessage("up: itlim=%d col=%d [%f,%f] obj=%f status=%d iter=%ld.\n", itlim, col_index, EPSFLOOR(value + 1.0, eps), ub,
    1629 lpi->solver->GetObjectiveValue(), (int) lpi->solver->GetProblemStatus(), lpi->solver->GetNumberOfIterations());
    1630 }
    1631 else
    1632 {
    1633 SCIPerrorMessage("error during solve");
    1634 *up = 0.0;
    1635 *upvalid = FALSE;
    1636 }
    1637 }
    1638 else
    1639 {
    1640 if (lpi->linear_program->IsMaximizationProblem())
    1641 *up = lpi->parameters->objective_lower_limit();
    1642 else
    1643 *up = lpi->parameters->objective_upper_limit();
    1644 *upvalid = TRUE;
    1645 }
    1646
    1647 /* Restore bound. */
    1648 lpi->scaled_lp->SetVariableBounds(col, lb, ub);
    1649 if ( iter != NULL )
    1650 *iter = num_iterations;
    1651
    1652 return SCIP_OKAY;
    1653}
    1654
    1655/** performs strong branching iterations on one @b fractional candidate */
    1657 SCIP_LPI* lpi, /**< LP interface structure */
    1658 int col_index, /**< column to apply strong branching on */
    1659 SCIP_Real psol, /**< fractional current primal solution value of column */
    1660 int itlim, /**< iteration limit for strong branchings */
    1661 SCIP_Real* down, /**< stores dual bound after branching column down */
    1662 SCIP_Real* up, /**< stores dual bound after branching column up */
    1663 SCIP_Bool* downvalid, /**< stores whether the returned down value is a valid dual bound;
    1664 * otherwise, it can only be used as an estimate value */
    1665 SCIP_Bool* upvalid, /**< stores whether the returned up value is a valid dual bound;
    1666 * otherwise, it can only be used as an estimate value */
    1667 int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
    1668 )
    1669{
    1670 assert( lpi != NULL );
    1671 assert( lpi->scaled_lp != NULL );
    1672 assert( down != NULL );
    1673 assert( up != NULL );
    1674 assert( downvalid != NULL );
    1675 assert( upvalid != NULL );
    1676
    1677 SCIPdebugMessage("calling strongbranching on fractional variable %d (%d iterations)\n", col_index, itlim);
    1678
    1679 SCIP_CALL( strongbranch(lpi, col_index, psol, itlim, down, up, downvalid, upvalid, iter) );
    1680
    1681 return SCIP_OKAY;
    1682}
    1683
    1684/** performs strong branching iterations on given @b fractional candidates */
    1686 SCIP_LPI* lpi, /**< LP interface structure */
    1687 int* cols, /**< columns to apply strong branching on */
    1688 int ncols, /**< number of columns */
    1689 SCIP_Real* psols, /**< fractional current primal solution values of columns */
    1690 int itlim, /**< iteration limit for strong branchings */
    1691 SCIP_Real* down, /**< stores dual bounds after branching columns down */
    1692 SCIP_Real* up, /**< stores dual bounds after branching columns up */
    1693 SCIP_Bool* downvalid, /**< stores whether the returned down values are valid dual bounds;
    1694 * otherwise, they can only be used as an estimate values */
    1695 SCIP_Bool* upvalid, /**< stores whether the returned up values are a valid dual bounds;
    1696 * otherwise, they can only be used as an estimate values */
    1697 int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
    1698 )
    1699{
    1700 assert( lpi != NULL );
    1701 assert( lpi->linear_program != NULL );
    1702 assert( cols != NULL );
    1703 assert( psols != NULL );
    1704 assert( down != NULL) ;
    1705 assert( up != NULL );
    1706 assert( downvalid != NULL );
    1707 assert( upvalid != NULL );
    1708
    1709 SCIPerrorMessage("SCIPlpiStrongbranchesFrac - not implemented.\n");
    1710
    1711 return SCIP_NOTIMPLEMENTED;
    1712}
    1713
    1714/** performs strong branching iterations on one candidate with @b integral value */
    1716 SCIP_LPI* lpi, /**< LP interface structure */
    1717 int col, /**< column to apply strong branching on */
    1718 SCIP_Real psol, /**< current integral primal solution value of column */
    1719 int itlim, /**< iteration limit for strong branchings */
    1720 SCIP_Real* down, /**< stores dual bound after branching column down */
    1721 SCIP_Real* up, /**< stores dual bound after branching column up */
    1722 SCIP_Bool* downvalid, /**< stores whether the returned down value is a valid dual bound;
    1723 * otherwise, it can only be used as an estimate value */
    1724 SCIP_Bool* upvalid, /**< stores whether the returned up value is a valid dual bound;
    1725 * otherwise, it can only be used as an estimate value */
    1726 int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
    1727 )
    1728{
    1729 assert( lpi != NULL );
    1730 assert( lpi->linear_program != NULL );
    1731 assert( down != NULL );
    1732 assert( up != NULL );
    1733 assert( downvalid != NULL );
    1734 assert( upvalid != NULL );
    1735
    1736 SCIP_CALL( strongbranch(lpi, col, psol, itlim, down, up, downvalid, upvalid, iter) );
    1737
    1738 return SCIP_OKAY;
    1739}
    1740
    1741/** performs strong branching iterations on given candidates with @b integral values */
    1743 SCIP_LPI* lpi, /**< LP interface structure */
    1744 int* cols, /**< columns to apply strong branching on */
    1745 int ncols, /**< number of columns */
    1746 SCIP_Real* psols, /**< current integral primal solution values of columns */
    1747 int itlim, /**< iteration limit for strong branchings */
    1748 SCIP_Real* down, /**< stores dual bounds after branching columns down */
    1749 SCIP_Real* up, /**< stores dual bounds after branching columns up */
    1750 SCIP_Bool* downvalid, /**< stores whether the returned down values are valid dual bounds;
    1751 * otherwise, they can only be used as an estimate values */
    1752 SCIP_Bool* upvalid, /**< stores whether the returned up values are a valid dual bounds;
    1753 * otherwise, they can only be used as an estimate values */
    1754 int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
    1755 )
    1756{
    1757 assert( lpi != NULL );
    1758 assert( lpi->linear_program != NULL );
    1759 assert( cols != NULL );
    1760 assert( psols != NULL );
    1761 assert( down != NULL) ;
    1762 assert( up != NULL );
    1763 assert( downvalid != NULL );
    1764 assert( upvalid != NULL );
    1765
    1766 SCIPerrorMessage("SCIPlpiStrongbranchesInt - not implemented.\n");
    1767
    1768 return SCIP_NOTIMPLEMENTED;
    1769}
    1770/**@} */
    1771
    1772
    1773
    1774
    1775/*
    1776 * Solution Information Methods
    1777 */
    1778
    1779/**@name Solution Information Methods */
    1780/**@{ */
    1781
    1782/** returns whether a solve method was called after the last modification of the LP */
    1784 SCIP_LPI* lpi /**< LP interface structure */
    1785 )
    1786{
    1787 assert( lpi != NULL );
    1788
    1789 /* @todo Track this to avoid uneeded resolving. */
    1790 return ( ! lpi->lp_modified_since_last_solve );
    1791}
    1792
    1793/** gets information about primal and dual feasibility of the current LP solution
    1794 *
    1795 * The feasibility information is with respect to the last solving call and it is only relevant if SCIPlpiWasSolved()
    1796 * returns true. If the LP is changed, this information might be invalidated.
    1797 *
    1798 * Note that @a primalfeasible and @a dualfeasible should only return true if the solver has proved the respective LP to
    1799 * be feasible. Thus, the return values should be equal to the values of SCIPlpiIsPrimalFeasible() and
    1800 * SCIPlpiIsDualFeasible(), respectively. Note that if feasibility cannot be proved, they should return false (even if
    1801 * the problem might actually be feasible).
    1802 */
    1804 SCIP_LPI* lpi, /**< LP interface structure */
    1805 SCIP_Bool* primalfeasible, /**< pointer to store primal feasibility status */
    1806 SCIP_Bool* dualfeasible /**< pointer to store dual feasibility status */
    1807 )
    1808{
    1809 assert( lpi != NULL );
    1810 assert( lpi->solver != NULL );
    1811 assert( primalfeasible != NULL );
    1812 assert( dualfeasible != NULL );
    1813
    1814 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1815 *primalfeasible = (status == ProblemStatus::OPTIMAL || status == ProblemStatus::PRIMAL_FEASIBLE);
    1816 *dualfeasible = (status == ProblemStatus::OPTIMAL || status == ProblemStatus::DUAL_FEASIBLE);
    1817
    1818 SCIPdebugMessage("SCIPlpiGetSolFeasibility primal:%u dual:%u\n", *primalfeasible, *dualfeasible);
    1819
    1820 return SCIP_OKAY;
    1821}
    1822
    1823/** returns TRUE iff LP is proven to have a primal unbounded ray (but not necessary a primal feasible point);
    1824 * this does not necessarily mean, that the solver knows and can return the primal ray
    1825 */
    1827 SCIP_LPI* lpi /**< LP interface structure */
    1828 )
    1829{
    1830 assert( lpi != NULL );
    1831 assert( lpi->solver != NULL );
    1832
    1833 return lpi->solver->GetProblemStatus() == ProblemStatus::PRIMAL_UNBOUNDED;
    1834}
    1835
    1836/** returns TRUE iff LP is proven to have a primal unbounded ray (but not necessary a primal feasible point),
    1837 * and the solver knows and can return the primal ray
    1838 */
    1840 SCIP_LPI* lpi /**< LP interface structure */
    1841 )
    1842{
    1843 assert( lpi != NULL );
    1844 assert( lpi->solver != NULL );
    1845
    1846 return lpi->solver->GetProblemStatus() == ProblemStatus::PRIMAL_UNBOUNDED;
    1847}
    1848
    1849/** returns TRUE iff LP is proven to be primal unbounded */
    1851 SCIP_LPI* lpi /**< LP interface structure */
    1852 )
    1853{
    1854 assert( lpi != NULL );
    1855 assert( lpi->solver != NULL );
    1856
    1857 return lpi->solver->GetProblemStatus() == ProblemStatus::PRIMAL_UNBOUNDED;
    1858}
    1859
    1860/** returns TRUE iff LP is proven to be primal infeasible */
    1862 SCIP_LPI* lpi /**< LP interface structure */
    1863 )
    1864{
    1865 assert( lpi != NULL );
    1866 assert( lpi->solver != NULL );
    1867
    1868 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1869
    1870 return status == ProblemStatus::DUAL_UNBOUNDED || status == ProblemStatus::PRIMAL_INFEASIBLE;
    1871}
    1872
    1873/** returns TRUE iff LP is proven to be primal feasible */
    1875 SCIP_LPI* lpi /**< LP interface structure */
    1876 )
    1877{
    1878 assert( lpi != NULL );
    1879 assert( lpi->solver != NULL );
    1880
    1881 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1882
    1883 return status == ProblemStatus::PRIMAL_FEASIBLE || status == ProblemStatus::OPTIMAL;
    1884}
    1885
    1886/** returns TRUE iff LP is proven to have a dual unbounded ray (but not necessary a dual feasible point);
    1887 * this does not necessarily mean, that the solver knows and can return the dual ray
    1888 */
    1890 SCIP_LPI* lpi /**< LP interface structure */
    1891 )
    1892{
    1893 assert( lpi != NULL );
    1894 assert( lpi->solver != NULL );
    1895
    1896 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1897
    1898 return status == ProblemStatus::DUAL_UNBOUNDED;
    1899}
    1900
    1901/** returns TRUE iff LP is proven to have a dual unbounded ray (but not necessary a dual feasible point),
    1902 * and the solver knows and can return the dual ray
    1903 */
    1905 SCIP_LPI* lpi /**< LP interface structure */
    1906 )
    1907{
    1908 assert( lpi != NULL );
    1909 assert( lpi->solver != NULL );
    1910
    1911 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1912
    1913 return status == ProblemStatus::DUAL_UNBOUNDED;
    1914}
    1915
    1916/** returns TRUE iff LP is proven to be dual unbounded */
    1918 SCIP_LPI* lpi /**< LP interface structure */
    1919 )
    1920{
    1921 assert( lpi != NULL );
    1922 assert( lpi->solver != NULL );
    1923
    1924 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1925 return status == ProblemStatus::DUAL_UNBOUNDED;
    1926}
    1927
    1928/** returns TRUE iff LP is proven to be dual infeasible */
    1930 SCIP_LPI* lpi /**< LP interface structure */
    1931 )
    1932{
    1933 assert( lpi != NULL );
    1934 assert( lpi->solver != NULL );
    1935
    1936 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1937 return status == ProblemStatus::PRIMAL_UNBOUNDED || status == ProblemStatus::DUAL_INFEASIBLE;
    1938}
    1939
    1940/** returns TRUE iff LP is proven to be dual feasible */
    1942 SCIP_LPI* lpi /**< LP interface structure */
    1943 )
    1944{
    1945 assert( lpi != NULL );
    1946 assert( lpi->solver != NULL );
    1947
    1948 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1949
    1950 return status == ProblemStatus::DUAL_FEASIBLE || status == ProblemStatus::OPTIMAL;
    1951}
    1952
    1953/** returns TRUE iff LP was solved to optimality */
    1955 SCIP_LPI* lpi /**< LP interface structure */
    1956 )
    1957{
    1958 assert( lpi != NULL );
    1959 assert( lpi->solver != NULL );
    1960
    1961 return lpi->solver->GetProblemStatus() == ProblemStatus::OPTIMAL;
    1962}
    1963
    1964/** returns TRUE iff current LP solution is stable
    1965 *
    1966 * This function should return true if the solution is reliable, i.e., feasible and optimal (or proven
    1967 * infeasible/unbounded) with respect to the original problem. The optimality status might be with respect to a scaled
    1968 * version of the problem, but the solution might not be feasible to the unscaled original problem; in this case,
    1969 * SCIPlpiIsStable() should return false.
    1970 */
    1972 SCIP_LPI* lpi /**< LP interface structure */
    1973 )
    1974{
    1975 assert( lpi != NULL );
    1976 assert( lpi->solver != NULL );
    1977
    1978 /* For correctness, we need to report "unstable" if Glop was not able to prove optimality because of numerical
    1979 * issues. Currently Glop still reports primal/dual feasible if at the end, one status is within the tolerance but not
    1980 * the other. */
    1981 const ProblemStatus status = lpi->solver->GetProblemStatus();
    1982 if ( (status == ProblemStatus::PRIMAL_FEASIBLE || status == ProblemStatus::DUAL_FEASIBLE) &&
    1984 {
    1985 SCIPdebugMessage("OPTIMAL not reached and no limit: unstable.\n");
    1986 return FALSE;
    1987 }
    1988
    1989 if ( status == ProblemStatus::ABNORMAL || status == ProblemStatus::INVALID_PROBLEM || status == ProblemStatus::IMPRECISE )
    1990 return FALSE;
    1991
    1992 /* If we have a regular basis and the condition limit is set, we compute (an upper bound on) the condition number of
    1993 * the basis; everything above the specified threshold is then counted as instable. */
    1994 if ( lpi->checkcondition && (SCIPlpiIsOptimal(lpi) || SCIPlpiIsObjlimExc(lpi)) )
    1995 {
    1996 SCIP_RETCODE retcode;
    1997 SCIP_Real kappa;
    1998
    2000 if ( retcode != SCIP_OKAY )
    2001 SCIPABORT();
    2002 assert( kappa != SCIP_INVALID ); /*lint !e777*/
    2003
    2004 if ( kappa > lpi->conditionlimit )
    2005 return FALSE;
    2006 }
    2007
    2008 return TRUE;
    2009}
    2010
    2011/** returns TRUE iff the objective limit was reached */
    2013 SCIP_LPI* lpi /**< LP interface structure */
    2014 )
    2015{
    2016 assert( lpi != NULL );
    2017 assert( lpi->solver != NULL );
    2018
    2019 return lpi->solver->objective_limit_reached();
    2020}
    2021
    2022/** returns TRUE iff the iteration limit was reached */
    2024 SCIP_LPI* lpi /**< LP interface structure */
    2025 )
    2026{
    2027 assert( lpi != NULL );
    2028 assert( lpi->solver != NULL );
    2029 assert( lpi->niterations >= (int) lpi->solver->GetNumberOfIterations() );
    2030
    2031 int maxiter = (int) lpi->parameters->max_number_of_iterations();
    2032 return maxiter >= 0 && lpi->niterations >= maxiter;
    2033}
    2034
    2035/** returns TRUE iff the time limit was reached */
    2037 SCIP_LPI* lpi /**< LP interface structure */
    2038 )
    2039{
    2040 assert( lpi != NULL );
    2041 assert( lpi->solver != NULL );
    2042
    2043 return lpi->lp_time_limit_was_reached;
    2044}
    2045
    2046/** returns the internal solution status of the solver */
    2048 SCIP_LPI* lpi /**< LP interface structure */
    2049 )
    2050{
    2051 assert( lpi != NULL );
    2052 assert( lpi->solver != NULL );
    2053
    2054 return static_cast<int>(lpi->solver->GetProblemStatus());
    2055}
    2056
    2057/** tries to reset the internal status of the LP solver in order to ignore an instability of the last solving call */
    2059 SCIP_LPI* lpi, /**< LP interface structure */
    2060 SCIP_Bool* success /**< pointer to store, whether the instability could be ignored */
    2061 )
    2062{
    2063 assert( lpi != NULL );
    2064 assert( lpi->solver != NULL );
    2065 assert( success != NULL );
    2066
    2067 *success = FALSE;
    2068
    2069 return SCIP_OKAY;
    2070}
    2071
    2072/** gets objective value of solution */
    2074 SCIP_LPI* lpi, /**< LP interface structure */
    2075 SCIP_Real* objval /**< stores the objective value */
    2076 )
    2077{
    2078 assert( lpi != NULL );
    2079 assert( lpi->solver != NULL );
    2080 assert( objval != NULL );
    2081
    2082 *objval = lpi->solver->GetObjectiveValue();
    2083
    2084 return SCIP_OKAY;
    2085}
    2086
    2087/** gets primal and dual solution vectors for feasible LPs
    2088 *
    2089 * Before calling this function, the caller must ensure that the LP has been solved to optimality, i.e., that
    2090 * SCIPlpiIsOptimal() returns true.
    2091 */
    2093 SCIP_LPI* lpi, /**< LP interface structure */
    2094 SCIP_Real* objval, /**< stores the objective value, may be NULL if not needed */
    2095 SCIP_Real* primsol, /**< primal solution vector, may be NULL if not needed */
    2096 SCIP_Real* dualsol, /**< dual solution vector, may be NULL if not needed */
    2097 SCIP_Real* activity, /**< row activity vector, may be NULL if not needed */
    2098 SCIP_Real* redcost /**< reduced cost vector, may be NULL if not needed */
    2099 )
    2100{
    2101 assert( lpi != NULL );
    2102 assert( lpi->solver != NULL );
    2103
    2104 SCIPdebugMessage("SCIPlpiGetSol\n");
    2105 if ( objval != NULL )
    2106 *objval = lpi->solver->GetObjectiveValue();
    2107
    2108 const ColIndex num_cols = lpi->linear_program->num_variables();
    2109 for (ColIndex col(0); col < num_cols; ++col)
    2110 {
    2111 int i = col.value();
    2112
    2113 if ( primsol != NULL )
    2114 primsol[i] = lpi->scaler->UnscaleVariableValue(col, lpi->solver->GetVariableValue(col));
    2115
    2116 if ( redcost != NULL )
    2117 redcost[i] = lpi->scaler->UnscaleReducedCost(col, lpi->solver->GetReducedCost(col));
    2118 }
    2119
    2120 const RowIndex num_rows = lpi->linear_program->num_constraints();
    2121 for (RowIndex row(0); row < num_rows; ++row)
    2122 {
    2123 int j = row.value();
    2124
    2125 if ( dualsol != NULL )
    2126 dualsol[j] = lpi->scaler->UnscaleDualValue(row, lpi->solver->GetDualValue(row));
    2127
    2128 if ( activity != NULL )
    2129 activity[j] = lpi->scaler->UnscaleConstraintActivity(row, lpi->solver->GetConstraintActivity(row));
    2130 }
    2131
    2132 return SCIP_OKAY;
    2133}
    2134
    2135/** gets primal ray for unbounded LPs */
    2137 SCIP_LPI* lpi, /**< LP interface structure */
    2138 SCIP_Real* ray /**< primal ray */
    2139 )
    2140{
    2141 assert( lpi != NULL );
    2142 assert( lpi->solver != NULL );
    2143 assert( ray != NULL );
    2144
    2145 SCIPdebugMessage("SCIPlpiGetPrimalRay\n");
    2146
    2147 const ColIndex num_cols = lpi->linear_program->num_variables();
    2148 const DenseRow& primal_ray = lpi->solver->GetPrimalRay();
    2149 for (ColIndex col(0); col < num_cols; ++col)
    2150 ray[col.value()] = lpi->scaler->UnscaleVariableValue(col, primal_ray[col]);
    2151
    2152 return SCIP_OKAY;
    2153}
    2154
    2155/** gets dual Farkas proof for infeasibility */
    2157 SCIP_LPI* lpi, /**< LP interface structure */
    2158 SCIP_Real* dualfarkas /**< dual Farkas row multipliers */
    2159 )
    2160{
    2161 assert( lpi != NULL );
    2162 assert( lpi->solver != NULL );
    2163 assert( dualfarkas != NULL );
    2164
    2165 SCIPdebugMessage("SCIPlpiGetDualfarkas\n");
    2166
    2167 const RowIndex num_rows = lpi->linear_program->num_constraints();
    2168 const DenseColumn& dual_ray = lpi->solver->GetDualRay();
    2169 for (RowIndex row(0); row < num_rows; ++row)
    2170 dualfarkas[row.value()] = -lpi->scaler->UnscaleDualValue(row, dual_ray[row]); /* reverse sign */
    2171
    2172 return SCIP_OKAY;
    2173}
    2174
    2175/** gets the number of LP iterations of the last solve call */
    2177 SCIP_LPI* lpi, /**< LP interface structure */
    2178 int* iterations /**< pointer to store the number of iterations of the last solve call */
    2179 )
    2180{
    2181 assert( lpi != NULL );
    2182 assert( lpi->solver != NULL );
    2183 assert( iterations != NULL );
    2184
    2185 *iterations = (int) lpi->niterations;
    2186
    2187 return SCIP_OKAY;
    2188}
    2189
    2190/** gets information about the quality of an LP solution
    2191 *
    2192 * Such information is usually only available, if also a (maybe not optimal) solution is available.
    2193 * The LPI should return SCIP_INVALID for @p quality, if the requested quantity is not available.
    2194 */
    2196 SCIP_LPI* lpi, /**< LP interface structure */
    2197 SCIP_LPSOLQUALITY qualityindicator, /**< indicates which quality should be returned */
    2198 SCIP_Real* quality /**< pointer to store quality number */
    2199 )
    2200{
    2201 assert( lpi != NULL );
    2202 assert( lpi->solver != NULL );
    2203 assert( quality != NULL );
    2204
    2205 SCIPdebugMessage("Requesting solution quality: quality %d\n", qualityindicator);
    2206
    2207 switch ( qualityindicator )
    2208 {
    2210 *quality = lpi->solver->GetBasisFactorization().ComputeInfinityNormConditionNumber();
    2211 break;
    2212
    2214 *quality = lpi->solver->GetBasisFactorization().ComputeInfinityNormConditionNumberUpperBound();
    2215 break;
    2216
    2217 default:
    2218 SCIPerrorMessage("Solution quality %d unknown.\n", qualityindicator);
    2219 return SCIP_INVALIDDATA;
    2220 }
    2221
    2222 return SCIP_OKAY;
    2223}
    2224
    2225/**@} */
    2226
    2227
    2228
    2229
    2230/*
    2231 * LP Basis Methods
    2232 */
    2233
    2234/**@name LP Basis Methods */
    2235/**@{ */
    2236
    2237/** convert Glop variable basis status to SCIP status */
    2238static
    2240 VariableStatus status, /**< variable status */
    2241 Fractional rc /**< reduced cost of variable */
    2242 )
    2243{
    2244 switch ( status )
    2245 {
    2246 case VariableStatus::BASIC:
    2247 return SCIP_BASESTAT_BASIC;
    2248 case VariableStatus::AT_UPPER_BOUND:
    2249 return SCIP_BASESTAT_UPPER;
    2250 case VariableStatus::AT_LOWER_BOUND:
    2251 return SCIP_BASESTAT_LOWER;
    2252 case VariableStatus::FREE:
    2253 return SCIP_BASESTAT_ZERO;
    2254 case VariableStatus::FIXED_VALUE:
    2255 return rc > 0.0 ? SCIP_BASESTAT_LOWER : SCIP_BASESTAT_UPPER;
    2256 default:
    2257 SCIPerrorMessage("invalid Glop basis status.\n");
    2258 abort();
    2259 }
    2260}
    2261
    2262/** convert Glop constraint basis status to SCIP status */
    2263static
    2265 ConstraintStatus status, /**< constraint status */
    2266 Fractional dual /**< dual variable value */
    2267 )
    2268{
    2269 switch ( status )
    2270 {
    2271 case ConstraintStatus::BASIC:
    2272 return SCIP_BASESTAT_BASIC;
    2273 case ConstraintStatus::AT_UPPER_BOUND:
    2274 return SCIP_BASESTAT_UPPER;
    2275 case ConstraintStatus::AT_LOWER_BOUND:
    2276 return SCIP_BASESTAT_LOWER;
    2277 case ConstraintStatus::FREE:
    2278 return SCIP_BASESTAT_ZERO;
    2279 case ConstraintStatus::FIXED_VALUE:
    2280 return dual > 0.0 ? SCIP_BASESTAT_LOWER : SCIP_BASESTAT_UPPER;
    2281 default:
    2282 SCIPerrorMessage("invalid Glop basis status.\n");
    2283 abort();
    2284 }
    2285}
    2286
    2287/** Convert SCIP variable status to Glop status */
    2288static
    2290 int status /**< SCIP variable status */
    2291 )
    2292{
    2293 switch ( status )
    2294 {
    2296 return VariableStatus::BASIC;
    2298 return VariableStatus::AT_UPPER_BOUND;
    2300 return VariableStatus::AT_LOWER_BOUND;
    2301 case SCIP_BASESTAT_ZERO:
    2302 return VariableStatus::FREE;
    2303 default:
    2304 SCIPerrorMessage("invalid SCIP basis status.\n");
    2305 abort();
    2306 }
    2307}
    2308
    2309/** Convert a SCIP constraint status to its corresponding Glop slack VariableStatus.
    2310 *
    2311 * Note that we swap the upper/lower bounds.
    2312 */
    2313static
    2315 int status /**< SCIP constraint status */
    2316 )
    2317{
    2318 switch ( status )
    2319 {
    2321 return VariableStatus::BASIC;
    2323 return VariableStatus::AT_LOWER_BOUND;
    2325 return VariableStatus::AT_UPPER_BOUND;
    2326 case SCIP_BASESTAT_ZERO:
    2327 return VariableStatus::FREE;
    2328 default:
    2329 SCIPerrorMessage("invalid SCIP basis status.\n");
    2330 abort();
    2331 }
    2332}
    2333
    2334/** gets current basis status for columns and rows; arrays must be large enough to store the basis status */
    2336 SCIP_LPI* lpi, /**< LP interface structure */
    2337 int* cstat, /**< array to store column basis status, or NULL */
    2338 int* rstat /**< array to store row basis status, or NULL */
    2339 )
    2340{
    2341 SCIPdebugMessage("SCIPlpiGetBase\n");
    2342
    2343 assert( lpi->solver->GetProblemStatus() == ProblemStatus::OPTIMAL );
    2344
    2345 if ( cstat != NULL )
    2346 {
    2347 const ColIndex num_cols = lpi->linear_program->num_variables();
    2348 for (ColIndex col(0); col < num_cols; ++col)
    2349 {
    2350 int i = col.value();
    2351 cstat[i] = (int) ConvertGlopVariableStatus(lpi->solver->GetVariableStatus(col), lpi->solver->GetReducedCost(col));
    2352 }
    2353 }
    2354
    2355 if ( rstat != NULL )
    2356 {
    2357 const RowIndex num_rows = lpi->linear_program->num_constraints();
    2358 for (RowIndex row(0); row < num_rows; ++row)
    2359 {
    2360 int i = row.value();
    2361 rstat[i] = (int) ConvertGlopConstraintStatus(lpi->solver->GetConstraintStatus(row), lpi->solver->GetDualValue(row));
    2362 }
    2363 }
    2364
    2365 return SCIP_OKAY;
    2366}
    2367
    2368/** sets current basis status for columns and rows */
    2370 SCIP_LPI* lpi, /**< LP interface structure */
    2371 const int* cstat, /**< array with column basis status */
    2372 const int* rstat /**< array with row basis status */
    2373 )
    2374{
    2375 assert( lpi != NULL );
    2376 assert( lpi->linear_program != NULL );
    2377 assert( lpi->solver != NULL );
    2378
    2379 const ColIndex num_cols = lpi->linear_program->num_variables();
    2380 const RowIndex num_rows = lpi->linear_program->num_constraints();
    2381
    2382 assert( cstat != NULL || num_cols == 0 );
    2383 assert( rstat != NULL || num_rows == 0 );
    2384
    2385 SCIPdebugMessage("SCIPlpiSetBase\n");
    2386
    2387 BasisState state;
    2388 state.statuses.reserve(ColIndex(num_cols.value() + num_rows.value()));
    2389
    2390 for (ColIndex col(0); col < num_cols; ++col)
    2391 state.statuses[col] = ConvertSCIPVariableStatus(cstat[col.value()]);
    2392
    2393 for (RowIndex row(0); row < num_rows; ++row)
    2394 state.statuses[num_cols + RowToColIndex(row)] = ConvertSCIPConstraintStatusToSlackStatus(cstat[row.value()]);
    2395
    2396 lpi->solver->LoadStateForNextSolve(state);
    2397
    2398 return SCIP_OKAY;
    2399}
    2400
    2401/** returns the indices of the basic columns and rows; basic column n gives value n, basic row m gives value -1-m */
    2403 SCIP_LPI* lpi, /**< LP interface structure */
    2404 int* bind /**< pointer to store basis indices ready to keep number of rows entries */
    2405 )
    2406{
    2407 assert( lpi != NULL );
    2408 assert( lpi->linear_program != NULL );
    2409 assert( lpi->solver != NULL );
    2410 assert( bind != NULL );
    2411
    2412 SCIPdebugMessage("SCIPlpiGetBasisInd\n");
    2413
    2414 /* the order is important! */
    2415 const ColIndex num_cols = lpi->linear_program->num_variables();
    2416 const RowIndex num_rows = lpi->linear_program->num_constraints();
    2417 for (RowIndex row(0); row < num_rows; ++row)
    2418 {
    2419 const ColIndex col = lpi->solver->GetBasis(row);
    2420 if (col < num_cols)
    2421 bind[row.value()] = col.value();
    2422 else
    2423 {
    2424 assert( col < num_cols.value() + num_rows.value() );
    2425 bind[row.value()] = -1 - (col - num_cols).value();
    2426 }
    2427 }
    2428
    2429 return SCIP_OKAY;
    2430}
    2431
    2432/** get row of inverse basis matrix B^-1
    2433 *
    2434 * @note The LP interface defines slack variables to have coefficient +1. This means that if, internally, the LP solver
    2435 * uses a -1 coefficient, then rows associated with slacks variables whose coefficient is -1, should be negated;
    2436 * see also the explanation in lpi.h.
    2437 */
    2439 SCIP_LPI* lpi, /**< LP interface structure */
    2440 int r, /**< row number */
    2441 SCIP_Real* coef, /**< pointer to store the coefficients of the row */
    2442 int* inds, /**< array to store the non-zero indices, or NULL */
    2443 int* ninds /**< pointer to store the number of non-zero indices, or NULL
    2444 * (-1: if we do not store sparsity information) */
    2445 )
    2446{
    2447 assert( lpi != NULL );
    2448 assert( lpi->linear_program != NULL );
    2449 assert( lpi->solver != NULL );
    2450 assert( lpi->scaler != NULL );
    2451 assert( coef != NULL );
    2452
    2453 lpi->solver->GetBasisFactorization().LeftSolveForUnitRow(ColIndex(r), lpi->tmp_row);
    2454 lpi->scaler->UnscaleUnitRowLeftSolve(lpi->solver->GetBasis(RowIndex(r)), lpi->tmp_row);
    2455
    2456 const ColIndex size = lpi->tmp_row->values.size();
    2457 assert( size.value() == lpi->linear_program->num_constraints() );
    2458
    2459 /* if we want a sparse vector */
    2460 if ( ninds != NULL && inds != NULL )
    2461 {
    2462 *ninds = 0;
    2463 /* Vectors in Glop might be stored in dense or sparse format depending on the values. If non_zeros are given, we
    2464 * can directly loop over the non_zeros, otherwise we have to collect the nonzeros. */
    2465 if ( ! lpi->tmp_row->non_zeros.empty() )
    2466 {
    2467 ScatteredRowIterator end = lpi->tmp_row->end();
    2468 for (ScatteredRowIterator iter = lpi->tmp_row->begin(); iter != end; ++iter)
    2469 {
    2470 int idx = (*iter).column().value();
    2471 assert( 0 <= idx && idx < lpi->linear_program->num_constraints() );
    2472 coef[idx] = (*iter).coefficient();
    2473 inds[(*ninds)++] = idx;
    2474 }
    2475 }
    2476 else
    2477 {
    2478 /* use dense access to tmp_row */
    2479 const Fractional eps = lpi->parameters->primal_feasibility_tolerance();
    2480 for (ColIndex col(0); col < size; ++col)
    2481 {
    2482 SCIP_Real val = (*lpi->tmp_row)[col];
    2483 if ( fabs(val) >= eps )
    2484 {
    2485 coef[col.value()] = val;
    2486 inds[(*ninds)++] = col.value();
    2487 }
    2488 }
    2489 }
    2490 return SCIP_OKAY;
    2491 }
    2492
    2493 /* dense version */
    2494 for (ColIndex col(0); col < size; ++col)
    2495 coef[col.value()] = (*lpi->tmp_row)[col];
    2496
    2497 if ( ninds != NULL )
    2498 *ninds = -1;
    2499
    2500 return SCIP_OKAY;
    2501}
    2502
    2503/** get column of inverse basis matrix B^-1
    2504 *
    2505 * @note The LP interface defines slack variables to have coefficient +1. This means that if, internally, the LP solver
    2506 * uses a -1 coefficient, then rows associated with slacks variables whose coefficient is -1, should be negated;
    2507 * see also the explanation in lpi.h.
    2508 */
    2510 SCIP_LPI* lpi, /**< LP interface structure */
    2511 int c, /**< column number of B^-1; this is NOT the number of the column in the LP;
    2512 * you have to call SCIPlpiGetBasisInd() to get the array which links the
    2513 * B^-1 column numbers to the row and column numbers of the LP!
    2514 * c must be between 0 and nrows-1, since the basis has the size
    2515 * nrows * nrows */
    2516 SCIP_Real* coef, /**< pointer to store the coefficients of the column */
    2517 int* inds, /**< array to store the non-zero indices, or NULL */
    2518 int* ninds /**< pointer to store the number of non-zero indices, or NULL
    2519 * (-1: if we do not store sparsity information) */
    2520 )
    2521{
    2522 assert( lpi != NULL );
    2523 assert( lpi->linear_program != NULL );
    2524 assert( lpi->solver != NULL );
    2525 assert( lpi->scaler != NULL );
    2526 assert( coef != NULL );
    2527
    2528 /* we need to loop through the rows to extract the values for column c */
    2529 const ColIndex col(c);
    2530 const RowIndex num_rows = lpi->linear_program->num_constraints();
    2531
    2532 /* if we want a sparse vector */
    2533 if ( ninds != NULL && inds != NULL )
    2534 {
    2535 const Fractional eps = lpi->parameters->primal_feasibility_tolerance();
    2536
    2537 *ninds = 0;
    2538 for (int row = 0; row < num_rows; ++row)
    2539 {
    2540 lpi->solver->GetBasisFactorization().LeftSolveForUnitRow(ColIndex(row), lpi->tmp_row);
    2541 lpi->scaler->UnscaleUnitRowLeftSolve(lpi->solver->GetBasis(RowIndex(row)), lpi->tmp_row);
    2542
    2543 SCIP_Real val = (*lpi->tmp_row)[col];
    2544 if ( fabs(val) >= eps )
    2545 {
    2546 coef[row] = val;
    2547 inds[(*ninds)++] = row;
    2548 }
    2549 }
    2550 return SCIP_OKAY;
    2551 }
    2552
    2553 /* dense version */
    2554 for (int row = 0; row < num_rows; ++row)
    2555 {
    2556 lpi->solver->GetBasisFactorization().LeftSolveForUnitRow(ColIndex(row), lpi->tmp_row);
    2557 lpi->scaler->UnscaleUnitRowLeftSolve(lpi->solver->GetBasis(RowIndex(row)), lpi->tmp_row);
    2558 coef[row] = (*lpi->tmp_row)[col];
    2559 }
    2560
    2561 if ( ninds != NULL )
    2562 *ninds = -1;
    2563
    2564 return SCIP_OKAY;
    2565}
    2566
    2567/** get row of inverse basis matrix times constraint matrix B^-1 * A
    2568 *
    2569 * @note The LP interface defines slack variables to have coefficient +1. This means that if, internally, the LP solver
    2570 * uses a -1 coefficient, then rows associated with slacks variables whose coefficient is -1, should be negated;
    2571 * see also the explanation in lpi.h.
    2572 */
    2574 SCIP_LPI* lpi, /**< LP interface structure */
    2575 int r, /**< row number */
    2576 const SCIP_Real* binvrow, /**< row in (A_B)^-1 from prior call to SCIPlpiGetBInvRow(), or NULL */
    2577 SCIP_Real* coef, /**< vector to return coefficients of the row */
    2578 int* inds, /**< array to store the non-zero indices, or NULL */
    2579 int* ninds /**< pointer to store the number of non-zero indices, or NULL
    2580 * (-1: if we do not store sparsity information) */
    2581 )
    2582{
    2583 assert( lpi != NULL );
    2584 assert( lpi->linear_program != NULL );
    2585 assert( lpi->solver != NULL );
    2586 assert( lpi->scaler != NULL );
    2587 assert( coef != NULL );
    2588
    2589 /* get row of basis inverse, loop through columns and muliply with matrix */
    2590 lpi->solver->GetBasisFactorization().LeftSolveForUnitRow(ColIndex(r), lpi->tmp_row);
    2591 lpi->scaler->UnscaleUnitRowLeftSolve(lpi->solver->GetBasis(RowIndex(r)), lpi->tmp_row);
    2592
    2593 const ColIndex num_cols = lpi->linear_program->num_variables();
    2594
    2595 /* if we want a sparse vector */
    2596 if ( ninds != NULL && inds != NULL )
    2597 {
    2598 const Fractional eps = lpi->parameters->primal_feasibility_tolerance();
    2599
    2600 *ninds = 0;
    2601 for (ColIndex col(0); col < num_cols; ++col)
    2602 {
    2603 SCIP_Real val = operations_research::glop::ScalarProduct(lpi->tmp_row->values, lpi->linear_program->GetSparseColumn(col));
    2604 if ( fabs(val) >= eps )
    2605 {
    2606 coef[col.value()] = val;
    2607 inds[(*ninds)++] = col.value();
    2608 }
    2609 }
    2610 return SCIP_OKAY;
    2611 }
    2612
    2613 /* dense version */
    2614 for (ColIndex col(0); col < num_cols; ++col)
    2615 coef[col.value()] = operations_research::glop::ScalarProduct(lpi->tmp_row->values, lpi->linear_program->GetSparseColumn(col));
    2616
    2617 if ( ninds != NULL )
    2618 *ninds = -1;
    2619
    2620 return SCIP_OKAY;
    2621}
    2622
    2623/** get column of inverse basis matrix times constraint matrix B^-1 * A
    2624 *
    2625 * @note The LP interface defines slack variables to have coefficient +1. This means that if, internally, the LP solver
    2626 * uses a -1 coefficient, then rows associated with slacks variables whose coefficient is -1, should be negated;
    2627 * see also the explanation in lpi.h.
    2628 */
    2630 SCIP_LPI* lpi, /**< LP interface structure */
    2631 int c, /**< column number */
    2632 SCIP_Real* coef, /**< vector to return coefficients of the column */
    2633 int* inds, /**< array to store the non-zero indices, or NULL */
    2634 int* ninds /**< pointer to store the number of non-zero indices, or NULL
    2635 * (-1: if we do not store sparsity information) */
    2636 )
    2637{
    2638 assert( lpi != NULL );
    2639 assert( lpi->linear_program != NULL );
    2640 assert( lpi->solver != NULL );
    2641 assert( lpi->scaler != NULL );
    2642 assert( coef != NULL );
    2643
    2644 lpi->solver->GetBasisFactorization().RightSolveForProblemColumn(ColIndex(c), lpi->tmp_column);
    2645 lpi->scaler->UnscaleColumnRightSolve(lpi->solver->GetBasisVector(), ColIndex(c), lpi->tmp_column);
    2646
    2647 const RowIndex num_rows = lpi->tmp_column->values.size();
    2648
    2649 /* if we want a sparse vector */
    2650 if ( ninds != NULL && inds != NULL )
    2651 {
    2652 *ninds = 0;
    2653 /* Vectors in Glop might be stored in dense or sparse format depending on the values. If non_zeros are given, we
    2654 * can directly loop over the non_zeros, otherwise we have to collect the nonzeros. */
    2655 if ( ! lpi->tmp_column->non_zeros.empty() )
    2656 {
    2657 ScatteredColumnIterator end = lpi->tmp_column->end();
    2658 for (ScatteredColumnIterator iter = lpi->tmp_column->begin(); iter != end; ++iter)
    2659 {
    2660 int idx = (*iter).row().value();
    2661 assert( 0 <= idx && idx < num_rows );
    2662 coef[idx] = (*iter).coefficient();
    2663 inds[(*ninds)++] = idx;
    2664 }
    2665 }
    2666 else
    2667 {
    2668 /* use dense access to tmp_column */
    2669 const Fractional eps = lpi->parameters->primal_feasibility_tolerance();
    2670 for (RowIndex row(0); row < num_rows; ++row)
    2671 {
    2672 SCIP_Real val = (*lpi->tmp_column)[row];
    2673 if ( fabs(val) > eps )
    2674 {
    2675 coef[row.value()] = val;
    2676 inds[(*ninds)++] = row.value();
    2677 }
    2678 }
    2679 }
    2680 return SCIP_OKAY;
    2681 }
    2682
    2683 /* dense version */
    2684 for (RowIndex row(0); row < num_rows; ++row)
    2685 coef[row.value()] = (*lpi->tmp_column)[row];
    2686
    2687 if ( ninds != NULL )
    2688 *ninds = -1;
    2689
    2690 return SCIP_OKAY;
    2691}
    2692
    2693/**@} */
    2694
    2695
    2696
    2697
    2698/*
    2699 * LP State Methods
    2700 */
    2701
    2702/**@name LP State Methods */
    2703/**@{ */
    2704
    2705/* SCIP_LPiState stores basis information and is implemented by the glop BasisState class.
    2706 * However, because in type_lpi.h there is
    2707 * typedef struct SCIP_LPiState SCIP_LPISTATE;
    2708 * We cannot just use a typedef here and SCIP_LPiState needs to be a struct.
    2709 */
    2710struct SCIP_LPiState : BasisState {};
    2711
    2712/** stores LPi state (like basis information) into lpistate object */
    2714 SCIP_LPI* lpi, /**< LP interface structure */
    2715 BMS_BLKMEM* blkmem, /**< block memory */
    2716 SCIP_LPISTATE** lpistate /**< pointer to LPi state information (like basis information) */
    2717 )
    2718{
    2719 assert( lpi != NULL );
    2720 assert( lpi->solver != NULL );
    2721 assert( lpistate != NULL );
    2722
    2723 *lpistate = static_cast<SCIP_LPISTATE*>(new BasisState(lpi->solver->GetState()));
    2724
    2725 return SCIP_OKAY;
    2726}
    2727
    2728/** loads LPi state (like basis information) into solver; note that the LP might have been extended with additional
    2729 * columns and rows since the state was stored with SCIPlpiGetState()
    2730 */
    2732 SCIP_LPI* lpi, /**< LP interface structure */
    2733 BMS_BLKMEM* blkmem, /**< block memory */
    2734 const SCIP_LPISTATE* lpistate /**< LPi state information (like basis information), or NULL */
    2735 )
    2736{
    2737 assert( lpi != NULL );
    2738 assert( lpi->solver != NULL );
    2739 assert( lpistate != NULL );
    2740
    2741 lpi->solver->LoadStateForNextSolve(*lpistate);
    2742
    2743 return SCIP_OKAY;
    2744}
    2745
    2746/** clears current LPi state (like basis information) of the solver */
    2748 SCIP_LPI* lpi /**< LP interface structure */
    2749 )
    2750{
    2751 assert( lpi != NULL );
    2752 assert( lpi->solver != NULL );
    2753
    2754 lpi->solver->ClearStateForNextSolve();
    2755
    2756 return SCIP_OKAY;
    2757}
    2758
    2759/** frees LPi state information */
    2761 SCIP_LPI* lpi, /**< LP interface structure */
    2762 BMS_BLKMEM* blkmem, /**< block memory */
    2763 SCIP_LPISTATE** lpistate /**< pointer to LPi state information (like basis information) */
    2764 )
    2765{
    2766 assert( lpi != NULL );
    2767 assert( lpi->solver != NULL );
    2768 assert( lpistate != NULL );
    2769
    2770 delete *lpistate;
    2771 *lpistate = NULL;
    2772
    2773 return SCIP_OKAY;
    2774}
    2775
    2776/** checks whether the given LP state contains simplex basis information */
    2778 SCIP_LPI* lpi, /**< LP interface structure */
    2779 SCIP_LPISTATE* lpistate /**< LP state information (like basis information), or NULL */
    2780 )
    2781{
    2782 assert( lpi != NULL );
    2783 assert( lpi->solver != NULL );
    2784
    2785 return lpistate != NULL;
    2786}
    2787
    2788/** reads LP state (like basis information from a file */
    2790 SCIP_LPI* lpi, /**< LP interface structure */
    2791 const char* fname /**< file name */
    2792 )
    2793{
    2794 assert( lpi != NULL );
    2795 assert( lpi->solver != NULL );
    2796
    2797 SCIPerrorMessage("SCIPlpiReadState - not implemented.\n");
    2798
    2799 return SCIP_NOTIMPLEMENTED;
    2800}
    2801
    2802/** writes LPi state (i.e. basis information) to a file */
    2804 SCIP_LPI* lpi, /**< LP interface structure */
    2805 const char* fname /**< file name */
    2806 )
    2807{
    2808 assert( lpi != NULL );
    2809 assert( lpi->solver != NULL );
    2810
    2811 SCIPerrorMessage("SCIPlpiWriteState - not implemented.\n");
    2812
    2813 return SCIP_NOTIMPLEMENTED;
    2814}
    2815
    2816/**@} */
    2817
    2818
    2819
    2820
    2821/*
    2822 * LP Pricing Norms Methods
    2823 */
    2824
    2825/**@name LP Pricing Norms Methods */
    2826/**@{ */
    2827
    2828/* SCIP_LPiNorms stores norm information so they are not recomputed from one state to the next. */
    2829/* @todo Implement this. */
    2830struct SCIP_LPiNorms {};
    2831
    2832/** stores LPi pricing norms information
    2833 *
    2834 * @todo store primal norms as well?
    2835 */
    2837 SCIP_LPI* lpi, /**< LP interface structure */
    2838 BMS_BLKMEM* blkmem, /**< block memory */
    2839 SCIP_LPINORMS** lpinorms /**< pointer to LPi pricing norms information */
    2840 )
    2841{
    2842 assert( lpi != NULL );
    2843 assert( blkmem != NULL );
    2844 assert( lpi->solver != NULL );
    2845 assert( lpinorms != NULL );
    2846
    2847 return SCIP_OKAY;
    2848}
    2849
    2850/** loads LPi pricing norms into solver; note that the LP might have been extended with additional
    2851 * columns and rows since the state was stored with SCIPlpiGetNorms()
    2852 */
    2854 SCIP_LPI* lpi, /**< LP interface structure */
    2855 BMS_BLKMEM* blkmem, /**< block memory */
    2856 const SCIP_LPINORMS* lpinorms /**< LPi pricing norms information, or NULL */
    2857 )
    2858{
    2859 assert( lpi != NULL );
    2860 assert( blkmem != NULL );
    2861 assert( lpi->solver != NULL );
    2862 assert( lpinorms != NULL );
    2863
    2864 return SCIP_OKAY;
    2865}
    2866
    2867/** frees pricing norms information */
    2869 SCIP_LPI* lpi, /**< LP interface structure */
    2870 BMS_BLKMEM* blkmem, /**< block memory */
    2871 SCIP_LPINORMS** lpinorms /**< pointer to LPi pricing norms information, or NULL */
    2872 )
    2873{
    2874 assert( lpi != NULL );
    2875 assert( blkmem != NULL );
    2876 assert( lpi->solver != NULL );
    2877 assert( lpinorms != NULL );
    2878
    2879 return SCIP_OKAY;
    2880}
    2881
    2882/**@} */
    2883
    2884
    2885
    2886
    2887/*
    2888 * Parameter Methods
    2889 */
    2890
    2891/**@name Parameter Methods */
    2892/**@{ */
    2893
    2894/** gets integer parameter of LP */
    2896 SCIP_LPI* lpi, /**< LP interface structure */
    2897 SCIP_LPPARAM type, /**< parameter number */
    2898 int* ival /**< buffer to store the parameter value */
    2899 )
    2900{
    2901 assert( lpi != NULL );
    2902 assert( lpi->parameters != NULL );
    2903
    2904 /* Not (yet) supported by Glop: SCIP_LPPAR_FASTMIP, SCIP_LPPAR_POLISHING, SCIP_LPPAR_REFACTOR */
    2905 switch ( type )
    2906 {
    2908 *ival = (int) lpi->from_scratch;
    2909 SCIPdebugMessage("SCIPlpiGetIntpar: SCIP_LPPAR_FROMSCRATCH = %d.\n", *ival);
    2910 break;
    2911 case SCIP_LPPAR_LPINFO:
    2912 *ival = (int) lpi->lp_info;
    2913 SCIPdebugMessage("SCIPlpiGetIntpar: SCIP_LPPAR_LPINFO = %d.\n", *ival);
    2914 break;
    2915 case SCIP_LPPAR_LPITLIM:
    2916 *ival = (int) lpi->parameters->max_number_of_iterations();
    2917 SCIPdebugMessage("SCIPlpiGetIntpar: SCIP_LPPAR_LPITLIM = %d.\n", *ival);
    2918 break;
    2920 *ival = lpi->parameters->use_preprocessing();
    2921 SCIPdebugMessage("SCIPlpiGetIntpar: SCIP_LPPAR_PRESOLVING = %d.\n", *ival);
    2922 break;
    2923 case SCIP_LPPAR_PRICING:
    2924 *ival = lpi->pricing;
    2925 SCIPdebugMessage("SCIPlpiGetIntpar: SCIP_LPPAR_PRICING = %d.\n", *ival);
    2926 break;
    2927#ifndef NOSCALING
    2928 case SCIP_LPPAR_SCALING:
    2929 *ival = lpi->parameters->use_scaling();
    2930 SCIPdebugMessage("SCIPlpiGetIntpar: SCIP_LPPAR_SCALING = %d.\n", *ival);
    2931 break;
    2932#endif
    2933 case SCIP_LPPAR_THREADS:
    2934 *ival = lpi->numthreads;
    2935 SCIPdebugMessage("SCIPlpiGetIntpar: SCIP_LPPAR_THREADS = %d.\n", *ival);
    2936 break;
    2937 case SCIP_LPPAR_TIMING:
    2938 *ival = lpi->timing;
    2939 SCIPdebugMessage("SCIPlpiGetIntpar: SCIP_LPPAR_TIMING = %d.\n", *ival);
    2940 break;
    2942 *ival = (int) lpi->parameters->random_seed();
    2943 SCIPdebugMessage("SCIPlpiGetIntpar: SCIP_LPPAR_RANDOMSEED = %d.\n", *ival);
    2944 break;
    2945 default:
    2946 return SCIP_PARAMETERUNKNOWN;
    2947 }
    2948
    2949 return SCIP_OKAY;
    2950}
    2951
    2952/** sets integer parameter of LP */
    2954 SCIP_LPI* lpi, /**< LP interface structure */
    2955 SCIP_LPPARAM type, /**< parameter number */
    2956 int ival /**< parameter value */
    2957 )
    2958{
    2959 assert( lpi != NULL );
    2960 assert( lpi->parameters != NULL );
    2961
    2962 switch ( type )
    2963 {
    2965 SCIPdebugMessage("SCIPlpiSetIntpar: SCIP_LPPAR_FROMSCRATCH -> %d.\n", ival);
    2966 lpi->from_scratch = (bool) ival;
    2967 break;
    2968 case SCIP_LPPAR_LPINFO:
    2969 SCIPdebugMessage("SCIPlpiSetIntpar: SCIP_LPPAR_LPINFO -> %d.\n", ival);
    2970 if ( ival == 0 )
    2971 {
    2972 (void) google::SetVLOGLevel("*", google::GLOG_INFO);
    2973 lpi->lp_info = false;
    2974 }
    2975 else
    2976 {
    2977 (void) google::SetVLOGLevel("*", google::GLOG_ERROR);
    2978 lpi->lp_info = true;
    2979 }
    2980 break;
    2981 case SCIP_LPPAR_LPITLIM:
    2982 SCIPdebugMessage("SCIPlpiSetIntpar: SCIP_LPPAR_LPITLIM -> %d.\n", ival);
    2983 lpi->parameters->set_max_number_of_iterations(ival);
    2984 break;
    2986 SCIPdebugMessage("SCIPlpiSetIntpar: SCIP_LPPAR_PRESOLVING -> %d.\n", ival);
    2987 lpi->parameters->set_use_preprocessing(ival);
    2988 break;
    2989 case SCIP_LPPAR_PRICING:
    2990 SCIPdebugMessage("SCIPlpiSetIntpar: SCIP_LPPAR_PRICING -> %d.\n", ival);
    2991 lpi->pricing = (SCIP_Pricing)ival;
    2992 switch ( lpi->pricing )
    2993 {
    2995 case SCIP_PRICING_AUTO:
    2997 case SCIP_PRICING_STEEP:
    2999 lpi->parameters->set_feasibility_rule(operations_research::glop::GlopParameters_PricingRule_STEEPEST_EDGE);
    3000 break;
    3001 case SCIP_PRICING_FULL:
    3002 /* Dantzig does not really fit, but use it anyway */
    3003 lpi->parameters->set_feasibility_rule(operations_research::glop::GlopParameters_PricingRule_DANTZIG);
    3004 break;
    3005 case SCIP_PRICING_DEVEX:
    3006 lpi->parameters->set_feasibility_rule(operations_research::glop::GlopParameters_PricingRule_DEVEX);
    3007 break;
    3008 default:
    3009 return SCIP_PARAMETERUNKNOWN;
    3010 }
    3011 break;
    3012#ifndef NOSCALING
    3013 case SCIP_LPPAR_SCALING:
    3014 SCIPdebugMessage("SCIPlpiSetIntpar: SCIP_LPPAR_SCALING -> %d.\n", ival);
    3015 lpi->parameters->set_use_scaling(ival);
    3016 break;
    3017#endif
    3018 case SCIP_LPPAR_THREADS:
    3019 SCIPdebugMessage("SCIPlpiSetIntpar: SCIP_LPPAR_THREADS -> %d.\n", ival);
    3020 assert( ival >= 0 );
    3021 lpi->numthreads = ival;
    3022 if ( ival == 0 )
    3023 lpi->parameters->set_num_omp_threads(1);
    3024 else
    3025 lpi->parameters->set_num_omp_threads(ival);
    3026 break;
    3027 case SCIP_LPPAR_TIMING:
    3028 SCIPdebugMessage("SCIPlpiSetIntpar: SCIP_LPPAR_TIMING -> %d.\n", ival);
    3029 assert( 0 <= ival && ival <= 2 );
    3030 lpi->timing = ival;
    3031 if ( ival == 1 )
    3032 absl::SetFlag(&FLAGS_time_limit_use_usertime, true);
    3033 else
    3034 absl::SetFlag(&FLAGS_time_limit_use_usertime, false);
    3035 break;
    3037 SCIPdebugMessage("SCIPlpiSetIntpar: SCIP_LPPAR_RANDOMSEED -> %d.\n", ival);
    3038 assert( ival >= 0 );
    3039 lpi->parameters->set_random_seed(ival);
    3040 break;
    3041 default:
    3042 return SCIP_PARAMETERUNKNOWN;
    3043 }
    3044
    3045 return SCIP_OKAY;
    3046}
    3047
    3048/** gets floating point parameter of LP */
    3050 SCIP_LPI* lpi, /**< LP interface structure */
    3051 SCIP_LPPARAM type, /**< parameter number */
    3052 SCIP_Real* dval /**< buffer to store the parameter value */
    3053 )
    3054{
    3055 assert( lpi != NULL );
    3056 assert( lpi->parameters != NULL );
    3057
    3058 /* Not (yet) supported by Glop: SCIP_LPPAR_ROWREPSWITCH, SCIP_LPPAR_BARRIERCONVTOL */
    3059 switch ( type )
    3060 {
    3061 case SCIP_LPPAR_FEASTOL:
    3062 *dval = lpi->parameters->primal_feasibility_tolerance();
    3063 SCIPdebugMessage("SCIPlpiGetRealpar: SCIP_LPPAR_FEASTOL = %g.\n", *dval);
    3064 break;
    3066 *dval = lpi->parameters->dual_feasibility_tolerance();
    3067 SCIPdebugMessage("SCIPlpiGetRealpar: SCIP_LPPAR_DUALFEASTOL = %g.\n", *dval);
    3068 break;
    3069 case SCIP_LPPAR_OBJLIM:
    3070 if (lpi->linear_program->IsMaximizationProblem())
    3071 *dval = lpi->parameters->objective_lower_limit();
    3072 else
    3073 *dval = lpi->parameters->objective_upper_limit();
    3074 SCIPdebugMessage("SCIPlpiGetRealpar: SCIP_LPPAR_OBJLIM = %f.\n", *dval);
    3075 break;
    3076 case SCIP_LPPAR_LPTILIM:
    3077 if ( absl::GetFlag(FLAGS_time_limit_use_usertime) )
    3078 *dval = lpi->parameters->max_time_in_seconds();
    3079 else
    3080 *dval = lpi->parameters->max_deterministic_time();
    3081 SCIPdebugMessage("SCIPlpiGetRealpar: SCIP_LPPAR_LPTILIM = %f.\n", *dval);
    3082 break;
    3084 *dval = lpi->conditionlimit;
    3085 break;
    3086#ifdef SCIP_DISABLED_CODE
    3087 /* currently do not apply Markowitz parameter, since the default value does not seem suitable for Glop */
    3089 *dval = lpi->parameters->markowitz_singularity_threshold();
    3090 SCIPdebugMessage("SCIPlpiGetRealpar: SCIP_LPPAR_MARKOWITZ = %f.\n", *dval);
    3091 break;
    3092#endif
    3093 default:
    3094 return SCIP_PARAMETERUNKNOWN;
    3095 }
    3096
    3097 return SCIP_OKAY;
    3098}
    3099
    3100/** sets floating point parameter of LP */
    3102 SCIP_LPI* lpi, /**< LP interface structure */
    3103 SCIP_LPPARAM type, /**< parameter number */
    3104 SCIP_Real dval /**< parameter value */
    3105 )
    3106{
    3107 assert( lpi != NULL );
    3108 assert( lpi->parameters != NULL );
    3109
    3110 switch( type )
    3111 {
    3112 case SCIP_LPPAR_FEASTOL:
    3113 SCIPdebugMessage("SCIPlpiSetRealpar: SCIP_LPPAR_FEASTOL -> %g.\n", dval);
    3114 lpi->parameters->set_primal_feasibility_tolerance(dval);
    3115 break;
    3117 SCIPdebugMessage("SCIPlpiSetRealpar: SCIP_LPPAR_DUALFEASTOL -> %g.\n", dval);
    3118 lpi->parameters->set_dual_feasibility_tolerance(dval);
    3119 break;
    3120 case SCIP_LPPAR_OBJLIM:
    3121 SCIPdebugMessage("SCIPlpiSetRealpar: SCIP_LPPAR_OBJLIM -> %f.\n", dval);
    3122 if (lpi->linear_program->IsMaximizationProblem())
    3123 lpi->parameters->set_objective_lower_limit(dval);
    3124 else
    3125 lpi->parameters->set_objective_upper_limit(dval);
    3126 break;
    3127 case SCIP_LPPAR_LPTILIM:
    3128 SCIPdebugMessage("SCIPlpiSetRealpar: SCIP_LPPAR_LPTILIM -> %f.\n", dval);
    3129 if ( absl::GetFlag(FLAGS_time_limit_use_usertime) )
    3130 lpi->parameters->set_max_time_in_seconds(dval);
    3131 else
    3132 lpi->parameters->set_max_deterministic_time(dval);
    3133 break;
    3135 lpi->conditionlimit = dval;
    3136 lpi->checkcondition = (dval >= 0.0);
    3137 break;
    3138#ifdef SCIP_DISABLED_CODE
    3139 /* currently do not apply Markowitz parameter, since the default value does not seem suitable for Glop */
    3141 SCIPdebugMessage("SCIPlpiSetRealpar: SCIP_LPPAR_MARKOWITZ -> %f.\n", dval);
    3142 lpi->parameters->set_markowitz_singularity_threshold(dval);
    3143 break;
    3144#endif
    3145 default:
    3146 return SCIP_PARAMETERUNKNOWN;
    3147 }
    3148
    3149 return SCIP_OKAY;
    3150}
    3151
    3152/** interrupts the currently ongoing lp solve or disables the interrupt */
    3154 SCIP_LPI* lpi, /**< LP interface structure */
    3155 SCIP_Bool interrupt /**< TRUE if interrupt should be set, FALSE if it should be disabled */
    3156 )
    3157{
    3158 /*lint --e{715}*/
    3159 assert(lpi != NULL);
    3160
    3161 return SCIP_OKAY;
    3162}
    3163
    3164/**@} */
    3165
    3166
    3167
    3168
    3169/*
    3170 * Numerical Methods
    3171 */
    3172
    3173/**@name Numerical Methods */
    3174/**@{ */
    3175
    3176/** returns value treated as infinity in the LP solver */
    3178 SCIP_LPI* lpi /**< LP interface structure */
    3179 )
    3180{
    3181 assert( lpi != NULL );
    3183}
    3184
    3185/** checks if given value is treated as infinity in the LP solver */
    3187 SCIP_LPI* lpi, /**< LP interface structure */
    3188 SCIP_Real val /**< value to be checked for infinity */
    3189 )
    3190{
    3191 assert( lpi != NULL );
    3192
    3194}
    3195
    3196/**@} */
    3197
    3198
    3199
    3200
    3201/*
    3202 * File Interface Methods
    3203 */
    3204
    3205/**@name File Interface Methods */
    3206/**@{ */
    3207
    3208/** reads LP from a file */
    3210 SCIP_LPI* lpi, /**< LP interface structure */
    3211 const char* fname /**< file name */
    3212 )
    3213{
    3214 assert( lpi != NULL );
    3215 assert( lpi->linear_program != NULL );
    3216 assert( fname != NULL );
    3217
    3218 const std::string filespec(fname);
    3219 MPModelProto proto;
    3220 if ( ! ReadFileToProto(filespec, &proto) )
    3221 {
    3222 SCIPerrorMessage("Could not read <%s>\n", fname);
    3223 return SCIP_READERROR;
    3224 }
    3225 lpi->linear_program->Clear();
    3226 MPModelProtoToLinearProgram(proto, lpi->linear_program);
    3227
    3228 return SCIP_OKAY;
    3229}
    3230
    3231/** writes LP to a file */
    3233 SCIP_LPI* lpi, /**< LP interface structure */
    3234 const char* fname /**< file name */
    3235 )
    3236{
    3237 assert( lpi != NULL );
    3238 assert( lpi->linear_program != NULL );
    3239 assert( fname != NULL );
    3240
    3241 MPModelProto proto;
    3242 LinearProgramToMPModelProto(*lpi->linear_program, &proto);
    3243 const std::string filespec(fname);
    3244 if ( ! WriteProtoToFile(filespec, proto, operations_research::ProtoWriteFormat::kProtoText, true) )
    3245 {
    3246 SCIPerrorMessage("Could not write <%s>\n", fname);
    3247 return SCIP_READERROR;
    3248 }
    3249
    3250 return SCIP_OKAY;
    3251}
    3252
    3253/**@} */
    SCIP_Real * r
    Definition: circlepacking.c:59
    #define NULL
    Definition: def.h:248
    #define SCIP_Longint
    Definition: def.h:141
    #define SCIP_INVALID
    Definition: def.h:178
    #define SCIP_Bool
    Definition: def.h:91
    #define SCIP_ALLOC(x)
    Definition: def.h:366
    #define SCIP_Real
    Definition: def.h:156
    #define TRUE
    Definition: def.h:93
    #define FALSE
    Definition: def.h:94
    #define EPSCEIL(x, eps)
    Definition: def.h:192
    #define SCIPABORT()
    Definition: def.h:327
    #define EPSFLOOR(x, eps)
    Definition: def.h:191
    #define SCIP_CALL(x)
    Definition: def.h:355
    #define infinity
    Definition: gastrans.c:80
    SCIP_RETCODE SCIPlpiChgSides(SCIP_LPI *lpi, int nrows, const int *ind, const SCIP_Real *lhs, const SCIP_Real *rhs)
    Definition: lpi_glop.cpp:706
    SCIP_RETCODE SCIPlpiSetState(SCIP_LPI *lpi, BMS_BLKMEM *blkmem, const SCIP_LPISTATE *lpistate)
    Definition: lpi_glop.cpp:2731
    SCIP_RETCODE SCIPlpiGetBInvACol(SCIP_LPI *lpi, int c, SCIP_Real *coef, int *inds, int *ninds)
    Definition: lpi_glop.cpp:2629
    SCIP_RETCODE SCIPlpiGetRealpar(SCIP_LPI *lpi, SCIP_LPPARAM type, SCIP_Real *dval)
    Definition: lpi_glop.cpp:3049
    SCIP_Real SCIPlpiInfinity(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:3177
    SCIP_Bool SCIPlpiIsObjlimExc(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:2012
    SCIP_RETCODE SCIPlpiChgObjsen(SCIP_LPI *lpi, SCIP_OBJSEN objsen)
    Definition: lpi_glop.cpp:754
    SCIP_Bool SCIPlpiIsInfinity(SCIP_LPI *lpi, SCIP_Real val)
    Definition: lpi_glop.cpp:3186
    SCIP_RETCODE SCIPlpiClear(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:651
    SCIP_RETCODE SCIPlpiClearState(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:2747
    SCIP_Bool SCIPlpiExistsDualRay(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1889
    SCIP_Bool SCIPlpiExistsPrimalRay(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1826
    SCIP_RETCODE SCIPlpiGetBase(SCIP_LPI *lpi, int *cstat, int *rstat)
    Definition: lpi_glop.cpp:2335
    SCIP_RETCODE SCIPlpiReadState(SCIP_LPI *lpi, const char *fname)
    Definition: lpi_glop.cpp:2789
    SCIP_RETCODE SCIPlpiAddRows(SCIP_LPI *lpi, int nrows, const SCIP_Real *lhs, const SCIP_Real *rhs, char **rownames, int nnonz, const int *beg, const int *ind, const SCIP_Real *val)
    Definition: lpi_glop.cpp:489
    SCIP_RETCODE SCIPlpiGetPrimalRay(SCIP_LPI *lpi, SCIP_Real *ray)
    Definition: lpi_glop.cpp:2136
    SCIP_RETCODE SCIPlpiGetIntpar(SCIP_LPI *lpi, SCIP_LPPARAM type, int *ival)
    Definition: lpi_glop.cpp:2895
    SCIP_RETCODE SCIPlpiWriteLP(SCIP_LPI *lpi, const char *fname)
    Definition: lpi_glop.cpp:3232
    SCIP_RETCODE SCIPlpiSetIntegralityInformation(SCIP_LPI *lpi, int ncols, int *intInfo)
    Definition: lpi_glop.cpp:183
    SCIP_Bool SCIPlpiIsDualInfeasible(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1929
    static void deleteRowsAndUpdateCurrentBasis(SCIP_LPI *lpi, const DenseBooleanColumn &rows_to_delete)
    Definition: lpi_glop.cpp:558
    SCIP_RETCODE SCIPlpiSetRealpar(SCIP_LPI *lpi, SCIP_LPPARAM type, SCIP_Real dval)
    Definition: lpi_glop.cpp:3101
    SCIP_RETCODE SCIPlpiStrongbranchFrac(SCIP_LPI *lpi, int col_index, SCIP_Real psol, int itlim, SCIP_Real *down, SCIP_Real *up, SCIP_Bool *downvalid, SCIP_Bool *upvalid, int *iter)
    Definition: lpi_glop.cpp:1656
    SCIP_RETCODE SCIPlpiSetNorms(SCIP_LPI *lpi, BMS_BLKMEM *blkmem, const SCIP_LPINORMS *lpinorms)
    Definition: lpi_glop.cpp:2853
    SCIP_RETCODE SCIPlpiGetNNonz(SCIP_LPI *lpi, int *nnonz)
    Definition: lpi_glop.cpp:994
    SCIP_Bool SCIPlpiHasPrimalSolve(void)
    Definition: lpi_glop.cpp:209
    SCIP_RETCODE SCIPlpiStrongbranchInt(SCIP_LPI *lpi, int col, SCIP_Real psol, int itlim, SCIP_Real *down, SCIP_Real *up, SCIP_Bool *downvalid, SCIP_Bool *upvalid, int *iter)
    Definition: lpi_glop.cpp:1715
    SCIP_RETCODE SCIPlpiGetBounds(SCIP_LPI *lpi, int firstcol, int lastcol, SCIP_Real *lbs, SCIP_Real *ubs)
    Definition: lpi_glop.cpp:1228
    SCIP_Bool SCIPlpiHasBarrierSolve(void)
    Definition: lpi_glop.cpp:225
    SCIP_RETCODE SCIPlpiGetDualfarkas(SCIP_LPI *lpi, SCIP_Real *dualfarkas)
    Definition: lpi_glop.cpp:2156
    SCIP_RETCODE SCIPlpiGetObjval(SCIP_LPI *lpi, SCIP_Real *objval)
    Definition: lpi_glop.cpp:2073
    SCIP_RETCODE SCIPlpiScaleCol(SCIP_LPI *lpi, int col, SCIP_Real scaleval)
    Definition: lpi_glop.cpp:865
    int SCIPlpiGetInternalStatus(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:2047
    SCIP_RETCODE SCIPlpiStartStrongbranch(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1508
    SCIP_RETCODE SCIPlpiGetSolFeasibility(SCIP_LPI *lpi, SCIP_Bool *primalfeasible, SCIP_Bool *dualfeasible)
    Definition: lpi_glop.cpp:1803
    SCIP_RETCODE SCIPlpiFreeNorms(SCIP_LPI *lpi, BMS_BLKMEM *blkmem, SCIP_LPINORMS **lpinorms)
    Definition: lpi_glop.cpp:2868
    SCIP_Bool SCIPlpiIsIterlimExc(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:2023
    SCIP_RETCODE SCIPlpiChgBounds(SCIP_LPI *lpi, int ncols, const int *ind, const SCIP_Real *lb, const SCIP_Real *ub)
    Definition: lpi_glop.cpp:667
    SCIP_Bool SCIPlpiIsPrimalUnbounded(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1850
    SCIP_RETCODE SCIPlpiIgnoreInstability(SCIP_LPI *lpi, SCIP_Bool *success)
    Definition: lpi_glop.cpp:2058
    SCIP_RETCODE SCIPlpiWriteState(SCIP_LPI *lpi, const char *fname)
    Definition: lpi_glop.cpp:2803
    SCIP_RETCODE SCIPlpiFree(SCIP_LPI **lpi)
    Definition: lpi_glop.cpp:283
    SCIP_RETCODE SCIPlpiStrongbranchesFrac(SCIP_LPI *lpi, int *cols, int ncols, SCIP_Real *psols, int itlim, SCIP_Real *down, SCIP_Real *up, SCIP_Bool *downvalid, SCIP_Bool *upvalid, int *iter)
    Definition: lpi_glop.cpp:1685
    SCIP_RETCODE SCIPlpiGetCoef(SCIP_LPI *lpi, int row, int col, SCIP_Real *val)
    Definition: lpi_glop.cpp:1292
    SCIP_Bool SCIPlpiIsPrimalFeasible(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1874
    SCIP_RETCODE SCIPlpiReadLP(SCIP_LPI *lpi, const char *fname)
    Definition: lpi_glop.cpp:3209
    SCIP_RETCODE SCIPlpiGetRealSolQuality(SCIP_LPI *lpi, SCIP_LPSOLQUALITY qualityindicator, SCIP_Real *quality)
    Definition: lpi_glop.cpp:2195
    SCIP_Bool SCIPlpiIsDualFeasible(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1941
    SCIP_RETCODE SCIPlpiGetNorms(SCIP_LPI *lpi, BMS_BLKMEM *blkmem, SCIP_LPINORMS **lpinorms)
    Definition: lpi_glop.cpp:2836
    SCIP_Bool SCIPlpiIsTimelimExc(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:2036
    SCIP_Bool SCIPlpiHasStateBasis(SCIP_LPI *lpi, SCIP_LPISTATE *lpistate)
    Definition: lpi_glop.cpp:2777
    SCIP_RETCODE SCIPlpiSetIntpar(SCIP_LPI *lpi, SCIP_LPPARAM type, int ival)
    Definition: lpi_glop.cpp:2953
    const char * SCIPlpiGetSolverName(void)
    Definition: lpi_glop.cpp:155
    SCIP_RETCODE SCIPlpiSetBase(SCIP_LPI *lpi, const int *cstat, const int *rstat)
    Definition: lpi_glop.cpp:2369
    SCIP_Bool SCIPlpiHasPrimalRay(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1839
    SCIP_RETCODE SCIPlpiGetBInvRow(SCIP_LPI *lpi, int r, SCIP_Real *coef, int *inds, int *ninds)
    Definition: lpi_glop.cpp:2438
    SCIP_RETCODE SCIPlpiDelRows(SCIP_LPI *lpi, int firstrow, int lastrow)
    Definition: lpi_glop.cpp:588
    SCIP_RETCODE SCIPlpiGetCols(SCIP_LPI *lpi, int firstcol, int lastcol, SCIP_Real *lb, SCIP_Real *ub, int *nnonz, int *beg, int *ind, SCIP_Real *val)
    Definition: lpi_glop.cpp:1014
    SCIP_RETCODE SCIPlpiGetBInvCol(SCIP_LPI *lpi, int c, SCIP_Real *coef, int *inds, int *ninds)
    Definition: lpi_glop.cpp:2509
    SCIP_RETCODE SCIPlpiGetColNames(SCIP_LPI *lpi, int firstcol, int lastcol, char **colnames, char *namestorage, int namestoragesize, int *storageleft)
    Definition: lpi_glop.cpp:1149
    SCIP_RETCODE SCIPlpiGetBInvARow(SCIP_LPI *lpi, int r, const SCIP_Real *binvrow, SCIP_Real *coef, int *inds, int *ninds)
    Definition: lpi_glop.cpp:2573
    SCIP_RETCODE SCIPlpiGetRows(SCIP_LPI *lpi, int firstrow, int lastrow, SCIP_Real *lhs, SCIP_Real *rhs, int *nnonz, int *beg, int *ind, SCIP_Real *val)
    Definition: lpi_glop.cpp:1082
    SCIP_Bool SCIPlpiWasSolved(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1783
    const char * SCIPlpiGetSolverDesc(void)
    Definition: lpi_glop.cpp:163
    SCIP_RETCODE SCIPlpiSolveBarrier(SCIP_LPI *lpi, SCIP_Bool crossover)
    Definition: lpi_glop.cpp:1492
    SCIP_Bool SCIPlpiIsOptimal(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1954
    SCIP_RETCODE SCIPlpiGetRowNames(SCIP_LPI *lpi, int firstrow, int lastrow, char **rownames, char *namestorage, int namestoragesize, int *storageleft)
    Definition: lpi_glop.cpp:1175
    SCIP_Bool SCIPlpiHasDualSolve(void)
    Definition: lpi_glop.cpp:217
    SCIP_RETCODE SCIPlpiEndStrongbranch(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1523
    SCIP_RETCODE SCIPlpiGetSides(SCIP_LPI *lpi, int firstrow, int lastrow, SCIP_Real *lhss, SCIP_Real *rhss)
    Definition: lpi_glop.cpp:1260
    SCIP_RETCODE SCIPlpiStrongbranchesInt(SCIP_LPI *lpi, int *cols, int ncols, SCIP_Real *psols, int itlim, SCIP_Real *down, SCIP_Real *up, SCIP_Bool *downvalid, SCIP_Bool *upvalid, int *iter)
    Definition: lpi_glop.cpp:1742
    SCIP_RETCODE SCIPlpiGetSol(SCIP_LPI *lpi, SCIP_Real *objval, SCIP_Real *primsol, SCIP_Real *dualsol, SCIP_Real *activity, SCIP_Real *redcost)
    Definition: lpi_glop.cpp:2092
    SCIP_Bool SCIPlpiHasDualRay(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1904
    SCIP_RETCODE SCIPlpiDelColset(SCIP_LPI *lpi, int *dstat)
    Definition: lpi_glop.cpp:454
    SCIP_RETCODE SCIPlpiGetObj(SCIP_LPI *lpi, int firstcol, int lastcol, SCIP_Real *vals)
    Definition: lpi_glop.cpp:1201
    SCIP_RETCODE SCIPlpiFreeState(SCIP_LPI *lpi, BMS_BLKMEM *blkmem, SCIP_LPISTATE **lpistate)
    Definition: lpi_glop.cpp:2760
    SCIP_Bool SCIPlpiIsPrimalInfeasible(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1861
    SCIP_RETCODE SCIPlpiSolveDual(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1474
    SCIP_RETCODE SCIPlpiAddCols(SCIP_LPI *lpi, int ncols, const SCIP_Real *obj, const SCIP_Real *lb, const SCIP_Real *ub, char **colnames, int nnonz, const int *beg, const int *ind, const SCIP_Real *val)
    Definition: lpi_glop.cpp:352
    SCIP_RETCODE SCIPlpiSolvePrimal(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1456
    SCIP_RETCODE SCIPlpiLoadColLP(SCIP_LPI *lpi, SCIP_OBJSEN objsen, int ncols, const SCIP_Real *obj, const SCIP_Real *lb, const SCIP_Real *ub, char **colnames, int nrows, const SCIP_Real *lhs, const SCIP_Real *rhs, char **rownames, int nnonz, const int *beg, const int *ind, const SCIP_Real *val)
    Definition: lpi_glop.cpp:316
    SCIP_Bool SCIPlpiIsDualUnbounded(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1917
    SCIP_RETCODE SCIPlpiGetIterations(SCIP_LPI *lpi, int *iterations)
    Definition: lpi_glop.cpp:2176
    SCIP_RETCODE SCIPlpiGetBasisInd(SCIP_LPI *lpi, int *bind)
    Definition: lpi_glop.cpp:2402
    SCIP_RETCODE SCIPlpiCreate(SCIP_LPI **lpi, SCIP_MESSAGEHDLR *messagehdlr, const char *name, SCIP_OBJSEN objsen)
    Definition: lpi_glop.cpp:241
    void * SCIPlpiGetSolverPointer(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:171
    SCIP_RETCODE SCIPlpiChgObj(SCIP_LPI *lpi, int ncols, const int *ind, const SCIP_Real *obj)
    Definition: lpi_glop.cpp:779
    SCIP_RETCODE SCIPlpiGetObjsen(SCIP_LPI *lpi, SCIP_OBJSEN *objsen)
    Definition: lpi_glop.cpp:977
    SCIP_Bool SCIPlpiIsStable(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1971
    SCIP_RETCODE SCIPlpiGetNCols(SCIP_LPI *lpi, int *ncols)
    Definition: lpi_glop.cpp:960
    SCIP_RETCODE SCIPlpiInterrupt(SCIP_LPI *lpi, SCIP_Bool interrupt)
    Definition: lpi_glop.cpp:3153
    SCIP_RETCODE SCIPlpiDelCols(SCIP_LPI *lpi, int firstcol, int lastcol)
    Definition: lpi_glop.cpp:424
    SCIP_RETCODE SCIPlpiDelRowset(SCIP_LPI *lpi, int *dstat)
    Definition: lpi_glop.cpp:617
    SCIP_RETCODE SCIPlpiScaleRow(SCIP_LPI *lpi, int row, SCIP_Real scaleval)
    Definition: lpi_glop.cpp:802
    SCIP_RETCODE SCIPlpiGetNRows(SCIP_LPI *lpi, int *nrows)
    Definition: lpi_glop.cpp:943
    SCIP_RETCODE SCIPlpiGetState(SCIP_LPI *lpi, BMS_BLKMEM *blkmem, SCIP_LPISTATE **lpistate)
    Definition: lpi_glop.cpp:2713
    SCIP_RETCODE SCIPlpiChgCoef(SCIP_LPI *lpi, int row, int col, SCIP_Real newval)
    Definition: lpi_glop.cpp:731
    interface methods for specific LP solvers
    static SCIP_BASESTAT ConvertGlopConstraintStatus(ConstraintStatus status, Fractional dual)
    Definition: lpi_glop.cpp:2264
    static char * glopname
    Definition: lpi_glop.cpp:145
    static SCIP_BASESTAT ConvertGlopVariableStatus(VariableStatus status, Fractional rc)
    Definition: lpi_glop.cpp:2239
    static VariableStatus ConvertSCIPConstraintStatusToSlackStatus(int status)
    Definition: lpi_glop.cpp:2314
    static bool IsDualBoundValid(ProblemStatus status)
    Definition: lpi_glop.cpp:1537
    static SCIP_RETCODE strongbranch(SCIP_LPI *lpi, int col_index, SCIP_Real psol, int itlim, SCIP_Real *down, SCIP_Real *up, SCIP_Bool *downvalid, SCIP_Bool *upvalid, int *iter)
    Definition: lpi_glop.cpp:1546
    static SCIP_RETCODE SolveInternal(SCIP_LPI *lpi, bool recursive, std::unique_ptr< TimeLimit > &time_limit)
    Definition: lpi_glop.cpp:1398
    char * initGlopName()
    Definition: lpi_glop.cpp:147
    static void updateScaledLP(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1324
    static bool checkUnscaledPrimalFeasibility(SCIP_LPI *lpi)
    Definition: lpi_glop.cpp:1344
    static VariableStatus ConvertSCIPVariableStatus(int status)
    Definition: lpi_glop.cpp:2289
    #define BMSfreeMemory(ptr)
    Definition: memory.h:145
    #define BMSallocMemoryArray(ptr, num)
    Definition: memory.h:123
    #define BMSfreeMemoryArray(ptr)
    Definition: memory.h:147
    struct BMS_BlkMem BMS_BLKMEM
    Definition: memory.h:437
    #define BMSallocMemory(ptr)
    Definition: memory.h:118
    real eps
    public methods for message output
    #define SCIPerrorMessage
    Definition: pub_message.h:64
    #define SCIPdebugMessage
    Definition: pub_message.h:96
    ScatteredColumn * tmp_column
    Definition: lpi_glop.cpp:125
    int numthreads
    Definition: lpi_glop.cpp:112
    operations_research::glop::LinearProgram * linear_program
    Definition: lpi_glop.cpp:98
    ScatteredRow * tmp_row
    Definition: lpi_glop.cpp:124
    operations_research::glop::LinearProgram * scaled_lp
    Definition: lpi_glop.cpp:99
    operations_research::glop::GlopParameters * parameters
    Definition: lpi_glop.cpp:101
    bool lp_info
    Definition: lpi_glop.cpp:109
    SCIP_Real conditionlimit
    Definition: lpi_cpx.c:174
    int timing
    Definition: lpi_glop.cpp:115
    bool lp_time_limit_was_reached
    Definition: lpi_glop.cpp:106
    bool lp_modified_since_last_solve
    Definition: lpi_glop.cpp:105
    operations_research::glop::LpScalingHelper * scaler
    Definition: lpi_glop.cpp:102
    SCIP_Longint niterations
    Definition: lpi_glop.cpp:118
    SCIP_PRICING pricing
    Definition: lpi_clp.cpp:112
    operations_research::glop::RevisedSimplex * solver
    Definition: lpi_glop.cpp:100
    bool checkcondition
    Definition: lpi_glop.cpp:114
    bool from_scratch
    Definition: lpi_glop.cpp:111
    SCIP_Bool checkcondition
    Definition: lpi_cpx.c:175
    SCIP_Pricing
    Definition: type_lpi.h:77
    @ SCIP_PRICING_STEEPQSTART
    Definition: type_lpi.h:83
    @ SCIP_PRICING_AUTO
    Definition: type_lpi.h:79
    @ SCIP_PRICING_DEVEX
    Definition: type_lpi.h:84
    @ SCIP_PRICING_STEEP
    Definition: type_lpi.h:82
    @ SCIP_PRICING_FULL
    Definition: type_lpi.h:80
    @ SCIP_PRICING_LPIDEFAULT
    Definition: type_lpi.h:78
    @ SCIP_PRICING_PARTIAL
    Definition: type_lpi.h:81
    enum SCIP_Pricing SCIP_PRICING
    Definition: type_lpi.h:86
    enum SCIP_LPParam SCIP_LPPARAM
    Definition: type_lpi.h:73
    @ SCIP_LPSOLQUALITY_EXACTCONDITION
    Definition: type_lpi.h:102
    @ SCIP_LPSOLQUALITY_ESTIMCONDITION
    Definition: type_lpi.h:101
    @ SCIP_LPPAR_PRICING
    Definition: type_lpi.h:54
    @ SCIP_LPPAR_THREADS
    Definition: type_lpi.h:66
    @ SCIP_LPPAR_LPINFO
    Definition: type_lpi.h:55
    @ SCIP_LPPAR_SCALING
    Definition: type_lpi.h:52
    @ SCIP_LPPAR_TIMING
    Definition: type_lpi.h:68
    @ SCIP_LPPAR_LPTILIM
    Definition: type_lpi.h:61
    @ SCIP_LPPAR_PRESOLVING
    Definition: type_lpi.h:53
    @ SCIP_LPPAR_CONDITIONLIMIT
    Definition: type_lpi.h:67
    @ SCIP_LPPAR_RANDOMSEED
    Definition: type_lpi.h:69
    @ SCIP_LPPAR_DUALFEASTOL
    Definition: type_lpi.h:57
    @ SCIP_LPPAR_FROMSCRATCH
    Definition: type_lpi.h:50
    @ SCIP_LPPAR_MARKOWITZ
    Definition: type_lpi.h:62
    @ SCIP_LPPAR_FEASTOL
    Definition: type_lpi.h:56
    @ SCIP_LPPAR_LPITLIM
    Definition: type_lpi.h:60
    @ SCIP_LPPAR_OBJLIM
    Definition: type_lpi.h:59
    @ SCIP_BASESTAT_BASIC
    Definition: type_lpi.h:92
    @ SCIP_BASESTAT_UPPER
    Definition: type_lpi.h:93
    @ SCIP_BASESTAT_LOWER
    Definition: type_lpi.h:91
    @ SCIP_BASESTAT_ZERO
    Definition: type_lpi.h:94
    enum SCIP_LPSolQuality SCIP_LPSOLQUALITY
    Definition: type_lpi.h:104
    @ SCIP_OBJSEN_MAXIMIZE
    Definition: type_lpi.h:42
    @ SCIP_OBJSEN_MINIMIZE
    Definition: type_lpi.h:43
    enum SCIP_ObjSen SCIP_OBJSEN
    Definition: type_lpi.h:45
    enum SCIP_BaseStat SCIP_BASESTAT
    Definition: type_lpi.h:96
    @ SCIP_LPERROR
    Definition: type_retcode.h:49
    @ SCIP_READERROR
    Definition: type_retcode.h:45
    @ SCIP_INVALIDDATA
    Definition: type_retcode.h:52
    @ SCIP_PARAMETERUNKNOWN
    Definition: type_retcode.h:55
    @ SCIP_OKAY
    Definition: type_retcode.h:42
    @ SCIP_NOTIMPLEMENTED
    Definition: type_retcode.h:61
    enum SCIP_Retcode SCIP_RETCODE
    Definition: type_retcode.h:63