GNU Octave  8.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
pr-output.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 1993-2023 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 "ovl.h"
54 #include "oct-stream.h"
56 #include "pager.h"
57 #include "parse.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'.
67 static 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.
75 static bool Vsplit_long_rows = true;
76 
77 // TRUE means don't do any fancy formatting.
78 static bool free_format = false;
79 
80 // TRUE means print plus sign for nonzero, blank for zero.
81 static bool plus_format = false;
82 
83 // First char for > 0, second for < 0, third for == 0.
84 static std::string plus_format_chars = "+- ";
85 
86 // TRUE means always print in a rational approximation
87 static bool rat_format = false;
88 
89 // Used to force the length of the rational approximation string for Frats
90 static int rat_string_len = -1;
91 
92 // TRUE means always print like dollars and cents.
93 static bool bank_format = false;
94 
95 // TRUE means print data in hexadecimal format.
96 static int hex_format = 0;
97 
98 // TRUE means print data in binary-bit-pattern format.
99 static int bit_format = 0;
100 
101 // TRUE means don't put newlines around the column number headers.
102 bool Vcompact_format = false;
103 
104 // TRUE means use an e format.
105 static bool print_e = false;
106 
107 // TRUE means use a g format.
108 static bool print_g = false;
109 
110 // TRUE means print uppercase E in exponent field and A-F in hex format.
111 static bool uppercase_format = false;
112 
113 // TRUE means use an engineering format.
114 static bool print_eng = false;
115 
116 static int
117 calc_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 
139 template <typename T>
140 static inline int
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 
162 template <typename T>
163 static inline int
165 {
166  return 1 + (print_eng
168  : static_cast<int> (std::floor (log10 (x))));
169 }
170 
171 template <typename T>
172 int
174 {
175  return engineering_exponent (m_val);
176 }
177 
178 template <typename T>
179 T
181 {
182  return m_val / std::pow (static_cast<T> (10), exponent ());
183 }
184 
185 template <typename T>
186 std::ostream&
187 operator << (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 
226 template <typename T>
227 std::ostream&
228 operator << (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 
247 template <typename T>
248 std::ostream&
249 operator << (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 
298 template <typename T>
299 static inline T
301 {
302  // We expect a 2-d array.
303  error_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 
331 template <typename T>
332 static inline T
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 
361 template <typename>
363 {
364  static const int DIGITS10;
365  static const int MAX_FIELD_WIDTH;
366 };
367 
368 template <>
369 struct pr_output_traits<double>
370 {
371  static const int DIGITS10;
372  static const int MAX_FIELD_WIDTH;
373 };
374 
377 
378 template <>
379 struct pr_output_traits<float>
380 {
381  static const int DIGITS10;
382  static const int MAX_FIELD_WIDTH;
383 };
384 
387 
388 // FIXME: it would be nice to share more code among these functions,..
389 
390 // Works for double and float.
391 
392 template <typename T>
393 static inline float_display_format
394 make_real_format (int digits, bool inf_or_nan, bool int_only)
395 {
396  float_format fmt;
397 
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 
456  && (print_e || print_g || print_eng
457  || ld + rd > pr_output_traits<T>::DIGITS10
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 
505 template <typename T>
507 make_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 
523 template <>
525 make_format (const double& d)
526 {
527  return make_scalar_format (d);
528 }
529 
530 template <>
532 make_format (const float& f)
533 {
534  return make_scalar_format (f);
535 }
536 
537 template <typename T>
538 static inline float_display_format
539 make_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 
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 
638  && (print_e || print_eng || print_g
639  || (! Vfixed_point_format
640  && (ld + rd > pr_output_traits<T>::DIGITS10
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 
679 template <typename MT>
680 static inline float_display_format
682 {
683  error_unless (m.ndims () == 2);
684 
685  if (free_format)
686  return float_display_format ();
687 
688  bool inf_or_nan = m.any_element_is_inf_or_nan ();
689 
690  bool int_or_inf_or_nan = m.all_elements_are_int_or_inf_or_nan ();
691 
692  MT m_abs = m.abs ();
693 
694  typedef typename MT::element_type ELT_T;
695 
696  ELT_T max_abs = pr_max_internal (m_abs);
697  ELT_T min_abs = pr_min_internal (m_abs);
698 
699  int x_max = (max_abs == 0 ? 0 : num_digits (max_abs));
700 
701  int x_min = (min_abs == 0 ? 0 : num_digits (min_abs));
702 
703  return make_real_matrix_format<ELT_T> (x_max, x_min, inf_or_nan,
704  int_or_inf_or_nan);
705 }
706 
707 template <>
710 {
711  return make_matrix_format (m);
712 }
713 
714 template <>
717 {
718  return make_matrix_format (m);
719 }
720 
721 template <typename T>
722 static inline float_display_format
723 make_complex_format (int x_max, int x_min, int r_x,
724  bool inf_or_nan, int int_only)
725 {
726  float_format r_fmt;
727  float_format i_fmt;
728 
730 
731  int i_fw = 0, r_fw = 0, ld = 0, rd = 0;
732 
733  if (rat_format)
734  {
735  i_fw = 0;
736  r_fw = 0;
737  rd = 0;
738  }
739  else if (bank_format)
740  {
741  int digits = r_x;
742  i_fw = 0;
743  r_fw = (digits <= 0 ? 5 : digits + 4);
744  rd = 2;
745  }
746  else if (hex_format)
747  {
748  r_fw = 2 * sizeof (T);
749  i_fw = 2 * sizeof (T);
750  rd = 0;
751  }
752  else if (bit_format)
753  {
754  r_fw = 8 * sizeof (T);
755  i_fw = 8 * sizeof (T);
756  rd = 0;
757  }
758  else if (inf_or_nan || int_only)
759  {
760  int digits = (x_max > x_min ? x_max : x_min);
761  i_fw = (digits <= 0 ? 1 : digits);
762  r_fw = i_fw + 1;
763  if (inf_or_nan && i_fw < 3)
764  {
765  i_fw = 3;
766  r_fw = 4;
767  }
768 
769  if (int_only)
770  {
771  ld = digits;
772  rd = 0;
773  }
774  }
775  else
776  {
777  int ld_max, rd_max;
778  if (x_max > 0)
779  {
780  ld_max = x_max;
781  rd_max = (prec > x_max ? prec - x_max : prec);
782  x_max++;
783  }
784  else if (x_max < 0)
785  {
786  ld_max = 1;
787  rd_max = (prec > x_max ? prec - x_max : prec);
788  x_max = -x_max + 1;
789  }
790  else
791  {
792  ld_max = 1;
793  rd_max = (prec > 1 ? prec - 1 : prec);
794  x_max = 1;
795  }
796 
797  int ld_min, rd_min;
798  if (x_min > 0)
799  {
800  ld_min = x_min;
801  rd_min = (prec > x_min ? prec - x_min : prec);
802  x_min++;
803  }
804  else if (x_min < 0)
805  {
806  ld_min = 1;
807  rd_min = (prec > x_min ? prec - x_min : prec);
808  x_min = -x_min + 1;
809  }
810  else
811  {
812  ld_min = 1;
813  rd_min = (prec > 1 ? prec - 1 : prec);
814  x_min = 1;
815  }
816 
817  ld = (ld_max > ld_min ? ld_max : ld_min);
818  rd = (rd_max > rd_min ? rd_max : rd_min);
819 
820  i_fw = ld + 1 + rd;
821  r_fw = i_fw + 1;
822  }
823 
825  && (print_e || print_eng || print_g
826  || ld + rd > pr_output_traits<T>::DIGITS10
829  || ld + rd > (1.5 * prec)))
830  {
831  if (print_g)
832  {
833  int width = prec + 6;
834  r_fmt = float_format (width, prec);
835  i_fmt = float_format (width, prec);
836  }
837  else
838  {
839  int ex = 4;
840  if (x_max > 100 || x_min > 100)
841  ex++;
842 
843  if (print_eng)
844  {
845  i_fw = 3 + prec + ex;
846  r_fw = i_fw + 1;
847  if (inf_or_nan && i_fw < 5)
848  {
849  i_fw = 5;
850  r_fw = 6;
851  }
852  r_fmt = float_format (r_fw, ex, prec - 1, std::ios::fixed);
853  i_fmt = float_format (i_fw, ex, prec - 1, std::ios::fixed);
854  }
855  else
856  {
857  i_fw = 1 + prec + ex;
858  r_fw = i_fw + 1;
859  if (inf_or_nan && i_fw < 3)
860  {
861  i_fw = 3;
862  r_fw = 4;
863  }
864  r_fmt = float_format (r_fw, prec - 1, std::ios::scientific);
865  i_fmt = float_format (i_fw, prec - 1, std::ios::scientific);
866  }
867  }
868 
869  if (uppercase_format)
870  {
871  r_fmt.uppercase ();
872  i_fmt.uppercase ();
873  }
874  }
875  else if (! bank_format && (inf_or_nan || int_only))
876  {
877  r_fmt = float_format (r_fw, ld);
878  i_fmt = float_format (i_fw, ld);
879  }
880  else
881  {
882  r_fmt = float_format (r_fw, rd, std::ios::fixed);
883  i_fmt = float_format (i_fw, rd, std::ios::fixed);
884  }
885 
886  return float_display_format (r_fmt, i_fmt);
887 }
888 
889 template <typename T>
891 make_complex_scalar_format (const std::complex<T>& c)
892 {
893  if (free_format)
894  return float_display_format ();
895 
896  T rp = c.real ();
897  T ip = c.imag ();
898 
899  bool inf_or_nan = (octave::math::isinf (c) || octave::math::isnan (c));
900 
901  bool int_only = (octave::math::x_nint (rp) == rp
902  && octave::math::x_nint (ip) == ip);
903 
904  T r_abs = (rp < 0 ? -rp : rp);
905  T i_abs = (ip < 0 ? -ip : ip);
906 
907  int r_x = (! octave::math::isfinite (rp)
908  || r_abs == 0) ? 0 : num_digits (r_abs);
909 
910  int i_x = (! octave::math::isfinite (ip)
911  || 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 
929 template <>
931 make_format (const std::complex<double>& c)
932 {
933  return make_complex_scalar_format (c);
934 }
935 
936 template <>
938 make_format (const std::complex<float>& fc)
939 {
940  return make_complex_scalar_format (fc);
941 }
942 
943 template <typename T>
944 static inline float_display_format
945 make_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 
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
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 
1127 template <typename CMT>
1128 static inline float_display_format
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 ? 0 : num_digits (r_max_abs));
1154 
1155  int r_x_min = (r_min_abs == 0 ? 0 : num_digits (r_min_abs));
1156 
1157  int i_x_max = (i_max_abs == 0 ? 0 : num_digits (i_max_abs));
1158 
1159  int i_x_min = (i_min_abs == 0 ? 0 : num_digits (i_min_abs));
1160 
1161  int x_max = (r_x_max > i_x_max ? r_x_max : i_x_max);
1162  int x_min = (r_x_min > i_x_min ? r_x_min : i_x_min);
1163 
1164  return make_complex_matrix_format<ELT_T> (x_max, x_min, r_x_max, r_x_min,
1165  inf_or_nan, int_or_inf_or_nan);
1166 }
1167 
1168 template <>
1171 {
1172  return make_complex_matrix_format (cm);
1173 }
1174 
1175 template <>
1178 {
1179  return make_complex_matrix_format (cm);
1180 }
1181 
1182 template <>
1185 {
1186  return float_display_format (float_format (1, 1));
1187 }
1188 
1189 template <typename T>
1190 static inline float_display_format
1191 make_range_format (int x_max, int x_min, int all_ints)
1192 {
1193  double scale = ((x_max == 0 || all_ints)
1194  ? 1 : std::pow (10.0, calc_scale_exp (x_max - 1)));
1195 
1196  float_format fmt;
1197 
1199 
1200  int fw = 0, ld = 0, rd = 0;
1201 
1202  if (rat_format)
1203  {
1204  fw = 9;
1205  rd = 0;
1206  }
1207  else if (bank_format)
1208  {
1209  int digits = (x_max > x_min ? x_max : x_min);
1210  fw = (digits < 0 ? 5 : digits + 4);
1211  rd = 2;
1212  }
1213  else if (hex_format)
1214  {
1215  fw = 2 * sizeof (T);
1216  rd = 0;
1217  }
1218  else if (bit_format)
1219  {
1220  fw = 8 * sizeof (T);
1221  rd = 0;
1222  }
1223  else if (all_ints)
1224  {
1225  int digits = (x_max > x_min ? x_max : x_min);
1226  fw = digits + 1;
1227  rd = fw;
1228  }
1229  else if (Vfixed_point_format && ! print_g)
1230  {
1231  rd = prec - 1;
1232  fw = rd + 3;
1233  }
1234  else
1235  {
1236  int ld_max, rd_max;
1237  if (x_max > 0)
1238  {
1239  ld_max = x_max;
1240  rd_max = (prec > x_max ? prec - x_max : prec);
1241  x_max++;
1242  }
1243  else if (x_max < 0)
1244  {
1245  ld_max = 1;
1246  rd_max = (prec > x_max ? prec - x_max : prec);
1247  x_max = -x_max + 1;
1248  }
1249  else
1250  {
1251  ld_max = 1;
1252  rd_max = (prec > 1 ? prec - 1 : prec);
1253  x_max = 1;
1254  }
1255 
1256  int ld_min, rd_min;
1257  if (x_min > 0)
1258  {
1259  ld_min = x_min;
1260  rd_min = (prec > x_min ? prec - x_min : prec);
1261  x_min++;
1262  }
1263  else if (x_min < 0)
1264  {
1265  ld_min = 1;
1266  rd_min = (prec > x_min ? prec - x_min : prec);
1267  x_min = -x_min + 1;
1268  }
1269  else
1270  {
1271  ld_min = 1;
1272  rd_min = (prec > 1 ? prec - 1 : prec);
1273  x_min = 1;
1274  }
1275 
1276  ld = (ld_max > ld_min ? ld_max : ld_min);
1277  rd = (rd_max > rd_min ? rd_max : rd_min);
1278 
1279  fw = ld + rd + 3;
1280  }
1281 
1282  if (! (rat_format || bank_format || hex_format || bit_format)
1283  && (print_e || print_eng || print_g
1284  || (! Vfixed_point_format
1285  && (ld + rd > pr_output_traits<T>::DIGITS10
1287  || ld + rd > (1.5 * prec)))))
1288  {
1289  if (print_g)
1290  fmt = float_format (prec+6, prec);
1291  else
1292  {
1293  int ex = 4;
1294  if (x_max > 100 || x_min > 100)
1295  ex++;
1296 
1297  if (print_eng)
1298  {
1299  fw = 5 + prec + ex;
1300  fmt = float_format (fw, ex, prec - 1, std::ios::fixed);
1301  }
1302  else
1303  {
1304  fw = 3 + prec + ex;
1305  fmt = float_format (fw, prec - 1, std::ios::scientific);
1306  }
1307  }
1308  }
1309  else if (! bank_format && all_ints)
1310  fmt = float_format (fw, rd);
1311  else
1312  fmt = float_format (fw, rd, std::ios::fixed);
1313 
1314  if (uppercase_format)
1315  fmt.uppercase ();
1316 
1317  return float_display_format (scale, fmt);
1318 }
1319 
1320 template <>
1322 make_format (const octave::range<double>& r)
1323 {
1324  if (free_format)
1325  return float_display_format ();
1326 
1327  double r_min = r.base ();
1328  double r_max = r.limit ();
1329 
1330  if (r_max < r_min)
1331  {
1332  double tmp = r_max;
1333  r_max = r_min;
1334  r_min = tmp;
1335  }
1336 
1337  bool all_ints = r.all_elements_are_ints ();
1338 
1339  double max_abs = (r_max < 0 ? -r_max : r_max);
1340  double min_abs = (r_min < 0 ? -r_min : r_min);
1341 
1342  int x_max = (max_abs == 0 ? 0 : num_digits (max_abs));
1343 
1344  int x_min = (min_abs == 0 ? 0 : num_digits (min_abs));
1345 
1346  return make_range_format<double> (x_max, x_min, all_ints);
1347 }
1348 
1349 template <typename T>
1350 union equiv
1351 {
1352  T val;
1353  unsigned char i[sizeof (T)];
1354 };
1355 
1356 #define PRINT_CHAR_BITS(os, c) \
1357  do \
1358  { \
1359  unsigned char ctmp = c; \
1360  char stmp[9]; \
1361  stmp[0] = (ctmp & 0x80) ? '1' : '0'; \
1362  stmp[1] = (ctmp & 0x40) ? '1' : '0'; \
1363  stmp[2] = (ctmp & 0x20) ? '1' : '0'; \
1364  stmp[3] = (ctmp & 0x10) ? '1' : '0'; \
1365  stmp[4] = (ctmp & 0x08) ? '1' : '0'; \
1366  stmp[5] = (ctmp & 0x04) ? '1' : '0'; \
1367  stmp[6] = (ctmp & 0x02) ? '1' : '0'; \
1368  stmp[7] = (ctmp & 0x01) ? '1' : '0'; \
1369  stmp[8] = '\0'; \
1370  os << stmp; \
1371  } \
1372  while (0)
1373 
1374 #define PRINT_CHAR_BITS_SWAPPED(os, c) \
1375  do \
1376  { \
1377  unsigned char ctmp = c; \
1378  char stmp[9]; \
1379  stmp[0] = (ctmp & 0x01) ? '1' : '0'; \
1380  stmp[1] = (ctmp & 0x02) ? '1' : '0'; \
1381  stmp[2] = (ctmp & 0x04) ? '1' : '0'; \
1382  stmp[3] = (ctmp & 0x08) ? '1' : '0'; \
1383  stmp[4] = (ctmp & 0x10) ? '1' : '0'; \
1384  stmp[5] = (ctmp & 0x20) ? '1' : '0'; \
1385  stmp[6] = (ctmp & 0x40) ? '1' : '0'; \
1386  stmp[7] = (ctmp & 0x80) ? '1' : '0'; \
1387  stmp[8] = '\0'; \
1388  os << stmp; \
1389  } \
1390  while (0)
1391 
1392 template <typename T>
1393 static inline void
1394 pr_any_float (std::ostream& os, const float_format& fmt, T val)
1395 {
1396  // Unless explicitly asked for, always print in big-endian format
1397  // for hex and bit formats.
1398  //
1399  // {bit,hex}_format == 1: print big-endian
1400  // {bit,hex}_format == 2: print native
1401 
1402  int fw = fmt.width ();
1403 
1404  if (hex_format)
1405  {
1406  octave::preserve_stream_state stream_state (os);
1407 
1408  equiv<T> tmp;
1409  tmp.val = val;
1410 
1411  // Unless explicitly asked for, always print in big-endian format.
1412 
1413  // FIXME: Will bad things happen if we are interrupted before resetting
1414  // the format flags and fill character?
1415 
1418 
1419  os.fill ('0');
1420  if (uppercase_format)
1421  os.flags (std::ios::right | std::ios::hex | std::ios::uppercase);
1422  else
1423  os.flags (std::ios::right | std::ios::hex);
1424 
1425  if (hex_format > 1
1427  {
1428  for (std::size_t i = 0; i < sizeof (T); i++)
1429  os << std::setw (2) << static_cast<int> (tmp.i[i]);
1430  }
1431  else
1432  {
1433  for (int i = sizeof (T) - 1; i >= 0; i--)
1434  os << std::setw (2) << static_cast<int> (tmp.i[i]);
1435  }
1436  }
1437  else if (bit_format)
1438  {
1439  equiv<T> tmp;
1440  tmp.val = val;
1441 
1444 
1446  {
1447  for (std::size_t i = 0; i < sizeof (T); i++)
1448  PRINT_CHAR_BITS (os, tmp.i[i]);
1449  }
1450  else
1451  {
1452  if (bit_format > 1)
1453  {
1454  for (std::size_t i = 0; i < sizeof (T); i++)
1455  PRINT_CHAR_BITS_SWAPPED (os, tmp.i[i]);
1456  }
1457  else
1458  {
1459  for (int i = sizeof (T) - 1; i >= 0; i--)
1460  PRINT_CHAR_BITS (os, tmp.i[i]);
1461  }
1462  }
1463  }
1464  else if (val == 0)
1465  {
1466  octave::preserve_stream_state stream_state (os);
1467 
1468  if (fw > 0)
1469  os << std::setw (fw) << "0";
1470  else
1471  os << "0";
1472  }
1473  else if (octave::math::isna (val))
1474  {
1475  octave::preserve_stream_state stream_state (os);
1476 
1477  if (fw > 0)
1478  os << std::setw (fw) << "NA";
1479  else
1480  os << "NA";
1481  }
1482  else if (rat_format)
1483  os << pr_rational_float<T> (fmt, val);
1484  else if (octave::math::isinf (val))
1485  {
1486  octave::preserve_stream_state stream_state (os);
1487 
1488  const char *s;
1489  if (val < 0)
1490  s = "-Inf";
1491  else
1492  s = "Inf";
1493 
1494  if (fw > 0)
1495  os << std::setw (fw) << s;
1496  else
1497  os << s;
1498  }
1499  else if (octave::math::isnan (val))
1500  {
1501  octave::preserve_stream_state stream_state (os);
1502 
1503  if (fw > 0)
1504  os << std::setw (fw) << "NaN";
1505  else
1506  os << "NaN";
1507  }
1508  else if (print_eng)
1509  os << pr_engineering_float<T> (fmt, val);
1510  else
1511  os << pr_formatted_float<T> (fmt, val);
1512 }
1513 
1514 template <typename T>
1515 static inline void
1516 pr_float (std::ostream& os, const float_display_format& fmt, T val)
1517 {
1518  double scale = fmt.scale_factor ();
1519 
1520  if (Vfixed_point_format && ! (print_g || print_e) && scale != 1)
1521  val /= scale;
1522 
1523  pr_any_float (os, fmt.real_format (), val);
1524 }
1525 
1526 template <typename T>
1527 static inline void
1528 pr_imag_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.imag_format (), val);
1536 }
1537 
1538 template <typename T>
1539 static inline void
1540 pr_float (std::ostream& os, const float_display_format& fmt,
1541  const std::complex<T>& cval)
1542 {
1543  T r = cval.real ();
1544 
1545  pr_float (os, fmt, r);
1546 
1547  if (! bank_format)
1548  {
1549  T i = cval.imag ();
1550  if (! (hex_format || bit_format) && lo_ieee_signbit (i))
1551  {
1552  os << " - ";
1553  i = -i;
1554  pr_imag_float (os, fmt, i);
1555  }
1556  else
1557  {
1558  if (hex_format || bit_format)
1559  os << " ";
1560  else
1561  os << " + ";
1562 
1563  pr_imag_float (os, fmt, i);
1564  }
1565  os << 'i';
1566  }
1567 }
1568 
1569 static inline void
1571  bool pr_as_read_syntax)
1572 {
1573  error_unless (nr == 0 || nc == 0);
1574 
1575  if (pr_as_read_syntax)
1576  {
1577  if (nr == 0 && nc == 0)
1578  os << "[]";
1579  else
1580  os << "zeros (" << nr << ", " << nc << ')';
1581  }
1582  else
1583  {
1584  os << "[]";
1585 
1587  os << '(' << nr << 'x' << nc << ')';
1588  }
1589 }
1590 
1591 static inline void
1592 print_empty_nd_array (std::ostream& os, const dim_vector& dims,
1593  bool pr_as_read_syntax)
1594 {
1595  error_unless (dims.any_zero ());
1596 
1597  if (pr_as_read_syntax)
1598  os << "zeros (" << dims.str (',') << ')';
1599  else
1600  {
1601  os << "[]";
1602 
1604  os << '(' << dims.str () << ')';
1605  }
1606 }
1607 
1608 static inline void
1609 pr_scale_header (std::ostream& os, double scale)
1610 {
1611  if (Vfixed_point_format && ! (print_g || print_e) && scale != 1)
1612  {
1613  octave::preserve_stream_state stream_state (os);
1614 
1615  os << " "
1616  << std::setw (8) << std::setprecision (1)
1617  << std::setiosflags (std::ios::scientific | std::ios::left)
1618  << scale
1619  << "*\n";
1620 
1621  if (! Vcompact_format)
1622  os << "\n";
1623  }
1624 }
1625 
1626 static inline void
1627 pr_col_num_header (std::ostream& os, octave_idx_type total_width, int max_width,
1628  octave_idx_type lim, octave_idx_type col, int extra_indent)
1629 {
1630  if (total_width > max_width && Vsplit_long_rows)
1631  {
1632  octave::preserve_stream_state stream_state (os);
1633 
1634  if (col != 0)
1635  {
1636  if (Vcompact_format)
1637  os << "\n";
1638  else
1639  os << "\n\n";
1640  }
1641 
1642  octave_idx_type num_cols = lim - col;
1643 
1644  os << std::setw (extra_indent) << "";
1645 
1646  if (num_cols == 1)
1647  os << " Column " << col + 1 << ":\n";
1648  else if (num_cols == 2)
1649  os << " Columns " << col + 1 << " and " << lim << ":\n";
1650  else
1651  os << " Columns " << col + 1 << " through " << lim << ":\n";
1652 
1653  if (! Vcompact_format)
1654  os << "\n";
1655  }
1656 }
1657 
1658 template <typename T>
1659 static inline void
1660 pr_plus_format (std::ostream& os, const T& val)
1661 {
1662  if (val > T (0))
1663  os << plus_format_chars[0];
1664  else if (val < T (0))
1665  os << plus_format_chars[1];
1666  else
1667  os << plus_format_chars[2];
1668 }
1669 
1670 // FIXME: all this mess with abs is an attempt to avoid seeing
1671 //
1672 // warning: comparison of unsigned expression < 0 is always false
1673 //
1674 // from GCC. Isn't there a better way?
1675 
1676 template <typename T>
1677 static inline T
1678 abs (T x)
1679 {
1680  return x < 0 ? -x : x;
1681 }
1682 
1683 #define INSTANTIATE_ABS(T) \
1684  template T abs (T)
1685 
1690 
1691 #define SPECIALIZE_UABS(T) \
1692  template <> \
1693  inline T \
1694  abs (T x) \
1695  { \
1696  return x; \
1697  }
1698 
1703 
1704 #define MAKE_INT_MATRIX_FORMAT(TYPE) \
1705  template <> \
1706  float_display_format \
1707  make_format (const intNDArray<TYPE>& nda) \
1708  { \
1709  bool isneg = false; \
1710  int digits = 0; \
1711  \
1712  for (octave_idx_type i = 0; i < nda.numel (); i++) \
1713  { \
1714  int new_digits \
1715  = static_cast<int> \
1716  (std::floor (log10 (double (abs (nda(i).value ()))) + 1)); \
1717  \
1718  if (new_digits > digits) \
1719  digits = new_digits; \
1720  \
1721  if (! isneg) \
1722  isneg = (abs (nda(i).value ()) != nda(i).value ()); \
1723  } \
1724  \
1725  return float_display_format (float_format (digits + isneg, 0, 0)); \
1726  }
1727 
1736 
1737 #define MAKE_INT_SCALAR_FORMAT(TYPE) \
1738  template <> \
1739  float_display_format \
1740  make_format (const octave_int<TYPE>& val) \
1741  { \
1742  bool isneg = false; \
1743  int digits \
1744  = static_cast<int> \
1745  (std::floor (log10 (double (abs (val.value ()))) + 1)); \
1746  \
1747  isneg = (abs (val.value ()) != val.value ()); \
1748  \
1749  return float_display_format (float_format (digits + isneg, 0, 0)); \
1750  }
1751 
1760 
1761 void
1763  bool d, bool pr_as_read_syntax)
1764 {
1765  octave_print_internal (os, fmt, octave_uint8 (d), pr_as_read_syntax);
1766 }
1767 
1768 void
1769 octave_print_internal (std::ostream& os, bool d, bool pr_as_read_syntax)
1770 {
1771  octave_print_internal (os, octave_uint8 (d), pr_as_read_syntax);
1772 }
1773 
1774 void
1776  char, bool)
1777 {
1778  panic_impossible ();
1779 }
1780 
1781 void
1782 octave_print_internal (std::ostream& os, const float_display_format& fmt,
1783  double d, bool pr_as_read_syntax)
1784 {
1785  if (pr_as_read_syntax)
1786  os << d;
1787  else if (plus_format)
1788  pr_plus_format (os, d);
1789  else
1790  {
1791  if (free_format)
1792  os << d;
1793  else
1794  pr_float (os, fmt, d);
1795  }
1796 }
1797 
1798 void
1799 octave_print_internal (std::ostream& os, const float_display_format& fmt,
1800  float d, bool pr_as_read_syntax)
1801 {
1802  if (pr_as_read_syntax)
1803  os << d;
1804  else if (plus_format)
1805  pr_plus_format (os, d);
1806  else
1807  {
1808  if (free_format)
1809  os << d;
1810  else
1811  pr_float (os, fmt, d);
1812  }
1813 }
1814 
1815 template <typename MT>
1816 static inline void
1817 octave_print_free (std::ostream& os, const MT& m, bool pr_as_read_syntax)
1818 {
1819  octave_idx_type nr = m.rows ();
1820  octave_idx_type nc = m.columns ();
1821 
1822  if (pr_as_read_syntax)
1823  os << "[\n";
1824 
1825  for (octave_idx_type i = 0; i < nr; i++)
1826  {
1827  for (octave_idx_type j = 0; j < nc; j++)
1828  os << ' ' << m.elem (i, j);
1829 
1830  if (i < nr - 1)
1831  os << "\n";
1832  }
1833 
1834  if (pr_as_read_syntax)
1835  os << ']';
1836 }
1837 
1838 template <typename MT>
1839 static inline void
1840 pr_plus_format_matrix (std::ostream& os, const MT& m)
1841 {
1842  octave_idx_type nr = m.rows ();
1843  octave_idx_type nc = m.columns ();
1844 
1845  for (octave_idx_type i = 0; i < nr; i++)
1846  {
1847  for (octave_idx_type j = 0; j < nc; j++)
1848  {
1849  octave_quit ();
1850 
1851  pr_plus_format (os, m(i, j));
1852  }
1853 
1854  if (i < nr - 1)
1855  os << "\n";
1856  }
1857 }
1858 
1859 static inline int
1861 {
1862  int r_fw = fmt.real_format().width ();
1863  int i_fw = fmt.imag_format().width ();
1864 
1865  int retval = r_fw + i_fw + 2;
1866 
1867  if (i_fw && ! (rat_format || bank_format || hex_format || bit_format))
1868  retval += 5;
1869 
1870  return retval;
1871 }
1872 
1873 template <typename MT>
1874 static void
1875 octave_print_matrix_internal (std::ostream& os, const MT& m,
1876  bool pr_as_read_syntax, int extra_indent)
1877 {
1878  octave_idx_type nr = m.rows ();
1879  octave_idx_type nc = m.columns ();
1880 
1881  if (nr == 0 || nc == 0)
1882  print_empty_matrix (os, nr, nc, pr_as_read_syntax);
1883  else if (plus_format && ! pr_as_read_syntax)
1884  pr_plus_format_matrix (os, m);
1885  else
1886  {
1888  int column_width = get_column_width (fmt);
1889  octave_idx_type total_width = nc * column_width;
1890  octave_idx_type max_width = octave::command_editor::terminal_cols ();
1891 
1892  if (pr_as_read_syntax)
1893  max_width -= 4;
1894  else
1895  max_width -= extra_indent;
1896 
1897  if (max_width < 0)
1898  max_width = 0;
1899 
1900  if (free_format)
1901  {
1902  octave_print_free (os, m, pr_as_read_syntax);
1903  return;
1904  }
1905 
1906  octave_idx_type inc = nc;
1907  if (total_width > max_width && Vsplit_long_rows)
1908  {
1909  inc = max_width / column_width;
1910  if (inc == 0)
1911  inc++;
1912  }
1913 
1914  if (pr_as_read_syntax)
1915  {
1916  for (octave_idx_type i = 0; i < nr; i++)
1917  {
1918  octave_idx_type col = 0;
1919  while (col < nc)
1920  {
1921  octave_idx_type lim = (col + inc < nc ? col + inc : nc);
1922 
1923  for (octave_idx_type j = col; j < lim; j++)
1924  {
1925  octave_quit ();
1926 
1927  if (i == 0 && j == 0)
1928  os << "[ ";
1929  else
1930  {
1931  if (j > col)
1932  os << ", ";
1933  else
1934  os << " ";
1935  }
1936 
1937  pr_float (os, fmt, m(i, j));
1938  }
1939 
1940  col += inc;
1941 
1942  if (col >= nc)
1943  {
1944  if (i == nr - 1)
1945  os << " ]";
1946  else
1947  os << ";\n";
1948  }
1949  else
1950  os << " ...\n";
1951  }
1952  }
1953  }
1954  else
1955  {
1956  octave::preserve_stream_state stream_state (os);
1957 
1958  pr_scale_header (os, fmt.scale_factor ());
1959 
1960  for (octave_idx_type col = 0; col < nc; col += inc)
1961  {
1962  octave_idx_type lim = (col + inc < nc ? col + inc : nc);
1963 
1964  pr_col_num_header (os, total_width, max_width, lim, col,
1965  extra_indent);
1966 
1967  for (octave_idx_type i = 0; i < nr; i++)
1968  {
1969  os << std::setw (extra_indent) << "";
1970 
1971  for (octave_idx_type j = col; j < lim; j++)
1972  {
1973  octave_quit ();
1974 
1975  os << " ";
1976 
1977  pr_float (os, fmt, m(i, j));
1978  }
1979 
1980  if (i < nr - 1)
1981  os << "\n";
1982  }
1983  }
1984  }
1985  }
1986 }
1987 
1988 template <typename DMT>
1989 static void
1990 octave_print_diag_matrix_internal (std::ostream& os, const DMT& m,
1991  bool pr_as_read_syntax, int extra_indent)
1992 {
1993  octave_idx_type nr = m.rows ();
1994  octave_idx_type nc = m.columns ();
1995 
1996  if (nr == 0 || nc == 0)
1997  print_empty_matrix (os, nr, nc, pr_as_read_syntax);
1998  else if (plus_format && ! pr_as_read_syntax)
1999  pr_plus_format_matrix (os, m);
2000  else
2001  {
2003  = make_format (typename DMT::full_matrix_type (m.diag ()));
2004  int column_width = get_column_width (fmt);
2005  octave_idx_type total_width = nc * column_width;
2006  octave_idx_type max_width = octave::command_editor::terminal_cols ();
2007 
2008  if (pr_as_read_syntax)
2009  max_width -= 4;
2010  else
2011  max_width -= extra_indent;
2012 
2013  if (max_width < 0)
2014  max_width = 0;
2015 
2016  if (free_format)
2017  {
2018  octave_print_free (os, m, pr_as_read_syntax);
2019  return;
2020  }
2021 
2022  octave_idx_type inc = nc;
2023  if (total_width > max_width && Vsplit_long_rows)
2024  {
2025  inc = max_width / column_width;
2026  if (inc == 0)
2027  inc++;
2028  }
2029 
2030  if (pr_as_read_syntax)
2031  {
2032  os << "diag (";
2033 
2034  octave_idx_type col = 0;
2035  while (col < nc)
2036  {
2037  octave_idx_type lim = (col + inc < nc ? col + inc : nc);
2038 
2039  for (octave_idx_type j = col; j < lim; j++)
2040  {
2041  octave_quit ();
2042 
2043  if (j == 0)
2044  os << "[ ";
2045  else
2046  {
2047  if (j > col)
2048  os << ", ";
2049  else
2050  os << " ";
2051  }
2052 
2053  pr_float (os, fmt, m(j, j));
2054  }
2055 
2056  col += inc;
2057 
2058  if (col >= nc)
2059  os << " ]";
2060  else
2061  os << " ...\n";
2062  }
2063  os << ')';
2064  }
2065  else
2066  {
2067  octave::preserve_stream_state stream_state (os);
2068 
2069  os << "Diagonal Matrix\n";
2070  if (! Vcompact_format)
2071  os << "\n";
2072 
2073  pr_scale_header (os, fmt.scale_factor ());
2074 
2075  // kluge. Get the true width of a number.
2076  int zero_fw;
2077  {
2078  std::ostringstream tmp_oss;
2079  typename DMT::element_type zero = 0;
2080  pr_float (tmp_oss, fmt, zero);
2081  zero_fw = tmp_oss.str ().length ();
2082  }
2083 
2084  for (octave_idx_type col = 0; col < nc; col += inc)
2085  {
2086  octave_idx_type lim = (col + inc < nc ? col + inc : nc);
2087 
2088  pr_col_num_header (os, total_width, max_width, lim, col,
2089  extra_indent);
2090 
2091  for (octave_idx_type i = 0; i < nr; i++)
2092  {
2093  os << std::setw (extra_indent) << "";
2094 
2095  for (octave_idx_type j = col; j < lim; j++)
2096  {
2097  octave_quit ();
2098 
2099  os << " ";
2100 
2101  if (i == j)
2102  pr_float (os, fmt, m(i, j));
2103  else
2104  os << std::setw (zero_fw) << '0';
2105  }
2106 
2107  if (i < nr - 1)
2108  os << "\n";
2109  }
2110  }
2111  }
2112  }
2113 }
2114 
2115 template <typename NDA_T, typename ELT_T, typename MAT_T>
2116 void print_nd_array (std::ostream& os, const NDA_T& nda,
2117  bool pr_as_read_syntax)
2118 {
2119 
2120  if (nda.isempty ())
2121  print_empty_nd_array (os, nda.dims (), pr_as_read_syntax);
2122  else
2123  {
2124 
2125  int ndims = nda.ndims ();
2126 
2127  dim_vector dims = nda.dims ();
2128 
2129  Array<octave_idx_type> ra_idx (dim_vector (ndims, 1), 0);
2130 
2131  octave_idx_type m = 1;
2132 
2133  for (int i = 2; i < ndims; i++)
2134  m *= dims(i);
2135 
2136  octave_idx_type nr = dims(0);
2137  octave_idx_type nc = dims(1);
2138 
2139  for (octave_idx_type i = 0; i < m; i++)
2140  {
2141  octave_quit ();
2142 
2143  std::string nm = "ans";
2144 
2145  if (m > 1)
2146  {
2147  nm += "(:,:,";
2148 
2149  std::ostringstream buf;
2150 
2151  for (int k = 2; k < ndims; k++)
2152  {
2153  buf << ra_idx(k) + 1;
2154 
2155  if (k < ndims - 1)
2156  buf << ',';
2157  else
2158  buf << ')';
2159  }
2160 
2161  nm += buf.str ();
2162  }
2163 
2164  Array<octave::idx_vector> idx (dim_vector (ndims, 1));
2165 
2166  idx(0) = octave::idx_vector (':');
2167  idx(1) = octave::idx_vector (':');
2168 
2169  for (int k = 2; k < ndims; k++)
2170  idx(k) = octave::idx_vector (ra_idx(k));
2171 
2172  octave_value page
2173  = MAT_T (Array<ELT_T> (nda.index (idx), dim_vector (nr, nc)));
2174 
2175  if (i != m - 1)
2176  {
2177  page.print_with_name (os, nm);
2178  }
2179  else
2180  {
2181  page.print_name_tag (os, nm);
2182  page.print_raw (os);
2183  }
2184 
2185  NDA_T::increment_index (ra_idx, dims, 2);
2186  }
2187  }
2188 }
2189 
2190 void
2191 octave_print_internal (std::ostream& os, const NDArray& nda,
2192  bool pr_as_read_syntax, int extra_indent)
2193 {
2194  switch (nda.ndims ())
2195  {
2196  case 1:
2197  case 2:
2198  octave_print_internal (os, Matrix (nda),
2199  pr_as_read_syntax, extra_indent);
2200  break;
2201 
2202  default:
2203  print_nd_array <NDArray, double, Matrix> (os, nda, pr_as_read_syntax);
2204  break;
2205  }
2206 }
2207 
2208 void
2209 octave_print_internal (std::ostream& os, const FloatNDArray& nda,
2210  bool pr_as_read_syntax, int extra_indent)
2211 {
2212  switch (nda.ndims ())
2213  {
2214  case 1:
2215  case 2:
2217  pr_as_read_syntax, extra_indent);
2218  break;
2219 
2220  default:
2221  print_nd_array <FloatNDArray, float, FloatMatrix> (os, nda, pr_as_read_syntax);
2222  break;
2223  }
2224 }
2225 
2226 template <typename T>
2227 static inline void
2228 pr_plus_format (std::ostream& os, const std::complex<T>& c)
2229 {
2230  T rp = c.real ();
2231  T ip = c.imag ();
2232 
2233  if (rp == 0)
2234  {
2235  if (ip == 0)
2236  os << ' ';
2237  else
2238  os << 'i';
2239  }
2240  else if (ip == 0)
2241  pr_plus_format (os, rp);
2242  else
2243  os << 'c';
2244 }
2245 
2246 extern void
2247 octave_print_internal (std::ostream& os, const float_display_format& fmt,
2248  const Complex& c, bool pr_as_read_syntax)
2249 {
2250  if (pr_as_read_syntax)
2251  os << c;
2252  else if (plus_format)
2253  pr_plus_format (os, c);
2254  else
2255  {
2256  if (free_format)
2257  os << c;
2258  else
2259  pr_float (os, fmt, c);
2260  }
2261 }
2262 
2263 void
2264 octave_print_internal (std::ostream& os, const float_display_format& fmt,
2265  const FloatComplex& c, bool pr_as_read_syntax)
2266 {
2267  if (pr_as_read_syntax)
2268  os << c;
2269  else if (plus_format)
2270  pr_plus_format (os, c);
2271  else
2272  {
2273  if (free_format)
2274  os << c;
2275  else
2276  pr_float (os, fmt, c);
2277  }
2278 }
2279 
2280 void
2281 octave_print_internal (std::ostream& os, const PermMatrix& m,
2282  bool pr_as_read_syntax, int extra_indent)
2283 {
2284  octave_idx_type nr = m.rows ();
2285  octave_idx_type nc = m.columns ();
2286 
2287  if (nr == 0 || nc == 0)
2288  print_empty_matrix (os, nr, nc, pr_as_read_syntax);
2289  else if (plus_format && ! pr_as_read_syntax)
2290  pr_plus_format_matrix (os, m);
2291  else
2292  {
2293  int fw = 2;
2294  int column_width = fw + 2;
2295  octave_idx_type total_width = nc * column_width;
2296  octave_idx_type max_width = octave::command_editor::terminal_cols ();
2297 
2298  if (pr_as_read_syntax)
2299  max_width -= 4;
2300  else
2301  max_width -= extra_indent;
2302 
2303  if (max_width < 0)
2304  max_width = 0;
2305 
2306  if (free_format)
2307  {
2308  octave_print_free (os, m, pr_as_read_syntax);
2309  return;
2310  }
2311 
2312  octave_idx_type inc = nc;
2313  if (total_width > max_width && Vsplit_long_rows)
2314  {
2315  inc = max_width / column_width;
2316  if (inc == 0)
2317  inc++;
2318  }
2319 
2320  if (pr_as_read_syntax)
2321  {
2322  Array<octave_idx_type> pvec = m.col_perm_vec ();
2323 
2324  os << "eye (";
2325  os << ":, ";
2326 
2327  octave_idx_type col = 0;
2328  while (col < nc)
2329  {
2330  octave_idx_type lim = (col + inc < nc ? col + inc : nc);
2331 
2332  for (octave_idx_type j = col; j < lim; j++)
2333  {
2334  octave_quit ();
2335 
2336  if (j == 0)
2337  os << "[ ";
2338  else
2339  {
2340  if (j > col)
2341  os << ", ";
2342  else
2343  os << " ";
2344  }
2345 
2346  os << pvec (j);
2347  }
2348 
2349  col += inc;
2350 
2351  if (col >= nc)
2352  os << " ]";
2353  else
2354  os << " ...\n";
2355  }
2356  os << ')';
2357  }
2358  else
2359  {
2360  octave::preserve_stream_state stream_state (os);
2361 
2362  os << "Permutation Matrix\n";
2363  if (! Vcompact_format)
2364  os << "\n";
2365 
2366  for (octave_idx_type col = 0; col < nc; col += inc)
2367  {
2368  octave_idx_type lim = (col + inc < nc ? col + inc : nc);
2369 
2370  pr_col_num_header (os, total_width, max_width, lim, col,
2371  extra_indent);
2372 
2373  for (octave_idx_type i = 0; i < nr; i++)
2374  {
2375  os << std::setw (extra_indent) << "";
2376 
2377  for (octave_idx_type j = col; j < lim; j++)
2378  {
2379  octave_quit ();
2380 
2381  os << " ";
2382 
2383  os << std::setw (fw) << m(i, j);
2384  }
2385 
2386  if (i < nr - 1)
2387  os << "\n";
2388  }
2389  }
2390  }
2391  }
2392 }
2393 
2394 void
2395 octave_print_internal (std::ostream& os, const ComplexNDArray& nda,
2396  bool pr_as_read_syntax, int extra_indent)
2397 {
2398  switch (nda.ndims ())
2399  {
2400  case 1:
2401  case 2:
2403  pr_as_read_syntax, extra_indent);
2404  break;
2405 
2406  default:
2407  print_nd_array <ComplexNDArray, Complex, ComplexMatrix>
2408  (os, nda, pr_as_read_syntax);
2409  break;
2410  }
2411 }
2412 
2413 void
2414 octave_print_internal (std::ostream& os, const FloatComplexNDArray& nda,
2415  bool pr_as_read_syntax, int extra_indent)
2416 {
2417  switch (nda.ndims ())
2418  {
2419  case 1:
2420  case 2:
2422  pr_as_read_syntax, extra_indent);
2423  break;
2424 
2425  default:
2426  print_nd_array <FloatComplexNDArray, FloatComplex, FloatComplexMatrix>
2427  (os, nda, pr_as_read_syntax);
2428  break;
2429  }
2430 }
2431 
2432 // FIXME: write single precision versions of the printing functions.
2433 
2434 void
2435 octave_print_internal (std::ostream& os, const Matrix& m,
2436  bool pr_as_read_syntax, int extra_indent)
2437 {
2438  octave_print_matrix_internal (os, m, pr_as_read_syntax, extra_indent);
2439 }
2440 
2441 void
2442 octave_print_internal (std::ostream& os, const FloatMatrix& m,
2443  bool pr_as_read_syntax, int extra_indent)
2444 {
2445  octave_print_matrix_internal (os, m, pr_as_read_syntax, extra_indent);
2446 }
2447 
2448 void
2449 octave_print_internal (std::ostream& os, const DiagMatrix& m,
2450  bool pr_as_read_syntax, int extra_indent)
2451 {
2452  octave_print_diag_matrix_internal (os, m, pr_as_read_syntax, extra_indent);
2453 }
2454 
2455 void
2456 octave_print_internal (std::ostream& os, const FloatDiagMatrix& m,
2457  bool pr_as_read_syntax, int extra_indent)
2458 {
2459  octave_print_diag_matrix_internal (os, m, pr_as_read_syntax, extra_indent);
2460 }
2461 
2462 void
2463 octave_print_internal (std::ostream& os, const ComplexMatrix& cm,
2464  bool pr_as_read_syntax, int extra_indent)
2465 {
2466  octave_print_matrix_internal (os, cm, pr_as_read_syntax, extra_indent);
2467 }
2468 
2469 void
2470 octave_print_internal (std::ostream& os, const FloatComplexMatrix& cm,
2471  bool pr_as_read_syntax, int extra_indent)
2472 {
2473  octave_print_matrix_internal (os, cm, pr_as_read_syntax, extra_indent);
2474 }
2475 
2476 void
2477 octave_print_internal (std::ostream& os, const ComplexDiagMatrix& cm,
2478  bool pr_as_read_syntax, int extra_indent)
2479 {
2480  octave_print_diag_matrix_internal (os, cm, pr_as_read_syntax, extra_indent);
2481 }
2482 
2483 void
2484 octave_print_internal (std::ostream& os, const FloatComplexDiagMatrix& cm,
2485  bool pr_as_read_syntax, int extra_indent)
2486 {
2487  octave_print_diag_matrix_internal (os, cm, pr_as_read_syntax, extra_indent);
2488 }
2489 
2490 void
2491 octave_print_internal (std::ostream& os, const octave::range<double>& r,
2492  bool pr_as_read_syntax, int extra_indent)
2493 {
2494  double base = r.base ();
2495  double increment = r.increment ();
2496  double limit = r.limit ();
2497  double final_value = r.final_value ();
2498  octave_idx_type num_elem = r.numel ();
2499 
2500  if (plus_format && ! pr_as_read_syntax)
2501  pr_plus_format_matrix (os, r);
2502  else
2503  {
2505 
2506  if (pr_as_read_syntax)
2507  {
2508  if (free_format)
2509  {
2510  os << base << " : ";
2511  if (increment != 1)
2512  os << increment << " : ";
2513  os << limit;
2514  }
2515  else
2516  {
2517  pr_float (os, fmt, base);
2518  os << " : ";
2519  if (increment != 1)
2520  {
2521  pr_float (os, fmt, increment);
2522  os << " : ";
2523  }
2524  pr_float (os, fmt, limit);
2525  }
2526  }
2527  else
2528  {
2529  octave::preserve_stream_state stream_state (os);
2530 
2531  int column_width = get_column_width (fmt);
2532  octave_idx_type total_width = num_elem * column_width;
2533  octave_idx_type max_width = octave::command_editor::terminal_cols ();
2534 
2535  if (free_format)
2536  {
2537  os << ' ';
2538  for (octave_idx_type i = 0; i < num_elem; i++)
2539  os << ' ' << r.elem (i);
2540  return;
2541  }
2542 
2543  octave_idx_type inc = num_elem;
2544  if (total_width > max_width && Vsplit_long_rows)
2545  {
2546  inc = max_width / column_width;
2547  if (inc == 0)
2548  inc++;
2549  }
2550 
2551  max_width -= extra_indent;
2552 
2553  if (max_width < 0)
2554  max_width = 0;
2555 
2556  pr_scale_header (os, fmt.scale_factor ());
2557 
2558  octave_idx_type col = 0;
2559  while (col < num_elem)
2560  {
2561  octave_idx_type lim = (col + inc < num_elem ? col + inc
2562  : num_elem);
2563 
2564  pr_col_num_header (os, total_width, max_width, lim, col,
2565  extra_indent);
2566 
2567  os << std::setw (extra_indent) << "";
2568 
2569  for (octave_idx_type i = col; i < lim; i++)
2570  {
2571  octave_quit ();
2572 
2573  double val;
2574  if (i == 0)
2575  val = base;
2576  else
2577  val = base + i * increment;
2578 
2579  if (i == num_elem - 1)
2580  val = final_value;
2581 
2582  os << " ";
2583 
2584  pr_float (os, fmt, val);
2585  }
2586 
2587  col += inc;
2588  }
2589  }
2590  }
2591 }
2592 
2593 void
2594 octave_print_internal (std::ostream& os, const boolMatrix& bm,
2595  bool pr_as_read_syntax,
2596  int extra_indent)
2597 {
2598  uint8NDArray tmp (bm);
2599  octave_print_internal (os, tmp, pr_as_read_syntax, extra_indent);
2600 }
2601 
2602 void
2603 octave_print_internal (std::ostream& os, const boolNDArray& nda,
2604  bool pr_as_read_syntax,
2605  int extra_indent)
2606 {
2607  switch (nda.ndims ())
2608  {
2609  case 1:
2610  case 2:
2611  octave_print_internal (os, boolMatrix (nda),
2612  pr_as_read_syntax, extra_indent);
2613  break;
2614 
2615  default:
2617  boolMatrix> (os, nda, pr_as_read_syntax);
2618  break;
2619  }
2620 }
2621 
2622 void
2623 octave_print_internal (std::ostream& os, const charMatrix& chm,
2624  bool pr_as_read_syntax,
2625  int /* FIXME: extra_indent */,
2626  bool pr_as_string)
2627 {
2628  if (pr_as_string)
2629  {
2630  octave_idx_type nstr = chm.rows ();
2631 
2632  if (pr_as_read_syntax && nstr > 1)
2633  os << "[ ";
2634 
2635  if (nstr != 0)
2636  {
2637  for (octave_idx_type i = 0; i < nstr; i++)
2638  {
2639  octave_quit ();
2640 
2641  std::string row = chm.row_as_string (i);
2642 
2643  if (pr_as_read_syntax)
2644  {
2645  os << '"' << octave::undo_string_escapes (row) << '"';
2646 
2647  if (i < nstr - 1)
2648  os << "; ";
2649  }
2650  else
2651  {
2652  os << row;
2653 
2654  if (i < nstr - 1)
2655  os << "\n";
2656  }
2657  }
2658  }
2659 
2660  if (pr_as_read_syntax && nstr > 1)
2661  os << " ]";
2662  }
2663  else
2664  {
2665  os << "sorry, printing char matrices not implemented yet\n";
2666  }
2667 }
2668 
2669 void
2670 octave_print_internal (std::ostream& os, const charNDArray& nda,
2671  bool pr_as_read_syntax, int extra_indent,
2672  bool pr_as_string)
2673 {
2674  switch (nda.ndims ())
2675  {
2676  case 1:
2677  case 2:
2678  octave_print_internal (os, charMatrix (nda),
2679  pr_as_read_syntax, extra_indent, pr_as_string);
2680  break;
2681 
2682  default:
2683  print_nd_array <charNDArray, char, charMatrix> (os, nda,
2684  pr_as_read_syntax);
2685  break;
2686  }
2687 }
2688 
2689 void
2690 octave_print_internal (std::ostream& os, const std::string& s,
2691  bool pr_as_read_syntax, int extra_indent)
2692 {
2693  Array<std::string> nda (dim_vector (1, 1), s);
2694 
2695  octave_print_internal (os, nda, pr_as_read_syntax, extra_indent);
2696 }
2697 
2698 void
2699 octave_print_internal (std::ostream& os, const Array<std::string>& nda,
2700  bool pr_as_read_syntax, int /* extra_indent */)
2701 {
2702  // FIXME: this mostly duplicates the code in the print_nd_array<>
2703  // function. Can fix this with std::is_same from C++11.
2704 
2705  if (nda.isempty ())
2706  print_empty_nd_array (os, nda.dims (), pr_as_read_syntax);
2707  else if (nda.numel () == 1)
2708  {
2709  os << nda(0);
2710  }
2711  else
2712  {
2713  int ndims = nda.ndims ();
2714 
2715  dim_vector dims = nda.dims ();
2716 
2717  Array<octave_idx_type> ra_idx (dim_vector (ndims, 1), 0);
2718 
2719  octave_idx_type m = 1;
2720 
2721  for (int i = 2; i < ndims; i++)
2722  m *= dims(i);
2723 
2724  octave_idx_type nr = dims(0);
2725  octave_idx_type nc = dims(1);
2726 
2727  for (octave_idx_type i = 0; i < m; i++)
2728  {
2729  std::string nm = "ans";
2730 
2731  if (m > 1)
2732  {
2733  nm += "(:,:,";
2734 
2735  std::ostringstream buf;
2736 
2737  for (int k = 2; k < ndims; k++)
2738  {
2739  buf << ra_idx(k) + 1;
2740 
2741  if (k < ndims - 1)
2742  buf << ',';
2743  else
2744  buf << ')';
2745  }
2746 
2747  nm += buf.str ();
2748  }
2749 
2750  Array<octave::idx_vector> idx (dim_vector (ndims, 1));
2751 
2752  idx(0) = octave::idx_vector (':');
2753  idx(1) = octave::idx_vector (':');
2754 
2755  for (int k = 2; k < ndims; k++)
2756  idx(k) = octave::idx_vector (ra_idx(k));
2757 
2758  Array<std::string> page (nda.index (idx), dim_vector (nr, nc));
2759 
2760  // FIXME: need to do some more work to put these
2761  // in neatly aligned columns...
2762 
2763  octave_idx_type n_rows = page.rows ();
2764  octave_idx_type n_cols = page.cols ();
2765 
2766  os << nm << " =\n";
2767  if (! Vcompact_format)
2768  os << "\n";
2769 
2770  for (octave_idx_type ii = 0; ii < n_rows; ii++)
2771  {
2772  for (octave_idx_type jj = 0; jj < n_cols; jj++)
2773  os << " " << page(ii, jj);
2774 
2775  os << "\n";
2776  }
2777 
2778  if (i < m - 1)
2779  os << "\n";
2780 
2781  increment_index (ra_idx, dims, 2);
2782  }
2783  }
2784 }
2785 
2786 template <typename T>
2787 class
2789 {
2790 public:
2791  typedef T print_conv_type;
2792 };
2793 
2794 #define PRINT_CONV(T1, T2) \
2795  template <> \
2796  class \
2797  octave_print_conv<T1> \
2798  { \
2799  public: \
2800  typedef T2 print_conv_type; \
2801  }
2802 
2805 
2806 #undef PRINT_CONV
2807 
2808 template <typename T>
2809 static inline void
2810 pr_int (std::ostream& os, const T& d, int fw = 0)
2811 {
2812  std::size_t sz = d.byte_size ();
2813  const unsigned char *tmpi = d.iptr ();
2814 
2815  // Unless explicitly asked for, always print in big-endian
2816  // format for hex and bit formats.
2817  //
2818  // {bit,hex}_format == 1: print big-endian
2819  // {bit,hex}_format == 2: print native
2820 
2821  if (hex_format)
2822  {
2823  octave::preserve_stream_state stream_state (os);
2824 
2825  os.fill ('0');
2826  if (uppercase_format)
2827  os.flags (std::ios::right | std::ios::hex | std::ios::uppercase);
2828  else
2829  os.flags (std::ios::right | std::ios::hex);
2830 
2832  {
2833  for (std::size_t i = 0; i < sz; i++)
2834  os << std::setw (2) << static_cast<int> (tmpi[i]);
2835  }
2836  else
2837  {
2838  for (int i = sz - 1; i >= 0; i--)
2839  os << std::setw (2) << static_cast<int> (tmpi[i]);
2840  }
2841  }
2842  else if (bit_format)
2843  {
2845  {
2846  for (std::size_t i = 0; i < sz; i++)
2847  PRINT_CHAR_BITS (os, tmpi[i]);
2848  }
2849  else
2850  {
2851  if (bit_format > 1)
2852  {
2853  for (std::size_t i = 0; i < sz; i++)
2854  PRINT_CHAR_BITS_SWAPPED (os, tmpi[i]);
2855  }
2856  else
2857  {
2858  for (int i = sz - 1; i >= 0; i--)
2859  PRINT_CHAR_BITS (os, tmpi[i]);
2860  }
2861  }
2862  }
2863  else
2864  {
2865  octave::preserve_stream_state stream_state (os);
2866 
2867  os << std::setw (fw)
2869 
2870  if (bank_format)
2871  os << ".00";
2872  }
2873 }
2874 
2875 template void
2876 pr_int (std::ostream&, const octave_int8&, int);
2877 
2878 template void
2879 pr_int (std::ostream&, const octave_int16&, int);
2880 
2881 template void
2882 pr_int (std::ostream&, const octave_int32&, int);
2883 
2884 template void
2885 pr_int (std::ostream&, const octave_int64&, int);
2886 
2887 template void
2888 pr_int (std::ostream&, const octave_uint8&, int);
2889 
2890 template void
2891 pr_int (std::ostream&, const octave_uint16&, int);
2892 
2893 template void
2894 pr_int (std::ostream&, const octave_uint32&, int);
2895 
2896 template void
2897 pr_int (std::ostream&, const octave_uint64&, int);
2898 
2899 template <typename T>
2900 void
2902  const float_display_format& fmt,
2903  const octave_int<T>& val, bool)
2904 {
2905  if (plus_format)
2906  pr_plus_format (os, val);
2907  else
2908  {
2909  if (free_format)
2910  os << typename octave_print_conv<octave_int<T>>::print_conv_type (val);
2911  else
2912  {
2913  float_format r_fmt = fmt.real_format ();
2914 
2915  pr_int (os, val, r_fmt.width ());
2916  }
2917  }
2918 }
2919 
2920 #define PRINT_INT_SCALAR_INTERNAL(TYPE) \
2921  void \
2922  octave_print_internal (std::ostream& os, \
2923  const float_display_format& fmt, \
2924  const octave_int<TYPE>& val, bool dummy) \
2925  { \
2926  octave_print_internal_template (os, fmt, val, dummy); \
2927  }
2928 
2937 
2938 template <typename T>
2939 static inline void
2941  bool pr_as_read_syntax, int extra_indent)
2942 {
2943  // FIXME: this mostly duplicates the code in the print_nd_array<>
2944  // function. Can fix this with std::is_same from C++11.
2945 
2946  if (nda.isempty ())
2947  print_empty_nd_array (os, nda.dims (), pr_as_read_syntax);
2948  else if (nda.numel () == 1)
2950  pr_as_read_syntax);
2951  else if (plus_format && ! pr_as_read_syntax)
2952  {
2953  int ndims = nda.ndims ();
2954 
2955  Array<octave_idx_type> ra_idx (dim_vector (ndims, 1), 0);
2956 
2957  dim_vector dims = nda.dims ();
2958 
2959  octave_idx_type m = 1;
2960 
2961  for (int i = 2; i < ndims; i++)
2962  m *= dims(i);
2963 
2964  octave_idx_type nr = dims(0);
2965  octave_idx_type nc = dims(1);
2966 
2967  for (octave_idx_type i = 0; i < m; i++)
2968  {
2969  if (m > 1)
2970  {
2971  std::string nm = "ans(:,:,";
2972 
2973  std::ostringstream buf;
2974 
2975  for (int k = 2; k < ndims; k++)
2976  {
2977  buf << ra_idx(k) + 1;
2978 
2979  if (k < ndims - 1)
2980  buf << ',';
2981  else
2982  buf << ')';
2983  }
2984 
2985  nm += buf.str ();
2986 
2987  os << nm << " =\n";
2988  if (! Vcompact_format)
2989  os << "\n";
2990  }
2991 
2992  Array<octave::idx_vector> idx (dim_vector (ndims, 1));
2993 
2994  idx(0) = octave::idx_vector (':');
2995  idx(1) = octave::idx_vector (':');
2996 
2997  for (int k = 2; k < ndims; k++)
2998  idx(k) = octave::idx_vector (ra_idx(k));
2999 
3000  Array<T> page (nda.index (idx), dim_vector (nr, nc));
3001 
3002  for (octave_idx_type ii = 0; ii < nr; ii++)
3003  {
3004  for (octave_idx_type jj = 0; jj < nc; jj++)
3005  {
3006  octave_quit ();
3007 
3008  pr_plus_format (os, page(ii, jj));
3009  }
3010 
3011  if ((ii < nr - 1) || (i < m -1))
3012  os << "\n";
3013  }
3014 
3015  if (i < m - 1)
3016  {
3017  os << "\n";
3018  increment_index (ra_idx, dims, 2);
3019  }
3020  }
3021  }
3022  else
3023  {
3024  int ndims = nda.ndims ();
3025 
3026  dim_vector dims = nda.dims ();
3027 
3028  Array<octave_idx_type> ra_idx (dim_vector (ndims, 1), 0);
3029 
3030  octave_idx_type m = 1;
3031 
3032  for (int i = 2; i < ndims; i++)
3033  m *= dims(i);
3034 
3035  octave_idx_type nr = dims(0);
3036  octave_idx_type nc = dims(1);
3037 
3038  int fw = 0;
3039  if (hex_format)
3040  fw = 2 * nda(0).byte_size ();
3041  else if (bit_format)
3042  fw = nda(0).nbits ();
3043  else
3044  {
3045  bool isneg = false;
3046  int digits = 0;
3047 
3048  for (octave_idx_type i = 0; i < dims.numel (); i++)
3049  {
3050  int new_digits
3051  = static_cast<int>
3052  (std::floor (log10 (double (abs (nda(i).value ()))) + 1));
3053 
3054  if (new_digits > digits)
3055  digits = new_digits;
3056 
3057  if (! isneg)
3058  isneg = (abs (nda(i).value ()) != nda(i).value ());
3059  }
3060 
3061  fw = digits + isneg;
3062  }
3063 
3064  int column_width = fw + (rat_format ? 0 : (bank_format ? 5 : 2));
3065  octave_idx_type total_width = nc * column_width;
3066  int max_width = octave::command_editor::terminal_cols () - extra_indent;
3067  octave_idx_type inc = nc;
3068  if (total_width > max_width && Vsplit_long_rows)
3069  {
3070  inc = max_width / column_width;
3071  if (inc == 0)
3072  inc++;
3073  }
3074 
3075  for (octave_idx_type i = 0; i < m; i++)
3076  {
3077  if (m > 1)
3078  {
3079  std::string nm = "ans(:,:,";
3080 
3081  std::ostringstream buf;
3082 
3083  for (int k = 2; k < ndims; k++)
3084  {
3085  buf << ra_idx(k) + 1;
3086 
3087  if (k < ndims - 1)
3088  buf << ',';
3089  else
3090  buf << ')';
3091  }
3092 
3093  nm += buf.str ();
3094 
3095  os << nm << " =\n";
3096  if (! Vcompact_format)
3097  os << "\n";
3098  }
3099 
3100  Array<octave::idx_vector> idx (dim_vector (ndims, 1));
3101 
3102  idx(0) = octave::idx_vector (':');
3103  idx(1) = octave::idx_vector (':');
3104 
3105  for (int k = 2; k < ndims; k++)
3106  idx(k) = octave::idx_vector (ra_idx(k));
3107 
3108  Array<T> page (nda.index (idx), dim_vector (nr, nc));
3109 
3110  if (free_format)
3111  {
3112  if (pr_as_read_syntax)
3113  os << "[\n";
3114 
3115  for (octave_idx_type ii = 0; ii < nr; ii++)
3116  {
3117  for (octave_idx_type jj = 0; jj < nc; jj++)
3118  {
3119  octave_quit ();
3120  os << " ";
3121  os << typename octave_print_conv<T>::print_conv_type (page(ii, jj));
3122  }
3123  os << "\n";
3124  }
3125 
3126  if (pr_as_read_syntax)
3127  os << ']';
3128  }
3129  else
3130  {
3131  octave::preserve_stream_state stream_state (os);
3132 
3133  octave_idx_type n_rows = page.rows ();
3134  octave_idx_type n_cols = page.cols ();
3135 
3136  for (octave_idx_type col = 0; col < n_cols; col += inc)
3137  {
3138  octave_idx_type lim = (col + inc < n_cols ? col + inc
3139  : n_cols);
3140 
3141  pr_col_num_header (os, total_width, max_width, lim, col,
3142  extra_indent);
3143 
3144  for (octave_idx_type ii = 0; ii < n_rows; ii++)
3145  {
3146  os << std::setw (extra_indent) << "";
3147 
3148  for (octave_idx_type jj = col; jj < lim; jj++)
3149  {
3150  octave_quit ();
3151  os << " ";
3152  pr_int (os, page(ii, jj), fw);
3153  }
3154  if ((ii < n_rows - 1) || (i < m -1))
3155  os << "\n";
3156  }
3157  }
3158  }
3159 
3160  if (i < m - 1)
3161  {
3162  os << "\n";
3163  increment_index (ra_idx, dims, 2);
3164  }
3165  }
3166  }
3167 }
3168 
3169 #define PRINT_INT_ARRAY_INTERNAL(TYPE) \
3170  OCTINTERP_API void \
3171  octave_print_internal (std::ostream& os, const intNDArray<TYPE>& nda, \
3172  bool pr_as_read_syntax, int extra_indent) \
3173  { \
3174  octave_print_internal_template (os, nda, pr_as_read_syntax, extra_indent); \
3175  }
3176 
3185 
3186 void
3187 octave_print_internal (std::ostream&, const Cell&, bool, int, bool)
3188 {
3189  panic_impossible ();
3190 }
3191 
3192 void
3193 octave_print_internal (std::ostream&, const octave_value&, bool)
3194 {
3195  panic_impossible ();
3196 }
3197 
3199 
3200 DEFUN (rats, args, ,
3201  doc: /* -*- texinfo -*-
3202 @deftypefn {} {@var{s} =} rats (@var{x})
3203 @deftypefnx {} {@var{s} =} rats (@var{x}, @var{len})
3204 Convert @var{x} into a rational approximation represented as a string.
3205 
3206 A rational approximation to a floating point number is a simple fraction
3207 with numerator @var{N} and denominator @var{D} such that
3208 @code{@var{x} = @var{N}/@var{D}}.
3209 
3210 The optional second argument defines the maximum length of the string
3211 representing the elements of @var{x}. By default, @var{len} is 9.
3212 
3213 If the length of the smallest possible rational approximation exceeds
3214 @var{len}, an asterisk (*) padded with spaces will be returned instead.
3215 
3216 Example conversion from matrix to string, and back again.
3217 
3218 @example
3219 @group
3220 r = rats (hilb (4));
3221 x = str2num (r)
3222 @end group
3223 @end example
3224 
3225 @seealso{rat, format}
3226 @end deftypefn */)
3227 {
3228  int nargin = args.length ();
3229 
3230  if (nargin < 1 || nargin > 2)
3231  print_usage ();
3232 
3233  octave_value arg = args(0);
3234 
3235  if (! arg.isnumeric ())
3236  error ("rats: X must be numeric");
3237 
3238  if (arg.isempty ())
3239  return ovl (octave_value (""));
3240 
3241  // Convert to N-D arrays to 2-D arrays for Matlab compatibility
3242  if (arg.ndims () > 2)
3243  {
3244  dim_vector dv (arg.rows (), arg.numel () / arg.rows ());
3245  arg = arg.reshape (dv);
3246  }
3247 
3248  unwind_protect frame;
3249 
3250  frame.protect_var (rat_string_len);
3251 
3252  rat_string_len = 9;
3253  if (nargin == 2)
3254  rat_string_len = args(1).nint_value ();
3255 
3256  frame.protect_var (rat_format);
3257 
3258  rat_format = true;
3259 
3260  std::ostringstream buf;
3261  arg.print (buf);
3262  std::string s = buf.str ();
3263 
3264  std::list<std::string> lst;
3265 
3266  std::size_t n = 0;
3267  std::size_t s_len = s.length ();
3268 
3269  while (n < s_len)
3270  {
3271  std::size_t m = s.find ('\n', n);
3272 
3273  if (m == std::string::npos)
3274  {
3275  lst.push_back (s.substr (n));
3276  break;
3277  }
3278  else
3279  {
3280  lst.push_back (s.substr (n, m - n));
3281  n = m + 1;
3282  }
3283  }
3284 
3285  return ovl (string_vector (lst));
3286 }
3287 
3288 /*
3289 %!test <*56941>
3290 %! [old_fmt, old_spacing] = format ();
3291 %! unwind_protect
3292 %! format short;
3293 %! assert (rats (-2.0005, 10), "-4001/2000");
3294 %! assert (strtrim (rats (2.0005, 30)), "4001/2000");
3295 %! assert (pi - str2num (rats (pi, 30)), 0, 4 * eps);
3296 %! assert (e - str2num (rats (e, 30)), 0, 4 * eps);
3297 %! unwind_protect_cleanup
3298 %! format (old_fmt);
3299 %! format (old_spacing);
3300 %! end_unwind_protect
3301 
3302 %!test <*57003>
3303 %! x = ones (2,1,3);
3304 %! s = rats (x,4);
3305 %! assert (ndims (s) == 2);
3306 %! assert (rows (s) == 2);
3307 %! assert (columns (s) == 3 * 6);
3308 
3309 %!assert <*57004> (rats ([]), '')
3310 
3311 %!test <57704>
3312 %! [old_fmt, old_spacing] = format ();
3313 %! unwind_protect
3314 %! format short;
3315 %! assert (rats (2.0005, 9), "4001/2000");
3316 %! assert (rats (123, 2), " *");
3317 %! unwind_protect_cleanup
3318 %! format (old_fmt);
3319 %! format (old_spacing);
3320 %! end_unwind_protect
3321 
3322 */
3323 
3324 DEFUN (disp, args, nargout,
3325  classes: cell char double function_handle int8 int16 int32 int64 logical single struct uint8 uint16 uint32 uint64
3326  doc: /* -*- texinfo -*-
3327 @deftypefn {} {} disp (@var{x})
3328 @deftypefnx {} {@var{str} =} disp (@var{x})
3329 Display the value of @var{x}.
3330 
3331 For example:
3332 
3333 @example
3334 @group
3335 disp ("The value of pi is:"), disp (pi)
3336 
3337  @print{} the value of pi is:
3338  @print{} 3.1416
3339 @end group
3340 @end example
3341 
3342 @noindent
3343 Note that the output from @code{disp} always ends with a newline.
3344 
3345 If an output value is requested, @code{disp} prints nothing and returns the
3346 formatted output in a string.
3347 @seealso{fdisp}
3348 @end deftypefn */)
3349 {
3350  if (args.length () != 1)
3351  print_usage ();
3352 
3353  octave_value_list retval;
3354 
3355  octave_value arg = args(0);
3356 
3357  if (nargout == 0)
3358  arg.print (octave_stdout);
3359  else
3360  {
3361  std::ostringstream buf;
3362  arg.print (buf);
3363  retval = (octave_value (buf.str (), arg.is_dq_string () ? '"' : '\''));
3364  }
3365 
3366  return retval;
3367 }
3368 
3369 DEFMETHOD (fdisp, interp, args, ,
3370  classes: cell char double function_handle int8 int16 int32 int64 logical single struct uint8 uint16 uint32 uint64
3371  doc: /* -*- texinfo -*-
3372 @deftypefn {} {} fdisp (@var{fid}, @var{x})
3373 Display the value of @var{x} on the stream @var{fid}.
3374 
3375 For example:
3376 
3377 @example
3378 @group
3379 fdisp (stdout, "The value of pi is:"), fdisp (stdout, pi)
3380 
3381  @print{} the value of pi is:
3382  @print{} 3.1416
3383 @end group
3384 @end example
3385 
3386 @noindent
3387 Note that the output from @code{fdisp} always ends with a newline.
3388 @seealso{disp}
3389 @end deftypefn */)
3390 {
3391  if (args.length () != 2)
3392  print_usage ();
3393 
3394  stream_list& streams = interp.get_stream_list ();
3395 
3396  int fid = streams.get_file_number (args(0));
3397 
3398  stream os = streams.lookup (fid, "fdisp");
3399 
3400  std::ostream *osp = os.output_stream ();
3401 
3402  if (! osp)
3403  error ("fdisp: stream FID not open for writing");
3404 
3405  octave_value arg = args(1);
3406  arg.print (*osp);
3407 
3408  return ovl ();
3409 }
3410 
3411 /*
3412 ## FIXME: This test writes values to a file, but then never checks them.
3413 %!test
3414 %! [old_fmt, old_spacing] = format ();
3415 %! unwind_protect
3416 %! format short
3417 %! fd = tmpfile ();
3418 %! for r = [0, Inf -Inf, NaN]
3419 %! for i = [0, Inf -Inf, NaN]
3420 %! fdisp (fd, complex (r, i));
3421 %! endfor
3422 %! endfor
3423 %! fclose (fd);
3424 %! unwind_protect_cleanup
3425 %! format (old_fmt);
3426 %! format (old_spacing);
3427 %! end_unwind_protect
3428 
3429 %!test
3430 %! [old_fmt, old_spacing] = format ();
3431 %! unwind_protect
3432 %! foo.real = pi * ones (3,20,3);
3433 %! foo.complex = pi * ones (3,20,3) + 1i;
3434 %! foo.char = repmat ("- Hello World -", [3, 20]);
3435 %! foo.cell = {foo.real, foo.complex, foo.char};
3436 %! fields = fieldnames (foo);
3437 %! for f = 1:numel (fields)
3438 %! format loose;
3439 %! loose = disp (foo.(fields{f}));
3440 %! format compact;
3441 %! compact = disp (foo.(fields{f}));
3442 %! expected = strrep (loose, "\n\n", "\n");
3443 %! assert (expected, compact);
3444 %! endfor
3445 %! unwind_protect_cleanup
3446 %! format (old_fmt);
3447 %! format (old_spacing);
3448 %! end_unwind_protect
3449 */
3450 
3451 DEFUN (display, args, ,
3452  classes: cell char double function_handle int8 int16 int32 int64 logical single struct uint8 uint16 uint32 uint64
3453  doc: /* -*- texinfo -*-
3454 @deftypefn {} {} display (@var{obj})
3455 Display the contents of the object @var{obj} prepended by its name.
3456 
3457 The Octave interpreter calls the @code{display} function whenever it needs
3458 to present a class on-screen. Typically, this would be a statement which
3459 does not end in a semicolon to suppress output. For example:
3460 
3461 @example
3462 myclass (@dots{})
3463 @end example
3464 
3465 Or:
3466 
3467 @example
3468 myobj = myclass (@dots{})
3469 @end example
3470 
3471 In general, user-defined classes should overload the @code{disp} method to
3472 avoid the default output:
3473 
3474 @example
3475 @group
3476 myobj = myclass (@dots{})
3477  @result{} myobj =
3478 
3479  <class myclass>
3480 @end group
3481 @end example
3482 
3483 When overloading the @code{display} method instead, one has to take care
3484 of properly displaying the object's name. This can be done by using the
3485 @code{inputname} function.
3486 
3487 @seealso{disp, class, subsref, subsasgn}
3488 @end deftypefn */)
3489 {
3490  int nargin = args.length ();
3491 
3492  // Matlab apparently accepts two arguments with the second set to the
3493  // inputname of the first. This is undocumented, but we'll use it.
3494  // However, we never call display methods for classes with more than
3495  // one argument.
3496 
3497  if (nargin < 1 || nargin > 2)
3498  print_usage ();
3499 
3500  std::string name;
3501 
3502  if (nargin == 2)
3503  name = args(1).xstring_value ("NAME must be a string");
3504  else
3505  {
3506  string_vector names = args.name_tags ();
3507  name = names(0);
3508  }
3509 
3510  // We are here because there is no overloaded display method for this
3511  // object type.
3512 
3513  octave_value value = args(0);
3514 
3515  // If print_name_tag displays a newline, then also print one after
3516  // disp is done.
3517 
3518  bool print_newlines = false;
3519  if (valid_identifier (name))
3520  print_newlines = value.print_name_tag (octave_stdout, name);
3521 
3522  // Use feval so that dispatch will also work for disp.
3523 
3524  feval ("disp", ovl (value));
3525 
3526  if (print_newlines)
3527  octave_stdout << std::endl;
3528 
3529  return ovl ();
3530 }
3531 
3532 /*
3533 %!test
3534 %! [old_fmt, old_spacing] = format ();
3535 %! unwind_protect
3536 %! format short;
3537 %! str = evalc ("x = 1.1; display (x)");
3538 %! assert (str, "x = 1.1000\n");
3539 %! unwind_protect_cleanup
3540 %! format (old_fmt);
3541 %! format (old_spacing);
3542 %! end_unwind_protect
3543 
3544 %!test
3545 %! [old_fmt, old_spacing] = format ();
3546 %! unwind_protect
3547 %! format short;
3548 %! str = evalc ("display (1.1)");
3549 %! assert (str, "1.1000\n");
3550 %! unwind_protect_cleanup
3551 %! format (old_fmt);
3552 %! format (old_spacing);
3553 %! end_unwind_protect
3554 
3555 ## Test input validation
3556 %!error display ()
3557 %!error display (1,2)
3558 */
3559 
3560 static inline void
3562 {
3563  free_format = false;
3564  plus_format = false;
3565  rat_format = false;
3566  bank_format = false;
3567  hex_format = 0;
3568  bit_format = 0;
3569  print_e = false;
3570  print_g = false;
3571  print_eng = false;
3572 }
3573 
3574 static std::string format_string ("short");
3575 
3576 static inline void
3577 set_format_style (int argc, const string_vector& argv)
3578 {
3579  if (--argc == 0)
3580  {
3581  // Special case of no options, reset to default values
3582  init_format_state ();
3583  set_output_prec (5);
3584  format_string = "short";
3585  Vcompact_format = false;
3586  uppercase_format = false;
3587  return;
3588  }
3589 
3590  int idx = 1;
3591  std::string format;
3592 
3593  unwind_protect frame;
3594 
3595  frame.protect_var (bank_format);
3596  frame.protect_var (bit_format);
3597  frame.protect_var (free_format);
3598  frame.protect_var (hex_format);
3599  frame.protect_var (plus_format);
3601  frame.protect_var (rat_format);
3602  frame.protect_var (print_e);
3603  frame.protect_var (print_eng);
3604  frame.protect_var (print_g);
3605  frame.protect_var (format_string);
3606  frame.protect_var (Vcompact_format);
3607  frame.protect_var (uppercase_format);
3608  int prec = output_precision ();
3609  frame.add ([=] (void) { set_output_prec (prec); });
3610 
3611  format = format_string; // Initialize with existing value
3612  while (argc-- > 0)
3613  {
3614  std::string arg = argv[idx++];
3615  std::transform (arg.begin (), arg.end (), arg.begin (), tolower);
3616 
3617  if (arg == "default")
3618  {
3619  format = "short";
3620  init_format_state ();
3621  set_output_prec (5);
3622  Vcompact_format = false;
3623  uppercase_format = false;
3624  }
3625  else if (arg == "short")
3626  {
3627  format = arg;
3628  init_format_state ();
3629  if (argc > 0)
3630  {
3631  arg = argv[idx];
3632  if (arg == "e")
3633  {
3634  print_e = true;
3635  format.append (arg);
3636  argc--; idx++;
3637  }
3638  else if (arg == "g")
3639  {
3640  init_format_state ();
3641  print_g = true;
3642  format.append (arg);
3643  argc--; idx++;
3644  }
3645  else if (arg == "eng")
3646  {
3647  init_format_state ();
3648  print_eng = true;
3649  format.append (arg);
3650  argc--; idx++;
3651  }
3652  }
3653  set_output_prec (5);
3654  }
3655  else if (arg == "shorte")
3656  {
3657  format = arg;
3658  init_format_state ();
3659  print_e = true;
3660  set_output_prec (5);
3661  }
3662  else if (arg == "shortg")
3663  {
3664  format = arg;
3665  init_format_state ();
3666  print_g = true;
3667  set_output_prec (5);
3668  }
3669  else if (arg == "shorteng")
3670  {
3671  format = arg;
3672  init_format_state ();
3673  print_eng = true;
3674  set_output_prec (5);
3675  }
3676  else if (arg == "long")
3677  {
3678  format = arg;
3679  init_format_state ();
3680  if (argc > 0)
3681  {
3682  arg = argv[idx];
3683 
3684  if (arg == "e")
3685  {
3686  print_e = true;
3687  format.append (arg);
3688  argc--; idx++;
3689  }
3690  else if (arg == "g")
3691  {
3692  print_g = true;
3693  format.append (arg);
3694  argc--; idx++;
3695  }
3696  else if (arg == "eng")
3697  {
3698  print_eng = true;
3699  format.append (arg);
3700  argc--; idx++;
3701  }
3702  }
3703  set_output_prec (16);
3704  }
3705  else if (arg == "longe")
3706  {
3707  format = arg;
3708  init_format_state ();
3709  print_e = true;
3710  set_output_prec (16);
3711  }
3712  else if (arg == "longg")
3713  {
3714  format = arg;
3715  init_format_state ();
3716  print_g = true;
3717  set_output_prec (16);
3718  }
3719  else if (arg == "longeng")
3720  {
3721  format = arg;
3722  init_format_state ();
3723  print_eng = true;
3724  set_output_prec (16);
3725  }
3726  else if (arg == "hex")
3727  {
3728  format = arg;
3729  init_format_state ();
3730  hex_format = 1;
3731  }
3732  else if (arg == "native-hex")
3733  {
3734  format = arg;
3735  init_format_state ();
3736  hex_format = 2;
3737  }
3738  else if (arg == "bit")
3739  {
3740  format = arg;
3741  init_format_state ();
3742  bit_format = 1;
3743  }
3744  else if (arg == "native-bit")
3745  {
3746  format = arg;
3747  init_format_state ();
3748  bit_format = 2;
3749  }
3750  else if (arg == "+" || arg == "plus")
3751  {
3752  format = arg;
3753  init_format_state ();
3754  plus_format = true;
3755  if (argc > 0)
3756  {
3757  arg = argv[idx];
3758 
3759  if (arg.length () == 3)
3760  {
3761  plus_format_chars = arg;
3762  format.append (arg);
3763  argc--; idx++;
3764  }
3765  else
3766  plus_format_chars = "+- ";
3767  }
3768  else
3769  plus_format_chars = "+- ";
3770  }
3771  else if (arg == "rat")
3772  {
3773  format = arg;
3774  init_format_state ();
3775  rat_format = true;
3776  }
3777  else if (arg == "bank")
3778  {
3779  format = arg;
3780  init_format_state ();
3781  bank_format = true;
3782  }
3783  else if (arg == "free")
3784  {
3785  format = arg;
3786  init_format_state ();
3787  free_format = true;
3788  }
3789  else if (arg == "none")
3790  {
3791  format = arg;
3792  init_format_state ();
3793  free_format = true;
3794  }
3795  else if (arg == "compact")
3796  Vcompact_format = true;
3797  else if (arg == "loose")
3798  Vcompact_format = false;
3799  else if (arg == "lowercase")
3800  uppercase_format = false;
3801  else if (arg == "uppercase")
3802  uppercase_format = true;
3803  else
3804  error ("format: unrecognized format state '%s'", arg.c_str ());
3805  }
3806 
3808 
3809  // If successful, discard unwind state information
3810  frame.discard ();
3811 }
3812 
3813 DEFUN (format, args, nargout,
3814  doc: /* -*- texinfo -*-
3815 @deftypefn {} {} format
3816 @deftypefnx {} {} format options
3817 @deftypefnx {} {} format (@var{options})
3818 @deftypefnx {} {[@var{format}, @var{formatspacing}, @var{uppercase}] =} format
3819 Reset or specify the format of the output produced by @code{disp} and Octave's
3820 normal echoing mechanism.
3821 
3822 This command only affects the display of numbers, not how they are stored
3823 or computed. To change the internal representation from the default double use
3824 one of the conversion functions such as @code{single}, @code{uint8},
3825 @code{int64}, etc. Any @code{format} options that change the number of
3826 displayed significant digits will also be reflected by the
3827 @code{output_precision} function.
3828 
3829 By default, Octave displays 5 significant digits in a human readable form
3830 (option @samp{short}, option @samp{lowercase}, and option @samp{loose} format
3831 for matrices). If @code{format} is invoked without any options, or the option
3832 @samp{default} is specified, then this default format is restored.
3833 
3834 Valid format options for floating point numbers are listed in the following
3835 table.
3836 
3837 @table @code
3838 @item default
3839 Restore the default format state described above.
3840 
3841 @item short
3842 Fixed point format with 5 significant figures (default).
3843 
3844 @item long
3845 Fixed point format with 16 significant figures.
3846 
3847 As with the @samp{short} format, Octave will switch to an exponential @samp{e}
3848 format if it is unable to format a matrix properly using the current format.
3849 
3850 @item shorte
3851 @itemx longe
3852 Exponential format. The number to be represented is split between a mantissa
3853 and an exponent (power of 10). The mantissa has 5 significant digits in the
3854 short format. In the long format, double values are displayed with 16
3855 significant digits and single values are displayed with 8. For example,
3856 with the @samp{shorte} format, @code{pi} is displayed as @code{3.1416e+00}.
3857 Optionally, the trailing @samp{e} can be split into a second argument.
3858 
3859 @item shortg
3860 @itemx longg
3861 Optimally choose between fixed point and exponential format based on the
3862 magnitude of the number. For example, with the @samp{shortg} format,
3863 @code{pi .^ [2; 4; 8; 16; 32]} is displayed as
3864 
3865 @example
3866 @group
3867 ans =
3868 
3869  9.8696
3870  97.409
3871  9488.5
3872  9.0032e+07
3873  8.1058e+15
3874 @end group
3875 @end example
3876 
3877 Optionally, the trailing @samp{g} can be split into a second argument.
3878 
3879 @item shorteng
3880 @itemx longeng
3881 Identical to @samp{shorte} or @samp{longe} but displays the value using an
3882 engineering format, where the exponent is divisible by 3. For example, with
3883 the @samp{shorteng} format, @code{10 * pi} is displayed as @code{31.416e+00}.
3884 Optionally, the trailing @samp{eng} can be split into a second argument.
3885 
3886 @item free
3887 @itemx none
3888 Print output in free format, without trying to line up columns of matrices on
3889 the decimal point. This is a raw format equivalent to the C++ code
3890 @code{std::cout << @var{variable}}. In general, the result is a presentation
3891 with 6 significant digits where unnecessary precision (such as trailing zeros
3892 for integers) is suppressed. Complex numbers are formatted as numeric pairs
3893 like this @samp{(0.60419, 0.60709)} instead of like this
3894 @samp{0.60419 + 0.60709i}.
3895 @end table
3896 
3897 The following formats affect all numeric output (floating point and integer
3898 types).
3899 
3900 @table @asis
3901 @item @qcode{"+"}
3902 @itemx @qcode{"+"} @qcode{"@var{chars}"}
3903 @itemx @code{plus}
3904 @itemx @code{plus @var{chars}}
3905 Print a @samp{+} symbol for matrix elements greater than zero, a @samp{-}
3906 symbol for elements less than zero, and a space for zero matrix elements. This
3907 format can be useful for examining the sparsity structure of a large matrix.
3908 For very large matrices the function @code{spy} which plots the sparsity
3909 pattern will be clearer.
3910 
3911 The optional argument @var{chars} specifies a list of 3 characters to use for
3912 printing values greater than zero, less than zero, and equal to zero. For
3913 example, with the format @qcode{"+" "+-."}, the matrix
3914 @code{[1, 0, -1; -1, 0, 1]} is displayed as
3915 
3916 @example
3917 @group
3918 ans =
3919 
3920 +.-
3921 -.+
3922 @end group
3923 @end example
3924 
3925 @item bank
3926 Print variable in a format appropriate for a currency (fixed format with two
3927 digits to the right of the decimal point). Only the real part of a variable is
3928 displayed, as the imaginary part makes no sense for a currency.
3929 
3930 @item native-hex
3931 Print the hexadecimal representation of numbers as they are stored in memory.
3932 For example, on a workstation which stores 8 byte real values in IEEE format
3933 with the least significant byte first, the value of @code{pi} when printed in
3934 @code{native-hex} format is @code{400921fb54442d18}.
3935 
3936 @item hex
3937 The same as @code{native-hex}, but always print the most significant byte
3938 first.
3939 
3940 @item native-bit
3941 Print the bit representation of numbers as stored in memory. For example, the
3942 value of @code{pi} is
3943 
3944 @example
3945 @group
3946 01000000000010010010000111111011
3947 01010100010001000010110100011000
3948 @end group
3949 @end example
3950 
3951 (shown here in two 32 bit sections for typesetting purposes) when printed in
3952 native-bit format on a workstation which stores 8 byte real values in IEEE
3953 format with the least significant byte first.
3954 
3955 @item bit
3956 The same as @code{native-bit}, but always print the most significant bits
3957 first.
3958 
3959 @item rat
3960 Print a rational approximation, i.e., values are approximated as the ratio of
3961 small integers. For example, with the @samp{rat} format, @code{pi} is
3962 displayed as @code{355/113}.
3963 @end table
3964 
3965 The following two options affect the display of scientific and hex notations.
3966 
3967 @table @code
3968 @item lowercase (default)
3969 Use a lowercase @samp{e} for the exponent character in scientific notation and
3970 lowercase @samp{a-f} for the hex digits representing 10-15.
3971 
3972 @item uppercase
3973 Use an uppercase @samp{E} for the exponent character in scientific notation and
3974 uppercase @samp{A-F} for the hex digits representing 10-15.
3975 @end table
3976 
3977 The following two options affect the display of all matrices.
3978 
3979 @table @code
3980 @item compact
3981 Remove blank lines around column number labels and between matrices producing
3982 more compact output with more data per page.
3983 
3984 @item loose (default)
3985 Insert blank lines above and below column number labels and between matrices to
3986 produce a more readable output with less data per page.
3987 @end table
3988 
3989 If @code{format} is called with multiple competing options, the rightmost one
3990 is used, except for @samp{default} which will override all other options. In
3991 case of an error the format remains unchanged.
3992 
3993 If called with one to three output arguments, and no inputs, return the current
3994 format, format spacing, and uppercase preference. Specifying both outputs and
3995 inputs will produce an error.
3996 
3997 @seealso{fixed_point_format, output_precision, split_long_rows,
3998 print_empty_dimensions, rats}
3999 @end deftypefn */)
4000 {
4001  octave_value_list retval (std::min (nargout, 2));
4002 
4003  int nargin = args.length ();
4004 
4005  if (nargout == 0)
4006  {
4007  int argc = nargin + 1;
4008 
4009  string_vector argv = args.make_argv ("format");
4010 
4011  set_format_style (argc, argv);
4012  }
4013  else
4014  {
4015  if (nargin > 0)
4016  warning ("format: cannot query and set format at the same time, ignoring set operation");
4017 
4018  if (nargout >= 3)
4019  retval(2) = (uppercase_format ? "uppercase" : "lowercase");
4020 
4021  if (nargout >= 2)
4022  retval(1) = (Vcompact_format ? "compact" : "loose");
4023 
4024  retval(0) = format_string;
4025  }
4026 
4027  return retval;
4028 }
4029 
4030 /*
4031 %!test
4032 %! [old_fmt, old_spacing, old_uppercase] = format ();
4033 %! unwind_protect
4034 %! ## Test one of the formats
4035 %! format long e;
4036 %! format uppercase;
4037 %! str = disp (pi);
4038 %! assert (str, "3.141592653589793E+00\n");
4039 %! str = disp (single (pi));
4040 %! assert (str, "3.1415927E+00\n");
4041 %! new_fmt = format ();
4042 %! assert (new_fmt, "longe");
4043 %! ## Test resetting format (method #1)
4044 %! format compact;
4045 %! [~, new_spacing] = format ();
4046 %! assert (new_spacing, "compact");
4047 %! format;
4048 %! [new_fmt, new_spacing, new_case] = format ();
4049 %! assert (new_fmt, "short");
4050 %! assert (new_spacing, "loose");
4051 %! assert (new_case, "lowercase");
4052 %! ## Test resetting format (method #2)
4053 %! format compact uppercase long e;
4054 %! [new_fmt, new_spacing, new_case] = format ();
4055 %! assert (new_fmt, "longe");
4056 %! assert (new_spacing, "compact");
4057 %! assert (new_case, "uppercase");
4058 %! format ("default");
4059 %! [new_fmt, new_spacing, new_case] = format ();
4060 %! assert (new_fmt, "short");
4061 %! assert (new_spacing, "loose");
4062 %! assert (new_case, "lowercase");
4063 %! unwind_protect_cleanup
4064 %! format (old_fmt);
4065 %! format (old_spacing);
4066 %! format (old_uppercase);
4067 %! end_unwind_protect
4068 
4069 %!test <*53427>
4070 %! [old_fmt, old_spacing] = format ();
4071 %! unwind_protect
4072 %! format; # reset format to short and loose
4073 %! format compact; # set compact format
4074 %! format long; # set long format
4075 %! [fmt, spacing] = format ();
4076 %! assert (fmt, "long");
4077 %! assert (spacing, "compact");
4078 %! unwind_protect_cleanup
4079 %! format (old_fmt);
4080 %! format (old_spacing);
4081 %! end_unwind_protect
4082 
4083 ## Test input validation
4084 %!test
4085 %! fail ("fmt = format ('long')", "warning", "cannot query and set format");
4086 
4087 */
4088 
4089 DEFUN (fixed_point_format, args, nargout,
4090  doc: /* -*- texinfo -*-
4091 @deftypefn {} {@var{val} =} fixed_point_format ()
4092 @deftypefnx {} {@var{old_val} =} fixed_point_format (@var{new_val})
4093 @deftypefnx {} {@var{old_val} =} fixed_point_format (@var{new_val}, "local")
4094 Query or set the internal variable that controls whether Octave will
4095 use a scaled format to print matrix values.
4096 
4097 The scaled format prints a scaling factor on the first line of output chosen
4098 such that the largest matrix element can be written with a single leading
4099 digit. For example:
4100 
4101 @example
4102 @group
4103 fixed_point_format (true)
4104 logspace (1, 7, 5)'
4105 ans =
4106 
4107  1.0e+07 *
4108 
4109  0.00000
4110  0.00003
4111  0.00100
4112  0.03162
4113  1.00000
4114 @end group
4115 @end example
4116 
4117 @noindent
4118 Notice that the first value appears to be 0 when it is actually 1. Because
4119 of the possibility for confusion you should be careful about enabling
4120 @code{fixed_point_format}.
4121 
4122 When called from inside a function with the @qcode{"local"} option, the
4123 variable is changed locally for the function and any subroutines it calls.
4124 The original variable value is restored when exiting the function.
4125 @seealso{format, output_precision}
4126 @end deftypefn */)
4127 {
4128  return set_internal_variable (Vfixed_point_format, args, nargout,
4129  "fixed_point_format");
4130 }
4131 
4132 DEFUN (print_empty_dimensions, args, nargout,
4133  doc: /* -*- texinfo -*-
4134 @deftypefn {} {@var{val} =} print_empty_dimensions ()
4135 @deftypefnx {} {@var{old_val} =} print_empty_dimensions (@var{new_val})
4136 @deftypefnx {} {@var{old_val} =} print_empty_dimensions (@var{new_val}, "local")
4137 Query or set the internal variable that controls whether the dimensions of
4138 empty matrices are printed along with the empty matrix symbol, @samp{[]}.
4139 
4140 For example, the expression
4141 
4142 @example
4143 zeros (3, 0)
4144 @end example
4145 
4146 @noindent
4147 will print
4148 
4149 @example
4150 ans = [](3x0)
4151 @end example
4152 
4153 When called from inside a function with the @qcode{"local"} option, the
4154 variable is changed locally for the function and any subroutines it calls.
4155 The original variable value is restored when exiting the function.
4156 @seealso{format}
4157 @end deftypefn */)
4158 {
4159  return set_internal_variable (Vprint_empty_dimensions, args, nargout,
4160  "print_empty_dimensions");
4161 }
4162 
4163 DEFUN (split_long_rows, args, nargout,
4164  doc: /* -*- texinfo -*-
4165 @deftypefn {} {@var{val} =} split_long_rows ()
4166 @deftypefnx {} {@var{old_val} =} split_long_rows (@var{new_val})
4167 @deftypefnx {} {@var{old_val} =} split_long_rows (@var{new_val}, "local")
4168 Query or set the internal variable that controls whether rows of a matrix
4169 may be split when displayed to a terminal window.
4170 
4171 If the rows are split, Octave will display the matrix in a series of smaller
4172 pieces, each of which can fit within the limits of your terminal width and
4173 each set of rows is labeled so that you can easily see which columns are
4174 currently being displayed. For example:
4175 
4176 @example
4177 @group
4178 octave:13> rand (2,10)
4179 ans =
4180 
4181  Columns 1 through 6:
4182 
4183  0.75883 0.93290 0.40064 0.43818 0.94958 0.16467
4184  0.75697 0.51942 0.40031 0.61784 0.92309 0.40201
4185 
4186  Columns 7 through 10:
4187 
4188  0.90174 0.11854 0.72313 0.73326
4189  0.44672 0.94303 0.56564 0.82150
4190 @end group
4191 @end example
4192 
4193 When called from inside a function with the @qcode{"local"} option, the
4194 variable is changed locally for the function and any subroutines it calls.
4195 The original variable value is restored when exiting the function.
4196 @seealso{format}
4197 @end deftypefn */)
4198 {
4199  return set_internal_variable (Vsplit_long_rows, args, nargout,
4200  "split_long_rows");
4201 }
4202 
void increment_index(Array< octave_idx_type > &ra_idx, const dim_vector &dimensions, int start_dimension)
Definition: Array-util.cc:60
OCTAVE_END_NAMESPACE(octave)
charNDArray max(char d, const charNDArray &m)
Definition: chNDArray.cc:230
charNDArray min(char d, const charNDArray &m)
Definition: chNDArray.cc:207
OCTARRAY_API Array< T, Alloc > index(const octave::idx_vector &i) const
Indexing without resizing.
Definition: Array-base.cc:719
OCTARRAY_OVERRIDABLE_FUNC_API bool isempty(void) const
Size of the specified dimension.
Definition: Array.h:651
OCTARRAY_OVERRIDABLE_FUNC_API octave_idx_type numel(void) const
Number of elements in the array.
Definition: Array.h:414
OCTARRAY_OVERRIDABLE_FUNC_API const dim_vector & dims(void) const
Return a const-reference so that dims ()(i) works efficiently.
Definition: Array.h:503
OCTARRAY_OVERRIDABLE_FUNC_API octave_idx_type rows(void) const
Definition: Array.h:459
OCTARRAY_OVERRIDABLE_FUNC_API int ndims(void) const
Size of the specified dimension.
Definition: Array.h:677
OCTARRAY_OVERRIDABLE_FUNC_API octave_idx_type cols(void) const
Definition: Array.h:469
Definition: Cell.h:43
Template for N-dimensional array classes with like-type math operators.
Definition: MArray.h:63
Definition: dMatrix.h:42
void add(F &&fcn, Args &&... args)
void discard(std::size_t num)
void protect_var(T &var)
OCTAVE_API 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:94
OCTAVE_API 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:335
bool any_zero(void) const
Definition: dim-vector.h:316
float_format real_format(void) const
Definition: pr-flt-fmt.h:234
double scale_factor(void) const
Definition: pr-flt-fmt.h:232
float_format imag_format(void) const
Definition: pr-flt-fmt.h:236
float_format & exponent_width(int w)
Definition: pr-flt-fmt.h:107
std::ios::fmtflags format_flags(void) const
Definition: pr-flt-fmt.h:120
float_format & precision(int p)
Definition: pr-flt-fmt.h:95
float_format & width(int w)
Definition: pr-flt-fmt.h:101
float_format & uppercase(void)
Definition: pr-flt-fmt.h:83
void print_raw(std::ostream &os, bool pr_as_read_syntax=false) const
Definition: ov.h:1431
void print_with_name(std::ostream &os, const std::string &name) const
Definition: ov.h:1437
octave_idx_type rows(void) const
Definition: ov.h:590
bool isnumeric(void) const
Definition: ov.h:795
octave_idx_type numel(void) const
Definition: ov.h:604
bool is_dq_string(void) const
Definition: ov.h:688
int ndims(void) const
Definition: ov.h:596
octave_value reshape(const dim_vector &dv) const
Definition: ov.h:616
bool isempty(void) const
Definition: ov.h:646
bool print_name_tag(std::ostream &os, const std::string &name) const
Definition: ov.h:1434
void print(std::ostream &os, bool pr_as_read_syntax=false)
Definition: ov.h:1428
T mantissa(void) const
Definition: pr-output.cc:180
const float_format m_ff
Definition: pr-output.h:516
int exponent(void) const
Definition: pr-output.cc:173
const float_format m_ff
Definition: pr-output.h:537
const float_format m_ff
Definition: pr-output.h:554
OCTINTERP_API int get_file_number(const octave_value &fid) const
Definition: oct-stream.cc:7673
OCTINTERP_API stream lookup(int fid, const std::string &who="") const
Definition: oct-stream.cc:7460
std::ostream * output_stream(void)
Definition: oct-stream.h:457
ColumnVector real(const ComplexColumnVector &a)
Definition: dColVector.cc:137
ColumnVector imag(const ComplexColumnVector &a)
Definition: dColVector.cc:143
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
OCTINTERP_API void print_usage(void)
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:1054
void error(const char *fmt,...)
Definition: error.cc:979
#define panic_impossible()
Definition: error.h:508
void error_unless(bool cond)
Definition: error.h:549
void scale(Matrix &m, double x, double y, double z)
Definition: graphics.cc:5851
ColumnVector transform(const Matrix &m, double x, double y, double z)
Definition: graphics.cc:5819
octave::idx_vector idx_vector
Definition: idx-vector.h:1039
#define lo_ieee_signbit(x)
Definition: lo-ieee.h:120
bool isna(double x)
Definition: lo-mappers.cc:47
bool isfinite(double x)
Definition: lo-mappers.h:192
bool isinf(double x)
Definition: lo-mappers.h:203
bool isnan(bool)
Definition: lo-mappers.h:178
std::complex< T > floor(const std::complex< T > &x)
Definition: lo-mappers.h:130
T x_nint(T x)
Definition: lo-mappers.h:269
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
bool words_big_endian(void)
Definition: mach-info.cc:72
float_format native_float_format(void)
Definition: mach-info.cc:65
float_format
Definition: mach-info.h:38
@ flt_fmt_ieee_big_endian
Definition: mach-info.h:44
class OCTAVE_API boolNDArray
Definition: mx-fwd.h:42
T octave_idx_type m
Definition: mx-inlines.cc:773
octave_idx_type n
Definition: mx-inlines.cc:753
T * r
Definition: mx-inlines.cc:773
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)
octave_value_list feval(const char *name, const octave_value_list &args, int nargout)
Evaluate an Octave function (built-in or interpreted) and return the list of result values.
Definition: oct-parse.cc:10370
std::string rational_approx(T val, int len)
Definition: oct-string.cc:696
const octave_base_value const Array< octave_idx_type > & ra_idx
return octave_value(v1.char_array_value() . concat(v2.char_array_value(), ra_idx),((a1.is_sq_string()||a2.is_sq_string()) ? '\'' :'"'))
octave_value_list ovl(const OV_Args &... args)
Construct an octave_value_list with less typing.
Definition: ovl.h:211
#define octave_stdout
Definition: pager.h:314
int output_precision(void)
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
static float_display_format make_range_format(int x_max, int x_min, int all_ints)
Definition: pr-output.cc:1191
static void init_format_state(void)
Definition: pr-output.cc:3561
static T pr_max_internal(const MArray< T > &m)
Definition: pr-output.cc:300
static void pr_col_num_header(std::ostream &os, octave_idx_type total_width, int max_width, octave_idx_type lim, octave_idx_type col, int extra_indent)
Definition: pr-output.cc:1627
#define PRINT_INT_SCALAR_INTERNAL(TYPE)
Definition: pr-output.cc:2920
static bool print_e
Definition: pr-output.cc:105
#define MAKE_INT_SCALAR_FORMAT(TYPE)
Definition: pr-output.cc:1737
static void print_empty_matrix(std::ostream &os, octave_idx_type nr, octave_idx_type nc, bool pr_as_read_syntax)
Definition: pr-output.cc:1570
#define PRINT_CONV(T1, T2)
Definition: pr-output.cc:2794
void octave_print_internal_template(std::ostream &os, const float_display_format &fmt, const octave_int< T > &val, bool)
Definition: pr-output.cc:2901
#define PRINT_CHAR_BITS(os, c)
Definition: pr-output.cc:1356
#define SPECIALIZE_UABS(T)
Definition: pr-output.cc:1691
void octave_print_internal(std::ostream &os, const float_display_format &fmt, bool d, bool pr_as_read_syntax)
Definition: pr-output.cc:1762
static void pr_float(std::ostream &os, const float_display_format &fmt, T val)
Definition: pr-output.cc:1516
static T abs(T x)
Definition: pr-output.cc:1678
static int rat_string_len
Definition: pr-output.cc:90
static void pr_plus_format_matrix(std::ostream &os, const MT &m)
Definition: pr-output.cc:1840
static void pr_plus_format(std::ostream &os, const T &val)
Definition: pr-output.cc:1660
static float_display_format make_complex_format(int x_max, int x_min, int r_x, bool inf_or_nan, int int_only)
Definition: pr-output.cc:723
#define MAKE_INT_MATRIX_FORMAT(TYPE)
Definition: pr-output.cc:1704
static void pr_int(std::ostream &os, const T &d, int fw=0)
Definition: pr-output.cc:2810
float_display_format make_format(const double &d)
Definition: pr-output.cc:525
static void print_empty_nd_array(std::ostream &os, const dim_vector &dims, bool pr_as_read_syntax)
Definition: pr-output.cc:1592
static bool print_g
Definition: pr-output.cc:108
static void pr_any_float(std::ostream &os, const float_format &fmt, T val)
Definition: pr-output.cc:1394
static void octave_print_matrix_internal(std::ostream &os, const MT &m, bool pr_as_read_syntax, int extra_indent)
Definition: pr-output.cc:1875
static bool plus_format
Definition: pr-output.cc:81
static float_display_format make_matrix_format(const MT &m)
Definition: pr-output.cc:681
static bool free_format
Definition: pr-output.cc:78
static bool uppercase_format
Definition: pr-output.cc:111
static int get_column_width(const float_display_format &fmt)
Definition: pr-output.cc:1860
static float_display_format make_real_format(int digits, bool inf_or_nan, bool int_only)
Definition: pr-output.cc:394
static bool rat_format
Definition: pr-output.cc:87
bool Vprint_empty_dimensions
Definition: pr-output.cc:71
static float_display_format make_complex_matrix_format(int x_max, int x_min, int r_x_max, int r_x_min, bool inf_or_nan, int int_or_inf_or_nan)
Definition: pr-output.cc:945
static int hex_format
Definition: pr-output.cc:96
static void octave_print_free(std::ostream &os, const MT &m, bool pr_as_read_syntax)
Definition: pr-output.cc:1817
static void set_format_style(int argc, const string_vector &argv)
Definition: pr-output.cc:3577
static std::string format_string("short")
static int calc_scale_exp(const int &x)
Definition: pr-output.cc:117
#define INSTANTIATE_ABS(T)
Definition: pr-output.cc:1683
#define PRINT_INT_ARRAY_INTERNAL(TYPE)
Definition: pr-output.cc:3169
static bool Vfixed_point_format
Definition: pr-output.cc:67
static int num_digits(T x)
Definition: pr-output.cc:164
static T pr_min_internal(const MArray< T > &m)
Definition: pr-output.cc:333
float_display_format make_complex_scalar_format(const std::complex< T > &c)
Definition: pr-output.cc:891
static void pr_imag_float(std::ostream &os, const float_display_format &fmt, T val)
Definition: pr-output.cc:1528
static bool Vsplit_long_rows
Definition: pr-output.cc:75
static void octave_print_diag_matrix_internal(std::ostream &os, const DMT &m, bool pr_as_read_syntax, int extra_indent)
Definition: pr-output.cc:1990
std::ostream & operator<<(std::ostream &os, const pr_engineering_float< T > &pef)
Definition: pr-output.cc:187
#define PRINT_CHAR_BITS_SWAPPED(os, c)
Definition: pr-output.cc:1374
static int engineering_exponent(T x)
Definition: pr-output.cc:141
static bool print_eng
Definition: pr-output.cc:114
static float_display_format make_real_matrix_format(int x_max, int x_min, bool inf_or_nan, int int_or_inf_or_nan)
Definition: pr-output.cc:539
static bool bank_format
Definition: pr-output.cc:93
static std::string plus_format_chars
Definition: pr-output.cc:84
static int bit_format
Definition: pr-output.cc:99
static void pr_scale_header(std::ostream &os, double scale)
Definition: pr-output.cc:1609
void print_nd_array(std::ostream &os, const NDA_T &nda, bool pr_as_read_syntax)
Definition: pr-output.cc:2116
bool Vcompact_format
Definition: pr-output.cc:102
static int left
Definition: randmtzig.cc:194
static const int DIGITS10
Definition: pr-output.cc:371
static const int MAX_FIELD_WIDTH
Definition: pr-output.cc:372
static const int DIGITS10
Definition: pr-output.cc:381
static const int MAX_FIELD_WIDTH
Definition: pr-output.cc:382
static const int DIGITS10
Definition: pr-output.cc:364
static const int MAX_FIELD_WIDTH
Definition: pr-output.cc:365
unsigned char i[sizeof(T)]
Definition: pr-output.cc:1353
std::string undo_string_escapes(const std::string &s)
Definition: utils.cc:1033
bool valid_identifier(const char *s)
Definition: utils.cc:78
std::size_t format(std::ostream &os, const char *fmt,...)
Definition: utils.cc:1473
octave_value set_internal_variable(bool &var, const octave_value_list &args, int nargout, const char *nm)
Definition: variables.cc:584