1
|
1 // g-builtins.cc -*- C++ -*- |
|
2 /* |
|
3 |
|
4 Copyright (C) 1992, 1993 John W. Eaton |
|
5 |
|
6 This file is part of Octave. |
|
7 |
|
8 Octave is free software; you can redistribute it and/or modify it |
|
9 under the terms of the GNU General Public License as published by the |
|
10 Free Software Foundation; either version 2, or (at your option) any |
|
11 later version. |
|
12 |
|
13 Octave is distributed in the hope that it will be useful, but WITHOUT |
|
14 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
|
15 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License |
|
16 for more details. |
|
17 |
|
18 You should have received a copy of the GNU General Public License |
|
19 along with Octave; see the file COPYING. If not, write to the Free |
|
20 Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA. |
|
21 |
|
22 */ |
|
23 |
|
24 /* |
|
25 |
|
26 The function builtin_pwd adapted from a similar function from GNU |
|
27 Bash, the Bourne Again SHell, copyright (C) 1987, 1989, 1991 Free |
|
28 Software Foundation, Inc. |
|
29 |
|
30 */ |
|
31 |
|
32 #ifdef __GNUG__ |
|
33 #pragma implementation |
|
34 #endif |
|
35 |
|
36 #include <sys/types.h> |
|
37 #ifdef HAVE_UNISTD_H |
|
38 #include <unistd.h> |
|
39 #endif |
|
40 #include <strstream.h> |
|
41 #include <iostream.h> |
|
42 #include <fstream.h> |
|
43 #include <stdio.h> |
|
44 #include <sys/wait.h> |
|
45 #include <sys/param.h> |
|
46 #include <signal.h> |
|
47 #include <math.h> |
|
48 |
9
|
49 #include "f-colloc.h" |
|
50 #include "f-dassl.h" |
|
51 #include "f-det.h" |
|
52 #include "f-eig.h" |
|
53 #include "f-fft.h" |
|
54 #include "f-fsolve.h" |
|
55 #include "f-fsqp.h" |
|
56 #include "f-hess.h" |
|
57 #include "f-ifft.h" |
|
58 #include "f-inv.h" |
|
59 #include "f-lpsolve.h" |
|
60 #include "f-lsode.h" |
|
61 #include "f-lu.h" |
|
62 #include "f-npsol.h" |
|
63 #include "f-qpsol.h" |
|
64 #include "f-qr.h" |
|
65 #include "f-quad.h" |
|
66 #include "f-rand.h" |
|
67 #include "f-schur.h" |
|
68 #include "f-svd.h" |
|
69 |
1
|
70 #include "procstream.h" |
|
71 #include "error.h" |
|
72 #include "variables.h" |
|
73 #include "builtins.h" |
|
74 #include "g-builtins.h" |
|
75 #include "user-prefs.h" |
|
76 #include "utils.h" |
|
77 #include "tree.h" |
9
|
78 #include "tree-const.h" |
1
|
79 #include "input.h" |
|
80 #include "pager.h" |
|
81 #include "octave.h" |
|
82 #include "version.h" |
|
83 #include "file-io.h" |
|
84 |
|
85 extern "C" |
|
86 { |
|
87 #include <readline/readline.h> |
|
88 } |
|
89 |
|
90 #ifndef MAXPATHLEN |
|
91 #define MAXPATHLEN 1024 |
|
92 #endif |
|
93 |
|
94 #ifdef WITH_DLD |
|
95 #include "dynamic-ld.h" |
|
96 #define Q_STR(name) #name |
|
97 #define DLD_FCN(name) Q_STR (builtin_##name##_2) |
|
98 #define DLD_OBJ(name) Q_STR (tc-##name##.o) |
|
99 #define DLD_BUILTIN(args,n_in,n_out,name,code) \ |
|
100 return octave_dld_tc2_and_go (args, n_in, n_out, Q_STR (name), \ |
|
101 DLD_FCN (name), DLD_OBJ (name)); |
|
102 |
|
103 #else |
|
104 #define DLD_BUILTIN(name,args,n_in,n_out,code) code |
|
105 #endif |
|
106 |
|
107 // Non-zero means that pwd always give verbatim directory, regardless |
|
108 // of symbolic link following. |
|
109 static int verbatim_pwd = 1; |
|
110 |
|
111 // Signal handler return type. |
|
112 #ifndef RETSIGTYPE |
|
113 #define RETSIGTYPE void |
|
114 #endif |
|
115 #ifndef BADSIG |
|
116 #define BADSIG (RETSIGTYPE (*)())-1 |
|
117 #endif |
|
118 |
|
119 typedef RETSIGTYPE sig_handler (...); |
|
120 |
|
121 /* |
|
122 * Are all elements of a constant nonzero? |
|
123 */ |
|
124 tree_constant * |
|
125 builtin_all (tree_constant *args, int nargin, int nargout) |
|
126 { |
|
127 tree_constant *retval = NULL_TREE_CONST; |
|
128 if (nargin != 2) |
|
129 usage ("all (M)"); |
|
130 else |
|
131 { |
|
132 if (args != NULL_TREE_CONST && args[1].is_defined ()) |
|
133 { |
|
134 retval = new tree_constant [2]; |
|
135 retval[0] = args[1].all (); |
|
136 } |
|
137 } |
|
138 return retval; |
|
139 } |
|
140 |
|
141 /* |
|
142 * Are any elements of a constant nonzero? |
|
143 */ |
|
144 tree_constant * |
|
145 builtin_any (tree_constant *args, int nargin, int nargout) |
|
146 { |
|
147 tree_constant *retval = NULL_TREE_CONST; |
|
148 if (nargin != 2) |
|
149 usage ("any (M)"); |
|
150 else |
|
151 { |
|
152 if (args != NULL_TREE_CONST && args[1].is_defined ()) |
|
153 { |
|
154 retval = new tree_constant [2]; |
|
155 retval[0] = args[1].any (); |
|
156 } |
|
157 } |
|
158 return retval; |
|
159 } |
|
160 |
|
161 /* |
18
|
162 * Balancing for eigenvalue problems |
|
163 */ |
|
164 tree_constant * |
|
165 builtin_balance (tree_constant *args, int nargin, int nargout) |
|
166 { |
|
167 tree_constant *retval = NULL_TREE_CONST; |
|
168 if (nargin <= 1 || nargin > 4 || nargout < 1 || nargout > 4) |
|
169 usage ("[aa {,dd}] = balance (a, {opt}) or \n\ |
|
170 [aa, bb {,cc, dd}] = balance (a, b {,opt}), opt = 'P' or 'S'"); |
|
171 else |
|
172 { |
|
173 DLD_BUILTIN (args, nargin, nargout, balance, |
|
174 retval = balance (args, nargin, nargout)); |
|
175 } |
|
176 return retval; |
|
177 } |
|
178 |
|
179 /* |
1
|
180 * Clear the screen? |
|
181 */ |
|
182 tree_constant * |
|
183 builtin_clc (tree_constant *args, int nargin, int nargout) |
|
184 { |
|
185 rl_beg_of_line (); |
|
186 rl_kill_line (1); |
|
187 rl_clear_screen (); |
|
188 return NULL_TREE_CONST; |
|
189 } |
|
190 |
|
191 /* |
|
192 * Time in a vector. |
|
193 */ |
|
194 tree_constant * |
|
195 builtin_clock (tree_constant *args, int nargin, int nargout) |
|
196 { |
|
197 tree_constant *retval = NULL_TREE_CONST; |
|
198 |
|
199 time_t now; |
|
200 struct tm *tm; |
|
201 |
|
202 time (&now); |
|
203 tm = localtime (&now); |
|
204 |
|
205 Matrix m (1, 6); |
|
206 m.elem (0, 0) = tm->tm_year + 1900; |
|
207 m.elem (0, 1) = tm->tm_mon + 1; |
|
208 m.elem (0, 2) = tm->tm_mday; |
|
209 m.elem (0, 3) = tm->tm_hour; |
|
210 m.elem (0, 4) = tm->tm_min; |
|
211 m.elem (0, 5) = tm->tm_sec; |
|
212 |
|
213 retval = new tree_constant [2]; |
|
214 retval[0] = tree_constant (m); |
|
215 |
|
216 return retval; |
|
217 } |
|
218 |
|
219 /* |
|
220 * Close the stream to the plotter. |
|
221 */ |
|
222 tree_constant * |
|
223 builtin_closeplot (tree_constant *args, int nargin, int nargout) |
|
224 { |
|
225 tree_constant *retval = NULL_TREE_CONST; |
|
226 close_plot_stream (); |
|
227 return retval; |
|
228 } |
|
229 |
|
230 /* |
|
231 * Collocation roots and weights. |
|
232 */ |
|
233 tree_constant * |
|
234 builtin_colloc (tree_constant *args, int nargin, int nargout) |
|
235 { |
|
236 tree_constant *retval = NULL_TREE_CONST; |
|
237 |
|
238 if (nargin < 2 || nargin > 4) |
|
239 usage ("[r, A, B, q] = colloc (n [, \"left\"] [, \"right\"])"); |
|
240 else |
|
241 DLD_BUILTIN (args, nargin, nargout, colloc, |
|
242 retval = collocation_weights (args, nargin);) |
|
243 |
|
244 return retval; |
|
245 } |
|
246 |
|
247 /* |
|
248 * Cumulative sums and products. |
|
249 */ |
|
250 tree_constant * |
|
251 builtin_cumprod (tree_constant *args, int nargin, int nargout) |
|
252 { |
|
253 tree_constant *retval = NULL_TREE_CONST; |
|
254 if (nargin != 2) |
|
255 usage ("cumprod (M)"); |
|
256 else |
|
257 { |
|
258 if (args != NULL_TREE_CONST && args[1].is_defined ()) |
|
259 { |
|
260 retval = new tree_constant [2]; |
|
261 retval[0] = args[1].cumprod (); |
|
262 } |
|
263 } |
|
264 return retval; |
|
265 } |
|
266 |
|
267 tree_constant * |
|
268 builtin_cumsum (tree_constant *args, int nargin, int nargout) |
|
269 { |
|
270 tree_constant *retval = NULL_TREE_CONST; |
|
271 if (nargin != 2) |
|
272 usage ("cumsum (M)"); |
|
273 else |
|
274 { |
|
275 if (args != NULL_TREE_CONST && args[1].is_defined ()) |
|
276 { |
|
277 retval = new tree_constant [2]; |
|
278 retval[0] = args[1].cumsum (); |
|
279 } |
|
280 } |
|
281 return retval; |
|
282 } |
|
283 |
|
284 /* |
|
285 * DAEs. |
|
286 */ |
|
287 static void |
|
288 dassl_usage (void) |
|
289 { |
|
290 usage ("dassl (\"function_name\", x_0, xdot_0, t_out\n\ |
|
291 dassl (\"function_name\", x_0, xdot_0, t_out, t_crit)\n\ |
|
292 \n\ |
|
293 The first argument is the name of the function to call to\n\ |
|
294 compute the vector of residuals. It must have the form\n\ |
|
295 \n\ |
|
296 res = f (x, xdot, t)\n\ |
|
297 \n\ |
|
298 where x, xdot, and res are vectors, and t is a scalar."); |
|
299 } |
|
300 |
|
301 tree_constant * |
|
302 builtin_dassl (tree_constant *args, int nargin, int nargout) |
|
303 { |
|
304 tree_constant *retval = new tree_constant [2]; |
|
305 |
|
306 if ((nargin == 5 || nargin == 6) && nargout > 0) |
|
307 DLD_BUILTIN (args, nargin, nargout, dassl, |
|
308 retval = dassl (args, nargin, nargout);) |
|
309 else |
|
310 dassl_usage (); |
|
311 |
|
312 return retval; |
|
313 } |
|
314 |
|
315 /* |
|
316 * Time in a string. |
|
317 */ |
|
318 tree_constant * |
|
319 builtin_date (tree_constant *args, int nargin, int nargout) |
|
320 { |
|
321 tree_constant *retval = NULL_TREE_CONST; |
|
322 |
|
323 time_t now; |
|
324 struct tm *tm; |
|
325 |
|
326 time (&now); |
|
327 tm = localtime (&now); |
|
328 char date[32]; |
|
329 int len = strftime (date, 31, "%d-%b-%y", tm); |
|
330 if (len > 0) |
|
331 { |
|
332 retval = new tree_constant [2]; |
|
333 retval[0] = tree_constant (date); |
|
334 } |
|
335 |
|
336 return retval; |
|
337 } |
|
338 |
|
339 /* |
|
340 * Determinant of a matrix. |
|
341 */ |
|
342 tree_constant * |
|
343 builtin_det (tree_constant *args, int nargin, int nargout) |
|
344 { |
|
345 tree_constant *retval = NULL_TREE_CONST; |
|
346 |
|
347 if (nargin == 2) |
|
348 DLD_BUILTIN (args, nargin, nargout, det, |
|
349 { |
|
350 retval = new tree_constant [2]; |
|
351 retval[0] = determinant (args[1]); |
|
352 }) |
|
353 else |
|
354 usage ("det (a)"); |
|
355 |
|
356 return retval; |
|
357 } |
|
358 |
|
359 /* |
|
360 * Diagonal elements of a matrix. |
|
361 */ |
|
362 tree_constant * |
|
363 builtin_diag (tree_constant *args, int nargin, int nargout) |
|
364 { |
|
365 tree_constant *retval = NULL_TREE_CONST; |
|
366 |
|
367 if (nargin == 2) |
|
368 { |
|
369 retval = new tree_constant [2]; |
|
370 retval[0] = args[1].diag (); |
|
371 } |
|
372 else if (nargin == 3) |
|
373 { |
|
374 retval = new tree_constant [2]; |
|
375 retval[0] = args[1].diag (args[2]); |
|
376 } |
|
377 else |
|
378 usage ("diag (X [, k])"); |
|
379 |
|
380 return retval; |
|
381 } |
|
382 |
|
383 /* |
|
384 * Display value without trimmings. |
|
385 */ |
|
386 tree_constant * |
|
387 builtin_disp (tree_constant *args, int nargin, int nargout) |
|
388 { |
|
389 tree_constant *retval = NULL_TREE_CONST; |
|
390 |
|
391 if (nargin == 2) |
|
392 args[1].eval (1); |
|
393 else |
|
394 usage ("disp (X)"); |
|
395 |
|
396 return retval; |
|
397 } |
|
398 |
|
399 /* |
|
400 * Compute eigenvalues and eigenvectors. |
|
401 */ |
|
402 tree_constant * |
|
403 builtin_eig (tree_constant *args, int nargin, int nargout) |
|
404 { |
|
405 tree_constant *retval = NULL_TREE_CONST; |
|
406 |
|
407 if (nargin == 2 && (nargout == 1 || nargout == 2)) |
|
408 DLD_BUILTIN (args, nargin, nargout, eig, |
|
409 retval = eig (args, nargin, nargout);) |
|
410 else |
|
411 usage ("lambda = eig (A)\n\ |
|
412 [v, d] = eig (A); d == diag (lambda)"); |
|
413 |
|
414 return retval; |
|
415 } |
|
416 |
|
417 /* |
|
418 * Print error message and jump to top level. |
|
419 */ |
|
420 tree_constant * |
|
421 builtin_error (tree_constant *args, int nargin, int nargout) |
|
422 { |
|
423 tree_constant *retval = NULL_TREE_CONST; |
|
424 |
|
425 if (nargin == 2 && args != NULL_TREE_CONST && args[1].is_defined ()) |
|
426 args[1].print_if_string (cerr, 1); |
|
427 else |
|
428 message ((char *) NULL, "unspecified error, jumping to top level..."); |
|
429 |
|
430 jump_to_top_level (); |
|
431 |
|
432 return retval; |
|
433 } |
|
434 |
|
435 /* |
|
436 * Evaluate text argument as octave source. |
|
437 */ |
|
438 tree_constant * |
|
439 builtin_eval (tree_constant *args, int nargin, int nargout) |
|
440 { |
|
441 tree_constant *retval = NULL_TREE_CONST; |
|
442 if (nargin == 2) |
|
443 { |
|
444 int parse_status = 0; |
|
445 retval = new tree_constant [2]; |
|
446 retval[0] = eval_string (args[1], parse_status); |
|
447 } |
|
448 else |
|
449 usage ("eval (\"string\")"); |
|
450 return retval; |
|
451 } |
|
452 |
|
453 /* |
|
454 * Check if variable or file exists. |
|
455 */ |
|
456 tree_constant * |
|
457 builtin_exist (tree_constant *args, int nargin, int nargout) |
|
458 { |
|
459 tree_constant *retval = NULL_TREE_CONST; |
|
460 if (nargin == 2 && args[1].is_string_type ()) |
|
461 { |
|
462 int status = identifier_exists (args[1].string_value ()); |
|
463 retval = new tree_constant [2]; |
|
464 retval[0] = tree_constant ((double) status); |
|
465 } |
|
466 else |
|
467 usage ("exist (\"string\")"); |
|
468 return retval; |
|
469 } |
|
470 |
|
471 /* |
|
472 * Matrix exponential. |
|
473 */ |
|
474 tree_constant * |
|
475 builtin_expm (tree_constant *args, int nargin, int nargout) |
|
476 { |
|
477 tree_constant *retval = NULL_TREE_CONST; |
|
478 |
|
479 if (nargin == 2) |
|
480 retval = matrix_exp (args[1]); |
|
481 else |
|
482 usage ("expm (A)"); |
|
483 |
|
484 return retval; |
|
485 } |
|
486 |
|
487 /* |
|
488 * Identity matrix. |
|
489 */ |
|
490 tree_constant * |
|
491 builtin_eye (tree_constant *args, int nargin, int nargout) |
|
492 { |
|
493 tree_constant *retval = NULL_TREE_CONST; |
|
494 |
|
495 switch (nargin) |
|
496 { |
|
497 case 2: |
|
498 retval = new tree_constant [2]; |
|
499 retval[0] = identity_matrix (args[1]); |
|
500 break; |
|
501 case 3: |
|
502 retval = new tree_constant [2]; |
|
503 retval[0] = identity_matrix (args[1], args[2]); |
|
504 break; |
|
505 default: |
|
506 usage ("eye (n)\n eye (A)\n eye (n, m)"); |
|
507 break; |
|
508 } |
|
509 return retval; |
|
510 } |
|
511 |
|
512 |
|
513 /* |
|
514 * Closing a file |
|
515 */ |
|
516 tree_constant * |
|
517 builtin_fclose (tree_constant *args, int nargin, int nargout) |
|
518 { |
|
519 tree_constant *retval = NULL_TREE_CONST; |
|
520 if (nargin != 2) |
|
521 usage ("success = fclose (\"filename\" or filenum)"); |
|
522 else |
|
523 retval = fclose_internal (args); |
|
524 return retval; |
|
525 } |
|
526 |
|
527 /* |
|
528 * Evaluate first argument as a function. |
|
529 */ |
|
530 tree_constant * |
|
531 builtin_feval (tree_constant *args, int nargin, int nargout) |
|
532 { |
|
533 tree_constant *retval = NULL_TREE_CONST; |
|
534 if (nargin > 1) |
|
535 retval = feval (args, nargin, nargout); |
|
536 else |
|
537 usage ("feval (\"function_name\" [, ...])"); |
|
538 return retval; |
|
539 } |
|
540 |
|
541 /* |
|
542 * Flushing output to a file |
|
543 */ |
|
544 tree_constant * |
|
545 builtin_fflush (tree_constant *args, int nargin, int nargout) |
|
546 { |
|
547 tree_constant *retval = NULL_TREE_CONST; |
|
548 if (nargin != 2) |
|
549 usage ("success = fflush (\"filename\" or filenum)"); |
|
550 else |
|
551 retval = fflush_internal (args); |
|
552 return retval; |
|
553 } |
|
554 |
|
555 /* |
|
556 * Fast Fourier Transform |
|
557 */ |
|
558 tree_constant * |
|
559 builtin_fft (tree_constant *args, int nargin, int nargout) |
|
560 { |
|
561 tree_constant *retval = NULL_TREE_CONST; |
|
562 |
|
563 if (nargin == 2) |
|
564 DLD_BUILTIN (args, nargin, nargout, fft, |
|
565 { |
|
566 retval = new tree_constant [2]; |
|
567 retval[0] = fft (args[1]); |
|
568 }) |
|
569 else |
|
570 usage ("fft (a)"); |
|
571 |
|
572 return retval; |
|
573 } |
|
574 |
|
575 /* |
|
576 * get a string from a file |
|
577 */ |
|
578 tree_constant * |
|
579 builtin_fgets (tree_constant *args, int nargin, int nargout) |
|
580 { |
|
581 tree_constant *retval = NULL_TREE_CONST; |
|
582 if (nargin != 3 && nargout < 3) |
|
583 usage ("string = fgets (\"filename\" or filenum, length)"); |
|
584 else |
|
585 retval = fgets_internal (args, nargout); |
|
586 return retval; |
|
587 } |
|
588 |
|
589 /* |
|
590 * Find nonzero elements. This should probably only work if |
|
591 * do_fortran_indexing is true... |
|
592 */ |
|
593 tree_constant * |
|
594 builtin_find (tree_constant *args, int nargin, int nargout) |
|
595 { |
|
596 tree_constant *retval = NULL_TREE_CONST; |
|
597 if (nargin == 2) |
|
598 { |
|
599 retval = new tree_constant [2]; |
|
600 retval[0] = find_nonzero_elem_idx (args[1]); |
|
601 } |
|
602 else |
|
603 usage ("find (x)"); |
|
604 return retval; |
|
605 } |
|
606 |
|
607 /* |
|
608 * Don\'t really count floating point operations. |
|
609 */ |
|
610 tree_constant * |
|
611 builtin_flops (tree_constant *args, int nargin, int nargout) |
|
612 { |
|
613 tree_constant *retval = NULL_TREE_CONST; |
|
614 if (nargin > 2) |
|
615 usage ("flops\n flops (0)"); |
|
616 |
|
617 warning ("flops always returns zero"); |
|
618 retval = new tree_constant [2]; |
|
619 retval[0] = tree_constant (0.0); |
|
620 return retval; |
|
621 } |
|
622 |
|
623 /* |
|
624 * Opening a file. |
|
625 */ |
|
626 tree_constant * |
|
627 builtin_fopen (tree_constant *args, int nargin, int nargout) |
|
628 { |
|
629 tree_constant *retval = NULL_TREE_CONST; |
|
630 if (nargin != 3) |
|
631 { |
|
632 usage ("filenum = fopen (\"file\", \"mode\")\n\n\ |
|
633 Legal values for mode include:\n\n\ |
|
634 r : open text file for reading\n\ |
|
635 w : open text file for writing; discard previous contents if any\n\ |
|
636 a : append; open or create text file for writing at end of file\n\ |
|
637 r+ : open text file for update (i.e., reading and writing)\n\ |
|
638 w+ : create text file for update; discard previous contents if any\n\ |
|
639 a+ : append; open or create text file for update, writing at end\n\n\ |
|
640 Update mode permits reading from and writing to the same file.\n"); |
|
641 } |
|
642 else |
|
643 retval = fopen_internal (args); |
|
644 return retval; |
|
645 } |
|
646 |
|
647 /* |
|
648 * Formatted printing to a file. |
|
649 */ |
|
650 tree_constant * |
|
651 builtin_fprintf (tree_constant *args, int nargin, int nargout) |
|
652 { |
|
653 tree_constant *retval = NULL_TREE_CONST; |
|
654 if (nargin < 3) |
|
655 usage ("fprintf (\"filename\" or filenum, \"fmt\" [, ...])"); |
|
656 else |
|
657 retval = do_printf ("fprintf", args, nargin, nargout); |
|
658 return retval; |
|
659 } |
|
660 |
|
661 /* |
|
662 * rewind a file |
|
663 */ |
|
664 tree_constant * |
|
665 builtin_frewind (tree_constant *args, int nargin, int nargout) |
|
666 { |
|
667 tree_constant *retval = NULL_TREE_CONST; |
|
668 if (nargin != 2) |
|
669 usage ("success = frewind (\"filename\" or filenum)"); |
|
670 else |
|
671 retval = frewind_internal (args); |
|
672 return retval; |
|
673 } |
|
674 |
|
675 /* |
|
676 * report on open files |
|
677 */ |
|
678 tree_constant * |
|
679 builtin_freport (tree_constant *args, int nargin, int nargout) |
|
680 { |
|
681 tree_constant *retval = NULL_TREE_CONST; |
|
682 if (nargin > 1) |
|
683 warning ("replot: ignoring extra arguments"); |
|
684 retval = freport_internal (); |
|
685 return retval; |
|
686 } |
|
687 |
|
688 /* |
|
689 * Formatted reading from a file. |
|
690 */ |
|
691 tree_constant * |
|
692 builtin_fscanf (tree_constant *args, int nargin, int nargout) |
|
693 { |
|
694 tree_constant *retval = NULL_TREE_CONST; |
|
695 if (nargin != 2 && nargin != 3) |
|
696 usage ("[...] = fscanf (\"file\", \"fmt\")"); |
|
697 else |
|
698 retval = do_scanf ("fscanf", args, nargin, nargout); |
|
699 return retval; |
|
700 } |
|
701 |
|
702 /* |
|
703 * seek a point in a file for reading and/or writing |
|
704 */ |
|
705 tree_constant * |
|
706 builtin_fseek (tree_constant *args, int nargin, int nargout) |
|
707 { |
|
708 tree_constant *retval = NULL_TREE_CONST; |
|
709 if (nargin != 3 && nargin != 4) |
|
710 usage ("success = fseek (\"filename\" or filenum, offset [,origin])"); |
|
711 else |
|
712 retval = fseek_internal (args, nargin); |
|
713 return retval; |
|
714 } |
|
715 |
|
716 /* |
|
717 * Nonlinear algebraic equations. |
|
718 */ |
|
719 static void |
|
720 fsolve_usage (void) |
|
721 { |
|
722 // usage ("[x, status, path] = fsolve (\"f\", x0 [, opts] [, par] [, \"jac\"] [, scale])"); |
|
723 |
|
724 usage ("[x, info] = fsolve (\"f\", x0)"); |
|
725 } |
|
726 |
|
727 tree_constant * |
|
728 builtin_fsolve (tree_constant *args, int nargin, int nargout) |
|
729 { |
|
730 tree_constant *retval = NULL_TREE_CONST; |
|
731 |
|
732 if (nargin >= 3 && nargin <= 7 && nargout >= 1 && nargout <= 3) |
|
733 DLD_BUILTIN (args, nargin, nargout, fsolve, |
|
734 retval = fsolve (args, nargin, nargout);) |
|
735 else |
|
736 fsolve_usage (); |
|
737 |
|
738 return retval; |
|
739 } |
|
740 |
|
741 /* |
|
742 * NLPs. |
|
743 */ |
|
744 static void |
|
745 fsqp_usage (void) |
|
746 { |
|
747 #if defined (FSQP_MISSING) |
|
748 message ("fsqp", "this function requires FSQP, which is not freely\n\ |
|
749 redistributable. For more information, read the file\n\ |
|
750 libcruft/fsqp/README.MISSING in the source distribution."); |
|
751 #else |
|
752 usage ("[x, phi] = fsqp (x, \"phi\" [, lb, ub] [, lb, A, ub] [, lb, \"g\", ub])\n\n\ |
|
753 Groups of arguments surrounded in `[]' are optional, but\n\ |
|
754 must appear in the same relative order shown above."); |
|
755 #endif |
|
756 } |
|
757 |
|
758 tree_constant * |
|
759 builtin_fsqp (tree_constant *args, int nargin, int nargout) |
|
760 { |
|
761 tree_constant *retval = NULL_TREE_CONST; |
|
762 |
|
763 #if defined (FSQP_MISSING) |
|
764 fsqp_usage (); |
|
765 #else |
|
766 if ((nargin == 3 || nargin == 5 || nargin == 6 || nargin == 8 |
|
767 || nargin == 9 || nargin == 11) |
|
768 && (nargout >= 1 && nargout <= 3)) |
|
769 DLD_BUILTIN (args, nargin, nargout, fsqp, |
|
770 retval = fsqp (args, nargin, nargout);) |
|
771 else |
|
772 fsqp_usage (); |
|
773 #endif |
|
774 |
|
775 return retval; |
|
776 } |
|
777 |
|
778 /* |
|
779 * tell current position of file |
|
780 */ |
|
781 tree_constant * |
|
782 builtin_ftell (tree_constant *args, int nargin, int nargout) |
|
783 { |
|
784 tree_constant *retval = NULL_TREE_CONST; |
|
785 if (nargin != 2) |
|
786 usage ("position = ftell (\"filename\" or filenumber)"); |
|
787 else |
|
788 retval = ftell_internal (args); |
|
789 return retval; |
|
790 } |
|
791 |
|
792 /* |
|
793 * Get the value of an environment variable. |
|
794 */ |
|
795 tree_constant * |
|
796 builtin_getenv (tree_constant *args, int nargin, int nargout) |
|
797 { |
|
798 tree_constant *retval = NULL_TREE_CONST; |
|
799 if (nargin == 2 && args[1].is_string_type ()) |
|
800 { |
|
801 retval = new tree_constant [2]; |
|
802 char *value = getenv (args[1].string_value ()); |
|
803 if (value != (char *) NULL) |
|
804 retval[0] = tree_constant (value); |
|
805 else |
|
806 retval[0] = tree_constant (""); |
|
807 } |
|
808 else |
|
809 usage ("getenv (\"string\")"); |
|
810 return retval; |
|
811 } |
|
812 |
|
813 /* |
|
814 * Inverse Fast Fourier Transform |
|
815 */ |
|
816 tree_constant * |
|
817 builtin_ifft (tree_constant *args, int nargin, int nargout) |
|
818 { |
|
819 tree_constant *retval = NULL_TREE_CONST; |
|
820 |
|
821 if (nargin == 2) |
|
822 DLD_BUILTIN (args, nargin, nargout, ifft, |
|
823 { |
|
824 retval = new tree_constant [2]; |
|
825 retval[0] = ifft (args[1]); |
|
826 }) |
|
827 else |
|
828 usage ("ifft (a)"); |
|
829 |
|
830 return retval; |
|
831 } |
|
832 |
|
833 /* |
|
834 * Inverse of a square matrix. |
|
835 */ |
|
836 tree_constant * |
|
837 builtin_inv (tree_constant *args, int nargin, int nargout) |
|
838 { |
|
839 tree_constant *retval = NULL_TREE_CONST; |
|
840 |
|
841 if (nargin == 2) |
|
842 DLD_BUILTIN (args, nargin, nargout, inv, |
|
843 { |
|
844 retval = new tree_constant [2]; |
|
845 retval[0] = inverse (args[1]); |
|
846 }) |
|
847 else |
|
848 usage ("inv (A)"); |
|
849 |
|
850 return retval; |
|
851 } |
|
852 |
|
853 /* |
|
854 * Prompt user for input. |
|
855 */ |
|
856 tree_constant * |
|
857 builtin_input (tree_constant *args, int nargin, int nargout) |
|
858 { |
|
859 tree_constant *retval = NULL_TREE_CONST; |
|
860 |
|
861 if (nargin == 2 || nargin == 3) |
|
862 { |
|
863 retval = new tree_constant [2]; |
|
864 retval[0] = get_user_input (args, nargin, nargout); |
|
865 } |
|
866 else |
|
867 usage ("input (\"prompt\" [, \"s\"])"); |
|
868 |
|
869 return retval; |
|
870 } |
|
871 |
|
872 /* |
|
873 * Is the argument a string? |
|
874 */ |
|
875 tree_constant * |
|
876 builtin_isstr (tree_constant *args, int nargin, int nargout) |
|
877 { |
|
878 tree_constant *retval = NULL_TREE_CONST; |
|
879 if (nargin != 2) |
|
880 usage ("isstr (value)"); |
|
881 else |
|
882 { |
|
883 if (args != NULL_TREE_CONST && args[1].is_defined ()) |
|
884 { |
|
885 retval = new tree_constant [2]; |
|
886 retval[0] = args[1].isstr (); |
|
887 } |
|
888 } |
|
889 return retval; |
|
890 } |
|
891 |
|
892 /* |
|
893 * Maybe help in debugging. |
|
894 */ |
|
895 tree_constant * |
|
896 builtin_keyboard (tree_constant *args, int nargin, int nargout) |
|
897 { |
|
898 tree_constant *retval = NULL_TREE_CONST; |
|
899 |
|
900 if (nargin == 1 || nargin == 2) |
|
901 { |
|
902 retval = new tree_constant [2]; |
|
903 retval[0] = get_user_input (args, nargin, nargout, 1); |
|
904 } |
|
905 else |
|
906 usage ("keyboard (\"prompt\")"); |
|
907 |
|
908 return retval; |
|
909 } |
|
910 |
|
911 /* |
|
912 * Matrix logarithm. |
|
913 */ |
|
914 tree_constant * |
|
915 builtin_logm (tree_constant *args, int nargin, int nargout) |
|
916 { |
|
917 tree_constant *retval = NULL_TREE_CONST; |
|
918 |
|
919 if (nargin == 2) |
|
920 retval = matrix_log (args[1]); |
|
921 else |
|
922 usage ("logm (A)"); |
|
923 |
|
924 return retval; |
|
925 } |
|
926 |
|
927 /* |
|
928 * LPs. |
|
929 */ |
|
930 static void |
|
931 lpsolve_usage (void) |
|
932 { |
|
933 usage ("[x, obj, info] = lpsolve (XXX FIXME XXX)"); |
|
934 } |
|
935 |
|
936 tree_constant * |
|
937 builtin_lpsolve (tree_constant *args, int nargin, int nargout) |
|
938 { |
|
939 tree_constant *retval = NULL_TREE_CONST; |
|
940 |
|
941 // Force a bad value of inform, and empty matrices for x and phi. |
|
942 retval = new tree_constant [4]; |
|
943 Matrix m; |
|
944 retval[0] = tree_constant (m); |
|
945 retval[1] = tree_constant (m); |
|
946 retval[2] = tree_constant (-1.0); |
|
947 |
|
948 if (nargin == 0) |
|
949 DLD_BUILTIN (args, nargin, nargout, lpsolve, |
|
950 retval = lpsolve (args, nargin, nargout);) |
|
951 else |
|
952 lpsolve_usage (); |
|
953 |
|
954 return retval; |
|
955 } |
|
956 |
|
957 /* |
|
958 * ODEs. |
|
959 */ |
|
960 static void |
|
961 lsode_usage (void) |
|
962 { |
|
963 usage ("lsode (\"function_name\", x0, t_out\n\ |
|
964 lsode (\"function_name\", x0, t_out, t_crit)\n\ |
|
965 \n\ |
|
966 The first argument is the name of the function to call to\n\ |
|
967 compute the vector of right hand sides. It must have the form\n\ |
|
968 \n\ |
|
969 xdot = f (x, t)\n\ |
|
970 \n\ |
|
971 where xdot and x are vectors and t is a scalar."); |
|
972 } |
|
973 |
|
974 tree_constant * |
|
975 builtin_lsode (tree_constant *args, int nargin, int nargout) |
|
976 { |
|
977 tree_constant *retval = NULL_TREE_CONST; |
|
978 |
|
979 if ((nargin == 4 || nargin == 5) && nargout == 1) |
|
980 DLD_BUILTIN (args, nargin, nargout, lsode, |
|
981 retval = lsode (args, nargin, nargout);) |
|
982 else |
|
983 lsode_usage (); |
|
984 |
|
985 return retval; |
|
986 } |
|
987 |
|
988 /* |
|
989 * LU factorization. |
|
990 */ |
|
991 tree_constant * |
|
992 builtin_lu (tree_constant *args, int nargin, int nargout) |
|
993 { |
|
994 tree_constant *retval = NULL_TREE_CONST; |
|
995 |
|
996 if (nargin == 2 && nargout < 4) |
|
997 DLD_BUILTIN (args, nargin, nargout, lu, |
|
998 retval = lu (args[1], nargout);) |
|
999 else |
|
1000 usage ("[L, U, P] = lu (A)"); |
|
1001 |
|
1002 return retval; |
|
1003 } |
|
1004 |
|
1005 /* |
|
1006 * Max values. |
|
1007 */ |
|
1008 tree_constant * |
|
1009 builtin_max (tree_constant *args, int nargin, int nargout) |
|
1010 { |
|
1011 tree_constant *retval = NULL_TREE_CONST; |
|
1012 |
|
1013 if ((nargin == 2 && (nargout == 1 || nargout == 2)) |
|
1014 || (nargin == 3 && nargout == 1)) |
|
1015 retval = column_max (args, nargin, nargout); |
|
1016 else |
|
1017 usage ("[X, I] = max (A)\n X = max (A)\n X = max (A, B)"); |
|
1018 |
|
1019 return retval; |
|
1020 } |
|
1021 |
|
1022 /* |
|
1023 * Min values. |
|
1024 */ |
|
1025 tree_constant * |
|
1026 builtin_min (tree_constant *args, int nargin, int nargout) |
|
1027 { |
|
1028 tree_constant *retval = NULL_TREE_CONST; |
|
1029 |
|
1030 if ((nargin == 2 && (nargout == 1 || nargout == 2)) |
|
1031 || (nargin == 3 && nargout == 1)) |
|
1032 retval = column_min (args, nargin, nargout); |
|
1033 else |
|
1034 usage ("[X, I] = min (A)\n X = min (A)\n X = min (A, B)"); |
|
1035 |
|
1036 return retval; |
|
1037 } |
|
1038 |
|
1039 /* |
|
1040 * NLPs. |
|
1041 */ |
|
1042 static void |
|
1043 npsol_usage (void) |
|
1044 { |
|
1045 #if defined (NPSOL_MISSING) |
|
1046 message ("npsol", "this function requires NPSOL, which is not freely\n\ |
|
1047 redistributable. For more information, read the file\n\ |
|
1048 libcruft/npsol/README.MISSING in the source distribution."); |
|
1049 #else |
|
1050 usage ("\n\n\ |
|
1051 [x, obj, info, lambda] = npsol (x, \"phi\" [, lb, ub] [, lb, A, ub] [, lb, \"g\", ub])\n\n\ |
|
1052 Groups of arguments surrounded in `[]' are optional, but\n\ |
|
1053 must appear in the same relative order shown above.\n\ |
|
1054 \n\ |
|
1055 The second argument is a string containing the name of the objective\n\ |
|
1056 function to call. The objective function must be of the form\n\ |
|
1057 \n\ |
|
1058 y = phi (x)\n\ |
|
1059 \n\ |
|
1060 where x is a vector and y is a scalar."); |
|
1061 #endif |
|
1062 } |
|
1063 |
|
1064 tree_constant * |
|
1065 builtin_npsol (tree_constant *args, int nargin, int nargout) |
|
1066 { |
|
1067 tree_constant *retval = NULL_TREE_CONST; |
|
1068 |
|
1069 #if defined (NPSOL_MISSING) |
|
1070 // Force a bad value of inform, and empty matrices for x, phi, and lambda. |
|
1071 retval = new tree_constant [4]; |
|
1072 Matrix m; |
|
1073 retval[0] = tree_constant (m); |
|
1074 retval[1] = tree_constant (m); |
|
1075 retval[2] = tree_constant (-1.0); |
|
1076 retval[3] = tree_constant (m); |
|
1077 npsol_usage (); |
|
1078 #else |
|
1079 if ((nargin == 3 || nargin == 5 || nargin == 6 || nargin == 8 |
|
1080 || nargin == 9 || nargin == 11) |
|
1081 && (nargout >= 1 && nargout <= 4)) |
|
1082 DLD_BUILTIN (args, nargin, nargout, npsol, |
|
1083 retval = npsol (args, nargin, nargout);) |
|
1084 else |
|
1085 npsol_usage (); |
|
1086 #endif |
|
1087 |
|
1088 return retval; |
|
1089 } |
|
1090 |
|
1091 /* |
|
1092 * A matrix of ones. |
|
1093 */ |
|
1094 tree_constant * |
|
1095 builtin_ones (tree_constant *args, int nargin, int nargout) |
|
1096 { |
|
1097 tree_constant *retval = NULL_TREE_CONST; |
|
1098 |
|
1099 switch (nargin) |
|
1100 { |
|
1101 case 2: |
|
1102 retval = new tree_constant [2]; |
|
1103 retval[0] = fill_matrix (args[1], 1.0, "ones"); |
|
1104 break; |
|
1105 case 3: |
|
1106 retval = new tree_constant [2]; |
|
1107 retval[0] = fill_matrix (args[1], args[2], 1.0, "ones"); |
|
1108 break; |
|
1109 default: |
|
1110 usage ("ones (n)\n ones (A)\n ones (n, m)"); |
|
1111 break; |
|
1112 } |
|
1113 return retval; |
|
1114 } |
|
1115 |
|
1116 /* |
|
1117 * You guessed it. |
|
1118 */ |
|
1119 tree_constant * |
|
1120 builtin_pause (tree_constant *args, int nargin, int nargout) |
|
1121 { |
|
1122 if (! (nargin == 1 || nargin == 2)) |
|
1123 { |
|
1124 usage ("pause ([delay])"); |
|
1125 return NULL_TREE_CONST; |
|
1126 } |
|
1127 |
|
1128 if (interactive) |
|
1129 { |
|
1130 if (nargin == 2) |
|
1131 sleep (NINT (args[1].double_value ())); |
|
1132 else if (kbhit () == EOF) |
|
1133 clean_up_and_exit (0); |
|
1134 } |
|
1135 return NULL_TREE_CONST; |
|
1136 } |
|
1137 |
|
1138 /* |
|
1139 * Delete turds from /tmp. |
|
1140 */ |
|
1141 tree_constant * |
|
1142 builtin_purge_tmp_files (tree_constant *, int, int) |
|
1143 { |
|
1144 cleanup_tmp_files (); |
|
1145 return NULL_TREE_CONST; |
|
1146 } |
|
1147 |
|
1148 /* |
|
1149 * Formatted printing. |
|
1150 */ |
|
1151 tree_constant * |
|
1152 builtin_printf (tree_constant *args, int nargin, int nargout) |
|
1153 { |
|
1154 tree_constant *retval = NULL_TREE_CONST; |
|
1155 if (nargin < 2) |
|
1156 usage ("printf (\"fmt\" [, ...])"); |
|
1157 else |
|
1158 retval = do_printf ("printf", args, nargin, nargout); |
|
1159 return retval; |
|
1160 } |
|
1161 |
|
1162 /* |
|
1163 * Product. |
|
1164 */ |
|
1165 tree_constant * |
|
1166 builtin_prod (tree_constant *args, int nargin, int nargout) |
|
1167 { |
|
1168 tree_constant *retval = NULL_TREE_CONST; |
|
1169 if (nargin != 2) |
|
1170 usage ("prod (M)"); |
|
1171 else |
|
1172 { |
|
1173 if (args != NULL_TREE_CONST && args[1].is_defined ()) |
|
1174 { |
|
1175 retval = new tree_constant [2]; |
|
1176 retval[0] = args[1].prod (); |
|
1177 } |
|
1178 } |
|
1179 return retval; |
|
1180 } |
|
1181 |
|
1182 /* |
|
1183 * Print name of current working directory. |
|
1184 */ |
|
1185 tree_constant * |
|
1186 builtin_pwd (tree_constant *args, int nargin, int nargout) |
|
1187 { |
|
1188 tree_constant *retval = NULL_TREE_CONST; |
|
1189 char *directory; |
|
1190 |
|
1191 if (verbatim_pwd) |
|
1192 { |
|
1193 char *buffer = new char [MAXPATHLEN]; |
|
1194 directory = getcwd (buffer, MAXPATHLEN); |
|
1195 |
|
1196 if (!directory) |
|
1197 { |
|
1198 message ("pwd", "can't find working directory!"); |
|
1199 delete buffer; |
|
1200 } |
|
1201 } |
|
1202 else |
|
1203 { |
|
1204 directory = get_working_directory ("pwd"); |
|
1205 } |
|
1206 |
|
1207 if (directory) |
|
1208 { |
|
1209 char *s = strconcat (directory, "\n"); |
|
1210 retval = new tree_constant [2]; |
|
1211 retval[0] = tree_constant (s); |
|
1212 delete [] s; |
|
1213 } |
|
1214 return retval; |
|
1215 } |
|
1216 |
|
1217 /* |
|
1218 * QPs. |
|
1219 */ |
|
1220 static void |
|
1221 qpsol_usage (void) |
|
1222 { |
|
1223 #if defined (QPSOL_MISSING) |
|
1224 message ("qpsol", "this function requires QPSOL, which is not freely\n\ |
|
1225 redistributable. For more information, read the file\n\ |
|
1226 libcruft/qpsol/README.MISSING in the source distribution."); |
|
1227 #else |
|
1228 usage ("[x, obj, info, lambda] = qpsol (x, H, c [, lb, ub] [, lb, A, ub])\n\ |
|
1229 \n\ |
|
1230 Groups of arguments surrounded in `[]' are optional, but\n\ |
|
1231 must appear in the same relative order shown above."); |
|
1232 #endif |
|
1233 } |
|
1234 |
|
1235 tree_constant * |
|
1236 builtin_qpsol (tree_constant *args, int nargin, int nargout) |
|
1237 { |
|
1238 tree_constant *retval = NULL_TREE_CONST; |
|
1239 |
|
1240 #if defined (QPSOL_MISSING) |
|
1241 // Force a bad value of inform, and empty matrices for x, phi, and lambda. |
|
1242 retval = new tree_constant [5]; |
|
1243 Matrix m; |
|
1244 retval[0] = tree_constant (m); |
|
1245 retval[1] = tree_constant (m); |
|
1246 retval[2] = tree_constant (-1.0); |
|
1247 retval[3] = tree_constant (m); |
|
1248 qpsol_usage (); |
|
1249 #else |
|
1250 if ((nargin == 4 || nargin == 6 || nargin == 7 || nargin == 9) |
|
1251 && (nargout >= 1 && nargout <= 4)) |
|
1252 DLD_BUILTIN (args, nargin, nargout, qpsol, |
|
1253 retval = qpsol (args, nargin, nargout);) |
|
1254 else |
|
1255 qpsol_usage (); |
|
1256 #endif |
|
1257 |
|
1258 return retval; |
|
1259 } |
|
1260 |
|
1261 /* |
|
1262 * QR factorization. |
|
1263 */ |
|
1264 tree_constant * |
|
1265 builtin_qr (tree_constant *args, int nargin, int nargout) |
|
1266 { |
|
1267 tree_constant *retval = NULL_TREE_CONST; |
|
1268 |
|
1269 if (nargin == 2 && nargout < 3) |
|
1270 DLD_BUILTIN (args, nargin, nargout, qr, |
|
1271 retval = qr (args[1], nargout);) |
|
1272 else |
|
1273 usage ("[Q, R] = qr (A)"); |
|
1274 |
|
1275 return retval; |
|
1276 } |
|
1277 |
|
1278 /* |
|
1279 * Random numbers. |
|
1280 */ |
|
1281 tree_constant * |
|
1282 builtin_quad (tree_constant *args, int nargin, int nargout) |
|
1283 { |
|
1284 tree_constant *retval = NULL_TREE_CONST; |
|
1285 |
|
1286 if ((nargin > 3 && nargin < 7) && (nargout > 0 && nargout < 5)) |
|
1287 DLD_BUILTIN (args, nargin, nargout, quad, |
|
1288 retval = do_quad (args, nargin, nargout);) |
|
1289 else |
|
1290 usage ("[v, ier, nfun, err] = quad (\"f\", a, b)\n\ |
|
1291 = quad (\"f\", a, b, tol)\n\ |
|
1292 = quad (\"f\", a, b, tol, sing)"); |
|
1293 |
|
1294 return retval; |
|
1295 } |
|
1296 |
|
1297 /* |
|
1298 * I'm outta here. |
|
1299 */ |
|
1300 tree_constant * |
|
1301 builtin_quit (tree_constant *args, int nargin, int nargout) |
|
1302 { |
|
1303 quitting_gracefully = 1; |
|
1304 clean_up_and_exit (0); |
|
1305 return NULL_TREE_CONST; |
|
1306 } |
|
1307 |
|
1308 /* |
|
1309 * Random numbers. |
|
1310 */ |
|
1311 tree_constant * |
|
1312 builtin_rand (tree_constant *args, int nargin, int nargout) |
|
1313 { |
|
1314 tree_constant *retval = NULL_TREE_CONST; |
|
1315 |
|
1316 if ((nargin > 0 && nargin < 4) && nargout == 1) |
|
1317 DLD_BUILTIN (args, nargin, nargout, rand, |
|
1318 retval = rand_internal (args, nargin, nargout);) |
|
1319 else |
|
1320 usage ("rand -- generate a random value\n\ |
|
1321 rand (n) -- generate N x N matrix\n\ |
|
1322 rand (A) -- generate matrix the size of A\n\ |
|
1323 rand (n, m) -- generate N x M matrix\n\ |
|
1324 rand (\"dist\") -- get current distribution\n\ |
|
1325 rand (\"distribution\") -- set distribution\n\ |
|
1326 rand (\"seed\") -- get current seed\n\ |
|
1327 rand (\"seed\", n) -- set seed"); |
|
1328 |
|
1329 return retval; |
|
1330 } |
|
1331 |
|
1332 /* |
|
1333 * Replot current plot. |
|
1334 */ |
|
1335 tree_constant * |
|
1336 builtin_replot (tree_constant *args, int nargin, int nargout) |
|
1337 { |
|
1338 tree_constant *retval = NULL_TREE_CONST; |
|
1339 |
|
1340 if (nargin > 1) |
|
1341 warning ("replot: ignoring extra arguments"); |
|
1342 |
|
1343 send_to_plot_stream ("replot\n"); |
|
1344 |
|
1345 return retval; |
|
1346 } |
|
1347 |
|
1348 /* |
|
1349 * Formatted reading. |
|
1350 */ |
|
1351 tree_constant * |
|
1352 builtin_scanf (tree_constant *args, int nargin, int nargout) |
|
1353 { |
|
1354 tree_constant *retval = NULL_TREE_CONST; |
|
1355 if (nargin != 2) |
|
1356 usage ("[...] = scanf (\"fmt\")"); |
|
1357 else |
|
1358 retval = do_scanf ("scanf", args, nargin, nargout); |
|
1359 return retval; |
|
1360 } |
|
1361 |
|
1362 /* |
|
1363 * Convert a vector to a string. |
|
1364 */ |
|
1365 tree_constant * |
|
1366 builtin_setstr (tree_constant *args, int nargin, int nargout) |
|
1367 { |
|
1368 tree_constant *retval = NULL_TREE_CONST; |
|
1369 |
|
1370 if (nargin == 2) |
|
1371 { |
|
1372 retval = new tree_constant [2]; |
|
1373 retval[0] = args[1].convert_to_str (); |
|
1374 } |
|
1375 else |
|
1376 usage ("setstr (v)"); |
|
1377 |
|
1378 return retval; |
|
1379 } |
|
1380 |
|
1381 /* |
|
1382 * Execute a shell command. |
|
1383 */ |
|
1384 tree_constant * |
|
1385 builtin_shell_command (tree_constant *args, int nargin, int nargout) |
|
1386 { |
|
1387 tree_constant *retval = NULL_TREE_CONST; |
|
1388 |
|
1389 if (nargin == 2 || nargin == 3) |
|
1390 { |
|
1391 if (args[1].is_string_type ()) |
|
1392 { |
|
1393 iprocstream cmd (args[1].string_value ()); |
|
1394 char ch; |
|
1395 ostrstream output_buf; |
|
1396 while (cmd.get (ch)) |
|
1397 output_buf.put (ch); |
|
1398 |
|
1399 output_buf << ends; |
|
1400 if (nargin == 2) |
|
1401 { |
|
1402 maybe_page_output (output_buf); |
|
1403 } |
|
1404 else |
|
1405 { |
|
1406 retval = new tree_constant [2]; |
|
1407 retval[0] = tree_constant (output_buf.str ()); |
|
1408 } |
|
1409 } |
|
1410 else |
|
1411 error ("shell_cmd: first argument must be a string"); |
|
1412 } |
|
1413 else |
|
1414 usage ("shell_cmd (string [, return_output])"); |
|
1415 |
|
1416 return retval; |
|
1417 } |
|
1418 |
|
1419 /* |
|
1420 * Report rows and columns. |
|
1421 */ |
|
1422 tree_constant * |
|
1423 builtin_size (tree_constant *args, int nargin, int nargout) |
|
1424 { |
|
1425 tree_constant *retval = NULL_TREE_CONST; |
|
1426 |
|
1427 if (nargin != 2) |
|
1428 usage ("size (x)"); |
|
1429 else |
|
1430 { |
|
1431 if (args != NULL_TREE_CONST && args[1].is_defined ()) |
|
1432 { |
|
1433 int nr = args[1].rows (); |
|
1434 int nc = args[1].columns (); |
|
1435 if (nargout == 1) |
|
1436 { |
|
1437 Matrix m (1, 2); |
|
1438 m.elem (0, 0) = nr; |
|
1439 m.elem (0, 1) = nc; |
|
1440 retval = new tree_constant [2]; |
|
1441 retval[0] = tree_constant (m); |
|
1442 } |
|
1443 else if (nargout == 2) |
|
1444 { |
|
1445 retval = new tree_constant [3]; |
|
1446 retval[0] = tree_constant ((double) nr); |
|
1447 retval[1] = tree_constant ((double) nc); |
|
1448 } |
|
1449 else |
|
1450 usage ("[n, m] = size (A)\n size (A)"); |
|
1451 } |
|
1452 } |
|
1453 return retval; |
|
1454 } |
|
1455 |
|
1456 /* |
|
1457 * Sort columns. |
|
1458 */ |
|
1459 tree_constant * |
|
1460 builtin_sort (tree_constant *args, int nargin, int nargout) |
|
1461 { |
|
1462 tree_constant *retval = NULL_TREE_CONST; |
|
1463 |
|
1464 if (nargin == 2) |
|
1465 retval = sort (args, nargin, nargout); |
|
1466 else |
|
1467 usage ("[s, i] = sort (x)"); |
|
1468 |
|
1469 return retval; |
|
1470 } |
|
1471 |
|
1472 /* |
|
1473 * Formatted printing to a string. |
|
1474 */ |
|
1475 tree_constant * |
|
1476 builtin_sprintf (tree_constant *args, int nargin, int nargout) |
|
1477 { |
|
1478 tree_constant *retval = NULL_TREE_CONST; |
|
1479 if (nargin < 2) |
|
1480 usage ("string = sprintf (\"fmt\" [, ...])"); |
|
1481 else |
|
1482 retval = do_printf ("sprintf", args, nargin, nargout); |
|
1483 return retval; |
|
1484 } |
|
1485 |
|
1486 /* |
|
1487 * Matrix sqrt. |
|
1488 */ |
|
1489 tree_constant * |
|
1490 builtin_sqrtm (tree_constant *args, int nargin, int nargout) |
|
1491 { |
|
1492 tree_constant *retval = NULL_TREE_CONST; |
|
1493 |
|
1494 if (nargin == 2) |
|
1495 retval = matrix_sqrt (args[1]); |
|
1496 else |
|
1497 usage ("sqrtm (A)"); |
|
1498 |
|
1499 return retval; |
|
1500 } |
|
1501 |
|
1502 /* |
|
1503 * Formatted reading from a string. |
|
1504 */ |
|
1505 tree_constant * |
|
1506 builtin_sscanf (tree_constant *args, int nargin, int nargout) |
|
1507 { |
|
1508 tree_constant *retval = NULL_TREE_CONST; |
|
1509 if (nargin != 3) |
|
1510 usage ("[...] = sscanf (string, \"fmt\")"); |
|
1511 else |
|
1512 retval = do_scanf ("sscanf", args, nargin, nargout); |
|
1513 return retval; |
|
1514 } |
|
1515 |
|
1516 /* |
|
1517 * Sum. |
|
1518 */ |
|
1519 tree_constant * |
|
1520 builtin_sum (tree_constant *args, int nargin, int nargout) |
|
1521 { |
|
1522 tree_constant *retval = NULL_TREE_CONST; |
|
1523 if (nargin != 2) |
|
1524 usage ("sum (M)"); |
|
1525 else |
|
1526 { |
|
1527 if (args != NULL_TREE_CONST && args[1].is_defined ()) |
|
1528 { |
|
1529 retval = new tree_constant [2]; |
|
1530 retval[0] = args[1].sum (); |
|
1531 } |
|
1532 } |
|
1533 return retval; |
|
1534 } |
|
1535 |
|
1536 /* |
|
1537 * Sum of squares. |
|
1538 */ |
|
1539 tree_constant * |
|
1540 builtin_sumsq (tree_constant *args, int nargin, int nargout) |
|
1541 { |
|
1542 tree_constant *retval = NULL_TREE_CONST; |
|
1543 if (nargin != 2) |
|
1544 usage ("sumsq (M)"); |
|
1545 else |
|
1546 { |
|
1547 if (args != NULL_TREE_CONST && args[1].is_defined ()) |
|
1548 { |
|
1549 retval = new tree_constant [2]; |
|
1550 retval[0] = args[1].sumsq (); |
|
1551 } |
|
1552 } |
|
1553 return retval; |
|
1554 } |
|
1555 |
|
1556 /* |
|
1557 * Singluar value decomposition. |
|
1558 */ |
|
1559 tree_constant * |
|
1560 builtin_svd (tree_constant *args, int nargin, int nargout) |
|
1561 { |
|
1562 tree_constant *retval = NULL_TREE_CONST; |
|
1563 |
|
1564 if (nargin == 2 && (nargout == 1 || nargout == 3)) |
|
1565 DLD_BUILTIN (args, nargin, nargout, svd, |
|
1566 retval = svd (args, nargin, nargout);) |
|
1567 else |
|
1568 usage ("[U, S, V] = svd (A)\n S = svd (A)"); |
|
1569 |
|
1570 return retval; |
|
1571 } |
|
1572 |
|
1573 /* |
|
1574 * Schur Decomposition |
|
1575 */ |
|
1576 tree_constant * |
|
1577 builtin_schur (tree_constant *args, int nargin, int nargout) |
|
1578 { |
|
1579 tree_constant *retval = NULL_TREE_CONST; |
|
1580 |
|
1581 if ((nargin == 3 || nargin == 2) && (nargout == 1 || nargout == 2)) |
|
1582 DLD_BUILTIN (args, nargin, nargout, hess, |
|
1583 retval = schur (args, nargin, nargout);) |
|
1584 else |
|
1585 usage ("[U, S] = schur (A)\n\ |
|
1586 S = schur (A)\n\n\ |
|
1587 or, for ordered Schur:\n\n\ |
|
1588 [U, S] = schur (A, \"A, D, or U\")\n\ |
|
1589 S = schur (A, \"A, D, or U\")\n\ |
|
1590 where:\n\n\ |
|
1591 A = continuous time poles\n\ |
|
1592 D = discrete time poles\n\ |
|
1593 U = unordered schur (default)"); |
|
1594 |
|
1595 return retval; |
|
1596 } |
|
1597 |
|
1598 /* |
|
1599 * Hessenburg Decomposition |
|
1600 */ |
|
1601 tree_constant * |
|
1602 builtin_hess (tree_constant *args, int nargin, int nargout) |
|
1603 { |
|
1604 tree_constant *retval = NULL_TREE_CONST; |
|
1605 |
|
1606 if (nargin == 2 && (nargout == 1 || nargout == 2)) |
|
1607 DLD_BUILTIN (args, nargin, nargout, hess, |
|
1608 retval = hess (args, nargin, nargout);) |
|
1609 else |
|
1610 usage ("[P, H] = hess (A)\n H = hess (A)"); |
|
1611 |
|
1612 return retval; |
|
1613 } |
|
1614 |
|
1615 /* |
|
1616 * Copying information. |
|
1617 */ |
|
1618 tree_constant * |
|
1619 builtin_warranty (tree_constant *args, int nargin, int nargout) |
|
1620 { |
|
1621 ostrstream output_buf; |
|
1622 output_buf << "\n Octave, version " << version_string |
|
1623 << ". Copyright (C) 1992, 1993, John W. Eaton\n" |
|
1624 << "\n\ |
|
1625 This program is free software; you can redistribute it and/or modify\n\ |
|
1626 it under the terms of the GNU General Public License as published by\n\ |
|
1627 the Free Software Foundation; either version 2 of the License, or\n\ |
|
1628 (at your option) any later version.\n\n\ |
|
1629 This program is distributed in the hope that it will be useful,\n\ |
|
1630 but WITHOUT ANY WARRANTY; without even the implied warranty of\n\ |
|
1631 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the\n\ |
|
1632 GNU General Public License for more details.\n\n\ |
|
1633 You should have received a copy of the GNU General Public License\n\ |
|
1634 along with this program. If not, write to the Free Software\n\ |
|
1635 Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.\n\n"; |
|
1636 |
|
1637 output_buf << ends; |
|
1638 maybe_page_output (output_buf); |
|
1639 |
|
1640 return NULL_TREE_CONST; |
|
1641 } |
|
1642 |
|
1643 /* |
|
1644 * A matrix of zeros. |
|
1645 */ |
|
1646 tree_constant * |
|
1647 builtin_zeros (tree_constant *args, int nargin, int nargout) |
|
1648 { |
|
1649 tree_constant *retval = NULL_TREE_CONST; |
|
1650 |
|
1651 switch (nargin) |
|
1652 { |
|
1653 case 2: |
|
1654 retval = new tree_constant [2]; |
|
1655 retval[0] = fill_matrix (args[1], 0.0, "zeros"); |
|
1656 break; |
|
1657 case 3: |
|
1658 retval = new tree_constant [2]; |
|
1659 retval[0] = fill_matrix (args[1], args[2], 0.0, "zeros"); |
|
1660 break; |
|
1661 default: |
|
1662 usage ("zeros (n)\n zeros (A)\n zeros (n, m)"); |
|
1663 break; |
|
1664 } |
|
1665 return retval; |
|
1666 } |
|
1667 |
|
1668 /* |
|
1669 ;;; Local Variables: *** |
|
1670 ;;; mode: C++ *** |
|
1671 ;;; page-delimiter: "^/\\*" *** |
|
1672 ;;; End: *** |
|
1673 */ |