GNU Octave 10.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
pr-output.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 1993-2025 The Octave Project Developers
4//
5// See the file COPYRIGHT.md in the top-level directory of this
6// distribution or <https://octave.org/copyright/>.
7//
8// This file is part of Octave.
9//
10// Octave is free software: you can redistribute it and/or modify it
11// under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// Octave is distributed in the hope that it will be useful, but
16// WITHOUT ANY WARRANTY; without even the implied warranty of
17// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18// GNU General Public License for more details.
19//
20// You should have received a copy of the GNU General Public License
21// along with Octave; see the file COPYING. If not, see
22// <https://www.gnu.org/licenses/>.
23//
24////////////////////////////////////////////////////////////////////////
25
26#if defined (HAVE_CONFIG_H)
27# include "config.h"
28#endif
29
30#include <cmath>
31
32#include <iomanip>
33#include <limits>
34#include <list>
35#include <sstream>
36#include <string>
37
38#include "Array-util.h"
39#include "CMatrix.h"
40#include "Range.h"
41#include "cmd-edit.h"
42#include "dMatrix.h"
43#include "lo-mappers.h"
44#include "mach-info.h"
45#include "oct-cmplx.h"
46#include "oct-string.h"
47#include "quit.h"
48
49#include "Cell.h"
50#include "defun.h"
51#include "error.h"
52#include "errwarn.h"
53#include "interpreter.h"
54#include "ovl.h"
55#include "oct-stream.h"
57#include "pager.h"
58#include "pr-flt-fmt.h"
59#include "pr-output.h"
60#include "sysdep.h"
61#include "unwind-prot.h"
62#include "utils.h"
63#include "variables.h"
64
65// TRUE means use a scaled fixed point format for 'format long' and
66// 'format short'.
67static bool Vfixed_point_format = false;
68
69// TRUE means that the dimensions of empty objects should be printed
70// like this: x = [](2x0).
72
73// TRUE means that the rows of big matrices should be split into
74// smaller slices that fit on the screen.
75static bool Vsplit_long_rows = true;
76
77// TRUE means don't do any fancy formatting.
78static bool free_format = false;
79
80// TRUE means print plus sign for nonzero, blank for zero.
81static bool plus_format = false;
82
83// First char for > 0, second for < 0, third for == 0.
84static std::string plus_format_chars = "+- ";
85
86// TRUE means always print in a rational approximation
87static bool rat_format = false;
88
89// Used to force the length of the rational approximation string for Frats
90static int rat_string_len = -1;
91
92// TRUE means always print like dollars and cents.
93static bool bank_format = false;
94
95// TRUE means print data in hexadecimal format.
96static int hex_format = 0;
97
98// TRUE means print data in binary-bit-pattern format.
99static int bit_format = 0;
100
101// TRUE means don't put newlines around the column number headers.
102bool Vcompact_format = false;
103
104// TRUE means use an e format.
105static bool print_e = false;
106
107// TRUE means use a g format.
108static bool print_g = false;
109
110// TRUE means print uppercase E in exponent field and A-F in hex format.
111static bool uppercase_format = false;
112
113// TRUE means use an engineering format.
114static bool print_eng = false;
115
116static int
117calc_scale_exp (const int& x)
118{
119 if (! print_eng)
120 return x;
121 else
122 return x - 3*static_cast<int> (x/3);
123
124 // The expression above is equivalent to x - (x % 3).
125
126 // According to the ISO specification for C++ the modulo operator is
127 // compiler dependent if any of the arguments are negative. Since
128 // this function will need to work on negative arguments, and we want
129 // to avoid portability issues, we re-implement the modulo function to
130 // the desired behavior (truncation). There may be a gnulib replacement.
131
132 // ISO/IEC 14882:2003 : Programming languages -- C++. 5.6.4: ISO,
133 // IEC. 2003 . "the binary % operator yields the remainder from the
134 // division of the first expression by the second. .... If both
135 // operands are nonnegative then the remainder is nonnegative; if not,
136 // the sign of the remainder is implementation-defined".
137}
138
139template <typename T>
140static inline int
141engineering_exponent (T x)
142{
143 int ex = 0;
144
145 if (x != 0)
146 {
147 T absval = (x < 0 ? -x : x);
148 int logabsval = static_cast<int> (std::floor (log10 (absval)));
149
150 // Avoid using modulo function with negative arguments for
151 // portability. See extended comment at calc_scale_exp
152
153 if (logabsval < 0)
154 ex = logabsval - 2 + ((-logabsval + 2) % 3);
155 else
156 ex = logabsval - (logabsval % 3);
157 }
158
159 return ex;
160}
161
162template <typename T>
163static inline int
164num_digits (T x)
165{
166 return 1 + (print_eng
167 ? engineering_exponent (x)
168 : static_cast<int> (std::floor (log10 (x))));
169}
170
171template <typename T>
172int
175 return engineering_exponent (m_val);
176}
177
178template <typename T>
179T
181{
182 return m_val / std::pow (static_cast<T> (10), exponent ());
183}
184
185template <typename T>
186std::ostream&
187operator << (std::ostream& os, const pr_engineering_float<T>& pef)
188{
189 octave::preserve_stream_state stream_state (os);
190
191 float_format real_fmt = pef.m_ff;
192
193 if (real_fmt.width () >= 0)
194 os << std::setw (real_fmt.width () - real_fmt.exponent_width ());
195
196 if (real_fmt.precision () >= 0)
197 os << std::setprecision (real_fmt.precision ());
198
199 os.flags (real_fmt.format_flags ());
200
201 os << pef.mantissa ();
202
203 int ex = pef.exponent ();
204 if (ex < 0)
205 {
206 if (uppercase_format)
207 os << std::setw (0) << "E-";
208 else
209 os << std::setw (0) << "e-";
210 ex = -ex;
211 }
212 else
213 {
214 if (uppercase_format)
215 os << std::setw (0) << "E+";
216 else
217 os << std::setw (0) << "e+";
218 }
219
220 os << std::setw (real_fmt.exponent_width () - 2)
221 << std::setfill ('0') << ex;
222
223 return os;
224}
225
226template <typename T>
227std::ostream&
228operator << (std::ostream& os, const pr_formatted_float<T>& pff)
229{
230 octave::preserve_stream_state stream_state (os);
231
232 float_format real_fmt = pff.m_ff;
233
234 if (real_fmt.width () >= 0)
235 os << std::setw (real_fmt.width ());
236
237 if (real_fmt.precision () >= 0)
238 os << std::setprecision (real_fmt.precision ());
239
240 os.flags (real_fmt.format_flags ());
241
242 os << pff.m_val;
243
244 return os;
245}
246
247template <typename T>
248std::ostream&
249operator << (std::ostream& os, const pr_rational_float<T>& prf)
250{
251 octave::preserve_stream_state stream_state (os);
252
253 float_format real_fmt = prf.m_ff;
254 bool have_neg_sign = prf.m_val < 0;
255
256 int fw = (rat_string_len > 0 ? rat_string_len : real_fmt.width ());
257 std::string s;
258
259 if (have_neg_sign)
260 s = rational_approx (prf.m_val, fw);
261 else
262 s = rational_approx (prf.m_val, fw-1);
263
264 if (fw >= 0)
265 os << std::setw (fw);
266
267 os.flags (real_fmt.format_flags ());
268
269 if (s == "0")
270 s = '*';
271 else if (fw > 0)
272 {
273 if (s.find ('/') != std::string::npos)
274 {
275 if (s.length () > (static_cast<unsigned int> (fw)))
276 s = '*';
277 }
278 else
279 {
280 if (have_neg_sign)
281 {
282 if (s.length () > (static_cast<unsigned int> (fw) - 2))
283 s = '*';
284 }
285 else
286 {
287 if (s.length () > (static_cast<unsigned int> (fw) - 3))
288 s = '*';
289 }
290 }
291 }
292
293 os << s;
294
295 return os;
296}
297
298template <typename T>
299static inline T
300pr_max_internal (const MArray<T>& m)
301{
302 // We expect a 2-d array.
303 panic_unless (m.ndims () == 2);
304
305 octave_idx_type nr = m.rows ();
306 octave_idx_type nc = m.columns ();
307
308 T result = std::numeric_limits<T>::lowest ();
309
310 bool all_inf_or_nan = true;
311
312 for (octave_idx_type j = 0; j < nc; j++)
313 for (octave_idx_type i = 0; i < nr; i++)
314 {
315 T val = m(i, j);
316 if (! octave::math::isfinite (val))
317 continue;
318
319 all_inf_or_nan = false;
320
321 if (val > result)
322 result = val;
323 }
324
325 if (all_inf_or_nan)
326 result = 0;
327
328 return result;
329}
330
331template <typename T>
332static inline T
333pr_min_internal (const MArray<T>& m)
334{
335 octave_idx_type nr = m.rows ();
336 octave_idx_type nc = m.columns ();
337
338 T result = std::numeric_limits<T>::max ();
339
340 bool all_inf_or_nan = true;
341
342 for (octave_idx_type j = 0; j < nc; j++)
343 for (octave_idx_type i = 0; i < nr; i++)
344 {
345 T val = m(i, j);
346 if (! octave::math::isfinite (val))
347 continue;
348
349 all_inf_or_nan = false;
350
351 if (val < result)
352 result = val;
353 }
354
355 if (all_inf_or_nan)
356 result = 0;
357
358 return result;
359}
360
361template <typename>
362struct pr_output_traits
363{
364 static const int DIGITS10;
365 static const int MAX_FIELD_WIDTH;
366};
367
368template <>
369struct pr_output_traits<double>
370{
371 static const int DIGITS10;
372 static const int MAX_FIELD_WIDTH;
373};
374
375const int pr_output_traits<double>::DIGITS10 = 16;
376const int pr_output_traits<double>::MAX_FIELD_WIDTH = 21;
377
378template <>
379struct pr_output_traits<float>
380{
381 static const int DIGITS10;
382 static const int MAX_FIELD_WIDTH;
383};
384
385const int pr_output_traits<float>::DIGITS10 = 8;
386const int pr_output_traits<float>::MAX_FIELD_WIDTH = 13;
387
388// FIXME: it would be nice to share more code among these functions,..
389
390// Works for double and float.
391
392template <typename T>
393static inline float_display_format
394make_real_format (int digits, bool inf_or_nan, bool int_only)
395{
396 float_format fmt;
397
398 int prec = std::min (output_precision (), pr_output_traits<T>::DIGITS10);
399
400 int fw = 0, ld = 0, rd = 0;
401
402 if (rat_format)
403 {
404 fw = 0;
405 rd = 0;
406 }
407 else if (bank_format)
408 {
409 fw = (digits < 0 ? 4 : digits + 3);
410 if (inf_or_nan)
411 fw = 3;
412 rd = 2;
413 }
414 else if (hex_format)
415 {
416 fw = 2 * sizeof (T);
417 rd = 0;
418 }
419 else if (bit_format)
420 {
421 fw = 8 * sizeof (T);
422 rd = 0;
423 }
424 else if (inf_or_nan)
425 {
426 fw = 3;
427 }
428 else if (int_only)
429 {
430 fw = digits;
431 ld = digits;
432 rd = 0;
433 }
434 else
435 {
436 if (digits > 0)
437 {
438 ld = digits;
439 rd = (prec > digits ? prec - digits : prec);
440 }
441 else if (digits < 0)
442 {
443 ld = 1;
444 rd = (prec > digits ? prec - digits : prec);
445 }
446 else
447 {
448 ld = 1;
449 rd = (prec > digits ? prec - 1 : prec);
450 }
451
452 fw = ld + 1 + rd;
453 }
454
455 if (! (rat_format || bank_format || hex_format || bit_format)
456 && (print_e || print_g || print_eng
457 || ld + rd > pr_output_traits<T>::DIGITS10
458 || fw > pr_output_traits<T>::MAX_FIELD_WIDTH
459 || ld + rd > (1.5 * prec)))
460 {
461 if (print_g)
462 fmt = float_format (prec, prec);
463 else
464 {
465 // e+ddd
466 int ex = 5;
467
468 if (print_eng)
469 {
470 // -ddd.
471 fw = 1 + prec + ex;
472 if (inf_or_nan)
473 {
474 fw = 3;
475 ex = 0;
476 }
477 fmt = float_format (fw, ex, prec - 1, std::ios::fixed);
478 }
479 else
480 {
481 // -d.
482 fw = prec + ex;
483 if (inf_or_nan)
484 {
485 fw = 3;
486 ex = 0;
487 }
488 fmt = float_format (fw, ex, prec - 1, std::ios::scientific);
489 }
490 }
491 }
492 else if (! bank_format && (inf_or_nan || int_only))
493 fmt = float_format (fw, ld);
494 else
495 fmt = float_format (fw, rd, std::ios::fixed);
496
497 if (uppercase_format)
498 fmt.uppercase ();
499
500 return float_display_format (fmt);
501}
502
503// Works for double and float.
504
505template <typename T>
507make_scalar_format (const T& val)
508{
509 if (free_format)
510 return float_display_format ();
511
512 bool inf_or_nan = (octave::math::isinf (val) || octave::math::isnan (val));
513
514 bool int_only = (! inf_or_nan && octave::math::x_nint (val) == val);
515
516 T val_abs = (val < 0 ? -val : val);
517
518 int digits = (inf_or_nan || val_abs == 0) ? 0 : num_digits (val_abs);
519
520 return make_real_format<T> (digits, inf_or_nan, int_only);
521}
522
523template <>
525make_format (const double& d)
526{
527 return make_scalar_format (d);
528}
529
530template <>
532make_format (const float& f)
533{
534 return make_scalar_format (f);
535}
536
537template <typename T>
538static inline float_display_format
539make_real_matrix_format (int x_max, int x_min, bool inf_or_nan,
540 int int_or_inf_or_nan)
541{
542 T scale = ((x_max == 0 || int_or_inf_or_nan)
543 ? 1 : std::pow (10.0, calc_scale_exp (x_max - 1)));
544
545 float_format fmt;
546
547 int prec = std::min (output_precision (), pr_output_traits<T>::DIGITS10);
548
549 int fw = 0, ld = 0, rd = 0;
550
551 if (rat_format)
552 {
553 fw = 9;
554 rd = 0;
555 }
556 else if (bank_format)
557 {
558 int digits = (x_max > x_min ? x_max : x_min);
559 fw = (digits <= 0 ? 5 : digits + 4);
560 rd = 2;
561 }
562 else if (hex_format)
563 {
564 fw = 2 * sizeof (T);
565 rd = 0;
566 }
567 else if (bit_format)
568 {
569 fw = 8 * sizeof (T);
570 rd = 0;
571 }
572 else if (Vfixed_point_format && ! print_g)
573 {
574 rd = prec - 1;
575 fw = rd + 3;
576 if (inf_or_nan && fw < 4)
577 fw = 4;
578 }
579 else if (int_or_inf_or_nan)
580 {
581 int digits = (x_max > x_min ? x_max : x_min);
582 fw = (digits <= 0 ? 2 : digits + 1);
583 if (inf_or_nan && fw < 4)
584 fw = 4;
585 rd = fw;
586 }
587 else
588 {
589 int ld_max, rd_max;
590 if (x_max > 0)
591 {
592 ld_max = x_max;
593 rd_max = (prec > x_max ? prec - x_max : prec);
594 x_max++;
595 }
596 else if (x_max < 0)
597 {
598 ld_max = 1;
599 rd_max = (prec > x_max ? prec - x_max : prec);
600 x_max = -x_max + 1;
601 }
602 else
603 {
604 ld_max = 1;
605 rd_max = (prec > 1 ? prec - 1 : prec);
606 x_max = 1;
607 }
608
609 int ld_min, rd_min;
610 if (x_min > 0)
611 {
612 ld_min = x_min;
613 rd_min = (prec > x_min ? prec - x_min : prec);
614 x_min++;
615 }
616 else if (x_min < 0)
617 {
618 ld_min = 1;
619 rd_min = (prec > x_min ? prec - x_min : prec);
620 x_min = -x_min + 1;
621 }
622 else
623 {
624 ld_min = 1;
625 rd_min = (prec > 1 ? prec - 1 : prec);
626 x_min = 1;
627 }
628
629 ld = (ld_max > ld_min ? ld_max : ld_min);
630 rd = (rd_max > rd_min ? rd_max : rd_min);
631
632 fw = 1 + ld + 1 + rd;
633 if (inf_or_nan && fw < 4)
634 fw = 4;
635 }
636
637 if (! (rat_format || bank_format || hex_format || bit_format)
638 && (print_e || print_eng || print_g
639 || (! Vfixed_point_format
640 && (ld + rd > pr_output_traits<T>::DIGITS10
641 || fw > pr_output_traits<T>::MAX_FIELD_WIDTH
642 || ld + rd > (1.5 * prec)))))
643 {
644 if (print_g)
645 fmt = float_format (prec+6, prec);
646 else
647 {
648 int ex = 4;
649 if (x_max > 100 || x_min > 100)
650 ex++;
651
652 if (print_eng)
653 {
654 fw = 4 + prec + ex;
655 if (inf_or_nan && fw < 6)
656 fw = 6;
657 fmt = float_format (fw, ex, prec - 1, std::ios::fixed);
658 }
659 else
660 {
661 fw = 2 + prec + ex;
662 if (inf_or_nan && fw < 4)
663 fw = 4;
664 fmt = float_format (fw, prec - 1, std::ios::scientific);
665 }
666 }
667 }
668 else if (! bank_format && int_or_inf_or_nan)
669 fmt = float_format (fw, rd);
670 else
671 fmt = float_format (fw, rd, std::ios::fixed);
672
673 if (uppercase_format)
674 fmt.uppercase ();
675
676 return float_display_format (scale, fmt);
677}
678
679template <typename MT>
680static inline float_display_format
681make_matrix_format (const MT& m)
682{
683 // We expect a 2-d array.
684 panic_unless (m.ndims () == 2);
685
686 if (free_format)
687 return float_display_format ();
688
689 bool inf_or_nan = m.any_element_is_inf_or_nan ();
690
691 bool int_or_inf_or_nan = m.all_elements_are_int_or_inf_or_nan ();
692
693 MT m_abs = m.abs ();
694
695 typedef typename MT::element_type ELT_T;
696
697 ELT_T max_abs = pr_max_internal (m_abs);
698 ELT_T min_abs = pr_min_internal (m_abs);
699
700 int x_max = ((inf_or_nan || max_abs == 0) ? 0 : num_digits (max_abs));
701
702 int x_min = ((inf_or_nan || min_abs == 0) ? 0 : num_digits (min_abs));
703
704 return make_real_matrix_format<ELT_T> (x_max, x_min, inf_or_nan,
705 int_or_inf_or_nan);
706}
707
708template <>
711{
712 return make_matrix_format (m);
713}
714
715template <>
718{
719 return make_matrix_format (m);
720}
721
722template <typename T>
723static inline float_display_format
724make_complex_format (int x_max, int x_min, int r_x,
725 bool inf_or_nan, int int_only)
726{
727 float_format r_fmt;
728 float_format i_fmt;
729
730 int prec = std::min (output_precision (), pr_output_traits<T>::DIGITS10);
731
732 int i_fw = 0, r_fw = 0, ld = 0, rd = 0;
733
734 if (rat_format)
735 {
736 i_fw = 0;
737 r_fw = 0;
738 rd = 0;
739 }
740 else if (bank_format)
741 {
742 int digits = r_x;
743 i_fw = 0;
744 r_fw = (digits <= 0 ? 5 : digits + 4);
745 rd = 2;
746 }
747 else if (hex_format)
748 {
749 r_fw = 2 * sizeof (T);
750 i_fw = 2 * sizeof (T);
751 rd = 0;
752 }
753 else if (bit_format)
754 {
755 r_fw = 8 * sizeof (T);
756 i_fw = 8 * sizeof (T);
757 rd = 0;
758 }
759 else if (int_only)
760 {
761 int digits = (x_max > x_min ? x_max : x_min);
762 i_fw = (digits <= 0 ? 1 : digits);
763 r_fw = i_fw + 1;
764 if (inf_or_nan && i_fw < 3)
765 {
766 i_fw = 3;
767 r_fw = 4;
768 }
769 ld = r_fw;
770 }
771 else // ordinary case of floating point numeric values
772 {
773 int ld_max, rd_max;
774 if (x_max > 0)
775 {
776 ld_max = x_max;
777 rd_max = (prec > x_max ? prec - x_max : prec);
778 x_max++;
779 }
780 else if (x_max < 0)
781 {
782 ld_max = 1;
783 rd_max = (prec > x_max ? prec - x_max : prec);
784 x_max = -x_max + 1;
785 }
786 else
787 {
788 ld_max = 1;
789 rd_max = (prec > 1 ? prec - 1 : prec);
790 x_max = 1;
791 }
792
793 int ld_min, rd_min;
794 if (x_min > 0)
795 {
796 ld_min = x_min;
797 rd_min = (prec > x_min ? prec - x_min : prec);
798 x_min++;
799 }
800 else if (x_min < 0)
801 {
802 ld_min = 1;
803 rd_min = (prec > x_min ? prec - x_min : prec);
804 x_min = -x_min + 1;
805 }
806 else
807 {
808 ld_min = 1;
809 rd_min = (prec > 1 ? prec - 1 : prec);
810 x_min = 1;
811 }
812
813 ld = (ld_max > ld_min ? ld_max : ld_min);
814 rd = (rd_max > rd_min ? rd_max : rd_min);
815
816 i_fw = ld + 1 + rd;
817 r_fw = i_fw + 1;
818 if (inf_or_nan && i_fw < 3)
819 {
820 i_fw = 3;
821 r_fw = 4;
822 }
823 }
824
825 if (! (rat_format || bank_format || hex_format || bit_format)
826 && (print_e || print_eng || print_g
827 || ld + rd > pr_output_traits<T>::DIGITS10
828 || r_fw > pr_output_traits<T>::MAX_FIELD_WIDTH
829 || i_fw > pr_output_traits<T>::MAX_FIELD_WIDTH
830 || ld + rd > (1.5 * prec)))
831 {
832 if (print_g)
833 {
834 int width = prec + 6;
835 r_fmt = float_format (width, prec);
836 i_fmt = float_format (width, prec);
837 }
838 else
839 {
840 int ex = 4;
841 if (x_max > 100 || x_min > 100)
842 ex++;
843
844 if (print_eng)
845 {
846 i_fw = 3 + prec + ex;
847 r_fw = i_fw + 1;
848 if (inf_or_nan && i_fw < 5)
849 {
850 i_fw = 5;
851 r_fw = 6;
852 }
853 r_fmt = float_format (r_fw, ex, prec - 1, std::ios::fixed);
854 i_fmt = float_format (i_fw, ex, prec - 1, std::ios::fixed);
855 }
856 else
857 {
858 i_fw = 1 + prec + ex;
859 r_fw = i_fw + 1;
860 if (inf_or_nan && i_fw < 3)
861 {
862 i_fw = 3;
863 r_fw = 4;
864 }
865 r_fmt = float_format (r_fw, prec - 1, std::ios::scientific);
866 i_fmt = float_format (i_fw, prec - 1, std::ios::scientific);
867 }
868 }
869
870 if (uppercase_format)
871 {
872 r_fmt.uppercase ();
873 i_fmt.uppercase ();
874 }
875 }
876 else if (! bank_format && int_only)
877 {
878 r_fmt = float_format (r_fw, ld);
879 i_fmt = float_format (i_fw, ld);
880 }
881 else
882 {
883 r_fmt = float_format (r_fw, rd, std::ios::fixed);
884 i_fmt = float_format (i_fw, rd, std::ios::fixed);
885 }
886
887 return float_display_format (r_fmt, i_fmt);
888}
889
890template <typename T>
892make_complex_scalar_format (const std::complex<T>& c)
893{
894 if (free_format)
895 return float_display_format ();
896
897 T rp = c.real ();
898 T ip = c.imag ();
899
900 bool r_inf_or_nan = (octave::math::isinf (rp) || octave::math::isnan (rp));
901 bool i_inf_or_nan = (octave::math::isinf (ip) || octave::math::isnan (ip));
902 bool inf_or_nan = r_inf_or_nan || i_inf_or_nan;
903
904 bool int_only = ((r_inf_or_nan || octave::math::x_nint (rp) == rp)
905 && (i_inf_or_nan || octave::math::x_nint (ip) == ip));
906
907 T r_abs = (rp < 0 ? -rp : rp);
908 T i_abs = (ip < 0 ? -ip : ip);
909
910 int r_x = ((r_inf_or_nan || r_abs == 0) ? 0 : num_digits (r_abs));
911 int i_x = ((i_inf_or_nan || i_abs == 0) ? 0 : num_digits (i_abs));
912
913 int x_max, x_min;
914
915 if (r_x > i_x)
916 {
917 x_max = r_x;
918 x_min = i_x;
919 }
920 else
921 {
922 x_max = i_x;
923 x_min = r_x;
924 }
925
926 return make_complex_format<T> (x_max, x_min, r_x, inf_or_nan, int_only);
927}
928
929template <>
931make_format (const std::complex<double>& c)
932{
934}
935
936template <>
938make_format (const std::complex<float>& fc)
939{
940 return make_complex_scalar_format (fc);
941}
942
943template <typename T>
944static inline float_display_format
945make_complex_matrix_format (int x_max, int x_min, int r_x_max,
946 int r_x_min, bool inf_or_nan,
947 int int_or_inf_or_nan)
948{
949 T scale = ((x_max == 0 || int_or_inf_or_nan)
950 ? 1 : std::pow (10.0, calc_scale_exp (x_max - 1)));
951
952 float_format r_fmt;
953 float_format i_fmt;
954
955 int prec = std::min (output_precision (), pr_output_traits<T>::DIGITS10);
956
957 int i_fw = 0, r_fw = 0, ld = 0, rd = 0;
958
959 if (rat_format)
960 {
961 i_fw = 9;
962 r_fw = 9;
963 rd = 0;
964 }
965 else if (bank_format)
966 {
967 int digits = (r_x_max > r_x_min ? r_x_max : r_x_min);
968 i_fw = 0;
969 r_fw = (digits <= 0 ? 5 : digits + 4);
970 rd = 2;
971 }
972 else if (hex_format)
973 {
974 r_fw = 2 * sizeof (T);
975 i_fw = 2 * sizeof (T);
976 rd = 0;
977 }
978 else if (bit_format)
979 {
980 r_fw = 8 * sizeof (T);
981 i_fw = 8 * sizeof (T);
982 rd = 0;
983 }
984 else if (Vfixed_point_format && ! print_g)
985 {
986 rd = prec - 1;
987 i_fw = rd + 1;
988 r_fw = i_fw + 2;
989 if (inf_or_nan && i_fw < 3)
990 {
991 i_fw = 3;
992 r_fw = 4;
993 }
994 }
995 else if (int_or_inf_or_nan)
996 {
997 int digits = (x_max > x_min ? x_max : x_min);
998 i_fw = (digits <= 0 ? 1 : digits);
999 r_fw = i_fw + 1;
1000 if (inf_or_nan && i_fw < 3)
1001 {
1002 i_fw = 3;
1003 r_fw = 4;
1004 }
1005 rd = r_fw;
1006 }
1007 else
1008 {
1009 int ld_max, rd_max;
1010 if (x_max > 0)
1011 {
1012 ld_max = x_max;
1013 rd_max = (prec > x_max ? prec - x_max : prec);
1014 x_max++;
1015 }
1016 else if (x_max < 0)
1017 {
1018 ld_max = 1;
1019 rd_max = (prec > x_max ? prec - x_max : prec);
1020 x_max = -x_max + 1;
1021 }
1022 else
1023 {
1024 ld_max = 1;
1025 rd_max = (prec > 1 ? prec - 1 : prec);
1026 x_max = 1;
1027 }
1028
1029 int ld_min, rd_min;
1030 if (x_min > 0)
1031 {
1032 ld_min = x_min;
1033 rd_min = (prec > x_min ? prec - x_min : prec);
1034 x_min++;
1035 }
1036 else if (x_min < 0)
1037 {
1038 ld_min = 1;
1039 rd_min = (prec > x_min ? prec - x_min : prec);
1040 x_min = -x_min + 1;
1041 }
1042 else
1043 {
1044 ld_min = 1;
1045 rd_min = (prec > 1 ? prec - 1 : prec);
1046 x_min = 1;
1047 }
1048
1049 ld = (ld_max > ld_min ? ld_max : ld_min);
1050 rd = (rd_max > rd_min ? rd_max : rd_min);
1051
1052 i_fw = ld + 1 + rd;
1053 r_fw = i_fw + 1;
1054 if (inf_or_nan && i_fw < 3)
1055 {
1056 i_fw = 3;
1057 r_fw = 4;
1058 }
1059 }
1060
1061 if (! (rat_format || bank_format || hex_format || bit_format)
1062 && (print_e || print_eng || print_g
1063 || (! Vfixed_point_format
1064 && (ld + rd > pr_output_traits<T>::DIGITS10
1065 || r_fw > pr_output_traits<T>::MAX_FIELD_WIDTH
1066 || i_fw > pr_output_traits<T>::MAX_FIELD_WIDTH
1067 || ld + rd > (1.5 * prec)))))
1068 {
1069 if (print_g)
1070 {
1071 int width = prec + 6;
1072 r_fmt = float_format (width, prec);
1073 i_fmt = float_format (width, prec);
1074 }
1075 else
1076 {
1077 int ex = 4;
1078 if (x_max > 100 || x_min > 100)
1079 ex++;
1080
1081 if (print_eng)
1082 {
1083 i_fw = 3 + prec + ex;
1084 r_fw = i_fw + 1;
1085 if (inf_or_nan && i_fw < 5)
1086 {
1087 i_fw = 5;
1088 r_fw = 6;
1089 }
1090 r_fmt = float_format (r_fw, ex, prec - 1, std::ios::fixed);
1091 i_fmt = float_format (i_fw, ex, prec - 1, std::ios::fixed);
1092 }
1093 else
1094 {
1095 i_fw = 1 + prec + ex;
1096 r_fw = i_fw + 1;
1097 if (inf_or_nan && i_fw < 3)
1098 {
1099 i_fw = 3;
1100 r_fw = 4;
1101 }
1102 r_fmt = float_format (r_fw, prec - 1, std::ios::scientific);
1103 i_fmt = float_format (i_fw, prec - 1, std::ios::scientific);
1104 }
1105 }
1106
1107 if (uppercase_format)
1108 {
1109 r_fmt.uppercase ();
1110 i_fmt.uppercase ();
1111 }
1112 }
1113 else if (! bank_format && int_or_inf_or_nan)
1114 {
1115 r_fmt = float_format (r_fw, rd);
1116 i_fmt = float_format (i_fw, rd);
1117 }
1118 else
1119 {
1120 r_fmt = float_format (r_fw, rd, std::ios::fixed);
1121 i_fmt = float_format (i_fw, rd, std::ios::fixed);
1122 }
1123
1124 return float_display_format (scale, r_fmt, i_fmt);
1125}
1126
1127template <typename CMT>
1128static inline float_display_format
1129make_complex_matrix_format (const CMT& cm)
1130{
1131 if (free_format)
1132 return float_display_format ();
1133
1134 typedef typename CMT::real_matrix_type RMT;
1135 typedef typename CMT::real_elt_type ELT_T;
1136
1137 RMT rp = real (cm);
1138 RMT ip = imag (cm);
1139
1140 bool inf_or_nan = cm.any_element_is_inf_or_nan ();
1141
1142 bool int_or_inf_or_nan = (rp.all_elements_are_int_or_inf_or_nan ()
1143 && ip.all_elements_are_int_or_inf_or_nan ());
1144
1145 RMT r_m_abs = rp.abs ();
1146 ELT_T r_max_abs = pr_max_internal (r_m_abs);
1147 ELT_T r_min_abs = pr_min_internal (r_m_abs);
1148
1149 RMT i_m_abs = ip.abs ();
1150 ELT_T i_max_abs = pr_max_internal (i_m_abs);
1151 ELT_T i_min_abs = pr_min_internal (i_m_abs);
1152
1153 int r_x_max = ((r_max_abs == 0 || octave::math::isinf (r_max_abs)
1154 || octave::math::isnan (r_max_abs))
1155 ? 0 : num_digits (r_max_abs));
1156
1157 int r_x_min = ((r_min_abs == 0 || octave::math::isinf (r_min_abs)
1158 || octave::math::isnan (r_min_abs))
1159 ? 0 : num_digits (r_min_abs));
1160
1161 int i_x_max = ((i_max_abs == 0 || octave::math::isinf (i_max_abs)
1162 || octave::math::isnan (i_max_abs))
1163 ? 0 : num_digits (i_max_abs));
1164
1165 int i_x_min = ((i_min_abs == 0 || octave::math::isinf (i_min_abs)
1166 || octave::math::isnan (i_min_abs))
1167 ? 0 : num_digits (i_min_abs));
1168
1169 int x_max = (r_x_max > i_x_max ? r_x_max : i_x_max);
1170 int x_min = (r_x_min > i_x_min ? r_x_min : i_x_min);
1171
1172 return make_complex_matrix_format<ELT_T> (x_max, x_min, r_x_max, r_x_min,
1173 inf_or_nan, int_or_inf_or_nan);
1174}
1175
1176template <>
1179{
1180 return make_complex_matrix_format (cm);
1181}
1182
1183template <>
1186{
1187 return make_complex_matrix_format (cm);
1188}
1189
1190template <>
1193{
1194 return float_display_format (float_format (1, 1));
1195}
1196
1197template <typename T>
1198static inline float_display_format
1199make_range_format (int x_max, int x_min, int all_ints)
1200{
1201 double scale = ((x_max == 0 || all_ints)
1202 ? 1 : std::pow (10.0, calc_scale_exp (x_max - 1)));
1203
1204 float_format fmt;
1205
1206 int prec = std::min (output_precision (), pr_output_traits<T>::DIGITS10);
1207
1208 int fw = 0, ld = 0, rd = 0;
1209
1210 if (rat_format)
1211 {
1212 fw = 9;
1213 rd = 0;
1214 }
1215 else if (bank_format)
1216 {
1217 int digits = (x_max > x_min ? x_max : x_min);
1218 fw = (digits < 0 ? 5 : digits + 4);
1219 rd = 2;
1220 }
1221 else if (hex_format)
1222 {
1223 fw = 2 * sizeof (T);
1224 rd = 0;
1225 }
1226 else if (bit_format)
1227 {
1228 fw = 8 * sizeof (T);
1229 rd = 0;
1230 }
1231 else if (all_ints)
1232 {
1233 int digits = (x_max > x_min ? x_max : x_min);
1234 fw = digits + 1;
1235 rd = fw;
1236 }
1237 else if (Vfixed_point_format && ! print_g)
1238 {
1239 rd = prec - 1;
1240 fw = rd + 3;
1241 }
1242 else
1243 {
1244 int ld_max, rd_max;
1245 if (x_max > 0)
1246 {
1247 ld_max = x_max;
1248 rd_max = (prec > x_max ? prec - x_max : prec);
1249 x_max++;
1250 }
1251 else if (x_max < 0)
1252 {
1253 ld_max = 1;
1254 rd_max = (prec > x_max ? prec - x_max : prec);
1255 x_max = -x_max + 1;
1256 }
1257 else
1258 {
1259 ld_max = 1;
1260 rd_max = (prec > 1 ? prec - 1 : prec);
1261 x_max = 1;
1262 }
1263
1264 int ld_min, rd_min;
1265 if (x_min > 0)
1266 {
1267 ld_min = x_min;
1268 rd_min = (prec > x_min ? prec - x_min : prec);
1269 x_min++;
1270 }
1271 else if (x_min < 0)
1272 {
1273 ld_min = 1;
1274 rd_min = (prec > x_min ? prec - x_min : prec);
1275 x_min = -x_min + 1;
1276 }
1277 else
1278 {
1279 ld_min = 1;
1280 rd_min = (prec > 1 ? prec - 1 : prec);
1281 x_min = 1;
1282 }
1283
1284 ld = (ld_max > ld_min ? ld_max : ld_min);
1285 rd = (rd_max > rd_min ? rd_max : rd_min);
1286
1287 fw = ld + rd + 3;
1288 }
1289
1290 if (! (rat_format || bank_format || hex_format || bit_format)
1291 && (print_e || print_eng || print_g
1292 || (! Vfixed_point_format
1293 && (ld + rd > pr_output_traits<T>::DIGITS10
1294 || fw > pr_output_traits<T>::MAX_FIELD_WIDTH
1295 || ld + rd > (1.5 * prec)))))
1296 {
1297 if (print_g)
1298 fmt = float_format (prec+6, prec);
1299 else
1300 {
1301 int ex = 4;
1302 if (x_max > 100 || x_min > 100)
1303 ex++;
1304
1305 if (print_eng)
1306 {
1307 fw = 5 + prec + ex;
1308 fmt = float_format (fw, ex, prec - 1, std::ios::fixed);
1309 }
1310 else
1311 {
1312 fw = 3 + prec + ex;
1313 fmt = float_format (fw, prec - 1, std::ios::scientific);
1314 }
1315 }
1316 }
1317 else if (! bank_format && all_ints)
1318 fmt = float_format (fw, rd);
1319 else
1320 fmt = float_format (fw, rd, std::ios::fixed);
1321
1322 if (uppercase_format)
1323 fmt.uppercase ();
1324
1325 return float_display_format (scale, fmt);
1326}
1327
1328template <>
1330make_format (const octave::range<double>& r)
1331{
1332 if (free_format)
1333 return float_display_format ();
1334
1335 double r_min = r.base ();
1336 double r_max = r.limit ();
1337
1338 if (r_max < r_min)
1339 {
1340 double tmp = r_max;
1341 r_max = r_min;
1342 r_min = tmp;
1343 }
1344
1345 bool all_ints = r.all_elements_are_ints ();
1346
1347 double max_abs = (r_max < 0 ? -r_max : r_max);
1348 double min_abs = (r_min < 0 ? -r_min : r_min);
1349
1350 int x_max = ((max_abs == 0 || octave::math::isinf (max_abs)
1351 || octave::math::isnan (max_abs))
1352 ? 0 : num_digits (max_abs));
1353
1354 int x_min = ((min_abs == 0 || octave::math::isinf (min_abs)
1355 || octave::math::isnan (min_abs))
1356 ? 0 : num_digits (min_abs));
1357
1358 return make_range_format<double> (x_max, x_min, all_ints);
1359}
1360
1361template <typename T>
1362union equiv
1363{
1364 T val;
1365 unsigned char i[sizeof (T)];
1366};
1367
1368#define PRINT_CHAR_BITS(os, c) \
1369 do \
1370 { \
1371 unsigned char ctmp = c; \
1372 char stmp[9]; \
1373 stmp[0] = (ctmp & 0x80) ? '1' : '0'; \
1374 stmp[1] = (ctmp & 0x40) ? '1' : '0'; \
1375 stmp[2] = (ctmp & 0x20) ? '1' : '0'; \
1376 stmp[3] = (ctmp & 0x10) ? '1' : '0'; \
1377 stmp[4] = (ctmp & 0x08) ? '1' : '0'; \
1378 stmp[5] = (ctmp & 0x04) ? '1' : '0'; \
1379 stmp[6] = (ctmp & 0x02) ? '1' : '0'; \
1380 stmp[7] = (ctmp & 0x01) ? '1' : '0'; \
1381 stmp[8] = '\0'; \
1382 os << stmp; \
1383 } \
1384 while (0)
1385
1386#define PRINT_CHAR_BITS_SWAPPED(os, c) \
1387 do \
1388 { \
1389 unsigned char ctmp = c; \
1390 char stmp[9]; \
1391 stmp[0] = (ctmp & 0x01) ? '1' : '0'; \
1392 stmp[1] = (ctmp & 0x02) ? '1' : '0'; \
1393 stmp[2] = (ctmp & 0x04) ? '1' : '0'; \
1394 stmp[3] = (ctmp & 0x08) ? '1' : '0'; \
1395 stmp[4] = (ctmp & 0x10) ? '1' : '0'; \
1396 stmp[5] = (ctmp & 0x20) ? '1' : '0'; \
1397 stmp[6] = (ctmp & 0x40) ? '1' : '0'; \
1398 stmp[7] = (ctmp & 0x80) ? '1' : '0'; \
1399 stmp[8] = '\0'; \
1400 os << stmp; \
1401 } \
1402 while (0)
1403
1404template <typename T>
1405static inline void
1406pr_any_float (std::ostream& os, const float_format& fmt, T val)
1407{
1408 // Unless explicitly asked for, always print in big-endian format
1409 // for hex and bit formats.
1410 //
1411 // {bit,hex}_format == 1: print big-endian
1412 // {bit,hex}_format == 2: print native
1413
1414 int fw = fmt.width ();
1415
1416 if (hex_format)
1417 {
1418 octave::preserve_stream_state stream_state (os);
1419
1420 equiv<T> tmp;
1421 tmp.val = val;
1422
1423 // Unless explicitly asked for, always print in big-endian format.
1424
1425 // FIXME: Will bad things happen if we are interrupted before resetting
1426 // the format flags and fill character?
1427
1428 octave::mach_info::float_format flt_fmt
1429 = octave::mach_info::native_float_format ();
1430
1431 os.fill ('0');
1432 if (uppercase_format)
1433 os.flags (std::ios::right | std::ios::hex | std::ios::uppercase);
1434 else
1435 os.flags (std::ios::right | std::ios::hex);
1436
1437 if (hex_format > 1
1438 || flt_fmt == octave::mach_info::flt_fmt_ieee_big_endian)
1439 {
1440 for (std::size_t i = 0; i < sizeof (T); i++)
1441 os << std::setw (2) << static_cast<int> (tmp.i[i]);
1442 }
1443 else
1444 {
1445 for (int i = sizeof (T) - 1; i >= 0; i--)
1446 os << std::setw (2) << static_cast<int> (tmp.i[i]);
1447 }
1448 }
1449 else if (bit_format)
1450 {
1451 equiv<T> tmp;
1452 tmp.val = val;
1453
1454 octave::mach_info::float_format flt_fmt
1455 = octave::mach_info::native_float_format ();
1456
1457 if (flt_fmt == octave::mach_info::flt_fmt_ieee_big_endian)
1458 {
1459 for (std::size_t i = 0; i < sizeof (T); i++)
1460 PRINT_CHAR_BITS (os, tmp.i[i]);
1461 }
1462 else
1463 {
1464 if (bit_format > 1)
1465 {
1466 for (std::size_t i = 0; i < sizeof (T); i++)
1467 PRINT_CHAR_BITS (os, tmp.i[i]);
1468 }
1469 else
1470 {
1471 for (int i = sizeof (T) - 1; i >= 0; i--)
1472 PRINT_CHAR_BITS (os, tmp.i[i]);
1473 }
1474 }
1475 }
1476 else if (val == 0)
1477 {
1478 octave::preserve_stream_state stream_state (os);
1479
1480 if (fw > 0)
1481 os << std::setw (fw) << "0";
1482 else
1483 os << "0";
1484 }
1485 else if (octave::math::isna (val))
1486 {
1487 octave::preserve_stream_state stream_state (os);
1488
1489 if (fw > 0)
1490 os << std::setw (fw) << "NA";
1491 else
1492 os << "NA";
1493 }
1494 else if (rat_format)
1495 os << pr_rational_float<T> (fmt, val);
1496 else if (octave::math::isinf (val))
1497 {
1498 octave::preserve_stream_state stream_state (os);
1499
1500 const char *s;
1501 if (val < 0)
1502 s = "-Inf";
1503 else
1504 s = "Inf";
1505
1506 if (fw > 0)
1507 os << std::setw (fw) << s;
1508 else
1509 os << s;
1510 }
1511 else if (octave::math::isnan (val))
1512 {
1513 octave::preserve_stream_state stream_state (os);
1514
1515 if (fw > 0)
1516 os << std::setw (fw) << "NaN";
1517 else
1518 os << "NaN";
1519 }
1520 else if (print_eng)
1521 os << pr_engineering_float<T> (fmt, val);
1522 else
1523 os << pr_formatted_float<T> (fmt, val);
1524}
1525
1526template <typename T>
1527static inline void
1528pr_float (std::ostream& os, const float_display_format& fmt, T val)
1529{
1530 double scale = fmt.scale_factor ();
1531
1532 if (Vfixed_point_format && ! (print_g || print_e) && scale != 1)
1533 val /= scale;
1534
1535 pr_any_float (os, fmt.real_format (), val);
1536}
1537
1538template <typename T>
1539static inline void
1540pr_imag_float (std::ostream& os, const float_display_format& fmt, T val)
1541{
1542 double scale = fmt.scale_factor ();
1543
1544 if (Vfixed_point_format && ! (print_g || print_e) && scale != 1)
1545 val /= scale;
1546
1547 pr_any_float (os, fmt.imag_format (), val);
1548}
1549
1550template <typename T>
1551static inline void
1552pr_float (std::ostream& os, const float_display_format& fmt,
1553 const std::complex<T>& cval)
1554{
1555 T r = cval.real ();
1556
1557 pr_float (os, fmt, r);
1558
1559 if (! bank_format)
1560 {
1561 T i = cval.imag ();
1562 if (! (hex_format || bit_format) && std::signbit (i))
1563 {
1564 os << " - ";
1565 i = -i;
1566 pr_imag_float (os, fmt, i);
1567 }
1568 else
1569 {
1570 if (hex_format || bit_format)
1571 os << " ";
1572 else
1573 os << " + ";
1574
1575 pr_imag_float (os, fmt, i);
1576 }
1577 os << 'i';
1578 }
1579}
1580
1581static inline void
1582print_empty_matrix (std::ostream& os, octave_idx_type nr, octave_idx_type nc,
1583 bool pr_as_read_syntax)
1584{
1585 panic_unless (nr == 0 || nc == 0);
1586
1587 if (pr_as_read_syntax)
1588 {
1589 if (nr == 0 && nc == 0)
1590 os << "[]";
1591 else
1592 os << "zeros (" << nr << ", " << nc << ')';
1593 }
1594 else
1595 {
1596 os << "[]";
1597
1599 os << '(' << nr << 'x' << nc << ')';
1600 }
1601}
1602
1603static inline void
1604print_empty_nd_array (std::ostream& os, const dim_vector& dims,
1605 bool pr_as_read_syntax)
1606{
1607 panic_unless (dims.any_zero ());
1608
1609 if (pr_as_read_syntax)
1610 os << "zeros (" << dims.str (',') << ')';
1611 else
1612 {
1613 os << "[]";
1614
1616 os << '(' << dims.str () << ')';
1617 }
1618}
1619
1620static inline void
1621pr_scale_header (std::ostream& os, double scale)
1622{
1623 if (Vfixed_point_format && ! (print_g || print_e) && scale != 1)
1624 {
1625 octave::preserve_stream_state stream_state (os);
1626
1627 os << " "
1628 << std::setw (8) << std::setprecision (1)
1629 << std::setiosflags (std::ios::scientific | std::ios::left)
1630 << scale
1631 << "*\n";
1632
1633 if (! Vcompact_format)
1634 os << "\n";
1635 }
1636}
1637
1638static inline void
1639pr_col_num_header (std::ostream& os, octave_idx_type total_width, int max_width,
1640 octave_idx_type lim, octave_idx_type col, int extra_indent)
1641{
1642 if (total_width > max_width && Vsplit_long_rows)
1643 {
1644 octave::preserve_stream_state stream_state (os);
1645
1646 if (col != 0)
1647 {
1648 if (Vcompact_format)
1649 os << "\n";
1650 else
1651 os << "\n\n";
1652 }
1653
1654 octave_idx_type num_cols = lim - col;
1655
1656 os << std::setw (extra_indent) << "";
1657
1658 if (num_cols == 1)
1659 os << " Column " << col + 1 << ":\n";
1660 else if (num_cols == 2)
1661 os << " Columns " << col + 1 << " and " << lim << ":\n";
1662 else
1663 os << " Columns " << col + 1 << " through " << lim << ":\n";
1664
1665 if (! Vcompact_format)
1666 os << "\n";
1667 }
1668}
1669
1670template <typename T>
1671static inline void
1672pr_plus_format (std::ostream& os, const T& val)
1673{
1674 if (val > T (0))
1675 os << plus_format_chars[0];
1676 else if (val < T (0))
1677 os << plus_format_chars[1];
1678 else
1679 os << plus_format_chars[2];
1680}
1681
1682// FIXME: all this mess with abs is an attempt to avoid seeing
1683//
1684// warning: comparison of unsigned expression < 0 is always false
1685//
1686// from GCC. Isn't there a better way?
1687
1688template <typename T>
1689static inline T
1690abs (T x)
1691{
1692 return x < 0 ? -x : x;
1693}
1694
1695#define INSTANTIATE_ABS(T) \
1696 template T abs (T)
1697
1702
1703#define SPECIALIZE_UABS(T) \
1704 template <> \
1705 inline T \
1706 abs (T x) \
1707 { \
1708 return x; \
1709 }
1710
1715
1716#define MAKE_INT_MATRIX_FORMAT(TYPE) \
1717 template <> \
1718 float_display_format \
1719 make_format (const intNDArray<TYPE>& nda) \
1720 { \
1721 bool isneg = false; \
1722 int digits = 0; \
1723 \
1724 for (octave_idx_type i = 0; i < nda.numel (); i++) \
1725 { \
1726 int new_digits \
1727 = static_cast<int> \
1728 (std::floor (log10 (double (abs (nda(i).value ()))) + 1)); \
1729 \
1730 if (new_digits > digits) \
1731 digits = new_digits; \
1732 \
1733 if (! isneg) \
1734 isneg = (abs (nda(i).value ()) != nda(i).value ()); \
1735 } \
1736 \
1737 return float_display_format (float_format (digits + isneg, 0, 0)); \
1738 }
1739
1748
1749#define MAKE_INT_SCALAR_FORMAT(TYPE) \
1750 template <> \
1751 float_display_format \
1752 make_format (const octave_int<TYPE>& val) \
1753 { \
1754 bool isneg = false; \
1755 int digits \
1756 = static_cast<int> \
1757 (std::floor (log10 (double (abs (val.value ()))) + 1)); \
1758 \
1759 isneg = (abs (val.value ()) != val.value ()); \
1760 \
1761 return float_display_format (float_format (digits + isneg, 0, 0)); \
1762 }
1763
1772
1773void
1775 bool d, bool pr_as_read_syntax)
1776{
1777 octave_print_internal (os, fmt, octave_uint8 (d), pr_as_read_syntax);
1778}
1779
1780void
1781octave_print_internal (std::ostream& os, bool d, bool pr_as_read_syntax)
1782{
1783 octave_print_internal (os, octave_uint8 (d), pr_as_read_syntax);
1784}
1785
1786void
1788 char, bool)
1789{
1790 error ("unexpected call to 'octave_print_internal (std::ostream&, const float_display_format&, char, bool)' - please report this bug");
1791}
1792
1793void
1794octave_print_internal (std::ostream& os, const float_display_format& fmt,
1795 double d, bool pr_as_read_syntax)
1796{
1797 if (pr_as_read_syntax)
1798 os << d;
1799 else if (plus_format)
1800 pr_plus_format (os, d);
1801 else
1802 {
1803 if (free_format)
1804 os << d;
1805 else
1806 pr_float (os, fmt, d);
1807 }
1808}
1809
1810void
1811octave_print_internal (std::ostream& os, const float_display_format& fmt,
1812 float d, bool pr_as_read_syntax)
1813{
1814 if (pr_as_read_syntax)
1815 os << d;
1816 else if (plus_format)
1817 pr_plus_format (os, d);
1818 else
1819 {
1820 if (free_format)
1821 os << d;
1822 else
1823 pr_float (os, fmt, d);
1824 }
1825}
1826
1827template <typename MT>
1828static inline void
1829octave_print_free (std::ostream& os, const MT& m, bool pr_as_read_syntax)
1830{
1831 octave_idx_type nr = m.rows ();
1832 octave_idx_type nc = m.columns ();
1833
1834 if (pr_as_read_syntax)
1835 os << "[\n";
1836
1837 for (octave_idx_type i = 0; i < nr; i++)
1838 {
1839 for (octave_idx_type j = 0; j < nc; j++)
1840 os << ' ' << m.elem (i, j);
1841
1842 if (i < nr - 1)
1843 os << "\n";
1844 }
1845
1846 if (pr_as_read_syntax)
1847 os << ']';
1848}
1849
1850template <typename MT>
1851static inline void
1852pr_plus_format_matrix (std::ostream& os, const MT& m)
1853{
1854 octave_idx_type nr = m.rows ();
1855 octave_idx_type nc = m.columns ();
1856
1857 for (octave_idx_type i = 0; i < nr; i++)
1858 {
1859 for (octave_idx_type j = 0; j < nc; j++)
1860 {
1861 octave_quit ();
1862
1863 pr_plus_format (os, m(i, j));
1864 }
1865
1866 if (i < nr - 1)
1867 os << "\n";
1868 }
1869}
1870
1871static inline int
1872get_column_width (const float_display_format& fmt)
1873{
1874 int r_fw = fmt.real_format().width ();
1875 int i_fw = fmt.imag_format().width ();
1876
1877 int retval = r_fw + i_fw + 2;
1878
1879 if (i_fw && ! (rat_format || bank_format || hex_format || bit_format))
1880 retval += 5;
1881
1882 return retval;
1883}
1884
1885template <typename MT>
1886static void
1887octave_print_matrix_internal (std::ostream& os, const MT& m,
1888 bool pr_as_read_syntax, int extra_indent)
1889{
1890 octave_idx_type nr = m.rows ();
1891 octave_idx_type nc = m.columns ();
1892
1893 if (nr == 0 || nc == 0)
1894 print_empty_matrix (os, nr, nc, pr_as_read_syntax);
1895 else if (plus_format && ! pr_as_read_syntax)
1896 pr_plus_format_matrix (os, m);
1897 else
1898 {
1900 int column_width = get_column_width (fmt);
1901 octave_idx_type total_width = nc * column_width;
1902 octave_idx_type max_width = octave::command_editor::terminal_cols ();
1903
1904 if (pr_as_read_syntax)
1905 max_width -= 4;
1906 else
1907 max_width -= extra_indent;
1908
1909 if (max_width < 0)
1910 max_width = 0;
1911
1912 if (free_format)
1913 {
1914 octave_print_free (os, m, pr_as_read_syntax);
1915 return;
1916 }
1917
1918 octave_idx_type inc = nc;
1919 if (total_width > max_width && Vsplit_long_rows)
1920 {
1921 inc = max_width / column_width;
1922 if (inc == 0)
1923 inc++;
1924 }
1925
1926 if (pr_as_read_syntax)
1927 {
1928 for (octave_idx_type i = 0; i < nr; i++)
1929 {
1930 octave_idx_type col = 0;
1931 while (col < nc)
1932 {
1933 octave_idx_type lim = (col + inc < nc ? col + inc : nc);
1934
1935 for (octave_idx_type j = col; j < lim; j++)
1936 {
1937 octave_quit ();
1938
1939 if (i == 0 && j == 0)
1940 os << "[ ";
1941 else
1942 {
1943 if (j > col)
1944 os << ", ";
1945 else
1946 os << " ";
1947 }
1948
1949 pr_float (os, fmt, m(i, j));
1950 }
1951
1952 col += inc;
1953
1954 if (col >= nc)
1955 {
1956 if (i == nr - 1)
1957 os << " ]";
1958 else
1959 os << ";\n";
1960 }
1961 else
1962 os << " ...\n";
1963 }
1964 }
1965 }
1966 else
1967 {
1968 octave::preserve_stream_state stream_state (os);
1969
1970 pr_scale_header (os, fmt.scale_factor ());
1971
1972 for (octave_idx_type col = 0; col < nc; col += inc)
1973 {
1974 octave_idx_type lim = (col + inc < nc ? col + inc : nc);
1975
1976 pr_col_num_header (os, total_width, max_width, lim, col,
1977 extra_indent);
1978
1979 for (octave_idx_type i = 0; i < nr; i++)
1980 {
1981 os << std::setw (extra_indent) << "";
1982
1983 for (octave_idx_type j = col; j < lim; j++)
1984 {
1985 octave_quit ();
1986
1987 os << " ";
1988
1989 pr_float (os, fmt, m(i, j));
1990 }
1991
1992 if (i < nr - 1)
1993 os << "\n";
1994 }
1995 }
1996 }
1997 }
1998}
1999
2000template <typename DMT>
2001static void
2002octave_print_diag_matrix_internal (std::ostream& os, const DMT& m,
2003 bool pr_as_read_syntax, int extra_indent)
2004{
2005 octave_idx_type nr = m.rows ();
2006 octave_idx_type nc = m.columns ();
2007
2008 if (nr == 0 || nc == 0)
2009 print_empty_matrix (os, nr, nc, pr_as_read_syntax);
2010 else if (plus_format && ! pr_as_read_syntax)
2011 pr_plus_format_matrix (os, m);
2012 else
2013 {
2015 = make_format (typename DMT::full_matrix_type (m.diag ()));
2016 int column_width = get_column_width (fmt);
2017 octave_idx_type total_width = nc * column_width;
2018 octave_idx_type max_width = octave::command_editor::terminal_cols ();
2019
2020 if (pr_as_read_syntax)
2021 max_width -= 4;
2022 else
2023 max_width -= extra_indent;
2024
2025 if (max_width < 0)
2026 max_width = 0;
2027
2028 if (free_format)
2029 {
2030 octave_print_free (os, m, pr_as_read_syntax);
2031 return;
2032 }
2033
2034 octave_idx_type inc = nc;
2035 if (total_width > max_width && Vsplit_long_rows)
2036 {
2037 inc = max_width / column_width;
2038 if (inc == 0)
2039 inc++;
2040 }
2041
2042 if (pr_as_read_syntax)
2043 {
2044 os << "diag (";
2045
2046 octave_idx_type col = 0;
2047 while (col < nc)
2048 {
2049 octave_idx_type lim = (col + inc < nc ? col + inc : nc);
2050
2051 for (octave_idx_type j = col; j < lim; j++)
2052 {
2053 octave_quit ();
2054
2055 if (j == 0)
2056 os << "[ ";
2057 else
2058 {
2059 if (j > col)
2060 os << ", ";
2061 else
2062 os << " ";
2063 }
2064
2065 pr_float (os, fmt, m(j, j));
2066 }
2067
2068 col += inc;
2069
2070 if (col >= nc)
2071 os << " ]";
2072 else
2073 os << " ...\n";
2074 }
2075 os << ')';
2076 }
2077 else
2078 {
2079 octave::preserve_stream_state stream_state (os);
2080
2081 os << "Diagonal Matrix\n";
2082 if (! Vcompact_format)
2083 os << "\n";
2084
2085 pr_scale_header (os, fmt.scale_factor ());
2086
2087 // kluge. Get the true width of a number.
2088 int zero_fw;
2089 {
2090 std::ostringstream tmp_oss;
2091 typename DMT::element_type zero = 0;
2092 pr_float (tmp_oss, fmt, zero);
2093 zero_fw = tmp_oss.str ().length ();
2094 }
2095
2096 for (octave_idx_type col = 0; col < nc; col += inc)
2097 {
2098 octave_idx_type lim = (col + inc < nc ? col + inc : nc);
2099
2100 pr_col_num_header (os, total_width, max_width, lim, col,
2101 extra_indent);
2102
2103 for (octave_idx_type i = 0; i < nr; i++)
2104 {
2105 os << std::setw (extra_indent) << "";
2106
2107 for (octave_idx_type j = col; j < lim; j++)
2108 {
2109 octave_quit ();
2110
2111 os << " ";
2112
2113 if (i == j)
2114 pr_float (os, fmt, m(i, j));
2115 else
2116 os << std::setw (zero_fw) << '0';
2117 }
2118
2119 if (i < nr - 1)
2120 os << "\n";
2121 }
2122 }
2123 }
2124 }
2125}
2126
2127template <typename NDA_T, typename ELT_T, typename MAT_T>
2128void
2129print_nd_array (std::ostream& os, const NDA_T& nda,
2130 bool pr_as_read_syntax)
2131{
2132
2133 if (nda.isempty ())
2134 print_empty_nd_array (os, nda.dims (), pr_as_read_syntax);
2135 else
2136 {
2137
2138 int ndims = nda.ndims ();
2139
2140 const dim_vector& dims = nda.dims ();
2141
2143
2144 octave_idx_type m = 1;
2145
2146 for (int i = 2; i < ndims; i++)
2147 m *= dims(i);
2148
2149 octave_idx_type nr = dims(0);
2150 octave_idx_type nc = dims(1);
2151
2152 for (octave_idx_type i = 0; i < m; i++)
2153 {
2154 octave_quit ();
2155
2156 std::string nm = "ans";
2157
2158 if (m > 1)
2159 {
2160 nm += "(:,:,";
2161
2162 std::ostringstream buf;
2163
2164 for (int k = 2; k < ndims; k++)
2165 {
2166 buf << ra_idx(k) + 1;
2167
2168 if (k < ndims - 1)
2169 buf << ',';
2170 else
2171 buf << ')';
2172 }
2173
2174 nm += buf.str ();
2175 }
2176
2177 Array<octave::idx_vector> idx (dim_vector (ndims, 1));
2178
2179 idx(0) = octave::idx_vector (':');
2180 idx(1) = octave::idx_vector (':');
2181
2182 for (int k = 2; k < ndims; k++)
2183 idx(k) = octave::idx_vector (ra_idx(k));
2184
2185 octave_value page
2186 = MAT_T (Array<ELT_T> (nda.index (idx), dim_vector (nr, nc)));
2187
2188 if (i != m - 1)
2189 {
2190 page.print_with_name (os, nm);
2191 }
2192 else
2193 {
2194 page.print_name_tag (os, nm);
2195 page.print_raw (os);
2196 }
2197
2198 NDA_T::increment_index (ra_idx, dims, 2);
2199 }
2200 }
2201}
2202
2203void
2204octave_print_internal (std::ostream& os, const NDArray& nda,
2205 bool pr_as_read_syntax, int extra_indent)
2206{
2207 switch (nda.ndims ())
2208 {
2209 case 1:
2210 case 2:
2211 octave_print_internal (os, Matrix (nda),
2212 pr_as_read_syntax, extra_indent);
2213 break;
2214
2215 default:
2216 print_nd_array <NDArray, double, Matrix> (os, nda, pr_as_read_syntax);
2217 break;
2218 }
2219}
2220
2221void
2222octave_print_internal (std::ostream& os, const FloatNDArray& nda,
2223 bool pr_as_read_syntax, int extra_indent)
2224{
2225 switch (nda.ndims ())
2226 {
2227 case 1:
2228 case 2:
2230 pr_as_read_syntax, extra_indent);
2231 break;
2232
2233 default:
2234 print_nd_array <FloatNDArray, float, FloatMatrix> (os, nda, pr_as_read_syntax);
2235 break;
2236 }
2237}
2238
2239template <typename T>
2240static inline void
2241pr_plus_format (std::ostream& os, const std::complex<T>& c)
2242{
2243 T rp = c.real ();
2244 T ip = c.imag ();
2245
2246 if (rp == 0)
2247 {
2248 if (ip == 0)
2249 os << ' ';
2250 else
2251 os << 'i';
2252 }
2253 else if (ip == 0)
2254 pr_plus_format (os, rp);
2255 else
2256 os << 'c';
2257}
2258
2259extern void
2260octave_print_internal (std::ostream& os, const float_display_format& fmt,
2261 const Complex& c, bool pr_as_read_syntax)
2262{
2263 if (pr_as_read_syntax)
2264 os << c;
2265 else if (plus_format)
2266 pr_plus_format (os, c);
2267 else
2268 {
2269 if (free_format)
2270 os << c;
2271 else
2272 pr_float (os, fmt, c);
2273 }
2274}
2275
2276void
2277octave_print_internal (std::ostream& os, const float_display_format& fmt,
2278 const FloatComplex& c, bool pr_as_read_syntax)
2279{
2280 if (pr_as_read_syntax)
2281 os << c;
2282 else if (plus_format)
2283 pr_plus_format (os, c);
2284 else
2285 {
2286 if (free_format)
2287 os << c;
2288 else
2289 pr_float (os, fmt, c);
2290 }
2291}
2292
2293void
2294octave_print_internal (std::ostream& os, const PermMatrix& m,
2295 bool pr_as_read_syntax, int extra_indent)
2296{
2297 octave_idx_type nr = m.rows ();
2298 octave_idx_type nc = m.columns ();
2299
2300 if (nr == 0 || nc == 0)
2301 print_empty_matrix (os, nr, nc, pr_as_read_syntax);
2302 else if (plus_format && ! pr_as_read_syntax)
2303 pr_plus_format_matrix (os, m);
2304 else
2305 {
2306 int fw = 2;
2307 int column_width = fw + 2;
2308 octave_idx_type total_width = nc * column_width;
2309 octave_idx_type max_width = octave::command_editor::terminal_cols ();
2310
2311 if (pr_as_read_syntax)
2312 max_width -= 4;
2313 else
2314 max_width -= extra_indent;
2315
2316 if (max_width < 0)
2317 max_width = 0;
2318
2319 if (free_format)
2320 {
2321 octave_print_free (os, m, pr_as_read_syntax);
2322 return;
2323 }
2324
2325 octave_idx_type inc = nc;
2326 if (total_width > max_width && Vsplit_long_rows)
2327 {
2328 inc = max_width / column_width;
2329 if (inc == 0)
2330 inc++;
2331 }
2332
2333 if (pr_as_read_syntax)
2334 {
2336
2337 os << "eye (";
2338 os << ":, ";
2339
2340 octave_idx_type col = 0;
2341 while (col < nc)
2342 {
2343 octave_idx_type lim = (col + inc < nc ? col + inc : nc);
2344
2345 for (octave_idx_type j = col; j < lim; j++)
2346 {
2347 octave_quit ();
2348
2349 if (j == 0)
2350 os << "[ ";
2351 else
2352 {
2353 if (j > col)
2354 os << ", ";
2355 else
2356 os << " ";
2357 }
2358
2359 os << pvec (j);
2360 }
2361
2362 col += inc;
2363
2364 if (col >= nc)
2365 os << " ]";
2366 else
2367 os << " ...\n";
2368 }
2369 os << ')';
2370 }
2371 else
2372 {
2373 octave::preserve_stream_state stream_state (os);
2374
2375 os << "Permutation Matrix\n";
2376 if (! Vcompact_format)
2377 os << "\n";
2378
2379 for (octave_idx_type col = 0; col < nc; col += inc)
2380 {
2381 octave_idx_type lim = (col + inc < nc ? col + inc : nc);
2382
2383 pr_col_num_header (os, total_width, max_width, lim, col,
2384 extra_indent);
2385
2386 for (octave_idx_type i = 0; i < nr; i++)
2387 {
2388 os << std::setw (extra_indent) << "";
2389
2390 for (octave_idx_type j = col; j < lim; j++)
2391 {
2392 octave_quit ();
2393
2394 os << " ";
2395
2396 os << std::setw (fw) << m(i, j);
2397 }
2398
2399 if (i < nr - 1)
2400 os << "\n";
2401 }
2402 }
2403 }
2404 }
2405}
2406
2407void
2408octave_print_internal (std::ostream& os, const ComplexNDArray& nda,
2409 bool pr_as_read_syntax, int extra_indent)
2410{
2411 switch (nda.ndims ())
2412 {
2413 case 1:
2414 case 2:
2416 pr_as_read_syntax, extra_indent);
2417 break;
2418
2419 default:
2420 print_nd_array <ComplexNDArray, Complex, ComplexMatrix>
2421 (os, nda, pr_as_read_syntax);
2422 break;
2423 }
2424}
2425
2426void
2427octave_print_internal (std::ostream& os, const FloatComplexNDArray& nda,
2428 bool pr_as_read_syntax, int extra_indent)
2429{
2430 switch (nda.ndims ())
2431 {
2432 case 1:
2433 case 2:
2435 pr_as_read_syntax, extra_indent);
2436 break;
2437
2438 default:
2439 print_nd_array <FloatComplexNDArray, FloatComplex, FloatComplexMatrix>
2440 (os, nda, pr_as_read_syntax);
2441 break;
2442 }
2443}
2444
2445// FIXME: write single precision versions of the printing functions.
2446
2447void
2448octave_print_internal (std::ostream& os, const Matrix& m,
2449 bool pr_as_read_syntax, int extra_indent)
2450{
2451 octave_print_matrix_internal (os, m, pr_as_read_syntax, extra_indent);
2452}
2453
2454void
2455octave_print_internal (std::ostream& os, const FloatMatrix& m,
2456 bool pr_as_read_syntax, int extra_indent)
2457{
2458 octave_print_matrix_internal (os, m, pr_as_read_syntax, extra_indent);
2459}
2460
2461void
2462octave_print_internal (std::ostream& os, const DiagMatrix& m,
2463 bool pr_as_read_syntax, int extra_indent)
2464{
2465 octave_print_diag_matrix_internal (os, m, pr_as_read_syntax, extra_indent);
2466}
2467
2468void
2469octave_print_internal (std::ostream& os, const FloatDiagMatrix& m,
2470 bool pr_as_read_syntax, int extra_indent)
2471{
2472 octave_print_diag_matrix_internal (os, m, pr_as_read_syntax, extra_indent);
2473}
2474
2475void
2476octave_print_internal (std::ostream& os, const ComplexMatrix& cm,
2477 bool pr_as_read_syntax, int extra_indent)
2478{
2479 octave_print_matrix_internal (os, cm, pr_as_read_syntax, extra_indent);
2480}
2481
2482void
2483octave_print_internal (std::ostream& os, const FloatComplexMatrix& cm,
2484 bool pr_as_read_syntax, int extra_indent)
2485{
2486 octave_print_matrix_internal (os, cm, pr_as_read_syntax, extra_indent);
2487}
2488
2489void
2490octave_print_internal (std::ostream& os, const ComplexDiagMatrix& cm,
2491 bool pr_as_read_syntax, int extra_indent)
2492{
2493 octave_print_diag_matrix_internal (os, cm, pr_as_read_syntax, extra_indent);
2494}
2495
2496void
2498 bool pr_as_read_syntax, int extra_indent)
2499{
2500 octave_print_diag_matrix_internal (os, cm, pr_as_read_syntax, extra_indent);
2501}
2502
2503void
2504octave_print_internal (std::ostream& os, const octave::range<double>& r,
2505 bool pr_as_read_syntax, int extra_indent)
2506{
2507 double base = r.base ();
2508 double increment = r.increment ();
2509 double limit = r.limit ();
2510 double final_value = r.final_value ();
2511 octave_idx_type num_elem = r.numel ();
2512
2513 if (plus_format && ! pr_as_read_syntax)
2514 pr_plus_format_matrix (os, r);
2515 else
2516 {
2518
2519 if (pr_as_read_syntax)
2520 {
2521 if (free_format)
2522 {
2523 os << base << " : ";
2524 if (increment != 1)
2525 os << increment << " : ";
2526 os << limit;
2527 }
2528 else
2529 {
2530 pr_float (os, fmt, base);
2531 os << " : ";
2532 if (increment != 1)
2533 {
2534 pr_float (os, fmt, increment);
2535 os << " : ";
2536 }
2537 pr_float (os, fmt, limit);
2538 }
2539 }
2540 else
2541 {
2542 octave::preserve_stream_state stream_state (os);
2543
2544 int column_width = get_column_width (fmt);
2545 octave_idx_type total_width = num_elem * column_width;
2546 octave_idx_type max_width = octave::command_editor::terminal_cols ();
2547
2548 if (free_format)
2549 {
2550 os << ' ';
2551 for (octave_idx_type i = 0; i < num_elem; i++)
2552 os << ' ' << r.elem (i);
2553 return;
2554 }
2555
2556 octave_idx_type inc = num_elem;
2557 if (total_width > max_width && Vsplit_long_rows)
2558 {
2559 inc = max_width / column_width;
2560 if (inc == 0)
2561 inc++;
2562 }
2563
2564 max_width -= extra_indent;
2565
2566 if (max_width < 0)
2567 max_width = 0;
2568
2569 pr_scale_header (os, fmt.scale_factor ());
2570
2571 octave_idx_type col = 0;
2572 while (col < num_elem)
2573 {
2574 octave_idx_type lim = (col + inc < num_elem ? col + inc
2575 : num_elem);
2576
2577 pr_col_num_header (os, total_width, max_width, lim, col,
2578 extra_indent);
2579
2580 os << std::setw (extra_indent) << "";
2581
2582 for (octave_idx_type i = col; i < lim; i++)
2583 {
2584 octave_quit ();
2585
2586 double val;
2587 if (i == 0)
2588 val = base;
2589 else
2590 val = base + i * increment;
2591
2592 if (i == num_elem - 1)
2593 val = final_value;
2594
2595 os << " ";
2596
2597 pr_float (os, fmt, val);
2598 }
2599
2600 col += inc;
2601 }
2602 }
2603 }
2604}
2605
2606void
2607octave_print_internal (std::ostream& os, const boolMatrix& bm,
2608 bool pr_as_read_syntax,
2609 int extra_indent)
2610{
2611 uint8NDArray tmp (bm);
2612 octave_print_internal (os, tmp, pr_as_read_syntax, extra_indent);
2613}
2614
2615void
2616octave_print_internal (std::ostream& os, const boolNDArray& nda,
2617 bool pr_as_read_syntax,
2618 int extra_indent)
2619{
2620 switch (nda.ndims ())
2621 {
2622 case 1:
2623 case 2:
2625 pr_as_read_syntax, extra_indent);
2626 break;
2627
2628 default:
2630 boolMatrix> (os, nda, pr_as_read_syntax);
2631 break;
2632 }
2633}
2634
2635void
2636octave_print_internal (std::ostream& os, const charMatrix& chm,
2637 bool pr_as_read_syntax,
2638 int /* FIXME: extra_indent */,
2639 bool pr_as_string)
2640{
2641 if (pr_as_string)
2642 {
2643 octave_idx_type nstr = chm.rows ();
2644
2645 if (pr_as_read_syntax && nstr > 1)
2646 os << "[ ";
2647
2648 if (nstr != 0)
2649 {
2650 for (octave_idx_type i = 0; i < nstr; i++)
2651 {
2652 octave_quit ();
2653
2654 std::string row = chm.row_as_string (i);
2655
2656 if (pr_as_read_syntax)
2657 {
2658 os << '"' << octave::undo_string_escapes (row) << '"';
2659
2660 if (i < nstr - 1)
2661 os << "; ";
2662 }
2663 else
2664 {
2665 os << row;
2666
2667 if (i < nstr - 1)
2668 os << "\n";
2669 }
2670 }
2671 }
2672
2673 if (pr_as_read_syntax && nstr > 1)
2674 os << " ]";
2675 }
2676 else
2677 {
2678 os << "sorry, printing char matrices not implemented yet\n";
2679 }
2680}
2681
2682void
2683octave_print_internal (std::ostream& os, const charNDArray& nda,
2684 bool pr_as_read_syntax, int extra_indent,
2685 bool pr_as_string)
2686{
2687 switch (nda.ndims ())
2688 {
2689 case 1:
2690 case 2:
2692 pr_as_read_syntax, extra_indent, pr_as_string);
2693 break;
2694
2695 default:
2696 print_nd_array <charNDArray, char, charMatrix> (os, nda,
2697 pr_as_read_syntax);
2698 break;
2699 }
2700}
2701
2702void
2703octave_print_internal (std::ostream& os, const std::string& s,
2704 bool pr_as_read_syntax, int extra_indent)
2705{
2706 Array<std::string> nda (dim_vector (1, 1), s);
2707
2708 octave_print_internal (os, nda, pr_as_read_syntax, extra_indent);
2709}
2710
2711void
2712octave_print_internal (std::ostream& os, const Array<std::string>& nda,
2713 bool pr_as_read_syntax, int /* extra_indent */)
2714{
2715 // FIXME: this mostly duplicates the code in the print_nd_array<>
2716 // function. Can fix this with std::is_same from C++11.
2717
2718 if (nda.isempty ())
2719 print_empty_nd_array (os, nda.dims (), pr_as_read_syntax);
2720 else if (nda.numel () == 1)
2721 {
2722 os << nda(0);
2723 }
2724 else
2725 {
2726 int ndims = nda.ndims ();
2727
2728 const dim_vector& dims = nda.dims ();
2729
2731
2732 octave_idx_type m = 1;
2733
2734 for (int i = 2; i < ndims; i++)
2735 m *= dims(i);
2736
2737 octave_idx_type nr = dims(0);
2738 octave_idx_type nc = dims(1);
2739
2740 for (octave_idx_type i = 0; i < m; i++)
2741 {
2742 std::string nm = "ans";
2743
2744 if (m > 1)
2745 {
2746 nm += "(:,:,";
2747
2748 std::ostringstream buf;
2749
2750 for (int k = 2; k < ndims; k++)
2751 {
2752 buf << ra_idx(k) + 1;
2753
2754 if (k < ndims - 1)
2755 buf << ',';
2756 else
2757 buf << ')';
2758 }
2759
2760 nm += buf.str ();
2761 }
2762
2763 Array<octave::idx_vector> idx (dim_vector (ndims, 1));
2764
2765 idx(0) = octave::idx_vector (':');
2766 idx(1) = octave::idx_vector (':');
2767
2768 for (int k = 2; k < ndims; k++)
2769 idx(k) = octave::idx_vector (ra_idx(k));
2770
2771 Array<std::string> page (nda.index (idx), dim_vector (nr, nc));
2772
2773 // FIXME: need to do some more work to put these
2774 // in neatly aligned columns...
2775
2776 octave_idx_type n_rows = page.rows ();
2777 octave_idx_type n_cols = page.cols ();
2778
2779 os << nm << " =\n";
2780 if (! Vcompact_format)
2781 os << "\n";
2782
2783 for (octave_idx_type ii = 0; ii < n_rows; ii++)
2784 {
2785 for (octave_idx_type jj = 0; jj < n_cols; jj++)
2786 os << " " << page(ii, jj);
2787
2788 os << "\n";
2789 }
2790
2791 if (i < m - 1)
2792 os << "\n";
2793
2794 increment_index (ra_idx, dims, 2);
2795 }
2796 }
2797}
2798
2799template <typename T>
2800class octave_print_conv
2801{
2802public:
2803 typedef T print_conv_type;
2804};
2805
2806#define PRINT_CONV(T1, T2) \
2807 template <> \
2808 class \
2809 octave_print_conv<T1> \
2810 { \
2811 public: \
2812 typedef T2 print_conv_type; \
2813 }
2814
2817
2818#undef PRINT_CONV
2819
2820template <typename T>
2821static inline void
2822pr_int (std::ostream& os, const T& d, int fw = 0)
2823{
2824 std::size_t sz = d.byte_size ();
2825 const unsigned char *tmpi = d.iptr ();
2826
2827 // Unless explicitly asked for, always print in big-endian
2828 // format for hex and bit formats.
2829 //
2830 // {bit,hex}_format == 1: print big-endian
2831 // {bit,hex}_format == 2: print native
2832
2833 if (hex_format)
2834 {
2835 octave::preserve_stream_state stream_state (os);
2836
2837 os.fill ('0');
2838 if (uppercase_format)
2839 os.flags (std::ios::right | std::ios::hex | std::ios::uppercase);
2840 else
2841 os.flags (std::ios::right | std::ios::hex);
2842
2843 if (hex_format > 1 || octave::mach_info::words_big_endian ())
2844 {
2845 for (std::size_t i = 0; i < sz; i++)
2846 os << std::setw (2) << static_cast<int> (tmpi[i]);
2847 }
2848 else
2849 {
2850 for (int i = sz - 1; i >= 0; i--)
2851 os << std::setw (2) << static_cast<int> (tmpi[i]);
2852 }
2853 }
2854 else if (bit_format)
2855 {
2856 if (octave::mach_info::words_big_endian ())
2857 {
2858 for (std::size_t i = 0; i < sz; i++)
2859 PRINT_CHAR_BITS (os, tmpi[i]);
2860 }
2861 else
2862 {
2863 if (bit_format > 1)
2864 {
2865 for (std::size_t i = 0; i < sz; i++)
2866 PRINT_CHAR_BITS (os, tmpi[i]);
2867 }
2868 else
2869 {
2870 for (int i = sz - 1; i >= 0; i--)
2871 PRINT_CHAR_BITS (os, tmpi[i]);
2872 }
2873 }
2874 }
2875 else
2876 {
2877 octave::preserve_stream_state stream_state (os);
2878
2879 os << std::setw (fw)
2880 << typename octave_print_conv<T>::print_conv_type (d);
2881
2882 if (bank_format)
2883 os << ".00";
2884 }
2885}
2886
2887template void
2888pr_int (std::ostream&, const octave_int8&, int);
2889
2890template void
2891pr_int (std::ostream&, const octave_int16&, int);
2892
2893template void
2894pr_int (std::ostream&, const octave_int32&, int);
2895
2896template void
2897pr_int (std::ostream&, const octave_int64&, int);
2898
2899template void
2900pr_int (std::ostream&, const octave_uint8&, int);
2901
2902template void
2903pr_int (std::ostream&, const octave_uint16&, int);
2904
2905template void
2906pr_int (std::ostream&, const octave_uint32&, int);
2907
2908template void
2909pr_int (std::ostream&, const octave_uint64&, int);
2910
2911template <typename T>
2912void
2914 const float_display_format& fmt,
2915 const octave_int<T>& val, bool)
2916{
2917 if (plus_format)
2918 pr_plus_format (os, val);
2919 else
2920 {
2921 if (free_format)
2922 os << typename octave_print_conv<octave_int<T>>::print_conv_type (val);
2923 else
2924 {
2925 float_format r_fmt = fmt.real_format ();
2926
2927 pr_int (os, val, r_fmt.width ());
2928 }
2929 }
2930}
2931
2932#define PRINT_INT_SCALAR_INTERNAL(TYPE) \
2933 void \
2934 octave_print_internal (std::ostream& os, \
2935 const float_display_format& fmt, \
2936 const octave_int<TYPE>& val, bool dummy) \
2937 { \
2938 octave_print_internal_template (os, fmt, val, dummy); \
2939 }
2940
2949
2950template <typename T>
2951static inline void
2952octave_print_internal_template (std::ostream& os, const intNDArray<T>& nda,
2953 bool pr_as_read_syntax, int extra_indent)
2954{
2955 // FIXME: this mostly duplicates the code in the print_nd_array<>
2956 // function. Can fix this with std::is_same from C++11.
2957
2958 if (nda.isempty ())
2959 print_empty_nd_array (os, nda.dims (), pr_as_read_syntax);
2960 else if (nda.numel () == 1)
2962 pr_as_read_syntax);
2963 else if (plus_format && ! pr_as_read_syntax)
2964 {
2965 int ndims = nda.ndims ();
2966
2968
2969 const dim_vector& dims = nda.dims ();
2970
2971 octave_idx_type m = 1;
2972
2973 for (int i = 2; i < ndims; i++)
2974 m *= dims(i);
2975
2976 octave_idx_type nr = dims(0);
2977 octave_idx_type nc = dims(1);
2978
2979 for (octave_idx_type i = 0; i < m; i++)
2980 {
2981 if (m > 1)
2982 {
2983 std::string nm = "ans(:,:,";
2984
2985 std::ostringstream buf;
2986
2987 for (int k = 2; k < ndims; k++)
2988 {
2989 buf << ra_idx(k) + 1;
2990
2991 if (k < ndims - 1)
2992 buf << ',';
2993 else
2994 buf << ')';
2995 }
2996
2997 nm += buf.str ();
2998
2999 os << nm << " =\n";
3000 if (! Vcompact_format)
3001 os << "\n";
3002 }
3003
3004 Array<octave::idx_vector> idx (dim_vector (ndims, 1));
3005
3006 idx(0) = octave::idx_vector (':');
3007 idx(1) = octave::idx_vector (':');
3008
3009 for (int k = 2; k < ndims; k++)
3010 idx(k) = octave::idx_vector (ra_idx(k));
3011
3012 Array<T> page (nda.index (idx), dim_vector (nr, nc));
3013
3014 for (octave_idx_type ii = 0; ii < nr; ii++)
3015 {
3016 for (octave_idx_type jj = 0; jj < nc; jj++)
3017 {
3018 octave_quit ();
3019
3020 pr_plus_format (os, page(ii, jj));
3021 }
3022
3023 if ((ii < nr - 1) || (i < m -1))
3024 os << "\n";
3025 }
3026
3027 if (i < m - 1)
3028 {
3029 os << "\n";
3030 increment_index (ra_idx, dims, 2);
3031 }
3032 }
3033 }
3034 else
3035 {
3036 int ndims = nda.ndims ();
3037
3038 const dim_vector& dims = nda.dims ();
3039
3041
3042 octave_idx_type m = 1;
3043
3044 for (int i = 2; i < ndims; i++)
3045 m *= dims(i);
3046
3047 octave_idx_type nr = dims(0);
3048 octave_idx_type nc = dims(1);
3049
3050 int fw = 0;
3051 if (hex_format)
3052 fw = 2 * nda(0).byte_size ();
3053 else if (bit_format)
3054 fw = nda(0).nbits ();
3055 else
3056 {
3057 bool isneg = false;
3058 int digits = 0;
3059
3060 for (octave_idx_type i = 0; i < dims.numel (); i++)
3061 {
3062 int new_digits
3063 = static_cast<int>
3064 (std::floor (log10 (double (abs (nda(i).value ()))) + 1));
3065
3066 if (new_digits > digits)
3067 digits = new_digits;
3068
3069 if (! isneg)
3070 isneg = (abs (nda(i).value ()) != nda(i).value ());
3071 }
3072
3073 fw = digits + isneg;
3074 }
3075
3076 int column_width = fw + (rat_format ? 0 : (bank_format ? 5 : 2));
3077 octave_idx_type total_width = nc * column_width;
3078 int max_width = octave::command_editor::terminal_cols () - extra_indent;
3079 octave_idx_type inc = nc;
3080 if (total_width > max_width && Vsplit_long_rows)
3081 {
3082 inc = max_width / column_width;
3083 if (inc == 0)
3084 inc++;
3085 }
3086
3087 for (octave_idx_type i = 0; i < m; i++)
3088 {
3089 if (m > 1)
3090 {
3091 std::string nm = "ans(:,:,";
3092
3093 std::ostringstream buf;
3094
3095 for (int k = 2; k < ndims; k++)
3096 {
3097 buf << ra_idx(k) + 1;
3098
3099 if (k < ndims - 1)
3100 buf << ',';
3101 else
3102 buf << ')';
3103 }
3104
3105 nm += buf.str ();
3106
3107 os << nm << " =\n";
3108 if (! Vcompact_format)
3109 os << "\n";
3110 }
3111
3112 Array<octave::idx_vector> idx (dim_vector (ndims, 1));
3113
3114 idx(0) = octave::idx_vector (':');
3115 idx(1) = octave::idx_vector (':');
3116
3117 for (int k = 2; k < ndims; k++)
3118 idx(k) = octave::idx_vector (ra_idx(k));
3119
3120 Array<T> page (nda.index (idx), dim_vector (nr, nc));
3121
3122 if (free_format)
3123 {
3124 if (pr_as_read_syntax)
3125 os << "[\n";
3126
3127 for (octave_idx_type ii = 0; ii < nr; ii++)
3128 {
3129 for (octave_idx_type jj = 0; jj < nc; jj++)
3130 {
3131 octave_quit ();
3132 os << " ";
3133 os << typename octave_print_conv<T>::print_conv_type (page(ii, jj));
3134 }
3135 os << "\n";
3136 }
3137
3138 if (pr_as_read_syntax)
3139 os << ']';
3140 }
3141 else
3142 {
3143 octave::preserve_stream_state stream_state (os);
3144
3145 octave_idx_type n_rows = page.rows ();
3146 octave_idx_type n_cols = page.cols ();
3147
3148 for (octave_idx_type col = 0; col < n_cols; col += inc)
3149 {
3150 octave_idx_type lim = (col + inc < n_cols ? col + inc
3151 : n_cols);
3152
3153 pr_col_num_header (os, total_width, max_width, lim, col,
3154 extra_indent);
3155
3156 for (octave_idx_type ii = 0; ii < n_rows; ii++)
3157 {
3158 os << std::setw (extra_indent) << "";
3159
3160 for (octave_idx_type jj = col; jj < lim; jj++)
3161 {
3162 octave_quit ();
3163 os << " ";
3164 pr_int (os, page(ii, jj), fw);
3165 }
3166 if ((ii < n_rows - 1) || (i < m -1))
3167 os << "\n";
3168 }
3169 }
3170 }
3171
3172 if (i < m - 1)
3173 {
3174 os << "\n";
3175 increment_index (ra_idx, dims, 2);
3176 }
3177 }
3178 }
3179}
3180
3181#define PRINT_INT_ARRAY_INTERNAL(TYPE) \
3182 OCTINTERP_API void \
3183 octave_print_internal (std::ostream& os, const intNDArray<TYPE>& nda, \
3184 bool pr_as_read_syntax, int extra_indent) \
3185 { \
3186 octave_print_internal_template (os, nda, pr_as_read_syntax, extra_indent); \
3187 }
3188
3197
3198void
3199octave_print_internal (std::ostream&, const Cell&, bool, int, bool)
3200{
3201 error ("unexpected call to 'octave_print_internal (std::ostream&, const Cell&, bool, int, bool)' - please report this bug");
3202}
3203
3204void
3205octave_print_internal (std::ostream&, const octave_value&, bool)
3206{
3207 error ("unexpected call to 'octave_print_internal (std::ostream&, const octave_value&, bool)' - please report this bug");
3208}
3209
3211
3212DEFUN (rats, args, ,
3213 doc: /* -*- texinfo -*-
3214@deftypefn {} {@var{s} =} rats (@var{x})
3215@deftypefnx {} {@var{s} =} rats (@var{x}, @var{len})
3216Convert @var{x} into a rational approximation represented as a string.
3217
3218A rational approximation to a floating point number is a simple fraction
3219with numerator @var{N} and denominator @var{D} such that
3220@code{@var{x} = @var{N}/@var{D}}.
3221
3222The optional second argument defines the maximum length of the string
3223representing the elements of @var{x}. By default, @var{len} is 13.
3224
3225If the length of the smallest possible rational approximation exceeds
3226@var{len}, an asterisk (*) padded with spaces will be returned instead.
3227
3228Example conversion from matrix to string, and back again.
3229
3230@example
3231@group
3232r = rats (hilb (4));
3233x = str2num (r)
3234@end group
3235@end example
3236
3237@seealso{rat, format}
3238@end deftypefn */)
3239{
3240 int nargin = args.length ();
3241
3242 if (nargin < 1 || nargin > 2)
3243 print_usage ();
3244
3245 octave_value arg = args(0);
3246
3247 if (! arg.isnumeric ())
3248 error ("rats: X must be numeric");
3249
3250 if (arg.isempty ())
3251 return ovl (octave_value (""));
3252
3253 // Convert to N-D arrays to 2-D arrays for Matlab compatibility
3254 if (arg.ndims () > 2)
3255 {
3256 dim_vector dv (arg.rows (), arg.numel () / arg.rows ());
3257 arg = arg.reshape (dv);
3258 }
3259
3260 unwind_protect frame;
3261
3262 frame.protect_var (rat_string_len);
3263
3264 rat_string_len = 13;
3265 if (nargin == 2)
3266 rat_string_len = args(1).strict_int_value ();
3267
3268 frame.protect_var (rat_format);
3269
3270 rat_format = true;
3271
3272 std::ostringstream buf;
3273 arg.print (buf);
3274 std::string s = buf.str ();
3275
3276 std::list<std::string> lst;
3277
3278 std::size_t n = 0;
3279 std::size_t s_len = s.length ();
3280
3281 while (n < s_len)
3282 {
3283 std::size_t m = s.find ('\n', n);
3284
3285 if (m == std::string::npos)
3286 {
3287 lst.push_back (s.substr (n));
3288 break;
3289 }
3290 else
3291 {
3292 lst.push_back (s.substr (n, m - n));
3293 n = m + 1;
3294 }
3295 }
3296
3297 return ovl (string_vector (lst));
3298}
3299
3300/*
3301%!test <*56941>
3302%! [old_fmt, old_spacing] = format ();
3303%! unwind_protect
3304%! format short;
3305%! assert (rats (-2.0005, 10), "-4001/2000");
3306%! assert (strtrim (rats (2.0005, 30)), "4001/2000");
3307%! assert (pi - str2num (rats (pi, 30)), 0, 4 * eps);
3308%! assert (e - str2num (rats (e, 30)), 0, 4 * eps);
3309%! unwind_protect_cleanup
3310%! format (old_fmt);
3311%! format (old_spacing);
3312%! end_unwind_protect
3313
3314%!test <*57003>
3315%! x = ones (2,1,3);
3316%! s = rats (x,4);
3317%! assert (ndims (s) == 2);
3318%! assert (rows (s) == 2);
3319%! assert (columns (s) == 3 * 6);
3320
3321%!assert <*57004> (rats ([]), '')
3322
3323%!test <57704>
3324%! [old_fmt, old_spacing] = format ();
3325%! unwind_protect
3326%! format short;
3327%! assert (rats (2.0005, 9), "4001/2000");
3328%! assert (rats (123, 2), " *");
3329%! unwind_protect_cleanup
3330%! format (old_fmt);
3331%! format (old_spacing);
3332%! end_unwind_protect
3333
3334*/
3335
3336DEFUN (disp, args, nargout,
3337 classes: cell char double function_handle int8 int16 int32 int64 logical single struct uint8 uint16 uint32 uint64
3338 doc: /* -*- texinfo -*-
3339@deftypefn {} {} disp (@var{x})
3340@deftypefnx {} {@var{str} =} disp (@var{x})
3341Display the value of @var{x}.
3342
3343For example:
3344
3345@example
3346@group
3347disp ("The value of pi is:"), disp (pi)
3348
3349 @print{} the value of pi is:
3350 @print{} 3.1416
3351@end group
3352@end example
3353
3354@noindent
3355Note that the output from @code{disp} always ends with a newline.
3356
3357If an output value is requested, @code{disp} prints nothing and returns the
3358formatted output in a string.
3359@seealso{fdisp}
3360@end deftypefn */)
3361{
3362 if (args.length () != 1)
3363 print_usage ();
3364
3365 octave_value_list retval;
3366
3367 octave_value arg = args(0);
3368
3369 if (nargout == 0)
3370 arg.print (octave_stdout);
3371 else
3372 {
3373 std::ostringstream buf;
3374 arg.print (buf);
3375 retval = (octave_value (buf.str (), arg.is_dq_string () ? '"' : '\''));
3376 }
3377
3378 return retval;
3379}
3380
3381DEFMETHOD (fdisp, interp, args, ,
3382 classes: cell char double function_handle int8 int16 int32 int64 logical single struct uint8 uint16 uint32 uint64
3383 doc: /* -*- texinfo -*-
3384@deftypefn {} {} fdisp (@var{fid}, @var{x})
3385Display the value of @var{x} on the stream @var{fid}.
3386
3387For example:
3388
3389@example
3390@group
3391fdisp (stdout, "The value of pi is:"), fdisp (stdout, pi)
3392
3393 @print{} the value of pi is:
3394 @print{} 3.1416
3395@end group
3396@end example
3397
3398@noindent
3399Note that the output from @code{fdisp} always ends with a newline.
3400@seealso{disp}
3401@end deftypefn */)
3402{
3403 if (args.length () != 2)
3404 print_usage ();
3405
3406 stream_list& streams = interp.get_stream_list ();
3407
3408 int fid = streams.get_file_number (args(0));
3409
3410 stream os = streams.lookup (fid, "fdisp");
3411
3412 std::ostream *osp = os.preferred_output_stream ();
3413
3414 if (! osp)
3415 error ("fdisp: stream FID not open for writing");
3416
3417 octave_value arg = args(1);
3418 arg.print (*osp);
3419
3420 return ovl ();
3421}
3422
3423/*
3424## FIXME: This test writes values to a file, but then never checks them.
3425%!test
3426%! [old_fmt, old_spacing] = format ();
3427%! unwind_protect
3428%! format short
3429%! fd = tmpfile ();
3430%! for r = [0, Inf -Inf, NaN]
3431%! for i = [0, Inf -Inf, NaN]
3432%! fdisp (fd, complex (r, i));
3433%! endfor
3434%! endfor
3435%! fclose (fd);
3436%! unwind_protect_cleanup
3437%! format (old_fmt);
3438%! format (old_spacing);
3439%! end_unwind_protect
3440
3441%!test
3442%! [old_fmt, old_spacing] = format ();
3443%! unwind_protect
3444%! foo.real = pi * ones (3,20,3);
3445%! foo.complex = pi * ones (3,20,3) + 1i;
3446%! foo.char = repmat ("- Hello World -", [3, 20]);
3447%! foo.cell = {foo.real, foo.complex, foo.char};
3448%! fields = fieldnames (foo);
3449%! for f = 1:numel (fields)
3450%! format loose;
3451%! loose = disp (foo.(fields{f}));
3452%! format compact;
3453%! compact = disp (foo.(fields{f}));
3454%! expected = strrep (loose, "\n\n", "\n");
3455%! assert (expected, compact);
3456%! endfor
3457%! unwind_protect_cleanup
3458%! format (old_fmt);
3459%! format (old_spacing);
3460%! end_unwind_protect
3461*/
3462
3463DEFMETHOD (display, interp, args, ,
3464 classes: cell char double function_handle int8 int16 int32 int64 logical single struct uint8 uint16 uint32 uint64
3465 doc: /* -*- texinfo -*-
3466@deftypefn {} {} display (@var{obj})
3467Display the contents of the object @var{obj} prepended by its name.
3468
3469The Octave interpreter calls the @code{display} function whenever it needs
3470to present a class on-screen. Typically, this would be a statement which
3471does not end in a semicolon to suppress output. For example:
3472
3473@example
3474myclass (@dots{})
3475@end example
3476
3477Or:
3478
3479@example
3480myobj = myclass (@dots{})
3481@end example
3482
3483In general, user-defined classes should overload the @code{disp} method to
3484avoid the default output:
3485
3486@example
3487@group
3488myobj = myclass (@dots{})
3489 @result{} myobj =
3490
3491 <class myclass>
3492@end group
3493@end example
3494
3495When overloading the @code{display} method instead, one has to take care
3496of properly displaying the object's name. This can be done by using the
3497@code{inputname} function.
3498
3499@seealso{disp, class, subsref, subsasgn}
3500@end deftypefn */)
3501{
3502 int nargin = args.length ();
3503
3504 // Matlab apparently accepts two arguments with the second set to the
3505 // inputname of the first. This is undocumented, but we'll use it.
3506 // However, we never call display methods for classes with more than
3507 // one argument.
3508
3509 if (nargin < 1 || nargin > 2)
3510 print_usage ();
3511
3512 std::string name;
3513
3514 if (nargin == 2)
3515 name = args(1).xstring_value ("NAME must be a string");
3516 else
3517 {
3518 string_vector names = args.name_tags ();
3519 name = names(0);
3520 }
3521
3522 // We are here because there is no overloaded display method for this
3523 // object type.
3524
3525 octave_value value = args(0);
3526
3527 // If print_name_tag displays a newline, then also print one after
3528 // disp is done.
3529
3530 bool print_newlines = false;
3531 if (valid_identifier (name))
3532 print_newlines = value.print_name_tag (octave_stdout, name);
3533
3534 // Use feval so that dispatch will also work for disp.
3535
3536 interp.feval ("disp", ovl (value));
3537
3538 if (print_newlines)
3539 octave_stdout << std::endl;
3540
3541 return ovl ();
3542}
3543
3544/*
3545%!test
3546%! [old_fmt, old_spacing] = format ();
3547%! unwind_protect
3548%! format short;
3549%! str = evalc ("x = 1.1; display (x)");
3550%! assert (str, "x = 1.1000\n");
3551%! unwind_protect_cleanup
3552%! format (old_fmt);
3553%! format (old_spacing);
3554%! end_unwind_protect
3555
3556%!test
3557%! [old_fmt, old_spacing] = format ();
3558%! unwind_protect
3559%! format short;
3560%! str = evalc ("display (1.1)");
3561%! assert (str, "1.1000\n");
3562%! unwind_protect_cleanup
3563%! format (old_fmt);
3564%! format (old_spacing);
3565%! end_unwind_protect
3566
3567## Test input validation
3568%!error display ()
3569%!error display (1,2)
3570*/
3571
3572static inline void
3573init_format_state ()
3574{
3575 free_format = false;
3576 plus_format = false;
3577 rat_format = false;
3578 bank_format = false;
3579 hex_format = 0;
3580 bit_format = 0;
3581 print_e = false;
3582 print_g = false;
3583 print_eng = false;
3584}
3585
3586static std::string format_string ("short");
3587
3588static inline void
3589set_format_style (int argc, const string_vector& argv)
3590{
3591 if (--argc == 0)
3592 {
3593 // Special case of no options, reset to default values
3594 init_format_state ();
3595 set_output_prec (5);
3596 format_string = "short";
3597 Vcompact_format = false;
3598 uppercase_format = false;
3599 return;
3600 }
3601
3602 int idx = 1;
3603 std::string format;
3604
3605 unwind_protect frame;
3606
3607 frame.protect_var (bank_format);
3608 frame.protect_var (bit_format);
3609 frame.protect_var (free_format);
3610 frame.protect_var (hex_format);
3611 frame.protect_var (plus_format);
3612 frame.protect_var (plus_format_chars);
3613 frame.protect_var (rat_format);
3614 frame.protect_var (print_e);
3615 frame.protect_var (print_eng);
3616 frame.protect_var (print_g);
3617 frame.protect_var (format_string);
3619 frame.protect_var (uppercase_format);
3620 int prec = output_precision ();
3621 frame.add ([prec] () { set_output_prec (prec); });
3622
3623 format = format_string; // Initialize with existing value
3624 while (argc-- > 0)
3625 {
3626 std::string arg = argv[idx++];
3627 std::transform (arg.begin (), arg.end (), arg.begin (), tolower);
3628
3629 if (arg == "default")
3630 {
3631 format = "short";
3632 init_format_state ();
3633 set_output_prec (5);
3634 Vcompact_format = false;
3635 uppercase_format = false;
3636 }
3637 else if (arg == "short")
3638 {
3639 format = arg;
3640 init_format_state ();
3641 if (argc > 0)
3642 {
3643 arg = argv[idx];
3644 if (arg == "e")
3645 {
3646 print_e = true;
3647 format.append (arg);
3648 argc--; idx++;
3649 }
3650 else if (arg == "g")
3651 {
3652 init_format_state ();
3653 print_g = true;
3654 format.append (arg);
3655 argc--; idx++;
3656 }
3657 else if (arg == "eng")
3658 {
3659 init_format_state ();
3660 print_eng = true;
3661 format.append (arg);
3662 argc--; idx++;
3663 }
3664 }
3665 set_output_prec (5);
3666 }
3667 else if (arg == "shorte")
3668 {
3669 format = arg;
3670 init_format_state ();
3671 print_e = true;
3672 set_output_prec (5);
3673 }
3674 else if (arg == "shortg")
3675 {
3676 format = arg;
3677 init_format_state ();
3678 print_g = true;
3679 set_output_prec (5);
3680 }
3681 else if (arg == "shorteng")
3682 {
3683 format = arg;
3684 init_format_state ();
3685 print_eng = true;
3686 set_output_prec (5);
3687 }
3688 else if (arg == "long")
3689 {
3690 format = arg;
3691 init_format_state ();
3692 if (argc > 0)
3693 {
3694 arg = argv[idx];
3695
3696 if (arg == "e")
3697 {
3698 print_e = true;
3699 format.append (arg);
3700 argc--; idx++;
3701 }
3702 else if (arg == "g")
3703 {
3704 print_g = true;
3705 format.append (arg);
3706 argc--; idx++;
3707 }
3708 else if (arg == "eng")
3709 {
3710 print_eng = true;
3711 format.append (arg);
3712 argc--; idx++;
3713 }
3714 }
3715 set_output_prec (16);
3716 }
3717 else if (arg == "longe")
3718 {
3719 format = arg;
3720 init_format_state ();
3721 print_e = true;
3722 set_output_prec (16);
3723 }
3724 else if (arg == "longg")
3725 {
3726 format = arg;
3727 init_format_state ();
3728 print_g = true;
3729 set_output_prec (16);
3730 }
3731 else if (arg == "longeng")
3732 {
3733 format = arg;
3734 init_format_state ();
3735 print_eng = true;
3736 set_output_prec (16);
3737 }
3738 else if (arg == "hex")
3739 {
3740 format = arg;
3741 init_format_state ();
3742 hex_format = 1;
3743 }
3744 else if (arg == "native-hex")
3745 {
3746 format = arg;
3747 init_format_state ();
3748 hex_format = 2;
3749 }
3750 else if (arg == "bit")
3751 {
3752 format = arg;
3753 init_format_state ();
3754 bit_format = 1;
3755 }
3756 else if (arg == "native-bit")
3757 {
3758 format = arg;
3759 init_format_state ();
3760 bit_format = 2;
3761 }
3762 else if (arg == "+" || arg == "plus")
3763 {
3764 format = arg;
3765 init_format_state ();
3766 plus_format = true;
3767 if (argc > 0)
3768 {
3769 arg = argv[idx];
3770
3771 if (arg.length () == 3)
3772 {
3773 plus_format_chars = arg;
3774 format.append (arg);
3775 argc--; idx++;
3776 }
3777 else
3778 plus_format_chars = "+- ";
3779 }
3780 else
3781 plus_format_chars = "+- ";
3782 }
3783 else if (arg == "rat")
3784 {
3785 format = arg;
3786 init_format_state ();
3787 rat_format = true;
3788 }
3789 else if (arg == "bank")
3790 {
3791 format = arg;
3792 init_format_state ();
3793 bank_format = true;
3794 }
3795 else if (arg == "free")
3796 {
3797 format = arg;
3798 init_format_state ();
3799 free_format = true;
3800 }
3801 else if (arg == "none")
3802 {
3803 format = arg;
3804 init_format_state ();
3805 free_format = true;
3806 }
3807 else if (arg == "compact")
3808 Vcompact_format = true;
3809 else if (arg == "loose")
3810 Vcompact_format = false;
3811 else if (arg == "lowercase")
3812 uppercase_format = false;
3813 else if (arg == "uppercase")
3814 uppercase_format = true;
3815 else
3816 error ("format: unrecognized format state '%s'", arg.c_str ());
3817 }
3818
3819 format_string = format;
3820
3821 // If successful, discard unwind state information
3822 frame.discard ();
3823}
3824
3825DEFUN (format, args, nargout,
3826 doc: /* -*- texinfo -*-
3827@deftypefn {} {} format
3828@deftypefnx {} {} format options
3829@deftypefnx {} {} format (@var{options})
3830@deftypefnx {} {[@var{format}, @var{formatspacing}, @var{uppercase}] =} format
3831Reset or specify the format of the output produced by @code{disp} and Octave's
3832normal echoing mechanism.
3833
3834This command only affects the display of numbers, not how they are stored
3835or computed. To change the internal representation from the default double use
3836one of the conversion functions such as @code{single}, @code{uint8},
3837@code{int64}, etc. Any @code{format} options that change the number of
3838displayed significant digits will also be reflected by the
3839@code{output_precision} function.
3840
3841By default, Octave displays 5 significant digits in a human readable form
3842(option @samp{short}, option @samp{lowercase}, and option @samp{loose} format
3843for matrices). If @code{format} is invoked without any options, or the option
3844@samp{default} is specified, then this default format is restored.
3845
3846Valid format options for floating point numbers are listed in the following
3847table.
3848
3849@table @code
3850@item default
3851Restore the default format state described above.
3852
3853@item short
3854Fixed point format with 5 significant figures (default).
3855
3856@item long
3857Fixed point format with 16 significant figures.
3858
3859As with the @samp{short} format, Octave will switch to an exponential @samp{e}
3860format if it is unable to format a matrix properly using the current format.
3861
3862@item shorte
3863@itemx longe
3864Exponential format. The number to be represented is split between a mantissa
3865and an exponent (power of 10). The mantissa has 5 significant digits in the
3866short format. In the long format, double values are displayed with 16
3867significant digits and single values are displayed with 8. For example,
3868with the @samp{shorte} format, @code{pi} is displayed as @code{3.1416e+00}.
3869Optionally, the trailing @samp{e} can be split into a second argument.
3870
3871@item shortg
3872@itemx longg
3873Optimally choose between fixed point and exponential format based on the
3874magnitude of the number. For example, with the @samp{shortg} format,
3875@w{@code{pi .^ [2; 4; 8; 16; 32]}} is displayed as
3876
3877@example
3878@group
3879ans =
3880
3881 9.8696
3882 97.409
3883 9488.5
3884 9.0032e+07
3885 8.1058e+15
3886@end group
3887@end example
3888
3889Optionally, the trailing @samp{g} can be split into a second argument.
3890
3891@item shorteng
3892@itemx longeng
3893Identical to @samp{shorte} or @samp{longe} but displays the value using an
3894engineering format, where the exponent is divisible by 3. For example, with
3895the @samp{shorteng} format, @code{10 * pi} is displayed as @code{31.416e+00}.
3896Optionally, the trailing @samp{eng} can be split into a second argument.
3897
3898@item free
3899@itemx none
3900Print output in free format, without trying to line up columns of matrices on
3901the decimal point. This is a raw format equivalent to the C++ code
3902@code{std::cout << @var{variable}}. In general, the result is a presentation
3903with 6 significant digits where unnecessary precision (such as trailing zeros
3904for integers) is suppressed. Complex numbers are formatted as numeric pairs
3905like this @code{(0.60419, 0.60709)} instead of like this
3906@code{0.60419 + 0.60709i}.
3907@end table
3908
3909The following formats affect all numeric output (floating point and integer
3910types).
3911
3912@table @asis
3913@item @qcode{"+"}
3914@itemx @qcode{"+"} @qcode{"@var{chars}"}
3915@itemx @code{plus}
3916@itemx @code{plus @var{chars}}
3917Print a @samp{+} symbol for matrix elements greater than zero, a @samp{-}
3918symbol for elements less than zero, and a space for zero matrix elements. This
3919format can be useful for examining the sparsity structure of a large matrix.
3920For very large matrices the function @code{spy} which plots the sparsity
3921pattern will be clearer.
3922
3923The optional argument @var{chars} specifies a list of 3 characters to use for
3924printing values greater than zero, less than zero, and equal to zero. For
3925example, with the format @qcode{"+" "+-."}, the matrix
3926@code{[1, 0, -1; -1, 0, 1]} is displayed as
3927
3928@example
3929@group
3930ans =
3931
3932+.-
3933-.+
3934@end group
3935@end example
3936@end table
3937
3938@table @code
3939@item bank
3940Print variable in a format appropriate for a currency (fixed format with two
3941digits to the right of the decimal point). Only the real part of a variable is
3942displayed, as the imaginary part makes no sense for a currency.
3943
3944@item bit
3945Print the bit representation of numbers in memory, always with the
3946most significant bit first. For example, @code{pi} is printed like this:
3947
3948@iftex
3949@w{@code{0 10000000000 1001001000011111101101010100010001000010110100011000}}
3950@end iftex
3951@ifnottex
3952@example
39530 10000000000 1001001000011111101101010100010001000010110100011000
3954@end example
3955@end ifnottex
3956
3957@noindent
3958where spaces have been added for clarity to show the sign bit, the 11-bit
3959exponent, and the 52-bit mantissa, in that order. Together they represent
3960@code{pi} as an IEEE 754 double precision floating point number in the normal
3961form. Single precision floating point numbers are analogous.
3962
3963@item native-bit
3964Print the bit representation of numbers as stored in memory. For big-endian
3965machines, this is identical to the @code{format bit} layout seen above.
3966For little-endian machines, it will print the bytes in the opposite order,
3967though bits within a byte will still be presented with the most significant
3968bit on the left.
3969
3970For example, the value of @code{pi} in this format on x86-64 is:
3971
3972@iftex
3973@w{@code{00011000 00101101 01000100 01010100 11111011 00100001 00001001 01000000}}
3974@end iftex
3975@ifnottex
3976@example
397700011000 00101101 01000100 01010100 11111011 00100001 00001001 01000000
3978@end example
3979@end ifnottex
3980
3981@noindent
3982shown here with spaces added for clarity. Compare with the previous bit string
3983from @code{format bit} to see the grouping into bytes and their ordering.
3984
3985@item hex
3986The same as @code{format bit} above, except that bits are grouped four at a
3987time into hexadecimal digits for brevity. Thus @code{pi} is represented as:
3988
3989@example
3990400921fb54442d18
3991@end example
3992
3993@item native-hex
3994The same as @code{format native-bit} above, except that bits are grouped four
3995at a time into hexadecimal digits for brevity. Thus @code{pi} is represented
3996on an x86-64 as:
3997
3998@example
3999182d4454fb210940
4000@end example
4001
4002@item rat
4003Print a rational approximation, i.e., values are approximated as the ratio of
4004small integers. For example, with the @samp{rat} format, @code{pi} is
4005displayed as @code{355/113}.
4006@end table
4007
4008The following two options affect the display of scientific and hex notations.
4009
4010@table @code
4011@item lowercase (default)
4012Use a lowercase @samp{e} for the exponent character in scientific notation and
4013lowercase @samp{a-f} for the hex digits representing 10-15.
4014
4015@item uppercase
4016Use an uppercase @samp{E} for the exponent character in scientific notation and
4017uppercase @samp{A-F} for the hex digits representing 10-15.
4018@end table
4019
4020The following two options affect the display of all matrices.
4021
4022@table @code
4023@item compact
4024Remove blank lines around column number labels and between matrices producing
4025more compact output with more data per page.
4026
4027@item loose (default)
4028Insert blank lines above and below column number labels and between matrices to
4029produce a more readable output with less data per page.
4030@end table
4031
4032If @code{format} is called with multiple competing options, the rightmost one
4033is used, except for @samp{default} which will override all other options. In
4034case of an error the format remains unchanged.
4035
4036If called with one to three output arguments, and no inputs, return the current
4037format, format spacing, and uppercase preference. Specifying both outputs and
4038inputs will produce an error.
4039
4040@seealso{fixed_point_format, output_precision, split_long_rows,
4041print_empty_dimensions, rats}
4042@end deftypefn */)
4043{
4044 octave_value_list retval (std::min (nargout, 2));
4045
4046 int nargin = args.length ();
4047
4048 if (nargout == 0)
4049 {
4050 int argc = nargin + 1;
4051
4052 string_vector argv = args.make_argv ("format");
4053
4054 set_format_style (argc, argv);
4055 }
4056 else
4057 {
4058 if (nargin > 0)
4059 warning ("format: cannot query and set format at the same time, ignoring set operation");
4060
4061 if (nargout >= 3)
4062 retval(2) = (uppercase_format ? "uppercase" : "lowercase");
4063
4064 if (nargout >= 2)
4065 retval(1) = (Vcompact_format ? "compact" : "loose");
4066
4067 retval(0) = format_string;
4068 }
4069
4070 return retval;
4071}
4072
4073/*
4074%!test
4075%! [old_fmt, old_spacing, old_uppercase] = format ();
4076%! unwind_protect
4077%! ## Test one of the formats
4078%! format long e;
4079%! format uppercase;
4080%! str = disp (pi);
4081%! assert (str, "3.141592653589793E+00\n");
4082%! str = disp (single (pi));
4083%! assert (str, "3.1415927E+00\n");
4084%! new_fmt = format ();
4085%! assert (new_fmt, "longe");
4086%! ## Test resetting format (method #1)
4087%! format compact;
4088%! [~, new_spacing] = format ();
4089%! assert (new_spacing, "compact");
4090%! format;
4091%! [new_fmt, new_spacing, new_case] = format ();
4092%! assert (new_fmt, "short");
4093%! assert (new_spacing, "loose");
4094%! assert (new_case, "lowercase");
4095%! ## Test resetting format (method #2)
4096%! format compact uppercase long e;
4097%! [new_fmt, new_spacing, new_case] = format ();
4098%! assert (new_fmt, "longe");
4099%! assert (new_spacing, "compact");
4100%! assert (new_case, "uppercase");
4101%! format ("default");
4102%! [new_fmt, new_spacing, new_case] = format ();
4103%! assert (new_fmt, "short");
4104%! assert (new_spacing, "loose");
4105%! assert (new_case, "lowercase");
4106%! unwind_protect_cleanup
4107%! format (old_fmt);
4108%! format (old_spacing);
4109%! format (old_uppercase);
4110%! end_unwind_protect
4111
4112%!test <*53427>
4113%! [old_fmt, old_spacing] = format ();
4114%! unwind_protect
4115%! format; # reset format to short and loose
4116%! format compact; # set compact format
4117%! format long; # set long format
4118%! [fmt, spacing] = format ();
4119%! assert (fmt, "long");
4120%! assert (spacing, "compact");
4121%! unwind_protect_cleanup
4122%! format (old_fmt);
4123%! format (old_spacing);
4124%! end_unwind_protect
4125
4126## Test input validation
4127%!test
4128%! fail ("fmt = format ('long')", "warning", "cannot query and set format");
4129
4130*/
4131
4132DEFUN (fixed_point_format, args, nargout,
4133 doc: /* -*- texinfo -*-
4134@deftypefn {} {@var{val} =} fixed_point_format ()
4135@deftypefnx {} {@var{old_val} =} fixed_point_format (@var{new_val})
4136@deftypefnx {} {@var{old_val} =} fixed_point_format (@var{new_val}, "local")
4137Query or set the internal variable that controls whether Octave will
4138use a scaled format to print matrix values.
4139
4140The scaled format prints a scaling factor on the first line of output chosen
4141such that the largest matrix element can be written with a single leading
4142digit. For example:
4143
4144@example
4145@group
4146fixed_point_format (true)
4147logspace (1, 7, 5)'
4148ans =
4149
4150 1.0e+07 *
4151
4152 0.00000
4153 0.00003
4154 0.00100
4155 0.03162
4156 1.00000
4157@end group
4158@end example
4159
4160@noindent
4161Notice that the first value appears to be 0 when it is actually 1. Because
4162of the possibility for confusion you should be careful about enabling
4163@code{fixed_point_format}.
4164
4165When called from inside a function with the @qcode{"local"} option, the
4166variable is changed locally for the function and any subroutines it calls.
4167The original variable value is restored when exiting the function.
4168@seealso{format, output_precision}
4169@end deftypefn */)
4170{
4171 return set_internal_variable (Vfixed_point_format, args, nargout,
4172 "fixed_point_format");
4173}
4174
4175DEFUN (print_empty_dimensions, args, nargout,
4176 doc: /* -*- texinfo -*-
4177@deftypefn {} {@var{val} =} print_empty_dimensions ()
4178@deftypefnx {} {@var{old_val} =} print_empty_dimensions (@var{new_val})
4179@deftypefnx {} {@var{old_val} =} print_empty_dimensions (@var{new_val}, "local")
4180Query or set the internal variable that controls whether the dimensions of
4181empty matrices are printed along with the empty matrix symbol, @samp{[]}.
4182
4183For example, the expression
4184
4185@example
4186zeros (3, 0)
4187@end example
4188
4189@noindent
4190will print
4191
4192@example
4193ans = [](3x0)
4194@end example
4195
4196When called from inside a function with the @qcode{"local"} option, the
4197variable is changed locally for the function and any subroutines it calls.
4198The original variable value is restored when exiting the function.
4199@seealso{format}
4200@end deftypefn */)
4201{
4202 return set_internal_variable (Vprint_empty_dimensions, args, nargout,
4203 "print_empty_dimensions");
4204}
4205
4206DEFUN (split_long_rows, args, nargout,
4207 doc: /* -*- texinfo -*-
4208@deftypefn {} {@var{val} =} split_long_rows ()
4209@deftypefnx {} {@var{old_val} =} split_long_rows (@var{new_val})
4210@deftypefnx {} {@var{old_val} =} split_long_rows (@var{new_val}, "local")
4211Query or set the internal variable that controls whether rows of a matrix
4212may be split when displayed to a terminal window.
4213
4214If the rows are split, Octave will display the matrix in a series of smaller
4215pieces, each of which can fit within the limits of your terminal width and
4216each set of rows is labeled so that you can easily see which columns are
4217currently being displayed. For example:
4218
4219@example
4220@group
4221octave:13> rand (2,10)
4222ans =
4223
4224 Columns 1 through 6:
4225
4226 0.75883 0.93290 0.40064 0.43818 0.94958 0.16467
4227 0.75697 0.51942 0.40031 0.61784 0.92309 0.40201
4228
4229 Columns 7 through 10:
4230
4231 0.90174 0.11854 0.72313 0.73326
4232 0.44672 0.94303 0.56564 0.82150
4233@end group
4234@end example
4235
4236When called from inside a function with the @qcode{"local"} option, the
4237variable is changed locally for the function and any subroutines it calls.
4238The original variable value is restored when exiting the function.
4239@seealso{format}
4240@end deftypefn */)
4241{
4242 return set_internal_variable (Vsplit_long_rows, args, nargout,
4243 "split_long_rows");
4244}
4245
4246OCTAVE_END_NAMESPACE(octave)
void increment_index(Array< octave_idx_type > &ra_idx, const dim_vector &dimensions, int start_dimension)
Definition Array-util.cc:60
N Dimensional Array with copy-on-write semantics.
Definition Array.h:130
const dim_vector & dims() const
Return a const-reference so that dims ()(i) works efficiently.
Definition Array.h:507
Array< T, Alloc > index(const octave::idx_vector &i) const
Indexing without resizing.
int ndims() const
Size of the specified dimension.
Definition Array.h:679
octave_idx_type rows() const
Definition Array.h:463
octave_idx_type columns() const
Definition Array.h:475
octave_idx_type cols() const
Definition Array.h:473
bool isempty() const
Size of the specified dimension.
Definition Array.h:652
octave_idx_type numel() const
Number of elements in the array.
Definition Array.h:418
Definition Cell.h:41
Template for N-dimensional array classes with like-type math operators.
Definition MArray.h:61
octave_idx_type rows() const
Definition PermMatrix.h:62
octave_idx_type columns() const
Definition PermMatrix.h:64
const Array< octave_idx_type > & col_perm_vec() const
Definition PermMatrix.h:83
void add(F &&fcn, Args &&... args)
void discard(std::size_t num)
void protect_var(T &var)
std::string row_as_string(octave_idx_type, bool strip_ws=false) const
Definition chMatrix.cc:82
Vector representing the dimensions (size) of an Array.
Definition dim-vector.h:90
std::string str(char sep='x') const
Definition dim-vector.cc:68
octave_idx_type numel(int n=0) const
Number of elements that a matrix with this dimensions would have.
Definition dim-vector.h:331
bool any_zero() const
Definition dim-vector.h:312
float_format imag_format() const
Definition pr-flt-fmt.h:232
float_format real_format() const
Definition pr-flt-fmt.h:230
double scale_factor() const
Definition pr-flt-fmt.h:228
float_format & exponent_width(int w)
Definition pr-flt-fmt.h:105
float_format & precision(int p)
Definition pr-flt-fmt.h:93
float_format & uppercase()
Definition pr-flt-fmt.h:81
float_format & width(int w)
Definition pr-flt-fmt.h:99
std::ios::fmtflags format_flags() const
Definition pr-flt-fmt.h:118
void print_raw(std::ostream &os, bool pr_as_read_syntax=false) const
Definition ov.h:1339
void print_with_name(std::ostream &os, const std::string &name) const
Definition ov.h:1345
bool is_dq_string() const
Definition ov.h:643
int ndims() const
Definition ov.h:551
octave_idx_type rows() const
Definition ov.h:545
bool isnumeric() const
Definition ov.h:750
bool isempty() const
Definition ov.h:601
octave_value reshape(const dim_vector &dv) const
Definition ov.h:571
octave_idx_type numel() const
Definition ov.h:559
bool print_name_tag(std::ostream &os, const std::string &name) const
Definition ov.h:1342
void print(std::ostream &os, bool pr_as_read_syntax=false)
Definition ov.h:1336
const float_format m_ff
Definition pr-output.h:515
int exponent() const
Definition pr-output.cc:173
const float_format m_ff
Definition pr-output.h:535
const float_format m_ff
Definition pr-output.h:551
int get_file_number(const octave_value &fid) const
stream lookup(int fid, const std::string &who="") const
std::ostream * preferred_output_stream()
Definition oct-stream.h:454
equiv
Definition cmach-info.c:41
ColumnVector real(const ComplexColumnVector &a)
ColumnVector imag(const ComplexColumnVector &a)
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
void print_usage()
Definition defun-int.h:72
#define DEFMETHOD(name, interp_name, args_name, nargout_name, doc)
Macro to define a builtin method.
Definition defun.h:111
#define DEFUN(name, args_name, nargout_name, doc)
Macro to define a builtin function.
Definition defun.h:56
void warning(const char *fmt,...)
Definition error.cc:1078
void error(const char *fmt,...)
Definition error.cc:1003
void scale(Matrix &m, double x, double y, double z)
Definition graphics.cc:5731
std::complex< T > floor(const std::complex< T > &x)
Definition lo-mappers.h:130
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
F77_RET_T const F77_DBLE * x
F77_RET_T const F77_DBLE const F77_DBLE * f
std::complex< double > Complex
Definition oct-cmplx.h:33
std::complex< float > FloatComplex
Definition oct-cmplx.h:34
octave_int< T > pow(const octave_int< T > &a, const octave_int< T > &b)
std::string rational_approx(T val, int len)
const octave_base_value const Array< octave_idx_type > & ra_idx
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
Definition ovl.h:217
#define octave_stdout
Definition pager.h:301
#define panic_unless(cond)
Definition panic.h:59
int output_precision()
Definition pr-flt-fmt.cc:40
void set_output_prec(int prec)
Definition pr-flt-fmt.cc:46
float_display_format make_scalar_format(const T &val)
Definition pr-output.cc:507
#define PRINT_INT_SCALAR_INTERNAL(TYPE)
#define MAKE_INT_SCALAR_FORMAT(TYPE)
#define PRINT_CONV(T1, T2)
void octave_print_internal_template(std::ostream &os, const float_display_format &fmt, const octave_int< T > &val, bool)
#define PRINT_CHAR_BITS(os, c)
#define SPECIALIZE_UABS(T)
void octave_print_internal(std::ostream &os, const float_display_format &fmt, bool d, bool pr_as_read_syntax)
std::ostream & operator<<(std::ostream &os, const pr_engineering_float< T > &pef)
Definition pr-output.cc:187
#define MAKE_INT_MATRIX_FORMAT(TYPE)
float_display_format make_format(const double &d)
Definition pr-output.cc:525
bool Vprint_empty_dimensions
Definition pr-output.cc:71
#define INSTANTIATE_ABS(T)
#define PRINT_INT_ARRAY_INTERNAL(TYPE)
float_display_format make_complex_scalar_format(const std::complex< T > &c)
Definition pr-output.cc:892
void print_nd_array(std::ostream &os, const NDA_T &nda, bool pr_as_read_syntax)
bool Vcompact_format
Definition pr-output.cc:102
bool valid_identifier(const char *s)
Definition utils.cc:79
std::size_t format(std::ostream &os, const char *fmt,...)
Definition utils.cc:1514
octave_value set_internal_variable(bool &var, const octave_value_list &args, int nargout, const char *nm)
Definition variables.cc:583