All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
lpi_cpx.c
Go to the documentation of this file.
30 /* CPLEX supports FASTMIP which fastens the lp solving process but therefor it might happen that there will be a loss in
34 /*--+----1----+----2----+----3----+----4----+----5----+----6----+----7----+----8----+----9----+----0----+----1----+----2*/
56 /* this macro is only called in functions returning SCIP_Bool; thus, we return FALSE if there is an error in optimized mode */
68 /* At several places we need to guarantee to have a factorization of an optimal basis and call the simplex to produce
69 * it. In a numerical perfect world, this should need no iterations. However, due to numerical inaccuracies after
70 * refactorization, it might be necessary to do a few extra pivot steps, in particular if FASTMIP is used. */
71 #define CPX_REFACTORMAXITERS 50 /* maximal number of iterations allowed for producing a refactorization of the basis */
73 /* CPLEX seems to ignore bounds with absolute value less than 1e-10. There is no interface define for this constant yet,
77 typedef SCIP_DUALPACKET COLPACKET; /* each column needs two bits of information (basic/on_lower/on_upper) */
79 typedef SCIP_DUALPACKET ROWPACKET; /* each row needs two bit of information (basic/on_lower/on_upper) */
160 SCIP_Real conditionlimit; /**< maximum condition number of LP basis counted as stable (-1.0: no limit) */
166 #if (CPX_VERSION == 1100 || (CPX_VERSION == 1220 && (CPX_SUBVERSION == 0 || CPX_SUBVERSION == 2)))
167 int pseudonthreads; /**< number of threads that SCIP set for the LP solver, but due to CPLEX bug,
358 /* because the basis status values are equally defined in SCIP and CPLEX, they don't need to be transformed */
483 CHECK_ZERO( lpi->messagehdlr, CPXgetintparam(lpi->cpxenv, intparam[i], &(cpxparam->intparval[i])) );
487 CHECK_ZERO( lpi->messagehdlr, CPXgetdblparam(lpi->cpxenv, dblparam[i], &(cpxparam->dblparval[i])) );
539 CHECK_ZERO( lpi->messagehdlr, CPXsetintparam(lpi->cpxenv, intparam[i], lpi->curparam.intparval[i]) );
549 CHECK_ZERO( lpi->messagehdlr, CPXsetdblparam(lpi->cpxenv, dblparam[i], lpi->curparam.dblparval[i]) );
756 * -> To keep SCIP's meaning of the rhs value, we would like to use negative range values: rng := lhs - rhs,
764 * -> Because of this bug, we have to use an additional rhsarray[] for the converted right hand sides and
765 * use rhsarray[i] = lhs[i] and rngarray[i] = rhs[i] - lhs[i] for ranged rows to keep the range values
834 /** converts CPLEX's sen/rhs/rng triplets into SCIP's lhs/rhs pairs, only storing the left hand side */
882 /** converts CPLEX's sen/rhs/rng triplets into SCIP's lhs/rhs pairs, only storing the right hand side */
948 /** after restoring the old lp data in CPLEX we need to resolve the lp to be able to retrieve correct information */
956 /* modifying the LP, restoring the old LP, and loading the old basis is not enough for CPLEX to be able to return the
959 * this may happen after manual strong branching on an integral variable, or after conflict analysis on a strong
960 * branching conflict created a constraint that is not able to modify the LP but trigger the additional call of the
963 * In a numerical perfect world, CPX_REFACTORMAXITERS below should be zero. However, due to numerical inaccuracies
969 SCIPmessagePrintWarning(lpi->messagehdlr, "CPLEX needed %d phase 1 iterations to restore optimal basis.\n", CPXgetphase1cnt(lpi->cpxenv, lpi->cpxlp));
971 SCIPmessagePrintWarning(lpi->messagehdlr, "CPLEX needed %d iterations to restore optimal basis.\n", CPXgetitcnt(lpi->cpxenv, lpi->cpxlp));
998 sprintf(cpxname, "CPLEX %d.%d.%d.%d", CPX_VERSION_VERSION, CPX_VERSION_RELEASE, CPX_VERSION_MODIFICATION, CPX_VERSION_FIX);
1000 sprintf(cpxname, "CPLEX %d.%d.%d.%d", CPX_VERSION/100, (CPX_VERSION%100)/10, CPX_VERSION%10, CPX_SUBVERSION);
1045 assert(sizeof(SCIP_Real) == sizeof(double)); /* CPLEX only works with doubles as floating points */
1057 #if (CPX_VERSION == 1100 || (CPX_VERSION == 1220 && (CPX_SUBVERSION == 0 || CPX_SUBVERSION == 2)))
1058 /* manually set number of threads to 1 to avoid huge system load due to CPLEX bug (version 1100) or segmentation fault (version 1220) */
1062 #if 0 /* turning presolve off seems to be faster than turning it off on demand (if presolve detects infeasibility) */
1206 CHECK_ZERO( lpi->messagehdlr, CPXcopylpwnames(lpi->cpxenv, lpi->cpxlp, ncols, nrows, cpxObjsen(objsen), obj,
1207 lpi->rhsarray, lpi->senarray, beg, cnt, ind, val, lb, ub, lpi->rngarray, colnames, rownames) );
1228 const int* beg, /**< start index of each column in ind- and val-array, or NULL if nnonz == 0 */
1243 CHECK_ZERO( lpi->messagehdlr, CPXaddcols(lpi->cpxenv, lpi->cpxlp, ncols, nnonz, obj, beg, ind, val, lb, ub, colnames) );
1247 CHECK_ZERO( lpi->messagehdlr, CPXnewcols(lpi->cpxenv, lpi->cpxlp, ncols, obj, lb, ub, NULL, colnames) );
1263 assert(0 <= firstcol && firstcol <= lastcol && lastcol < CPXgetnumcols(lpi->cpxenv, lpi->cpxlp));
1274 /** deletes columns from SCIP_LP; the new position of a column must not be greater that its old position */
1326 CHECK_ZERO( lpi->messagehdlr, CPXaddrows(lpi->cpxenv, lpi->cpxlp, 0, nrows, nnonz, lpi->rhsarray, lpi->senarray, beg, ind, val, NULL,
1331 CHECK_ZERO( lpi->messagehdlr, CPXnewrows(lpi->cpxenv, lpi->cpxlp, nrows, lpi->rhsarray, lpi->senarray, NULL, rownames) );
1342 CHECK_ZERO( lpi->messagehdlr, CPXchgrngval(lpi->cpxenv, lpi->cpxlp, rngcount, lpi->rngindarray, lpi->rngarray) );
1358 assert(0 <= firstrow && firstrow <= lastrow && lastrow < CPXgetnumrows(lpi->cpxenv, lpi->cpxlp));
1369 /** deletes rows from SCIP_LP; the new position of a row must not be greater that its old position */
1446 CHECK_ZERO( lpi->messagehdlr, CPXchgbds(lpi->cpxenv, lpi->cpxlp, ncols, ind, lpi->larray, (SCIP_Real*)lb) );
1447 CHECK_ZERO( lpi->messagehdlr, CPXchgbds(lpi->cpxenv, lpi->cpxlp, ncols, ind, lpi->uarray, (SCIP_Real*)ub) );
1496 CHECK_ZERO( lpi->messagehdlr, CPXchgsense(lpi->cpxenv, lpi->cpxlp, nrows, ind, lpi->senarray) );
1509 CHECK_ZERO( lpi->messagehdlr, CPXchgrngval(lpi->cpxenv, lpi->cpxlp, rngcount, lpi->rngindarray, lpi->rngarray) );
1574 /** multiplies a row with a non-zero scalar; for negative scalars, the row's sense is switched accordingly */
1599 SCIP_CALL( SCIPlpiGetRows(lpi, row, row, &lhs, &rhs, &nnonz, &beg, lpi->indarray, lpi->valarray) );
1628 /** multiplies a column with a non-zero scalar; the objective value is multiplied with the scalar, and the bounds
1656 SCIP_CALL( SCIPlpiGetCols(lpi, col, col, &lb, &ub, &nnonz, &beg, lpi->indarray, lpi->valarray) );
1755 /** gets columns from LP problem object; the arrays have to be large enough to store all values
1774 assert(0 <= firstcol && firstcol <= lastcol && lastcol < CPXgetnumcols(lpi->cpxenv, lpi->cpxlp));
1832 assert(0 <= firstrow && firstrow <= lastrow && lastrow < CPXgetnumrows(lpi->cpxenv, lpi->cpxlp));
1840 CHECK_ZERO( lpi->messagehdlr, CPXgetsense(lpi->cpxenv, lpi->cpxlp, lpi->senarray, firstrow, lastrow) );
1841 CHECK_ZERO( lpi->messagehdlr, CPXgetrhs(lpi->cpxenv, lpi->cpxlp, lpi->rhsarray, firstrow, lastrow) );
1884 int namestoragesize, /**< size of namestorage (if 0, storageleft returns the storage needed) */
1897 assert(0 <= firstcol && firstcol <= lastcol && lastcol < CPXgetnumcols(lpi->cpxenv, lpi->cpxlp));
1901 retcode = CPXgetcolname(lpi->cpxenv, lpi->cpxlp, colnames, namestorage, namestoragesize, storageleft, firstcol, lastcol);
1918 int namestoragesize, /**< size of namestorage (if 0, -storageleft returns the storage needed) */
1931 assert(0 <= firstrow && firstrow <= lastrow && lastrow < CPXgetnumrows(lpi->cpxenv, lpi->cpxlp));
1935 retcode = CPXgetrowname(lpi->cpxenv, lpi->cpxlp, rownames, namestorage, namestoragesize, storageleft, firstrow, lastrow);
1953 assert(CPXgetobjsen(lpi->cpxenv, lpi->cpxlp) == CPX_MIN || CPXgetobjsen(lpi->cpxenv, lpi->cpxlp) == CPX_MAX);
1957 *objsen = (CPXgetobjsen(lpi->cpxenv, lpi->cpxlp) == CPX_MIN) ? SCIP_OBJSEN_MINIMIZE : SCIP_OBJSEN_MAXIMIZE;
2032 CHECK_ZERO( lpi->messagehdlr, CPXgetsense(lpi->cpxenv, lpi->cpxlp, lpi->senarray, firstrow, lastrow) );
2033 CHECK_ZERO( lpi->messagehdlr, CPXgetrhs(lpi->cpxenv, lpi->cpxlp, lpi->rhsarray, firstrow, lastrow) );
2105 lpi->iterations = CPXgetphase1cnt(lpi->cpxenv, lpi->cpxlp) + CPXgetitcnt(lpi->cpxenv, lpi->cpxlp);
2118 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, NULL, &primalfeasible, &dualfeasible) );
2128 /* maybe the preprocessor solved the problem; but we need a solution, so solve again without preprocessing */
2129 SCIPdebugMessage("presolver may have solved the problem -> calling CPLEX primal simplex again without presolve\n");
2146 lpi->iterations += CPXgetphase1cnt(lpi->cpxenv, lpi->cpxlp) + CPXgetitcnt(lpi->cpxenv, lpi->cpxlp);
2149 SCIPdebugMessage(" -> CPLEX returned solstat=%d (%d iterations)\n", lpi->solstat, lpi->iterations);
2158 SCIPerrorMessage("CPLEX primal simplex returned CPX_STAT_INForUNBD after presolving was turned off\n");
2162 /* check whether the solution is basic: if Cplex, e.g., hits a time limit in data setup, this might not be the case,
2170 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, &solntype, NULL, NULL) );
2177 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, &solntype, NULL, NULL) );
2210 lpi->iterations = CPXgetphase1cnt(lpi->cpxenv, lpi->cpxlp) + CPXgetitcnt(lpi->cpxenv, lpi->cpxlp);
2221 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, &solntype, NULL, NULL) );
2225 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, NULL, &primalfeasible, &dualfeasible) );
2235 /* maybe the preprocessor solved the problem; but we need a solution, so solve again without preprocessing */
2236 SCIPdebugMessage("presolver may have solved the problem -> calling CPLEX dual simplex again without presolve\n");
2253 lpi->iterations += CPXgetphase1cnt(lpi->cpxenv, lpi->cpxlp) + CPXgetitcnt(lpi->cpxenv, lpi->cpxlp);
2256 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, NULL, &primalfeasible, &dualfeasible) );
2257 SCIPdebugMessage(" -> CPLEX returned solstat=%d (%d iterations)\n", lpi->solstat, lpi->iterations);
2266 SCIPerrorMessage("CPLEX dual simplex returned CPX_STAT_INForUNBD after presolving was turned off\n");
2270 /* check whether the solution is basic: if Cplex, e.g., hits a time limit in data setup, this might not be the case,
2278 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, &solntype, NULL, NULL) );
2285 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, &solntype, NULL, NULL) );
2290 /* this fixes the strange behavior of CPLEX, that in case of the objective limit exceedance, it returns the
2292 * (using this "wrong" dual solution can cause column generation algorithms to fail to find an improving column)
2314 SCIPdebugMessage("dual solution %g does not exceed objective limit [%g,%g] (%d iterations) -> calling CPLEX dual simplex again for one iteration\n",
2336 lpi->iterations += CPXgetphase1cnt(lpi->cpxenv, lpi->cpxlp) + CPXgetitcnt(lpi->cpxenv, lpi->cpxlp);
2358 lpi->iterations += CPXgetphase1cnt(lpi->cpxenv, lpi->cpxlp) + CPXgetitcnt(lpi->cpxenv, lpi->cpxlp);
2361 SCIPdebugMessage(" -> CPLEX returned solstat=%d (%d iterations)\n", lpi->solstat, lpi->iterations);
2369 /** calls barrier or interior point algorithm to solve the LP with crossover to simplex basis */
2405 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, &solntype, NULL, NULL) );
2410 SCIPdebugMessage(" -> CPLEX returned solstat=%d (%d iterations)\n", lpi->solstat, lpi->iterations);
2414 /* maybe the preprocessor solved the problem; but we need a solution, so solve again without preprocessing */
2415 SCIPdebugMessage("CPLEX returned INForUNBD -> calling CPLEX barrier again without presolve\n");
2432 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, &solntype, NULL, NULL) );
2443 SCIPerrorMessage("CPLEX barrier returned CPX_STAT_INForUNBD after presolving was turned off\n");
2506 *down = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJULIM) : getDblParam(lpi, CPX_PARAM_OBJLLIM);
2512 *down = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJLLIM) : getDblParam(lpi, CPX_PARAM_OBJULIM);
2524 *down = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJULIM) : getDblParam(lpi, CPX_PARAM_OBJLLIM);
2533 *up = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJULIM) : getDblParam(lpi, CPX_PARAM_OBJLLIM);
2539 *up = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJLLIM) : getDblParam(lpi, CPX_PARAM_OBJULIM);
2551 *up = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJLLIM) : getDblParam(lpi, CPX_PARAM_OBJULIM);
2642 *down = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJULIM) : getDblParam(lpi, CPX_PARAM_OBJLLIM);
2648 *down = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJLLIM) : getDblParam(lpi, CPX_PARAM_OBJULIM);
2661 *down = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJULIM) : getDblParam(lpi, CPX_PARAM_OBJLLIM);
2670 *up = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJULIM) : getDblParam(lpi, CPX_PARAM_OBJLLIM);
2676 *up = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJLLIM) : getDblParam(lpi, CPX_PARAM_OBJULIM);
2689 *up = objsen == CPX_MIN ? getDblParam(lpi, CPX_PARAM_OBJLLIM) : getDblParam(lpi, CPX_PARAM_OBJULIM);
2710 SCIPdebugMessage("calling CPLEX strongbranching on fractional variable %d (%d iterations)\n", col, itlim);
2778 SCIPdebugMessage("calling CPLEX strongbranching on %d fractional variables (%d iterations)\n", ncols, itlim);
2837 SCIPdebugMessage("calling CPLEX strongbranching on variable %d with integral value (%d iterations)\n", col, itlim);
2844 SCIP_CALL( lpiStrongbranchIntegral(lpi, col, psol, itlim, down, up, downvalid, upvalid, iter) );
2876 SCIPdebugMessage("calling CPLEX strongbranching on %d variables with integer values (%d iterations)\n", ncols, itlim);
2885 SCIP_CALL( lpiStrongbranchIntegral(lpi, cols[j], psols[j], itlim, &(down[j]), &(up[j]), &(downvalid[j]), &(upvalid[j]), iter) );
2930 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, NULL, &pfeas, &dfeas) );
2937 /** returns TRUE iff LP is proven to have a primal unbounded ray (but not necessary a primal feasible point);
2948 return (lpi->solstat == CPX_STAT_UNBOUNDED || lpi->solstat == CPX_STAT_OPTIMAL_FACE_UNBOUNDED);
2951 /** returns TRUE iff LP is proven to have a primal unbounded ray (but not necessary a primal feasible point),
2963 return (lpi->solstat == CPX_STAT_UNBOUNDED && CPXgetmethod(lpi->cpxenv, lpi->cpxlp) == CPX_ALG_PRIMAL);
2982 /* If the solution status of CPLEX is CPX_STAT_UNBOUNDED, it only means, there is an unbounded ray,
2983 * but not necessarily a feasible primal solution. If primalfeasible == FALSE, we cannot conclude,
2986 return ((primalfeasible && (lpi->solstat == CPX_STAT_UNBOUNDED || lpi->solstat == CPX_STAT_INForUNBD))
3006 return (lpi->solstat == CPX_STAT_INFEASIBLE || (lpi->solstat == CPX_STAT_INForUNBD && dualfeasible));
3028 /** returns TRUE iff LP is proven to have a dual unbounded ray (but not necessary a dual feasible point);
3041 /** returns TRUE iff LP is proven to have a dual unbounded ray (but not necessary a dual feasible point),
3053 return (lpi->solstat == CPX_STAT_INFEASIBLE && CPXgetmethod(lpi->cpxenv, lpi->cpxlp) == CPX_ALG_DUAL);
3072 return (dualfeasible && (lpi->solstat == CPX_STAT_INFEASIBLE || lpi->solstat == CPX_STAT_INForUNBD));
3138 /* If the solution status of CPLEX is CPX_STAT_UNBOUNDED, it only means, there is an unbounded ray,
3139 * but not necessarily a feasible primal solution. If primalfeasible == FALSE, we interpret this
3152 /* If the condition number of the basis should be checked, everything above the specified threshold is counted
3161 /* if the kappa could not be computed (e.g., because we do not have a basis), we cannot check the condition */
3215 /** tries to reset the internal status of the LP solver in order to ignore an instability of the last solving call */
3286 CHECK_ZERO( lpi->messagehdlr, CPXsolution(lpi->cpxenv, lpi->cpxlp, &dummy, objval, primsol, dualsol, NULL, redcost) );
3291 CHECK_ZERO( lpi->messagehdlr, CPXgetax(lpi->cpxenv, lpi->cpxlp, activity, 0, CPXgetnumrows(lpi->cpxenv, lpi->cpxlp)-1) );
3352 * Such information is usually only available, if also a (maybe not optimal) solution is available.
3353 * The LPI should return SCIP_INVALID for @p quality, if the requested quantity is not available.
3386 CHECK_ZERO( lpi->messagehdlr, CPXsolninfo(lpi->cpxenv, lpi->cpxlp, NULL, &solntype, NULL, NULL) );
3408 /** gets current basis status for columns and rows; arrays must be large enough to store the basis status */
3427 /* correct rstat values for "<=" constraints: Here CPX_AT_LOWER bound means that the slack is 0, i.e., the upper bound is tight */
3439 /* because the basis status values are equally defined in SCIP and CPLEX, they don't need to be transformed */
3469 /* because the basis status values are equally defined in SCIP and CPLEX, they don't need to be transformed */
3475 /* correct rstat values for ">=" constraints: Here CPX_AT_LOWER bound means that the slack is 0, i.e., the upper bound is tight */
3492 /** returns the indices of the basic columns and rows; basic column n gives value n, basic row m gives value -1-m */
3506 /* this might be turned off if the user as called SCIPlpiClearState() or set SCIP_LPPAR_FROMSCRATCH to TRUE */
3511 if( retval == CPXERR_NO_SOLN || retval == CPXERR_NO_LU_FACTOR || retval == CPXERR_NO_BASIC_SOLN || retval == CPXERR_NO_BASIS )
3536 /* this might be turned off if the user as called SCIPlpiClearState() or set SCIP_LPPAR_FROMSCRATCH to TRUE */
3541 if( retval == CPXERR_NO_SOLN || retval == CPXERR_NO_LU_FACTOR || retval == CPXERR_NO_BASIC_SOLN || retval == CPXERR_NO_BASIS )
3570 /* this might be turned off if the user as called SCIPlpiClearState() or set SCIP_LPPAR_FROMSCRATCH to TRUE */
3575 if( retval == CPXERR_NO_SOLN || retval == CPXERR_NO_LU_FACTOR || retval == CPXERR_NO_BASIC_SOLN || retval == CPXERR_NO_BASIS )
3589 const SCIP_Real* binvrow, /**< row in (A_B)^-1 from prior call to SCIPlpiGetBInvRow(), or NULL */
3601 /* this might be turned off if the user as called SCIPlpiClearState() or set SCIP_LPPAR_FROMSCRATCH to TRUE */
3606 if( retval == CPXERR_NO_SOLN || retval == CPXERR_NO_LU_FACTOR || retval == CPXERR_NO_BASIC_SOLN || retval == CPXERR_NO_BASIS )
3631 /* this might be turned off if the user as called SCIPlpiClearState() or set SCIP_LPPAR_FROMSCRATCH to TRUE */
3636 if( retval == CPXERR_NO_SOLN || retval == CPXERR_NO_LU_FACTOR || retval == CPXERR_NO_BASIC_SOLN || retval == CPXERR_NO_BASIS )
3674 /* if there is no basis information available (e.g. after barrier without crossover), or no state can be saved; if
3691 SCIPdebugMessage("storing CPLEX LPI state in %p (%d cols, %d rows)\n", (void *) *lpistate, ncols, nrows);
3704 /** loads LPi state (like basis information) into solver; note that the LP might have been extended with additional
3731 SCIPdebugMessage("loading LPI state %p (%d cols, %d rows) into CPLEX LP with %d cols and %d rows\n",
3876 /* if there is no basis information available (e.g. after barrier without crossover), norms cannot be saved; if
3894 SCIPdebugMessage("storing CPLEX LPI pricing norms in %p (%d rows)\n", (void *) *lpinorms, nrows);
3897 retval = CPXgetdnorms(lpi->cpxenv, lpi->cpxlp, (*lpinorms)->norm, (*lpinorms)->head, &((*lpinorms)->normlen));
3899 /* if CPLEX used the primal simplex in the last optimization call, we do not have dual norms (error 1264) */
3917 /** loads LPi pricing norms into solver; note that the LP might have been extended with additional
3947 CHECK_ZERO( lpi->messagehdlr, CPXcopydnorms(lpi->cpxenv, lpi->cpxlp, lpinorms->norm, lpinorms->head, lpinorms->normlen) );
4053 #if (CPX_VERSION == 1100 || (CPX_VERSION == 1220 && (CPX_SUBVERSION == 0 || CPX_SUBVERSION == 2)))
4054 /* Due to CPLEX bug, we always set the thread count to 1. In order to fulfill an assert in lp.c, we have to
4152 #if (CPX_VERSION == 1100 || (CPX_VERSION == 1220 && (CPX_SUBVERSION == 0 || CPX_SUBVERSION == 2)))
|