Scippy

    SCIP

    Solving Constraint Integer Programs

    lpiexact_spx.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 lpiexact_spx.cpp
    26 * @ingroup LPIEXACTS
    27 * @brief exact LP interface for SoPlex
    28 * @author Leon Eifler
    29 *
    30 * This is an implementation of SCIP's LP interface for SoPlex.
    31 *
    32 * For debugging purposes, the SoPlex results can be double checked with CPLEX if WITH_LPSCHECK is defined. This may
    33 * yield false positives, since the LP is dumped to a file for transferring it to CPLEX, hence, precision may be lost.
    34 */
    35
    36/*--+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
    37#include "scip/def.h"
    39#include "scip/type_retcode.h"
    41#include "lpiexact/lpiexact.h"
    42#include "scip/rational.h"
    43#include "scip/rationalgmp.h"
    44
    45/* check the return value of setParam methods */
    46#define CHECK_SOPLEX_PARAM(x) \
    47 if( !x ) \
    48 { \
    49 SCIPmessagePrintWarning(_messagehdlr, "SoPlex: unsupported parameter value\n"); \
    50 }
    51
    52/* remember the original value of the SCIP_DEBUG define and undefine it */
    53#ifdef SCIP_DEBUG
    54#define ___DEBUG
    55#undef SCIP_DEBUG
    56#endif
    57
    58/* disable -Wclass-memaccess warnings due to dubious memcpy/realloc calls in SoPlex headers, e.g.,
    59 * dataarray.h:314:16: warning: ‘void* memcpy(void*, const void*, size_t)’ writing to an object of type ‘struct soplex::SPxParMultPR::SPxParMultPr_Tmp’ with no trivial copy-assignment; use copy-assignment or copy-initialization instead [-Wclass-memaccess]
    60 */
    61#ifdef __GNUC__
    62#if __GNUC__ >= 8
    63#pragma GCC diagnostic ignored "-Wclass-memaccess"
    64#endif
    65#endif
    66
    67/* include SoPlex solver */
    68#include "soplex.h"
    69
    70/* define subversion for versions <= 1.5.0.1 */
    71#ifndef SOPLEX_SUBVERSION
    72#define SOPLEX_SUBVERSION 0
    73#endif
    74/* define API version for versions <= 3.0.0 */
    75#ifndef SOPLEX_APIVERSION
    76#define SOPLEX_APIVERSION 0
    77#endif
    78
    79/* check version */
    80#if (SOPLEX_APIVERSION < 11)
    81#error "This interface requires SoPlex with API version 11 or higher, which is available from SoPlex 5.0.0."
    82#endif
    83
    84#if (SOPLEX_APIVERSION <= 5)
    85#include "spxgithash.h"
    86#endif
    87
    88#ifndef SOPLEX_WITH_GMP
    89#error "SOPLEX_WITH_GMP not defined: SoPlex must be build with GMP to be used with LPSEXACT=spx"
    90#endif
    91
    92#ifndef SOPLEX_WITH_MPFR
    93#error "SOPLEX_WITH_MPFR not defined: SoPlex must be build with MPFR to be used with LPSEXACT=spx"
    94#endif
    95
    96#ifndef SOPLEX_WITH_BOOST
    97#error "SOPLEX_WITH_BOOST not defined: SoPlex must be build with BOOST to be used with LPSEXACT=spx"
    98#endif
    99
    100/* reset the SCIP_DEBUG define to its original SCIP value */
    101#undef SCIP_DEBUG
    102#ifdef ___DEBUG
    103#define SCIP_DEBUG
    104#undef ___DEBUG
    105#endif
    106
    107/* define snprintf when using a too old MSVC version */
    108#if defined(_MSC_VER) && _MSC_VER < 1900
    109#ifndef snprintf
    110#define snprintf _snprintf
    111#endif
    112#endif
    113
    114#define SOPLEX_VERBLEVEL 5 /**< verbosity level for LPINFO */
    115
    116#include "scip/pub_message.h"
    117#include "scip/struct_rational.h"
    118
    119/********************************************************************/
    120/*----------------------------- C++ --------------------------------*/
    121/********************************************************************/
    122
    123/* in C++ we have to use "0" instead of "(void*)0" */
    124#undef NULL
    125#define NULL 0
    126
    127#include <cassert>
    128using namespace soplex;
    129
    130
    131/** Macro for a single SoPlex call for which exceptions have to be catched - return an LP error. We
    132 * make no distinction between different exception types, e.g., between memory allocation and other
    133 * exceptions.
    134 */
    135#ifndef NDEBUG
    136#define SOPLEX_TRY(messagehdlr, x) do \
    137 { \
    138 try \
    139 { \
    140 (x); \
    141 } \
    142 catch( const SPxMemoryException& E ) \
    143 { \
    144 std::string s = E.what(); \
    145 SCIPerrorMessage("SoPlex threw a memory exception: %s\n", s.c_str()); \
    146 return SCIP_ERROR; \
    147 } \
    148 catch( const SPxException& E ) \
    149 { \
    150 std::string s = E.what(); \
    151 SCIPmessagePrintWarning((messagehdlr), "SoPlex threw an exception: %s\n", s.c_str()); \
    152 return SCIP_LPERROR; \
    153 } \
    154 } \
    155 while( FALSE )
    156
    157#else
    158#define SOPLEX_TRY(messagehdlr, x) do \
    159 { \
    160 try \
    161 { \
    162 (x); \
    163 } \
    164 catch( const SPxMemoryException& E ) \
    165 { \
    166 std::string s = E.what(); \
    167 SCIPerrorMessage("SoPlex threw a memory exception: %s\n", s.c_str()); \
    168 return SCIP_ERROR; \
    169 } \
    170 catch( const SPxException& ) \
    171 { \
    172 return SCIP_LPERROR; \
    173 } \
    174 } \
    175 while( FALSE )
    176#endif
    177
    178/* Macro for a single SoPlex call for which exceptions have to be catched - abort if they
    179 * arise. SCIP_ABORT() is not accessible here.
    180 */
    181#define SOPLEX_TRY_ABORT(x) do \
    182 { \
    183 try \
    184 { \
    185 (x); \
    186 } \
    187 catch( const SPxException& E ) \
    188 { \
    189 std::string s = E.what(); \
    190 SCIPerrorMessage("SoPlex threw an exception: %s\n", s.c_str()); \
    191 abort(); \
    192 } \
    193 } \
    194 while( FALSE )
    195
    196/* Set the value of a SCIP_RATIONAL* from a SoPlex Rational */
    197static void RsetSpxR(
    198 SCIP_LPIEXACT* lpi, /**< exact lpi*/
    199 SCIP_RATIONAL* r, /**< scip rational */
    200 const soplex::Rational& spxr /**< soplex rational */
    201 )
    202{
    203 if( SCIPlpiExactIsInfinity(lpi, double(spxr)) )
    205 else if( SCIPlpiExactIsInfinity(lpi, -double(spxr)) )
    207 else
    208 {
    209#if defined(SCIP_WITH_BOOST) && defined(SCIP_WITH_GMP)
    210 r->val = spxr;
    211#else
    212 r->val = (double) spxr;
    213#endif
    214 r->isinf = FALSE;
    215 r->isfprepresentable = SCIP_ISFPREPRESENTABLE_UNKNOWN;
    216 }
    217}
    218
    219/* Set the value of a SoPlex Rational vector from a SCIP_RATIONAL** */
    220static void RsetSpxVector(
    221 SCIP_LPIEXACT* lpi, /**< exact LPI */
    222 SCIP_RATIONAL** r, /**< SCIP_RATIONAL array */
    223 VectorRational src /**< SoPlex rational vector */
    224 )
    225{
    226 for( int i = 0; i < src.dim(); ++i )
    227 {
    228 assert(r[i] != NULL);
    229
    230 RsetSpxR(lpi, r[i], src[i]);
    231 }
    232}
    233
    234/** set the value of a SoPlex rational from a SCIP_RATIONAL */
    235static void SpxRSetRat(
    236 SCIP_LPIEXACT* lpi, /**< exact LPI */
    237 soplex::Rational& spxr, /**< SoPlex Rational*/
    238 SCIP_RATIONAL* src /**< SCIP_RATIONAL */
    239 )
    240{
    242 {
    243 if( SCIPrationalIsPositive(src) )
    244 spxr = soplex::Rational(SCIPlpiExactInfinity(lpi));
    245 else
    246 spxr = soplex::Rational(-SCIPlpiExactInfinity(lpi));
    247 }
    248 else
    249 {
    250#if defined(SOPLEX_WITH_GMP) && defined(SCIP_WITH_BOOST) && defined(SCIP_WITH_GMP)
    251 spxr = soplex::Rational(*SCIPrationalGetGMP(src));
    252#else
    253 spxr = soplex::Rational(SCIPrationalGetReal(src));
    254#endif
    255 }
    256}
    257
    258
    259
    260/** SCIP's SoPlex class */
    261class SPxexSCIP : public SoPlex
    262{/*lint !e1790*/
    263 bool _lpinfo;
    264 bool _fromscratch;
    265 char* _probname;
    266 DataArray<SPxSolver::VarStatus> _colStat; /**< column basis status used for strong branching */
    267 DataArray<SPxSolver::VarStatus> _rowStat; /**< row basis status used for strong branching */
    268 SCIP_MESSAGEHDLR* _messagehdlr; /**< messagehdlr handler for printing messages, or NULL */
    269
    270public:
    272 SCIP_MESSAGEHDLR* messagehdlr = NULL, /**< message handler */
    273 const char* probname = NULL /**< name of problem */
    274 )
    275 : _lpinfo(false),
    276 _fromscratch(false),
    277 _probname(NULL),
    278 _colStat(0),
    279 _rowStat(0),
    280 _messagehdlr(messagehdlr)
    281 {
    282 if ( probname != NULL )
    283 SOPLEX_TRY_ABORT( setProbname(probname) );
    284
    285#if SOPLEX_APIVERSION >= 2
    286 (void) setBoolParam(SoPlex::ENSURERAY, true);
    287#endif
    288 }
    289
    290 virtual ~SPxexSCIP()
    291 {
    292 if( _probname != NULL )
    293 spx_free(_probname); /*lint !e1551*/
    294
    295 freePreStrongbranchingBasis(); /*lint !e1551*/
    296 }/*lint -e1579*/
    297
    298 /** get objective limit according to objective sense */
    299 Real getObjLimit() const
    300 {
    301 return (intParam(SoPlex::OBJSENSE) == SoPlex::OBJSENSE_MINIMIZE)
    302 ? realParam(SoPlex::OBJLIMIT_UPPER)
    303 : realParam(SoPlex::OBJLIMIT_LOWER);
    304 }
    305
    306 // @todo realize this with a member variable as before
    307 bool getFromScratch() const
    308 {
    309 return _fromscratch;
    310 }
    311
    312 void setFromScratch(bool fs)
    313 {
    314 _fromscratch = fs;
    315 }
    316
    317 // @todo member variable?
    318 bool getLpInfo() const
    319 {
    320 return _lpinfo;
    321 }
    322
    323 void setLpInfo(bool lpinfo)
    324 {
    325 _lpinfo = lpinfo;
    326 }
    327
    328 // @todo member variable?
    329 void setProbname(const char* probname)
    330 {
    331 size_t len;
    332
    333 assert(probname != NULL);
    334 if( _probname != NULL )
    335 spx_free(_probname);
    336
    337 len = strlen(probname);
    338 spx_alloc(_probname, (int) len + 1);
    339 memcpy(_probname, probname, len + 1);
    340 }
    341
    342 void setRep(SPxSolver::Representation p_rep)
    343 {
    344 if( p_rep == SPxSolver::COLUMN && intParam(SoPlex::REPRESENTATION) == SoPlex::REPRESENTATION_ROW )
    345 {
    346 SCIPdebugMessage("switching to column representation of the basis\n");
    347 CHECK_SOPLEX_PARAM(setIntParam(SoPlex::REPRESENTATION, (int) SoPlex::REPRESENTATION_COLUMN));
    348 }
    349 else if( (p_rep == SPxSolver::ROW && intParam(SoPlex::REPRESENTATION) == SoPlex::REPRESENTATION_COLUMN) )
    350 {
    351 SCIPdebugMessage("switching to row representation of the basis\n");
    352 CHECK_SOPLEX_PARAM(setIntParam(SoPlex::REPRESENTATION, (int) SoPlex::REPRESENTATION_ROW));
    353 }
    354 }
    355
    356#ifndef NDEBUG
    358 {
    359 for( int i = 0; i < numColsRational(); ++i )
    360 {
    361 if( lowerRational(i) > upperRational(i) )
    362 {
    363 SCIPerrorMessage("inconsistent bounds on column %d: lower=%s, upper=%s\n",
    364 i, lowerRational(i).str().c_str(), upperRational(i).str().c_str());
    365 return false;
    366 }
    367 }
    368
    369 return true;
    370 }
    371
    373 {
    374 for( int i = 0; i < numRowsRational(); ++i )
    375 {
    376 if( lhsRational(i) > rhsRational(i) )
    377 {
    378 SCIPerrorMessage("inconsistent sides on row %d: lhs=%s, rhs=%s\n",
    379 i, lhsRational(i).str().c_str(), rhsRational(i).str().c_str());
    380 return false;
    381 }
    382 }
    383
    384 return true;
    385 }
    386#endif
    387
    388 void trySolve(bool printwarning = true)
    389 {
    390 Real timespent;
    391 Real timelimit;
    392
    393 try
    394 {
    395 (void) optimize();
    396 }
    397 catch(const SPxException& x)
    398 {
    399 std::string s = x.what();
    400 if( printwarning )
    401 {
    402 SCIPmessagePrintWarning(_messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    403 }
    404
    405 /* since it is not clear if the status in SoPlex are set correctly
    406 * we want to make sure that if an error is thrown the status is
    407 * not OPTIMAL anymore.
    408 */
    409 assert(status() != SPxSolver::OPTIMAL);
    410 }
    411
    412 assert(intParam(ITERLIMIT) < 0 || numIterations() <= intParam(ITERLIMIT));
    413
    414 /* update time limit */
    415 timespent = solveTime();
    416 if( timespent > 0 )
    417 {
    418 /* get current time limit */
    419 timelimit = realParam(TIMELIMIT);
    420 if( timelimit > timespent )
    421 timelimit -= timespent;
    422 else
    423 timelimit = 0;
    424 /* set new time limit */
    425 assert(timelimit >= 0);
    426 CHECK_SOPLEX_PARAM(setRealParam(TIMELIMIT, timelimit));
    427 }
    428 }
    429
    430 SPxSolver::Status doSolve(bool printwarning = true)
    431 {
    432 SPxOut::Verbosity verbosity;
    433
    434 SPxSolver::Status spxStatus;
    435
    436 /* store and set verbosity */
    437 verbosity = spxout.getVerbosity();
    438 spxout.setVerbosity((SPxOut::Verbosity)(getLpInfo() ? SOPLEX_VERBLEVEL : 0));
    439
    440 assert(checkConsistentBounds());
    441 assert(checkConsistentSides());
    442
    443 trySolve(printwarning);
    444 spxStatus = status();
    445
    446 /* restore verbosity */
    447 spxout.setVerbosity(verbosity);
    448
    449 return spxStatus;
    450 }
    451
    452 /** save the current basis */
    454 {
    455 _rowStat.reSize(numRowsRational());
    456 _colStat.reSize(numColsRational());
    457
    458 try
    459 {
    460 getBasis(_rowStat.get_ptr(), _colStat.get_ptr());
    461 }
    462#ifndef NDEBUG
    463 catch(const SPxException& x)
    464 {
    465 std::string s = x.what();
    466 SCIPmessagePrintWarning(_messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    467
    468 /* since it is not clear if the status in SoPlex are set correctly
    469 * we want to make sure that if an error is thrown the status is
    470 * not OPTIMAL anymore.
    471 */
    472 assert(status() != SPxSolver::OPTIMAL);
    473 }
    474#else
    475 catch(const SPxException&)
    476 { }
    477#endif
    478 }
    479
    480 /** restore basis */
    482 {
    483 assert(_rowStat.size() == numRowsRational());
    484 assert(_colStat.size() == numColsRational());
    485
    486 try
    487 {
    488 setBasis(_rowStat.get_ptr(), _colStat.get_ptr());
    489 }
    490#ifndef NDEBUG
    491 catch(const SPxException& x)
    492 {
    493 std::string s = x.what();
    494 SCIPmessagePrintWarning(_messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    495#else
    496 catch(const SPxException&)
    497 {
    498#endif
    499 /* since it is not clear if the status in SoPlex are set correctly
    500 * we want to make sure that if an error is thrown the status is
    501 * not OPTIMAL anymore.
    502 */
    503 assert(status() != SPxSolver::OPTIMAL);
    504 }
    505 }
    506
    507 /** if basis is in store, delete it without restoring it */
    509 {
    510 _rowStat.clear();
    511 _colStat.clear();
    512 }
    513
    514 /** is pre-strong-branching basis freed? */
    516 {
    517 return ((_rowStat.size() == 0 ) && (_colStat.size() == 0));
    518 }
    519
    520 /** provides access for temporary storage of basis status of rows */
    521 DataArray<SPxSolver::VarStatus>& rowStat()
    522 {
    523 return _rowStat; /*lint !e1536*/
    524 }
    525
    526 /** provides access for temporary storage of basis status or columns */
    527 DataArray<SPxSolver::VarStatus>& colStat()
    528 {
    529 return _colStat; /*lint !e1536*/
    530 }
    531
    532}; /*lint !e1748*/
    533
    534
    535
    536
    537/********************************************************************/
    538/*----------------------------- C --------------------------------*/
    539/********************************************************************/
    540
    541#include "lpiexact/lpiexact.h"
    542#include "scip/bitencode.h"
    543
    544typedef SCIP_DUALPACKET COLPACKET; /* each column needs two bits of information (basic/on_lower/on_upper) */
    545#define COLS_PER_PACKET SCIP_DUALPACKETSIZE
    546typedef SCIP_DUALPACKET ROWPACKET; /* each row needs two bit of information (basic/on_lower/on_upper) */
    547#define ROWS_PER_PACKET SCIP_DUALPACKETSIZE
    548
    549
    550
    551/** LP interface */
    552struct SCIP_LPiExact
    553{
    554 SPxexSCIP* spx; /**< our SoPlex implementation */
    555 int* cstat; /**< array for storing column basis status */
    556 int* rstat; /**< array for storing row basis status */
    557 int cstatsize; /**< size of cstat array */
    558 int rstatsize; /**< size of rstat array */
    559 SCIP_PRICING pricing; /**< current pricing strategy */
    560 SCIP_Bool solved; /**< was the current LP solved? */
    561 SCIP_Real conditionlimit; /**< maximum condition number of LP basis counted as stable (-1.0: no limit) */
    562 SCIP_Bool checkcondition; /**< should condition number of LP basis be checked for stability? */
    563 SCIP_MESSAGEHDLR* messagehdlr; /**< messagehdlr handler to printing messages, or NULL */
    564};
    565
    566/** LPi state stores basis information */
    567struct SCIP_LPiState
    568{
    569 int ncols; /**< number of LP columns */
    570 int nrows; /**< number of LP rows */
    571 COLPACKET* packcstat; /**< column basis status in compressed form */
    572 ROWPACKET* packrstat; /**< row basis status in compressed form */
    573};
    574
    575
    576/*
    577 * dynamic memory arrays
    578 */
    579
    580/** resizes cstat array to have at least num entries */
    581static
    583 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    584 int num /**< minimal number of entries in array */
    585 )
    586{
    587 assert(lpi != NULL);
    588
    589 if( num > lpi->cstatsize )
    590 {
    591 int newsize;
    592
    593 newsize = MAX(2*lpi->cstatsize, num);
    594 SCIP_ALLOC( BMSreallocMemoryArray(&lpi->cstat, newsize) );
    595 lpi->cstatsize = newsize;
    596 }
    597 assert(num <= lpi->cstatsize);
    598
    599 return SCIP_OKAY;
    600}
    601
    602/** resizes rstat array to have at least num entries */
    603static
    605 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    606 int num /**< minimal number of entries in array */
    607 )
    608{
    609 assert(lpi != NULL);
    610
    611 if( num > lpi->rstatsize )
    612 {
    613 int newsize;
    614
    615 newsize = MAX(2*lpi->rstatsize, num);
    616 SCIP_ALLOC( BMSreallocMemoryArray(&lpi->rstat, newsize) );
    617 lpi->rstatsize = newsize;
    618 }
    619 assert(num <= lpi->rstatsize);
    620
    621 return SCIP_OKAY;
    622}
    623
    624
    625
    626
    627/*
    628 * LPi state methods
    629 */
    630
    631/** returns the number of packets needed to store column packet information */
    632static
    634 int ncols /**< number of columns to store */
    635 )
    636{
    637 return (ncols+(int)COLS_PER_PACKET-1)/(int)COLS_PER_PACKET;
    638}
    639
    640/** returns the number of packets needed to store row packet information */
    641static
    643 int nrows /**< number of rows to store */
    644 )
    645{
    646 return (nrows+(int)ROWS_PER_PACKET-1)/(int)ROWS_PER_PACKET;
    647}
    648
    649/** store row and column basis status in a packed LPi state object */
    650static
    652 SCIP_LPISTATE* lpistate, /**< pointer to LPi state data */
    653 const int* cstat, /**< basis status of columns in unpacked format */
    654 const int* rstat /**< basis status of rows in unpacked format */
    655 )
    656{
    657 assert(lpistate != NULL);
    658 assert(lpistate->packcstat != NULL);
    659 assert(lpistate->packrstat != NULL);
    660
    661 SCIPencodeDualBit(cstat, lpistate->packcstat, lpistate->ncols);
    662 SCIPencodeDualBit(rstat, lpistate->packrstat, lpistate->nrows);
    663}
    664
    665/** unpacks row and column basis status from a packed LPi state object */
    666static
    668 const SCIP_LPISTATE* lpistate, /**< pointer to LPi state data */
    669 int* cstat, /**< buffer for storing basis status of columns in unpacked format */
    670 int* rstat /**< buffer for storing basis status of rows in unpacked format */
    671 )
    672{
    673 assert(lpistate != NULL);
    674 assert(lpistate->packcstat != NULL);
    675 assert(lpistate->packrstat != NULL);
    676
    677 SCIPdecodeDualBit(lpistate->packcstat, cstat, lpistate->ncols);
    678 SCIPdecodeDualBit(lpistate->packrstat, rstat, lpistate->nrows);
    679}
    680
    681/** creates LPi state information object */
    682static
    684 SCIP_LPISTATE** lpistate, /**< pointer to LPi state */
    685 BMS_BLKMEM* blkmem, /**< block memory */
    686 int ncols, /**< number of columns to store */
    687 int nrows /**< number of rows to store */
    688 )
    689{
    690 assert(lpistate != NULL);
    691 assert(blkmem != NULL);
    692 assert(ncols >= 0);
    693 assert(nrows >= 0);
    694
    695 int nColPackets = colpacketNum(ncols);
    696 int nRowPackets = rowpacketNum(nrows);
    697
    698 SCIP_ALLOC( BMSallocBlockMemory(blkmem, lpistate) );
    699 SCIP_ALLOC( BMSallocBlockMemoryArray(blkmem, &(*lpistate)->packcstat, nColPackets) );
    700 SCIP_ALLOC( BMSallocBlockMemoryArray(blkmem, &(*lpistate)->packrstat, nRowPackets) );
    701
    702 return SCIP_OKAY;
    703}
    704
    705/** frees LPi state information */
    706static
    708 SCIP_LPISTATE** lpistate, /**< pointer to LPi state information (like basis information) */
    709 BMS_BLKMEM* blkmem /**< block memory */
    710 )
    711{
    712 assert(blkmem != NULL);
    713 assert(lpistate != NULL);
    714 assert(*lpistate != NULL);
    715
    716 int nColPackets = colpacketNum((*lpistate)->ncols);
    717 int nRowPackets = rowpacketNum((*lpistate)->nrows);
    718
    719 BMSfreeBlockMemoryArray(blkmem, &(*lpistate)->packcstat, nColPackets);
    720 BMSfreeBlockMemoryArray(blkmem, &(*lpistate)->packrstat, nRowPackets);
    721 BMSfreeBlockMemory(blkmem, lpistate);
    722}
    723
    724
    725
    726
    727/*
    728 * local methods
    729 */
    730
    731
    732/** marks the current LP to be unsolved */
    733static
    735{
    736 assert(lpi != NULL);
    737 lpi->solved = FALSE;
    738}
    739
    740
    741
    742/*
    743 * LP Interface Methods
    744 */
    745
    746
    747/*
    748 * Miscellaneous Methods
    749 */
    750
    751#if (SOPLEX_SUBVERSION > 0)
    752 const static char spxname[20] = {'S', 'o', 'P', 'l', 'e', 'x', ' ', SOPLEX_VERSION/100 + '0', '.', (SOPLEX_VERSION % 100)/10 + '0', '.', SOPLEX_VERSION % 10 + '0', '.', SOPLEX_SUBVERSION + '0', '\0'};
    753#else
    754 const static char spxname[20] = {'S', 'o', 'P', 'l', 'e', 'x', ' ', SOPLEX_VERSION/100 + '0', '.', (SOPLEX_VERSION % 100)/10 + '0', '.', SOPLEX_VERSION % 10 + '0', '\0'};
    755#endif
    756const static char spxdesc[200] = {'L', 'i', 'n', 'e', 'a', 'r', ' ', 'p', 'r', 'o', 'g', 'r', 'a', 'm', 'm', 'i', 'n', 'g',
    757 ' ', 's', 'o', 'l', 'v', 'e', 'r', ' ' , 'd', 'e', 'v', 'e', 'l', 'o', 'p', 'e', 'd',
    758 ' ', 'a', 't', ' ', 'Z', 'u', 's', 'e', ' ', 'I', 'n', 's', 't', 'i', 't', 'u', 't', 'e',
    759 ' ', 'B', 'e', 'r', 'l', 'i', 'n', ' ', '(', 's', 'o', 'p', 'l', 'e', 'x', '.', 'z', 'i', 'b', '.', 'd', 'e', ')',
    760 ' ', '[', 'G', 'i', 't', 'H', 'a', 's', 'h', ':', ' ',
    761 getGitHash()[0], getGitHash()[1], getGitHash()[2], getGitHash()[3],
    762 getGitHash()[4], getGitHash()[5], getGitHash()[6], getGitHash()[7],
    763 ']', '\0'};
    764
    765/**@name Miscellaneous Methods */
    766/**@{ */
    767
    768/** gets name and version of LP solver */
    770 void
    771 )
    772{
    773 SCIPdebugMessage("calling SCIPlpiExactGetSolverName()\n");
    774 return spxname;
    775}
    776
    777/** gets description of LP solver (developer, webpage, ...) */
    779 void
    780 )
    781{
    782 return spxdesc;
    783}
    784
    785/** gets pointer for LP solver - use only with great care */
    787 SCIP_LPIEXACT* lpi /**< pointer to an LP interface structure */
    788 )
    789{
    790 return (void*) lpi->spx;
    791}
    792
    793/** pass integrality information about variables to the solver */
    795 SCIP_LPIEXACT* lpi, /**< pointer to an LP interface structure */
    796 int ncols, /**< length of integrality array */
    797 int* intInfo /**< integrality array (0: continuous, 1: integer). May be NULL iff ncols is 0. */
    798 )
    799{
    800 assert(ncols == lpi->spx->numColsRational() || (ncols == 0 && intInfo == NULL));
    801 lpi->spx->setIntegralityInformation(ncols, intInfo);
    802 return SCIP_OKAY;
    803}
    804
    805/** informs about availability of a primal simplex solving method */
    807 void
    808 )
    809{
    810 return TRUE;
    811}
    812
    813/** informs about availability of a dual simplex solving method */
    815 void
    816 )
    817{
    818 return TRUE;
    819}
    820
    821/** informs about availability of a barrier solving method */
    823 void
    824 )
    825{
    826 return FALSE;
    827}
    828
    829/**@} */
    830
    831
    832/*
    833 * LPI Creation and Destruction Methods
    834 */
    835
    836/**@name LPI Creation and Destruction Methods */
    837/**@{ */
    838
    839/** creates an LP problem object */
    841 SCIP_LPIEXACT** lpi, /**< pointer to an LP interface structure */
    842 SCIP_MESSAGEHDLR* messagehdlr, /**< message handler to use for printing messages, or NULL */
    843 const char* name, /**< problem name */
    844 SCIP_OBJSEN objsen /**< objective sense */
    845 )
    846{
    847 assert(lpi != NULL);
    848 assert(name != NULL);
    849
    850 /* create SoPlex object */
    852
    853 /* we use this construction to allocate the memory for the SoPlex class also via the blockmemshell */
    854 (*lpi)->spx = static_cast<SPxexSCIP*>(BMSallocMemoryCPP(sizeof(SPxexSCIP)));
    855 SOPLEX_TRY( messagehdlr, (*lpi)->spx = new ((*lpi)->spx) SPxexSCIP(messagehdlr, name) );
    856 (void) (*lpi)->spx->setIntParam(SoPlex::SYNCMODE, (int) SoPlex::SYNCMODE_AUTO);
    857 (void) (*lpi)->spx->setIntParam(SoPlex::SOLVEMODE, (int) SoPlex::SOLVEMODE_RATIONAL);
    858 (void) (*lpi)->spx->setRealParam(SoPlex::FEASTOL, 0.0, true);
    859 (void) (*lpi)->spx->setRealParam(SoPlex::OPTTOL, 0.0, true);
    860 (void) (*lpi)->spx->setBoolParam(SoPlex::FORCEBASIC, true);
    861
    862 (*lpi)->cstat = NULL;
    863 (*lpi)->rstat = NULL;
    864 (*lpi)->cstatsize = 0;
    865 (*lpi)->rstatsize = 0;
    866 (*lpi)->pricing = SCIP_PRICING_LPIDEFAULT;
    867 (*lpi)->conditionlimit = -1.0;
    868 (*lpi)->checkcondition = FALSE;
    869 (*lpi)->messagehdlr = messagehdlr;
    870
    871 invalidateSolution(*lpi);
    872
    873 /* set objective sense */
    874 SCIP_CALL( SCIPlpiExactChgObjsen(*lpi, objsen) );
    875
    876 /* set default pricing */
    877 SCIP_CALL( SCIPlpiExactSetIntpar(*lpi, SCIP_LPPAR_PRICING, (int)(*lpi)->pricing) );
    878
    879 {
    880 SPxOut::Verbosity verbosity = (*lpi)->spx->spxout.getVerbosity();
    881 (*lpi)->spx->spxout.setVerbosity((SPxOut::Verbosity)((*lpi)->spx->getLpInfo() ? SOPLEX_VERBLEVEL : 0));
    882 (*lpi)->spx->printVersion();
    883 (*lpi)->spx->spxout.setVerbosity(verbosity);
    884 }
    885
    886 return SCIP_OKAY;
    887}
    888
    889/** deletes an LP problem object */
    891 SCIP_LPIEXACT** lpi /**< pointer to an LP interface structure */
    892 )
    893{
    894 assert(lpi != NULL);
    895 assert(*lpi != NULL);
    896 assert((*lpi)->spx != NULL);
    897
    898 /* free LP using destructor and free memory via blockmemshell */
    899 (*lpi)->spx->~SPxexSCIP();
    900 BMSfreeMemory(&((*lpi)->spx));
    901
    902 /* free memory */
    903 BMSfreeMemoryArrayNull(&(*lpi)->cstat);
    904 BMSfreeMemoryArrayNull(&(*lpi)->rstat);
    905 BMSfreeMemory(lpi);
    906
    907 return SCIP_OKAY;
    908}
    909
    910/**@} */
    911
    912
    913/*
    914 * Modification Methods
    915 */
    916
    917/**@name Modification Methods */
    918/**@{ */
    919
    920/** copies LP data with column matrix into LP solver */
    922 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    923 SCIP_OBJSEN objsen, /**< objective sense */
    924 int ncols, /**< number of columns */
    925 SCIP_RATIONAL** obj, /**< objective function values of columns */
    926 SCIP_RATIONAL** lb, /**< lower bounds of columns */
    927 SCIP_RATIONAL** ub, /**< upper bounds of columns */
    928 char** colnames, /**< column names, or NULL */
    929 int nrows, /**< number of rows */
    930 SCIP_RATIONAL** lhs, /**< left hand sides of rows */
    931 SCIP_RATIONAL** rhs, /**< right hand sides of rows */
    932 char** /*rownames*/, /**< row names, or NULL */
    933 int nnonz, /**< number of nonzero elements in the constraint matrix */
    934 int* beg, /**< start index of each column in ind- and val-array */
    935 int* ind, /**< row indices of constraint matrix entries */
    936 SCIP_RATIONAL** val /**< values of constraint matrix entries */
    937 )
    938{
    939#ifndef NDEBUG
    940 {
    941 int j;
    942 for( j = 0; j < nnonz; j++ )
    943 {
    944 assert(val[j] != NULL);
    945 assert(!SCIPrationalIsZero(val[j]));
    946 }
    947 }
    948#endif
    949
    950 SCIPdebugMessage("calling SCIPlpiExactLoadColLP()\n");
    951
    952 assert(lpi != NULL);
    953 assert(lpi->spx != NULL);
    954 assert(lhs != NULL);
    955 assert(rhs != NULL);
    956 assert(obj != NULL);
    957 assert(lb != NULL);
    958 assert(ub != NULL);
    959 assert(beg != NULL);
    960 assert(ind != NULL);
    961 assert(val != NULL);
    962
    964 assert(lpi->spx->preStrongbranchingBasisFreed());
    965
    966 try
    967 {
    968 SPxexSCIP* spx = lpi->spx;
    969 LPRowSetRational rows(nrows);
    970 DSVectorRational emptyVector(0);
    971 int i;
    972
    973 spx->clearLPRational();
    974
    975 /* set objective sense */
    976 (void) spx->setIntParam(SoPlex::OBJSENSE, (int) (objsen == SCIP_OBJSEN_MINIMIZE ? SoPlex::OBJSENSE_MINIMIZE : SoPlex::OBJSENSE_MAXIMIZE));
    977
    978 /* create empty rows with given sides */
    979 for( i = 0; i < nrows; ++i )
    980 {
    981 soplex::Rational spxlhs;
    982 soplex::Rational spxrhs;
    983 SpxRSetRat(lpi, spxlhs, lhs[i]);
    984 SpxRSetRat(lpi, spxlhs, rhs[i]);
    985 rows.add(spxlhs, emptyVector, spxrhs);
    986 }
    987 spx->addRowsRational(rows);
    988
    989 /* create column vectors with coefficients and bounds */
    990 SCIP_CALL( SCIPlpiExactAddCols(lpi, ncols, obj, lb, ub, colnames, nnonz, beg, ind, val) );
    991 //spx->syncLPReal();
    992 }
    993#ifndef NDEBUG
    994 catch( const SPxException& x )
    995 {
    996 std::string s = x.what();
    997 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    998#else
    999 catch( const SPxException& )
    1000 {
    1001#endif
    1002 return SCIP_LPERROR;
    1003 }
    1004
    1005 return SCIP_OKAY;
    1006}
    1007
    1008/** adds columns to the LP */
    1010 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1011 int ncols, /**< number of columns to be added */
    1012 SCIP_RATIONAL** obj, /**< objective function values of new columns */
    1013 SCIP_RATIONAL** lb, /**< lower bounds of new columns */
    1014 SCIP_RATIONAL** ub, /**< upper bounds of new columns */
    1015 char** colnames, /**< column names, or NULL */
    1016 int nnonz, /**< number of nonzero elements to be added to the constraint matrix */
    1017 int* beg, /**< start index of each column in ind- and val-array, or NULL if nnonz == 0 */
    1018 int* ind, /**< row indices of constraint matrix entries, or NULL if nnonz == 0 */
    1019 SCIP_RATIONAL** val /**< values of constraint matrix entries, or NULL if nnonz == 0 */
    1020 )
    1021{ /*lint --e{715}*/
    1022 SCIPdebugMessage("calling SCIPlpiAddCols()\n");
    1023
    1024 assert(lpi != NULL);
    1025 assert(lpi->spx != NULL);
    1026 assert(obj != NULL);
    1027 assert(lb != NULL);
    1028 assert(ub != NULL);
    1029 assert(nnonz == 0 || beg != NULL);
    1030 assert(nnonz == 0 || ind != NULL);
    1031 assert(nnonz == 0 || val != NULL);
    1032 assert(nnonz >= 0);
    1033 assert(ncols >= 0);
    1034
    1035 invalidateSolution(lpi);
    1036
    1037 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1038
    1039#ifndef NDEBUG
    1040 if ( nnonz > 0 )
    1041 {
    1042 /* perform check that no new rows are added - this is likely to be a mistake */
    1043 int nrows = lpi->spx->numRowsRational();
    1044 for (int j = 0; j < nnonz; ++j)
    1045 {
    1046 assert( 0 <= ind[j] && ind[j] < nrows );
    1047 assert( val[j] != NULL );
    1048 assert( !SCIPrationalIsZero(val[j]) );
    1049 }
    1050 }
    1051#endif
    1052
    1053 SPxexSCIP* spx = lpi->spx;
    1054 try
    1055 {
    1056 LPColSetRational cols(ncols);
    1057 DSVectorRational colVector(ncols);
    1058 int start;
    1059 int last;
    1060 int i;
    1061
    1062 /* create column vectors with coefficients and bounds */
    1063 for( i = 0; i < ncols; ++i )
    1064 {
    1065 int j;
    1066 soplex::Rational spxlb;
    1067 soplex::Rational spxub;
    1068 soplex::Rational spxobj;
    1069 std::vector<soplex::Rational> spxval;
    1070
    1071 SpxRSetRat(lpi, spxlb, lb[i]);
    1072 SpxRSetRat(lpi, spxub, ub[i]);
    1073 SpxRSetRat(lpi, spxobj, obj[i]);
    1074
    1075 colVector.clear();
    1076 if( nnonz > 0 )
    1077 {
    1078 soplex::Rational tmp;
    1079
    1080 start = beg[i];
    1081 last = (i == ncols-1 ? nnonz : beg[i+1]);
    1082 spxval.reserve((size_t) (last - start));
    1083 for( j = start; j < last; ++j )
    1084 {
    1085 spxval.push_back(tmp);
    1086 SpxRSetRat(lpi, spxval[(size_t) (j - start)], val[j]);
    1087 }
    1088 colVector.add(last - start, ind + start, spxval.data());
    1089 }
    1090 cols.add(spxobj, spxlb, colVector, spxub);
    1091 }
    1092 spx->addColsRational(cols);
    1093 }
    1094#ifndef NDEBUG
    1095 catch( const SPxException& x )
    1096 {
    1097 std::string s = x.what();
    1098 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    1099#else
    1100 catch( const SPxException& )
    1101 {
    1102#endif
    1103 return SCIP_LPERROR;
    1104 }
    1105
    1106 return SCIP_OKAY;
    1107}
    1108
    1109/** deletes all columns in the given range from LP */
    1111 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1112 int firstcol, /**< first column to be deleted */
    1113 int lastcol /**< last column to be deleted */
    1114 )
    1115{
    1116 SCIPdebugMessage("calling SCIPlpiDelCols()\n");
    1117
    1118 assert(lpi != NULL);
    1119 assert(lpi->spx != NULL);
    1120 assert(0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsRational());
    1121
    1122 invalidateSolution(lpi);
    1123
    1124 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1125
    1126 SOPLEX_TRY( lpi->messagehdlr, lpi->spx->removeColRangeRational(firstcol, lastcol) );
    1127
    1128 return SCIP_OKAY;
    1129}
    1130
    1131/** deletes columns from SCIP_LP; the new position of a column must not be greater that its old position */
    1133 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1134 int* dstat /**< deletion status of columns
    1135 * input: 1 if column should be deleted, 0 if not
    1136 * output: new position of column, -1 if column was deleted */
    1137 )
    1138{
    1139 int ncols;
    1140 int i;
    1141
    1142 SCIPdebugMessage("calling SCIPlpiDelColset()\n");
    1143
    1144 assert(lpi != NULL);
    1145 assert(lpi->spx != NULL);
    1146 assert(dstat != NULL);
    1147
    1148 invalidateSolution(lpi);
    1149
    1150 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1151
    1152 ncols = lpi->spx->numColsRational();
    1153
    1154 /* SoPlex removeCols() method deletes the columns with dstat[i] < 0, so we have to negate the values */
    1155 for( i = 0; i < ncols; ++i )
    1156 dstat[i] *= -1;
    1157
    1158 SOPLEX_TRY( lpi->messagehdlr, lpi->spx->removeColsRational(dstat) );
    1159
    1160 return SCIP_OKAY;
    1161}
    1162
    1163/** adds rows to the LP */
    1165 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1166 int nrows, /**< number of rows to be added */
    1167 SCIP_RATIONAL** lhs, /**< left hand sides of new rows */
    1168 SCIP_RATIONAL** rhs, /**< right hand sides of new rows */
    1169 char** rownames, /**< row names, or NULL */
    1170 int nnonz, /**< number of nonzero elements to be added to the constraint matrix */
    1171 int* beg, /**< start index of each row in ind- and val-array, or NULL if nnonz == 0 */
    1172 int* ind, /**< column indices of constraint matrix entries, or NULL if nnonz == 0 */
    1173 SCIP_RATIONAL** val /**< values of constraint matrix entries, or NULL if nnonz == 0 */
    1174 )
    1175{ /*lint --e{715}*/
    1176 SCIPdebugMessage("calling SCIPlpiAddRows()\n");
    1177
    1178 assert(lpi != NULL);
    1179 assert(lpi->spx != NULL);
    1180 assert(lhs != NULL);
    1181 assert(rhs != NULL);
    1182 assert(nnonz == 0 || beg != NULL);
    1183 assert(nnonz == 0 || ind != NULL);
    1184 assert(nnonz == 0 || val != NULL);
    1185
    1186 invalidateSolution(lpi);
    1187
    1188 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1189
    1190#ifndef NDEBUG
    1191 if ( nnonz > 0 )
    1192 {
    1193 /* perform check that no new columns are added - this is likely to be a mistake */
    1194 int ncols = lpi->spx->numColsRational();
    1195 for (int j = 0; j < nnonz; ++j)
    1196 {
    1197 assert( !SCIPrationalIsZero(val[j]) );
    1198 assert( 0 <= ind[j] && ind[j] < ncols );
    1199 }
    1200 }
    1201#endif
    1202
    1203 try
    1204 {
    1205 SPxexSCIP* spx = lpi->spx;
    1206 LPRowSetRational rows(nrows);
    1207 DSVectorRational rowVector;
    1208 int start;
    1209 int last;
    1210 int i;
    1211
    1212 /* create row vectors with given sides */
    1213 for( i = 0; i < nrows; ++i )
    1214 {
    1215 soplex::Rational spxlhs;
    1216 soplex::Rational spxrhs;
    1217 std::vector<soplex::Rational> spxval;
    1218
    1219 SpxRSetRat(lpi, spxlhs, lhs[i]);
    1220 SpxRSetRat(lpi, spxrhs, rhs[i]);
    1221
    1222 rowVector.clear();
    1223 if( nnonz > 0 )
    1224 {
    1225 start = beg[i];
    1226 last = (i == nrows-1 ? nnonz : beg[i+1]);
    1227 spxval.reserve((size_t) (last - start));
    1228
    1229 for( int j = start; j < last; ++j )
    1230 {
    1231 soplex::Rational tmp;
    1232 spxval.push_back(tmp);
    1233 SpxRSetRat(lpi, spxval[(size_t) (j - start)], val[j]);
    1234 }
    1235 rowVector.add(last - start, ind + start, spxval.data());
    1236 }
    1237 rows.add(spxlhs, rowVector, spxrhs);
    1238 }
    1239 spx->addRowsRational(rows);
    1240 }
    1241#ifndef NDEBUG
    1242 catch( const SPxException& x )
    1243 {
    1244 std::string s = x.what();
    1245 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    1246#else
    1247 catch( const SPxException& )
    1248 {
    1249#endif
    1250 return SCIP_LPERROR;
    1251 }
    1252
    1253 return SCIP_OKAY;
    1254}
    1255
    1256/** deletes all rows in the given range from LP */
    1258 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1259 int firstrow, /**< first row to be deleted */
    1260 int lastrow /**< last row to be deleted */
    1261 )
    1262{
    1263 SCIPdebugMessage("calling SCIPlpiDelRows()\n");
    1264
    1265 assert(lpi != NULL);
    1266 assert(lpi->spx != NULL);
    1267 assert(0 <= firstrow && firstrow <= lastrow && lastrow < lpi->spx->numRowsRational());
    1268
    1269 invalidateSolution(lpi);
    1270
    1271 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1272
    1273 SOPLEX_TRY( lpi->messagehdlr, lpi->spx->removeRowRangeRational(firstrow, lastrow) );
    1274
    1275 return SCIP_OKAY;
    1276}
    1277
    1278/** deletes rows from SCIP_LP; the new position of a row must not be greater that its old position */
    1280 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1281 int* dstat /**< deletion status of rows
    1282 * input: 1 if row should be deleted, 0 if not
    1283 * output: new position of row, -1 if row was deleted */
    1284 )
    1285{
    1286 int nrows;
    1287 int i;
    1288
    1289 SCIPdebugMessage("calling SCIPlpiDelRowset()\n");
    1290
    1291 assert(lpi != NULL);
    1292 assert(lpi->spx != NULL);
    1293
    1294 invalidateSolution(lpi);
    1295
    1296 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1297
    1298 nrows = lpi->spx->numRowsRational();
    1299
    1300 /* SoPlex removeRows() method deletes the rows with dstat[i] < 0, so we have to negate the values */
    1301 for( i = 0; i < nrows; ++i )
    1302 dstat[i] *= -1;
    1303
    1304 SOPLEX_TRY( lpi->messagehdlr, lpi->spx->removeRowsRational(dstat) );
    1305
    1306 return SCIP_OKAY;
    1307}
    1308
    1309/** clears the whole LP */
    1311 SCIP_LPIEXACT* lpi /**< LP interface structure */
    1312 )
    1313{
    1314 SCIPdebugMessage("calling SCIPlpiExactClear()\n");
    1315
    1316 assert(lpi != NULL);
    1317 assert(lpi->spx != NULL);
    1318
    1319 invalidateSolution(lpi);
    1320
    1321 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1322 SOPLEX_TRY( lpi->messagehdlr, lpi->spx->clearLPRational() );
    1323
    1324 return SCIP_OKAY;
    1325}
    1326
    1327/** changes lower and upper bounds of columns */
    1329 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1330 int ncols, /**< number of columns to change bounds for */
    1331 int* ind, /**< column indices or NULL if ncols is zero */
    1332 SCIP_RATIONAL** lb, /**< values for the new lower bounds or NULL if ncols is zero */
    1333 SCIP_RATIONAL** ub /**< values for the new upper bounds or NULL if ncols is zero */
    1334 )
    1335{
    1336 int i;
    1337
    1338 SCIPdebugMessage("calling SCIPlpiChgBounds()\n");
    1339
    1340 assert(lpi != NULL);
    1341 assert(lpi->spx != NULL);
    1342 assert(ncols == 0 || (ind != NULL && lb != NULL && ub != NULL));
    1343 if( ncols <= 0 )
    1344 return SCIP_OKAY;
    1345
    1346 invalidateSolution(lpi);
    1347
    1348 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1349
    1350 try
    1351 {
    1352 soplex::Rational spxlb;
    1353 soplex::Rational spxub;
    1354
    1355 for( i = 0; i < ncols; ++i )
    1356 {
    1357 assert(0 <= ind[i] && ind[i] < lpi->spx->numColsRational());
    1358 assert(lb[i] != NULL && ub[i] != NULL);
    1359
    1360 if( SCIPrationalIsInfinity(lb[i]) )
    1361 {
    1362 SCIPerrorMessage("LP Error: fixing lower bound for variable %d to infinity.\n", ind[i]);
    1363 return SCIP_LPERROR;
    1364 }
    1365 if( SCIPrationalIsNegInfinity(ub[i]) )
    1366 {
    1367 SCIPerrorMessage("LP Error: fixing upper bound for variable %d to -infinity.\n", ind[i]);
    1368 return SCIP_LPERROR;
    1369 }
    1370
    1371 SpxRSetRat(lpi, spxlb, lb[i]);
    1372 SpxRSetRat(lpi, spxub, ub[i]);
    1373
    1374 lpi->spx->changeBoundsRational(ind[i], spxlb, spxub);
    1375 assert(lpi->spx->lowerRational(ind[i]) <= lpi->spx->upperRational(ind[i]));
    1376 }
    1377 }
    1378#ifndef NDEBUG
    1379 catch( const SPxException& x )
    1380 {
    1381 std::string s = x.what();
    1382 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    1383#else
    1384 catch( const SPxException& )
    1385 {
    1386#endif
    1387 return SCIP_LPERROR;
    1388 }
    1389
    1390 return SCIP_OKAY;
    1391}
    1392
    1393/** changes left and right hand sides of rows */
    1395 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1396 int nrows, /**< number of rows to change sides for */
    1397 int* ind, /**< row indices */
    1398 SCIP_RATIONAL** lhs, /**< new values for left hand sides */
    1399 SCIP_RATIONAL** rhs /**< new values for right hand sides */
    1400 )
    1401{
    1402 int i;
    1403
    1404 SCIPdebugMessage("calling SCIPlpiChgSides()\n");
    1405
    1406 assert(lpi != NULL);
    1407 assert(lpi->spx != NULL);
    1408 assert(ind != NULL);
    1409 assert(lhs != NULL);
    1410 assert(rhs != NULL);
    1411 if( nrows <= 0 )
    1412 return SCIP_OKAY;
    1413
    1414 invalidateSolution(lpi);
    1415
    1416 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1417
    1418 try
    1419 {
    1420 soplex::Rational spxlhs;
    1421 soplex::Rational spxrhs;
    1422
    1423 for( i = 0; i < nrows; ++i )
    1424 {
    1425 assert(0 <= ind[i] && ind[i] < lpi->spx->numRowsRational());
    1426
    1427 SpxRSetRat(lpi, spxlhs, lhs[i]);
    1428 SpxRSetRat(lpi, spxrhs, rhs[i]);
    1429
    1430 lpi->spx->changeRangeRational(ind[i], spxlhs, spxrhs);
    1431 assert(lpi->spx->lhsRational(ind[i]) <= lpi->spx->rhsRational(ind[i]));
    1432 }
    1433 }
    1434#ifndef NDEBUG
    1435 catch( const SPxException& x )
    1436 {
    1437 std::string s = x.what();
    1438 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    1439#else
    1440 catch( const SPxException& )
    1441 {
    1442#endif
    1443 return SCIP_LPERROR;
    1444 }
    1445
    1446 return SCIP_OKAY;
    1447}
    1448
    1449/** changes a single coefficient */
    1451 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1452 int row, /**< row number of coefficient to change */
    1453 int col, /**< column number of coefficient to change */
    1454 SCIP_RATIONAL* newval /**< new value of coefficient */
    1455 )
    1456{
    1457 soplex::Rational spxval;
    1458
    1459 SCIPdebugMessage("calling SCIPlpiChgCoef()\n");
    1460
    1461 assert(lpi != NULL);
    1462 assert(lpi->spx != NULL);
    1463 assert(0 <= row && row < lpi->spx->numRowsRational());
    1464 assert(0 <= col && col < lpi->spx->numColsRational());
    1465
    1466 invalidateSolution(lpi);
    1467
    1468 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1469
    1470 SpxRSetRat(lpi, spxval, newval);
    1471 SOPLEX_TRY( lpi->messagehdlr, lpi->spx->changeElementRational(row, col, spxval) );
    1472
    1473 return SCIP_OKAY;
    1474}
    1475
    1476/** changes the objective sense */
    1478 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1479 SCIP_OBJSEN objsen /**< new objective sense */
    1480 )
    1481{
    1482 SCIPdebugMessage("calling SCIPlpiChgObjsen()\n");
    1483
    1484 assert(lpi != NULL);
    1485 assert(lpi->spx != NULL);
    1486
    1487 invalidateSolution(lpi);
    1488
    1489 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1490
    1491 SOPLEX_TRY( lpi->messagehdlr, (void) lpi->spx->setIntParam(SoPlex::OBJSENSE, (int) (objsen == SCIP_OBJSEN_MINIMIZE ? SoPlex::OBJSENSE_MINIMIZE : SoPlex::OBJSENSE_MAXIMIZE) ) );
    1492
    1493 return SCIP_OKAY;
    1494}
    1495
    1496/** changes objective values of columns in the LP */
    1498 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1499 int ncols, /**< number of columns to change objective value for */
    1500 int* ind, /**< column indices to change objective value for */
    1501 SCIP_RATIONAL** obj /**< new objective values for columns */
    1502 )
    1503{
    1504 int i;
    1505
    1506 SCIPdebugMessage("calling SCIPlpiChgObj()\n");
    1507
    1508 assert(lpi != NULL);
    1509 assert(lpi->spx != NULL);
    1510 assert(ind != NULL);
    1511 assert(obj != NULL);
    1512
    1513 invalidateSolution(lpi);
    1514
    1515 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1516
    1517 try
    1518 {
    1519 soplex::Rational spxobj;
    1520
    1521 for( i = 0; i < ncols; ++i )
    1522 {
    1523 assert(obj[i] != NULL);
    1524 assert(0 <= ind[i] && ind[i] < lpi->spx->numColsRational());
    1525
    1526 SpxRSetRat(lpi, spxobj, obj[i]);
    1527 lpi->spx->changeObjRational(ind[i], spxobj);
    1528 }
    1529 }
    1530#ifndef NDEBUG
    1531 catch( const SPxException& x )
    1532 {
    1533 std::string s = x.what();
    1534 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    1535#else
    1536 catch( const SPxException& )
    1537 {
    1538#endif
    1539 return SCIP_LPERROR;
    1540 }
    1541
    1542 return SCIP_OKAY;
    1543}
    1544
    1545/**@} */
    1546
    1547
    1548/*
    1549 * Data Accessing Methods
    1550 */
    1551
    1552/**@name Data Accessing Methods */
    1553/**@{ */
    1554
    1555/** gets the number of rows in the LP */
    1557 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1558 int* nrows /**< pointer to store the number of rows */
    1559 )
    1560{
    1561 SCIPdebugMessage("calling SCIPlpiExactGetNRows()\n");
    1562
    1563 assert(lpi != NULL);
    1564 assert(lpi->spx != NULL);
    1565 assert(nrows != NULL);
    1566
    1567 *nrows = lpi->spx->numRowsRational();
    1568
    1569 return SCIP_OKAY;
    1570}
    1571
    1572/** gets the number of columns in the LP */
    1574 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1575 int* ncols /**< pointer to store the number of cols */
    1576 )
    1577{
    1578 SCIPdebugMessage("calling SCIPlpiExactGetNCols()\n");
    1579
    1580 assert(lpi != NULL);
    1581 assert(lpi->spx != NULL);
    1582 assert(ncols != NULL);
    1583
    1584 *ncols = lpi->spx->numColsRational();
    1585
    1586 return SCIP_OKAY;
    1587}
    1588
    1589/** gets the number of nonzero elements in the LP constraint matrix */
    1591 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1592 int* nnonz /**< pointer to store the number of nonzeros */
    1593 )
    1594{
    1595 int i;
    1596
    1597 SCIPdebugMessage("calling SCIPlpiExactGetNNonz()\n");
    1598
    1599 assert(lpi != NULL);
    1600 assert(lpi->spx != NULL);
    1601 assert(nnonz != NULL);
    1602
    1603 /* SoPlex has no direct method to return the number of nonzeros, so we have to count them manually */
    1604 *nnonz = 0;
    1605 if( lpi->spx->numRowsRational() < lpi->spx->numColsRational() )
    1606 {
    1607 for( i = 0; i < lpi->spx->numRowsRational(); ++i )
    1608 (*nnonz) += lpi->spx->rowVectorRational(i).size();
    1609 }
    1610 else
    1611 {
    1612 for( i = 0; i < lpi->spx->numColsRational(); ++i )
    1613 (*nnonz) += lpi->spx->colVectorRational(i).size();
    1614 }
    1615
    1616 return SCIP_OKAY;
    1617}
    1618
    1619/** gets columns from LP problem object; the arrays have to be large enough to store all values
    1620 * Either both, lb and ub, have to be NULL, or both have to be non-NULL,
    1621 * either nnonz, beg, ind, and val have to be NULL, or all of them have to be non-NULL.
    1622 */
    1624 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1625 int firstcol, /**< first column to get from LP */
    1626 int lastcol, /**< last column to get from LP */
    1627 SCIP_RATIONAL** lb, /**< buffer to store the lower bound vector, or NULL */
    1628 SCIP_RATIONAL** ub, /**< buffer to store the upper bound vector, or NULL */
    1629 int* nnonz, /**< pointer to store the number of nonzero elements returned, or NULL */
    1630 int* beg, /**< buffer to store start index of each column in ind- and val-array, or NULL */
    1631 int* ind, /**< buffer to store row indices of constraint matrix entries, or NULL */
    1632 SCIP_RATIONAL** val /**< buffer to store values of constraint matrix entries, or NULL */
    1633 )
    1634{
    1635 int i;
    1636 int j;
    1637
    1638 SCIPdebugMessage("calling SCIPlpiExactGetCols()\n");
    1639
    1640 assert(lpi != NULL);
    1641 assert(lpi->spx != NULL);
    1642 assert(0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsRational());
    1643 assert((lb != NULL && ub != NULL) || (lb == NULL && ub == NULL));
    1644 assert((nnonz != NULL && beg != NULL && ind != NULL && val != NULL) || (nnonz == NULL && beg == NULL && ind == NULL && val == NULL));
    1645
    1646 if( lb != NULL )
    1647 {
    1648 if( lpi->spx->boolParam(SoPlex::PERSISTENTSCALING) )
    1649 {
    1650 const VectorRational& lbvec = lpi->spx->lowerRational();
    1651 const VectorRational& ubvec = lpi->spx->upperRational();
    1652
    1653 for( i = firstcol; i <= lastcol; ++i )
    1654 {
    1655 RsetSpxR(lpi, lb[i-firstcol], lbvec[i]);
    1656 RsetSpxR(lpi, ub[i-firstcol], ubvec[i]);
    1657 }
    1658 }
    1659 else
    1660 {
    1661 const VectorRational& lbvec = lpi->spx->lowerRational();
    1662 const VectorRational& ubvec = lpi->spx->upperRational();
    1663 for( i = firstcol; i <= lastcol; ++i )
    1664 {
    1665 RsetSpxR(lpi, lb[i-firstcol], lbvec[i]);
    1666 RsetSpxR(lpi, ub[i-firstcol], ubvec[i]);
    1667 }
    1668 }
    1669 }
    1670
    1671 if( nnonz != NULL )
    1672 {
    1673 *nnonz = 0;
    1674 for( i = firstcol; i <= lastcol; ++i )
    1675 {
    1676 const SVectorRational& cvec = lpi->spx->colVectorRational(i);
    1677 beg[i-firstcol] = *nnonz;
    1678
    1679 for( j = 0; j < cvec.size(); ++j )
    1680 {
    1681 ind[*nnonz] = cvec.index(j);
    1682 RsetSpxR(lpi, val[*nnonz], cvec.value(j));
    1683 (*nnonz)++;
    1684 }
    1685 }
    1686 }
    1687
    1688 return SCIP_OKAY;
    1689}
    1690
    1691/** gets rows from LP problem object; the arrays have to be large enough to store all values.
    1692 * Either both, lhs and rhs, have to be NULL, or both have to be non-NULL,
    1693 * either nnonz, beg, ind, and val have to be NULL, or all of them have to be non-NULL.
    1694 */
    1696 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1697 int firstrow, /**< first row to get from LP */
    1698 int lastrow, /**< last row to get from LP */
    1699 SCIP_RATIONAL** lhs, /**< buffer to store left hand side vector, or NULL */
    1700 SCIP_RATIONAL** rhs, /**< buffer to store right hand side vector, or NULL */
    1701 int* nnonz, /**< pointer to store the number of nonzero elements returned, or NULL */
    1702 int* beg, /**< buffer to store start index of each row in ind- and val-array, or NULL */
    1703 int* ind, /**< buffer to store column indices of constraint matrix entries, or NULL */
    1704 SCIP_RATIONAL** val /**< buffer to store values of constraint matrix entries, or NULL */
    1705 )
    1706{
    1707 int i;
    1708 int j;
    1709
    1710 SCIPdebugMessage("calling SCIPlpiExactGetRows()\n");
    1711
    1712 assert(lpi != NULL);
    1713 assert(lpi->spx != NULL);
    1714 assert(0 <= firstrow && firstrow <= lastrow && lastrow < lpi->spx->numRowsRational());
    1715 assert((lhs != NULL && rhs != NULL) || (lhs == NULL && rhs == NULL));
    1716 assert((nnonz != NULL && beg != NULL && ind != NULL && val != NULL) || (nnonz == NULL && beg == NULL && ind == NULL && val == NULL));
    1717
    1718 if( lhs != NULL )
    1719 {
    1720 if( lpi->spx->boolParam(SoPlex::PERSISTENTSCALING) && FALSE )
    1721 {
    1722 }
    1723 else
    1724 {
    1725 const VectorRational& lhsvec = lpi->spx->lhsRational();
    1726 const VectorRational& rhsvec = lpi->spx->rhsRational();
    1727 for( i = firstrow; i <= lastrow; ++i )
    1728 {
    1729 RsetSpxR(lpi, lhs[i-firstrow], lhsvec[i]);
    1730 RsetSpxR(lpi, rhs[i-firstrow], rhsvec[i]);
    1731 }
    1732 }
    1733 }
    1734
    1735 if( nnonz != NULL )
    1736 {
    1737 *nnonz = 0;
    1738 for( i = firstrow; i <= lastrow; ++i )
    1739 {
    1740 beg[i-firstrow] = *nnonz;
    1741
    1742 if( lpi->spx->boolParam(SoPlex::PERSISTENTSCALING) )
    1743 {
    1744 }
    1745 else
    1746 {
    1747 const SVectorRational& rvec = lpi->spx->rowVectorRational(i);
    1748 for( j = 0; j < rvec.size(); ++j )
    1749 {
    1750 ind[*nnonz] = rvec.index(j);
    1751 RsetSpxR(lpi, val[*nnonz], rvec.value(j));
    1752 (*nnonz)++;
    1753 }
    1754 }
    1755 }
    1756 }
    1757
    1758 return SCIP_OKAY;
    1759}
    1760
    1761/** gets column names */
    1763 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1764 int firstcol, /**< first column to get name from LP */
    1765 int lastcol, /**< last column to get name from LP */
    1766 char** colnames, /**< pointers to column names (of size at least lastcol-firstcol+1) or NULL if namestoragesize is zero */
    1767 char* namestorage, /**< storage for col names or NULL if namestoragesize is zero */
    1768 int namestoragesize, /**< size of namestorage (if 0, storageleft returns the storage needed) */
    1769 int* storageleft /**< amount of storage left (if < 0 the namestorage was not big enough) or NULL if namestoragesize is zero */
    1770 )
    1771{
    1772 assert( lpi != NULL );
    1773 assert( lpi->spx != NULL );
    1774 assert( colnames != NULL || namestoragesize == 0 );
    1775 assert( namestorage != NULL || namestoragesize == 0 );
    1776 assert( namestoragesize >= 0 );
    1777 assert( storageleft != NULL );
    1778 assert( 0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsRational() );
    1779
    1780 SCIPdebugMessage("getting column names %d to %d\n", firstcol, lastcol);
    1781
    1782 SCIPerrorMessage("SCIPlpiExactGetColNames() not implemented\n");
    1783
    1784 return SCIP_NOTIMPLEMENTED;
    1785}
    1786
    1787/** gets row names */
    1789 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1790 int firstrow, /**< first row to get name from LP */
    1791 int lastrow, /**< last row to get name from LP */
    1792 char** rownames, /**< pointers to row names (of size at least lastrow-firstrow+1) or NULL if namestoragesize is zero */
    1793 char* namestorage, /**< storage for row names or NULL if namestoragesize is zero */
    1794 int namestoragesize, /**< size of namestorage (if 0, -storageleft returns the storage needed) */
    1795 int* storageleft /**< amount of storage left (if < 0 the namestorage was not big enough) or NULL if namestoragesize is zero */
    1796 )
    1797{
    1798 assert( lpi != NULL );
    1799 assert( lpi->spx != NULL );
    1800 assert( rownames != NULL || namestoragesize == 0 );
    1801 assert( namestorage != NULL || namestoragesize == 0 );
    1802 assert( namestoragesize >= 0 );
    1803 assert( storageleft != NULL );
    1804 assert( 0 <= firstrow && firstrow <= lastrow && lastrow < lpi->spx->numRowsRational() );
    1805
    1806 SCIPdebugMessage("getting row names %d to %d\n", firstrow, lastrow);
    1807
    1808 SCIPerrorMessage("SCIPlpiExactGetRowNames() not implemented\n");
    1809
    1810 return SCIP_NOTIMPLEMENTED;
    1811}
    1812
    1813/** gets objective sense of the LP */
    1815 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1816 SCIP_OBJSEN* objsen /**< pointer to store objective sense */
    1817 )
    1818{
    1819 SCIPdebugMessage("calling SCIPlpiExactGetObjsen()\n");
    1820
    1821 assert(lpi != NULL);
    1822 assert(lpi->spx != NULL);
    1823 assert(objsen != NULL);
    1824
    1825 *objsen = (lpi->spx->intParam(SoPlex::OBJSENSE) == SoPlex::OBJSENSE_MINIMIZE) ? SCIP_OBJSEN_MINIMIZE : SCIP_OBJSEN_MAXIMIZE;
    1826
    1827 return SCIP_OKAY;
    1828}
    1829
    1830/** gets objective coefficients from LP problem object */
    1832 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1833 int firstcol, /**< first column to get objective coefficient for */
    1834 int lastcol, /**< last column to get objective coefficient for */
    1835 SCIP_RATIONAL** vals /**< array to store objective coefficients */
    1836 )
    1837{
    1838 int i;
    1839
    1840 SCIPdebugMessage("calling SCIPlpiExactGetObj()\n");
    1841
    1842 assert(lpi != NULL);
    1843 assert(lpi->spx != NULL);
    1844 assert(0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsRational());
    1845 assert(vals != NULL);
    1846
    1847 for( i = firstcol; i <= lastcol; ++i )
    1848 {
    1849 assert(vals[i-firstcol] != NULL);
    1850 RsetSpxR(lpi, vals[i-firstcol], lpi->spx->objRational(i));
    1851 }
    1852
    1853 return SCIP_OKAY;
    1854}
    1855
    1856/** gets current bounds from LP problem object */
    1858 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1859 int firstcol, /**< first column to get objective value for */
    1860 int lastcol, /**< last column to get objective value for */
    1861 SCIP_RATIONAL** lbs, /**< array to store lower bound values, or NULL */
    1862 SCIP_RATIONAL** ubs /**< array to store upper bound values, or NULL */
    1863 )
    1864{
    1865 int i;
    1866
    1867 SCIPdebugMessage("calling SCIPlpiExactGetBounds()\n");
    1868
    1869 assert(lpi != NULL);
    1870 assert(lpi->spx != NULL);
    1871 assert(0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsRational());
    1872
    1873 for( i = firstcol; i <= lastcol; ++i )
    1874 {
    1875 if( lbs != NULL )
    1876 {
    1877 assert(lbs[i-firstcol] != NULL);
    1878 RsetSpxR(lpi, lbs[i-firstcol], lpi->spx->lowerRational(i));
    1879 }
    1880 if( ubs != NULL )
    1881 {
    1882 assert(ubs[i-firstcol] != NULL);
    1883 RsetSpxR(lpi, ubs[i-firstcol], lpi->spx->upperRational(i));
    1884 }
    1885 }
    1886
    1887 return SCIP_OKAY;
    1888}
    1889
    1890/** gets current row sides from LP problem object */
    1892 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1893 int firstrow, /**< first row to get sides for */
    1894 int lastrow, /**< last row to get sides for */
    1895 SCIP_RATIONAL** lhss, /**< array to store left hand side values, or NULL */
    1896 SCIP_RATIONAL** rhss /**< array to store right hand side values, or NULL */
    1897 )
    1898{
    1899 int i;
    1900
    1901 SCIPdebugMessage("calling SCIPlpiExactGetSides()\n");
    1902
    1903 assert(lpi != NULL);
    1904 assert(lpi->spx != NULL);
    1905 assert(0 <= firstrow && firstrow <= lastrow && lastrow < lpi->spx->numRowsRational());
    1906
    1907 for( i = firstrow; i <= lastrow; ++i )
    1908 {
    1909 if( lhss != NULL )
    1910 {
    1911 assert(lhss[i-firstrow] != NULL);
    1912 RsetSpxR(lpi, lhss[i-firstrow], lpi->spx->lhsRational(i));
    1913 }
    1914 if( rhss != NULL )
    1915 {
    1916 assert(rhss[i-firstrow] != NULL);
    1917 RsetSpxR(lpi, rhss[i-firstrow], lpi->spx->rhsReal(i));
    1918 }
    1919 }
    1920
    1921 return SCIP_OKAY;
    1922}
    1923
    1924/** gets a single coefficient */
    1926 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    1927 int row, /**< row number of coefficient */
    1928 int col, /**< column number of coefficient */
    1929 SCIP_RATIONAL* val /**< pointer to store the value of the coefficient */
    1930 )
    1931{
    1932 SCIPdebugMessage("calling SCIPlpiExactGetCoef()\n");
    1933
    1934 assert(lpi != NULL);
    1935 assert(lpi->spx != NULL);
    1936 assert(0 <= col && col < lpi->spx->numColsRational());
    1937 assert(0 <= row && row < lpi->spx->numRowsRational());
    1938 assert(val != NULL);
    1939
    1940 RsetSpxR(lpi, val, lpi->spx->colVectorRational(col)[row]);
    1941
    1942 return SCIP_OKAY;
    1943}
    1944
    1945/**@} */
    1946
    1947
    1948/*
    1949 * Solving Methods
    1950 */
    1951
    1952/**@name Solving Methods */
    1953/**@{ */
    1954
    1955/** solves LP -- used for both, primal and dual simplex, because SoPlex doesn't distinct the two cases */
    1956static
    1958 SCIP_LPIEXACT* lpi /**< LP interface structure */
    1959 )
    1960{
    1961 assert( lpi != NULL );
    1962 assert( lpi->spx != NULL );
    1963
    1964 SPxOut::Verbosity verbosity;
    1965 /* store and set verbosity */
    1966 verbosity = lpi->spx->spxout.getVerbosity();
    1967 lpi->spx->spxout.setVerbosity((SPxOut::Verbosity)(lpi->spx->getLpInfo() ? SOPLEX_VERBLEVEL : 0));
    1968
    1969 SCIPdebugMessage("calling exact SoPlex solve(): %d cols, %d rows\n", lpi->spx->numColsRational(), lpi->spx->numRowsRational());
    1970
    1971 invalidateSolution(lpi);
    1972
    1973 assert( lpi->spx->preStrongbranchingBasisFreed() );
    1974
    1975#ifdef WITH_LPSCHECK
    1976 lpi->spx->setDoubleCheck(CHECK_SPXSOLVE);
    1977#endif
    1978
    1979 /* delete starting basis if solving from scratch */
    1980 if( lpi->spx->getFromScratch() )
    1981 {
    1982 try
    1983 {
    1984 lpi->spx->clearBasis();
    1985 }
    1986#ifndef NDEBUG
    1987 catch(const SPxException& x)
    1988 {
    1989 std::string s = x.what();
    1990 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    1991#else
    1992 catch(const SPxException&)
    1993 {
    1994#endif
    1995 assert( lpi->spx->status() != SPxSolver::OPTIMAL );
    1996 return SCIP_LPERROR;
    1997 }
    1998 }
    1999 assert(!lpi->spx->getFromScratch() || lpi->spx->status() == SPxSolver::NO_PROBLEM);
    2000
    2001 SPxSolver::Status status = lpi->spx->doSolve();
    2002 SCIPdebugMessage(" -> SoPlex status: %d, basis status: %d\n", lpi->spx->status(), lpi->spx->basisStatus());
    2003 lpi->solved = TRUE;
    2004
    2005 /* restore verbosity */
    2006 lpi->spx->spxout.setVerbosity(verbosity);
    2007
    2008 switch( status )
    2009 {
    2010 case SPxSolver::ABORT_TIME:
    2011 case SPxSolver::ABORT_ITER:
    2012 case SPxSolver::ABORT_VALUE:
    2013 case SPxSolver::SINGULAR:
    2014 case SPxSolver::REGULAR:
    2015 case SPxSolver::OPTIMAL:
    2016#if SOPLEX_APIVERSION >= 3
    2017 case SPxSolver::OPTIMAL_UNSCALED_VIOLATIONS:
    2018#endif
    2019 case SPxSolver::UNBOUNDED:
    2020 case SPxSolver::INFEASIBLE:
    2021 return SCIP_OKAY;
    2022 case SPxSolver::UNKNOWN:
    2023 default:
    2024 return SCIP_LPERROR;
    2025 } /*lint !e788*/
    2026}
    2027
    2028/** calls primal simplex to solve the LP */
    2030 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2031 )
    2032{
    2033 SCIPdebugMessage("calling SCIPlpiSolvePrimal()\n");
    2034
    2035 assert(lpi != NULL);
    2036 assert(lpi->spx != NULL);
    2037
    2038 (void) lpi->spx->setIntParam(SoPlex::ALGORITHM, (int) SoPlex::ALGORITHM_PRIMAL);
    2039 return spxSolve(lpi);
    2040}
    2041
    2042/** calls dual simplex to solve the LP */
    2044 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2045 )
    2046{
    2047 SCIPdebugMessage("calling SCIPlpiSolveDual()\n");
    2048
    2049 assert(lpi != NULL);
    2050 assert(lpi->spx != NULL);
    2051
    2052 (void) lpi->spx->setIntParam(SoPlex::ALGORITHM, (int) SoPlex::ALGORITHM_DUAL);
    2053 return spxSolve(lpi);
    2054}
    2055
    2056/** calls barrier or interior point algorithm to solve the LP with crossover to simplex basis */
    2058 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2059 SCIP_Bool crossover /**< perform crossover */
    2060 )
    2061{ /*lint --e{715}*/
    2062 assert(lpi != NULL);
    2063 assert(lpi->spx != NULL);
    2064
    2065 SCIPdebugMessage("calling SCIPlpiSolveBarrier()\n");
    2066
    2067 /* Since SoPlex does not support barrier we switch to DUAL */
    2068 return SCIPlpiExactSolveDual(lpi);
    2069}
    2070
    2071/** start strong branching - call before any strongbranching */
    2073 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2074 )
    2075{
    2076 assert(lpi != NULL);
    2077 assert(lpi->spx != NULL);
    2078
    2079 assert( lpi->spx->preStrongbranchingBasisFreed() );
    2081
    2082 return SCIP_OKAY;
    2083}
    2084
    2085/** end strong branching - call after any strongbranching */
    2087 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2088 )
    2089{
    2090 assert(lpi != NULL);
    2091 assert(lpi->spx != NULL);
    2092
    2093 assert( ! lpi->spx->preStrongbranchingBasisFreed() );
    2096
    2097 return SCIP_OKAY;
    2098}
    2099
    2100/**@} */
    2101
    2102
    2103/*
    2104 * Solution Information Methods
    2105 */
    2106
    2107/**@name Solution Information Methods */
    2108/**@{ */
    2109
    2110/** returns whether a solve method was called after the last modification of the LP */
    2112 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2113 )
    2114{
    2115 assert(lpi != NULL);
    2116
    2117 return lpi->solved;
    2118}
    2119
    2120/** gets information about primal and dual feasibility of the current LP solution
    2121 *
    2122 * The feasibility information is with respect to the last solving call and it is only relevant if SCIPlpiWasSolved()
    2123 * returns true. If the LP is changed, this information might be invalidated.
    2124 *
    2125 * Note that @param primalfeasible and @param dualfeasible should only return true if the solver has proved the respective LP to
    2126 * be feasible. Thus, the return values should be equal to the values of SCIPlpiIsPrimalFeasible() and
    2127 * SCIPlpiIsDualFeasible(), respectively. Note that if feasibility cannot be proved, they should return false (even if
    2128 * the problem might actually be feasible).
    2129 */
    2131 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2132 SCIP_Bool* primalfeasible, /**< pointer to store primal feasibility status */
    2133 SCIP_Bool* dualfeasible /**< pointer to store dual feasibility status */
    2134 )
    2135{
    2136 SCIPdebugMessage("calling SCIPlpiExactGetSolFeasibility()\n");
    2137
    2138 assert(lpi != NULL);
    2139 assert(primalfeasible != NULL);
    2140 assert(dualfeasible != NULL);
    2141
    2142 *primalfeasible = SCIPlpiExactIsPrimalFeasible(lpi);
    2143 *dualfeasible = SCIPlpiExactIsDualFeasible(lpi);
    2144
    2145 return SCIP_OKAY;
    2146}
    2147
    2148/** returns TRUE iff LP is proven to have a primal unbounded ray (but not necessary a primal feasible point);
    2149 * this does not necessarily mean, that the solver knows and can return the primal ray
    2150 */
    2152 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2153 )
    2154{
    2155 SCIPdebugMessage("calling SCIPlpiExactExistsPrimalRay()\n");
    2156
    2157 assert(lpi != NULL);
    2158 assert(lpi->spx != NULL);
    2159
    2160 return (lpi->spx->status() == SPxSolver::UNBOUNDED);
    2161}
    2162
    2163/** returns TRUE iff LP is proven to have a primal unbounded ray (but not necessary a primal feasible point),
    2164 * and the solver knows and can return the primal ray
    2165 */
    2167 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2168 )
    2169{
    2170 SCIPdebugMessage("calling SCIPlpiHasPrimalRay()\n");
    2171
    2172 assert(lpi != NULL);
    2173 assert(lpi->spx != NULL);
    2174
    2175 return lpi->spx->hasPrimalRay();
    2176}
    2177
    2178/** returns TRUE iff LP is proven to be primal unbounded */
    2180 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2181 )
    2182{
    2183 SCIPdebugMessage("calling SCIPlpiIsPrimalUnbounded()\n");
    2184
    2185 assert(lpi != NULL);
    2186 assert(lpi->spx != NULL);
    2187
    2188 assert(lpi->spx->status() != SPxSolver::UNBOUNDED || lpi->spx->basisStatus() == SPxBasis::UNBOUNDED);
    2189
    2190 /* if SoPlex returns unbounded, this may only mean that an unbounded ray is available, not necessarily a primal
    2191 * feasible point; hence we have to check the perturbation
    2192 */
    2193 return lpi->spx->status() == SPxSolver::UNBOUNDED;
    2194}
    2195
    2196/** returns TRUE iff LP is proven to be primal infeasible */
    2198 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2199 )
    2200{
    2201 SCIPdebugMessage("calling SCIPlpiIsPrimalInfeasible()\n");
    2202
    2203 assert(lpi != NULL);
    2204 assert(lpi->spx != NULL);
    2205
    2206 return (lpi->spx->status() == SPxSolver::INFEASIBLE);
    2207}
    2208
    2209/** returns TRUE iff LP is proven to be primal feasible */
    2211 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2212 )
    2213{
    2214 SPxBasis::SPxStatus basestatus;
    2215
    2216 SCIPdebugMessage("calling SCIPlpiIsPrimalFeasible()\n");
    2217
    2218 assert(lpi != NULL);
    2219 assert(lpi->spx != NULL);
    2220
    2221 basestatus = lpi->spx->basisStatus();
    2222
    2223 /* note that the solver status may be ABORT_VALUE and the basis status optimal; if we are optimal, isPerturbed() may
    2224 * still return true as long as perturbation plus violation is within tolerances
    2225 */
    2226 assert(basestatus == SPxBasis::OPTIMAL || lpi->spx->status() != SPxSolver::OPTIMAL);
    2227
    2228 return basestatus == SPxBasis::OPTIMAL || basestatus == SPxBasis::PRIMAL;
    2229}
    2230
    2231/** returns TRUE iff LP is proven to have a dual unbounded ray (but not necessary a dual feasible point);
    2232 * this does not necessarily mean, that the solver knows and can return the dual ray
    2233 */
    2235 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2236 )
    2237{
    2238 SCIPdebugMessage("calling SCIPlpiExactExistsDualRay()\n");
    2239
    2240 assert(lpi != NULL);
    2241 assert(lpi->spx != NULL);
    2242
    2243 return (lpi->spx->status() == SPxSolver::INFEASIBLE);
    2244}
    2245
    2246/** returns TRUE iff LP is proven to have a dual unbounded ray (but not necessary a dual feasible point),
    2247 * and the solver knows and can return the dual ray
    2248 */
    2250 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2251 )
    2252{
    2253 SCIPdebugMessage("calling SCIPlpiHasDualRay()\n");
    2254
    2255 assert(lpi != NULL);
    2256 assert(lpi->spx != NULL);
    2257
    2258 return lpi->spx->hasDualFarkas();
    2259}
    2260
    2261/** returns TRUE iff LP is dual unbounded */
    2263 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2264 )
    2265{
    2266 SCIPdebugMessage("calling SCIPlpiIsDualUnbounded()\n");
    2267
    2268 assert(lpi != NULL);
    2269 assert(lpi->spx != NULL);
    2270
    2271 return lpi->spx->status() == SPxSolver::INFEASIBLE && lpi->spx->basisStatus() == SPxBasis::DUAL;
    2272}
    2273
    2274/** returns TRUE iff LP is dual infeasible */
    2276 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2277 )
    2278{
    2279 SCIPdebugMessage("calling SCIPlpiIsDualInfeasible()\n");
    2280
    2281 assert(lpi != NULL);
    2282 assert(lpi->spx != NULL);
    2283
    2284 return (lpi->spx->status() == SPxSolver::UNBOUNDED);
    2285}
    2286
    2287/** returns TRUE iff LP is proven to be dual feasible */
    2289 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2290 )
    2291{
    2292 SCIPdebugMessage("calling SCIPlpiIsDualFeasible()\n");
    2293
    2294 assert(lpi != NULL);
    2295 assert(lpi->spx != NULL);
    2296
    2297 /* note that the solver status may be ABORT_VALUE and the basis status optimal; if we are optimal, isPerturbed() may
    2298 * still return true as long as perturbation plus violation is within tolerances
    2299 */
    2300 assert(lpi->spx->basisStatus() == SPxBasis::OPTIMAL || lpi->spx->status() != SPxSolver::OPTIMAL);
    2301
    2302 return (lpi->spx->basisStatus() == SPxBasis::OPTIMAL) || lpi->spx->basisStatus() == SPxBasis::DUAL;
    2303}
    2304
    2305/** returns TRUE iff LP was solved to optimality */
    2307 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2308 )
    2309{
    2310 SCIPdebugMessage("calling SCIPlpiIsOptimal()\n");
    2311
    2312 assert(lpi != NULL);
    2313 assert(lpi->spx != NULL);
    2314 assert((lpi->spx->basisStatus() == SPxBasis::OPTIMAL)
    2316
    2317 /* note that the solver status may be ABORT_VALUE and the basis status optimal; if we are optimal, isPerturbed() may
    2318 * still return true as long as perturbation plus violation is within tolerances
    2319 */
    2320 return (lpi->spx->basisStatus() == SPxBasis::OPTIMAL);
    2321}
    2322
    2323/** returns TRUE iff the objective limit was reached */
    2325 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2326 )
    2327{
    2328 SCIPdebugMessage("calling SCIPlpiIsObjlimExc()\n");
    2329
    2330 assert(lpi != NULL);
    2331 assert(lpi->spx != NULL);
    2332
    2333 return (lpi->spx->status() == SPxSolver::ABORT_VALUE);
    2334}
    2335
    2336/** returns TRUE iff the iteration limit was reached */
    2338 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2339 )
    2340{
    2341 SCIPdebugMessage("calling SCIPlpiIsIterlimExc()\n");
    2342
    2343 assert(lpi != NULL);
    2344 assert(lpi->spx != NULL);
    2345
    2346 return (lpi->spx->status() == SPxSolver::ABORT_ITER);
    2347}
    2348
    2349/** returns TRUE iff the time limit was reached */
    2351 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2352 )
    2353{
    2354 SCIPdebugMessage("calling SCIPlpiIsTimelimExc()\n");
    2355
    2356 assert(lpi != NULL);
    2357 assert(lpi->spx != NULL);
    2358
    2359 return (lpi->spx->status() == SPxSolver::ABORT_TIME);
    2360}
    2361
    2362/** returns the internal solution status of the solver */
    2364 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2365 )
    2366{
    2367 SCIPdebugMessage("calling SCIPlpiExactGetInternalStatus()\n");
    2368
    2369 assert(lpi != NULL);
    2370 assert(lpi->spx != NULL);
    2371
    2372 return static_cast<int>(lpi->spx->status());
    2373}
    2374
    2375/** tries to reset the internal status of the LP solver in order to ignore an instability of the last solving call */
    2377 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2378 SCIP_Bool* success /**< pointer to store, whether the instability could be ignored */
    2379 )
    2380{ /*lint --e{715}*/
    2381 SCIPdebugMessage("calling SCIPlpiIgnoreInstability()\n");
    2382
    2383 assert(lpi != NULL);
    2384 assert(lpi->spx != NULL);
    2385 assert(success != NULL);
    2386
    2387#if SOPLEX_APIVERSION >= 4
    2388 *success = lpi->spx->ignoreUnscaledViolations();
    2389#else
    2390 *success = FALSE;
    2391#endif
    2392
    2393 return SCIP_OKAY;
    2394}
    2395
    2396/** gets objective value of solution */
    2398 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2399 SCIP_RATIONAL* objval /**< stores the objective value */
    2400 )
    2401{
    2402 SCIPdebugMessage("calling SCIPlpiExactGetObjval()\n");
    2403
    2404 assert(lpi != NULL);
    2405 assert(lpi->spx != NULL);
    2406 assert(objval != NULL);
    2407
    2408 RsetSpxR(lpi, objval, lpi->spx->objValueRational());
    2409
    2410 return SCIP_OKAY;
    2411}
    2412
    2413
    2414/** gets primal and dual solution vectors for feasible LPs
    2415 *
    2416 * Before calling this function, the caller must ensure that the LP has been solved to optimality, i.e., that
    2417 * SCIPlpiIsOptimal() returns true.
    2418 */
    2420 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2421 SCIP_RATIONAL* objval, /**< stores the objective value, may be NULL if not needed */
    2422 SCIP_RATIONAL** primsol, /**< primal solution vector, may be NULL if not needed */
    2423 SCIP_RATIONAL** dualsol, /**< dual solution vector, may be NULL if not needed */
    2424 SCIP_RATIONAL** activity, /**< row activity vector, may be NULL if not needed */
    2425 SCIP_RATIONAL** redcost /**< reduced cost vector, may be NULL if not needed */
    2426 )
    2427{
    2428 DVectorRational* tmpvec;
    2429 SCIPdebugMessage("calling SCIPlpiExactGetSol()\n");
    2430
    2431 assert(lpi != NULL);
    2432 assert(lpi->spx != NULL);
    2433
    2434 tmpvec = new DVectorRational(MAX(lpi->spx->numColsRational(), lpi->spx->numRowsRational()));
    2435 tmpvec->clear();
    2436
    2437 if( objval != NULL )
    2438 RsetSpxR(lpi, objval, lpi->spx->objValueRational());
    2439
    2440 try
    2441 {
    2442 if( primsol != NULL )
    2443 {
    2444 VectorRational tmp(lpi->spx->numColsRational(), tmpvec->get_ptr());
    2445 (void)lpi->spx->getPrimalRational(tmp);
    2446 RsetSpxVector(lpi, primsol, tmp);
    2447 }
    2448 if( dualsol != NULL )
    2449 {
    2450 VectorRational tmp(lpi->spx->numRowsRational(), tmpvec->get_ptr());
    2451 (void)lpi->spx->getDualRational(tmp);
    2452 RsetSpxVector(lpi, dualsol, tmp);
    2453 }
    2454 if( activity != NULL )
    2455 {
    2456 VectorRational tmp(lpi->spx->numRowsRational(), tmpvec->get_ptr());
    2457 (void)lpi->spx->getSlacksRational(tmp); /* in SoPlex, the activities are called "slacks" */
    2458 RsetSpxVector(lpi, activity, tmp);
    2459 }
    2460 if( redcost != NULL )
    2461 {
    2462 VectorRational tmp(lpi->spx->numColsRational(), tmpvec->get_ptr());
    2463 (void)lpi->spx->getRedCostRational(tmp);
    2464 RsetSpxVector(lpi, redcost, tmp);
    2465 }
    2466
    2467 delete tmpvec;
    2468 }
    2469#ifndef NDEBUG
    2470 catch( const SPxException& x )
    2471 {
    2472 std::string s = x.what();
    2473 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    2474#else
    2475 catch( const SPxException& )
    2476 {
    2477#endif
    2478 return SCIP_LPERROR;
    2479 }
    2480
    2481 return SCIP_OKAY;
    2482}
    2483
    2484
    2485/** gets primal ray for unbounded LPs */
    2487 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2488 SCIP_RATIONAL** ray /**< primal ray */
    2489 )
    2490{ /*lint --e{715}*/
    2491 DVectorRational* tmpvec;
    2492 SCIPdebugMessage("calling SCIPlpiExactGetPrimalRay()\n");
    2493
    2494 assert(lpi != NULL);
    2495 assert(lpi->spx != NULL);
    2496 assert(lpi->spx->hasPrimalRay());
    2497 assert(ray != NULL);
    2498
    2499 tmpvec = new DVectorRational(MAX(lpi->spx->numColsRational(), lpi->spx->numRowsRational()));
    2500
    2501 try
    2502 {
    2503 VectorRational tmp(lpi->spx->numColsRational(), tmpvec->get_ptr());
    2504 (void)lpi->spx->getPrimalRayRational(tmp);
    2505 RsetSpxVector(lpi, ray, tmp);
    2506 }
    2507#ifndef NDEBUG
    2508 catch( const SPxException& x )
    2509 {
    2510 std::string s = x.what();
    2511 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    2512#else
    2513 catch( const SPxException& )
    2514 {
    2515#endif
    2516 delete tmpvec;
    2517 return SCIP_LPERROR;
    2518 }
    2519
    2520 delete tmpvec;
    2521
    2522 return SCIP_OKAY;
    2523}
    2524
    2525/** gets dual farkas proof for infeasibility */
    2527 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2528 SCIP_RATIONAL** dualfarkas /**< dual farkas row multipliers */
    2529 )
    2530{
    2531 DVectorRational* tmpvec;
    2532 SCIPdebugMessage("calling SCIPlpiExactGetDualfarkas()\n");
    2533
    2534 assert(lpi != NULL);
    2535 assert(lpi->spx != NULL);
    2536 assert(lpi->spx->hasDualFarkas());
    2537 assert(dualfarkas != NULL);
    2538
    2539 tmpvec = new DVectorRational(MAX(lpi->spx->numColsRational(), lpi->spx->numRowsRational()));
    2540
    2541 try
    2542 {
    2543 VectorRational tmp(lpi->spx->numRowsRational(), tmpvec->get_ptr());
    2544 (void)lpi->spx->getDualFarkasRational(tmp);
    2545 RsetSpxVector(lpi, dualfarkas, tmp);
    2546 }
    2547#ifndef NDEBUG
    2548 catch( const SPxException& x )
    2549 {
    2550 std::string s = x.what();
    2551 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    2552#else
    2553 catch( const SPxException& )
    2554 {
    2555#endif
    2556 return SCIP_LPERROR;
    2557 }
    2558 delete tmpvec;
    2559
    2560 return SCIP_OKAY;
    2561}
    2562
    2563/** gets the number of LP iterations of the last solve call */
    2565 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2566 int* iterations /**< pointer to store the number of iterations of the last solve call */
    2567 )
    2568{
    2569 SCIPdebugMessage("calling SCIPlpiExactGetIterations()\n");
    2570
    2571 assert(lpi != NULL);
    2572 assert(lpi->spx != NULL);
    2573 assert(iterations != NULL);
    2574
    2575 *iterations = lpi->spx->numIterations();
    2576
    2577 return SCIP_OKAY;
    2578}
    2579
    2580/**@} */
    2581
    2582
    2583/*
    2584 * LP Basis Methods
    2585 */
    2586
    2587/**@name LP Basis Methods */
    2588/**@{ */
    2589
    2590
    2591/** gets current basis status for columns and rows; arrays must be large enough to store the basis status */
    2593 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2594 int* cstat, /**< array to store column basis status, or NULL */
    2595 int* rstat /**< array to store row basis status, or NULL */
    2596 )
    2597{
    2598 int i;
    2599
    2600 SCIPdebugMessage("calling SCIPlpiExactGetBase()\n");
    2601
    2602 assert(lpi != NULL);
    2603 assert(lpi->spx != NULL);
    2604
    2605 assert( lpi->spx->preStrongbranchingBasisFreed() );
    2606
    2607 if( rstat != NULL )
    2608 {
    2609 for( i = 0; i < lpi->spx->numRowsRational(); ++i )
    2610 {
    2611 switch( lpi->spx->basisRowStatus(i) )
    2612 {
    2613 case SPxSolver::BASIC:
    2614 rstat[i] = SCIP_BASESTAT_BASIC; /*lint !e641*/
    2615 break;
    2616 case SPxSolver::FIXED:
    2617 case SPxSolver::ON_LOWER:
    2618 rstat[i] = SCIP_BASESTAT_LOWER; /*lint !e641*/
    2619 break;
    2620 case SPxSolver::ON_UPPER:
    2621 rstat[i] = SCIP_BASESTAT_UPPER; /*lint !e641*/
    2622 break;
    2623 case SPxSolver::ZERO:
    2624 SCIPerrorMessage("slack variable has basis status ZERO (should not occur)\n");
    2625 return SCIP_LPERROR;
    2626 case SPxSolver::UNDEFINED:
    2627 default:
    2628 SCIPerrorMessage("invalid basis status\n");
    2629 SCIPABORT();
    2630 return SCIP_INVALIDDATA; /*lint !e527*/
    2631 }
    2632 }
    2633 }
    2634
    2635 if( cstat != NULL )
    2636 {
    2637 for( i = 0; i < lpi->spx->numColsRational(); ++i )
    2638 {
    2639 switch( lpi->spx->basisColStatus(i) )
    2640 {
    2641 case SPxSolver::BASIC:
    2642 cstat[i] = SCIP_BASESTAT_BASIC; /*lint !e641*/
    2643 break;
    2644 case SPxSolver::FIXED:
    2645 /* Get reduced cost estimation. If the estimation is not correct this should not hurt:
    2646 * If the basis is loaded into SoPlex again, the status is converted to FIXED again; in
    2647 * this case there is no problem at all. If the basis is saved and/or used in some other
    2648 * solver, it usually is very cheap to perform the pivots necessary to get an optimal
    2649 * basis.
    2650 * @todo implement getRedCostEst()
    2651 * */
    2652// SCIP_CALL( getRedCostEst(lpi->spx, i, &val) );
    2653// if( val < 0.0 ) /* reduced costs < 0 => UPPER else => LOWER */
    2654// cstat[i] = SCIP_BASESTAT_UPPER; /*lint !e641*/
    2655// else
    2656 cstat[i] = SCIP_BASESTAT_LOWER; /*lint !e641*/
    2657 break;
    2658 case SPxSolver::ON_LOWER:
    2659 cstat[i] = SCIP_BASESTAT_LOWER; /*lint !e641*/
    2660 break;
    2661 case SPxSolver::ON_UPPER:
    2662 cstat[i] = SCIP_BASESTAT_UPPER; /*lint !e641*/
    2663 break;
    2664 case SPxSolver::ZERO:
    2665 cstat[i] = SCIP_BASESTAT_ZERO; /*lint !e641*/
    2666 break;
    2667 case SPxSolver::UNDEFINED:
    2668 default:
    2669 SCIPerrorMessage("invalid basis status\n");
    2670 SCIPABORT();
    2671 return SCIP_INVALIDDATA; /*lint !e527*/
    2672 }
    2673 }
    2674 }
    2675
    2676 return SCIP_OKAY;
    2677}
    2678
    2679/** sets current basis status for columns and rows */
    2681 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2682 int* cstat, /**< array with column basis status */
    2683 int* rstat /**< array with row basis status */
    2684 )
    2685{
    2686 int i;
    2687 int ncols;
    2688 int nrows;
    2689
    2690 SCIPdebugMessage("calling SCIPlpiExactSetBase()\n");
    2691
    2692 assert(lpi != NULL);
    2693 assert(lpi->spx != NULL);
    2694
    2695 SCIP_CALL( SCIPlpiExactGetNRows(lpi, &nrows) );
    2696 SCIP_CALL( SCIPlpiExactGetNCols(lpi, &ncols) );
    2697
    2698 assert(cstat != NULL || ncols == 0);
    2699 assert(rstat != NULL || nrows == 0);
    2700
    2701 assert( lpi->spx->preStrongbranchingBasisFreed() );
    2702 invalidateSolution(lpi);
    2703
    2704 DataArray<SPxSolver::VarStatus>& _colstat = lpi->spx->colStat();
    2705 DataArray<SPxSolver::VarStatus>& _rowstat = lpi->spx->rowStat();
    2706
    2707 _colstat.reSize(ncols);
    2708 _rowstat.reSize(nrows);
    2709
    2710 for( i = 0; i < nrows; ++i )
    2711 {
    2712 switch( rstat[i] ) /*lint !e613*/
    2713 {
    2715 _rowstat[i] = SPxSolver::ON_LOWER;
    2716 break;
    2718 _rowstat[i] = SPxSolver::BASIC;
    2719 break;
    2721 _rowstat[i] = SPxSolver::ON_UPPER;
    2722 break;
    2723 case SCIP_BASESTAT_ZERO:
    2724 SCIPerrorMessage("slack variable has basis status ZERO (should not occur)\n");
    2725 return SCIP_LPERROR; /*lint !e429*/
    2726 default:
    2727 SCIPerrorMessage("invalid basis status\n");
    2728 SCIPABORT();
    2729 return SCIP_INVALIDDATA; /*lint !e527*/
    2730 }
    2731 }
    2732
    2733 for( i = 0; i < ncols; ++i )
    2734 {
    2735 switch( cstat[i] ) /*lint !e613*/
    2736 {
    2738 _colstat[i] = SPxSolver::ON_LOWER;
    2739 break;
    2741 _colstat[i] = SPxSolver::BASIC;
    2742 break;
    2744 _colstat[i] = SPxSolver::ON_UPPER;
    2745 break;
    2746 case SCIP_BASESTAT_ZERO:
    2747 _colstat[i] = SPxSolver::ZERO;
    2748 break;
    2749 default:
    2750 SCIPerrorMessage("invalid basis status\n");
    2751 SCIPABORT();
    2752 return SCIP_INVALIDDATA; /*lint !e527*/
    2753 }
    2754 }
    2755
    2756 SOPLEX_TRY( lpi->messagehdlr, lpi->spx->setBasis(_rowstat.get_ptr(), _colstat.get_ptr()) );
    2758
    2759 return SCIP_OKAY;
    2760}
    2761
    2762/** returns the indices of the basic columns and rows; basic column n gives value n, basic row m gives value -1-m */
    2764 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2765 int* bind /**< pointer to store basis indices ready to keep number of rows entries */
    2766 )
    2767{
    2768 SCIPdebugMessage("calling SCIPlpiExactGetBasisInd()\n");
    2769
    2770 assert(lpi != NULL);
    2771 assert(lpi->spx != NULL);
    2772 assert(bind != NULL);
    2773
    2774 assert(lpi->spx->preStrongbranchingBasisFreed());
    2775
    2776 lpi->spx->getBasisInd(bind);
    2777
    2778 return SCIP_OKAY;
    2779}
    2780
    2781/** get row of inverse basis matrix B^-1
    2782 *
    2783 * @note The LP interface defines slack variables to have coefficient +1. This means that if, internally, the LP solver
    2784 * uses a -1 coefficient, then rows associated with slacks variables whose coefficient is -1, should be negated;
    2785 * see also the explanation in lpi.h.
    2786 */
    2788 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2789 int r, /**< row number */
    2790 SCIP_RATIONAL** coef, /**< pointer to store the coefficients of the row */
    2791 int* inds, /**< array to store the non-zero indices, or NULL */
    2792 int* ninds /**< pointer to store the number of non-zero indices, or NULL
    2793 * (-1: if we do not store sparsity information) */
    2794 )
    2795{
    2796 int i;
    2797 SSVectorBase<soplex::Rational> tmpvec(0);
    2798 SCIPdebugMessage("calling SCIPlpiExactGetBInvRow()\n");
    2799
    2800 assert(lpi != NULL);
    2801 assert(lpi->spx != NULL);
    2802 assert(lpi->spx->preStrongbranchingBasisFreed());
    2803 assert(coef != NULL);
    2804
    2805 assert(r >= 0);
    2806 assert(r < lpi->spx->numRowsRational());
    2807
    2808 if( !lpi->spx->getBasisInverseRowRational(r, tmpvec) )
    2809 return SCIP_LPERROR;
    2810
    2811 for( i = 0; i < tmpvec.size(); ++i )
    2812 {
    2813 inds[i] = tmpvec.index(i);
    2814 RsetSpxR(lpi, coef[i], tmpvec.value(i));
    2815 }
    2816
    2817 *ninds = tmpvec.size();
    2818
    2819 return SCIP_OKAY;
    2820}
    2821
    2822/** get column of inverse basis matrix B^-1
    2823 *
    2824 * @note The LP interface defines slack variables to have coefficient +1. This means that if, internally, the LP solver
    2825 * uses a -1 coefficient, then rows associated with slacks variables whose coefficient is -1, should be negated;
    2826 * see also the explanation in lpi.h.
    2827 */
    2829 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2830 int c, /**< column number of B^-1; this is NOT the number of the column in the LP;
    2831 * you have to call SCIPlpiExactGetBasisInd() to get the array which links the
    2832 * B^-1 column numbers to the row and column numbers of the LP!
    2833 * c must be between 0 and nrows-1, since the basis has the size
    2834 * nrows * nrows */
    2835 SCIP_RATIONAL** coef, /**< pointer to store the coefficients of the column */
    2836 int* inds, /**< array to store the non-zero indices, or NULL */
    2837 int* ninds /**< pointer to store the number of non-zero indices, or NULL
    2838 * (-1: if we do not store sparsity information) */
    2839 )
    2840{
    2841 int i;
    2842 SSVectorRational tmpvec(0);
    2843
    2844 SCIPdebugMessage("calling SCIPlpiExactGetBInvCol()\n");
    2845
    2846 assert( lpi != NULL );
    2847 assert( lpi->spx != NULL );
    2848 assert( lpi->spx->preStrongbranchingBasisFreed() );
    2849 assert(coef != NULL);
    2850
    2851 if( ! lpi->spx->getBasisInverseColRational(c, tmpvec) )
    2852 return SCIP_LPERROR;
    2853
    2854 for( i = 0; i < tmpvec.size(); ++i )
    2855 {
    2856 inds[i] = tmpvec.index(i);
    2857 RsetSpxR(lpi, coef[i], tmpvec.value(i));
    2858 }
    2859
    2860 *ninds = tmpvec.size();
    2861
    2862
    2863 return SCIP_OKAY;
    2864}
    2865
    2866/**@} */
    2867
    2868
    2869/*
    2870 * LP State Methods
    2871 */
    2872
    2873/**@name LP State Methods */
    2874/**@{ */
    2875
    2876/** stores LPi state (like basis information) into lpistate object */
    2878 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2879 BMS_BLKMEM* blkmem, /**< block memory */
    2880 SCIP_LPISTATE** lpistate /**< pointer to LPi state information (like basis information) */
    2881 )
    2882{
    2883 int ncols;
    2884 int nrows;
    2885
    2886 SCIPdebugMessage("calling SCIPlpiExactGetState()\n");
    2887
    2888 assert(blkmem != NULL);
    2889 assert(lpi != NULL);
    2890 assert(lpi->spx != NULL);
    2891 assert(lpistate != NULL);
    2892
    2893 assert( lpi->spx->preStrongbranchingBasisFreed() );
    2894
    2895 ncols = lpi->spx->numColsRational();
    2896 nrows = lpi->spx->numRowsRational();
    2897 assert(ncols >= 0);
    2898 assert(nrows >= 0);
    2899
    2900 /* allocate lpistate data */
    2901 SCIP_CALL( lpistateCreate(lpistate, blkmem, ncols, nrows) );
    2902
    2903 /* allocate enough memory for storing uncompressed basis information */
    2904 SCIP_CALL( ensureCstatMem(lpi, ncols) );
    2905 SCIP_CALL( ensureRstatMem(lpi, nrows) );
    2906
    2907 /* get unpacked basis information */
    2908 SCIP_CALL( SCIPlpiExactGetBase(lpi, lpi->cstat, lpi->rstat) );
    2909
    2910 /* pack LPi state data */
    2911 (*lpistate)->ncols = ncols;
    2912 (*lpistate)->nrows = nrows;
    2913 lpistatePack(*lpistate, lpi->cstat, lpi->rstat);
    2914
    2915 return SCIP_OKAY;
    2916}
    2917
    2918/** loads LPi state (like basis information) into solver; note that the LP might have been extended with additional
    2919 * columns and rows since the state was stored with SCIPlpiExactGetState()
    2920 */
    2922 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    2923 BMS_BLKMEM* blkmem, /**< block memory */
    2924 SCIP_LPISTATE* lpistate /**< LPi state information (like basis information), or NULL */
    2925 )
    2926{ /*lint --e{715}*/
    2927 int lpncols;
    2928 int lpnrows;
    2929 int i;
    2930
    2931 SCIPdebugMessage("calling SCIPlpiExactSetState()\n");
    2932
    2933 assert(lpi != NULL);
    2934 assert(lpi->spx != NULL);
    2935 assert(lpistate != NULL);
    2936
    2937 assert( lpi->spx->preStrongbranchingBasisFreed() );
    2938
    2939 lpncols = lpi->spx->numColsRational();
    2940 lpnrows = lpi->spx->numRowsRational();
    2941 assert(lpistate->ncols <= lpncols);
    2942 assert(lpistate->nrows <= lpnrows);
    2943
    2944 /* allocate enough memory for storing uncompressed basis information */
    2945 SCIP_CALL( ensureCstatMem(lpi, lpncols) );
    2946 SCIP_CALL( ensureRstatMem(lpi, lpnrows) );
    2947
    2948 /* unpack LPi state data */
    2949 lpistateUnpack(lpistate, lpi->cstat, lpi->rstat);
    2950
    2951 /* extend the basis to the current LP beyond the previously existing columns */
    2952 for( i = lpistate->ncols; i < lpncols; ++i )
    2953 {
    2954 SCIP_Real bnd = lpi->spx->lowerReal(i);
    2955 if ( SCIPlpiExactIsInfinity(lpi, REALABS(bnd)) )
    2956 {
    2957 /* if lower bound is +/- infinity -> try upper bound */
    2958 bnd = lpi->spx->lowerReal(i);
    2959 if ( SCIPlpiExactIsInfinity(lpi, REALABS(bnd)) )
    2960 /* variable is free */
    2961 lpi->cstat[i] = SCIP_BASESTAT_ZERO; /*lint !e641*/
    2962 else
    2963 /* use finite upper bound */
    2964 lpi->cstat[i] = SCIP_BASESTAT_UPPER; /*lint !e641*/
    2965 }
    2966 else
    2967 /* use finite lower bound */
    2968 lpi->cstat[i] = SCIP_BASESTAT_LOWER; /*lint !e641*/
    2969 }
    2970 for( i = lpistate->nrows; i < lpnrows; ++i )
    2971 lpi->rstat[i] = SCIP_BASESTAT_BASIC; /*lint !e641*/
    2972
    2973 /* load basis information */
    2974 SCIP_CALL( SCIPlpiExactSetBase(lpi, lpi->cstat, lpi->rstat) );
    2975
    2976 return SCIP_OKAY;
    2977}
    2978
    2979/** clears current LPi state (like basis information) of the solver */
    2981 SCIP_LPIEXACT* lpi /**< LP interface structure */
    2982 )
    2983{ /*lint --e{715}*/
    2984 SCIPdebugMessage("calling SCIPlpiClearState()\n");
    2985
    2986 assert(lpi != NULL);
    2987 assert(lpi->spx != NULL);
    2988
    2989 try
    2990 {
    2991 lpi->spx->clearBasis();
    2992 }
    2993#ifndef NDEBUG
    2994 catch( const SPxException& x )
    2995 {
    2996 std::string s = x.what();
    2997 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    2998#else
    2999 catch( const SPxException& )
    3000 {
    3001#endif
    3002 assert( lpi->spx->status() != SPxSolver::OPTIMAL );
    3003 return SCIP_LPERROR;
    3004 }
    3005
    3006 return SCIP_OKAY;
    3007}
    3008
    3009/** frees LPi state information */
    3011 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3012 BMS_BLKMEM* blkmem, /**< block memory */
    3013 SCIP_LPISTATE** lpistate /**< pointer to LPi state information (like basis information) */
    3014 )
    3015{ /*lint --e{715}*/
    3016 SCIPdebugMessage("calling SCIPlpiFreeState()\n");
    3017
    3018 assert(lpi != NULL);
    3019 assert(lpistate != NULL);
    3020 assert(blkmem != NULL);
    3021
    3022 if ( *lpistate != NULL )
    3023 lpistateFree(lpistate, blkmem);
    3024
    3025 return SCIP_OKAY;
    3026}
    3027
    3028/** checks, whether the given LP state contains simplex basis information */
    3030 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3031 SCIP_LPISTATE* lpistate /**< LP state information (like basis information), or NULL */
    3032 )
    3033{ /*lint --e{715}*/
    3034 assert(lpi != NULL);
    3035 return TRUE;
    3036}
    3037
    3038/** reads LP state (like basis information from a file */
    3040 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3041 const char* fname /**< file name */
    3042 )
    3043{
    3044 SCIPdebugMessage("calling SCIPlpiReadState()\n");
    3045 assert(lpi != NULL);
    3046 assert(lpi->spx != NULL);
    3047 assert(fname != NULL);
    3048
    3049 assert( lpi->spx->preStrongbranchingBasisFreed() );
    3050
    3051 bool success;
    3052 SOPLEX_TRY( lpi->messagehdlr, success = lpi->spx->readBasisFile(fname, 0, 0) );
    3053
    3054 return success ? SCIP_OKAY : SCIP_LPERROR;
    3055}
    3056
    3057/** writes LPi state (i.e. basis information) to a file */
    3059 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3060 const char* fname /**< file name */
    3061 )
    3062{
    3063 assert(lpi != NULL);
    3064 assert(lpi->spx != NULL);
    3065 assert(fname != NULL);
    3066 SCIPdebugMessage("calling SCIPlpiWriteState()\n");
    3067
    3068 assert( lpi->spx->preStrongbranchingBasisFreed() );
    3069
    3070 bool res;
    3071 SOPLEX_TRY( lpi->messagehdlr, res = lpi->spx->writeBasisFile(fname, 0, 0) );
    3072
    3073 if ( ! res )
    3074 return SCIP_LPERROR;
    3075
    3076 return SCIP_OKAY;
    3077}
    3078
    3079/**@} */
    3080
    3081/*
    3082 * Parameter Methods
    3083 */
    3084
    3085/**@name Parameter Methods */
    3086/**@{ */
    3087
    3088/** gets integer parameter of LP */
    3090 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3091 SCIP_LPPARAM type, /**< parameter number */
    3092 int* ival /**< buffer to store the parameter value */
    3093 )
    3094{
    3095 int scaleparam;
    3096
    3097 SCIPdebugMessage("calling SCIPlpiExactGetIntpar()\n");
    3098
    3099 assert(lpi != NULL);
    3100 assert(lpi->spx != NULL);
    3101 assert(ival != NULL);
    3102
    3103 switch( type )
    3104 {
    3106 *ival = (int) lpi->spx->getFromScratch();
    3107 break;
    3108 case SCIP_LPPAR_LPINFO:
    3109 *ival = (int) lpi->spx->getLpInfo();
    3110 break;
    3111 case SCIP_LPPAR_LPITLIM:
    3112 *ival = lpi->spx->intParam(SoPlex::ITERLIMIT);
    3113 if( *ival == -1 )
    3114 *ival = INT_MAX;
    3115 break;
    3117 *ival = lpi->spx->intParam(SoPlex::SIMPLIFIER) == (int) SoPlex::SIMPLIFIER_AUTO;
    3118 break;
    3119 case SCIP_LPPAR_PRICING:
    3120 *ival = (int) lpi->pricing;
    3121 break;
    3122 case SCIP_LPPAR_SCALING:
    3123 scaleparam = lpi->spx->intParam(SoPlex::SCALER);
    3124
    3125 if( scaleparam == SoPlex::SCALER_OFF )
    3126 *ival = 0;
    3127 else if( scaleparam == SoPlex::SCALER_BIEQUI )
    3128 *ival = 1;
    3129#if SOPLEX_VERSION > 221 || (SOPLEX_VERSION == 221 && SOPLEX_SUBVERSION >= 2)
    3130 else
    3131 {
    3132 assert(scaleparam == SoPlex::SCALER_LEASTSQ);
    3133 *ival = 2;
    3134 }
    3135#else
    3136 else
    3137 {
    3138 assert(scaleparam == SoPlex::SCALER_GEO8);
    3139 *ival = 2;
    3140 }
    3141#endif
    3142 break;
    3143#if SOPLEX_VERSION >= 201
    3144 case SCIP_LPPAR_TIMING:
    3145 *ival = (int) (lpi->spx->intParam(SoPlex::TIMER));
    3146 break;
    3147#endif
    3148#if SOPLEX_VERSION >= 230 || (SOPLEX_VERSION == 220 && SOPLEX_SUBVERSION >= 3)
    3150 *ival = (int) lpi->spx->randomSeed();
    3151 break;
    3152#endif
    3153#if SOPLEX_APIVERSION >= 1
    3155 *ival = (int) lpi->spx->intParam(SoPlex::FACTOR_UPDATE_MAX);
    3156 break;
    3157#endif
    3158 default:
    3159 return SCIP_PARAMETERUNKNOWN;
    3160 } /*lint !e788*/
    3161
    3162 return SCIP_OKAY;
    3163}
    3164
    3165/** sets integer parameter of LP */
    3167 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3168 SCIP_LPPARAM type, /**< parameter number */
    3169 int ival /**< parameter value */
    3170 )
    3171{
    3172 SCIPdebugMessage("calling SCIPlpiExactSetIntpar()\n");
    3173
    3174 assert(lpi != NULL);
    3175 assert(lpi->spx != NULL);
    3176
    3177 switch( type )
    3178 {
    3180 assert(ival == TRUE || ival == FALSE);
    3181 lpi->spx->setFromScratch(bool(ival));
    3182 break;
    3183 case SCIP_LPPAR_LPINFO:
    3184 assert(ival == TRUE || ival == FALSE);
    3185 lpi->spx->setLpInfo(bool(ival));
    3186 break;
    3187 case SCIP_LPPAR_LPITLIM:
    3188 assert( ival >= 0 );
    3189 /* -1 <= ival, -1 meaning no time limit, 0 stopping immediately */
    3190 if( ival >= INT_MAX )
    3191 ival = -1;
    3192 (void) lpi->spx->setIntParam(SoPlex::ITERLIMIT, ival);
    3193 break;
    3195 assert(ival == TRUE || ival == FALSE);
    3196 (void) lpi->spx->setIntParam(SoPlex::SIMPLIFIER, (int) (ival ? SoPlex::SIMPLIFIER_AUTO : SoPlex::SIMPLIFIER_OFF));
    3197 break;
    3198 case SCIP_LPPAR_PRICING:
    3199 lpi->pricing = (SCIP_PRICING)ival;
    3200 switch( lpi->pricing )
    3201 {
    3203 case SCIP_PRICING_AUTO:
    3204 (void) lpi->spx->setIntParam(SoPlex::PRICER, (int) SoPlex::PRICER_AUTO);
    3205 break;
    3206 case SCIP_PRICING_FULL:
    3207 (void) lpi->spx->setIntParam(SoPlex::PRICER, (int) SoPlex::PRICER_STEEP);
    3208 break;
    3210 (void) lpi->spx->setIntParam(SoPlex::PRICER, (int) SoPlex::PRICER_PARMULT);
    3211 break;
    3212 case SCIP_PRICING_STEEP:
    3213 (void) lpi->spx->setIntParam(SoPlex::PRICER, (int) SoPlex::PRICER_STEEP);
    3214 break;
    3216 (void) lpi->spx->setIntParam(SoPlex::PRICER, (int) SoPlex::PRICER_QUICKSTEEP);
    3217 break;
    3218 case SCIP_PRICING_DEVEX:
    3219 (void) lpi->spx->setIntParam(SoPlex::PRICER, (int) SoPlex::PRICER_DEVEX);
    3220 break;
    3221 default:
    3222 return SCIP_LPERROR;
    3223 }
    3224 break;
    3225 case SCIP_LPPAR_SCALING:
    3226 assert(ival >= 0 && ival <= 2);
    3227 if( ival == 0 )
    3228 (void) lpi->spx->setIntParam(SoPlex::SCALER, (int) SoPlex::SCALER_OFF);
    3229 else if( ival == 1 )
    3230 (void) lpi->spx->setIntParam(SoPlex::SCALER, (int) SoPlex::SCALER_BIEQUI);
    3231 else
    3232#if SOPLEX_VERSION > 221 || (SOPLEX_VERSION == 221 && SOPLEX_SUBVERSION >= 2)
    3233 (void) lpi->spx->setIntParam(SoPlex::SCALER, (int) SoPlex::SCALER_LEASTSQ);
    3234#else
    3235 (void) lpi->spx->setIntParam(SoPlex::SCALER, (int) SoPlex::SCALER_GEO8);
    3236#endif
    3237
    3238 break;
    3239#if SOPLEX_VERSION >= 201
    3240 case SCIP_LPPAR_TIMING:
    3241 assert(ival >= 0 && ival < 3);
    3242 (void) lpi->spx->setIntParam(SoPlex::TIMER, ival);
    3243 break;
    3244#endif
    3245#if SOPLEX_VERSION > 221 || (SOPLEX_VERSION == 221 && SOPLEX_SUBVERSION >= 3)
    3247 lpi->spx->setRandomSeed((unsigned long)(long)ival);
    3248 break;
    3249#endif
    3250#if SOPLEX_VERSION > 221 || (SOPLEX_VERSION >= 221 && SOPLEX_SUBVERSION >= 3)
    3252 assert(ival >= 0 && ival < 3);
    3253 (void) lpi->spx->setIntParam(SoPlex::SOLUTION_POLISHING, ival);
    3254 break;
    3255#endif
    3256#if SOPLEX_APIVERSION >= 1
    3258 assert(ival >= 0);
    3259 (void) lpi->spx->setIntParam(SoPlex::FACTOR_UPDATE_MAX, ival);
    3260 break;
    3261#endif
    3262
    3263 default:
    3264 return SCIP_PARAMETERUNKNOWN;
    3265 } /*lint !e788*/
    3266
    3267 return SCIP_OKAY;
    3268}
    3269
    3270/** gets floating point parameter of LP */
    3272 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3273 SCIP_LPPARAM type, /**< parameter number */
    3274 SCIP_Real* dval /**< buffer to store the parameter value */
    3275 )
    3276{
    3277 SCIPdebugMessage("calling SCIPlpiExactGetRealpar()\n");
    3278
    3279 assert(lpi != NULL);
    3280 assert(lpi->spx != NULL);
    3281 assert(dval != NULL);
    3282
    3283 switch( type )
    3284 {
    3285 case SCIP_LPPAR_FEASTOL:
    3286 *dval = 0.0; //lpi->spx->feastol();
    3287 break;
    3289 *dval = 0.0; //lpi->spx->opttol();
    3290 break;
    3291 case SCIP_LPPAR_OBJLIM:
    3292 if ( lpi->spx->intParam(SoPlex::OBJSENSE) == SoPlex::OBJSENSE_MINIMIZE )
    3293 *dval = lpi->spx->realParam(SoPlex::OBJLIMIT_UPPER);
    3294 else
    3295 *dval = lpi->spx->realParam(SoPlex::OBJLIMIT_LOWER);
    3296 break;
    3297 case SCIP_LPPAR_LPTILIM:
    3298 *dval = lpi->spx->realParam(SoPlex::TIMELIMIT);
    3299 break;
    3301 *dval = lpi->spx->realParam(SoPlex::REPRESENTATION_SWITCH);
    3302 if( *dval >= SCIPlpiExactInfinity(lpi) )
    3303 *dval = -1.0;
    3304 break;
    3306 *dval = lpi->conditionlimit;
    3307 break;
    3308 default:
    3309 return SCIP_PARAMETERUNKNOWN;
    3310 } /*lint !e788*/
    3311
    3312 return SCIP_OKAY;
    3313}
    3314
    3315/** sets floating point parameter of LP */
    3317 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3318 SCIP_LPPARAM type, /**< parameter number */
    3319 SCIP_Real dval /**< parameter value */
    3320 )
    3321{
    3322 SCIPdebugMessage("calling SCIPlpiExactSetRealpar()\n");
    3323
    3324 assert(lpi != NULL);
    3325 assert(lpi->spx != NULL);
    3326
    3327 switch( type )
    3328 {
    3329 case SCIP_LPPAR_OBJLIM:
    3330 /* no restrictions on dval */
    3331 if ( lpi->spx->intParam(SoPlex::OBJSENSE) == SoPlex::OBJSENSE_MINIMIZE )
    3332 (void) lpi->spx->setRealParam(SoPlex::OBJLIMIT_UPPER, dval);
    3333 else
    3334 (void) lpi->spx->setRealParam(SoPlex::OBJLIMIT_LOWER, dval);
    3335 break;
    3336 case SCIP_LPPAR_LPTILIM:
    3337 assert( dval > 0.0 );
    3338 /* soplex requires 0 < dval < DEFAULT_INFINITY (= 1e100), -1 means unlimited */
    3339 (void) lpi->spx->setRealParam(SoPlex::TIMELIMIT, dval);
    3340 break;
    3342 /* 0 <= dval <= inf */
    3343 assert( dval >= 0.0 || dval == -1.0 );
    3344 if( dval == -1 )
    3345 (void) lpi->spx->setRealParam(SoPlex::REPRESENTATION_SWITCH, SCIPlpiExactInfinity(lpi));
    3346 else
    3347 (void) lpi->spx->setRealParam(SoPlex::REPRESENTATION_SWITCH, dval);
    3348 break;
    3350 lpi->conditionlimit = dval;
    3351 lpi->checkcondition = (dval >= 0.0);
    3352 break;
    3353 default:
    3354 return SCIP_PARAMETERUNKNOWN;
    3355 } /*lint !e788*/
    3356
    3357 return SCIP_OKAY;
    3358}
    3359
    3360/**@} */
    3361
    3362
    3363/*
    3364 * Numerical Methods
    3365 */
    3366
    3367/**@name Numerical Methods */
    3368/**@{ */
    3369
    3370/** returns value treated as infinity in the LP solver */
    3372 SCIP_LPIEXACT* lpi /**< LP interface structure */
    3373 )
    3374{
    3375 assert(lpi != NULL);
    3376 SCIPdebugMessage("calling SCIPlpiInfinity()\n");
    3377
    3378 return lpi->spx->realParam(SoPlex::INFTY);
    3379}
    3380
    3381/** checks if given value is treated as infinity in the LP solver */
    3383 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3384 SCIP_Real val /**< the value */
    3385 )
    3386{
    3387 assert(lpi != NULL);
    3388 SCIPdebugMessage("calling SCIPlpiExactIsInfinity()\n");
    3389
    3390 return (val >= lpi->spx->realParam(SoPlex::INFTY));
    3391}
    3392
    3393/**@} */
    3394
    3395
    3396/*
    3397 * File Interface Methods
    3398 */
    3399
    3400/**@name File Interface Methods */
    3401/**@{ */
    3402
    3403/** returns, whether the given file exists */
    3404static
    3406 const char* filename /**< file name */
    3407 )
    3408{
    3409 FILE* f;
    3410
    3411 f = fopen(filename, "r");
    3412 if( f == NULL )
    3413 return FALSE;
    3414
    3415 fclose(f);
    3416
    3417 return TRUE;
    3418}
    3419
    3420/** reads LP from a file */
    3422 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3423 const char* fname /**< file name */
    3424 )
    3425{
    3426 SCIPdebugMessage("calling SCIPlpiReadLP()\n");
    3427
    3428 assert(lpi != NULL);
    3429 assert(lpi->spx != NULL);
    3430 assert(fname != NULL);
    3431
    3432 assert( lpi->spx->preStrongbranchingBasisFreed() );
    3433
    3434 if( !fileExists(fname) )
    3435 return SCIP_NOFILE;
    3436
    3437 try
    3438 {
    3439 assert(lpi->spx->intParam(SoPlex::READMODE) == SoPlex::READMODE_RATIONAL);
    3440 if( !lpi->spx->readFile(fname) )
    3441 return SCIP_READERROR;
    3442 }
    3443#ifndef NDEBUG
    3444 catch( const SPxException& x )
    3445 {
    3446 std::string s = x.what();
    3447 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    3448#else
    3449 catch( const SPxException& )
    3450 {
    3451#endif
    3452 return SCIP_READERROR;
    3453 }
    3454
    3455 return SCIP_OKAY;
    3456}
    3457
    3458/** writes LP to a file */
    3460 SCIP_LPIEXACT* lpi, /**< LP interface structure */
    3461 const char* fname /**< file name */
    3462 )
    3463{
    3464 SCIPdebugMessage("calling SCIPlpiWriteLP()\n");
    3465
    3466 assert(lpi != NULL);
    3467 assert(lpi->spx != NULL);
    3468 assert(fname != NULL);
    3469
    3470 try
    3471 {
    3472 (void) lpi->spx->writeFileRational(fname);
    3473 }
    3474#ifndef NDEBUG
    3475 catch( const SPxException& x )
    3476 {
    3477 std::string s = x.what();
    3478 SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
    3479#else
    3480 catch( const SPxException& )
    3481 {
    3482#endif
    3483 return SCIP_WRITEERROR;
    3484 }
    3485
    3486 return SCIP_OKAY;
    3487}
    3488
    3489/**@} */
    void SCIPdecodeDualBit(const SCIP_DUALPACKET *inp, int *out, int count)
    Definition: bitencode.c:308
    void SCIPencodeDualBit(const int *inp, SCIP_DUALPACKET *out, int count)
    Definition: bitencode.c:238
    packing single and dual bit values
    unsigned int SCIP_DUALPACKET
    Definition: bitencode.h:42
    SCIP_Real * r
    Definition: circlepacking.c:59
    SCIP_VAR ** x
    Definition: circlepacking.c:63
    void setProbname(const char *probname)
    DataArray< SPxSolver::VarStatus > & colStat()
    bool checkConsistentSides() const
    bool getFromScratch() const
    SPxSolver::Status doSolve(bool printwarning=true)
    void savePreStrongbranchingBasis()
    bool preStrongbranchingBasisFreed() const
    SPxexSCIP(SCIP_MESSAGEHDLR *messagehdlr=NULL, const char *probname=NULL)
    void freePreStrongbranchingBasis()
    Real getObjLimit() const
    void setRep(SPxSolver::Representation p_rep)
    void trySolve(bool printwarning=true)
    bool getLpInfo() const
    bool checkConsistentBounds() const
    void setFromScratch(bool fs)
    void setLpInfo(bool lpinfo)
    DataArray< SPxSolver::VarStatus > & rowStat()
    virtual ~SPxexSCIP()
    void restorePreStrongbranchingBasis()
    common defines and data types used in all packages of SCIP
    #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 MAX(x, y)
    Definition: def.h:220
    #define SCIPABORT()
    Definition: def.h:327
    #define REALABS(x)
    Definition: def.h:182
    #define SCIP_CALL(x)
    Definition: def.h:355
    void * SCIPlpiExactGetSolverPointer(SCIP_LPIEXACT *lpi)
    SCIP_Bool SCIPlpiExactHasDualRay(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactSetRealpar(SCIP_LPIEXACT *lpi, SCIP_LPPARAM type, SCIP_Real dval)
    SCIP_RETCODE SCIPlpiExactSetBase(SCIP_LPIEXACT *lpi, int *cstat, int *rstat)
    SCIP_RETCODE SCIPlpiExactReadState(SCIP_LPIEXACT *lpi, const char *fname)
    SCIP_Bool SCIPlpiExactHasStateBasis(SCIP_LPIEXACT *lpi, SCIP_LPISTATE *lpistate)
    SCIP_RETCODE SCIPlpiExactGetObj(SCIP_LPIEXACT *lpi, int firstcol, int lastcol, SCIP_RATIONAL **vals)
    SCIP_RETCODE SCIPlpiExactStartStrongbranch(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactIgnoreInstability(SCIP_LPIEXACT *lpi, SCIP_Bool *success)
    SCIP_RETCODE SCIPlpiExactGetObjval(SCIP_LPIEXACT *lpi, SCIP_RATIONAL *objval)
    SCIP_Bool SCIPlpiExactIsDualUnbounded(SCIP_LPIEXACT *lpi)
    SCIP_Bool SCIPlpiExactHasBarrierSolve(void)
    SCIP_RETCODE SCIPlpiExactWriteLP(SCIP_LPIEXACT *lpi, const char *fname)
    SCIP_Bool SCIPlpiExactHasPrimalRay(SCIP_LPIEXACT *lpi)
    SCIP_Bool SCIPlpiExactExistsDualRay(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactGetDualfarkas(SCIP_LPIEXACT *lpi, SCIP_RATIONAL **dualfarkas)
    SCIP_RETCODE SCIPlpiExactEndStrongbranch(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactSolveDual(SCIP_LPIEXACT *lpi)
    const char * SCIPlpiExactGetSolverDesc(void)
    SCIP_RETCODE SCIPlpiExactChgObjsen(SCIP_LPIEXACT *lpi, SCIP_OBJSEN objsen)
    SCIP_RETCODE SCIPlpiExactSetIntpar(SCIP_LPIEXACT *lpi, SCIP_LPPARAM type, int ival)
    SCIP_RETCODE SCIPlpiExactWriteState(SCIP_LPIEXACT *lpi, const char *fname)
    SCIP_RETCODE SCIPlpiExactChgCoef(SCIP_LPIEXACT *lpi, int row, int col, SCIP_RATIONAL *newval)
    SCIP_Bool SCIPlpiExactWasSolved(SCIP_LPIEXACT *lpi)
    SCIP_Bool SCIPlpiExactIsPrimalUnbounded(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactGetSides(SCIP_LPIEXACT *lpi, int firstrow, int lastrow, SCIP_RATIONAL **lhss, SCIP_RATIONAL **rhss)
    SCIP_Bool SCIPlpiExactHasDualSolve(void)
    const char * SCIPlpiExactGetSolverName(void)
    SCIP_RETCODE SCIPlpiExactGetCols(SCIP_LPIEXACT *lpi, int firstcol, int lastcol, SCIP_RATIONAL **lb, SCIP_RATIONAL **ub, int *nnonz, int *beg, int *ind, SCIP_RATIONAL **val)
    SCIP_Bool SCIPlpiExactIsPrimalInfeasible(SCIP_LPIEXACT *lpi)
    SCIP_Real SCIPlpiExactInfinity(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactGetRows(SCIP_LPIEXACT *lpi, int firstrow, int lastrow, SCIP_RATIONAL **lhs, SCIP_RATIONAL **rhs, int *nnonz, int *beg, int *ind, SCIP_RATIONAL **val)
    SCIP_RETCODE SCIPlpiExactGetSolFeasibility(SCIP_LPIEXACT *lpi, SCIP_Bool *primalfeasible, SCIP_Bool *dualfeasible)
    SCIP_RETCODE SCIPlpiExactGetBInvCol(SCIP_LPIEXACT *lpi, int c, SCIP_RATIONAL **coef, int *inds, int *ninds)
    SCIP_Bool SCIPlpiExactExistsPrimalRay(SCIP_LPIEXACT *lpi)
    SCIP_Bool SCIPlpiExactIsInfinity(SCIP_LPIEXACT *lpi, SCIP_Real val)
    SCIP_Bool SCIPlpiExactHasPrimalSolve(void)
    SCIP_RETCODE SCIPlpiExactAddRows(SCIP_LPIEXACT *lpi, int nrows, SCIP_RATIONAL **lhs, SCIP_RATIONAL **rhs, char **rownames, int nnonz, int *beg, int *ind, SCIP_RATIONAL **val)
    SCIP_RETCODE SCIPlpiExactCreate(SCIP_LPIEXACT **lpi, SCIP_MESSAGEHDLR *messagehdlr, const char *name, SCIP_OBJSEN objsen)
    SCIP_Bool SCIPlpiExactIsIterlimExc(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactLoadColLP(SCIP_LPIEXACT *lpi, SCIP_OBJSEN objsen, int ncols, SCIP_RATIONAL **obj, SCIP_RATIONAL **lb, SCIP_RATIONAL **ub, char **colnames, int nrows, SCIP_RATIONAL **lhs, SCIP_RATIONAL **rhs, char **, int nnonz, int *beg, int *ind, SCIP_RATIONAL **val)
    SCIP_Bool SCIPlpiExactIsPrimalFeasible(SCIP_LPIEXACT *lpi)
    SCIP_Bool SCIPlpiExactIsDualFeasible(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactReadLP(SCIP_LPIEXACT *lpi, const char *fname)
    SCIP_RETCODE SCIPlpiExactGetBasisInd(SCIP_LPIEXACT *lpi, int *bind)
    SCIP_RETCODE SCIPlpiExactGetObjsen(SCIP_LPIEXACT *lpi, SCIP_OBJSEN *objsen)
    SCIP_Bool SCIPlpiExactIsOptimal(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactGetRowNames(SCIP_LPIEXACT *lpi, int firstrow, int lastrow, char **rownames, char *namestorage, int namestoragesize, int *storageleft)
    SCIP_RETCODE SCIPlpiExactAddCols(SCIP_LPIEXACT *lpi, int ncols, SCIP_RATIONAL **obj, SCIP_RATIONAL **lb, SCIP_RATIONAL **ub, char **colnames, int nnonz, int *beg, int *ind, SCIP_RATIONAL **val)
    SCIP_RETCODE SCIPlpiExactGetBase(SCIP_LPIEXACT *lpi, int *cstat, int *rstat)
    SCIP_RETCODE SCIPlpiExactDelCols(SCIP_LPIEXACT *lpi, int firstcol, int lastcol)
    SCIP_RETCODE SCIPlpiExactChgObj(SCIP_LPIEXACT *lpi, int ncols, int *ind, SCIP_RATIONAL **obj)
    SCIP_RETCODE SCIPlpiExactGetPrimalRay(SCIP_LPIEXACT *lpi, SCIP_RATIONAL **ray)
    SCIP_RETCODE SCIPlpiExactGetBounds(SCIP_LPIEXACT *lpi, int firstcol, int lastcol, SCIP_RATIONAL **lbs, SCIP_RATIONAL **ubs)
    int SCIPlpiExactGetInternalStatus(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiDelColset(SCIP_LPIEXACT *lpi, int *dstat)
    SCIP_RETCODE SCIPlpiExactGetCoef(SCIP_LPIEXACT *lpi, int row, int col, SCIP_RATIONAL *val)
    SCIP_RETCODE SCIPlpiExactGetColNames(SCIP_LPIEXACT *lpi, int firstcol, int lastcol, char **colnames, char *namestorage, int namestoragesize, int *storageleft)
    SCIP_RETCODE SCIPlpiExactClear(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactGetNNonz(SCIP_LPIEXACT *lpi, int *nnonz)
    SCIP_RETCODE SCIPlpiExactFree(SCIP_LPIEXACT **lpi)
    SCIP_Bool SCIPlpiExactIsDualInfeasible(SCIP_LPIEXACT *lpi)
    SCIP_Bool SCIPlpiExactIsObjlimExc(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactClearState(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactGetNCols(SCIP_LPIEXACT *lpi, int *ncols)
    SCIP_RETCODE SCIPlpiExactSolveBarrier(SCIP_LPIEXACT *lpi, SCIP_Bool crossover)
    SCIP_RETCODE SCIPlpiExactGetBInvRow(SCIP_LPIEXACT *lpi, int r, SCIP_RATIONAL **coef, int *inds, int *ninds)
    SCIP_RETCODE SCIPlpiExactChgSides(SCIP_LPIEXACT *lpi, int nrows, int *ind, SCIP_RATIONAL **lhs, SCIP_RATIONAL **rhs)
    SCIP_RETCODE SCIPlpiExactGetState(SCIP_LPIEXACT *lpi, BMS_BLKMEM *blkmem, SCIP_LPISTATE **lpistate)
    SCIP_RETCODE SCIPlpiExactDelRows(SCIP_LPIEXACT *lpi, int firstrow, int lastrow)
    SCIP_RETCODE SCIPlpiExactGetNRows(SCIP_LPIEXACT *lpi, int *nrows)
    SCIP_Bool SCIPlpiExactIsTimelimExc(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactGetRealpar(SCIP_LPIEXACT *lpi, SCIP_LPPARAM type, SCIP_Real *dval)
    SCIP_RETCODE SCIPlpiExactChgBounds(SCIP_LPIEXACT *lpi, int ncols, int *ind, SCIP_RATIONAL **lb, SCIP_RATIONAL **ub)
    SCIP_RETCODE SCIPlpiExactFreeState(SCIP_LPIEXACT *lpi, BMS_BLKMEM *blkmem, SCIP_LPISTATE **lpistate)
    SCIP_RETCODE SCIPlpiExactDelRowset(SCIP_LPIEXACT *lpi, int *dstat)
    SCIP_RETCODE SCIPlpiExactGetSol(SCIP_LPIEXACT *lpi, SCIP_RATIONAL *objval, SCIP_RATIONAL **primsol, SCIP_RATIONAL **dualsol, SCIP_RATIONAL **activity, SCIP_RATIONAL **redcost)
    SCIP_RETCODE SCIPlpiExactSolvePrimal(SCIP_LPIEXACT *lpi)
    SCIP_RETCODE SCIPlpiExactGetIntpar(SCIP_LPIEXACT *lpi, SCIP_LPPARAM type, int *ival)
    SCIP_RETCODE SCIPlpiExactSetState(SCIP_LPIEXACT *lpi, BMS_BLKMEM *blkmem, SCIP_LPISTATE *lpistate)
    SCIP_RETCODE SCIPlpiExactGetIterations(SCIP_LPIEXACT *lpi, int *iterations)
    SCIP_RETCODE SCIPlpiExactSetIntegralityInformation(SCIP_LPIEXACT *lpi, int ncols, int *intInfo)
    void SCIPrationalSetInfinity(SCIP_RATIONAL *res)
    Definition: rational.cpp:618
    SCIP_Real SCIPrationalGetReal(SCIP_RATIONAL *rational)
    Definition: rational.cpp:2085
    SCIP_Bool SCIPrationalIsAbsInfinity(SCIP_RATIONAL *rational)
    Definition: rational.cpp:1680
    SCIP_Bool SCIPrationalIsPositive(SCIP_RATIONAL *rational)
    Definition: rational.cpp:1640
    SCIP_Bool SCIPrationalIsZero(SCIP_RATIONAL *rational)
    Definition: rational.cpp:1624
    void SCIPrationalSetNegInfinity(SCIP_RATIONAL *res)
    Definition: rational.cpp:630
    SCIP_Bool SCIPrationalIsInfinity(SCIP_RATIONAL *rational)
    Definition: rational.cpp:1660
    SCIP_Bool SCIPrationalIsNegInfinity(SCIP_RATIONAL *rational)
    Definition: rational.cpp:1670
    static SCIP_RETCODE optimize(SCIP *scip, SCIP_SOL *worksol, SCIP_VAR **vars, int *blockstart, int *blockend, int nblocks, OPTTYPE opttype, SCIP_Real *activities, int nrows, SCIP_Bool *improvement, SCIP_Bool *varboundserr, SCIP_HEURDATA *heurdata)
    Definition: heur_twoopt.c:967
    SCIP_DUALPACKET ROWPACKET
    Definition: lpi_clp.cpp:128
    SCIP_DUALPACKET COLPACKET
    Definition: lpi_clp.cpp:126
    static void setIntParam(SCIP_LPI *lpi, int const param, int const parval)
    Definition: lpi_cpx.c:654
    #define FEASTOL
    Definition: lpi_qso.c:99
    interface methods for specific exact LP solvers
    static void RsetSpxR(SCIP_LPIEXACT *lpi, SCIP_RATIONAL *r, const soplex::Rational &spxr)
    static void lpistatePack(SCIP_LPISTATE *lpistate, const int *cstat, const int *rstat)
    #define NULL
    static void lpistateUnpack(const SCIP_LPISTATE *lpistate, int *cstat, int *rstat)
    static SCIP_Bool fileExists(const char *filename)
    #define SOPLEX_TRY(messagehdlr, x)
    static int rowpacketNum(int nrows)
    SCIP_DUALPACKET ROWPACKET
    #define SOPLEX_VERBLEVEL
    static const char spxdesc[200]
    #define COLS_PER_PACKET
    static void invalidateSolution(SCIP_LPIEXACT *lpi)
    static void lpistateFree(SCIP_LPISTATE **lpistate, BMS_BLKMEM *blkmem)
    SCIP_DUALPACKET COLPACKET
    static SCIP_RETCODE ensureCstatMem(SCIP_LPIEXACT *lpi, int num)
    static SCIP_RETCODE spxSolve(SCIP_LPIEXACT *lpi)
    static int colpacketNum(int ncols)
    static SCIP_RETCODE ensureRstatMem(SCIP_LPIEXACT *lpi, int num)
    #define CHECK_SOPLEX_PARAM(x)
    #define SOPLEX_SUBVERSION
    static void SpxRSetRat(SCIP_LPIEXACT *lpi, soplex::Rational &spxr, SCIP_RATIONAL *src)
    static SCIP_RETCODE lpistateCreate(SCIP_LPISTATE **lpistate, BMS_BLKMEM *blkmem, int ncols, int nrows)
    #define SOPLEX_TRY_ABORT(x)
    #define ROWS_PER_PACKET
    static const char spxname[20]
    static void RsetSpxVector(SCIP_LPIEXACT *lpi, SCIP_RATIONAL **r, VectorRational src)
    memory allocation routines
    #define BMSfreeMemory(ptr)
    Definition: memory.h:145
    #define BMSfreeBlockMemory(mem, ptr)
    Definition: memory.h:465
    #define BMSallocBlockMemory(mem, ptr)
    Definition: memory.h:451
    #define BMSreallocMemoryArray(ptr, num)
    Definition: memory.h:127
    #define BMSallocMemoryCPP(size)
    Definition: memory.h:121
    #define BMSallocBlockMemoryArray(mem, ptr, num)
    Definition: memory.h:454
    #define BMSfreeBlockMemoryArray(mem, ptr, num)
    Definition: memory.h:467
    struct BMS_BlkMem BMS_BLKMEM
    Definition: memory.h:437
    #define BMSfreeMemoryArrayNull(ptr)
    Definition: memory.h:148
    #define BMSallocMemory(ptr)
    Definition: memory.h:118
    void SCIPmessagePrintWarning(SCIP_MESSAGEHDLR *messagehdlr, const char *formatstr,...)
    Definition: message.c:427
    static int COI_CALLCONV Status(int MODSTA, int SOLSTA, int ITER, double OBJVAL, void *USRMEM)
    Definition: nlpi_conopt.c:551
    public methods for message output
    #define SCIPerrorMessage
    Definition: pub_message.h:64
    #define SCIPdebugMessage
    Definition: pub_message.h:96
    wrapper for rational number arithmetic
    wrapper for rational number arithmetic that interacts with GMP
    #define UNKNOWN
    Definition: sepa_mcf.c:4110
    SCIP_PRICING pricing
    SCIP_Bool solved
    SCIP_MESSAGEHDLR * messagehdlr
    SPxexSCIP * spx
    SCIP_Bool checkcondition
    SCIP_Real conditionlimit
    COLPACKET * packcstat
    Definition: lpi_clp.cpp:136
    ROWPACKET * packrstat
    Definition: lpi_clp.cpp:137
    definition of wrapper class for rational numbers
    @ 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_LPPAR_PRICING
    Definition: type_lpi.h:54
    @ SCIP_LPPAR_REFACTOR
    Definition: type_lpi.h:71
    @ SCIP_LPPAR_LPINFO
    Definition: type_lpi.h:55
    @ SCIP_LPPAR_POLISHING
    Definition: type_lpi.h:70
    @ 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_FEASTOL
    Definition: type_lpi.h:56
    @ SCIP_LPPAR_LPITLIM
    Definition: type_lpi.h:60
    @ SCIP_LPPAR_ROWREPSWITCH
    Definition: type_lpi.h:63
    @ 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
    @ 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
    type definitions for specific exact LP solvers interface
    @ SCIP_ISFPREPRESENTABLE_UNKNOWN
    Definition: type_rational.h:48
    type definitions for return codes for SCIP methods
    @ SCIP_LPERROR
    Definition: type_retcode.h:49
    @ SCIP_NOFILE
    Definition: type_retcode.h:47
    @ SCIP_READERROR
    Definition: type_retcode.h:45
    @ SCIP_INVALIDDATA
    Definition: type_retcode.h:52
    @ SCIP_PARAMETERUNKNOWN
    Definition: type_retcode.h:55
    @ SCIP_WRITEERROR
    Definition: type_retcode.h:46
    @ SCIP_OKAY
    Definition: type_retcode.h:42
    @ SCIP_NOTIMPLEMENTED
    Definition: type_retcode.h:61
    enum SCIP_Retcode SCIP_RETCODE
    Definition: type_retcode.h:63