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