Scippy

SCIP

Solving Constraint Integer Programs

lpi_spx2.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-2014 Konrad-Zuse-Zentrum */
7 /* fuer Informationstechnik Berlin */
8 /* */
9 /* SCIP is distributed under the terms of the ZIB Academic License. */
10 /* */
11 /* You should have received a copy of the ZIB Academic License */
12 /* along with SCIP; see the file COPYING. If not email to scip@zib.de. */
13 /* */
14 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
15 
16 /**@file lpi_spx.cpp
17  * @ingroup LPIS
18  * @brief LP interface for SoPlex version 1.4 and higher
19  * @author Tobias Achterberg
20  * @author Timo Berthold
21  * @author Ambros Gleixner
22  * @author Marc Pfetsch
23  *
24  * This is an implementation of SCIP's LP interface for SoPlex. While the ratio test is fixed to SoPlex's standard,
25  * different pricing methods can be chosen and an autopricing strategy (start with devex and switch to steepest edge
26  * after too many iterations) is implemented directly. Scaler and simplifier may be applied if solving from scratch.
27  *
28  * For debugging purposes, the SoPlex results can be double checked with CPLEX if WITH_LPSCHECK is defined. This may
29  * yield false positives, since the LP is dumped to a file for transfering it to CPLEX, hence, precision may be lost.
30  */
31 
32 /*--+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
33 
34 #define STRONGBRANCH_RESTOREBASIS /**< if defined then in SCIPlpiStrongbranch() we restore the basis after the
35  * down branch and after the up branch; if false only after the end of a
36  * strong branching phase, which however seems to mostly increase strong
37  * branching time and iterations */
38 
39 /* in this case the SoPlex results are double checked using CPLEX */
40 #ifdef WITH_LPSCHECK
41 #include <cplex.h>
42 
43 #define CHECK_SPXSOLVE true /**< shall the SoPlex results in spxSolve() be double checked using CPLEX? */
44 #define CHECK_SPXSTRONGBRANCH true /**< shall the SoPlex results in SCIPlpStrongbranch() be double checked using CPLEX? */
45 #define CHECK_START 0 /**< skip first CHECK_START number of checks */
46 #define EXIT_AT_WRONG_RESULT false/**< shall program be exited if CPLEX returns different result than SoPlex? */
47 #define EXIT_AT_CPXERROR false/**< shall program be exited if CPLEX returns an error? */
48 
49 #define CPX_CALL(x) do \
50  { \
51  int _cpxstat_; \
52  if( (_cpxstat_ = (x)) != 0 ) \
53  { \
54  SCIPmessagePrintWarning(_messagehdlr, "CPLEX error <%d>; SoPlex result unchecked\n", _cpxstat_); \
55  if( EXIT_AT_CPXERROR ) \
56  { \
57  exit(1); \
58  } \
59  else \
60  { \
61  goto ENDCHECK; \
62  } \
63  } \
64  } \
65  while( false )
66 #endif
67 
68 /* remember the original value of the SCIP_DEBUG define and undefine it */
69 #ifdef SCIP_DEBUG
70 #define ___DEBUG
71 #undef SCIP_DEBUG
72 #endif
73 
74 /* include SoPlex solver */
75 #include "soplex.h"
76 
77 /* define subversion for versions <= 1.5.0.1 */
78 #ifndef SOPLEX_SUBVERSION
79 #define SOPLEX_SUBVERSION 0
80 #endif
81 
82 /**@todo update this check to version 2.0 */
83 /* check version */
84 #if (SOPLEX_VERSION < 172 || (SOPLEX_VERSION == 172 && SOPLEX_SUBVERSION < 8))
85 #error "This interface is not compatible with SoPlex versions prior to 1.7.2.8"
86 #endif
87 
88 #include "spxgithash.h"
89 
90 /* reset the SCIP_DEBUG define to its original SCIP value */
91 #undef SCIP_DEBUG
92 #ifdef ___DEBUG
93 #define SCIP_DEBUG
94 #undef ___DEBUG
95 #endif
96 
97 #define SOPLEX_VERBLEVEL 5 /**< verbosity level for LPINFO */
98 
99 #include "scip/pub_message.h"
100 
101 /********************************************************************/
102 /*----------------------------- C++ --------------------------------*/
103 /********************************************************************/
104 
105 /* in C++ we have to use "0" instead of "(void*)0" */
106 #undef NULL
107 #define NULL 0
108 
109 #include <cassert>
110 using namespace soplex;
111 
112 
113 /** Macro for a single SoPlex call for which exceptions have to be catched - return an LP error. We
114  * make no distinction between different exception types, e.g., between memory allocation and other
115  * exceptions.
116  */
117 #ifndef NDEBUG
118 #define SOPLEX_TRY(messagehdlr, x) do \
119  { \
120  try \
121  { \
122  (x); \
123  } \
124  catch(SPxException E) \
125  { \
126  std::string s = E.what(); \
127  SCIPmessagePrintWarning((messagehdlr), "SoPlex threw an exception: %s\n", s.c_str()); \
128  return SCIP_LPERROR; \
129  } \
130  } \
131  while( FALSE )
132 
133 #else
134 #define SOPLEX_TRY(messagehdlr, x) do \
135  { \
136  try \
137  { \
138  (x); \
139  } \
140  catch(SPxException E) \
141  { \
142  return SCIP_LPERROR; \
143  } \
144  } \
145  while( FALSE )
146 #endif
147 
148 #define SOPLEX_TRYLPI(x) SOPLEX_TRY(lpi->messagehdlr, x)
149 #define SOPLEX_TRYLPIPTR(x) SOPLEX_TRY((*lpi)->messagehdlr, x)
150 
151 /* Macro for a single SoPlex call for which exceptions have to be catched - abort if they
152  * arise. SCIP_ABORT() is not accessible here.
153  */
154 #define SOPLEX_TRY_ABORT(x) do \
155  { \
156  try \
157  { \
158  (x); \
159  } \
160  catch(SPxException E) \
161  { \
162  std::string s = E.what(); \
163  SCIPerrorMessage("SoPlex threw an exception: %s\n", s.c_str()); \
164  abort(); \
165  } \
166  } \
167  while( FALSE )
168 
169 
170 
171 /** SCIP's SoPlex class */
172 class SPxSCIP : public SoPlex
173 {
174  bool _lpinfo;
175  bool _fromscratch;
176  char* _probname;
177  SPxSolver::VarStatus* _colStat; /**< column basis status used for strong branching */
178  SPxSolver::VarStatus* _rowStat; /**< row basis status used for strong branching */
179 #ifdef WITH_LPSCHECK
180  int _checknum;
181  bool _doublecheck;
182  CPXENVptr _cpxenv; /**< CPLEX memory environment */
183  CPXLPptr _cpxlp; /**< CPLEX lp structure */
184 #endif
185  SCIP_MESSAGEHDLR* _messagehdlr; /**< messagehdlr handler for printing messages, or NULL */
186 
187 public:
188  SPxSCIP(
189  SCIP_MESSAGEHDLR* messagehdlr = NULL, /**< message handler */
190  const char* probname = NULL /**< name of problem */
191  )
192  : _lpinfo(false),
193  _fromscratch(false),
194  _probname(NULL),
195  _colStat(NULL),
196  _rowStat(NULL),
197  _messagehdlr(messagehdlr)
198  {
199  if ( probname != NULL )
200  SOPLEX_TRY_ABORT( setProbname(probname) );
201 
202 #ifdef WITH_LPSCHECK
203  int cpxstat;
204  _checknum = 0;
205  _doublecheck = false;
206  _cpxenv = CPXopenCPLEX(&cpxstat);
207  assert(_cpxenv != NULL);
208  _cpxlp = CPXcreateprob(_cpxenv, &cpxstat, probname != NULL ? probname : "spxcheck");
209  (void) CPXsetintparam(_cpxenv, CPX_PARAM_SCRIND, 0);
210 #endif
211  }
212 
213  virtual ~SPxSCIP()
214  {
215  if( _probname != NULL )
216  spx_free(_probname); /*lint !e1551*/
217 
218  freePreStrongbranchingBasis();
219 
220 #ifdef WITH_LPSCHECK
221  (void) CPXfreeprob(_cpxenv, &_cpxlp);
222  (void) CPXcloseCPLEX(&_cpxenv);
223 #endif
224  }
225 
226  // we might need these methods to return the original values SCIP provided, even if they could not be set
227  /** return feastol set by SCIPlpiSetRealpar(), which might be tighter than what SoPlex accepted */
228  Real feastol()
229  {
230  return realParam(FEASTOL);
231  }
232 
233  /** set feastol and store value in case SoPlex only accepts a larger tolerance */
234  void setFeastol(
235  const Real d
236  )
237  {
238  setRealParam(FEASTOL, d);
239  }
240 
241  /** return opttol set by SCIPlpiSetRealpar(), which might be tighter than what SoPlex accepted */
242  Real opttol()
243  {
244  return realParam(OPTTOL);
245  }
246 
247  /** set opttol and store value in case SoPlex only accepts a larger tolerance */
248  void setOpttol(
249  const Real d
250  )
251  {
252  setRealParam(OPTTOL, d);
253  }
254 
255  /** get objective limit according to objective sense */
256  Real getObjLimit() const
257  {
258  return (intParam(SoPlex::OBJSENSE) == SoPlex::OBJSENSE_MINIMIZE)
259  ? realParam(SoPlex::OBJLIMIT_UPPER)
260  : realParam(SoPlex::OBJLIMIT_LOWER);
261  }
262 
263  // @todo realize this with a member variable as before
264  bool getFromScratch() const
265  {
266  return _fromscratch;
267  }
268 
269  void setFromScratch(bool fs)
270  {
271  _fromscratch = fs;
272  }
273 
274  // @todo member variable?
275  bool getLpInfo() const
276  {
277  return _lpinfo;
278  }
279 
280  void setLpInfo(bool lpinfo)
281  {
282  _lpinfo = lpinfo;
283  }
284 
285  // @todo member variable?
286  void setProbname(const char* probname)
287  {
288  int len;
289 
290  assert(probname != NULL);
291  if( _probname != NULL )
292  spx_free(_probname);
293  len = (int)strlen(probname);
294  spx_alloc(_probname, len + 1);
295  strncpy(_probname, probname, len);
296  _probname[len] = '\0';
297  }
298 
299  void setRep(SPxSolver::Representation p_rep)
300  {
301  if( p_rep == SPxSolver::COLUMN && intParam(REPRESENTATION) == REPRESENTATION_ROW )
302  {
303  SCIPdebugMessage("switching to column representation of the basis\n");
304  setIntParam(REPRESENTATION, REPRESENTATION_COLUMN);
305  }
306  else if( (p_rep == SPxSolver::ROW && intParam(REPRESENTATION) == REPRESENTATION_COLUMN) )
307  {
308  SCIPdebugMessage("switching to row representation of the basis\n");
309  setIntParam(REPRESENTATION, REPRESENTATION_ROW);
310  }
311  }
312 
313 #ifdef WITH_LPSCHECK
314  bool getDoubleCheck()
315  {
316  _checknum++;
317  return _doublecheck && _checknum + 1 >= CHECK_START;
318  }
319 
320  void setDoubleCheck(bool dc)
321  {
322  _doublecheck = dc;
323  }
324 
325  const char* spxStatusString(const SPxSolver::Status stat)
326  {
327  switch( stat )
328  {
329  case SPxSolver::ABORT_TIME:
330  return "ABORT_TIME";
331  case SPxSolver::ABORT_ITER:
332  return "ABORT_ITER";
333  case SPxSolver::ABORT_VALUE:
334  return "ABORT_VALUE";
335  case SPxSolver::SINGULAR:
336  return "SINGULAR";
337  case SPxSolver::REGULAR:
338  return "REGULAR";
339  case SPxSolver::UNKNOWN:
340  return "UNKNOWN";
341  case SPxSolver::OPTIMAL:
342  return "OPTIMAL";
343  case SPxSolver::UNBOUNDED:
344  return "UNBOUNDED";
345  case SPxSolver::INFEASIBLE:
346  return "INFEASIBLE";
347  default:
348  return "UNKNOWN";
349  } /*lint !e788*/
350 
351  return "UNKNOWN";
352  }
353 
354  const char* cpxStatusString(const int stat)
355  {
356  switch( stat )
357  {
358  case CPX_STAT_ABORT_TIME_LIM:
359  return "ABORT_TIME";
360  case CPX_STAT_ABORT_IT_LIM:
361  return "ABORT_ITER";
362  case CPX_STAT_ABORT_OBJ_LIM:
363  return "ABORT_VALUE";
364  case CPX_STAT_OPTIMAL:
365  return "OPTIMAL";
366  case CPX_STAT_OPTIMAL_INFEAS:
367  return "CPX_STAT_OPTIMAL_INFEAS: OPT SOL INFEASIBLE AFTER UNSCALING";
368  case CPX_STAT_UNBOUNDED:
369  return "UNBOUNDED";
370  case CPX_STAT_INFEASIBLE:
371  return "INFEASIBLE";
372  case CPX_STAT_INForUNBD:
373  return "INFEASIBLE or UNBOUNDED";
374  case CPX_STAT_NUM_BEST:
375  return "CPX_STAT_NUM_BEST: SOL AVAILABLE BUT NOT PROVEN OPTIMAL DUE TO NUM TROUBLE";
376  default:
377  return "UNKNOWN";
378  } /*lint !e788*/
379 
380  return "UNKNOWN";
381  }
382 #endif
383 
384 #ifndef NDEBUG
385  bool checkConsistentBounds()
386  {
387  for( int i = 0; i < numColsReal(); ++i )
388  {
389  if( lowerReal(i) > upperReal(i) )
390  {
391  SCIPerrorMessage("inconsistent bounds on column %d: lower=%.17g, upper=%.17g\n",
392  i, lowerReal(i), upperReal(i));
393  return false;
394  }
395  }
396 
397  return true;
398  }
399 
400  bool checkConsistentSides()
401  {
402  for( int i = 0; i < numRowsReal(); ++i )
403  {
404  if( lhsReal(i) > rhsReal(i) )
405  {
406  SCIPerrorMessage("inconsistent sides on row %d: lhs=%.17g, rhs=%.17g\n",
407  i, lhsReal(i), rhsReal(i));
408  return false;
409  }
410  }
411 
412  return true;
413  }
414 #endif
415 
416  void trySolve(bool printwarning = true)
417  {
418  Real timespent;
419  Real timelimit;
420 
421  try
422  {
423  solve();
424  }
425  catch(SPxException x)
426  {
427  std::string s = x.what();
428  if( printwarning )
429  {
430  SCIPmessagePrintWarning(_messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
431  }
432 
433  /* since it is not clear if the status in SoPlex are set correctly
434  * we want to make sure that if an error is thrown the status is
435  * not OPTIMAL anymore.
436  */
437  assert(status() != SPxSolver::OPTIMAL);
438  }
439 
440  assert(intParam(ITERLIMIT) < 0 || numIterations() <= intParam(ITERLIMIT));
441 
442  /* update time limit */
443  timespent = solveTime();
444  if( timespent > 0 )
445  {
446  /* get current time limit */
447  timelimit = realParam(TIMELIMIT);
448  if( timelimit > timespent )
449  timelimit -= timespent;
450  else
451  timelimit = 0;
452  /* set new time limit */
453  assert(timelimit >= 0);
454  setRealParam(TIMELIMIT, timelimit);
455  }
456  }
457 
458  virtual SPxSolver::Status doSolve(bool printwarning = true)
459  {
460  int verbosity;
461 
462  SPxSolver::Status spxStatus;
463 
464  /* store and set verbosity */
465  verbosity = Param::verbose();
466  Param::setVerbose(getLpInfo() ? SOPLEX_VERBLEVEL : 0);
467 
468  assert(checkConsistentBounds());
469  assert(checkConsistentSides());
470 
471 #ifdef WITH_LPSCHECK
472  /* dump LP with current basis and settings saved in SoPlex */
473  if( getDoubleCheck() )
474  writeStateReal("spxcheck", NULL, NULL);
475 #endif
476 
477  trySolve(printwarning);
478  spxStatus = status();
479 
480  /* for safety reset iteration limit */
481 // setTerminationIter(_itlim);
482 
483 #ifdef WITH_LPSCHECK
484  bool minimize = intParam(OBJSENSE) == OBJSENSE_MINIMIZE;
485  Real objLimitUpper = realParam(OBJLIMIT_UPPER);
486  Real objLimitLower = realParam(OBJLIMIT_LOWER);
487 
488  /* if SoPlex gave a definite answer, we double check if it is consistent with CPLEX's answer */
489  if( getDoubleCheck() && (spxStatus == SPxSolver::OPTIMAL || spxStatus == SPxSolver::UNBOUNDED || spxStatus == SPxSolver::INFEASIBLE || spxStatus == SPxSolver::ABORT_VALUE) )
490  {
491  SCIP_Real cpxobj;
492  int cpxstat;
493 
494  /* read LP with basis */
495  CPX_CALL( CPXreadcopyprob(_cpxenv, _cpxlp, "spxcheck.mps", NULL) );
496  CPX_CALL( CPXreadcopybase(_cpxenv, _cpxlp, "spxcheck.bas") );
497 
498  /* set tolerances */
499  CPX_CALL( CPXsetdblparam(_cpxenv, CPX_PARAM_EPOPT, MAX(opttol(), 1e-9)) );
500  CPX_CALL( CPXsetdblparam(_cpxenv, CPX_PARAM_EPRHS, MAX(feastol(), 1e-9)) );
501 
502  /* solve LP */
503  CPX_CALL( CPXlpopt(_cpxenv, _cpxlp) );
504 
505  /* get solution status and objective value */
506  CPX_CALL( CPXsolution(_cpxenv, _cpxlp, &cpxstat, &cpxobj, NULL, NULL, NULL, NULL) );
507  if( !minimize )
508  cpxobj *= -1.0;
509 
510  /* check for inconsistent statuses */
511  if( cpxstat == CPX_STAT_OPTIMAL_INFEAS )
512  {
513  SCIPerrorMessage("In %s: SoPlex status=%d (%s) while CPLEX status=%d (%s)\n",
514  _probname, spxStatus, spxStatusString(spxStatus), cpxstat, cpxStatusString(cpxstat));
515  if( EXIT_AT_CPXERROR )
516  exit(1);
517  }
518  else if( (spxStatus == SPxSolver::OPTIMAL && cpxstat != CPX_STAT_OPTIMAL)
519  || (spxStatus == SPxSolver::UNBOUNDED && cpxstat != CPX_STAT_UNBOUNDED)
520  || (spxStatus == SPxSolver::INFEASIBLE && cpxstat != CPX_STAT_INFEASIBLE) )
521  {
522  SCIPerrorMessage("In %s: SoPlex status=%d (%s) while CPLEX status=%d (%s) (checknum=%d)\n",
523  _probname, spxStatus, spxStatusString(spxStatus), cpxstat, cpxStatusString(cpxstat), _checknum);
524  if( EXIT_AT_WRONG_RESULT )
525  exit(1);
526  }
527  else if( spxStatus == SPxSolver::ABORT_VALUE )
528  {
529  switch( cpxstat )
530  {
531  case CPX_STAT_OPTIMAL:
532  if( (minimize && LTrel(cpxobj, objLimitUpper, 2*opttol()))
533  || (!minimize && GTrel(cpxobj, objLimitLower, 2*opttol())) )
534  {
535  SCIPerrorMessage("In %s: SoPlex returned status=%d (%s) while CPLEX claims obj=%.10f %s %.10f=obj.limit (%s) (checknum=%d)\n",
536  _probname, spxStatus, spxStatusString(spxStatus), cpxobj, minimize ? "<" : ">",
537  minimize ? objLimitUpper : objLimitLower, cpxStatusString(cpxstat), _checknum);
538  if( EXIT_AT_WRONG_RESULT )
539  exit(1);
540  }
541  else if( (minimize && cpxobj < objLimitUpper) || (!minimize && cpxobj > objLimitLower) )
542  {
543  SCIPerrorMessage("In %s: SoPlex returned status=%d (%s) while CPLEX claims obj=%.10f %s %.10f=obj.limit (%s) (checknum=%d)\n",
544  _probname, spxStatus, spxStatusString(spxStatus), cpxobj, minimize? "<" : ">",
545  minimize ? objLimitUpper : objLimitLower, cpxStatusString(cpxstat), _checknum);
546  }
547  break;
548  case CPX_STAT_OPTIMAL_INFEAS:
549  case CPX_STAT_NUM_BEST:
550  if( (minimize && cpxobj < objLimitUpper) || (!minimize && cpxobj > objLimitLower) )
551  {
552  SCIPerrorMessage("In %s: SoPlex returned status=%d (%s) while CPLEX claims obj=%.10f %s %.10f=obj.limit (%s) (checknum=%d)\n",
553  _probname, spxStatus, spxStatusString(spxStatus), cpxobj, minimize ? "<" : ">",
554  minimize ? objLimitUpper : objLimitLower, cpxStatusString(cpxstat), _checknum);
555  }
556  break;
557  case CPX_STAT_INFEASIBLE:
558  break;
559  case CPX_STAT_UNBOUNDED:
560  SCIPerrorMessage("In %s: SoPlex status=%d (%s) while CPLEX status=%d (%s) (checknum=%d)\n",
561  _probname, spxStatus, spxStatusString(spxStatus), cpxstat, cpxStatusString(cpxstat), _checknum);
562  if( EXIT_AT_WRONG_RESULT )
563  exit(1);
564  break;
565  case CPX_STAT_INForUNBD:
566  default:
567  SCIPerrorMessage("In %s: SoPlex status=%d (%s) while CPLEX status=%d (%s) (checknum=%d)\n",
568  _probname, spxStatus, spxStatusString(spxStatus), cpxstat, cpxStatusString(cpxstat), _checknum);
569  break;
570  } /*lint !e788*/
571  }
572  /* check for same objective values */
573  else if( spxStatus == SPxSolver::OPTIMAL )
574  {
575  if( (minimize && LTrel(objValueReal(), cpxobj, 2*opttol()))
576  || (!minimize && GTrel(objValueReal(), cpxobj, 2*opttol())) )
577  {
578  /* SCIPerrorMessage("In %s: LP optimal; SoPlex value=%.10f %s CPLEX value=%.10f too good (checknum=%d)\n", value(),
579  _probname, getSense() == SPxSolver::MINIMIZE ? "<" : ">", cpxobj, _checknum); */
580  }
581  else if( (minimize && GTrel(objValueReal(), cpxobj, 2*opttol()))
582  || (!minimize && LTrel(objValueReal(), cpxobj, 2*opttol())) )
583  {
584  SCIPerrorMessage("In %s: LP optimal; SoPlex value=%.10f %s CPLEX value=%.10f suboptimal (checknum=%d)\n", objValueReal(),
585  _probname, minimize ? ">" : "<", cpxobj, _checknum);
586  if( EXIT_AT_WRONG_RESULT )
587  exit(1);
588  }
589  }
590  }
591 
592  ENDCHECK:
593 #endif
594 
595  /* restore verbosity */
596  Param::setVerbose(verbosity);
597 
598  return spxStatus;
599  }
600 
601  /** save the current basis */
602  void savePreStrongbranchingBasis()
603  {
604  assert(_rowStat == NULL);
605  assert(_colStat == NULL);
606 
607  spx_alloc(_rowStat, numRowsReal());
608  spx_alloc(_colStat, numColsReal());
609 
610  try
611  {
612  getBasis(_rowStat, _colStat);
613  }
614  catch(SPxException x)
615  {
616 #ifndef NDEBUG
617  std::string s = x.what();
618  SCIPmessagePrintWarning(_messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
619 
620  /* since it is not clear if the status in SoPlex are set correctly
621  * we want to make sure that if an error is thrown the status is
622  * not OPTIMAL anymore.
623  */
624  assert(status() != SPxSolver::OPTIMAL);
625 #endif
626  }
627  }
628 
629  /** restore basis */
630  void restorePreStrongbranchingBasis()
631  {
632  assert(_rowStat != NULL);
633  assert(_colStat != NULL);
634 
635  try
636  {
637  setBasis(_rowStat, _colStat);
638  }
639  catch(SPxException x)
640  {
641 #ifndef NDEBUG
642  std::string s = x.what();
643  SCIPmessagePrintWarning(_messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
644 #endif
645  /* since it is not clear if the status in SoPlex are set correctly
646  * we want to make sure that if an error is thrown the status is
647  * not OPTIMAL anymore.
648  */
649  assert(status() != SPxSolver::OPTIMAL);
650  }
651  }
652 
653  /** if basis is in store, delete it without restoring it */
654  void freePreStrongbranchingBasis()
655  {
656  if( _rowStat != NULL )
657  spx_free(_rowStat);
658  if( _colStat != NULL )
659  spx_free(_colStat);
660  }
661 
662  /** is pre-strong-branching basis freed? */
663  bool preStrongbranchingBasisFreed()
664  {
665  return ((_rowStat == NULL ) && (_colStat == NULL));
666  }
667 
668 }; /*lint !e1748*/
669 
670 
671 
672 
673 /********************************************************************/
674 /*----------------------------- C --------------------------------*/
675 /********************************************************************/
676 
677 #include "lpi/lpi.h"
678 #include "scip/bitencode.h"
679 
680 typedef SCIP_DUALPACKET COLPACKET; /* each column needs two bits of information (basic/on_lower/on_upper) */
681 #define COLS_PER_PACKET SCIP_DUALPACKETSIZE
682 typedef SCIP_DUALPACKET ROWPACKET; /* each row needs two bit of information (basic/on_lower/on_upper) */
683 #define ROWS_PER_PACKET SCIP_DUALPACKETSIZE
684 
685 
686 
687 /** LP interface */
688 struct SCIP_LPi
689 {
690  SPxSCIP* spx; /**< our SoPlex implementation */
691  int* cstat; /**< array for storing column basis status */
692  int* rstat; /**< array for storing row basis status */
693  int cstatsize; /**< size of cstat array */
694  int rstatsize; /**< size of rstat array */
695  SCIP_PRICING pricing; /**< current pricing strategy */
696  SCIP_Bool solved; /**< was the current LP solved? */
697  SCIP_Real rowrepswitch; /**< use row representation if number of rows divided by number of columns exceeds this value */
698  SCIP_Real conditionlimit; /**< maximum condition number of LP basis counted as stable (-1.0: no limit) */
699  SCIP_Bool checkcondition; /**< should condition number of LP basis be checked for stability? */
700  SCIP_MESSAGEHDLR* messagehdlr; /**< messagehdlr handler to printing messages, or NULL */
701 };
702 
703 /** LPi state stores basis information */
704 struct SCIP_LPiState
705 {
706  int ncols; /**< number of LP columns */
707  int nrows; /**< number of LP rows */
708  COLPACKET* packcstat; /**< column basis status in compressed form */
709  ROWPACKET* packrstat; /**< row basis status in compressed form */
710 };
711 
712 
713 
714 
715 /*
716  * dynamic memory arrays
717  */
718 
719 /** resizes cstat array to have at least num entries */
720 static
721 SCIP_RETCODE ensureCstatMem(
722  SCIP_LPI* lpi, /**< LP interface structure */
723  int num /**< minimal number of entries in array */
724  )
725 {
726  assert(lpi != NULL);
727 
728  if( num > lpi->cstatsize )
729  {
730  int newsize;
731 
732  newsize = MAX(2*lpi->cstatsize, num);
733  SCIP_ALLOC( BMSreallocMemoryArray(&lpi->cstat, newsize) );
734  lpi->cstatsize = newsize;
735  }
736  assert(num <= lpi->cstatsize);
737 
738  return SCIP_OKAY;
739 }
740 
741 /** resizes rstat array to have at least num entries */
742 static
743 SCIP_RETCODE ensureRstatMem(
744  SCIP_LPI* lpi, /**< LP interface structure */
745  int num /**< minimal number of entries in array */
746  )
747 {
748  assert(lpi != NULL);
749 
750  if( num > lpi->rstatsize )
751  {
752  int newsize;
753 
754  newsize = MAX(2*lpi->rstatsize, num);
755  SCIP_ALLOC( BMSreallocMemoryArray(&lpi->rstat, newsize) );
756  lpi->rstatsize = newsize;
757  }
758  assert(num <= lpi->rstatsize);
759 
760  return SCIP_OKAY;
761 }
762 
763 
764 
765 
766 /*
767  * LPi state methods
768  */
769 
770 /** returns the number of packets needed to store column packet information */
771 static
772 int colpacketNum(
773  int ncols /**< number of columns to store */
774  )
775 {
776  return (ncols+(int)COLS_PER_PACKET-1)/(int)COLS_PER_PACKET;
777 }
778 
779 /** returns the number of packets needed to store row packet information */
780 static
781 int rowpacketNum(
782  int nrows /**< number of rows to store */
783  )
784 {
785  return (nrows+(int)ROWS_PER_PACKET-1)/(int)ROWS_PER_PACKET;
786 }
787 
788 /** store row and column basis status in a packed LPi state object */
789 static
790 void lpistatePack(
791  SCIP_LPISTATE* lpistate, /**< pointer to LPi state data */
792  const int* cstat, /**< basis status of columns in unpacked format */
793  const int* rstat /**< basis status of rows in unpacked format */
794  )
795 {
796  assert(lpistate != NULL);
797  assert(lpistate->packcstat != NULL);
798  assert(lpistate->packrstat != NULL);
799 
800  SCIPencodeDualBit(cstat, lpistate->packcstat, lpistate->ncols);
801  SCIPencodeDualBit(rstat, lpistate->packrstat, lpistate->nrows);
802 }
803 
804 /** unpacks row and column basis status from a packed LPi state object */
805 static
806 void lpistateUnpack(
807  const SCIP_LPISTATE* lpistate, /**< pointer to LPi state data */
808  int* cstat, /**< buffer for storing basis status of columns in unpacked format */
809  int* rstat /**< buffer for storing basis status of rows in unpacked format */
810  )
811 {
812  assert(lpistate != NULL);
813  assert(lpistate->packcstat != NULL);
814  assert(lpistate->packrstat != NULL);
815 
816  SCIPdecodeDualBit(lpistate->packcstat, cstat, lpistate->ncols);
817  SCIPdecodeDualBit(lpistate->packrstat, rstat, lpistate->nrows);
818 }
819 
820 /** creates LPi state information object */
821 static
822 SCIP_RETCODE lpistateCreate(
823  SCIP_LPISTATE** lpistate, /**< pointer to LPi state */
824  BMS_BLKMEM* blkmem, /**< block memory */
825  int ncols, /**< number of columns to store */
826  int nrows /**< number of rows to store */
827  )
828 {
829  assert(lpistate != NULL);
830  assert(blkmem != NULL);
831  assert(ncols >= 0);
832  assert(nrows >= 0);
833 
834  SCIP_ALLOC( BMSallocBlockMemory(blkmem, lpistate) );
835  SCIP_ALLOC( BMSallocBlockMemoryArray(blkmem, &(*lpistate)->packcstat, colpacketNum(ncols)) );
836  SCIP_ALLOC( BMSallocBlockMemoryArray(blkmem, &(*lpistate)->packrstat, rowpacketNum(nrows)) );
837 
838  return SCIP_OKAY;
839 }
840 
841 /** frees LPi state information */
842 static
843 void lpistateFree(
844  SCIP_LPISTATE** lpistate, /**< pointer to LPi state information (like basis information) */
845  BMS_BLKMEM* blkmem /**< block memory */
846  )
847 {
848  assert(blkmem != NULL);
849  assert(lpistate != NULL);
850  assert(*lpistate != NULL);
851 
852  BMSfreeBlockMemoryArray(blkmem, &(*lpistate)->packcstat, colpacketNum((*lpistate)->ncols));
853  BMSfreeBlockMemoryArray(blkmem, &(*lpistate)->packrstat, rowpacketNum((*lpistate)->nrows));
854  BMSfreeBlockMemory(blkmem, lpistate);
855 }
856 
857 
858 
859 
860 /*
861  * local methods
862  */
863 
864 
865 /** marks the current LP to be unsolved */
866 static
867 void invalidateSolution(SCIP_LPI* lpi)
868 {
869  assert(lpi != NULL);
870  lpi->solved = FALSE;
871 }
872 
873 
874 
875 /*
876  * LP Interface Methods
877  */
878 
879 
880 /*
881  * Miscellaneous Methods
882  */
883 
884 static char spxname[100];
885 static char spxdesc[200];
886 
887 /**@name Miscellaneous Methods */
888 /**@{ */
889 
890 /** gets name and version of LP solver */
892  void
893  )
894 {
895  SCIPdebugMessage("calling SCIPlpiGetSolverName()\n");
896 
897 #if (SOPLEX_SUBVERSION > 0)
898  sprintf(spxname, "SoPlex2 %d.%d.%d.%d", SOPLEX_VERSION/100, (SOPLEX_VERSION % 100)/10, SOPLEX_VERSION % 10, SOPLEX_SUBVERSION);
899 #else
900  sprintf(spxname, "SoPlex2 %d.%d.%d", SOPLEX_VERSION/100, (SOPLEX_VERSION % 100)/10, SOPLEX_VERSION % 10);
901 #endif
902  return spxname;
903 }
904 
905 /** gets description of LP solver (developer, webpage, ...) */
907  void
908  )
909 {
910  sprintf(spxdesc, "%s", "Linear Programming Solver developed at Zuse Institute Berlin (soplex.zib.de)");
911  sprintf(spxdesc, "%s [GitHash: %s]", spxdesc, getGitHash());
912 #ifdef WITH_LPSCHECK
913  sprintf(spxdesc, "%s %s", spxdesc, "- including CPLEX double check");
914 #endif
915  return spxdesc;
916 }
917 
918 /** gets pointer for LP solver - use only with great care */
920  SCIP_LPI* lpi /**< pointer to an LP interface structure */
921  )
922 {
923  return (void*) lpi->spx;
924 }
925 /**@} */
926 
927 
928 
929 
930 /*
931  * LPI Creation and Destruction Methods
932  */
933 
934 /**@name LPI Creation and Destruction Methods */
935 /**@{ */
936 
937 /** creates an LP problem object */
939  SCIP_LPI** lpi, /**< pointer to an LP interface structure */
940  SCIP_MESSAGEHDLR* messagehdlr, /**< message handler to use for printing messages, or NULL */
941  const char* name, /**< problem name */
942  SCIP_OBJSEN objsen /**< objective sense */
943  )
944 {
945  assert(lpi != NULL);
946 
947  /* create SoPlex object */
948  SCIP_ALLOC( BMSallocMemory(lpi) );
949 
950  /* we use this construction to allocate the memory for the SoPlex class also via the blockmemshell */
951  (*lpi)->spx = static_cast<SPxSCIP*>(BMSallocMemoryCPP(sizeof(SPxSCIP)));
952  SOPLEX_TRY( messagehdlr, (*lpi)->spx = new ((*lpi)->spx) SPxSCIP(messagehdlr, name) );
953  (*lpi)->cstat = NULL;
954  (*lpi)->rstat = NULL;
955  (*lpi)->cstatsize = 0;
956  (*lpi)->rstatsize = 0;
957  (*lpi)->pricing = SCIP_PRICING_LPIDEFAULT;
958  (*lpi)->rowrepswitch = SCIPlpiInfinity(*lpi);
959  (*lpi)->conditionlimit = -1.0;
960  (*lpi)->checkcondition = FALSE;
961  (*lpi)->messagehdlr = messagehdlr;
962 
963  invalidateSolution(*lpi);
964 
965  /* set objective sense */
966  SCIP_CALL( SCIPlpiChgObjsen(*lpi, objsen) );
967 
968  /* set default pricing */
969  SCIP_CALL( SCIPlpiSetIntpar(*lpi, SCIP_LPPAR_PRICING, (int)(*lpi)->pricing) );
970 
971  {
972  int verbosity = Param::verbose();
973  Param::setVerbose((*lpi)->spx->getLpInfo() ? SOPLEX_VERBLEVEL : 0);
974  (*lpi)->spx->printVersion();
975  Param::setVerbose(verbosity);
976  }
977 
978  return SCIP_OKAY;
979 }
980 
981 /** deletes an LP problem object */
983  SCIP_LPI** lpi /**< pointer to an LP interface structure */
984  )
985 {
986  assert(lpi != NULL);
987  assert(*lpi != NULL);
988  assert((*lpi)->spx != NULL);
989 
990  /* free LP using destructor and free memory via blockmemshell */
991  (*lpi)->spx->~SPxSCIP();
992  BMSfreeMemory(&((*lpi)->spx));
993 
994  /* free memory */
995  BMSfreeMemoryArrayNull(&(*lpi)->cstat);
996  BMSfreeMemoryArrayNull(&(*lpi)->rstat);
997  BMSfreeMemory(lpi);
998 
999  return SCIP_OKAY;
1000 }
1001 
1002 /**@} */
1003 
1004 
1005 
1006 
1007 /*
1008  * Modification Methods
1009  */
1010 
1011 /**@name Modification Methods */
1012 /**@{ */
1013 
1014 /** copies LP data with column matrix into LP solver */
1016  SCIP_LPI* lpi, /**< LP interface structure */
1017  SCIP_OBJSEN objsen, /**< objective sense */
1018  int ncols, /**< number of columns */
1019  const SCIP_Real* obj, /**< objective function values of columns */
1020  const SCIP_Real* lb, /**< lower bounds of columns */
1021  const SCIP_Real* ub, /**< upper bounds of columns */
1022  char** colnames, /**< column names, or NULL */
1023  int nrows, /**< number of rows */
1024  const SCIP_Real* lhs, /**< left hand sides of rows */
1025  const SCIP_Real* rhs, /**< right hand sides of rows */
1026  char** /*rownames*/, /**< row names, or NULL */
1027  int nnonz, /**< number of nonzero elements in the constraint matrix */
1028  const int* beg, /**< start index of each column in ind- and val-array */
1029  const int* ind, /**< row indices of constraint matrix entries */
1030  const SCIP_Real* val /**< values of constraint matrix entries */
1031  )
1032 {
1033  SCIPdebugMessage("calling SCIPlpiLoadColLP()\n");
1034 
1035  assert(lpi != NULL);
1036  assert(lpi->spx != NULL);
1037  assert(lhs != NULL);
1038  assert(rhs != NULL);
1039 
1040  invalidateSolution(lpi);
1041  assert(lpi->spx->preStrongbranchingBasisFreed());
1042 
1043  try
1044  {
1045  SPxSCIP* spx = lpi->spx;
1046  LPRowSet rows(nrows);
1047  DSVector emptyVector(0);
1048  int i;
1049 
1050  spx->clearLPReal();
1051 
1052  /* set objective sense */
1053  spx->setIntParam(SoPlex::OBJSENSE, (objsen == SCIP_OBJSEN_MINIMIZE ? SoPlex::OBJSENSE_MINIMIZE : SoPlex::OBJSENSE_MAXIMIZE));
1054 
1055  /* create empty rows with given sides */
1056  for( i = 0; i < nrows; ++i )
1057  rows.add(lhs[i], emptyVector, rhs[i]);
1058  spx->addRowsReal(rows);
1059 
1060  /* create column vectors with coefficients and bounds */
1061  SCIP_CALL( SCIPlpiAddCols(lpi, ncols, obj, lb, ub, colnames, nnonz, beg, ind, val) );
1062  }
1063  catch(SPxException x)
1064  {
1065 #ifndef NDEBUG
1066  std::string s = x.what();
1067  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
1068 #endif
1069  return SCIP_LPERROR;
1070  }
1071 
1072  return SCIP_OKAY;
1073 }
1074 
1075 /** adds columns to the LP */
1077  SCIP_LPI* lpi, /**< LP interface structure */
1078  int ncols, /**< number of columns to be added */
1079  const SCIP_Real* obj, /**< objective function values of new columns */
1080  const SCIP_Real* lb, /**< lower bounds of new columns */
1081  const SCIP_Real* ub, /**< upper bounds of new columns */
1082  char** /*colnames*/, /**< column names, or NULL */
1083  int nnonz, /**< number of nonzero elements to be added to the constraint matrix */
1084  const int* beg, /**< start index of each column in ind- and val-array, or NULL if nnonz == 0 */
1085  const int* ind, /**< row indices of constraint matrix entries, or NULL if nnonz == 0 */
1086  const SCIP_Real* val /**< values of constraint matrix entries, or NULL if nnonz == 0 */
1087  )
1088 {
1089  SCIPdebugMessage("calling SCIPlpiAddCols()\n");
1090 
1091  assert(lpi != NULL);
1092  assert(lpi->spx != NULL);
1093  assert(obj != NULL);
1094  assert(lb != NULL);
1095  assert(ub != NULL);
1096  assert(nnonz == 0 || beg != NULL);
1097  assert(nnonz == 0 || ind != NULL);
1098  assert(nnonz == 0 || val != NULL);
1099 
1100  invalidateSolution(lpi);
1101 
1102  assert( lpi->spx->preStrongbranchingBasisFreed() );
1103 
1104  SPxSCIP* spx = lpi->spx;
1105  try
1106  {
1107  LPColSet cols(ncols);
1108  DSVector colVector(ncols);
1109  int start;
1110  int last;
1111  int i;
1112 
1113  /* create column vectors with coefficients and bounds */
1114  for( i = 0; i < ncols; ++i )
1115  {
1116  colVector.clear();
1117  if( nnonz > 0 )
1118  {
1119  start = beg[i];
1120  last = (i == ncols-1 ? nnonz : beg[i+1]);
1121  colVector.add( last-start, &ind[start], &val[start] );
1122  }
1123  cols.add(obj[i], lb[i], colVector, ub[i]);
1124  }
1125  spx->addColsReal(cols);
1126  }
1127  catch(SPxException x)
1128  {
1129 #ifndef NDEBUG
1130  std::string s = x.what();
1131  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
1132 #endif
1133  return SCIP_LPERROR;
1134  }
1135 
1136  return SCIP_OKAY;
1137 }
1138 
1139 /** deletes all columns in the given range from LP */
1141  SCIP_LPI* lpi, /**< LP interface structure */
1142  int firstcol, /**< first column to be deleted */
1143  int lastcol /**< last column to be deleted */
1144  )
1145 {
1146  SCIPdebugMessage("calling SCIPlpiDelCols()\n");
1147 
1148  assert(lpi != NULL);
1149  assert(lpi->spx != NULL);
1150  assert(0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsReal());
1151 
1152  invalidateSolution(lpi);
1153 
1154  assert( lpi->spx->preStrongbranchingBasisFreed() );
1155 
1156  SOPLEX_TRY( lpi->messagehdlr, lpi->spx->removeColRangeReal(firstcol, lastcol) );
1157 
1158  return SCIP_OKAY;
1159 }
1160 
1161 /** deletes columns from SCIP_LP; the new position of a column must not be greater that its old position */
1163  SCIP_LPI* lpi, /**< LP interface structure */
1164  int* dstat /**< deletion status of columns
1165  * input: 1 if column should be deleted, 0 if not
1166  * output: new position of column, -1 if column was deleted */
1167  )
1168 {
1169  int ncols;
1170  int i;
1171 
1172  SCIPdebugMessage("calling SCIPlpiDelColset()\n");
1173 
1174  assert(lpi != NULL);
1175  assert(lpi->spx != NULL);
1176 
1177  invalidateSolution(lpi);
1178 
1179  assert( lpi->spx->preStrongbranchingBasisFreed() );
1180 
1181  ncols = lpi->spx->numColsReal();
1182 
1183  /* SoPlex removeCols() method deletes the columns with dstat[i] < 0, so we have to negate the values */
1184  for( i = 0; i < ncols; ++i )
1185  dstat[i] *= -1;
1186 
1187  SOPLEX_TRY( lpi->messagehdlr, lpi->spx->removeColsReal(dstat) );
1188 
1189  return SCIP_OKAY;
1190 }
1191 
1192 /** adds rows to the LP */
1194  SCIP_LPI* lpi, /**< LP interface structure */
1195  int nrows, /**< number of rows to be added */
1196  const SCIP_Real* lhs, /**< left hand sides of new rows */
1197  const SCIP_Real* rhs, /**< right hand sides of new rows */
1198  char** /*rownames*/, /**< row names, or NULL */
1199  int nnonz, /**< number of nonzero elements to be added to the constraint matrix */
1200  const int* beg, /**< start index of each row in ind- and val-array, or NULL if nnonz == 0 */
1201  const int* ind, /**< column indices of constraint matrix entries, or NULL if nnonz == 0 */
1202  const SCIP_Real* val /**< values of constraint matrix entries, or NULL if nnonz == 0 */
1203  )
1204 {
1205  SCIPdebugMessage("calling SCIPlpiAddRows()\n");
1206 
1207  assert(lpi != NULL);
1208  assert(lpi->spx != NULL);
1209  assert(lhs != NULL);
1210  assert(rhs != NULL);
1211  assert(nnonz == 0 || beg != NULL);
1212  assert(nnonz == 0 || ind != NULL);
1213  assert(nnonz == 0 || val != NULL);
1214 
1215  invalidateSolution(lpi);
1216 
1217  assert( lpi->spx->preStrongbranchingBasisFreed() );
1218 
1219  try
1220  {
1221  SPxSCIP* spx = lpi->spx;
1222  LPRowSet rows(nrows);
1223  DSVector rowVector;
1224  int start;
1225  int last;
1226  int i;
1227 
1228  /* create row vectors with given sides */
1229  for( i = 0; i < nrows; ++i )
1230  {
1231  rowVector.clear();
1232  if( nnonz > 0 )
1233  {
1234  start = beg[i];
1235  last = (i == nrows-1 ? nnonz : beg[i+1]);
1236  rowVector.add( last-start, &ind[start], &val[start] );
1237  }
1238  rows.add(lhs[i], rowVector, rhs[i]);
1239  }
1240  spx->addRowsReal(rows);
1241  }
1242  catch(SPxException x)
1243  {
1244 #ifndef NDEBUG
1245  std::string s = x.what();
1246  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
1247 #endif
1248  return SCIP_LPERROR;
1249  }
1250 
1251  return SCIP_OKAY;
1252 }
1253 
1254 /** deletes all rows in the given range from LP */
1256  SCIP_LPI* lpi, /**< LP interface structure */
1257  int firstrow, /**< first row to be deleted */
1258  int lastrow /**< last row to be deleted */
1259  )
1260 {
1261  SCIPdebugMessage("calling SCIPlpiDelRows()\n");
1262 
1263  assert(lpi != NULL);
1264  assert(lpi->spx != NULL);
1265  assert(0 <= firstrow && firstrow <= lastrow && lastrow < lpi->spx->numRowsReal());
1266 
1267  invalidateSolution(lpi);
1268 
1269  assert( lpi->spx->preStrongbranchingBasisFreed() );
1270 
1271  SOPLEX_TRY( lpi->messagehdlr, lpi->spx->removeRowRangeReal(firstrow, lastrow) );
1272 
1273  return SCIP_OKAY;
1274 }
1275 
1276 /** deletes rows from SCIP_LP; the new position of a row must not be greater that its old position */
1278  SCIP_LPI* lpi, /**< LP interface structure */
1279  int* dstat /**< deletion status of rows
1280  * input: 1 if row should be deleted, 0 if not
1281  * output: new position of row, -1 if row was deleted */
1282  )
1283 {
1284  int nrows;
1285  int i;
1286 
1287  SCIPdebugMessage("calling SCIPlpiDelRowset()\n");
1288 
1289  assert(lpi != NULL);
1290  assert(lpi->spx != NULL);
1291 
1292  invalidateSolution(lpi);
1293 
1294  assert( lpi->spx->preStrongbranchingBasisFreed() );
1295 
1296  nrows = lpi->spx->numRowsReal();
1297 
1298  /* SoPlex removeRows() method deletes the rows with dstat[i] < 0, so we have to negate the values */
1299  for( i = 0; i < nrows; ++i )
1300  dstat[i] *= -1;
1301 
1302  SOPLEX_TRY( lpi->messagehdlr, lpi->spx->removeRowsReal(dstat) );
1303 
1304  return SCIP_OKAY;
1305 }
1306 
1307 /** clears the whole LP */
1309  SCIP_LPI* lpi /**< LP interface structure */
1310  )
1311 {
1312  SCIPdebugMessage("calling SCIPlpiClear()\n");
1313 
1314  assert(lpi != NULL);
1315  assert(lpi->spx != NULL);
1316 
1317  invalidateSolution(lpi);
1318 
1319  assert( lpi->spx->preStrongbranchingBasisFreed() );
1320  SOPLEX_TRY( lpi->messagehdlr, lpi->spx->clearLPReal() );
1321 
1322  return SCIP_OKAY;
1323 }
1324 
1325 /** changes lower and upper bounds of columns */
1327  SCIP_LPI* lpi, /**< LP interface structure */
1328  int ncols, /**< number of columns to change bounds for */
1329  const int* ind, /**< column indices */
1330  const SCIP_Real* lb, /**< values for the new lower bounds */
1331  const SCIP_Real* ub /**< values for the new upper bounds */
1332  )
1333 {
1334  int i;
1335 
1336  SCIPdebugMessage("calling SCIPlpiChgBounds()\n");
1337 
1338  assert(lpi != NULL);
1339  assert(lpi->spx != NULL);
1340  assert(ind != NULL);
1341  assert(lb != NULL);
1342  assert(ub != NULL);
1343 
1344  invalidateSolution(lpi);
1345 
1346  assert( lpi->spx->preStrongbranchingBasisFreed() );
1347 
1348  try
1349  {
1350  for( i = 0; i < ncols; ++i )
1351  {
1352  assert(0 <= ind[i] && ind[i] < lpi->spx->numColsReal());
1353  lpi->spx->changeBoundsReal(ind[i], lb[i], ub[i]);
1354  assert(lpi->spx->lowerReal(ind[i]) <= lpi->spx->upperReal(ind[i]));
1355  }
1356  }
1357  catch(SPxException x)
1358  {
1359 #ifndef NDEBUG
1360  std::string s = x.what();
1361  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
1362 #endif
1363  return SCIP_LPERROR;
1364  }
1365 
1366  return SCIP_OKAY;
1367 }
1368 
1369 /** changes left and right hand sides of rows */
1371  SCIP_LPI* lpi, /**< LP interface structure */
1372  int nrows, /**< number of rows to change sides for */
1373  const int* ind, /**< row indices */
1374  const SCIP_Real* lhs, /**< new values for left hand sides */
1375  const SCIP_Real* rhs /**< new values for right hand sides */
1376  )
1377 {
1378  int i;
1379 
1380  SCIPdebugMessage("calling SCIPlpiChgSides()\n");
1381 
1382  assert(lpi != NULL);
1383  assert(lpi->spx != NULL);
1384  assert(ind != NULL);
1385  assert(lhs != NULL);
1386  assert(rhs != NULL);
1387 
1388  invalidateSolution(lpi);
1389 
1390  assert( lpi->spx->preStrongbranchingBasisFreed() );
1391 
1392  try
1393  {
1394  for( i = 0; i < nrows; ++i )
1395  {
1396  assert(0 <= ind[i] && ind[i] < lpi->spx->numRowsReal());
1397  lpi->spx->changeRangeReal(ind[i], lhs[i], rhs[i]);
1398  assert(lpi->spx->lhsReal(ind[i]) <= lpi->spx->rhsReal(ind[i]));
1399  }
1400  }
1401  catch(SPxException x)
1402  {
1403 #ifndef NDEBUG
1404  std::string s = x.what();
1405  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
1406 #endif
1407  return SCIP_LPERROR;
1408  }
1409 
1410  return SCIP_OKAY;
1411 }
1412 
1413 /** changes a single coefficient */
1415  SCIP_LPI* lpi, /**< LP interface structure */
1416  int row, /**< row number of coefficient to change */
1417  int col, /**< column number of coefficient to change */
1418  SCIP_Real newval /**< new value of coefficient */
1419  )
1420 {
1421  SCIPdebugMessage("calling SCIPlpiChgCoef()\n");
1422 
1423  assert(lpi != NULL);
1424  assert(lpi->spx != NULL);
1425  assert(0 <= row && row < lpi->spx->numRowsReal());
1426  assert(0 <= col && col < lpi->spx->numColsReal());
1427 
1428  invalidateSolution(lpi);
1429 
1430  assert( lpi->spx->preStrongbranchingBasisFreed() );
1431 
1432  SOPLEX_TRY( lpi->messagehdlr, lpi->spx->changeElementReal(row, col, newval) );
1433 
1434  return SCIP_OKAY;
1435 }
1436 
1437 /** changes the objective sense */
1439  SCIP_LPI* lpi, /**< LP interface structure */
1440  SCIP_OBJSEN objsen /**< new objective sense */
1441  )
1442 {
1443  SCIPdebugMessage("calling SCIPlpiChgObjsen()\n");
1444 
1445  assert(lpi != NULL);
1446  assert(lpi->spx != NULL);
1447 
1448  invalidateSolution(lpi);
1449 
1450  assert( lpi->spx->preStrongbranchingBasisFreed() );
1451 
1452  SOPLEX_TRY( lpi->messagehdlr, lpi->spx->setIntParam(SoPlex::OBJSENSE, objsen == SCIP_OBJSEN_MINIMIZE ? SoPlex::OBJSENSE_MINIMIZE : SoPlex::OBJSENSE_MAXIMIZE ) );
1453 
1454  return SCIP_OKAY;
1455 }
1456 
1457 /** changes objective values of columns in the LP */
1459  SCIP_LPI* lpi, /**< LP interface structure */
1460  int ncols, /**< number of columns to change objective value for */
1461  int* ind, /**< column indices to change objective value for */
1462  SCIP_Real* obj /**< new objective values for columns */
1463  )
1464 {
1465  int i;
1466 
1467  SCIPdebugMessage("calling SCIPlpiChgObj()\n");
1468 
1469  assert(lpi != NULL);
1470  assert(lpi->spx != NULL);
1471  assert(ind != NULL);
1472  assert(obj != NULL);
1473 
1474  invalidateSolution(lpi);
1475 
1476  assert( lpi->spx->preStrongbranchingBasisFreed() );
1477 
1478  try
1479  {
1480  for( i = 0; i < ncols; ++i )
1481  {
1482  assert(0 <= ind[i] && ind[i] < lpi->spx->numColsReal());
1483  lpi->spx->changeObjReal(ind[i], obj[i]);
1484  }
1485  }
1486  catch(SPxException x)
1487  {
1488 #ifndef NDEBUG
1489  std::string s = x.what();
1490  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
1491 #endif
1492  return SCIP_LPERROR;
1493  }
1494 
1495  return SCIP_OKAY;
1496 }
1497 
1498 /** multiplies a row with a non-zero scalar; for negative scalars, the row's sense is switched accordingly */
1500  SCIP_LPI* lpi, /**< LP interface structure */
1501  int row, /**< row number to scale */
1502  SCIP_Real scaleval /**< scaling multiplier */
1503  )
1504 {
1505  SCIP_Real lhs;
1506  SCIP_Real rhs;
1507 
1508  SCIPdebugMessage("calling SCIPlpiScaleRow()\n");
1509 
1510  assert(lpi != NULL);
1511  assert(lpi->spx != NULL);
1512  assert(scaleval != 0.0);
1513 
1514  try
1515  {
1516  invalidateSolution(lpi);
1517 
1518  assert( lpi->spx->preStrongbranchingBasisFreed() );
1519 
1520  /* get the row vector and the row's sides */
1521  SVector rowvec = lpi->spx->rowVectorReal(row);
1522  lhs = lpi->spx->lhsReal(row);
1523  rhs = lpi->spx->rhsReal(row);
1524 
1525  /* scale the row vector */
1526  rowvec *= scaleval;
1527 
1528  /* adjust the sides */
1529  if( lhs > -lpi->spx->realParam(SoPlex::INFTY) )
1530  lhs *= scaleval;
1531  else if( scaleval < 0.0 )
1532  lhs = lpi->spx->realParam(SoPlex::INFTY);
1533  if( rhs < lpi->spx->realParam(SoPlex::INFTY) )
1534  rhs *= scaleval;
1535  else if( scaleval < 0.0 )
1536  rhs = -lpi->spx->realParam(SoPlex::INFTY);
1537  if( scaleval < 0.0 )
1538  {
1539  SCIP_Real oldlhs = lhs;
1540  lhs = rhs;
1541  rhs = oldlhs;
1542  }
1543 
1544  /* create the new row */
1545  LPRow lprow(lhs, rowvec, rhs);
1546 
1547  /* change the row in the LP */
1548  lpi->spx->changeRowReal(row, lprow);
1549  assert(lpi->spx->lhsReal(row) <= lpi->spx->rhsReal(row));
1550  }
1551  catch(SPxException x)
1552  {
1553 #ifndef NDEBUG
1554  std::string s = x.what();
1555  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
1556 #endif
1557  return SCIP_LPERROR;
1558  }
1559 
1560  return SCIP_OKAY;
1561 }
1562 
1563 /** multiplies a column with a non-zero scalar; the objective value is multiplied with the scalar, and the bounds
1564  * are divided by the scalar; for negative scalars, the column's bounds are switched
1565  */
1567  SCIP_LPI* lpi, /**< LP interface structure */
1568  int col, /**< column number to scale */
1569  SCIP_Real scaleval /**< scaling multiplier */
1570  )
1571 {
1572  SCIP_Real obj;
1573  SCIP_Real lb;
1574  SCIP_Real ub;
1575 
1576  SCIPdebugMessage("calling SCIPlpiScaleCol()\n");
1577 
1578  assert(lpi != NULL);
1579  assert(lpi->spx != NULL);
1580  assert(scaleval != 0.0);
1581 
1582  try
1583  {
1584  invalidateSolution(lpi);
1585 
1586  assert( lpi->spx->preStrongbranchingBasisFreed() );
1587 
1588  /* get the col vector and the col's bounds and objective value */
1589  SVector colvec = lpi->spx->colVectorReal(col);
1590  obj = lpi->spx->objReal(col);
1591  lb = lpi->spx->lowerReal(col);
1592  ub = lpi->spx->upperReal(col);
1593 
1594  /* scale the col vector */
1595  colvec *= scaleval;
1596 
1597  /* scale the objective value */
1598  obj *= scaleval;
1599 
1600  /* adjust the bounds */
1601  if( lb > -lpi->spx->realParam(SoPlex::INFTY) )
1602  lb /= scaleval;
1603  else if( scaleval < 0.0 )
1604  lb = lpi->spx->realParam(SoPlex::INFTY);
1605  if( ub < lpi->spx->realParam(SoPlex::INFTY) )
1606  ub /= scaleval;
1607  else if( scaleval < 0.0 )
1608  ub = -lpi->spx->realParam(SoPlex::INFTY);
1609  if( scaleval < 0.0 )
1610  {
1611  SCIP_Real oldlb = lb;
1612  lb = ub;
1613  ub = oldlb;
1614  }
1615 
1616  /* create the new col (in LPCol's constructor, the upper bound is given first!) */
1617  LPCol lpcol(obj, colvec, ub, lb);
1618 
1619  /* change the col in the LP */
1620  lpi->spx->changeColReal(col, lpcol);
1621  assert(lpi->spx->lowerReal(col) <= lpi->spx->upperReal(col));
1622  }
1623  catch(SPxException x)
1624  {
1625 #ifndef NDEBUG
1626  std::string s = x.what();
1627  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
1628 #endif
1629  return SCIP_LPERROR;
1630  }
1631 
1632  return SCIP_OKAY;
1633 }
1634 
1635 /**@} */
1636 
1637 
1638 
1639 
1640 /*
1641  * Data Accessing Methods
1642  */
1643 
1644 /**@name Data Accessing Methods */
1645 /**@{ */
1646 
1647 /** gets the number of rows in the LP */
1649  SCIP_LPI* lpi, /**< LP interface structure */
1650  int* nrows /**< pointer to store the number of rows */
1651  )
1652 {
1653  SCIPdebugMessage("calling SCIPlpiGetNRows()\n");
1654 
1655  assert(lpi != NULL);
1656  assert(lpi->spx != NULL);
1657  assert(nrows != NULL);
1658 
1659  *nrows = lpi->spx->numRowsReal();
1660 
1661  return SCIP_OKAY;
1662 }
1663 
1664 /** gets the number of columns in the LP */
1666  SCIP_LPI* lpi, /**< LP interface structure */
1667  int* ncols /**< pointer to store the number of cols */
1668  )
1669 {
1670  SCIPdebugMessage("calling SCIPlpiGetNCols()\n");
1671 
1672  assert(lpi != NULL);
1673  assert(lpi->spx != NULL);
1674  assert(ncols != NULL);
1675 
1676  *ncols = lpi->spx->numColsReal();
1677 
1678  return SCIP_OKAY;
1679 }
1680 
1681 /** gets the number of nonzero elements in the LP constraint matrix */
1683  SCIP_LPI* lpi, /**< LP interface structure */
1684  int* nnonz /**< pointer to store the number of nonzeros */
1685  )
1686 {
1687  int i;
1688 
1689  SCIPdebugMessage("calling SCIPlpiGetNNonz()\n");
1690 
1691  assert(lpi != NULL);
1692  assert(lpi->spx != NULL);
1693  assert(nnonz != NULL);
1694 
1695  /* SoPlex has no direct method to return the number of nonzeros, so we have to count them manually */
1696  *nnonz = 0;
1697  if( lpi->spx->numRowsReal() < lpi->spx->numColsReal() )
1698  {
1699  for( i = 0; i < lpi->spx->numRowsReal(); ++i )
1700  (*nnonz) += lpi->spx->rowVectorReal(i).size();
1701  }
1702  else
1703  {
1704  for( i = 0; i < lpi->spx->numColsReal(); ++i )
1705  (*nnonz) += lpi->spx->colVectorReal(i).size();
1706  }
1707 
1708  return SCIP_OKAY;
1709 }
1710 
1711 /** gets columns from LP problem object; the arrays have to be large enough to store all values
1712  * Either both, lb and ub, have to be NULL, or both have to be non-NULL,
1713  * either nnonz, beg, ind, and val have to be NULL, or all of them have to be non-NULL.
1714  */
1716  SCIP_LPI* lpi, /**< LP interface structure */
1717  int firstcol, /**< first column to get from LP */
1718  int lastcol, /**< last column to get from LP */
1719  SCIP_Real* lb, /**< buffer to store the lower bound vector, or NULL */
1720  SCIP_Real* ub, /**< buffer to store the upper bound vector, or NULL */
1721  int* nnonz, /**< pointer to store the number of nonzero elements returned, or NULL */
1722  int* beg, /**< buffer to store start index of each column in ind- and val-array, or NULL */
1723  int* ind, /**< buffer to store column indices of constraint matrix entries, or NULL */
1724  SCIP_Real* val /**< buffer to store values of constraint matrix entries, or NULL */
1725  )
1726 {
1727  int i;
1728  int j;
1729 
1730  SCIPdebugMessage("calling SCIPlpiGetCols()\n");
1731 
1732  assert(lpi != NULL);
1733  assert(lpi->spx != NULL);
1734  assert(0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsReal());
1735 
1736  if( lb != NULL )
1737  {
1738  assert(ub != NULL);
1739 
1740  const Vector& lbvec = lpi->spx->lowerReal();
1741  const Vector& ubvec = lpi->spx->upperReal();
1742  for( i = firstcol; i <= lastcol; ++i )
1743  {
1744  lb[i-firstcol] = lbvec[i];
1745  ub[i-firstcol] = ubvec[i];
1746  }
1747  }
1748  else
1749  assert(ub == NULL);
1750 
1751  if( nnonz != NULL )
1752  {
1753  *nnonz = 0;
1754  for( i = firstcol; i <= lastcol; ++i )
1755  {
1756  beg[i-firstcol] = *nnonz;
1757  const SVector& cvec = lpi->spx->colVectorReal(i);
1758  for( j = 0; j < cvec.size(); ++j )
1759  {
1760  ind[*nnonz] = cvec.index(j);
1761  val[*nnonz] = cvec.value(j);
1762  (*nnonz)++;
1763  }
1764  }
1765  }
1766  else
1767  {
1768  assert(beg == NULL);
1769  assert(ind == NULL);
1770  assert(val == NULL);
1771  }
1772 
1773  return SCIP_OKAY;
1774 }
1775 
1776 /** gets rows from LP problem object; the arrays have to be large enough to store all values.
1777  * Either both, lhs and rhs, have to be NULL, or both have to be non-NULL,
1778  * either nnonz, beg, ind, and val have to be NULL, or all of them have to be non-NULL.
1779  */
1781  SCIP_LPI* lpi, /**< LP interface structure */
1782  int firstrow, /**< first row to get from LP */
1783  int lastrow, /**< last row to get from LP */
1784  SCIP_Real* lhs, /**< buffer to store left hand side vector, or NULL */
1785  SCIP_Real* rhs, /**< buffer to store right hand side vector, or NULL */
1786  int* nnonz, /**< pointer to store the number of nonzero elements returned, or NULL */
1787  int* beg, /**< buffer to store start index of each row in ind- and val-array, or NULL */
1788  int* ind, /**< buffer to store row indices of constraint matrix entries, or NULL */
1789  SCIP_Real* val /**< buffer to store values of constraint matrix entries, or NULL */
1790  )
1791 {
1792  int i;
1793  int j;
1794 
1795  SCIPdebugMessage("calling SCIPlpiGetRows()\n");
1796 
1797  assert(lpi != NULL);
1798  assert(lpi->spx != NULL);
1799  assert(0 <= firstrow && firstrow <= lastrow && lastrow < lpi->spx->numRowsReal());
1800 
1801  if( lhs != NULL )
1802  {
1803  assert(rhs != NULL);
1804 
1805  const Vector& lhsvec = lpi->spx->lhsReal();
1806  const Vector& rhsvec = lpi->spx->rhsReal();
1807  for( i = firstrow; i <= lastrow; ++i )
1808  {
1809  lhs[i-firstrow] = lhsvec[i];
1810  rhs[i-firstrow] = rhsvec[i];
1811  }
1812  }
1813  else
1814  assert(rhs == NULL);
1815 
1816  if( nnonz != NULL )
1817  {
1818  *nnonz = 0;
1819  for( i = firstrow; i <= lastrow; ++i )
1820  {
1821  beg[i-firstrow] = *nnonz;
1822  const SVector& rvec = lpi->spx->rowVectorReal(i);
1823  for( j = 0; j < rvec.size(); ++j )
1824  {
1825  ind[*nnonz] = rvec.index(j);
1826  val[*nnonz] = rvec.value(j);
1827  (*nnonz)++;
1828  }
1829  }
1830  }
1831  else
1832  {
1833  assert(beg == NULL);
1834  assert(ind == NULL);
1835  assert(val == NULL);
1836  }
1837 
1838  return SCIP_OKAY;
1839 }
1840 
1841 /** gets column names */
1843  SCIP_LPI* lpi, /**< LP interface structure */
1844  int firstcol, /**< first column to get name from LP */
1845  int lastcol, /**< last column to get name from LP */
1846  char** colnames, /**< pointers to column names (of size at least lastcol-firstcol+1) */
1847  char* namestorage, /**< storage for col names */
1848  int namestoragesize, /**< size of namestorage (if 0, storageleft returns the storage needed) */
1849  int* storageleft /**< amount of storage left (if < 0 the namestorage was not big enough) */
1850  )
1851 {
1852  assert( lpi != NULL );
1853  assert( lpi->spx != NULL );
1854  assert( colnames != NULL || namestoragesize == 0 );
1855  assert( namestorage != NULL || namestoragesize == 0 );
1856  assert( namestoragesize >= 0 );
1857  assert( storageleft != NULL );
1858  assert( 0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsReal() );
1859 
1860  SCIPdebugMessage("getting column names %d to %d\n", firstcol, lastcol);
1861 
1862 // lpi->spx->getColNames(firstcol, lastcol, colnames, namestorage, namestoragesize, storageleft);
1863 
1864  return SCIP_OKAY;
1865 }
1866 
1867 /** gets row names */
1868 extern
1870  SCIP_LPI* lpi, /**< LP interface structure */
1871  int firstrow, /**< first row to get name from LP */
1872  int lastrow, /**< last row to get name from LP */
1873  char** rownames, /**< pointers to row names (of size at least lastrow-firstrow+1) */
1874  char* namestorage, /**< storage for row names */
1875  int namestoragesize, /**< size of namestorage (if 0, -storageleft returns the storage needed) */
1876  int* storageleft /**< amount of storage left (if < 0 the namestorage was not big enough) */
1877  )
1878 {
1879  assert( lpi != NULL );
1880  assert( lpi->spx != NULL );
1881  assert( rownames != NULL || namestoragesize == 0 );
1882  assert( namestorage != NULL || namestoragesize == 0 );
1883  assert( namestoragesize >= 0 );
1884  assert( storageleft != NULL );
1885  assert( 0 <= firstrow && firstrow <= lastrow && lastrow < lpi->spx->numRowsReal() );
1886 
1887  SCIPdebugMessage("getting row names %d to %d\n", firstrow, lastrow);
1888 
1889 // lpi->spx->getRowNames(firstrow, lastrow, rownames, namestorage, namestoragesize, storageleft);
1890 
1891  return SCIP_OKAY;
1892 }
1893 
1894 /** gets objective sense of the LP */
1896  SCIP_LPI* lpi, /**< LP interface structure */
1897  SCIP_OBJSEN* objsen /**< pointer to store objective sense */
1898  )
1899 {
1900  SCIPdebugMessage("calling SCIPlpiGetObjsen()\n");
1901 
1902  assert(lpi != NULL);
1903  assert(lpi->spx != NULL);
1904  assert(objsen != NULL);
1905 
1906  *objsen = (lpi->spx->intParam(SoPlex::OBJSENSE) == SoPlex::OBJSENSE_MINIMIZE) ? SCIP_OBJSEN_MINIMIZE : SCIP_OBJSEN_MAXIMIZE;
1907 
1908  return SCIP_OKAY;
1909 }
1910 
1911 /** gets objective coefficients from LP problem object */
1913  SCIP_LPI* lpi, /**< LP interface structure */
1914  int firstcol, /**< first column to get objective coefficient for */
1915  int lastcol, /**< last column to get objective coefficient for */
1916  SCIP_Real* vals /**< array to store objective coefficients */
1917  )
1918 {
1919  int i;
1920 
1921  SCIPdebugMessage("calling SCIPlpiGetObj()\n");
1922 
1923  assert(lpi != NULL);
1924  assert(lpi->spx != NULL);
1925  assert(0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsReal());
1926  assert(vals != NULL);
1927 
1928  for( i = firstcol; i <= lastcol; ++i )
1929  vals[i-firstcol] = lpi->spx->objReal(i);
1930 
1931  return SCIP_OKAY;
1932 }
1933 
1934 /** gets current bounds from LP problem object */
1936  SCIP_LPI* lpi, /**< LP interface structure */
1937  int firstcol, /**< first column to get objective value for */
1938  int lastcol, /**< last column to get objective value for */
1939  SCIP_Real* lbs, /**< array to store lower bound values, or NULL */
1940  SCIP_Real* ubs /**< array to store upper bound values, or NULL */
1941  )
1942 {
1943  int i;
1944 
1945  SCIPdebugMessage("calling SCIPlpiGetBounds()\n");
1946 
1947  assert(lpi != NULL);
1948  assert(lpi->spx != NULL);
1949  assert(0 <= firstcol && firstcol <= lastcol && lastcol < lpi->spx->numColsReal());
1950 
1951  for( i = firstcol; i <= lastcol; ++i )
1952  {
1953  if( lbs != NULL )
1954  lbs[i-firstcol] = lpi->spx->lowerReal(i);
1955  if( ubs != NULL )
1956  ubs[i-firstcol] = lpi->spx->upperReal(i);
1957  }
1958 
1959  return SCIP_OKAY;
1960 }
1961 
1962 /** gets current row sides from LP problem object */
1964  SCIP_LPI* lpi, /**< LP interface structure */
1965  int firstrow, /**< first row to get sides for */
1966  int lastrow, /**< last row to get sides for */
1967  SCIP_Real* lhss, /**< array to store left hand side values, or NULL */
1968  SCIP_Real* rhss /**< array to store right hand side values, or NULL */
1969  )
1970 {
1971  int i;
1972 
1973  SCIPdebugMessage("calling SCIPlpiGetSides()\n");
1974 
1975  assert(lpi != NULL);
1976  assert(lpi->spx != NULL);
1977  assert(0 <= firstrow && firstrow <= lastrow && lastrow < lpi->spx->numRowsReal());
1978 
1979  for( i = firstrow; i <= lastrow; ++i )
1980  {
1981  if( lhss != NULL )
1982  lhss[i-firstrow] = lpi->spx->lhsReal(i);
1983  if( rhss != NULL )
1984  rhss[i-firstrow] = lpi->spx->rhsReal(i);
1985  }
1986 
1987  return SCIP_OKAY;
1988 }
1989 
1990 /** gets a single coefficient */
1992  SCIP_LPI* lpi, /**< LP interface structure */
1993  int row, /**< row number of coefficient */
1994  int col, /**< column number of coefficient */
1995  SCIP_Real* val /**< pointer to store the value of the coefficient */
1996  )
1997 {
1998  SCIPdebugMessage("calling SCIPlpiGetCoef()\n");
1999 
2000  assert(lpi != NULL);
2001  assert(lpi->spx != NULL);
2002  assert(0 <= col && col < lpi->spx->numColsReal());
2003  assert(0 <= row && row < lpi->spx->numRowsReal());
2004  assert(val != NULL);
2005 
2006  *val = lpi->spx->colVectorReal(col)[row];
2007 
2008  return SCIP_OKAY;
2009 }
2010 
2011 /**@} */
2012 
2013 
2014 
2015 
2016 /*
2017  * Solving Methods
2018  */
2019 
2020 /**@name Solving Methods */
2021 /**@{ */
2022 
2023 /** solves LP -- used for both, primal and dual simplex, because SoPlex doesn't distinct the two cases */
2024 static
2025 SCIP_RETCODE spxSolve(
2026  SCIP_LPI* lpi, /**< LP interface structure */
2027  SPxSolver::Representation rep, /**< basis representation */
2028  SPxSolver::Type type /**< algorithm type */
2029  )
2030 {
2031  assert( lpi != NULL );
2032  assert( lpi->spx != NULL );
2033  assert( rep == SPxSolver::ROW || rep == SPxSolver::COLUMN );
2034  assert( type == SPxSolver::ENTER || type == SPxSolver::LEAVE );
2035 
2036  int verbosity;
2037  /* store and set verbosity */
2038  verbosity = Param::verbose();
2039  Param::setVerbose(lpi->spx->getLpInfo() ? SOPLEX_VERBLEVEL : 0);
2040 
2041  SCIPdebugMessage("calling SoPlex solve(): %d cols, %d rows\n", lpi->spx->numColsReal(), lpi->spx->numRowsReal());
2042 
2043  invalidateSolution(lpi);
2044 
2045  assert( lpi->spx->preStrongbranchingBasisFreed() );
2046 
2047  /* set basis representation and algorithm type */
2048  if( rep == SPxSolver::COLUMN )
2049  lpi->spx->setIntParam(SoPlex::REPRESENTATION, SoPlex::REPRESENTATION_COLUMN);
2050  else
2051  lpi->spx->setIntParam(SoPlex::REPRESENTATION, SoPlex::REPRESENTATION_ROW);
2052  if( type == SPxSolver::LEAVE )
2053  lpi->spx->setIntParam(SoPlex::ALGORITHM, SoPlex::ALGORITHM_LEAVE);
2054  else
2055  lpi->spx->setIntParam(SoPlex::ALGORITHM, SoPlex::ALGORITHM_ENTER);
2056 
2057 #ifdef WITH_LPSCHECK
2058  lpi->spx->setDoubleCheck(CHECK_SPXSOLVE);
2059 #endif
2060 
2061  SPxSolver::Status status = lpi->spx->doSolve();
2062  SCIPdebugMessage(" -> SoPlex status: %d, basis status: %d\n", lpi->spx->status(), lpi->spx->basisStatus());
2063  lpi->solved = TRUE;
2064 
2065  /* restore verbosity */
2066  Param::setVerbose(verbosity);
2067 
2068  switch( status )
2069  {
2070  case SPxSolver::ABORT_TIME:
2071  case SPxSolver::ABORT_ITER:
2072  case SPxSolver::ABORT_VALUE:
2073  case SPxSolver::SINGULAR:
2074  case SPxSolver::REGULAR:
2075  case SPxSolver::UNKNOWN:
2076  case SPxSolver::OPTIMAL:
2077  case SPxSolver::UNBOUNDED:
2078  case SPxSolver::INFEASIBLE:
2079  return SCIP_OKAY;
2080  default:
2081  return SCIP_LPERROR;
2082  } /*lint !e788*/
2083 }
2084 
2085 /** calls primal simplex to solve the LP */
2087  SCIP_LPI* lpi /**< LP interface structure */
2088  )
2089 {
2090  SCIPdebugMessage("calling SCIPlpiSolvePrimal()\n");
2091 
2092  SCIP_RETCODE retcode;
2093  SCIP_Bool rowrep;
2094 
2095  assert(lpi != NULL);
2096  assert(lpi->spx != NULL);
2097 
2098  /* first decide if we want to switch the basis representation; in order to avoid oscillatory behaviour, we add the
2099  factor 1.1 for switching back to column representation */
2100  if( lpi->rowrepswitch >= 0 )
2101  {
2102  rowrep = lpi->spx->intParam(SoPlex::REPRESENTATION) == SoPlex::REPRESENTATION_ROW;
2103 
2104  if( !rowrep )
2105  rowrep = lpi->spx->numRowsReal() > lpi->spx->numColsReal() * (lpi->rowrepswitch);
2106  else
2107  rowrep = lpi->spx->numRowsReal() * 1.1 > lpi->spx->numColsReal() * (lpi->rowrepswitch);
2108  }
2109  else
2110  rowrep = FALSE;
2111 
2112  /* SoPlex doesn't distinct between the primal and dual simplex; however
2113  * we can force SoPlex to start with the desired method:
2114  * If the representation is COLUMN:
2115  * - ENTER = PRIMAL
2116  * - LEAVE = DUAL
2117  *
2118  * If the representation is ROW:
2119  * - ENTER = DUAL
2120  * - LEAVE = PRIMAL
2121  */
2122  retcode = rowrep ? spxSolve(lpi, SPxSolver::ROW, SPxSolver::LEAVE) : spxSolve(lpi, SPxSolver::COLUMN, SPxSolver::ENTER);
2123  assert(!rowrep || lpi->spx->intParam(SoPlex::REPRESENTATION) == SoPlex::REPRESENTATION_ROW);
2124  assert(rowrep || lpi->spx->intParam(SoPlex::REPRESENTATION) == SoPlex::REPRESENTATION_COLUMN);
2125 
2126  return retcode;
2127 }
2128 
2129 /** calls dual simplex to solve the LP */
2131  SCIP_LPI* lpi /**< LP interface structure */
2132  )
2133 {
2134  SCIPdebugMessage("calling SCIPlpiSolveDual()\n");
2135 
2136  SCIP_RETCODE retcode;
2137  SCIP_Bool rowrep;
2138 
2139  assert(lpi != NULL);
2140  assert(lpi->spx != NULL);
2141 
2142  /* first decide if we want to switch the basis representation; in order to avoid oscillatory behaviour, we add the
2143  factor 1.1 for switching back to column representation */
2144  if( lpi->rowrepswitch >= 0 )
2145  {
2146  rowrep = lpi->spx->intParam(SoPlex::REPRESENTATION) == SoPlex::REPRESENTATION_ROW;
2147 
2148  if( !rowrep )
2149  rowrep = lpi->spx->numRowsReal() > lpi->spx->numColsReal() * (lpi->rowrepswitch);
2150  else
2151  rowrep = lpi->spx->numRowsReal() * 1.1 > lpi->spx->numColsReal() * (lpi->rowrepswitch);
2152  }
2153  else
2154  rowrep = FALSE;
2155 
2156  /* SoPlex doesn't distinct between the primal and dual simplex; however
2157  * we can force SoPlex to start with the desired method:
2158  * If the representation is COLUMN:
2159  * - ENTER = PRIMAL
2160  * - LEAVE = DUAL
2161  *
2162  * If the representation is ROW:
2163  * - ENTER = DUAL
2164  * - LEAVE = PRIMAL
2165  */
2166  retcode = rowrep ? spxSolve(lpi, SPxSolver::ROW, SPxSolver::ENTER) : spxSolve(lpi, SPxSolver::COLUMN, SPxSolver::LEAVE);
2167  assert(!rowrep || lpi->spx->intParam(SoPlex::REPRESENTATION) == SoPlex::REPRESENTATION_ROW);
2168  assert(rowrep || lpi->spx->intParam(SoPlex::REPRESENTATION) == SoPlex::REPRESENTATION_COLUMN);
2169 
2170  return retcode;
2171 }
2172 
2173 /** calls barrier or interior point algorithm to solve the LP with crossover to simplex basis */
2175  SCIP_LPI* lpi, /**< LP interface structure */
2176  SCIP_Bool crossover /**< perform crossover */
2177  )
2178 { /*lint --e{715}*/
2179  SCIPdebugMessage("calling SCIPlpiSolveBarrier()\n");
2180 
2181  /* Since SoPlex does not support barrier we switch to DUAL */
2182  return SCIPlpiSolveDual(lpi);
2183 }
2184 
2185 /** start strong branching - call before any strongbranching */
2187  SCIP_LPI* lpi /**< LP interface structure */
2188  )
2189 {
2190  assert( lpi->spx->preStrongbranchingBasisFreed() );
2191  lpi->spx->savePreStrongbranchingBasis();
2192 
2193  return SCIP_OKAY;
2194 }
2195 
2196 /** end strong branching - call after any strongbranching */
2198  SCIP_LPI* lpi /**< LP interface structure */
2199  )
2200 {
2201  assert( ! lpi->spx->preStrongbranchingBasisFreed() );
2202  lpi->spx->restorePreStrongbranchingBasis();
2203  lpi->spx->freePreStrongbranchingBasis();
2204 
2205  return SCIP_OKAY;
2206 }
2207 
2208 /** performs strong branching iterations on one arbitrary candidate */
2209 static
2210 SCIP_RETCODE lpiStrongbranch(
2211  SCIP_LPI* lpi, /**< LP interface structure */
2212  int col, /**< column to apply strong branching on */
2213  SCIP_Real psol, /**< current primal solution value of column */
2214  int itlim, /**< iteration limit for strong branchings */
2215  SCIP_Real* down, /**< stores dual bound after branching column down */
2216  SCIP_Real* up, /**< stores dual bound after branching column up */
2217  SCIP_Bool* downvalid, /**< stores whether the returned down value is a valid dual bound;
2218  * otherwise, it can only be used as an estimate value */
2219  SCIP_Bool* upvalid, /**< stores whether the returned up value is a valid dual bound;
2220  * otherwise, it can only be used as an estimate value */
2221  int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
2222  )
2223 {
2224  SPxSCIP* spx;
2225  SPxSolver::Status status;
2226  SCIP_Real oldlb;
2227  SCIP_Real oldub;
2228  SCIP_Real newlb;
2229  SCIP_Real newub;
2230  bool fromparentbasis;
2231  bool error;
2232  int oldItlim;
2233  int verbosity;
2234 
2235  /* store and set verbosity */
2236  verbosity = Param::verbose();
2237  Param::setVerbose(lpi->spx->getLpInfo() ? SOPLEX_VERBLEVEL : 0);
2238 
2239  SCIPdebugMessage("calling SCIPlpiStrongbranch() on variable %d (%d iterations)\n", col, itlim);
2240 
2241  assert(lpi != NULL);
2242  assert(lpi->spx != NULL);
2243  /* assert(down != NULL);
2244  * assert(up != NULL); temporary hack for cloud branching */
2245  assert(downvalid != NULL);
2246  assert(upvalid != NULL);
2247 
2248  spx = lpi->spx;
2249  status = SPxSolver::UNKNOWN;
2250  fromparentbasis = false;
2251  error = false;
2252  oldItlim = spx->intParam(SoPlex::ITERLIMIT);
2253 
2254  /* get current bounds of column */
2255  oldlb = spx->lowerReal(col);
2256  oldub = spx->upperReal(col);
2257 
2258  *downvalid = FALSE;
2259  *upvalid = FALSE;
2260 
2261  if( iter != NULL )
2262  *iter = 0;
2263 
2264  /* set the algorithm type to use dual simplex */
2265  spx->setIntParam(SoPlex::ALGORITHM, spx->intParam(SoPlex::REPRESENTATION) == SoPlex::REPRESENTATION_ROW
2266  ? SoPlex::ALGORITHM_ENTER : SoPlex::ALGORITHM_LEAVE);
2267 
2268  /* down branch */
2269  newub = EPSCEIL(psol-1.0, lpi->spx->feastol());
2270  if( newub >= oldlb - 0.5 && down != NULL )
2271  {
2272  SCIPdebugMessage("strong branching down on x%d (%g) with %d iterations\n", col, psol, itlim);
2273 
2274  spx->changeUpperReal(col, newub);
2275  assert(spx->lowerReal(col) <= spx->upperReal(col));
2276 
2277  spx->setIntParam(SoPlex::ITERLIMIT, itlim);
2278  do
2279  {
2280 #ifdef WITH_LPSCHECK
2281  spx->setDoubleCheck(CHECK_SPXSTRONGBRANCH);
2282 #endif
2283  status = spx->solve();
2284  SCIPdebugMessage(" --> Terminate with status %d\n", status);
2285  switch( status )
2286  {
2287  case SPxSolver::OPTIMAL:
2288  *down = spx->objValueReal();
2289  *downvalid = TRUE;
2290  SCIPdebugMessage(" --> Terminate with value %f\n", *down);
2291  break;
2292  case SPxSolver::ABORT_TIME: /* SoPlex does not return a proven dual bound, if it is aborted */
2293  case SPxSolver::ABORT_ITER:
2294  case SPxSolver::ABORT_CYCLING:
2295  *down = spx->objValueReal();
2296  break;
2297  case SPxSolver::ABORT_VALUE:
2298  case SPxSolver::INFEASIBLE:
2299  *down = spx->getObjLimit();
2300  *downvalid = TRUE;
2301  break;
2302  default:
2303  error = true;
2304  break;
2305  } /*lint !e788*/
2306  if( iter != NULL )
2307  (*iter) += spx->numIterations();
2308 
2309 #ifdef STRONGBRANCH_RESTOREBASIS
2310  /* we restore the pre-strong-branching basis by default (and don't solve again) */
2311  assert( ! spx->preStrongbranchingBasisFreed() );
2312  spx->restorePreStrongbranchingBasis();
2313  fromparentbasis = false;
2314 #else
2315  /* if cycling or singular basis occured and we started not from the pre-strong-branching basis, then we restore the
2316  * pre-strong-branching basis and try again with reduced iteration limit */
2317  if( (status == SPxSolver::ABORT_CYCLING || status == SPxSolver::SINGULAR) && !fromparentbasis && spx->numIterations() < itlim )
2318  {
2319  SCIPdebugMessage(" --> Repeat strong branching down with %d iterations after restoring basis\n",
2320  itlim - spx->numIterations());
2321  spx->setIntParam(SoPlex::ITERLIMIT, itlim - spx->numIterations());
2322  assert( ! spx->hasPreStrongbranchingBasis() );
2323  spx->restorePreStrongbranchingBasis();
2324  fromparentbasis = true;
2325  error = false;
2326  }
2327  /* otherwise don't solve again */
2328  else
2329  fromparentbasis = false;
2330 #endif
2331  }
2332  while( fromparentbasis );
2333 
2334  spx->changeUpperReal(col, oldub);
2335  assert(spx->lowerReal(col) <= spx->upperReal(col));
2336  }
2337  else if( down != NULL )
2338  {
2339  *down = spx->getObjLimit();
2340  *downvalid = TRUE;
2341  }
2342  else
2343  *downvalid = TRUE;
2344 
2345  /* up branch */
2346  if( !error )
2347  {
2348  newlb = EPSFLOOR(psol+1.0, lpi->spx->feastol());
2349  if( newlb <= oldub + 0.5 && up != NULL )
2350  {
2351  SCIPdebugMessage("strong branching up on x%d (%g) with %d iterations\n", col, psol, itlim);
2352 
2353  spx->changeLowerReal(col, newlb);
2354  assert(spx->lowerReal(col) <= spx->upperReal(col));
2355 
2356  spx->setIntParam(SoPlex::ITERLIMIT, itlim);
2357  do
2358  {
2359 #ifdef WITH_LPSCHECK
2360  spx->setDoubleCheck(CHECK_SPXSTRONGBRANCH);
2361 #endif
2362  status = spx->solve();
2363  SCIPdebugMessage(" --> Terminate with status %d\n", status);
2364  switch( status )
2365  {
2366  case SPxSolver::OPTIMAL:
2367  *up = spx->objValueReal();
2368  *upvalid = TRUE;
2369  SCIPdebugMessage(" --> Terminate with value %f\n", spx->objValueReal());
2370  break;
2371  case SPxSolver::ABORT_TIME: /* SoPlex does not return a proven dual bound, if it is aborted */
2372  case SPxSolver::ABORT_ITER:
2373  case SPxSolver::ABORT_CYCLING:
2374  *up = spx->objValueReal();
2375  break;
2376  case SPxSolver::ABORT_VALUE:
2377  case SPxSolver::INFEASIBLE:
2378  *up = spx->getObjLimit();
2379  *upvalid = TRUE;
2380  break;
2381  default:
2382  error = true;
2383  break;
2384  } /*lint !e788*/
2385  if( iter != NULL )
2386  (*iter) += spx->numIterations();
2387 
2388 #ifdef STRONGBRANCH_RESTOREBASIS
2389  /* we restore the pre-strong-branching basis by default (and don't solve again) */
2390  assert( ! spx->preStrongbranchingBasisFreed() );
2391  spx->restorePreStrongbranchingBasis();
2392  fromparentbasis = false;
2393 #else
2394  /* if cycling or singular basis occured and we started not from the pre-strong-branching basis, then we restore the
2395  * pre-strong-branching basis and try again with reduced iteration limit */
2396  else if( (status == SPxSolver::ABORT_CYCLING || status == SPxSolver::SINGULAR) && !fromparentbasis && spx->numIterations() < itlim )
2397  {
2398  SCIPdebugMessage(" --> Repeat strong branching up with %d iterations after restoring basis\n", itlim - spx->numIterations());
2399  assert( ! spx->hasPreStrongbranchingBasis() );
2400  spx->restorePreStrongbranchingBasis();
2401  spx->setIntParam(SoPlex::ITERLIMIT, itlim - spx->numIterations());
2402  error = false;
2403  fromparentbasis = true;
2404  }
2405  /* otherwise don't solve again */
2406  else
2407  fromparentbasis = false;
2408 #endif
2409  }
2410  while( fromparentbasis );
2411 
2412  spx->changeLowerReal(col, oldlb);
2413  assert(spx->lowerReal(col) <= spx->upperReal(col));
2414  }
2415  else if( up != NULL )
2416  {
2417  *up = spx->getObjLimit();
2418  *upvalid = TRUE;
2419  }
2420  else
2421  *upvalid = TRUE;
2422  }
2423 
2424  /* reset old iteration limit */
2425  spx->setIntParam(SoPlex::ITERLIMIT, oldItlim);
2426 
2427  /* restore verbosity */
2428  Param::setVerbose(verbosity);
2429 
2430  if( error )
2431  {
2432  SCIPdebugMessage("SCIPlpiStrongbranch() returned SoPlex status %d\n", int(status));
2433  return SCIP_LPERROR;
2434  }
2435 
2436  return SCIP_OKAY;
2437 }
2438 
2439 /** performs strong branching iterations on one @b fractional candidate */
2441  SCIP_LPI* lpi, /**< LP interface structure */
2442  int col, /**< column to apply strong branching on */
2443  SCIP_Real psol, /**< fractional current primal solution value of column */
2444  int itlim, /**< iteration limit for strong branchings */
2445  SCIP_Real* down, /**< stores dual bound after branching column down */
2446  SCIP_Real* up, /**< stores dual bound after branching column up */
2447  SCIP_Bool* downvalid, /**< stores whether the returned down value is a valid dual bound;
2448  * otherwise, it can only be used as an estimate value */
2449  SCIP_Bool* upvalid, /**< stores whether the returned up value is a valid dual bound;
2450  * otherwise, it can only be used as an estimate value */
2451  int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
2452  )
2453 {
2454  SCIP_RETCODE retcode;
2455 
2456  /* pass call on to lpiStrongbranch() */
2457  retcode = lpiStrongbranch(lpi, col, psol, itlim, down, up, downvalid, upvalid, iter);
2458 
2459  /* pass SCIP_LPERROR to SCIP without a back trace */
2460  if( retcode == SCIP_LPERROR )
2461  return SCIP_LPERROR;
2462 
2463  /* evaluate retcode */
2464  SCIP_CALL( retcode );
2465 
2466  return SCIP_OKAY;
2467 }
2468 
2469 /** performs strong branching iterations on given @b fractional candidates */
2471  SCIP_LPI* lpi, /**< LP interface structure */
2472  int* cols, /**< columns to apply strong branching on */
2473  int ncols, /**< number of columns */
2474  SCIP_Real* psols, /**< fractional current primal solution values of columns */
2475  int itlim, /**< iteration limit for strong branchings */
2476  SCIP_Real* down, /**< stores dual bounds after branching columns down */
2477  SCIP_Real* up, /**< stores dual bounds after branching columns up */
2478  SCIP_Bool* downvalid, /**< stores whether the returned down values are valid dual bounds;
2479  * otherwise, they can only be used as an estimate values */
2480  SCIP_Bool* upvalid, /**< stores whether the returned up values are a valid dual bounds;
2481  * otherwise, they can only be used as an estimate values */
2482  int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
2483  )
2484 {
2485  SCIP_RETCODE retcode;
2486 
2487  assert( iter != NULL );
2488  assert( cols != NULL );
2489  assert( psols != NULL );
2490  assert( down != NULL );
2491  assert( up != NULL );
2492  assert( downvalid != NULL );
2493  assert( upvalid != NULL );
2494  assert( down != NULL );
2495 
2496  if ( iter != NULL )
2497  *iter = 0;
2498 
2499  for (int j = 0; j < ncols; ++j)
2500  {
2501  /* pass call on to lpiStrongbranch() */
2502  retcode = lpiStrongbranch(lpi, cols[j], psols[j], itlim, &(down[j]), &(up[j]), &(downvalid[j]), &(upvalid[j]), iter);
2503 
2504  /* pass SCIP_LPERROR to SCIP without a back trace */
2505  if( retcode == SCIP_LPERROR )
2506  return SCIP_LPERROR;
2507 
2508  /* evaluate retcode */
2509  SCIP_CALL( retcode );
2510  }
2511  return SCIP_OKAY;
2512 }
2513 
2514 /** performs strong branching iterations on one candidate with @b integral value */
2516  SCIP_LPI* lpi, /**< LP interface structure */
2517  int col, /**< column to apply strong branching on */
2518  SCIP_Real psol, /**< current integral primal solution value of column */
2519  int itlim, /**< iteration limit for strong branchings */
2520  SCIP_Real* down, /**< stores dual bound after branching column down */
2521  SCIP_Real* up, /**< stores dual bound after branching column up */
2522  SCIP_Bool* downvalid, /**< stores whether the returned down value is a valid dual bound;
2523  * otherwise, it can only be used as an estimate value */
2524  SCIP_Bool* upvalid, /**< stores whether the returned up value is a valid dual bound;
2525  * otherwise, it can only be used as an estimate value */
2526  int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
2527  )
2528 {
2529  SCIP_RETCODE retcode;
2530 
2531  /* pass call on to lpiStrongbranch() */
2532  retcode = lpiStrongbranch(lpi, col, psol, itlim, down, up, downvalid, upvalid, iter);
2533 
2534  /* pass SCIP_LPERROR to SCIP without a back trace */
2535  if( retcode == SCIP_LPERROR )
2536  return SCIP_LPERROR;
2537 
2538  /* evaluate retcode */
2539  SCIP_CALL( retcode );
2540 
2541  return SCIP_OKAY;
2542 }
2543 
2544 /** performs strong branching iterations on given candidates with @b integral values */
2546  SCIP_LPI* lpi, /**< LP interface structure */
2547  int* cols, /**< columns to apply strong branching on */
2548  int ncols, /**< number of columns */
2549  SCIP_Real* psols, /**< current integral primal solution values of columns */
2550  int itlim, /**< iteration limit for strong branchings */
2551  SCIP_Real* down, /**< stores dual bounds after branching columns down */
2552  SCIP_Real* up, /**< stores dual bounds after branching columns up */
2553  SCIP_Bool* downvalid, /**< stores whether the returned down values are valid dual bounds;
2554  * otherwise, they can only be used as an estimate values */
2555  SCIP_Bool* upvalid, /**< stores whether the returned up values are a valid dual bounds;
2556  * otherwise, they can only be used as an estimate values */
2557  int* iter /**< stores total number of strong branching iterations, or -1; may be NULL */
2558  )
2559 {
2560  SCIP_RETCODE retcode;
2561 
2562  assert( iter != NULL );
2563  assert( cols != NULL );
2564  assert( psols != NULL );
2565  assert( down != NULL );
2566  assert( up != NULL );
2567  assert( downvalid != NULL );
2568  assert( upvalid != NULL );
2569  assert( down != NULL );
2570 
2571  if ( iter != NULL )
2572  *iter = 0;
2573 
2574  for (int j = 0; j < ncols; ++j)
2575  {
2576  /* pass call on to lpiStrongbranch() */
2577  retcode = lpiStrongbranch(lpi, cols[j], psols[j], itlim, &(down[j]), &(up[j]), &(downvalid[j]), &(upvalid[j]), iter);
2578 
2579  /* pass SCIP_LPERROR to SCIP without a back trace */
2580  if( retcode == SCIP_LPERROR )
2581  return SCIP_LPERROR;
2582 
2583  /* evaluate retcode */
2584  SCIP_CALL( retcode );
2585  }
2586 
2587  return SCIP_OKAY;
2588 }
2589 /**@} */
2590 
2591 
2592 
2593 
2594 /*
2595  * Solution Information Methods
2596  */
2597 
2598 /**@name Solution Information Methods */
2599 /**@{ */
2600 
2601 /** returns whether a solve method was called after the last modification of the LP */
2603  SCIP_LPI* lpi /**< LP interface structure */
2604  )
2605 {
2606  assert(lpi != NULL);
2607 
2608  return lpi->solved;
2609 }
2610 
2611 /** gets information about primal and dual feasibility of the current LP solution */
2613  SCIP_LPI* lpi, /**< LP interface structure */
2614  SCIP_Bool* primalfeasible, /**< stores primal feasibility status */
2615  SCIP_Bool* dualfeasible /**< stores dual feasibility status */
2616  )
2617 {
2618  SCIPdebugMessage("calling SCIPlpiGetSolFeasibility()\n");
2619 
2620  assert(lpi != NULL);
2621  assert(primalfeasible != NULL);
2622  assert(dualfeasible != NULL);
2623 
2624  *primalfeasible = SCIPlpiIsPrimalFeasible(lpi);
2625  *dualfeasible = SCIPlpiIsDualFeasible(lpi);
2626 
2627  return SCIP_OKAY;
2628 }
2629 
2630 /** returns TRUE iff LP is proven to have a primal unbounded ray (but not necessary a primal feasible point);
2631  * this does not necessarily mean, that the solver knows and can return the primal ray
2632  */
2634  SCIP_LPI* lpi /**< LP interface structure */
2635  )
2636 {
2637  SCIPdebugMessage("calling SCIPlpiExistsPrimalRay()\n");
2638 
2639  assert(lpi != NULL);
2640  assert(lpi->spx != NULL);
2641 
2642  return (lpi->spx->status() == SPxSolver::UNBOUNDED);
2643 }
2644 
2645 /** returns TRUE iff LP is proven to have a primal unbounded ray (but not necessary a primal feasible point),
2646  * and the solver knows and can return the primal ray
2647  */
2649  SCIP_LPI* lpi /**< LP interface structure */
2650  )
2651 {
2652  SCIPdebugMessage("calling SCIPlpiHasPrimalRay()\n");
2653 
2654  assert(lpi != NULL);
2655  assert(lpi->spx != NULL);
2656 
2657  return (lpi->spx->status() == SPxSolver::UNBOUNDED);
2658 }
2659 
2660 /** returns TRUE iff LP is proven to be primal unbounded */
2662  SCIP_LPI* lpi /**< LP interface structure */
2663  )
2664 {
2665  SCIPdebugMessage("calling SCIPlpiIsPrimalUnbounded()\n");
2666 
2667  assert(lpi != NULL);
2668  assert(lpi->spx != NULL);
2669 
2670  assert(lpi->spx->status() != SPxSolver::UNBOUNDED || lpi->spx->basisStatus() == SPxBasis::UNBOUNDED);
2671 
2672  /* if SoPlex returns unbounded, this may only mean that an unbounded ray is available, not necessarily a primal
2673  * feasible point; hence we have to check the perturbation
2674  */
2675  return lpi->spx->status() == SPxSolver::UNBOUNDED;
2676 }
2677 
2678 /** returns TRUE iff LP is proven to be primal infeasible */
2680  SCIP_LPI* lpi /**< LP interface structure */
2681  )
2682 {
2683  SCIPdebugMessage("calling SCIPlpiIsPrimalInfeasible()\n");
2684 
2685  assert(lpi != NULL);
2686  assert(lpi->spx != NULL);
2687 
2688  return (lpi->spx->status() == SPxSolver::INFEASIBLE);
2689 }
2690 
2691 /** returns TRUE iff LP is proven to be primal feasible */
2693  SCIP_LPI* lpi /**< LP interface structure */
2694  )
2695 {
2696  SPxBasis::SPxStatus basestatus;
2697 
2698  SCIPdebugMessage("calling SCIPlpiIsPrimalFeasible()\n");
2699 
2700  assert(lpi != NULL);
2701  assert(lpi->spx != NULL);
2702 
2703  basestatus = lpi->spx->basisStatus();
2704 
2705  /* note that the solver status may be ABORT_VALUE and the basis status optimal; if we are optimal, isPerturbed() may
2706  * still return true as long as perturbation plus violation is within tolerances
2707  */
2708  assert(basestatus == SPxBasis::OPTIMAL || lpi->spx->status() != SPxSolver::OPTIMAL);
2709 
2710  return basestatus == SPxBasis::OPTIMAL || basestatus == SPxBasis::PRIMAL;
2711 }
2712 
2713 /** returns TRUE iff LP is proven to have a dual unbounded ray (but not necessary a dual feasible point);
2714  * this does not necessarily mean, that the solver knows and can return the dual ray
2715  */
2717  SCIP_LPI* lpi /**< LP interface structure */
2718  )
2719 {
2720  SCIPdebugMessage("calling SCIPlpiExistsDualRay()\n");
2721 
2722  assert(lpi != NULL);
2723  assert(lpi->spx != NULL);
2724 
2725  return (lpi->spx->status() == SPxSolver::INFEASIBLE);
2726 }
2727 
2728 /** returns TRUE iff LP is proven to have a dual unbounded ray (but not necessary a dual feasible point),
2729  * and the solver knows and can return the dual ray
2730  */
2732  SCIP_LPI* lpi /**< LP interface structure */
2733  )
2734 {
2735  SCIPdebugMessage("calling SCIPlpiHasDualRay()\n");
2736 
2737  assert(lpi != NULL);
2738  assert(lpi->spx != NULL);
2739 
2740  return (lpi->spx->status() == SPxSolver::INFEASIBLE);
2741 }
2742 
2743 /** returns TRUE iff LP is dual unbounded */
2745  SCIP_LPI* lpi /**< LP interface structure */
2746  )
2747 {
2748  SCIPdebugMessage("calling SCIPlpiIsDualUnbounded()\n");
2749 
2750  assert(lpi != NULL);
2751  assert(lpi->spx != NULL);
2752 
2753  return lpi->spx->status() == SPxSolver::INFEASIBLE && lpi->spx->basisStatus() == SPxBasis::DUAL;
2754 }
2755 
2756 /** returns TRUE iff LP is dual infeasible */
2758  SCIP_LPI* lpi /**< LP interface structure */
2759  )
2760 {
2761  SCIPdebugMessage("calling SCIPlpiIsDualInfeasible()\n");
2762 
2763  assert(lpi != NULL);
2764  assert(lpi->spx != NULL);
2765 
2766  return (lpi->spx->status() == SPxSolver::UNBOUNDED);
2767 }
2768 
2769 /** returns TRUE iff LP is proven to be dual feasible */
2771  SCIP_LPI* lpi /**< LP interface structure */
2772  )
2773 {
2774  SCIPdebugMessage("calling SCIPlpiIsDualFeasible()\n");
2775 
2776  assert(lpi != NULL);
2777  assert(lpi->spx != NULL);
2778 
2779  /* note that the solver status may be ABORT_VALUE and the basis status optimal; if we are optimal, isPerturbed() may
2780  * still return true as long as perturbation plus violation is within tolerances
2781  */
2782  assert(lpi->spx->basisStatus() == SPxBasis::OPTIMAL || lpi->spx->status() != SPxSolver::OPTIMAL);
2783 
2784  return (lpi->spx->basisStatus() == SPxBasis::OPTIMAL) || lpi->spx->basisStatus() == SPxBasis::DUAL;
2785 }
2786 
2787 /** returns TRUE iff LP was solved to optimality */
2789  SCIP_LPI* lpi /**< LP interface structure */
2790  )
2791 {
2792  SCIPdebugMessage("calling SCIPlpiIsOptimal()\n");
2793 
2794  assert(lpi != NULL);
2795  assert((lpi->spx->basisStatus() == SPxBasis::OPTIMAL)
2797 
2798  /* note that the solver status may be ABORT_VALUE and the basis status optimal; if we are optimal, isPerturbed() may
2799  * still return true as long as perturbation plus violation is within tolerances
2800  */
2801  return (lpi->spx->basisStatus() == SPxBasis::OPTIMAL);
2802 }
2803 
2804 /** returns TRUE iff current LP basis is stable */
2806  SCIP_LPI* lpi /**< LP interface structure */
2807  )
2808 {
2809  SCIPdebugMessage("calling SCIPlpiIsStable()\n");
2810 
2811  assert(lpi != NULL);
2812  assert(lpi->spx != NULL);
2813 
2814  if( lpi->spx->status() == SPxSolver::ERROR || lpi->spx->status() == SPxSolver::SINGULAR )
2815  return FALSE;
2816 
2817  /* only if we have a regular basis and the condition limit is set, we compute the condition number of the basis;
2818  * everything above the specified threshold is then counted as instable
2819  */
2820  if( lpi->checkcondition && (SCIPlpiIsOptimal(lpi) || SCIPlpiIsObjlimExc(lpi)) )
2821  {
2822  SCIP_RETCODE retcode;
2823  SCIP_Real kappa;
2824 
2826  if( retcode != SCIP_OKAY )
2827  {
2828  SCIPABORT();
2829  }
2830  assert(kappa != SCIP_INVALID);
2831 
2832  if( kappa > lpi->conditionlimit )
2833  return FALSE;
2834  }
2835 
2836  return TRUE;
2837 }
2838 
2839 /** returns TRUE iff the objective limit was reached */
2841  SCIP_LPI* lpi /**< LP interface structure */
2842  )
2843 {
2844  SCIPdebugMessage("calling SCIPlpiIsObjlimExc()\n");
2845 
2846  assert(lpi != NULL);
2847  assert(lpi->spx != NULL);
2848 
2849  return (lpi->spx->status() == SPxSolver::ABORT_VALUE);
2850 }
2851 
2852 /** returns TRUE iff the iteration limit was reached */
2854  SCIP_LPI* lpi /**< LP interface structure */
2855  )
2856 {
2857  SCIPdebugMessage("calling SCIPlpiIsIterlimExc()\n");
2858 
2859  assert(lpi != NULL);
2860  assert(lpi->spx != NULL);
2861 
2862  return (lpi->spx->status() == SPxSolver::ABORT_ITER);
2863 }
2864 
2865 /** returns TRUE iff the time limit was reached */
2867  SCIP_LPI* lpi /**< LP interface structure */
2868  )
2869 {
2870  SCIPdebugMessage("calling SCIPlpiIsTimelimExc()\n");
2871 
2872  assert(lpi != NULL);
2873  assert(lpi->spx != NULL);
2874 
2875  return (lpi->spx->status() == SPxSolver::ABORT_TIME);
2876 }
2877 
2878 /** returns the internal solution status of the solver */
2880  SCIP_LPI* lpi /**< LP interface structure */
2881  )
2882 {
2883  SCIPdebugMessage("calling SCIPlpiIsTimelimExc()\n");
2884 
2885  assert(lpi != NULL);
2886  assert(lpi->spx != NULL);
2887 
2888  return static_cast<int>(lpi->spx->status());
2889 }
2890 
2891 /** tries to reset the internal status of the LP solver in order to ignore an instability of the last solving call */
2893  SCIP_LPI* lpi, /**< LP interface structure */
2894  SCIP_Bool* success /**< pointer to store, whether the instability could be ignored */
2895  )
2896 { /*lint --e{715}*/
2897  SCIPdebugMessage("calling SCIPlpiIgnoreInstability()\n");
2898 
2899  assert(lpi != NULL);
2900  assert(lpi->spx != NULL);
2901 
2902  /* instable situations cannot be ignored */
2903  *success = FALSE;
2904 
2905  return SCIP_OKAY;
2906 }
2907 
2908 /** gets objective value of solution */
2910  SCIP_LPI* lpi, /**< LP interface structure */
2911  SCIP_Real* objval /**< stores the objective value */
2912  )
2913 {
2914  SCIPdebugMessage("calling SCIPlpiGetObjval()\n");
2915 
2916  assert(lpi != NULL);
2917  assert(lpi->spx != NULL);
2918  assert(objval != NULL);
2919 
2920  *objval = lpi->spx->objValueReal();
2921 
2922  return SCIP_OKAY;
2923 }
2924 
2925 /** gets primal and dual solution vectors */
2927  SCIP_LPI* lpi, /**< LP interface structure */
2928  SCIP_Real* objval, /**< stores the objective value, may be NULL if not needed */
2929  SCIP_Real* primsol, /**< primal solution vector, may be NULL if not needed */
2930  SCIP_Real* dualsol, /**< dual solution vector, may be NULL if not needed */
2931  SCIP_Real* activity, /**< row activity vector, may be NULL if not needed */
2932  SCIP_Real* redcost /**< reduced cost vector, may be NULL if not needed */
2933  )
2934 {
2935  SCIPdebugMessage("calling SCIPlpiGetSol()\n");
2936 
2937  assert(lpi != NULL);
2938  assert(lpi->spx != NULL);
2939 
2940  if( objval != NULL )
2941  *objval = lpi->spx->objValueReal();
2942 
2943  try
2944  {
2945  if( primsol != NULL )
2946  {
2947  Vector tmp(lpi->spx->numColsReal(), primsol);
2948  (void)lpi->spx->getPrimalReal(tmp);
2949  }
2950  if( dualsol != NULL )
2951  {
2952  Vector tmp(lpi->spx->numRowsReal(), dualsol);
2953  (void)lpi->spx->getDualReal(tmp);
2954  }
2955  if( activity != NULL )
2956  {
2957  Vector tmp(lpi->spx->numRowsReal(), activity);
2958  (void)lpi->spx->getSlacksReal(tmp); /* in SoPlex, the activities are called "slacks" */
2959  }
2960  if( redcost != NULL )
2961  {
2962  Vector tmp(lpi->spx->numColsReal(), redcost);
2963  (void)lpi->spx->getRedCostReal(tmp);
2964  }
2965  }
2966  catch(SPxException x)
2967  {
2968 #ifndef NDEBUG
2969  std::string s = x.what();
2970  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
2971 #endif
2972  return SCIP_LPERROR;
2973  }
2974 
2975  return SCIP_OKAY;
2976 }
2977 
2978 /** gets primal ray for unbounded LPs */
2980  SCIP_LPI* lpi, /**< LP interface structure */
2981  SCIP_Real* ray /**< primal ray */
2982  )
2983 { /*lint --e{715}*/
2984  SCIPdebugMessage("calling SCIPlpiGetPrimalRay()\n");
2985 
2986  assert(lpi != NULL);
2987  assert(lpi->spx != NULL);
2988 
2989  try
2990  {
2991  Vector tmp(lpi->spx->numColsReal(), ray);
2992  (void)lpi->spx->getPrimalRayReal(tmp);
2993  }
2994  catch(SPxException x)
2995  {
2996 #ifndef NDEBUG
2997  std::string s = x.what();
2998  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
2999 #endif
3000  return SCIP_LPERROR;
3001  }
3002 
3003  return SCIP_OKAY;
3004 }
3005 
3006 /** gets dual farkas proof for infeasibility */
3008  SCIP_LPI* lpi, /**< LP interface structure */
3009  SCIP_Real* dualfarkas /**< dual farkas row multipliers */
3010  )
3011 {
3012  SCIPdebugMessage("calling SCIPlpiGetDualfarkas()\n");
3013 
3014  assert(lpi != NULL);
3015  assert(lpi->spx != NULL);
3016 
3017  try
3018  {
3019  Vector tmp(lpi->spx->numRowsReal(), dualfarkas);
3020  (void)lpi->spx->getDualFarkasReal(tmp);
3021  }
3022  catch(SPxException x)
3023  {
3024 #ifndef NDEBUG
3025  std::string s = x.what();
3026  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
3027 #endif
3028  return SCIP_LPERROR;
3029  }
3030 
3031  return SCIP_OKAY;
3032 }
3033 
3034 /** gets the number of LP iterations of the last solve call */
3036  SCIP_LPI* lpi, /**< LP interface structure */
3037  int* iterations /**< pointer to store the number of iterations of the last solve call */
3038  )
3039 {
3040  SCIPdebugMessage("calling SCIPlpiGetIterations()\n");
3041 
3042  assert(lpi != NULL);
3043  assert(lpi->spx != NULL);
3044 
3045  *iterations = lpi->spx->numIterations();
3046 
3047  return SCIP_OKAY;
3048 }
3049 
3050 /** gets information about the quality of an LP solution
3051  *
3052  * Such information is usually only available, if also a (maybe not optimal) solution is available.
3053  * The LPI should return SCIP_INVALID for @p quality, if the requested quantity is not available.
3054  */
3055 extern
3057  SCIP_LPI* lpi, /**< LP interface structure */
3058  SCIP_LPSOLQUALITY qualityindicator, /**< indicates which quality should be returned */
3059  SCIP_Real* quality /**< pointer to store quality number */
3060  )
3061 {
3062  SCIPdebugMessage("calling SCIPlpiGetRealSolQuality()\n");
3063 
3064  assert(lpi != NULL);
3065  assert(quality != NULL);
3066 
3067  bool success;
3068 
3069  assert(lpi != NULL);
3070  assert(quality != NULL);
3071 
3072  SCIPdebugMessage("requesting solution quality from SoPlex: quality %d\n", qualityindicator);
3073 
3074  switch( qualityindicator )
3075  {
3077  success = lpi->spx->getEstimatedCondition(*quality);
3078  break;
3079 
3081  success = lpi->spx->getExactCondition(*quality);
3082  break;
3083 
3084  default:
3085  SCIPerrorMessage("Solution quality %d unknown.\n", qualityindicator);
3086  return SCIP_INVALIDDATA;
3087  }
3088 
3089  if( !success )
3090  {
3091  SCIPdebugMessage("problem computing condition number\n");
3092  *quality = SCIP_INVALID;
3093  }
3094 
3095  return SCIP_OKAY;
3096 }
3097 
3098 /**@} */
3099 
3100 
3101 
3102 
3103 /*
3104  * LP Basis Methods
3105  */
3106 
3107 /**@name LP Basis Methods */
3108 /**@{ */
3109 
3110 
3111 /** gets current basis status for columns and rows; arrays must be large enough to store the basis status */
3113  SCIP_LPI* lpi, /**< LP interface structure */
3114  int* cstat, /**< array to store column basis status, or NULL */
3115  int* rstat /**< array to store row basis status, or NULL */
3116  )
3117 {
3118  int i;
3119 
3120  SCIPdebugMessage("calling SCIPlpiGetBase()\n");
3121 
3122  assert(lpi != NULL);
3123  assert(lpi->spx != NULL);
3124 
3125  assert( lpi->spx->preStrongbranchingBasisFreed() );
3126 
3127  if( rstat != NULL && cstat != NULL )
3128  {
3129  for( i = 0; i < lpi->spx->numRowsReal(); ++i )
3130  {
3131  switch( lpi->spx->basisRowStatus(i) )
3132  {
3133  case SPxSolver::BASIC:
3134  rstat[i] = SCIP_BASESTAT_BASIC; /*lint !e641*/
3135  break;
3136  case SPxSolver::FIXED:
3137  case SPxSolver::ON_LOWER:
3138  rstat[i] = SCIP_BASESTAT_LOWER; /*lint !e641*/
3139  break;
3140  case SPxSolver::ON_UPPER:
3141  rstat[i] = SCIP_BASESTAT_UPPER; /*lint !e641*/
3142  break;
3143  case SPxSolver::ZERO:
3144  SCIPerrorMessage("slack variable has basis status ZERO (should not occur)\n");
3145  return SCIP_LPERROR;
3146  default:
3147  SCIPerrorMessage("invalid basis status\n");
3148  SCIPABORT();
3149  return SCIP_INVALIDDATA; /*lint !e527*/
3150  }
3151  }
3152  }
3153 
3154  if( cstat != NULL )
3155  {
3156  for( i = 0; i < lpi->spx->numColsReal(); ++i )
3157  {
3158  SCIP_Real val = 0.0;
3159  switch( lpi->spx->basisColStatus(i) )
3160  {
3161  case SPxSolver::BASIC:
3162  cstat[i] = SCIP_BASESTAT_BASIC; /*lint !e641*/
3163  break;
3164  case SPxSolver::FIXED:
3165  /* Get reduced cost estimation. If the estimation is not correct this should not hurt:
3166  * If the basis is loaded into SoPlex again, the status is converted to FIXED again; in
3167  * this case there is no problem at all. If the basis is saved and/or used in some other
3168  * solver, it usually is very cheap to perform the pivots necessary to get an optimal
3169  * basis. */
3170 // SCIP_CALL( getRedCostEst(lpi->spx, i, &val) );
3171  if( val < 0.0 ) /* reduced costs < 0 => UPPER else => LOWER */
3172  cstat[i] = SCIP_BASESTAT_UPPER; /*lint !e641*/
3173  else
3174  cstat[i] = SCIP_BASESTAT_LOWER; /*lint !e641*/
3175  break;
3176  case SPxSolver::ON_LOWER:
3177  cstat[i] = SCIP_BASESTAT_LOWER; /*lint !e641*/
3178  break;
3179  case SPxSolver::ON_UPPER:
3180  cstat[i] = SCIP_BASESTAT_UPPER; /*lint !e641*/
3181  break;
3182  case SPxSolver::ZERO:
3183  cstat[i] = SCIP_BASESTAT_ZERO; /*lint !e641*/
3184  break;
3185  default:
3186  SCIPerrorMessage("invalid basis status\n");
3187  SCIPABORT();
3188  return SCIP_INVALIDDATA; /*lint !e527*/
3189  }
3190  }
3191  }
3192 
3193  return SCIP_OKAY;
3194 }
3195 
3196 /** sets current basis status for columns and rows */
3198  SCIP_LPI* lpi, /**< LP interface structure */
3199  int* cstat, /**< array with column basis status */
3200  int* rstat /**< array with row basis status */
3201  )
3202 {
3203  int i;
3204 
3205  SCIPdebugMessage("calling SCIPlpiSetBase()\n");
3206 
3207  assert(lpi != NULL);
3208  assert(lpi->spx != NULL);
3209  assert(cstat != NULL || lpi->spx->numColsReal() == 0);
3210  assert(rstat != NULL || lpi->spx->numRowsReal() == 0);
3211 
3212  assert( lpi->spx->preStrongbranchingBasisFreed() );
3213  invalidateSolution(lpi);
3214 
3215  SPxSolver::VarStatus* spxcstat = NULL;
3216  SPxSolver::VarStatus* spxrstat = NULL;
3217  SCIP_ALLOC( BMSallocMemoryArray(&spxcstat, lpi->spx->numColsReal()) );
3218  SCIP_ALLOC( BMSallocMemoryArray(&spxrstat, lpi->spx->numRowsReal()) );
3219 
3220  for( i = 0; i < lpi->spx->numRowsReal(); ++i )
3221  {
3222  switch( rstat[i] )
3223  {
3224  case SCIP_BASESTAT_LOWER:
3225  spxrstat[i] = SPxSolver::ON_LOWER;
3226  break;
3227  case SCIP_BASESTAT_BASIC:
3228  spxrstat[i] = SPxSolver::BASIC;
3229  break;
3230  case SCIP_BASESTAT_UPPER:
3231  spxrstat[i] = SPxSolver::ON_UPPER;
3232  break;
3233  case SCIP_BASESTAT_ZERO:
3234  SCIPerrorMessage("slack variable has basis status ZERO (should not occur)\n");
3235  BMSfreeMemoryArrayNull(&spxcstat);
3236  BMSfreeMemoryArrayNull(&spxrstat);
3237  return SCIP_LPERROR; /*lint !e429*/
3238  default:
3239  SCIPerrorMessage("invalid basis status\n");
3240  SCIPABORT();
3241  return SCIP_INVALIDDATA; /*lint !e527*/
3242  }
3243  }
3244 
3245  for( i = 0; i < lpi->spx->numColsReal(); ++i )
3246  {
3247  switch( cstat[i] )
3248  {
3249  case SCIP_BASESTAT_LOWER:
3250  spxcstat[i] = SPxSolver::ON_LOWER;
3251  break;
3252  case SCIP_BASESTAT_BASIC:
3253  spxcstat[i] = SPxSolver::BASIC;
3254  break;
3255  case SCIP_BASESTAT_UPPER:
3256  spxcstat[i] = SPxSolver::ON_UPPER;
3257  break;
3258  case SCIP_BASESTAT_ZERO:
3259  spxcstat[i] = SPxSolver::ZERO;
3260  break;
3261  default:
3262  SCIPerrorMessage("invalid basis status\n");
3263  SCIPABORT();
3264  return SCIP_INVALIDDATA; /*lint !e527*/
3265  }
3266  }
3267 
3268  SOPLEX_TRY( lpi->messagehdlr, lpi->spx->setBasis(spxrstat, spxcstat) );
3269  // do we still need this?
3270 // lpi->spx->updateStatus();
3271 
3272  BMSfreeMemoryArrayNull(&spxcstat);
3273  BMSfreeMemoryArrayNull(&spxrstat);
3274 
3275  return SCIP_OKAY;
3276 }
3277 
3278 /** returns the indices of the basic columns and rows; basic column n gives value n, basic row m gives value -1-m */
3280  SCIP_LPI* lpi, /**< LP interface structure */
3281  int* bind /**< pointer to store basis indices ready to keep number of rows entries */
3282  )
3283 {
3284  SCIPdebugMessage("calling SCIPlpiGetBasisInd()\n");
3285 
3286  assert(lpi != NULL);
3287  assert(lpi->spx != NULL);
3288 
3289  assert(lpi->spx->preStrongbranchingBasisFreed());
3290 
3291  lpi->spx->getBasisInd(bind);
3292 
3293  return SCIP_OKAY;
3294 }
3295 
3296 
3297 /** get dense row of inverse basis matrix B^-1 */
3299  SCIP_LPI* lpi, /**< LP interface structure */
3300  int r, /**< row number */
3301  SCIP_Real* coef /**< pointer to store the coefficients of the row */
3302  )
3303 {
3304  SCIPdebugMessage("calling SCIPlpiGetBInvRow()\n");
3305 
3306  assert( lpi != NULL);
3307  assert(lpi->spx != NULL);
3308  assert(lpi->spx->preStrongbranchingBasisFreed());
3309 
3310  assert(r >= 0);
3311  assert(r < lpi->spx->numRowsReal());
3312 
3313  if( ! lpi->spx->getBasisInverseRowReal(r, coef) )
3314  return SCIP_LPERROR;
3315 
3316  return SCIP_OKAY;
3317 }
3318 
3319 /** get dense column of inverse basis matrix B^-1 */
3321  SCIP_LPI* lpi, /**< LP interface structure */
3322  int c, /**< column number of B^-1; this is NOT the number of the column in the LP;
3323  * you have to call SCIPlpiGetBasisInd() to get the array which links the
3324  * B^-1 column numbers to the row and column numbers of the LP!
3325  * c must be between 0 and nrows-1, since the basis has the size
3326  * nrows * nrows */
3327  SCIP_Real* coef /**< pointer to store the coefficients of the column */
3328  )
3329 {
3330  SCIPdebugMessage("calling SCIPlpiGetBInvCol()\n");
3331 
3332  assert( lpi != NULL );
3333  assert( lpi->spx != NULL );
3334  assert( lpi->spx->preStrongbranchingBasisFreed() );
3335 
3336  if( ! lpi->spx->getBasisInverseColReal(c, coef) )
3337  return SCIP_LPERROR;
3338 
3339  return SCIP_OKAY;
3340 }
3341 
3342 /** get dense row of inverse basis matrix times constraint matrix B^-1 * A */
3344  SCIP_LPI* lpi, /**< LP interface structure */
3345  int r, /**< row number */
3346  const SCIP_Real* binvrow, /**< row in (A_B)^-1 from prior call to SCIPlpiGetBInvRow(), or NULL */
3347  SCIP_Real* coef /**< vector to return coefficients */
3348  )
3349 {
3350  SCIP_Real* buf;
3351  SCIP_Real* binv;
3352  int nrows;
3353  int ncols;
3354  int c;
3355 
3356  SCIPdebugMessage("calling SCIPlpiGetBInvARow()\n");
3357 
3358  assert(lpi != NULL);
3359  assert(lpi->spx != NULL);
3360  assert( lpi->spx->preStrongbranchingBasisFreed() );
3361 
3362  nrows = lpi->spx->numRowsReal();
3363  ncols = lpi->spx->numColsReal();
3364  buf = NULL;
3365 
3366  /* get (or calculate) the row in B^-1 */
3367  if( binvrow == NULL )
3368  {
3369  SCIP_ALLOC( BMSallocMemoryArray(&buf, nrows) );
3370  SCIP_CALL( SCIPlpiGetBInvRow(lpi, r, buf) );
3371  binv = buf;
3372  }
3373  else
3374  binv = const_cast<SCIP_Real*>(binvrow);
3375 
3376  assert(binv != NULL);
3377 
3378  // @todo exploit sparsity in binv by looping over nrows
3379  /* calculate the scalar product of the row in B^-1 and A */
3380  Vector binvvec(nrows, binv);
3381  for( c = 0; c < ncols; ++c )
3382  coef[c] = binvvec * lpi->spx->colVectorReal(c); /* scalar product */ /*lint !e1702*/
3383 
3384  /* free memory if it was temporarily allocated */
3385  BMSfreeMemoryArrayNull(&buf);
3386 
3387  return SCIP_OKAY;
3388 }
3389 
3390 /** get dense column of inverse basis matrix times constraint matrix B^-1 * A */
3392  SCIP_LPI* lpi, /**< LP interface structure */
3393  int c, /**< column number */
3394  SCIP_Real* coef /**< vector to return coefficients */
3395  )
3396 {
3397  DVector col(lpi->spx->numRowsReal());
3398 
3399  SCIPdebugMessage("calling SCIPlpiGetBInvACol()\n");
3400 
3401  assert( lpi != NULL );
3402  assert( lpi->spx != NULL );
3403  assert( lpi->spx->preStrongbranchingBasisFreed() );
3404 
3405  /* extract column c of A */
3406  assert(c >= 0);
3407  assert(c < lpi->spx->numColsReal());
3408 
3409  // @todo why is clear() needed?
3410  col.clear();
3411  col = lpi->spx->colVectorReal(c);
3412  // @todo col should be of this size already
3413  col.reDim(lpi->spx->numRowsReal());
3414 
3415  /* solve */
3416  if( ! lpi->spx->getBasisInverseTimesVecReal(col.get_ptr(), coef) )
3417  return SCIP_LPERROR;
3418 
3419  return SCIP_OKAY;
3420 }
3421 
3422 /**@} */
3423 
3424 
3425 
3426 
3427 /*
3428  * LP State Methods
3429  */
3430 
3431 /**@name LP State Methods */
3432 /**@{ */
3433 
3434 /** stores LPi state (like basis information) into lpistate object */
3436  SCIP_LPI* lpi, /**< LP interface structure */
3437  BMS_BLKMEM* blkmem, /**< block memory */
3438  SCIP_LPISTATE** lpistate /**< pointer to LPi state information (like basis information) */
3439  )
3440 {
3441  int ncols;
3442  int nrows;
3443 
3444  SCIPdebugMessage("calling SCIPlpiGetState()\n");
3445 
3446  assert(blkmem != NULL);
3447  assert(lpi != NULL);
3448  assert(lpi->spx != NULL);
3449  assert(lpistate != NULL);
3450 
3451  assert( lpi->spx->preStrongbranchingBasisFreed() );
3452 
3453  ncols = lpi->spx->numColsReal();
3454  nrows = lpi->spx->numRowsReal();
3455  assert(ncols >= 0);
3456  assert(nrows >= 0);
3457 
3458  /* allocate lpistate data */
3459  SCIP_CALL( lpistateCreate(lpistate, blkmem, ncols, nrows) );
3460 
3461  /* allocate enough memory for storing uncompressed basis information */
3462  SCIP_CALL( ensureCstatMem(lpi, ncols) );
3463  SCIP_CALL( ensureRstatMem(lpi, nrows) );
3464 
3465  /* get unpacked basis information */
3466  SCIP_CALL( SCIPlpiGetBase(lpi, lpi->cstat, lpi->rstat) );
3467 
3468  /* pack LPi state data */
3469  (*lpistate)->ncols = ncols;
3470  (*lpistate)->nrows = nrows;
3471  lpistatePack(*lpistate, lpi->cstat, lpi->rstat);
3472 
3473  return SCIP_OKAY;
3474 }
3475 
3476 /** loads LPi state (like basis information) into solver; note that the LP might have been extended with additional
3477  * columns and rows since the state was stored with SCIPlpiGetState()
3478  */
3480  SCIP_LPI* lpi, /**< LP interface structure */
3481  BMS_BLKMEM* /*blkmem*/, /**< block memory */
3482  SCIP_LPISTATE* lpistate /**< LPi state information (like basis information) */
3483  )
3484 {
3485  int lpncols;
3486  int lpnrows;
3487  int i;
3488 
3489  SCIPdebugMessage("calling SCIPlpiSetState()\n");
3490 
3491  assert(lpi != NULL);
3492  assert(lpi->spx != NULL);
3493  assert(lpistate != NULL);
3494 
3495  assert( lpi->spx->preStrongbranchingBasisFreed() );
3496 
3497  lpncols = lpi->spx->numColsReal();
3498  lpnrows = lpi->spx->numRowsReal();
3499  assert(lpistate->ncols <= lpncols);
3500  assert(lpistate->nrows <= lpnrows);
3501 
3502  /* allocate enough memory for storing uncompressed basis information */
3503  SCIP_CALL( ensureCstatMem(lpi, lpncols) );
3504  SCIP_CALL( ensureRstatMem(lpi, lpnrows) );
3505 
3506  /* unpack LPi state data */
3507  lpistateUnpack(lpistate, lpi->cstat, lpi->rstat);
3508 
3509  /* extend the basis to the current LP beyond the previously existing columns */
3510  for( i = lpistate->ncols; i < lpncols; ++i )
3511  {
3512  SCIP_Real bnd = lpi->spx->lowerReal(i);
3513  if ( SCIPlpiIsInfinity(lpi, REALABS(bnd)) )
3514  {
3515  /* if lower bound is +/- infinity -> try upper bound */
3516  bnd = lpi->spx->lowerReal(i);
3517  if ( SCIPlpiIsInfinity(lpi, REALABS(bnd)) )
3518  lpi->cstat[i] = SCIP_BASESTAT_ZERO; /* variable is free */
3519  else
3520  lpi->cstat[i] = SCIP_BASESTAT_UPPER; /* use finite upper bound */
3521  }
3522  else
3523  lpi->cstat[i] = SCIP_BASESTAT_LOWER; /* use finite lower bound */
3524  }
3525  for( i = lpistate->nrows; i < lpnrows; ++i )
3526  lpi->rstat[i] = SCIP_BASESTAT_BASIC; /*lint !e641*/
3527 
3528  /* load basis information */
3529  SCIP_CALL( SCIPlpiSetBase(lpi, lpi->cstat, lpi->rstat) );
3530 
3531  return SCIP_OKAY;
3532 }
3533 
3534 /** clears current LPi state (like basis information) of the solver */
3536  SCIP_LPI* lpi /**< LP interface structure */
3537  )
3538 { /*lint --e{715}*/
3539  SCIPdebugMessage("calling SCIPlpiClearState()\n");
3540 
3541  assert(lpi != NULL);
3542  assert(lpi->spx != NULL);
3543 
3544  try
3545  {
3546  lpi->spx->clearBasis();
3547  }
3548  catch(SPxException x)
3549  {
3550 #ifndef NDEBUG
3551  std::string s = x.what();
3552  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
3553 #endif
3554  assert( lpi->spx->status() != SPxSolver::OPTIMAL );
3555  return SCIP_LPERROR;
3556  }
3557 
3558  return SCIP_OKAY;
3559 }
3560 
3561 /** frees LPi state information */
3563  SCIP_LPI* lpi, /**< LP interface structure */
3564  BMS_BLKMEM* blkmem, /**< block memory */
3565  SCIP_LPISTATE** lpistate /**< pointer to LPi state information (like basis information) */
3566  )
3567 { /*lint --e{715}*/
3568  SCIPdebugMessage("calling SCIPlpiFreeState()\n");
3569 
3570  assert(lpi != NULL);
3571  assert(lpistate != NULL);
3572 
3573  if ( *lpistate != NULL )
3574  lpistateFree(lpistate, blkmem);
3575 
3576  return SCIP_OKAY;
3577 }
3578 
3579 /** checks, whether the given LP state contains simplex basis information */
3581  SCIP_LPI* lpi, /**< LP interface structure */
3582  SCIP_LPISTATE* lpistate /**< LP state information (like basis information) */
3583  )
3584 { /*lint --e{715}*/
3585  return TRUE;
3586 }
3587 
3588 /** reads LP state (like basis information from a file */
3590  SCIP_LPI* lpi, /**< LP interface structure */
3591  const char* fname /**< file name */
3592  )
3593 {
3594  SCIPdebugMessage("calling SCIPlpiReadState()\n");
3595 
3596  assert( lpi->spx->preStrongbranchingBasisFreed() );
3597 
3598  bool success;
3599  SOPLEX_TRY( lpi->messagehdlr, success = lpi->spx->readBasisFile(fname, 0, 0) );
3600 
3601  return success ? SCIP_OKAY : SCIP_LPERROR;
3602 }
3603 
3604 /** writes LP state (like basis information) to a file */
3606  SCIP_LPI* lpi, /**< LP interface structure */
3607  const char* fname /**< file name */
3608  )
3609 {
3610  SCIPdebugMessage("calling SCIPlpiWriteState()\n");
3611 
3612  assert( lpi->spx->preStrongbranchingBasisFreed() );
3613 
3614  bool res;
3615  SOPLEX_TRY( lpi->messagehdlr, res = lpi->spx->writeBasisFile(fname, 0, 0) );
3616 
3617  if ( ! res )
3618  return SCIP_LPERROR;
3619 
3620  return SCIP_OKAY;
3621 }
3622 
3623 /**@} */
3624 
3625 
3626 
3627 
3628 /*
3629  * LP Pricing Norms Methods
3630  */
3631 
3632 /**@name LP Pricing Norms Methods */
3633 /**@{ */
3634 
3635 /** stores LPi pricing norms information
3636  * @todo should we store norm information?
3637  */
3639  SCIP_LPI* lpi, /**< LP interface structure */
3640  BMS_BLKMEM* blkmem, /**< block memory */
3641  SCIP_LPINORMS** lpinorms /**< pointer to LPi pricing norms information */
3642  )
3643 {
3644  assert(lpinorms != NULL);
3645 
3646  (*lpinorms) = NULL;
3647 
3648  return SCIP_OKAY;
3649 }
3650 
3651 /** loads LPi pricing norms into solver; note that the LP might have been extended with additional
3652  * columns and rows since the state was stored with SCIPlpiGetNorms()
3653  */
3655  SCIP_LPI* lpi, /**< LP interface structure */
3656  BMS_BLKMEM* blkmem, /**< block memory */
3657  SCIP_LPINORMS* lpinorms /**< LPi pricing norms information */
3658  )
3659 {
3660  assert(lpinorms == NULL);
3661 
3662  /* no work necessary */
3663  return SCIP_OKAY;
3664 }
3665 
3666 /** frees pricing norms information */
3668  SCIP_LPI* lpi, /**< LP interface structure */
3669  BMS_BLKMEM* blkmem, /**< block memory */
3670  SCIP_LPINORMS** lpinorms /**< pointer to LPi pricing norms information */
3671  )
3672 {
3673  assert(lpinorms == NULL);
3674 
3675  /* no work necessary */
3676  return SCIP_OKAY;
3677 }
3678 
3679 /**@} */
3680 
3681 
3682 
3683 
3684 /*
3685  * Parameter Methods
3686  */
3687 
3688 /**@name Parameter Methods */
3689 /**@{ */
3690 
3691 /** gets integer parameter of LP */
3693  SCIP_LPI* lpi, /**< LP interface structure */
3694  SCIP_LPPARAM type, /**< parameter number */
3695  int* ival /**< buffer to store the parameter value */
3696  )
3697 {
3698  SCIPdebugMessage("calling SCIPlpiGetIntpar()\n");
3699 
3700  assert(lpi != NULL);
3701  assert(lpi->spx != NULL);
3702  assert(ival != NULL);
3703 
3704  switch( type )
3705  {
3707  *ival = lpi->spx->getFromScratch();
3708  break;
3709  case SCIP_LPPAR_LPINFO:
3710  *ival = lpi->spx->getLpInfo();
3711  break;
3712  case SCIP_LPPAR_LPITLIM:
3713  *ival = lpi->spx->intParam(SoPlex::ITERLIMIT);
3714  break;
3715  case SCIP_LPPAR_PRESOLVING:
3716  *ival = lpi->spx->intParam(SoPlex::SIMPLIFIER) == SoPlex::SIMPLIFIER_AUTO;
3717  break;
3718  case SCIP_LPPAR_PRICING:
3719  *ival = (int) lpi->pricing;
3720  break;
3721  case SCIP_LPPAR_SCALING:
3722  *ival = (int) (lpi->spx->intParam(SoPlex::SCALER) != SoPlex::SCALER_OFF);
3723  break;
3724  default:
3725  return SCIP_PARAMETERUNKNOWN;
3726  } /*lint !e788*/
3727 
3728  return SCIP_OKAY;
3729 }
3730 
3731 /** sets integer parameter of LP */
3733  SCIP_LPI* lpi, /**< LP interface structure */
3734  SCIP_LPPARAM type, /**< parameter number */
3735  int ival /**< parameter value */
3736  )
3737 {
3738  SCIPdebugMessage("calling SCIPlpiSetIntpar()\n");
3739 
3740  assert(lpi != NULL);
3741  assert(lpi->spx != NULL);
3742 
3743  switch( type )
3744  {
3746  assert(ival == TRUE || ival == FALSE);
3747  lpi->spx->setFromScratch(bool(ival));
3748  break;
3749  case SCIP_LPPAR_LPINFO:
3750  assert(ival == TRUE || ival == FALSE);
3751  lpi->spx->setLpInfo(bool(ival));
3752  break;
3753  case SCIP_LPPAR_LPITLIM:
3754  assert(ival >= -1);
3755  lpi->spx->setIntParam(SoPlex::ITERLIMIT, ival);
3756  break;
3757  case SCIP_LPPAR_PRESOLVING:
3758  assert(ival == TRUE || ival == FALSE);
3759  lpi->spx->setIntParam(SoPlex::SIMPLIFIER, (ival ? SoPlex::SIMPLIFIER_AUTO : SoPlex::SIMPLIFIER_OFF));
3760  break;
3761  case SCIP_LPPAR_PRICING:
3762  lpi->pricing = (SCIP_PRICING)ival;
3763  switch( lpi->pricing )
3764  {
3766  case SCIP_PRICING_AUTO:
3767  lpi->spx->setIntParam(SoPlex::PRICER, SoPlex::PRICER_AUTO);
3768  break;
3769  case SCIP_PRICING_FULL:
3770  lpi->spx->setIntParam(SoPlex::PRICER, SoPlex::PRICER_STEEP);
3771  break;
3772  case SCIP_PRICING_PARTIAL:
3773  lpi->spx->setIntParam(SoPlex::PRICER, SoPlex::PRICER_PARMULT);
3774  break;
3775  case SCIP_PRICING_STEEP:
3776  lpi->spx->setIntParam(SoPlex::PRICER, SoPlex::PRICER_STEEP);
3777  break;
3779  lpi->spx->setIntParam(SoPlex::PRICER, SoPlex::PRICER_QUICKSTEEP);
3780  break;
3781  case SCIP_PRICING_DEVEX:
3782  lpi->spx->setIntParam(SoPlex::PRICER, SoPlex::PRICER_DEVEX);
3783  break;
3784  default:
3785  return SCIP_LPERROR;
3786  }
3787  break;
3788  case SCIP_LPPAR_SCALING:
3789  assert(ival == TRUE || ival == FALSE);
3790  lpi->spx->setIntParam(SoPlex::SCALER, ( ival ? SoPlex::SCALER_BIEQUI : SoPlex::SCALER_OFF));
3791  break;
3792  default:
3793  return SCIP_PARAMETERUNKNOWN;
3794  } /*lint !e788*/
3795 
3796  return SCIP_OKAY;
3797 }
3798 
3799 /** gets floating point parameter of LP */
3801  SCIP_LPI* lpi, /**< LP interface structure */
3802  SCIP_LPPARAM type, /**< parameter number */
3803  SCIP_Real* dval /**< buffer to store the parameter value */
3804  )
3805 {
3806  SCIPdebugMessage("calling SCIPlpiGetRealpar()\n");
3807 
3808  assert(lpi != NULL);
3809  assert(lpi->spx != NULL);
3810  assert(dval != NULL);
3811 
3812  switch( type )
3813  {
3814  case SCIP_LPPAR_FEASTOL:
3815  *dval = lpi->spx->feastol();
3816  break;
3818  *dval = lpi->spx->opttol();
3819  break;
3820  case SCIP_LPPAR_LOBJLIM:
3821  *dval = lpi->spx->realParam(SoPlex::OBJLIMIT_LOWER);
3822  break;
3823  case SCIP_LPPAR_UOBJLIM:
3824  *dval = lpi->spx->realParam(SoPlex::OBJLIMIT_UPPER);
3825  break;
3826  case SCIP_LPPAR_LPTILIM:
3827  *dval = lpi->spx->realParam(SoPlex::TIMELIMIT);
3828  break;
3830  *dval = lpi->rowrepswitch;
3831  break;
3833  *dval = lpi->conditionlimit;
3834  default:
3835  return SCIP_PARAMETERUNKNOWN;
3836  } /*lint !e788*/
3837 
3838  return SCIP_OKAY;
3839 }
3840 
3841 /** sets floating point parameter of LP */
3843  SCIP_LPI* lpi, /**< LP interface structure */
3844  SCIP_LPPARAM type, /**< parameter number */
3845  SCIP_Real dval /**< parameter value */
3846  )
3847 {
3848  SCIPdebugMessage("calling SCIPlpiSetRealpar()\n");
3849 
3850  assert(lpi != NULL);
3851  assert(lpi->spx != NULL);
3852 
3853  switch( type )
3854  {
3855  case SCIP_LPPAR_FEASTOL:
3856  lpi->spx->setFeastol(dval);
3857  break;
3859  lpi->spx->setOpttol(dval);
3860  break;
3861  case SCIP_LPPAR_LOBJLIM:
3862  lpi->spx->setRealParam(SoPlex::OBJLIMIT_LOWER, dval);
3863  break;
3864  case SCIP_LPPAR_UOBJLIM:
3865  lpi->spx->setRealParam(SoPlex::OBJLIMIT_UPPER, dval);
3866  break;
3867  case SCIP_LPPAR_LPTILIM:
3868  lpi->spx->setRealParam(SoPlex::TIMELIMIT, dval);
3869  break;
3871  assert(dval >= -1.5);
3872  lpi->rowrepswitch = dval;
3873  break;
3875  lpi->conditionlimit = dval;
3876  lpi->checkcondition = (dval >= 0);
3877  break;
3878  default:
3879  return SCIP_PARAMETERUNKNOWN;
3880  } /*lint !e788*/
3881 
3882  return SCIP_OKAY;
3883 }
3884 
3885 /**@} */
3886 
3887 
3888 
3889 
3890 /*
3891  * Numerical Methods
3892  */
3893 
3894 /**@name Numerical Methods */
3895 /**@{ */
3896 
3897 /** returns value treated as infinity in the LP solver */
3899  SCIP_LPI* lpi /**< LP interface structure */
3900  )
3901 {
3902  SCIPdebugMessage("calling SCIPlpiInfinity()\n");
3903 
3904  return lpi->spx->realParam(SoPlex::INFTY);
3905 }
3906 
3907 /** checks if given value is treated as infinity in the LP solver */
3909  SCIP_LPI* lpi, /**< LP interface structure */
3910  SCIP_Real val
3911  )
3912 {
3913  SCIPdebugMessage("calling SCIPlpiIsInfinity()\n");
3914 
3915  return (val >= lpi->spx->realParam(SoPlex::INFTY));
3916 }
3917 
3918 /**@} */
3919 
3920 
3921 
3922 
3923 /*
3924  * File Interface Methods
3925  */
3926 
3927 /**@name File Interface Methods */
3928 /**@{ */
3929 
3930 /** returns, whether the given file exists */
3931 static
3932 SCIP_Bool fileExists(
3933  const char* filename /**< file name */
3934  )
3935 {
3936  FILE* f;
3937 
3938  f = fopen(filename, "r");
3939  if( f == NULL )
3940  return FALSE;
3941 
3942  fclose(f);
3943 
3944  return TRUE;
3945 }
3946 
3947 /** reads LP from a file */
3949  SCIP_LPI* lpi, /**< LP interface structure */
3950  const char* fname /**< file name */
3951  )
3952 {
3953  SCIPdebugMessage("calling SCIPlpiReadLP()\n");
3954 
3955  assert(lpi != NULL);
3956  assert(lpi->spx != NULL);
3957 
3958  assert( lpi->spx->preStrongbranchingBasisFreed() );
3959 
3960  if( !fileExists(fname) )
3961  return SCIP_NOFILE;
3962 
3963  try
3964  {
3965  assert(lpi->spx->intParam(SoPlex::READMODE) == SoPlex::READMODE_REAL);
3966  if( !lpi->spx->readFile(fname) )
3967  return SCIP_READERROR;
3968  }
3969  catch(SPxException x)
3970  {
3971 #ifndef NDEBUG
3972  std::string s = x.what();
3973  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
3974 #endif
3975  return SCIP_READERROR;
3976  }
3977 
3978  return SCIP_OKAY;
3979 }
3980 
3981 /** writes LP to a file */
3983  SCIP_LPI* lpi, /**< LP interface structure */
3984  const char* fname /**< file name */
3985  )
3986 {
3987  SCIPdebugMessage("calling SCIPlpiWriteLP()\n");
3988 
3989  assert(lpi != NULL);
3990  assert(lpi->spx != NULL);
3991 
3992  try
3993  {
3994  lpi->spx->writeFileReal(fname);
3995  }
3996  catch(SPxException x)
3997  {
3998 #ifndef NDEBUG
3999  std::string s = x.what();
4000  SCIPmessagePrintWarning(lpi->messagehdlr, "SoPlex threw an exception: %s\n", s.c_str());
4001 #endif
4002  return SCIP_WRITEERROR;
4003  }
4004 
4005  return SCIP_OKAY;
4006 }
4007 
4008 /**@} */
4009