5234
|
1 /* |
|
2 |
|
3 Copyright (C) 2005 Nicolo' Giorgetti |
|
4 |
|
5 This file is part of Octave. |
|
6 |
|
7 Octave is free software; you can redistribute it and/or modify it |
|
8 under the terms of the GNU General Public License as published by the |
|
9 Free Software Foundation; either version 2, or (at your option) any |
|
10 later version. |
|
11 |
|
12 Octave is distributed in the hope that it will be useful, but WITHOUT |
|
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
|
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
|
15 for more details. |
|
16 |
|
17 You should have received a copy of the GNU General Public License |
|
18 along with Octave; see the file COPYING. If not, write to the Free |
5307
|
19 Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA |
|
20 02110-1301, USA. |
5234
|
21 |
|
22 */ |
|
23 |
|
24 #ifdef HAVE_CONFIG_H |
|
25 #include <config.h> |
|
26 #endif |
5232
|
27 |
|
28 #include <cfloat> |
|
29 #include <csetjmp> |
|
30 #include <ctime> |
|
31 |
5234
|
32 #include "defun-dld.h" |
|
33 #include "error.h" |
5235
|
34 #include "gripes.h" |
|
35 #include "oct-map.h" |
5234
|
36 #include "oct-obj.h" |
5235
|
37 #include "pager.h" |
|
38 |
|
39 #if defined (HAVE_GLPK) |
5232
|
40 |
6333
|
41 extern "C" |
|
42 { |
5234
|
43 #include <glpk.h> |
6333
|
44 |
6381
|
45 #ifdef GLPK_PRE_4_14 |
6333
|
46 |
6287
|
47 #ifndef _GLPLIB_H |
|
48 #include <glplib.h> |
|
49 #endif |
|
50 #ifndef lib_set_fault_hook |
|
51 #define lib_set_fault_hook lib_fault_hook |
|
52 #endif |
|
53 #ifndef lib_set_print_hook |
|
54 #define lib_set_print_hook lib_print_hook |
|
55 #endif |
5232
|
56 |
6333
|
57 #else |
|
58 |
|
59 void _glp_lib_print_hook (int (*func)(void *info, char *buf), void *info); |
|
60 void _glp_lib_fault_hook (int (*func)(void *info, char *buf), void *info); |
|
61 |
|
62 #endif |
6472
|
63 } |
6333
|
64 |
5232
|
65 #define NIntP 17 |
|
66 #define NRealP 10 |
|
67 |
5234
|
68 int lpxIntParam[NIntP] = { |
5240
|
69 0, |
5234
|
70 1, |
|
71 0, |
|
72 1, |
|
73 0, |
|
74 -1, |
|
75 0, |
|
76 200, |
|
77 1, |
|
78 2, |
|
79 0, |
|
80 1, |
|
81 0, |
|
82 0, |
|
83 2, |
|
84 2, |
|
85 1 |
5232
|
86 }; |
|
87 |
5234
|
88 int IParam[NIntP] = { |
|
89 LPX_K_MSGLEV, |
|
90 LPX_K_SCALE, |
|
91 LPX_K_DUAL, |
|
92 LPX_K_PRICE, |
|
93 LPX_K_ROUND, |
|
94 LPX_K_ITLIM, |
|
95 LPX_K_ITCNT, |
|
96 LPX_K_OUTFRQ, |
|
97 LPX_K_MPSINFO, |
|
98 LPX_K_MPSOBJ, |
|
99 LPX_K_MPSORIG, |
|
100 LPX_K_MPSWIDE, |
|
101 LPX_K_MPSFREE, |
|
102 LPX_K_MPSSKIP, |
|
103 LPX_K_BRANCH, |
|
104 LPX_K_BTRACK, |
|
105 LPX_K_PRESOL |
5232
|
106 }; |
|
107 |
|
108 |
5234
|
109 double lpxRealParam[NRealP] = { |
|
110 0.07, |
|
111 1e-7, |
|
112 1e-7, |
|
113 1e-9, |
|
114 -DBL_MAX, |
|
115 DBL_MAX, |
|
116 -1.0, |
|
117 0.0, |
|
118 1e-6, |
|
119 1e-7 |
5232
|
120 }; |
|
121 |
5234
|
122 int RParam[NRealP] = { |
|
123 LPX_K_RELAX, |
|
124 LPX_K_TOLBND, |
|
125 LPX_K_TOLDJ, |
|
126 LPX_K_TOLPIV, |
|
127 LPX_K_OBJLL, |
|
128 LPX_K_OBJUL, |
|
129 LPX_K_TMLIM, |
|
130 LPX_K_OUTDLY, |
|
131 LPX_K_TOLINT, |
|
132 LPX_K_TOLOBJ |
5232
|
133 }; |
|
134 |
5241
|
135 static jmp_buf mark; //-- Address for long jump to jump to |
5232
|
136 |
5234
|
137 int |
5235
|
138 glpk_fault_hook (void * /* info */, char *msg) |
5234
|
139 { |
5240
|
140 error ("CRITICAL ERROR in GLPK: %s", msg); |
5234
|
141 longjmp (mark, -1); |
5232
|
142 } |
|
143 |
5234
|
144 int |
5235
|
145 glpk_print_hook (void * /* info */, char *msg) |
5232
|
146 { |
5241
|
147 message (0, "%s", msg); |
5234
|
148 return 1; |
5232
|
149 } |
|
150 |
5234
|
151 int |
|
152 glpk (int sense, int n, int m, double *c, int nz, int *rn, int *cn, |
|
153 double *a, double *b, char *ctype, int *freeLB, double *lb, |
|
154 int *freeUB, double *ub, int *vartype, int isMIP, int lpsolver, |
|
155 int save_pb, double *xmin, double *fmin, double *status, |
|
156 double *lambda, double *redcosts, double *time, double *mem) |
5232
|
157 { |
5240
|
158 int errnum; |
5234
|
159 int typx = 0; |
|
160 int method; |
|
161 |
|
162 clock_t t_start = clock(); |
|
163 |
6381
|
164 #ifdef GLPK_PRE_4_14 |
6333
|
165 lib_set_fault_hook (0, glpk_fault_hook); |
|
166 #else |
|
167 _glp_lib_fault_hook (glpk_fault_hook, 0); |
|
168 #endif |
5232
|
169 |
5234
|
170 if (lpxIntParam[0] > 1) |
6381
|
171 #ifdef GLPK_PRE_4_14 |
6333
|
172 lib_set_print_hook (0, glpk_print_hook); |
|
173 #else |
|
174 _glp_lib_print_hook (glpk_print_hook, 0); |
|
175 #endif |
5234
|
176 |
|
177 LPX *lp = lpx_create_prob (); |
|
178 |
5232
|
179 |
5234
|
180 //-- Set the sense of optimization |
|
181 if (sense == 1) |
|
182 lpx_set_obj_dir (lp, LPX_MIN); |
|
183 else |
|
184 lpx_set_obj_dir (lp, LPX_MAX); |
|
185 |
|
186 //-- If the problem has integer structural variables switch to MIP |
|
187 if (isMIP) |
|
188 lpx_set_class (lp, LPX_MIP); |
5232
|
189 |
5234
|
190 lpx_add_cols (lp, n); |
|
191 for (int i = 0; i < n; i++) |
|
192 { |
5232
|
193 //-- Define type of the structural variables |
5234
|
194 if (! freeLB[i] && ! freeUB[i]) |
|
195 lpx_set_col_bnds (lp, i+1, LPX_DB, lb[i], ub[i]); |
|
196 else |
|
197 { |
|
198 if (! freeLB[i] && freeUB[i]) |
|
199 lpx_set_col_bnds (lp, i+1, LPX_LO, lb[i], ub[i]); |
|
200 else |
|
201 { |
|
202 if (freeLB[i] && ! freeUB[i]) |
|
203 lpx_set_col_bnds (lp, i+1, LPX_UP, lb[i], ub[i]); |
|
204 else |
|
205 lpx_set_col_bnds (lp, i+1, LPX_FR, lb[i], ub[i]); |
|
206 } |
|
207 } |
|
208 |
|
209 // -- Set the objective coefficient of the corresponding |
5232
|
210 // -- structural variable. No constant term is assumed. |
5234
|
211 lpx_set_obj_coef(lp,i+1,c[i]); |
5232
|
212 |
5234
|
213 if (isMIP) |
|
214 lpx_set_col_kind (lp, i+1, vartype[i]); |
|
215 } |
|
216 |
|
217 lpx_add_rows (lp, m); |
|
218 |
|
219 for (int i = 0; i < m; i++) |
|
220 { |
5232
|
221 /* If the i-th row has no lower bound (types F,U), the |
|
222 corrispondent parameter will be ignored. |
|
223 If the i-th row has no upper bound (types F,L), the corrispondent |
|
224 parameter will be ignored. |
|
225 If the i-th row is of S type, the i-th LB is used, but |
|
226 the i-th UB is ignored. |
|
227 */ |
5234
|
228 |
|
229 switch (ctype[i]) |
|
230 { |
|
231 case 'F': |
|
232 typx = LPX_FR; |
|
233 break; |
|
234 |
|
235 case 'U': |
|
236 typx = LPX_UP; |
|
237 break; |
|
238 |
|
239 case 'L': |
|
240 typx = LPX_LO; |
|
241 break; |
5232
|
242 |
5234
|
243 case 'S': |
|
244 typx = LPX_FX; |
|
245 break; |
5232
|
246 |
5234
|
247 case 'D': |
|
248 typx = LPX_DB; |
|
249 break; |
|
250 } |
|
251 |
|
252 lpx_set_row_bnds (lp, i+1, typx, b[i], b[i]); |
|
253 |
|
254 } |
|
255 |
|
256 lpx_load_matrix (lp, nz, rn, cn, a); |
5232
|
257 |
5234
|
258 if (save_pb) |
|
259 { |
6484
|
260 static char tmp[] = "outpb.lp"; |
|
261 if (lpx_write_cpxlp (lp, tmp) != 0) |
5234
|
262 { |
5240
|
263 error ("__glpk__: unable to write problem"); |
5234
|
264 longjmp (mark, -1); |
|
265 } |
|
266 } |
|
267 |
|
268 //-- scale the problem data (if required) |
|
269 //-- if (scale && (!presol || method == 1)) lpx_scale_prob(lp); |
|
270 //-- LPX_K_SCALE=IParam[1] LPX_K_PRESOL=IParam[16] |
|
271 if (lpxIntParam[1] && (! lpxIntParam[16] || lpsolver != 1)) |
|
272 lpx_scale_prob (lp); |
|
273 |
|
274 //-- build advanced initial basis (if required) |
|
275 if (lpsolver == 1 && ! lpxIntParam[16]) |
|
276 lpx_adv_basis (lp); |
|
277 |
|
278 for(int i = 0; i < NIntP; i++) |
|
279 lpx_set_int_parm (lp, IParam[i], lpxIntParam[i]); |
5232
|
280 |
5234
|
281 for (int i = 0; i < NRealP; i++) |
|
282 lpx_set_real_parm (lp, RParam[i], lpxRealParam[i]); |
|
283 |
|
284 if (lpsolver == 1) |
|
285 method = 'S'; |
|
286 else |
|
287 method = 'T'; |
5232
|
288 |
5234
|
289 switch (method) |
|
290 { |
|
291 case 'S': |
|
292 { |
|
293 if (isMIP) |
|
294 { |
|
295 method = 'I'; |
5240
|
296 errnum = lpx_simplex (lp); |
|
297 errnum = lpx_integer (lp); |
5234
|
298 } |
|
299 else |
5240
|
300 errnum = lpx_simplex(lp); |
5234
|
301 } |
5232
|
302 break; |
5234
|
303 |
|
304 case 'T': |
5240
|
305 errnum = lpx_interior(lp); |
5234
|
306 break; |
|
307 |
|
308 default: |
6381
|
309 #ifdef GLPK_PRE_4_14 |
5234
|
310 insist (method != method); |
6333
|
311 #else |
6484
|
312 static char tmp[] = "method != method"; |
|
313 glpk_fault_hook (0, tmp); |
6333
|
314 #endif |
5234
|
315 } |
|
316 |
5240
|
317 /* errnum assumes the following results: |
|
318 errnum = 0 <=> No errors |
|
319 errnum = 1 <=> Iteration limit exceeded. |
|
320 errnum = 2 <=> Numerical problems with basis matrix. |
5234
|
321 */ |
5240
|
322 if (errnum == LPX_E_OK) |
5234
|
323 { |
|
324 if (isMIP) |
|
325 { |
5241
|
326 *status = lpx_mip_status (lp); |
5234
|
327 *fmin = lpx_mip_obj_val (lp); |
|
328 } |
|
329 else |
|
330 { |
|
331 if (lpsolver == 1) |
|
332 { |
5241
|
333 *status = lpx_get_status (lp); |
5234
|
334 *fmin = lpx_get_obj_val (lp); |
|
335 } |
|
336 else |
|
337 { |
5241
|
338 *status = lpx_ipt_status (lp); |
5234
|
339 *fmin = lpx_ipt_obj_val (lp); |
|
340 } |
|
341 } |
5232
|
342 |
5234
|
343 if (isMIP) |
|
344 { |
|
345 for (int i = 0; i < n; i++) |
|
346 xmin[i] = lpx_mip_col_val (lp, i+1); |
|
347 } |
|
348 else |
|
349 { |
|
350 /* Primal values */ |
|
351 for (int i = 0; i < n; i++) |
|
352 { |
|
353 if (lpsolver == 1) |
|
354 xmin[i] = lpx_get_col_prim (lp, i+1); |
|
355 else |
|
356 xmin[i] = lpx_ipt_col_prim (lp, i+1); |
|
357 } |
5232
|
358 |
5234
|
359 /* Dual values */ |
|
360 for (int i = 0; i < m; i++) |
|
361 { |
|
362 if (lpsolver == 1) |
|
363 lambda[i] = lpx_get_row_dual (lp, i+1); |
|
364 else |
|
365 lambda[i] = lpx_ipt_row_dual (lp, i+1); |
|
366 } |
|
367 |
|
368 /* Reduced costs */ |
|
369 for (int i = 0; i < lpx_get_num_cols (lp); i++) |
|
370 { |
|
371 if (lpsolver == 1) |
|
372 redcosts[i] = lpx_get_col_dual (lp, i+1); |
|
373 else |
|
374 redcosts[i] = lpx_ipt_col_dual (lp, i+1); |
|
375 } |
|
376 } |
|
377 |
5241
|
378 *time = (clock () - t_start) / CLOCKS_PER_SEC; |
6333
|
379 |
6381
|
380 #ifdef GLPK_PRE_4_14 |
5241
|
381 *mem = (lib_env_ptr () -> mem_tpeak); |
6333
|
382 #else |
|
383 *mem = 0; |
|
384 #endif |
5234
|
385 |
|
386 lpx_delete_prob (lp); |
|
387 return 0; |
|
388 } |
|
389 |
|
390 lpx_delete_prob (lp); |
|
391 |
5241
|
392 *status = errnum; |
5234
|
393 |
5240
|
394 return errnum; |
5232
|
395 } |
|
396 |
5235
|
397 #endif |
|
398 |
5240
|
399 #define OCTAVE_GLPK_GET_REAL_PARAM(NAME, IDX) \ |
|
400 do \ |
|
401 { \ |
|
402 if (PARAM.contains (NAME)) \ |
|
403 { \ |
|
404 Cell tmp = PARAM.contents (NAME); \ |
|
405 \ |
|
406 if (! tmp.is_empty ()) \ |
|
407 { \ |
|
408 lpxRealParam[IDX] = tmp(0).scalar_value (); \ |
|
409 \ |
|
410 if (error_state) \ |
|
411 { \ |
|
412 error ("glpk: invalid value in param." NAME); \ |
|
413 return retval; \ |
|
414 } \ |
|
415 } \ |
|
416 else \ |
|
417 { \ |
|
418 error ("glpk: invalid value in param." NAME); \ |
|
419 return retval; \ |
|
420 } \ |
|
421 } \ |
|
422 } \ |
|
423 while (0) |
|
424 |
|
425 #define OCTAVE_GLPK_GET_INT_PARAM(NAME, VAL) \ |
|
426 do \ |
|
427 { \ |
|
428 if (PARAM.contains (NAME)) \ |
|
429 { \ |
|
430 Cell tmp = PARAM.contents (NAME); \ |
|
431 \ |
|
432 if (! tmp.is_empty ()) \ |
|
433 { \ |
|
434 VAL = tmp(0).int_value (); \ |
|
435 \ |
|
436 if (error_state) \ |
|
437 { \ |
|
438 error ("glpk: invalid value in param." NAME); \ |
|
439 return retval; \ |
|
440 } \ |
|
441 } \ |
|
442 else \ |
|
443 { \ |
|
444 error ("glpk: invalid value in param." NAME); \ |
|
445 return retval; \ |
|
446 } \ |
|
447 } \ |
|
448 } \ |
|
449 while (0) |
|
450 |
5235
|
451 DEFUN_DLD (__glpk__, args, , |
5245
|
452 "-*- texinfo -*-\n\ |
|
453 @deftypefn {Loadable Function} {[@var{values}] =} __glpk__ (@var{args})\n\ |
|
454 Internal interface for the GNU GLPK library.\n\ |
|
455 You should be using using the @code{glpk} function instead.\n\ |
|
456 @end deftypefn") |
5232
|
457 { |
5234
|
458 // The list of values to return. See the declaration in oct-obj.h |
|
459 octave_value_list retval; |
5232
|
460 |
5235
|
461 #if defined (HAVE_GLPK) |
|
462 |
5234
|
463 int nrhs = args.length (); |
|
464 |
5240
|
465 if (nrhs != 9) |
5234
|
466 { |
5823
|
467 print_usage (); |
5232
|
468 return retval; |
5234
|
469 } |
|
470 |
5237
|
471 //-- 1nd Input. A column array containing the objective function |
|
472 //-- coefficients. |
5880
|
473 volatile int mrowsc = args(0).rows(); |
5234
|
474 |
5237
|
475 Matrix C (args(0).matrix_value ()); |
5240
|
476 |
|
477 if (error_state) |
|
478 { |
|
479 error ("__glpk__: invalid value of C"); |
|
480 return retval; |
|
481 } |
|
482 |
5234
|
483 double *c = C.fortran_vec (); |
5455
|
484 Array<int> rn; |
|
485 Array<int> cn; |
|
486 ColumnVector a; |
5515
|
487 volatile int mrowsA; |
5455
|
488 volatile int nz = 0; |
5232
|
489 |
5237
|
490 //-- 2nd Input. A matrix containing the constraints coefficients. |
5234
|
491 // If matrix A is NOT a sparse matrix |
5631
|
492 if (args(1).is_sparse_type ()) |
5455
|
493 { |
5603
|
494 SparseMatrix A = args(1).sparse_matrix_value (); // get the sparse matrix |
5455
|
495 |
|
496 if (error_state) |
|
497 { |
|
498 error ("__glpk__: invalid value of A"); |
|
499 return retval; |
|
500 } |
5234
|
501 |
5455
|
502 mrowsA = A.rows (); |
|
503 octave_idx_type Anc = A.cols (); |
5604
|
504 octave_idx_type Anz = A.nzmax (); |
5455
|
505 rn.resize (Anz+1); |
|
506 cn.resize (Anz+1); |
|
507 a.resize (Anz+1, 0.0); |
|
508 |
|
509 if (Anc != mrowsc) |
|
510 { |
|
511 error ("__glpk__: invalid value of A"); |
|
512 return retval; |
|
513 } |
|
514 |
|
515 for (octave_idx_type j = 0; j < Anc; j++) |
|
516 for (octave_idx_type i = A.cidx(j); i < A.cidx(j+1); i++) |
|
517 { |
|
518 nz++; |
|
519 rn(nz) = A.ridx(i) + 1; |
|
520 cn(nz) = j + 1; |
5603
|
521 a(nz) = A.data(i); |
5455
|
522 } |
|
523 } |
5631
|
524 else |
|
525 { |
|
526 Matrix A (args(1).matrix_value ()); // get the matrix |
|
527 |
|
528 if (error_state) |
|
529 { |
|
530 error ("__glpk__: invalid value of A"); |
|
531 return retval; |
|
532 } |
|
533 |
|
534 mrowsA = A.rows (); |
|
535 rn.resize (mrowsA*mrowsc+1); |
|
536 cn.resize (mrowsA*mrowsc+1); |
|
537 a.resize (mrowsA*mrowsc+1, 0.0); |
|
538 |
|
539 for (int i = 0; i < mrowsA; i++) |
|
540 { |
|
541 for (int j = 0; j < mrowsc; j++) |
|
542 { |
|
543 if (A(i,j) != 0) |
|
544 { |
|
545 nz++; |
|
546 rn(nz) = i + 1; |
|
547 cn(nz) = j + 1; |
|
548 a(nz) = A(i,j); |
|
549 } |
|
550 } |
|
551 } |
|
552 |
|
553 } |
5232
|
554 |
5237
|
555 //-- 3rd Input. A column array containing the right-hand side value |
5234
|
556 // for each constraint in the constraint matrix. |
5237
|
557 Matrix B (args(2).matrix_value ()); |
5240
|
558 |
|
559 if (error_state) |
|
560 { |
|
561 error ("__glpk__: invalid value of b"); |
|
562 return retval; |
|
563 } |
|
564 |
5234
|
565 double *b = B.fortran_vec (); |
|
566 |
5237
|
567 //-- 4th Input. An array of length mrowsc containing the lower |
5234
|
568 //-- bound on each of the variables. |
5237
|
569 Matrix LB (args(3).matrix_value ()); |
5240
|
570 |
|
571 if (error_state) |
|
572 { |
|
573 error ("__glpk__: invalid value of lb"); |
|
574 return retval; |
|
575 } |
|
576 |
5234
|
577 double *lb = LB.fortran_vec (); |
|
578 |
|
579 //-- LB argument, default: Free |
|
580 Array<int> freeLB (mrowsc); |
|
581 for (int i = 0; i < mrowsc; i++) |
|
582 { |
5455
|
583 if (xisinf (lb[i])) |
5234
|
584 { |
|
585 freeLB(i) = 1; |
|
586 lb[i] = -octave_Inf; |
|
587 } |
|
588 else |
|
589 freeLB(i) = 0; |
|
590 } |
|
591 |
5237
|
592 //-- 5th Input. An array of at least length numcols containing the upper |
5234
|
593 //-- bound on each of the variables. |
5237
|
594 Matrix UB (args(4).matrix_value ()); |
5234
|
595 |
5240
|
596 if (error_state) |
|
597 { |
|
598 error ("__glpk__: invalid value of ub"); |
|
599 return retval; |
|
600 } |
|
601 |
5234
|
602 double *ub = UB.fortran_vec (); |
5232
|
603 |
5234
|
604 Array<int> freeUB (mrowsc); |
|
605 for (int i = 0; i < mrowsc; i++) |
|
606 { |
5455
|
607 if (xisinf (ub[i])) |
5234
|
608 { |
|
609 freeUB(i) = 1; |
|
610 ub[i] = octave_Inf; |
|
611 } |
|
612 else |
|
613 freeUB(i) = 0; |
|
614 } |
|
615 |
5237
|
616 //-- 6th Input. A column array containing the sense of each constraint |
|
617 //-- in the constraint matrix. |
|
618 charMatrix CTYPE (args(5).char_matrix_value ()); |
5240
|
619 |
|
620 if (error_state) |
|
621 { |
|
622 error ("__glpk__: invalid value of ctype"); |
|
623 return retval; |
|
624 } |
|
625 |
5237
|
626 char *ctype = CTYPE.fortran_vec (); |
|
627 |
|
628 //-- 7th Input. A column array containing the types of the variables. |
|
629 charMatrix VTYPE (args(6).char_matrix_value ()); |
5232
|
630 |
5240
|
631 if (error_state) |
|
632 { |
|
633 error ("__glpk__: invalid value of vtype"); |
|
634 return retval; |
|
635 } |
|
636 |
5234
|
637 Array<int> vartype (mrowsc); |
5235
|
638 volatile int isMIP = 0; |
5234
|
639 for (int i = 0; i < mrowsc ; i++) |
|
640 { |
|
641 if (VTYPE(i,0) == 'I') |
|
642 { |
|
643 isMIP = 1; |
|
644 vartype(i) = LPX_IV; |
|
645 } |
|
646 else |
5235
|
647 vartype(i) = LPX_CV; |
5234
|
648 } |
|
649 |
5237
|
650 //-- 8th Input. Sense of optimization. |
|
651 volatile int sense; |
|
652 double SENSE = args(7).scalar_value (); |
5240
|
653 |
|
654 if (error_state) |
|
655 { |
|
656 error ("__glpk__: invalid value of sense"); |
|
657 return retval; |
|
658 } |
|
659 |
5237
|
660 if (SENSE >= 0) |
|
661 sense = 1; |
|
662 else |
|
663 sense = -1; |
|
664 |
5234
|
665 //-- 9th Input. A structure containing the control parameters. |
|
666 Octave_map PARAM = args(8).map_value (); |
|
667 |
5240
|
668 if (error_state) |
|
669 { |
|
670 error ("__glpk__: invalid value of param"); |
|
671 return retval; |
|
672 } |
|
673 |
5234
|
674 //-- ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|
675 //-- Integer parameters |
|
676 //-- ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|
677 |
|
678 //-- Level of messages output by the solver |
5240
|
679 OCTAVE_GLPK_GET_INT_PARAM ("msglev", lpxIntParam[0]); |
|
680 if (lpxIntParam[0] < 0 || lpxIntParam[0] > 3) |
5234
|
681 { |
5240
|
682 error ("__glpk__: param.msglev must be 0 (no output [default]) or 1 (error messages only) or 2 (normal output) or 3 (full output)"); |
|
683 return retval; |
|
684 } |
5232
|
685 |
5234
|
686 //-- scaling option |
5240
|
687 OCTAVE_GLPK_GET_INT_PARAM ("scale", lpxIntParam[1]); |
|
688 if (lpxIntParam[1] < 0 || lpxIntParam[1] > 2) |
5234
|
689 { |
5240
|
690 error ("__glpk__: param.scale must be 0 (no scaling) or 1 (equilibration scaling [default]) or 2 (geometric mean scaling)"); |
|
691 return retval; |
5234
|
692 } |
5232
|
693 |
5234
|
694 //-- Dual dimplex option |
5240
|
695 OCTAVE_GLPK_GET_INT_PARAM ("dual", lpxIntParam[2]); |
|
696 if (lpxIntParam[2] < 0 || lpxIntParam[2] > 1) |
5234
|
697 { |
5240
|
698 error ("__glpk__: param.dual must be 0 (do NOT use dual simplex [default]) or 1 (use dual simplex)"); |
|
699 return retval; |
5234
|
700 } |
5232
|
701 |
5234
|
702 //-- Pricing option |
5240
|
703 OCTAVE_GLPK_GET_INT_PARAM ("price", lpxIntParam[3]); |
|
704 if (lpxIntParam[3] < 0 || lpxIntParam[3] > 1) |
5234
|
705 { |
5240
|
706 error ("__glpk__: param.price must be 0 (textbook pricing) or 1 (steepest edge pricing [default])"); |
|
707 return retval; |
|
708 } |
5232
|
709 |
5234
|
710 //-- Solution rounding option |
5240
|
711 OCTAVE_GLPK_GET_INT_PARAM ("round", lpxIntParam[4]); |
|
712 if (lpxIntParam[4] < 0 || lpxIntParam[4] > 1) |
5234
|
713 { |
5240
|
714 error ("__glpk__: param.round must be 0 (report all primal and dual values [default]) or 1 (replace tiny primal and dual values by exact zero)"); |
|
715 return retval; |
5234
|
716 } |
|
717 |
|
718 //-- Simplex iterations limit |
5240
|
719 OCTAVE_GLPK_GET_INT_PARAM ("itlim", lpxIntParam[5]); |
5232
|
720 |
5234
|
721 //-- Simplex iterations count |
5240
|
722 OCTAVE_GLPK_GET_INT_PARAM ("itcnt", lpxIntParam[6]); |
5234
|
723 |
|
724 //-- Output frequency, in iterations |
5240
|
725 OCTAVE_GLPK_GET_INT_PARAM ("outfrq", lpxIntParam[7]); |
5234
|
726 |
|
727 //-- Branching heuristic option |
5240
|
728 OCTAVE_GLPK_GET_INT_PARAM ("branch", lpxIntParam[14]); |
|
729 if (lpxIntParam[14] < 0 || lpxIntParam[14] > 2) |
5234
|
730 { |
5240
|
731 error ("__glpk__: param.branch must be (MIP only) 0 (branch on first variable) or 1 (branch on last variable) or 2 (branch using a heuristic by Driebeck and Tomlin [default]"); |
|
732 return retval; |
|
733 } |
5232
|
734 |
5234
|
735 //-- Backtracking heuristic option |
5240
|
736 OCTAVE_GLPK_GET_INT_PARAM ("btrack", lpxIntParam[15]); |
|
737 if (lpxIntParam[15] < 0 || lpxIntParam[15] > 2) |
5234
|
738 { |
5240
|
739 error ("__glpk__: param.btrack must be (MIP only) 0 (depth first search) or 1 (breadth first search) or 2 (backtrack using the best projection heuristic [default]"); |
|
740 return retval; |
|
741 } |
5232
|
742 |
5234
|
743 //-- Presolver option |
5240
|
744 OCTAVE_GLPK_GET_INT_PARAM ("presol", lpxIntParam[16]); |
|
745 if (lpxIntParam[16] < 0 || lpxIntParam[16] > 1) |
5234
|
746 { |
5240
|
747 error ("__glpk__: param.presol must be 0 (do NOT use LP presolver) or 1 (use LP presolver [default])"); |
|
748 return retval; |
5234
|
749 } |
|
750 |
5237
|
751 //-- LPsolver option |
|
752 volatile int lpsolver = 1; |
5240
|
753 OCTAVE_GLPK_GET_INT_PARAM ("lpsolver", lpsolver); |
|
754 if (lpsolver < 1 || lpsolver > 2) |
5237
|
755 { |
5240
|
756 error ("__glpk__: param.lpsolver must be 1 (simplex method) or 2 (interior point method)"); |
|
757 return retval; |
5237
|
758 } |
|
759 |
|
760 //-- Save option |
|
761 volatile int save_pb = 0; |
5240
|
762 OCTAVE_GLPK_GET_INT_PARAM ("save", save_pb); |
|
763 save_pb = save_pb != 0; |
5237
|
764 |
5234
|
765 //-- ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|
766 //-- Real parameters |
|
767 //-- ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|
768 |
|
769 //-- Ratio test option |
5240
|
770 OCTAVE_GLPK_GET_REAL_PARAM ("relax", 0); |
5232
|
771 |
5234
|
772 //-- Relative tolerance used to check if the current basic solution |
|
773 //-- is primal feasible |
5240
|
774 OCTAVE_GLPK_GET_REAL_PARAM ("tolbnd", 1); |
5234
|
775 |
|
776 //-- Absolute tolerance used to check if the current basic solution |
|
777 //-- is dual feasible |
5240
|
778 OCTAVE_GLPK_GET_REAL_PARAM ("toldj", 2); |
5232
|
779 |
5234
|
780 //-- Relative tolerance used to choose eligible pivotal elements of |
|
781 //-- the simplex table in the ratio test |
5240
|
782 OCTAVE_GLPK_GET_REAL_PARAM ("tolpiv", 3); |
5234
|
783 |
5240
|
784 OCTAVE_GLPK_GET_REAL_PARAM ("objll", 4); |
5234
|
785 |
5240
|
786 OCTAVE_GLPK_GET_REAL_PARAM ("objul", 5); |
5232
|
787 |
5240
|
788 OCTAVE_GLPK_GET_REAL_PARAM ("tmlim", 6); |
5234
|
789 |
5240
|
790 OCTAVE_GLPK_GET_REAL_PARAM ("outdly", 7); |
5234
|
791 |
5240
|
792 OCTAVE_GLPK_GET_REAL_PARAM ("tolint", 8); |
5234
|
793 |
5240
|
794 OCTAVE_GLPK_GET_REAL_PARAM ("tolobj", 9); |
5234
|
795 |
|
796 //-- Assign pointers to the output parameters |
|
797 ColumnVector xmin (mrowsc); |
|
798 ColumnVector fmin (1); |
|
799 ColumnVector status (1); |
|
800 ColumnVector lambda (mrowsA); |
|
801 ColumnVector redcosts (mrowsc); |
|
802 ColumnVector time (1); |
|
803 ColumnVector mem (1); |
5232
|
804 |
5234
|
805 int jmpret = setjmp (mark); |
5235
|
806 |
5234
|
807 if (jmpret == 0) |
5235
|
808 glpk (sense, mrowsc, mrowsA, c, nz, rn.fortran_vec (), |
|
809 cn.fortran_vec (), a.fortran_vec (), b, ctype, |
|
810 freeLB.fortran_vec (), lb, freeUB.fortran_vec (), |
|
811 ub, vartype.fortran_vec (), isMIP, lpsolver, |
|
812 save_pb, xmin.fortran_vec (), fmin.fortran_vec (), |
|
813 status.fortran_vec (), lambda.fortran_vec (), |
|
814 redcosts.fortran_vec (), time.fortran_vec (), |
|
815 mem.fortran_vec ()); |
5232
|
816 |
5234
|
817 Octave_map extra; |
|
818 |
5238
|
819 if (! isMIP) |
|
820 { |
|
821 extra.assign ("lambda", octave_value (lambda)); |
|
822 extra.assign ("redcosts", octave_value (redcosts)); |
|
823 } |
|
824 |
5234
|
825 extra.assign ("time", octave_value (time)); |
|
826 extra.assign ("mem", octave_value (mem)); |
5232
|
827 |
5234
|
828 retval(3) = extra; |
|
829 retval(2) = octave_value(status); |
|
830 retval(1) = octave_value(fmin); |
|
831 retval(0) = octave_value(xmin); |
|
832 |
5235
|
833 #else |
|
834 |
|
835 gripe_not_supported ("glpk"); |
|
836 |
|
837 #endif |
|
838 |
5234
|
839 return retval; |
5232
|
840 } |