glpapi07.c 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452
  1. /* glpapi07.c (exact simplex solver) */
  2. /***********************************************************************
  3. * This code is part of GLPK (GNU Linear Programming Kit).
  4. *
  5. * Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008,
  6. * 2009, 2010 Andrew Makhorin, Department for Applied Informatics,
  7. * Moscow Aviation Institute, Moscow, Russia. All rights reserved.
  8. * E-mail: <mao@gnu.org>.
  9. *
  10. * GLPK is free software: you can redistribute it and/or modify it
  11. * under the terms of the GNU General Public License as published by
  12. * the Free Software Foundation, either version 3 of the License, or
  13. * (at your option) any later version.
  14. *
  15. * GLPK is distributed in the hope that it will be useful, but WITHOUT
  16. * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
  17. * or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public
  18. * License for more details.
  19. *
  20. * You should have received a copy of the GNU General Public License
  21. * along with GLPK. If not, see <http://www.gnu.org/licenses/>.
  22. ***********************************************************************/
  23. #include "glpapi.h"
  24. #include "glpssx.h"
  25. /***********************************************************************
  26. * NAME
  27. *
  28. * glp_exact - solve LP problem in exact arithmetic
  29. *
  30. * SYNOPSIS
  31. *
  32. * int glp_exact(glp_prob *lp, const glp_smcp *parm);
  33. *
  34. * DESCRIPTION
  35. *
  36. * The routine glp_exact is a tentative implementation of the primal
  37. * two-phase simplex method based on exact (rational) arithmetic. It is
  38. * similar to the routine glp_simplex, however, for all internal
  39. * computations it uses arithmetic of rational numbers, which is exact
  40. * in mathematical sense, i.e. free of round-off errors unlike floating
  41. * point arithmetic.
  42. *
  43. * Note that the routine glp_exact uses inly two control parameters
  44. * passed in the structure glp_smcp, namely, it_lim and tm_lim.
  45. *
  46. * RETURNS
  47. *
  48. * 0 The LP problem instance has been successfully solved. This code
  49. * does not necessarily mean that the solver has found optimal
  50. * solution. It only means that the solution process was successful.
  51. *
  52. * GLP_EBADB
  53. * Unable to start the search, because the initial basis specified
  54. * in the problem object is invalid--the number of basic (auxiliary
  55. * and structural) variables is not the same as the number of rows in
  56. * the problem object.
  57. *
  58. * GLP_ESING
  59. * Unable to start the search, because the basis matrix correspodning
  60. * to the initial basis is exactly singular.
  61. *
  62. * GLP_EBOUND
  63. * Unable to start the search, because some double-bounded variables
  64. * have incorrect bounds.
  65. *
  66. * GLP_EFAIL
  67. * The problem has no rows/columns.
  68. *
  69. * GLP_EITLIM
  70. * The search was prematurely terminated, because the simplex
  71. * iteration limit has been exceeded.
  72. *
  73. * GLP_ETMLIM
  74. * The search was prematurely terminated, because the time limit has
  75. * been exceeded. */
  76. static void set_d_eps(mpq_t x, double val)
  77. { /* convert double val to rational x obtaining a more adequate
  78. fraction than provided by mpq_set_d due to allowing a small
  79. approximation error specified by a given relative tolerance;
  80. for example, mpq_set_d would give the following
  81. 1/3 ~= 0.333333333333333314829616256247391... ->
  82. -> 6004799503160661/18014398509481984
  83. while this routine gives exactly 1/3 */
  84. int s, n, j;
  85. double f, p, q, eps = 1e-9;
  86. mpq_t temp;
  87. xassert(-DBL_MAX <= val && val <= +DBL_MAX);
  88. #if 1 /* 30/VII-2008 */
  89. if (val == floor(val))
  90. { /* if val is integral, do not approximate */
  91. mpq_set_d(x, val);
  92. goto done;
  93. }
  94. #endif
  95. if (val > 0.0)
  96. s = +1;
  97. else if (val < 0.0)
  98. s = -1;
  99. else
  100. { mpq_set_si(x, 0, 1);
  101. goto done;
  102. }
  103. f = frexp(fabs(val), &n);
  104. /* |val| = f * 2^n, where 0.5 <= f < 1.0 */
  105. fp2rat(f, 0.1 * eps, &p, &q);
  106. /* f ~= p / q, where p and q are integers */
  107. mpq_init(temp);
  108. mpq_set_d(x, p);
  109. mpq_set_d(temp, q);
  110. mpq_div(x, x, temp);
  111. mpq_set_si(temp, 1, 1);
  112. for (j = 1; j <= abs(n); j++)
  113. mpq_add(temp, temp, temp);
  114. if (n > 0)
  115. mpq_mul(x, x, temp);
  116. else if (n < 0)
  117. mpq_div(x, x, temp);
  118. mpq_clear(temp);
  119. if (s < 0) mpq_neg(x, x);
  120. /* check that the desired tolerance has been attained */
  121. xassert(fabs(val - mpq_get_d(x)) <= eps * (1.0 + fabs(val)));
  122. done: return;
  123. }
  124. static void load_data(SSX *ssx, LPX *lp)
  125. { /* load LP problem data into simplex solver workspace */
  126. int m = ssx->m;
  127. int n = ssx->n;
  128. int nnz = ssx->A_ptr[n+1]-1;
  129. int j, k, type, loc, len, *ind;
  130. double lb, ub, coef, *val;
  131. xassert(lpx_get_num_rows(lp) == m);
  132. xassert(lpx_get_num_cols(lp) == n);
  133. xassert(lpx_get_num_nz(lp) == nnz);
  134. /* types and bounds of rows and columns */
  135. for (k = 1; k <= m+n; k++)
  136. { if (k <= m)
  137. { type = lpx_get_row_type(lp, k);
  138. lb = lpx_get_row_lb(lp, k);
  139. ub = lpx_get_row_ub(lp, k);
  140. }
  141. else
  142. { type = lpx_get_col_type(lp, k-m);
  143. lb = lpx_get_col_lb(lp, k-m);
  144. ub = lpx_get_col_ub(lp, k-m);
  145. }
  146. switch (type)
  147. { case LPX_FR: type = SSX_FR; break;
  148. case LPX_LO: type = SSX_LO; break;
  149. case LPX_UP: type = SSX_UP; break;
  150. case LPX_DB: type = SSX_DB; break;
  151. case LPX_FX: type = SSX_FX; break;
  152. default: xassert(type != type);
  153. }
  154. ssx->type[k] = type;
  155. set_d_eps(ssx->lb[k], lb);
  156. set_d_eps(ssx->ub[k], ub);
  157. }
  158. /* optimization direction */
  159. switch (lpx_get_obj_dir(lp))
  160. { case LPX_MIN: ssx->dir = SSX_MIN; break;
  161. case LPX_MAX: ssx->dir = SSX_MAX; break;
  162. default: xassert(lp != lp);
  163. }
  164. /* objective coefficients */
  165. for (k = 0; k <= m+n; k++)
  166. { if (k == 0)
  167. coef = lpx_get_obj_coef(lp, 0);
  168. else if (k <= m)
  169. coef = 0.0;
  170. else
  171. coef = lpx_get_obj_coef(lp, k-m);
  172. set_d_eps(ssx->coef[k], coef);
  173. }
  174. /* constraint coefficients */
  175. ind = xcalloc(1+m, sizeof(int));
  176. val = xcalloc(1+m, sizeof(double));
  177. loc = 0;
  178. for (j = 1; j <= n; j++)
  179. { ssx->A_ptr[j] = loc+1;
  180. len = lpx_get_mat_col(lp, j, ind, val);
  181. for (k = 1; k <= len; k++)
  182. { loc++;
  183. ssx->A_ind[loc] = ind[k];
  184. set_d_eps(ssx->A_val[loc], val[k]);
  185. }
  186. }
  187. xassert(loc == nnz);
  188. xfree(ind);
  189. xfree(val);
  190. return;
  191. }
  192. static int load_basis(SSX *ssx, LPX *lp)
  193. { /* load current LP basis into simplex solver workspace */
  194. int m = ssx->m;
  195. int n = ssx->n;
  196. int *type = ssx->type;
  197. int *stat = ssx->stat;
  198. int *Q_row = ssx->Q_row;
  199. int *Q_col = ssx->Q_col;
  200. int i, j, k;
  201. xassert(lpx_get_num_rows(lp) == m);
  202. xassert(lpx_get_num_cols(lp) == n);
  203. /* statuses of rows and columns */
  204. for (k = 1; k <= m+n; k++)
  205. { if (k <= m)
  206. stat[k] = lpx_get_row_stat(lp, k);
  207. else
  208. stat[k] = lpx_get_col_stat(lp, k-m);
  209. switch (stat[k])
  210. { case LPX_BS:
  211. stat[k] = SSX_BS;
  212. break;
  213. case LPX_NL:
  214. stat[k] = SSX_NL;
  215. xassert(type[k] == SSX_LO || type[k] == SSX_DB);
  216. break;
  217. case LPX_NU:
  218. stat[k] = SSX_NU;
  219. xassert(type[k] == SSX_UP || type[k] == SSX_DB);
  220. break;
  221. case LPX_NF:
  222. stat[k] = SSX_NF;
  223. xassert(type[k] == SSX_FR);
  224. break;
  225. case LPX_NS:
  226. stat[k] = SSX_NS;
  227. xassert(type[k] == SSX_FX);
  228. break;
  229. default:
  230. xassert(stat != stat);
  231. }
  232. }
  233. /* build permutation matix Q */
  234. i = j = 0;
  235. for (k = 1; k <= m+n; k++)
  236. { if (stat[k] == SSX_BS)
  237. { i++;
  238. if (i > m) return 1;
  239. Q_row[k] = i, Q_col[i] = k;
  240. }
  241. else
  242. { j++;
  243. if (j > n) return 1;
  244. Q_row[k] = m+j, Q_col[m+j] = k;
  245. }
  246. }
  247. xassert(i == m && j == n);
  248. return 0;
  249. }
  250. int glp_exact(glp_prob *lp, const glp_smcp *parm)
  251. { glp_smcp _parm;
  252. SSX *ssx;
  253. int m = lpx_get_num_rows(lp);
  254. int n = lpx_get_num_cols(lp);
  255. int nnz = lpx_get_num_nz(lp);
  256. int i, j, k, type, pst, dst, ret, *stat;
  257. double lb, ub, *prim, *dual, sum;
  258. if (parm == NULL)
  259. parm = &_parm, glp_init_smcp((glp_smcp *)parm);
  260. /* check control parameters */
  261. if (parm->it_lim < 0)
  262. xerror("glp_exact: it_lim = %d; invalid parameter\n",
  263. parm->it_lim);
  264. if (parm->tm_lim < 0)
  265. xerror("glp_exact: tm_lim = %d; invalid parameter\n",
  266. parm->tm_lim);
  267. /* the problem must have at least one row and one column */
  268. if (!(m > 0 && n > 0))
  269. { xprintf("glp_exact: problem has no rows/columns\n");
  270. return GLP_EFAIL;
  271. }
  272. #if 1
  273. /* basic solution is currently undefined */
  274. lp->pbs_stat = lp->dbs_stat = GLP_UNDEF;
  275. lp->obj_val = 0.0;
  276. lp->some = 0;
  277. #endif
  278. /* check that all double-bounded variables have correct bounds */
  279. for (k = 1; k <= m+n; k++)
  280. { if (k <= m)
  281. { type = lpx_get_row_type(lp, k);
  282. lb = lpx_get_row_lb(lp, k);
  283. ub = lpx_get_row_ub(lp, k);
  284. }
  285. else
  286. { type = lpx_get_col_type(lp, k-m);
  287. lb = lpx_get_col_lb(lp, k-m);
  288. ub = lpx_get_col_ub(lp, k-m);
  289. }
  290. if (type == LPX_DB && lb >= ub)
  291. { xprintf("glp_exact: %s %d has invalid bounds\n",
  292. k <= m ? "row" : "column", k <= m ? k : k-m);
  293. return GLP_EBOUND;
  294. }
  295. }
  296. /* create the simplex solver workspace */
  297. xprintf("glp_exact: %d rows, %d columns, %d non-zeros\n",
  298. m, n, nnz);
  299. #ifdef HAVE_GMP
  300. xprintf("GNU MP bignum library is being used\n");
  301. #else
  302. xprintf("GLPK bignum module is being used\n");
  303. xprintf("(Consider installing GNU MP to attain a much better perf"
  304. "ormance.)\n");
  305. #endif
  306. ssx = ssx_create(m, n, nnz);
  307. /* load LP problem data into the workspace */
  308. load_data(ssx, lp);
  309. /* load current LP basis into the workspace */
  310. if (load_basis(ssx, lp))
  311. { xprintf("glp_exact: initial LP basis is invalid\n");
  312. ret = GLP_EBADB;
  313. goto done;
  314. }
  315. /* inherit some control parameters from the LP object */
  316. #if 0
  317. ssx->it_lim = lpx_get_int_parm(lp, LPX_K_ITLIM);
  318. ssx->it_cnt = lpx_get_int_parm(lp, LPX_K_ITCNT);
  319. ssx->tm_lim = lpx_get_real_parm(lp, LPX_K_TMLIM);
  320. #else
  321. ssx->it_lim = parm->it_lim;
  322. ssx->it_cnt = lp->it_cnt;
  323. ssx->tm_lim = (double)parm->tm_lim / 1000.0;
  324. #endif
  325. ssx->out_frq = 5.0;
  326. ssx->tm_beg = xtime();
  327. ssx->tm_lag = xlset(0);
  328. /* solve LP */
  329. ret = ssx_driver(ssx);
  330. /* copy back some statistics to the LP object */
  331. #if 0
  332. lpx_set_int_parm(lp, LPX_K_ITLIM, ssx->it_lim);
  333. lpx_set_int_parm(lp, LPX_K_ITCNT, ssx->it_cnt);
  334. lpx_set_real_parm(lp, LPX_K_TMLIM, ssx->tm_lim);
  335. #else
  336. lp->it_cnt = ssx->it_cnt;
  337. #endif
  338. /* analyze the return code */
  339. switch (ret)
  340. { case 0:
  341. /* optimal solution found */
  342. ret = 0;
  343. pst = LPX_P_FEAS, dst = LPX_D_FEAS;
  344. break;
  345. case 1:
  346. /* problem has no feasible solution */
  347. ret = 0;
  348. pst = LPX_P_NOFEAS, dst = LPX_D_INFEAS;
  349. break;
  350. case 2:
  351. /* problem has unbounded solution */
  352. ret = 0;
  353. pst = LPX_P_FEAS, dst = LPX_D_NOFEAS;
  354. #if 1
  355. xassert(1 <= ssx->q && ssx->q <= n);
  356. lp->some = ssx->Q_col[m + ssx->q];
  357. xassert(1 <= lp->some && lp->some <= m+n);
  358. #endif
  359. break;
  360. case 3:
  361. /* iteration limit exceeded (phase I) */
  362. ret = GLP_EITLIM;
  363. pst = LPX_P_INFEAS, dst = LPX_D_INFEAS;
  364. break;
  365. case 4:
  366. /* iteration limit exceeded (phase II) */
  367. ret = GLP_EITLIM;
  368. pst = LPX_P_FEAS, dst = LPX_D_INFEAS;
  369. break;
  370. case 5:
  371. /* time limit exceeded (phase I) */
  372. ret = GLP_ETMLIM;
  373. pst = LPX_P_INFEAS, dst = LPX_D_INFEAS;
  374. break;
  375. case 6:
  376. /* time limit exceeded (phase II) */
  377. ret = GLP_ETMLIM;
  378. pst = LPX_P_FEAS, dst = LPX_D_INFEAS;
  379. break;
  380. case 7:
  381. /* initial basis matrix is singular */
  382. ret = GLP_ESING;
  383. goto done;
  384. default:
  385. xassert(ret != ret);
  386. }
  387. /* obtain final basic solution components */
  388. stat = xcalloc(1+m+n, sizeof(int));
  389. prim = xcalloc(1+m+n, sizeof(double));
  390. dual = xcalloc(1+m+n, sizeof(double));
  391. for (k = 1; k <= m+n; k++)
  392. { if (ssx->stat[k] == SSX_BS)
  393. { i = ssx->Q_row[k]; /* x[k] = xB[i] */
  394. xassert(1 <= i && i <= m);
  395. stat[k] = LPX_BS;
  396. prim[k] = mpq_get_d(ssx->bbar[i]);
  397. dual[k] = 0.0;
  398. }
  399. else
  400. { j = ssx->Q_row[k] - m; /* x[k] = xN[j] */
  401. xassert(1 <= j && j <= n);
  402. switch (ssx->stat[k])
  403. { case SSX_NF:
  404. stat[k] = LPX_NF;
  405. prim[k] = 0.0;
  406. break;
  407. case SSX_NL:
  408. stat[k] = LPX_NL;
  409. prim[k] = mpq_get_d(ssx->lb[k]);
  410. break;
  411. case SSX_NU:
  412. stat[k] = LPX_NU;
  413. prim[k] = mpq_get_d(ssx->ub[k]);
  414. break;
  415. case SSX_NS:
  416. stat[k] = LPX_NS;
  417. prim[k] = mpq_get_d(ssx->lb[k]);
  418. break;
  419. default:
  420. xassert(ssx != ssx);
  421. }
  422. dual[k] = mpq_get_d(ssx->cbar[j]);
  423. }
  424. }
  425. /* and store them into the LP object */
  426. pst = pst - LPX_P_UNDEF + GLP_UNDEF;
  427. dst = dst - LPX_D_UNDEF + GLP_UNDEF;
  428. for (k = 1; k <= m+n; k++)
  429. stat[k] = stat[k] - LPX_BS + GLP_BS;
  430. sum = lpx_get_obj_coef(lp, 0);
  431. for (j = 1; j <= n; j++)
  432. sum += lpx_get_obj_coef(lp, j) * prim[m+j];
  433. lpx_put_solution(lp, 1, &pst, &dst, &sum,
  434. &stat[0], &prim[0], &dual[0], &stat[m], &prim[m], &dual[m]);
  435. xfree(stat);
  436. xfree(prim);
  437. xfree(dual);
  438. done: /* delete the simplex solver workspace */
  439. ssx_delete(ssx);
  440. /* return to the application program */
  441. return ret;
  442. }
  443. /* eof */