68 #ifdef SOPLEX_WITH_GMP 72 #ifdef SOPLEX_WITH_BOOST 73 #ifdef SOPLEX_WITH_MPFR 75 #include <boost/multiprecision/mpfr.hpp> 77 #ifdef SOPLEX_WITH_CPPMPF 78 #include <boost/multiprecision/cpp_dec_float.hpp> 85 #define DEFAULT_RANDOM_SEED 0 // used to suppress output when the seed was not changed 230 bool getDualNorms(
int& nnormsRow,
int& nnormsCol, R* norms)
const;
233 bool setDualNorms(
int nnormsRow,
int nnormsCol, R* norms);
428 #ifdef SOPLEX_WITH_GMP 430 void addRowRational(
const mpq_t* lhs,
const mpq_t* rowValues,
const int* rowIndices,
431 const int rowSize,
const mpq_t* rhs);
434 void addRowsRational(
const mpq_t* lhs,
const mpq_t* rowValues,
const int* rowIndices,
435 const int* rowStarts,
const int* rowLengths,
const int numRows,
const int numValues,
445 #ifdef SOPLEX_WITH_GMP 447 void addColRational(
const mpq_t* obj,
const mpq_t* lower,
const mpq_t* colValues,
448 const int* colIndices,
const int colSize,
const mpq_t* upper);
451 void addColsRational(
const mpq_t* obj,
const mpq_t* lower,
const mpq_t* colValues,
452 const int* colIndices,
const int* colStarts,
const int* colLengths,
const int numCols,
453 const int numValues,
const mpq_t* upper);
468 #ifdef SOPLEX_WITH_GMP 476 #ifdef SOPLEX_WITH_GMP 490 #ifdef SOPLEX_WITH_GMP 504 #ifdef SOPLEX_WITH_GMP 515 #ifdef SOPLEX_WITH_GMP 526 #ifdef SOPLEX_WITH_GMP 537 #ifdef SOPLEX_WITH_GMP 545 #ifdef SOPLEX_WITH_GMP 709 #ifdef SOPLEX_WITH_GMP 792 bool unscale =
true);
801 bool unscale =
true);
809 bool multBasis(R* vec,
bool unscale =
true);
892 const DIdxSet* intvars = 0,
const bool unscale =
true)
const;
899 const DIdxSet* intvars = 0,
const bool unscale =
true)
const;
914 const bool cpxFormat =
false)
const;
919 const bool cpxFormat =
false)
const;
924 const NameSet* colNames = 0,
const bool cpxFormat =
false)
const;
1421 #ifdef SOPLEX_WITH_RATIONALPARAM 1426 RATIONALPARAM_COUNT = 0
1478 #ifdef SOPLEX_WITH_RATIONALPARAM 1479 static struct RationalParam
1505 #ifdef SOPLEX_WITH_RATIONALPARAM 1531 #ifdef SOPLEX_WITH_RATIONALPARAM 1533 Rational rationalParam(
const RationalParam param)
const;
1548 #ifdef SOPLEX_WITH_RATIONALPARAM 1550 bool setRationalParam(
const RationalParam param,
const Rational value,
const bool init =
true);
1557 void resetSettings(
const bool quiet =
false,
const bool init =
true);
1563 bool saveSettingsFile(
const char* filename,
const bool onlyChanged =
false,
1564 int solvemode = 1)
const;
1609 bool areLPsInSync(
const bool checkVecVals =
true,
const bool checkMatVals =
false,
1610 const bool quiet =
false)
const;
1945 void _idxToPerm(
int* idx,
int idxSize,
int* perm,
int permSize)
const;
1948 void _rangeToPerm(
int start,
int end,
int* perm,
int permSize)
const;
2135 bool acceptUnbounded,
2136 bool acceptInfeasible,
2138 bool& primalFeasible,
2148 bool& stoppedIter,
bool& error);
2152 bool& stoppedIter,
bool& error);
2233 const bool forceNoSimplifier =
false);
2242 bool& stoppedIter,
bool& error,
bool& optimal);
2248 const Rational& denomBoundSquared);
2255 void _optimize(
volatile bool* interrupt = NULL);
2310 bool applyPreprocessing);
2318 int* nnonposind,
bool& stop);
2322 int* rowsforremoval,
int* colsforremoval,
2323 int nnonposind,
int* ncompatind,
bool formRedProb,
bool& stop);
2330 int* ncompatboundcons,
2331 int nnonposind,
bool& stop);
2335 int* nrowsforremoval,
2433 int& nNonBasicRows);
2448 #include "soplex/solvedbds.hpp" 2449 #include "soplex.hpp" 2450 #include "soplex/solverational.hpp" 2451 #include "soplex/testsoplex.hpp" 2452 #include "soplex/solvereal.hpp" 2454 #endif // _SOPLEX_H_ void getBasis(typename SPxSolverBase< R >::VarStatus rows[], typename SPxSolverBase< R >::VarStatus cols[]) const
gets current basis via arrays of statuses
SPxMainSM< R > _simplifierMainSM
SPxSteepPR< R > _pricerQuickSteep
SPxFastRT< R > _ratiotesterFast
decompStatus _currentProb
LPRowRational::Type rowTypeRational(int i) const
returns inequality type of row i
Fast shifting ratio test.
const VectorRational & upperRational() const
returns upper bound vector
Settings & operator=(const Settings &settings)
assignment operator
const char * getRatiotesterName()
name of currently loaded ratiotester
bool writeFile(const char *filename, const NameSet *rowNames=0, const NameSet *colNames=0, const DIdxSet *intvars=0, const bool unscale=true) const
Templated write function Real writes real LP to file; LP or MPS format is chosen from the extension i...
minimal Markowitz threshold to control sparsity/stability in LU factorization
void getNdualNorms(int &nnormsRow, int &nnormsCol) const
gets number of available dual norms
Rational _rationalNegInfty
void printSolutionStatistics(std::ostream &os)
prints solution statistics
void _writeOriginalProblemBasis(const char *filename, NameSet *rowNames, NameSet *colNames, bool cpxFormat)
function to build a basis for the original problem as given by the solution to the reduced problem ...
void _changeLhsReal(const VectorBase< R > &lhs)
changes left-hand side vector for constraints to lhs and adjusts basis
void _removeRowReal(int i)
removes row i and adjusts basis
void _ensureRealLPLoaded()
ensures that the real LP and the basis are loaded in the solver; performs no sync ...
steepest edge pricer with exact initialization of norms
SoPlex start basis generation base class.
virtual ~SoPlexBase()
destructor
bool ignoreUnscaledViolations()
sets the status to OPTIMAL in case the LP has been solved with unscaled violations ...
void clearLPReal()
clears the LP
void _changeRangeReal(const VectorBase< R > &lhs, const VectorBase< R > &rhs)
changes left- and right-hand side vectors and adjusts basis
working tolerance for optimality in floating-point solver during iterative refinement ...
const UnitVectorRational * _unitVectorRational(const int i)
returns pointer to a constant unit vector available until destruction of the SoPlexBase class ...
Bound flipping ratio test ("long step dual") for SoPlex.Class SPxBoundFlippingRT provides an implemen...
maximum number of updates without fresh factorization
std::string name[SoPlexBase< R >::BOOLPARAM_COUNT]
array of names for boolean parameters
Types of solution classes.
number< gmp_rational, et_off > Rational
DataArray< SPxRowId > _decompReducedProbColRowIDs
void _loadRealLP(bool initBasis)
load original LP and possibly setup a slack basis
void _getCompatibleBoundCons(LPRowSetBase< R > &boundcons, int *compatboundcons, int *nonposind, int *ncompatboundcons, int nnonposind, bool &stop)
computes the compatible bound constraints and adds them to the reduced problem
const char * getPricerName()
name of currently loaded pricer
VectorBase< R > _manualUpper
bool getExactCondition(R &condition)
computes the exact condition number for the current basis matrix using the power method; returns true...
void _updateDecompComplementaryDualProblem(bool origObj)
update the dual complementary problem with additional columns and rows
lower threshold in lifting (nonzero matrix coefficients with smaller absolute value will be reformula...
bool getBasisInverseRowRational(const int r, SSVectorRational &vec)
computes row r of basis inverse; performs rational factorization if not available; returns true on su...
int numNonzerosRational() const
bool _readFileRational(const char *filename, NameSet *rowNames=0, NameSet *colNames=0, DIdxSet *intVars=0)
reads rational LP in LP or MPS format from file and returns true on success; gets row names...
SPxLPBase< R > _manualRealLP
Implementation of Sparse Linear Solver with Rational precision.This class implements a SLinSolverRati...
lower and upper bound finite, but different
void _untransformFeasibility(SolRational &sol, bool infeasible)
undoes transformation to feasibility problem
const VectorRational & lowerRational() const
returns lower bound vector
SPxEquiliSC< R > _scalerUniequi
bool getDualFarkasReal(R *vector, int dim)
bool getDualNorms(int &nnormsRow, int &nnormsCol, R *norms) const
gets steepest edge norms and returns false if they are not available
the maximum number of rows that are added in each iteration of the decomposition based simplex ...
void _disableSimplifierAndScaler()
disables simplifier and scaler
void _removeColReal(int i)
removes column i
R maxAbsNonzeroReal() const
returns biggest non-zero element in absolute value
bool getBoundViolationRational(Rational &maxviol, Rational &sumviol)
R operator()(RowViolation i, RowViolation j) const
Devex pricer.The Devex Pricer for SoPlex implements an approximate steepest edge pricing, that does without solving an extra linear system and computing the scalar products.
void _computeReducedProbObjCoeff(bool &stop)
computes the reduced problem objective coefficients
void _optimize(volatile bool *interrupt=NULL)
solves the templated LP
R upperReal(int i) const
returns upper bound of column i
R coefReal(int row, int col) const
returns (unscaled) coefficient
int dlcmSizePrimalRational(const int base=2)
get size of least common multiple of denominators in primal solution
maximum number of conjugate gradient iterations in least square scaling
void addColsRational(const mpq_t *obj, const mpq_t *lower, const mpq_t *colValues, const int *colIndices, const int *colStarts, const int *colLengths, const int numCols, const int numValues, const mpq_t *upper)
adds a set of columns (GMP only method)
textbook ratio test without stabilization
Result
Result of the simplification.
void changeRhsRational(const VectorRational &rhs)
changes right-hand side vector to rhs
DataArray< SPxRowId > _decompCompPrimalRowIDs
Steepest edge pricer.Class SPxSteepExPR implements a steepest edge pricer to be used with SoPlex...
refactor threshold for fill-in in current factor update compared to fill-in in last factorization ...
SPxSolverBase< R >::Status status() const
returns the current solver status
Dense vector.Class VectorBase provides dense linear algebra vectors. Internally, VectorBase wraps std...
void _optimizeRational(volatile bool *interrupt=NULL)
temporary fix for Rational
Geometric mean row/column scaling.This SPxScaler implementation performs geometric mean scaling of th...
void _updateDecompComplementaryPrimalProblem(bool origObj)
update the primal complementary problem with additional columns and rows
automatic sync of real and rational LP
bool getDualRational(VectorRational &vector)
bool setBoolParam(const BoolParam param, const bool value, const bool init=true)
sets boolean parameter value; returns true on success
geometric mean scaling on rows and columns, max 1 round
crash basis from a greedy solution
DataArray< SPxColId > _decompFixedVarDualIDs
void _getRowsForRemovalComplementaryProblem(int *nonposind, int *bind, int *rowsforremoval, int *nrowsforremoval, int nnonposind)
computes the rows to remove from the complementary problem
const char * getSimplifierName()
name of simplifier
Abstract pricer base class.
VectorBase< R > _decompFeasVector
int numNonzeros() const
returns number of nonzeros
SPxParMultPR< R > _pricerParMult
threshold on number of rows vs. number of columns for switching from column to row representations in...
SPxSumST< R > _starterSum
void removeColsReal(int perm[])
removes all columns with an index i such that perm[i] < 0; upon completion, perm[i] >= 0 indicates th...
void addRowsRational(const mpq_t *lhs, const mpq_t *rowValues, const int *rowIndices, const int *rowStarts, const int *rowLengths, const int numRows, const int numValues, const mpq_t *rhs)
adds a set of rows (GMP only method)
void _completeRangeTypesRational()
completes range type arrays after adding columns and/or rows
Solution vector based start basis.
void getRowsRational(int start, int end, LPRowSetRational &lprowset) const
gets rows start, ..., end.
void addRowReal(const LPRowBase< R > &lprow)
adds a single row
void _findViolatedRows(R compObjValue, Array< RowViolation > &violatedrows, int &nviolatedrows)
builds the update rows with those violated in the complmentary problem
void removeRowRangeReal(int start, int end, int perm[]=0)
removes rows start to end including both; an array perm of size numRows() may be passed as buffer mem...
SPxDefaultRT< R > _ratiotesterTextbook
void changeUpperRational(const VectorRational &upper)
changes vector of upper bounds to upper
VectorRational _feasLower
Sequential object-oriented SimPlex.SPxSolverBase is an LP solver class using the revised Simplex algo...
primal simplex algorithm, i.e., entering for column and leaving for row representation ...
void _recomputeRangeTypesReal()
recomputes range types from scratch using real LP
working tolerance for feasibility in floating-point solver during iterative refinement ...
int dmaxSizePrimalRational(const int base=2)
get size of largest denominator in primal solution
bool getSlacksRational(VectorRational &vector)
gets the vector of slack values if available; returns true on success
type of algorithm, i.e., primal or dual
round scaling factors for iterative refinement to powers of two?
void _addColsReal(const LPColSetReal &lpcolset)
adds multiple columns to the real LP and adjusts basis
void _checkBasisScaling()
check correctness of (un)scaled basis matrix operations
LP geometric mean scaling.
int numRowsRational() const
bool getPrimalReal(R *p_vector, int size)
Simplex basis.Consider the linear program as provided from class SPxLP: where , and ...
bool getBasisMetric(R &metric, int type=0)
bool getPrimalRational(VectorRational &vector)
SPxBasisBase< R > _decompTransBasis
Abstract ratio test base class.
dual feasibility tolerance
SPxAutoPR< R > _pricerAuto
static struct soplex::SoPlexBase::Settings::BoolParam boolParam
Dynamic sparse vectors.Class DSVectorBase implements dynamic sparse vectors, i.e. SVectorBases with a...
void _getCompatibleColumns(VectorBase< R > feasVector, int *nonposind, int *compatind, int *rowsforremoval, int *colsforremoval, int nnonposind, int *ncompatind, bool formRedProb, bool &stop)
retrieves the compatible columns from the constraint matrix
const SVectorRational & colVectorRational(int i) const
returns vector of column i
Least squares scaling.This SPxScaler implementation performs least squares scaling as suggested by Cu...
int * _decompViolatedBounds
void _updateComplementaryDualFixedPrimalVars(int *currFixedVars)
updating the dual columns related to the fixed primal variables.
std::string description[SoPlexBase< R >::BOOLPARAM_COUNT]
array of descriptions for boolean parameters
Safe arrays of arbitrary types.Class Array provides safe arrays of arbitrary type. Array elements are accessed just like ordinary C++ array elements by means of the index operator[](). Safety is provided by.
DataArray< SPxColId > _decompCompPrimalColIDs
void setIntegralityInformation(int ncols, int *intInfo)
pass integrality information about the variables to the solver
SPxSolverBase< R > _solver
VectorBase< R > _manualLower
bool writeDualFileReal(const char *filename, const NameSet *rowNames=0, const NameSet *colNames=0, const DIdxSet *intvars=0) const
writes the dual of the real LP to file; LP or MPS format is chosen from the extension in filename; if...
Implementation of Sparse Linear Solver.
void _updateDecompReducedProblem(R objVal, VectorBase< R > dualVector, VectorBase< R > redcostVector, VectorBase< R > compPrimalVector, VectorBase< R > compDualVector)
update the reduced problem with additional columns and rows
bool computeBasisInverseRational()
compute rational basis inverse; returns true on success
Implementation of Sparse Linear Solver with Rational precision.
SoPlexBase< R > & operator=(const SoPlexBase< R > &rhs)
assignment operator
bool getBasisInverseTimesVecRational(const SVectorRational &rhs, SSVectorRational &sol)
computes solution of basis matrix B * sol = rhs; performs rational factorization if not available; re...
void changeColRational(int i, const LPColRational &lpcol)
replaces column i with lpcol
SoPlex start basis generation base class.SPxStarter is the virtual base class for classes generating ...
const VectorRational & lhsRational() const
returns left-hand side vector
Settings()
default constructor initializing default settings
number of boolean parameters
bool isDualFeasible() const
is stored dual solution feasible?
void _factorizeColumnRational(SolRational &sol, DataArray< typename SPxSolverBase< R >::VarStatus > &basisStatusRows, DataArray< typename SPxSolverBase< R >::VarStatus > &basisStatusCols, bool &stoppedTime, bool &stoppedIter, bool &error, bool &optimal)
factorizes rational basis matrix in column representation
bool getDual(VectorBase< R > &vector)
gets the dual solution vector if available; returns true on success
void syncLPRational()
synchronizes rational LP with real LP, i.e., copies real LP to rational LP, if sync mode is manual ...
DataArray< SPxRowId > _decompElimPrimalRowIDs
bool decompTerminate(R timeLimit)
function call to terminate the decomposition simplex
void removeColRangeRational(int start, int end, int perm[]=0)
removes columns start to end including both; an array perm of size numColsRational() may be passed as...
const char * getStarterName()
name of starter
DataArray< typename SPxSolverBase< R >::VarStatus > _storedBasisStatusCols
General methods in LP preprocessing.
SLUFactorRational _rationalLUSolver
R objValueReal()
returns the objective value if a primal solution is available
VectorRational _unboundedLhs
column representation Ax - s = 0, lower <= x <= upper, lhs <= s <= rhs
int totalSizePrimalRational(const int base=2)
get size of primal solution
DataArray< SPxColId > _decompPrimalColIDs
Ids for LP columns.Class SPxColId provides DataKeys for the column indices of an SPxLP.
void changeLowerRational(const VectorRational &lower)
changes vector of lower bounds to lower
void _updateComplementaryPrimalSlackColCoeff()
updating the slack column coefficients to adjust for equality constraints
void _removeColsReal(int perm[])
removes all columns with an index i such that perm[i] < 0; upon completion, perm[i] >= 0 indicates th...
void _changeBoundsReal(const VectorBase< R > &lower, const VectorBase< R > &upper)
changes vectors of column bounds to lower and upper and adjusts basis
decide according to READMODE
refactor threshold for nonzeros in last factorized basis matrix compared to updated basis matrix ...
bool writeBasisFile(const char *filename, const NameSet *rowNames=0, const NameSet *colNames=0, const bool cpxFormat=false) const
writes basis information to filename; if rowNames and colNames are NULL, default names are used; retu...
Auto pricer.This pricer switches between Devex and Steepest edge pricer based on the difficulty of th...
void _project(SolRational &sol)
undoes lifting
SPxVectorST< R > _starterVector
bool getRowViolationRational(Rational &maxviol, Rational &sumviol)
the verbosity of the decomposition based simplex
SPxSolverBase< R >::Status optimize(volatile bool *interrupt=NULL)
optimize the given LP
DataArray< SPxColId > _decompCompPrimalVarBoundIDs
bool getEstimatedCondition(R &condition)
computes an estimated condition number for the current basis matrix using the power method; returns t...
const VectorRational & maxObjRational() const
returns objective function vector after transformation to a maximization problem; since this is how i...
void changeObjRational(const VectorRational &obj)
changes objective function vector to obj
void _getZeroDualMultiplierIndices(VectorBase< R > feasVector, int *nonposind, int *colsforremoval, int *nnonposind, bool &stop)
identifies the columns of the row-form basis that correspond to rows with zero dual multipliers...
void _identifyComplementaryDualFixedPrimalVars(int *currFixedVars)
removing the dual columns related to the fixed variables
should dual infeasibility be tested in order to try to return a dual solution even if primal infeasib...
bool getDualReal(R *p_vector, int dim)
bool getPrimalRay(VectorBase< R > &vector)
gets the primal ray if available; returns true on success
const SVectorRational & rowVectorRational(int i) const
returns vector of row i
void _changeUpperReal(const VectorBase< R > &upper)
changes vector of upper bounds to upper and adjusts basis
void changeElementReal(int i, int j, const R &val)
changes matrix entry in row i and column j to val
bool hasSol() const
is a solution available (not neccessarily feasible)?
void addRowsReal(const LPRowSetBase< R > &lprowset)
adds multiple rows
void _resolveWithoutPreprocessing(typename SPxSimplifier< R >::Result simplificationStatus)
loads original problem into solver and solves again after it has been solved to optimality with prepr...
bool hasBasis() const
is an advanced starting basis available?
bool getRedCostViolation(R &maxviol, R &sumviol)
gets violation of reduced costs; returns true on success
void _syncLPRational(bool time=true)
synchronizes rational LP with real LP, i.e., copies real LP to rational LP, without looking at the sy...
LP simplification base class.
Saving LPs in a form suitable for SoPlex.
minimize number of basic slack variables, i.e. more variables between bounds
Statistics * _statistics
statistics since last call to solveReal() or solveRational()
R maxObjReal(int i) const
returns objective value of column i after transformation to a maximization problem; since this is how...
void getOriginalProblemBasisColStatus(int &nNonBasicCols)
function to retrieve the column status for the original problem basis from the reduced and complement...
should cycling solutions be accepted during iterative refinement?
mode for a posteriori feasibility checks
void changeRangeRational(const VectorRational &lhs, const VectorRational &rhs)
changes left- and right-hand side vectors
Partial multiple pricing.Class SPxParMultPr is an implementation class for SPxPricer implementing Dan...
void printStatus(std::ostream &os, typename SPxSolverBase< R >::Status status)
prints status
const VectorBase< R > & maxObjRealInternal() const
returns objective function vector after transformation to a maximization problem; since this is how i...
bool parseSettingsString(char *str)
parses one setting string and returns true on success; note that string is modified ...
equilibrium scaling on rows or columns
RangeType
type of bounds and sides
should lifting be used to reduce range of nonzero matrix coefficients?
const VectorRational & rhsRational() const
returns right-hand side vector
void _removeComplementaryDualFixedPrimalVars(int *currFixedVars)
removing the dual columns related to the fixed variables
bound flipping ratio test for long steps in the dual simplex
Rational maxAbsNonzeroRational() const
returns biggest non-zero element in absolute value
void printShortStatistics(std::ostream &os)
prints short statistics
void getRhsReal(VectorBase< R > &rhs) const
gets right-hand side vector
DataArray< SPxColId > _decompVarBoundDualIDs
void _idxToPerm(int *idx, int idxSize, int *perm, int permSize) const
creates a permutation for removing rows/columns from an array of indices
void changeRowRational(int i, const LPRowRational &lprow)
replaces row i with lprow
Real solveTime() const
time spent in last call to solve
Fast shifting ratio test.Class SPxFastRT is an implementation class of SPxRatioTester providing fast ...
bool _reapplyPersistentScaling() const
check whether persistent scaling is supposed to be reapplied again after unscaling ...
decide according to problem size
const RowViolation * entry
should the degeneracy be computed for each basis?
void getColsRational(int start, int end, LPColSetRational &lpcolset) const
gets columns start, ..., end
should a rational factorization be performed after iterative refinement?
apply standard floating-point algorithm
sparse pricing threshold (#violations < dimension * SPARSITY_THRESHOLD activates sparse pricing) ...
void _addRowReal(const LPRowBase< R > &lprow)
adds a single row to the real LP and adjusts basis
bool getRedCostReal(R *vector, int dim)
row representation (lower,lhs) <= (x,Ax) <= (upper,rhs)
DSVectorRational _primalDualDiff
should the dual of the complementary problem be used in the decomposition simplex?
R minAbsNonzeroReal() const
returns smallest non-zero element in absolute value
DataArray< SPxRowId > _decompDualRowIDs
DataArray< RangeType > _rowTypes
apply rational reconstruction after each iterative refinement?
bool getDualViolationRational(Rational &maxviol, Rational &sumviol)
void _checkScaling(SPxLPBase< R > *origLP) const
check scaling of LP
partial multiple pricer based on Dantzig pricing
bool getPrimalRayReal(R *vector, int dim)
void _removeRowRangeReal(int start, int end, int perm[])
removes rows start to end including both; an array perm of size numRows() may be passed as buffer mem...
bool getRedCostRational(VectorRational &vector)
VectorBase< R > _manualObj
SPxLPBase< R > * _decompLP
const char * getScalerName()
name of scaling method
void removeRowRational(int i)
removes row i
SPxSolverBase< R >::Status solve(volatile bool *interrupt=NULL)
void _verifySolutionReal()
verify computed solution and resolve if necessary
DataArray< SPxColId > _decompCompPrimalFixedVarIDs
SPxSolverBase< R >::Status _solveRealStable(bool acceptUnbounded, bool acceptInfeasible, VectorBase< R > &primal, VectorBase< R > &dual, DataArray< typename SPxSolverBase< R >::VarStatus > &basisStatusRows, DataArray< typename SPxSolverBase< R >::VarStatus > &basisStatusCols, const bool forceNoSimplifier=false)
solves real LP with recovery mechanism
DataArray< SPxRowId > _decompReducedProbRowIDs
SPxSolverBase< R >::VarStatus basisColStatus(int col) const
returns basis status for a single column
void getColVectorReal(int i, DSVectorBase< R > &col) const
gets vector of col i
SPxLeastSqSC< R > _scalerLeastsq
equilibrium scaling on rows and columns
SPxSolverBase< R >::VarStatus basisRowStatus(int row) const
returns basis status for a single row
void changeBoundsRational(const VectorRational &lower, const VectorRational &upper)
changes vectors of column bounds to lower and upper
LP simplification abstract base class.Instances of classes derived from SPxSimplifier may be loaded t...
void removeRowsReal(int perm[])
removes all rows with an index i such that perm[i] < 0; upon completion, perm[i] >= 0 indicates the n...
int dlcmSizeDualRational(const int base=2)
get size of least common multiple of denominators in dual solution
bool _lowerFinite(const RangeType &rangeType) const
checks whether RangeType corresponds to finite lower bound
void _solveDecompositionDualSimplex()
solves LP using the decomposition based dual simplex
void _addColReal(const LPColReal &lpcol)
adds a single column to the real LP and adjusts basis
void _rangeToPerm(int start, int end, int *perm, int permSize) const
creates a permutation for removing rows/columns from a range of indices
void changeLowerReal(const VectorBase< R > &lower)
changes vector of lower bounds to lower
void getColRational(int i, LPColRational &lpcol) const
gets column i
type of timer for statistics
int totalSizeDualRational(const int base=2)
get size of dual solution
R lowerReal(int i) const
returns lower bound of column i
bool getDualFarkas(VectorBase< R > &vector)
gets the Farkas proof if available; returns true on success
steepest edge pricer with initialization to unit norms
void _unscaleSolutionReal(SPxLPBase< R > &LP, bool persistent=true)
unscales stored solution to remove internal or external scaling of LP
Wrapper for several output streams. A verbosity level is used to decide which stream to use and wheth...
decide depending on tolerances whether to apply iterative refinement
iteration limit (-1 if unlimited)
bool getBoundViolation(R &maxviol, R &sumviol)
gets violation of bounds; returns true on success
void _deleteAndUpdateRowsComplementaryProblem(SPxRowId rangedRowIds[], int &naddedrows)
removing rows from the complementary problem.
type of computational form, i.e., column or row representation
void _updateComplementaryDualSlackColCoeff()
updating the slack column coefficients to adjust for equality constraints
primal feasibility tolerance
bool hasDual() const
deprecated: use hasSol() instead
void _removeColRangeReal(int start, int end, int perm[])
removes columns start to end including both; an array perm of size numColsReal() may be passed as buf...
void getObjRational(VectorRational &obj) const
gets objective function vector
void _decompSimplifyAndSolve(SPxSolverBase< R > &solver, SLUFactor< R > &sluFactor, bool fromScratch, bool applyPreprocessing)
simplifies the problem and solves
SPxLPRational * _rationalLP
void _untransformUnbounded(SolRational &sol, bool unbounded)
undoes transformation to unboundedness problem
void _recomputeRangeTypesRational()
recomputes range types from scratch using rational LP
void removeColsRational(int perm[])
removes all columns with an index i such that perm[i] < 0; upon completion, perm[i] >= 0 indicates th...
bool saveSettingsFile(const char *filename, const bool onlyChanged=false, int solvemode=1) const
writes settings file; returns true on success
const VectorBase< R > & lowerRealInternal() const
returns lower bound vector
void removeColRangeReal(int start, int end, int perm[]=0)
removes columns start to end including both; an array perm of size numColsReal() may be passed as buf...
bool * _decompReducedProbCols
bool getPrimal(VectorBase< R > &vector)
gets the primal solution vector if available; returns true on success
void _computeBasisInverseRational()
computes rational inverse of basis matrix as defined by _rationalLUSolverBind
LPRowBase< R >::Type rowTypeReal(int i) const
returns inequality type of row i
bool hasPrimalRay() const
is a primal unbounded ray available?
void _solveRealLPAndRecordStatistics(volatile bool *interrupt=NULL)
call floating-point solver and update statistics on iterations etc.
void changeRhsReal(const VectorBase< R > &rhs)
changes right-hand side vector to rhs
bool getDecompBoundViolation(R &maxviol, R &sumviol)
gets violation of bounds; returns true on success
bool readBasisFile(const char *filename, const NameSet *rowNames=0, const NameSet *colNames=0)
reads basis information from filename and returns true on success; if rowNames and colNames are NULL...
upper limit on objective value
should LP be transformed to equality form before a rational solve?
void _preprocessAndSolveReal(bool applyPreprocessing, volatile bool *interrupt=NULL)
solves real LP with/without preprocessing
void _changeElementReal(int i, int j, const R &val)
changes matrix entry in row i and column j to val and adjusts basis
print condition number during the solve
void syncLPReal()
synchronizes real LP with rational LP, i.e., copies (rounded) rational LP into real LP...
bool hasPrimal() const
deprecated: use hasSol() instead
SLUFactor< R > _slufactor
void _performFeasIRStable(SolRational &sol, bool &withDualFarkas, bool &stoppedTime, bool &stoppedIter, bool &error)
performs iterative refinement on the auxiliary problem for testing feasibility
bool multBasis(R *vec, bool unscale=true)
multiply with basis matrix; B * vec (inplace)
void _setComplementaryDualOriginalObjective()
updating the complementary dual problem with the original objective function
SPxSolverBase< R > _compSolver
void _addRowsReal(const LPRowSetBase< R > &lprowset)
adds multiple rows to the real LP and adjusts basis
Solution vector based start basis.This version of SPxWeightST can be used to construct a starting bas...
void _evaluateSolutionDecomp(SPxSolverBase< R > &solver, SLUFactor< R > &sluFactor, typename SPxSimplifier< R >::Result result)
evaluates the solution of the reduced problem for the DBDS
SPxDevexPR< R > _pricerDevex
BoolParam
boolean parameters
bool getDualFarkasRational(VectorRational &vector)
const Settings & settings() const
returns current parameter settings
int numCols() const
Templated function that returns number of columns.
Dynamic index set.Class DIdxSet provides dynamic IdxSet in the sense, that no restrictions are posed ...
void setRandomSeed(unsigned int seed)
set the random seeds of the solver instance
Rational _rationalPosInfty
IntParam
integer parameters
SPxBasisBase< R >::SPxStatus basisStatus() const
returns the current basis status
unsigned int randomSeed() const
returns the current random seed of the solver instance
DataArray< RangeType > _colTypes
SPxHarrisRT< R > _ratiotesterHarris
SPxEquiliSC< R > _scalerBiequi
void _invalidateSolution()
invalidates solution
bool getBasisIndRational(DataArray< int > &bind)
gets an array of indices for the columns of the rational basis matrix; bind[i] >= 0 means that the i-...
void _transformEquality()
introduces slack variables to transform inequality constraints into equations for both rational and r...
void printStatistics(std::ostream &os)
prints complete statistics
Simple heuristic SPxStarter.
void removeRowReal(int i)
removes row i
refactor threshold for memory growth in factorization since last refactorization
Class for storing a primal-dual solution with basis information.
only error and warning output
bool getRedCostViolationRational(Rational &maxviol, Rational &sumviol)
try to enforce that the optimal solution is a basic solutiong
SPxSolverBase< R >::Status _status
void writeStateReal(const char *filename, const NameSet *rowNames=0, const NameSet *colNames=0, const bool cpxFormat=false) const
writes internal LP, basis information, and parameter settings; if rowNames and colNames are NULL...
R getCompSlackVarCoeff(int primalRowNum)
gets the coefficient of the slack variable in the primal complementary problem
void removeColReal(int i)
removes column i
void _setComplementaryPrimalOriginalObjective()
updating the complementary primal problem with the original objective function
(In)equality for LPs.Class LPRowBase provides constraints for linear programs in the form where a is...
void removeRowsRational(int perm[])
removes all rows with an index i such that perm[i] < 0; upon completion, perm[i] >= 0 indicates the n...
int * _decompCompProbColIDsIdx
bool setRealParam(const RealParam param, const Real value, const bool init=true)
sets real parameter value; returns true on success
int * _decompViolatedRows
type of starter used to create crash basis
bool getSlacksReal(VectorBase< R > &vector)
gets the vector of slack values if available; returns true on success
void removeColRational(int i)
removes column i
SPxSimplifier< R > * _simplifier
Sparse vector .A UnitVectorBase is an SVectorBase that can take only one nonzero value with value 1 b...
Bound flipping ratio test (long step dual) for SoPlex.
void _updateDecompReducedProblemViol(bool allrows)
update the reduced problem with additional columns and rows based upon the violated original bounds a...
VectorRational _unboundedLower
bool _isConsistent() const
checks consistency
VectorBase< R > _manualLhs
const VectorBase< R > & upperRealInternal() const
returns upper bound vector
Debugging, floating point type and parameter definitions.
void _ensureDSVectorRationalMemory(DSVectorRational &vec, const int newmax) const
extends sparse vector to hold newmax entries if and only if it holds no more free entries ...
bool getBasisInverseColRational(const int c, SSVectorRational &vec)
computes column c of basis inverse; performs rational factorization if not available; returns true on...
Set of strings.Class NameSet implements a symbol or name table. It allows to store or remove names (i...
mode for solution polishing
bool _parseSettingsLine(char *line, const int lineNumber)
parses one line in a settings file and returns true on success; note that the string is modified ...
mode for hyper sparse pricing
RangeType _switchRangeType(const RangeType &rangeType) const
switches RANGETYPE_LOWER to RANGETYPE_UPPER and vice versa
bool getDualViolation(R &maxviol, R &sumviol)
gets violation of dual multipliers; returns true on success
pivot zero tolerance used in factorization
SPxSolverBase< R >::Status _solveRealForRational(bool fromscratch, VectorBase< R > &primal, VectorBase< R > &dual, DataArray< typename SPxSolverBase< R >::VarStatus > &basisStatusRows, DataArray< typename SPxSolverBase< R >::VarStatus > &basisStatusCols)
solves real LP during iterative refinement
void _evaluateSolutionReal(typename SPxSimplifier< R >::Result simplificationStatus)
checks result of the solving process and solves again without preprocessing if necessary ...
bool _upperFinite(const RangeType &rangeType) const
checks whether RangeType corresponds to finite upper bound
const VectorBase< R > & rhsRealInternal() const
returns right-hand side vector, ignoring scaling
LP least squares scaling.
geometric mean scaling (max 8 rounds) followed by equilibrium scaling (rows and columns) ...
static struct soplex::SoPlexBase::Settings::RealParam realParam
number of real parameters
upper bound is finite, lower bound is infinite
Collection of dense, sparse, and semi-sparse vectors.
RangeType _rangeTypeReal(const R &lower, const R &upper) const
determines RangeType from real bounds
RangeType _rangeTypeRational(const Rational &lower, const Rational &upper) const
determines RangeType from rational bounds
Implementation of Sparse Linear Solver.This class implements a SLinSolver interface by using the spar...
void getOriginalProblemStatistics()
stores the problem statistics of the original problem
bool defaultValue[SoPlexBase< R >::BOOLPARAM_COUNT]
array of default values for boolean parameters
bool _reconstructSolutionRational(SolRational &sol, DataArray< typename SPxSolverBase< R >::VarStatus > &basisStatusRows, DataArray< typename SPxSolverBase< R >::VarStatus > &basisStatusCols, const Rational &denomBoundSquared)
attempts rational reconstruction of primal-dual solution
SPxSteepExPR< R > _pricerSteep
void resetSettings(const bool quiet=false, const bool init=true)
resets default parameter settings
Everything should be within this namespace.
force iterative refinement
void removeRowRangeRational(int start, int end, int perm[]=0)
removes rows start to end including both; an array perm of size numRowsRational() may be passed as bu...
SLUFactor< R > _compSlufactor
lower bound is finite, upper bound is infinite
returns the current git hash of SoPlex
mode for iterative refinement strategy
maximum increase of scaling factors between refinements
Rational _rationalMaxscaleincr
maximize number of basic slack variables, i.e. more variables on bounds
Dantzig pricer.Class SPxDantzigPR is an implementation class of an SPxPricer implementing Dantzig's d...
Harris pricing with shifting.
void _storeSolutionRealFromPresol()
stores solution from the simplifier because problem vanished in presolving step
void printOriginalProblemStatistics(std::ostream &os)
stores the problem statistics of the original problem
zero tolerance used in update of the factorization
bool getRowViolation(R &maxviol, R &sumviol)
gets violation of constraints; returns true on success
void changeColReal(int i, const LPColReal &lpcol)
replaces column i with lpcol
R objReal(int i) const
returns objective value of column i
void changeLhsReal(const VectorBase< R > &lhs)
changes left-hand side vector for constraints to lhs
void _storeSolutionReal(bool verify=true)
stores solution of the real LP; before calling this, the real LP must be loaded in the solver and sol...
int numRows() const
returns number of rows
accuracy of conjugate gradient method in least squares scaling (higher value leads to more iterations...
Rational minAbsNonzeroRational() const
returns smallest non-zero element in absolute value
Saving LPs in a form suitable for SoPlex.Class SPxLPBase provides the data structures required for sa...
void _syncRealSolution()
synchronizes rational solution with real solution, i.e., copies (rounded) rational solution to real s...
void addRowRational(const LPRowRational &lprow)
adds a single row
time limit in seconds (INFTY if unlimited)
Steepest edge pricer with exact initialization of weights.
bool getDecompRowViolation(R &maxviol, R &sumviol)
gets violation of constraints; returns true on success
void changeElementRational(int i, int j, const Rational &val)
changes matrix entry in row i and column j to val
user sync of real and rational LP
DataArray< typename SPxSolverBase< R >::VarStatus > _basisStatusRows
Type
(In)Equality type of an LP row.
Forrest-Tomlin type update.
std::string statisticString() const
statistical information in form of a string
mode for reading LP files
use bound flipping also for row representation?
int getOrigVarFixedDirection(int colNum)
determining which bound the primal variables will be fixed to.
bool setSettings(const Settings &newSettings, const bool init=true)
sets parameter settings; returns true on success
void _updateComplementaryPrimalFixedPrimalVars(int *currFixedVars)
updating the dual columns related to the fixed primal variables.
DataArray< typename SPxSolverBase< R >::VarStatus > _storedBasisStatusRows
const SVectorBase< R > & colVectorRealInternal(int i) const
returns vector of col i, ignoring scaling
Simple heuristic SPxStarter.Testing version of an SPxVectorST using a very simplistic heuristic to bu...
geometric mean scaling on rows and columns, max 8 rounds
Weighted start basis.Class SPxWeightST is an implementation of a SPxStarter for generating a Simplex ...
SPxBoundFlippingRT< R > _ratiotesterBoundFlipping
bool readFile(const char *filename, NameSet *rowNames=0, NameSet *colNames=0, DIdxSet *intVars=0)
reads LP file in LP or MPS format according to READMODE parameter; gets row names, column names, and integer variables if desired; returns true on success
SPxGeometSC< R > _scalerGeoequi
re-optimize the original problem to get a proof (ray) of infeasibility/unboundedness?
should row and bound violations be computed explicitly in the update of reduced problem in the decomp...
SPxStarter< R > * _starter
LPRowSetBase< R > _transformedRows
Partial multiple pricing.
minimum number of stalling refinements since last pivot to trigger rational factorization ...
void getRowRational(int i, LPRowRational &lprow) const
gets row i
Verbosity
Verbosity level.
void addColReal(const LPColBase< R > &lpcol)
adds a single column
SoPlexBase< Real > SoPlex
lower bound equals upper bound
Real _realParamValues[SoPlexBase< R >::REALPARAM_COUNT]
array of current real parameter values
DataArray< SPxColId > _decompDualColIDs
stalling refinement limit (-1 if unlimited)
Textbook ratio test for SoPlex.Class SPxDefaultRT provides an implementation of the textbook ratio te...
void changeBoundsReal(const VectorBase< R > &lower, const VectorBase< R > &upper)
changes vectors of column bounds to lower and upper
void printVersion() const
prints version and compilation options
int _intParamValues[SoPlexBase< R >::INTPARAM_COUNT]
array of current integer parameter values
void _transformUnbounded()
transforms LP to unboundedness problem by moving the objective function to the constraints, changing right-hand side and bounds to zero, and adding an auxiliary variable for the decrease in the objective function
bool areLPsInSync(const bool checkVecVals=true, const bool checkMatVals=false, const bool quiet=false) const
checks if real LP and rational LP are in sync; dimensions will always be compared, vector and matrix values only if the respective parameter is set to true. If quiet is set to true the function will only display which vectors are different.
const SVectorBase< R > & rowVectorRealInternal(int i) const
returns vector of row i, ignoring scaling
void _checkOriginalProblemOptimality(VectorBase< R > primalVector, bool printViol)
checking the optimality of the original problem.
LP scaler abstract base class.Instances of classes derived from SPxScaler may be loaded to SoPlex in ...
bool getPrimalRayRational(VectorRational &vector)
Set of LP rows.Class LPRowSetBase implements a set of LPRowBase%s. Unless for memory limitations...
VectorRational _unboundedUpper
void printDecompDisplayLine(SPxSolverBase< R > &solver, const SPxOut::Verbosity origVerb, bool force, bool forceHead)
prints a display line of the flying table for the DBDS
DataArray< typename SPxSolverBase< R >::VarStatus > _basisStatusCols
void _storeBasis()
store basis
VectorRational _feasUpper
void _createDecompReducedAndComplementaryProblems()
creating copies of the original problem that will be manipulated to form the reduced and complementar...
int numIterations() const
number of iterations since last call to solve
void changeRowReal(int i, const LPRowBase< R > &lprow)
replaces row i with lprow
int dmaxSizeDualRational(const int base=2)
get size of largest denominator in dual solution
automatic choice according to number of rows and columns
void _changeLowerReal(const VectorBase< R > &lower)
changes vector of lower bounds to lower and adjusts basis
SoPlex chooses automatically (currently always "internal")
VectorBase< R > _transformedObj
only error, warning, and debug output
bool multBasisTranspose(R *vec, bool unscale=true)
multiply with transpose of basis matrix; vec * B^T (inplace)
dual simplex algorithm, i.e., leaving for column and entering for row representation ...
void _syncLPReal(bool time=true)
synchronizes real LP with rational LP, i.e., copies (rounded) rational LP into real LP...
LPColSetRational _slackCols
minimal modification threshold to apply presolve reductions
void setTimings(const Timer::TYPE ttype)
set statistic timers to a certain type
void writeStateRational(const char *filename, const NameSet *rowNames=0, const NameSet *colNames=0, const bool cpxFormat=false) const
writes internal LP, basis information, and parameter settings; if rowNames and colNames are NULL...
void changeUpperReal(const VectorBase< R > &upper)
changes vector of upper bounds to upper
void clearBasis()
clears starting basis
continue iterative refinement with exact basic solution if not optimal?
standard floating-point parsing
perturb the entire problem or only the relevant bounds of s single pivot?
void _ensureRationalLP()
ensures that the rational LP is available; performs no sync
void _restoreBasis()
restore basis
DataArray< int > _rationalLUSolverBind
DataArray< SPxColId > _decompReducedProbColIDs
using internal presolving methods
DualSign getExpectedDualVariableSign(int rowNumber)
returns the expected sign of the dual variables for the reduced problem
R rhsReal(int i) const
returns right-hand side of row i
SPxDantzigPR< R > _pricerDantzig
minimal reduction (sum of removed rows/cols) to continue simplification
geometric frequency at which to apply rational reconstruction
int numColsRational() const
bool isPrimalFeasible() const
is stored primal solution feasible?
void getOriginalProblemBasisRowStatus(DataArray< int > °enerateRowNums, DataArray< typename SPxSolverBase< R >::VarStatus > °enerateRowStatus, int &nDegenerateRows, int &nNonBasicRows)
function to retrieve the original problem row basis status from the reduced and complementary problem...
SPxGeometSC< R > _scalerGeo1
Class for storing a primal-dual solution with basis information.
void _changeRowReal(int i, const LPRowBase< R > &lprow)
replaces row i with lprow and adjusts basis
Sparse vectors.Class SVectorBase provides packed sparse vectors. Such are a sparse vectors...
bool _readFileReal(const char *filename, NameSet *rowNames=0, NameSet *colNames=0, DIdxSet *intVars=0)
reads real LP in LP or MPS format from file and returns true on success; gets row names...
VectorRational _unboundedRhs
SPxGeometSC< R > _scalerGeo8
number of integer parameters
should the decomposition based dual simplex be used to solve the LP? Setting this to true forces the ...
class of parameter settings
void _restoreLPReal()
restores objective, bounds, and sides of real LP
void _formDecompComplementaryProblem()
forms the complementary problem
Ids for LP rows.Class SPxRowId provides DataKeys for the row indices of an SPxLP. ...
static struct soplex::SoPlexBase::Settings::IntParam intParam
void addColRational(const LPColRational &lpcol)
adds a single column
void changeLhsRational(const VectorRational &lhs)
changes left-hand side vector for constraints to lhs
the iteration frequency at which the decomposition solve output is displayed.
void getUpperReal(VectorBase< R > &upper) const
gets upper bound vector
void _syncRationalSolution()
synchronizes real solution with rational solution, i.e., copies real solution to rational solution ...
refinement limit (-1 if unlimited)
void _transformFeasibility()
transforms LP to feasibility problem by removing the objective function, shifting variables...
void _solveDecompReducedProblem()
solves the reduced problem
SoPlexBase()
default constructor
void _enableSimplifierAndScaler()
enables simplifier and scaler according to current parameters
void getBasisInd(int *bind) const
gets the indices of the basic columns and rows; basic column n gives value n, basic row m gives value...
bool checkBasisDualFeasibility(VectorBase< R > feasVec)
checks the dual feasibility of the current basis
void clearLPRational()
clears the LP
Textbook ratio test for SoPlex.
bool _boolParamValues[SoPlexBase< R >::BOOLPARAM_COUNT]
array of current boolean parameter values
Array< UnitVectorRational *> _unitMatrixRational
upper threshold in lifting (nonzero matrix coefficients with larger absolute value will be reformulat...
void getLhsReal(VectorBase< R > &lhs) const
gets left-hand side vector
the number of iterations before the decomposition simplex initialisation is terminated.
DualSign getOrigProbDualVariableSign(int rowNumber)
returns the expected sign of the dual variables for the original problem
SPxRowId _compSlackDualRowId
standard Harris ratio test
const VectorBase< R > & lhsRealInternal() const
returns left-hand side vector, ignoring scaling
bool writeFileReal(const char *filename, const NameSet *rowNames=0, const NameSet *colNames=0, const DIdxSet *intvars=0, const bool unscale=true) const
bool setIntParam(const IntParam param, const int value, const bool init=true)
sets integer parameter value; returns true on success
mode for synchronizing real and rational LP
lower limit on objective value
using the presolve lib papilo
modified Harris ratio test
Rational objValueRational()
returns the objective value if a primal solution is available
void changeRangeReal(const VectorBase< R > &lhs, const VectorBase< R > &rhs)
changes left- and right-hand side vectors
void _computeInfeasBox(SolRational &sol, bool transformed)
void _untransformEquality(SolRational &sol)
undoes transformation to equality form
Rational _rationalFeastol
SPxWeightST< R > _starterWeight
void _decompResolveWithoutPreprocessing(SPxSolverBase< R > &solver, SLUFactor< R > &sluFactor, typename SPxSimplifier< R >::Result result)
loads original problem into solver and solves again after it has been solved to optimality with prepr...
Equilibrium row/column scaling.This SPxScaler implementation performs equilibrium scaling of the LPs ...
void _storeLPReal()
stores objective, bounds, and sides of real LP
bool * _decompReducedProbRows
void addColsReal(const LPColSetBase< R > &lpcolset)
adds multiple columns
bool getBasisInverseColReal(int c, R *coef, int *inds=NULL, int *ninds=NULL, bool unscale=true)
computes column c of basis inverse; returns true on success
void _performOptIRStable(SolRational &sol, bool acceptUnbounded, bool acceptInfeasible, int minRounds, bool &primalFeasible, bool &dualFeasible, bool &infeasible, bool &unbounded, bool &stoppedTime, bool &stoppedIter, bool &error)
solves current problem with iterative refinement and recovery mechanism
greedy crash basis weighted by objective, bounds, and sides
VectorBase< R > _manualRhs
void changeObjReal(const VectorBase< R > &obj)
changes objective function vector to obj
bool writeFileRational(const char *filename, const NameSet *rowNames=0, const NameSet *colNames=0, const DIdxSet *intvars=0) const
void getLowerReal(VectorBase< R > &lower) const
gets lower bound vector
zero tolerance used in factorization
R lhsReal(int i) const
returns left-hand side of row i
generic solution-based crash basis
bool getBasisInverseTimesVecReal(R *rhs, R *sol, bool unscale=true)
computes dense solution of basis matrix B * sol = rhs; returns true on success
void getObjReal(VectorBase< R > &obj) const
gets objective function vector
LP column.Class LPColBase provides a datatype for storing the column of an LP a the form similar to ...
void printSolvingStatistics(std::ostream &os)
prints statistics on solving process
bool getBasisInverseRowReal(int r, R *coef, int *inds=NULL, int *ninds=NULL, bool unscale=true)
computes row r of basis inverse; returns true on success
bool hasDualFarkas() const
is Farkas proof of infeasibility available?
DataArray< SPxRowId > _decompPrimalRowIDs
bool loadSettingsFile(const char *filename)
reads settings file; returns true on success
void _changeRhsReal(const VectorBase< R > &rhs)
changes right-hand side vector to rhs and adjusts basis
void _removeRowsReal(int perm[])
removes all rows with an index i such that perm[i] < 0; upon completion, perm[i] >= 0 indicates the n...
Steepest edge pricer.Class SPxSteepPR implements a steepest edge pricer to be used with SoPlex...
bool setDualNorms(int nnormsRow, int nnormsCol, R *norms)
sets steepest edge norms and returns false if that's not possible
void setBasis(const typename SPxSolverBase< R >::VarStatus rows[], const typename SPxSolverBase< R >::VarStatus cols[])
sets starting basis via arrays of statuses
void getRowVectorReal(int i, DSVectorBase< R > &row) const
gets vector of row i
Presol< R > _simplifierPaPILO
void _performUnboundedIRStable(SolRational &sol, bool &hasUnboundedRay, bool &stoppedTime, bool &stoppedIter, bool &error)
performs iterative refinement on the auxiliary problem for testing unboundedness
Harris pricing with shifting.Class SPxHarrisRT is a stable implementation of a SPxRatioTester class a...
Rational objRational(int i) const
returns objective value of column i
bool getRedCost(VectorBase< R > &vector)
gets the vector of reduced cost values if available; returns true on success
DSVectorRational _tauColVector
void _formDecompReducedProblem(bool &stop)
forms the reduced problem
LP simplifier for removing uneccessary row/columns.This SPxSimplifier is mainly based on the paper "P...
void _lift()
reduces matrix coefficient in absolute value by the lifting procedure of Thiele et al...
void _changeColReal(int i, const LPColReal &lpcol)
replaces column i with lpcol and adjusts basis
void _identifyComplementaryPrimalFixedPrimalVars(int *currFixedVars)
removing the dual columns related to the fixed variables
bool _isSolveStopped(bool &stoppedTime, bool &stoppedIter) const
should solving process be stopped?
Settings * _currentSettings
void printUserSettings()
print non-default parameter values