GNU Octave  6.2.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
max.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // Copyright (C) 1996-2021 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 "lo-ieee.h"
33 #include "lo-mappers.h"
34 #include "dNDArray.h"
35 #include "CNDArray.h"
36 #include "quit.h"
37 
38 #include "defun.h"
39 #include "error.h"
40 #include "errwarn.h"
41 #include "ovl.h"
42 
43 #include "ov-cx-mat.h"
44 #include "ov-re-sparse.h"
45 #include "ov-cx-sparse.h"
46 
47 template <typename ArrayType>
48 static octave_value_list
50  int nargout, int dim, bool ismin)
51 {
52  octave_value_list retval (nargout > 1 ? 2 : 1);
53  ArrayType array = octave_value_extract<ArrayType> (arg);
54 
55  if (nargout <= 1)
56  {
57  if (ismin)
58  retval(0) = array.min (dim);
59  else
60  retval(0) = array.max (dim);
61  }
62  else
63  {
65  if (ismin)
66  retval(0) = array.min (idx, dim);
67  else
68  retval(0) = array.max (idx, dim);
69 
70  retval(1) = octave_value (idx, true, true);
71  }
72 
73  return retval;
74 }
75 
76 // Matlab returns double arrays for min/max operations on character
77 // arrays, so we specialize here to get that behavior. Other possible
78 // solutions are to convert the argument to double here and call the
79 // code for double, but that could waste memory, or to have the
80 // underlying charNDArray::min/max functions return NDArray instead of
81 // charNDArray, but that is inconsistent with the way other min/max
82 // functions work.
83 
84 template <>
87  int nargout, int dim, bool ismin)
88 {
89  octave_value_list retval (nargout > 1 ? 2 : 1);
91 
92  if (nargout <= 1)
93  {
94  if (ismin)
95  retval(0) = NDArray (array.min (dim));
96  else
97  retval(0) = NDArray (array.max (dim));
98  }
99  else
100  {
102  if (ismin)
103  retval(0) = NDArray (array.min (idx, dim));
104  else
105  retval(0) = NDArray (array.max (idx, dim));
106 
107  retval(1) = octave_value (idx, true, true);
108  }
109 
110  return retval;
111 }
112 
113 // Specialization for bool arrays (dense or sparse).
114 template <>
117  int nargout, int dim, bool ismin)
118 {
120 
121  if (! arg.issparse ())
122  {
123  if (nargout <= 1)
124  {
125  // This case can be handled using any/all.
126  boolNDArray array = arg.bool_array_value ();
127 
128  if (array.isempty ())
129  retval(0) = array;
130  else if (ismin)
131  retval(0) = array.all (dim);
132  else
133  retval(0) = array.any (dim);
134  }
135  else
136  {
137  // any/all don't have indexed versions, so do it via a conversion.
138  retval = do_minmax_red_op<int8NDArray> (arg, nargout, dim, ismin);
139 
140  retval(0) = retval(0).bool_array_value ();
141  }
142  }
143  else
144  {
145  // Sparse: Don't use any/all trick, as full matrix could exceed memory.
146  // Instead, convert to double.
147  retval = do_minmax_red_op<SparseMatrix> (arg, nargout, dim, ismin);
148 
149  retval(0) = retval(0).sparse_bool_matrix_value ();
150  }
151 
152  return retval;
153 }
154 
155 template <typename ArrayType>
156 static octave_value
157 do_minmax_bin_op (const octave_value& argx, const octave_value& argy,
158  bool ismin)
159 {
160  typedef typename ArrayType::element_type ScalarType;
161 
163 
164  if (argx.is_scalar_type ())
165  {
166  ScalarType x = octave_value_extract<ScalarType> (argx);
167  ArrayType y = octave_value_extract<ArrayType> (argy);
168 
169  if (ismin)
170  retval = min (x, y);
171  else
172  retval = max (x, y);
173  }
174  else if (argy.is_scalar_type ())
175  {
176  ArrayType x = octave_value_extract<ArrayType> (argx);
177  ScalarType y = octave_value_extract<ScalarType> (argy);
178 
179  if (ismin)
180  retval = min (x, y);
181  else
182  retval = max (x, y);
183  }
184  else
185  {
186  ArrayType x = octave_value_extract<ArrayType> (argx);
187  ArrayType y = octave_value_extract<ArrayType> (argy);
188 
189  if (ismin)
190  retval = min (x, y);
191  else
192  retval = max (x, y);
193  }
194 
195  return retval;
196 }
197 
198 // Matlab returns double arrays for min/max operations on character
199 // arrays, so we specialize here to get that behavior. Other possible
200 // solutions are to convert the arguments to double here and call the
201 // code for double, but that could waste a lot of memory, or to have the
202 // underlying charNDArray::min/max functions return NDArray instead of
203 // charNDArray, but that is inconsistent with the way other min/max
204 // functions work.
205 
206 template <>
209  const octave_value& argy, bool ismin)
210 {
212 
215 
216  if (ismin)
217  {
218  if (x.numel () == 1)
219  retval = NDArray (min (x(0), y));
220  else if (y.numel () == 1)
221  retval = NDArray (min (x, y(0)));
222  else
223  retval = NDArray (min (x, y));
224  }
225  else
226  {
227  if (x.numel () == 1)
228  retval = NDArray (max (x(0), y));
229  else if (y.numel () == 1)
230  retval = NDArray (max (x, y(0)));
231  else
232  retval = NDArray (max (x, y));
233  }
234 
235  return retval;
236 }
237 
238 static octave_value_list
240  int nargout, bool ismin)
241 {
242  int nargin = args.length ();
243 
244  if (nargin < 1 || nargin > 3)
245  print_usage ();
246 
247  octave_value_list retval (nargout > 1 ? 2 : 1);
248 
249  const char *func = (ismin ? "min" : "max");
250 
251  if (nargin == 3 || nargin == 1)
252  {
253  octave_value arg = args(0);
254  int dim = -1;
255  if (nargin == 3)
256  {
257  dim = args(2).int_value (true) - 1;
258 
259  if (dim < 0)
260  error ("%s: DIM must be a valid dimension", func);
261 
262  if (! args(1).isempty ())
263  warning ("%s: second argument is ignored", func);
264  }
265 
266  switch (arg.builtin_type ())
267  {
268  case btyp_double:
269  {
270  if (arg.is_range () && (dim == -1 || dim == 1))
271  {
272  Range range = arg.range_value ();
273  if (range.numel () < 1)
274  {
275  retval(0) = arg;
276  if (nargout > 1)
277  retval(1) = arg;
278  }
279  else if (ismin)
280  {
281  retval(0) = range.min ();
282  if (nargout > 1)
283  retval(1) = static_cast<double>
284  (range.inc () < 0 ? range.numel () : 1);
285  }
286  else
287  {
288  retval(0) = range.max ();
289  if (nargout > 1)
290  retval(1) = static_cast<double>
291  (range.inc () >= 0 ? range.numel () : 1);
292  }
293  }
294  else if (arg.issparse ())
295  retval = do_minmax_red_op<SparseMatrix> (arg, nargout, dim,
296  ismin);
297  else
298  retval = do_minmax_red_op<NDArray> (arg, nargout, dim, ismin);
299 
300  }
301  break;
302 
303  case btyp_complex:
304  {
305  if (arg.issparse ())
306  retval = do_minmax_red_op<SparseComplexMatrix> (arg, nargout, dim,
307  ismin);
308  else
309  retval = do_minmax_red_op<ComplexNDArray> (arg, nargout, dim,
310  ismin);
311  }
312  break;
313 
314  case btyp_float:
315  retval = do_minmax_red_op<FloatNDArray> (arg, nargout, dim, ismin);
316  break;
317 
318  case btyp_float_complex:
319  retval = do_minmax_red_op<FloatComplexNDArray> (arg, nargout, dim,
320  ismin);
321  break;
322 
323  case btyp_char:
324  retval = do_minmax_red_op<charNDArray> (arg, nargout, dim, ismin);
325  break;
326 
327 #define MAKE_INT_BRANCH(X) \
328  case btyp_ ## X: \
329  retval = do_minmax_red_op<X ## NDArray> (arg, nargout, dim, ismin); \
330  break;
331 
332  MAKE_INT_BRANCH (int8);
333  MAKE_INT_BRANCH (int16);
334  MAKE_INT_BRANCH (int32);
335  MAKE_INT_BRANCH (int64);
336  MAKE_INT_BRANCH (uint8);
337  MAKE_INT_BRANCH (uint16);
338  MAKE_INT_BRANCH (uint32);
339  MAKE_INT_BRANCH (uint64);
340 
341 #undef MAKE_INT_BRANCH
342 
343  case btyp_bool:
344  retval = do_minmax_red_op<boolNDArray> (arg, nargout, dim, ismin);
345  break;
346 
347  default:
348  err_wrong_type_arg (func, arg);
349  }
350  }
351  else
352  {
353  octave_value argx = args(0);
354  octave_value argy = args(1);
355  builtin_type_t xtyp = argx.builtin_type ();
356  builtin_type_t ytyp = argy.builtin_type ();
357  builtin_type_t rtyp;
358  if (xtyp == btyp_char && ytyp == btyp_char)
359  rtyp = btyp_char;
360  // FIXME: This is what should happen when boolNDArray has max()
361  // else if (xtyp == btyp_bool && ytyp == btyp_bool)
362  // rtyp = btyp_bool;
363  else
364  rtyp = btyp_mixed_numeric (xtyp, ytyp);
365 
366  switch (rtyp)
367  {
368  case btyp_double:
369  {
370  if ((argx.issparse ()
371  && (argy.issparse () || argy.is_scalar_type ()))
372  || (argy.issparse () && argx.is_scalar_type ()))
373  retval = do_minmax_bin_op<SparseMatrix> (argx, argy, ismin);
374  else
375  retval = do_minmax_bin_op<NDArray> (argx, argy, ismin);
376  }
377  break;
378 
379  case btyp_complex:
380  {
381  if ((argx.issparse ()
382  && (argy.issparse () || argy.is_scalar_type ()))
383  || (argy.issparse () && argx.is_scalar_type ()))
384  retval = do_minmax_bin_op<SparseComplexMatrix> (argx, argy,
385  ismin);
386  else
387  retval = do_minmax_bin_op<ComplexNDArray> (argx, argy, ismin);
388  }
389  break;
390 
391  case btyp_float:
392  retval = do_minmax_bin_op<FloatNDArray> (argx, argy, ismin);
393  break;
394 
395  case btyp_float_complex:
396  retval = do_minmax_bin_op<FloatComplexNDArray> (argx, argy, ismin);
397  break;
398 
399  case btyp_char:
400  retval = do_minmax_bin_op<charNDArray> (argx, argy, ismin);
401  break;
402 
403 #define MAKE_INT_BRANCH(X) \
404  case btyp_ ## X: \
405  retval = do_minmax_bin_op<X ## NDArray> (argx, argy, ismin); \
406  break;
407 
408  MAKE_INT_BRANCH (int8);
409  MAKE_INT_BRANCH (int16);
410  MAKE_INT_BRANCH (int32);
411  MAKE_INT_BRANCH (int64);
412  MAKE_INT_BRANCH (uint8);
413  MAKE_INT_BRANCH (uint16);
414  MAKE_INT_BRANCH (uint32);
415  MAKE_INT_BRANCH (uint64);
416 
417 #undef MAKE_INT_BRANCH
418 
419  // FIXME: This is what should happen when boolNDArray has max()
420  // case btyp_bool:
421  // retval = do_minmax_bin_op<boolNDArray> (argx, argy, ismin);
422  // break;
423 
424  default:
425  error ("%s: cannot compute %s (%s, %s)", func, func,
426  argx.type_name ().c_str (), argy.type_name ().c_str ());
427  }
428 
429  // FIXME: Delete when boolNDArray has max()
430  if (xtyp == btyp_bool && ytyp == btyp_bool)
431  retval(0) = retval(0).bool_array_value ();
432 
433  }
434 
435  return retval;
436 }
437 
438 DEFUN (min, args, nargout,
439  doc: /* -*- texinfo -*-
440 @deftypefn {} {} min (@var{x})
441 @deftypefnx {} {} min (@var{x}, [], @var{dim})
442 @deftypefnx {} {[@var{w}, @var{iw}] =} min (@var{x})
443 @deftypefnx {} {} min (@var{x}, @var{y})
444 Find minimum values in the array @var{x}.
445 
446 For a vector argument, return the minimum value. For a matrix argument,
447 return a row vector with the minimum value of each column. For a
448 multi-dimensional array, @code{min} operates along the first non-singleton
449 dimension.
450 
451 If the optional third argument @var{dim} is present then operate along
452 this dimension. In this case the second argument is ignored and should be
453 set to the empty matrix.
454 
455 For two inputs (@var{x} and @var{y}), return the pairwise minimum according to
456 the rules for @ref{Broadcasting}.
457 
458 Thus,
459 
460 @example
461 min (min (@var{x}))
462 @end example
463 
464 @noindent
465 returns the smallest element of the 2-D matrix @var{x}, and
466 
467 @example
468 @group
469 min (2:5, pi)
470  @result{} 2.0000 3.0000 3.1416 3.1416
471 @end group
472 @end example
473 
474 @noindent
475 compares each element of the range @code{2:5} with @code{pi}, and returns a
476 row vector of the minimum values.
477 
478 For complex arguments, the magnitude of the elements are used for
479 comparison. If the magnitudes are identical, then the results are ordered
480 by phase angle in the range (-pi, pi]. Hence,
481 
482 @example
483 @group
484 min ([-1 i 1 -i])
485  @result{} -i
486 @end group
487 @end example
488 
489 @noindent
490 because all entries have magnitude 1, but -i has the smallest phase angle
491 with value -pi/2.
492 
493 If called with one input and two output arguments, @code{min} also returns
494 the first index of the minimum value(s). Thus,
495 
496 @example
497 @group
498 [x, ix] = min ([1, 3, 0, 2, 0])
499  @result{} x = 0
500  ix = 3
501 @end group
502 @end example
503 @seealso{max, cummin, cummax}
504 @end deftypefn */)
505 {
506  return do_minmax_body (args, nargout, true);
507 }
508 
509 /*
510 ## Test generic double class
511 %!assert (min ([1, 4, 2, 3]), 1)
512 %!assert (min ([1; -10; 5; -2]), -10)
513 %!assert (min ([4, 2i 4.999; -2, 2, 3+4i]), [-2, 2, 4.999])
514 ## Special routines for char arrays
515 %!assert (min (["abc", "ABC"]), 65)
516 %!assert (min (["abc"; "CBA"]), [67 66 65])
517 ## Special routines for logical arrays
518 %!assert (min (logical ([])), logical ([]))
519 %!assert (min (logical ([0 0 1 0])), false)
520 %!assert (min (logical ([0 0 1 0; 0 1 1 0])), logical ([0 0 1 0]))
521 ## Single values
522 %!assert (min (single ([1, 4, 2, 3])), single (1))
523 %!assert (min (single ([1; -10; 5; -2])), single (-10))
524 %!assert (min (single ([4, 2i 4.999; -2, 2, 3+4i])), single ([-2, 2, 4.999]))
525 ## Integer values
526 %!assert (min (uint8 ([1, 4, 2, 3])), uint8 (1))
527 %!assert (min (uint8 ([1; -10; 5; -2])), uint8 (-10))
528 %!assert (min (int8 ([1, 4, 2, 3])), int8 (1))
529 %!assert (min (int8 ([1; -10; 5; -2])), int8 (-10))
530 %!assert (min (uint16 ([1, 4, 2, 3])), uint16 (1))
531 %!assert (min (uint16 ([1; -10; 5; -2])), uint16 (-10))
532 %!assert (min (int16 ([1, 4, 2, 3])), int16 (1))
533 %!assert (min (int16 ([1; -10; 5; -2])), int16 (-10))
534 %!assert (min (uint32 ([1, 4, 2, 3])), uint32 (1))
535 %!assert (min (uint32 ([1; -10; 5; -2])), uint32 (-10))
536 %!assert (min (int32 ([1, 4, 2, 3])), int32 (1))
537 %!assert (min (int32 ([1; -10; 5; -2])), int32 (-10))
538 %!assert (min (uint64 ([1, 4, 2, 3])), uint64 (1))
539 %!assert (min (uint64 ([1; -10; 5; -2])), uint64 (-10))
540 %!assert (min (int64 ([1, 4, 2, 3])), int64 (1))
541 %!assert (min (int64 ([1; -10; 5; -2])), int64 (-10))
542 ## Sparse double values
543 %!assert (min (sparse ([1, 4, 2, 3])), sparse (1))
544 %!assert (min (sparse ([1; -10; 5; -2])), sparse(-10))
545 ## FIXME: sparse doesn't order complex values by phase angle
546 %!test <51307>
547 %! assert (min (sparse ([4, 2i 4.999; -2, 2, 3+4i])), sparse ([-2, 2, 4.999]));
548 
549 ## Test dimension argument
550 %!test
551 %! x = reshape (1:8, [2,2,2]);
552 %! assert (min (x, [], 1), reshape ([1, 3, 5, 7], [1,2,2]));
553 %! assert (min (x, [], 2), reshape ([1, 2, 5, 6], [2,1,2]));
554 %! [y, i] = min (x, [], 3);
555 %! assert (ndims (y), 2);
556 %! assert (y, [1, 3; 2, 4]);
557 %! assert (ndims (i), 2);
558 %! assert (i, [1, 1; 1, 1]);
559 
560 ## Test 2-output forms for various arg types
561 ## Special routines for char arrays
562 %!test
563 %! [y, i] = min (["abc", "ABC"]);
564 %! assert (y, 65);
565 %! assert (i, 4);
566 ## Special routines for logical arrays
567 %!test
568 %! x = logical ([0 0 1 0]);
569 %! [y, i] = min (x);
570 %! assert (y, false);
571 %! assert (i, 1);
572 ## Special handling of ranges
573 %!test
574 %! rng = 1:2:10;
575 %! [y, i] = min (rng);
576 %! assert (y, 1);
577 %! assert (i, 1);
578 %! rng = 10:-2:1;
579 %! [y, i] = min (rng);
580 %! assert (y, 2);
581 %! assert (i, 5);
582 
583 ## Test 2-input calling form for various arg types
584 ## Test generic double class
585 %!test
586 %! x = [1, 2, 3, 4]; y = fliplr (x);
587 %! assert (min (x, y), [1 2 2 1]);
588 %! assert (min (x, 3), [1 2 3 3]);
589 %! assert (min (2, x), [1 2 2 2]);
590 %! assert (min (x, 2.1i), [1 2 2.1i 2.1i]);
591 ## FIXME: Ordering of complex results with equal magnitude is not by phase
592 ## angle in the 2-input form. Instead, it is in the order in which it
593 ## appears in the argument list.
594 %!test <51307>
595 %! x = [1, 2, 3, 4]; y = fliplr (x);
596 %! assert (min (x, 2i), [2i 2i 3 4]);
597 ## Special routines for char arrays
598 %!assert (min ("abc", "b"), [97 98 98])
599 %!assert (min ("b", "cba"), [98 98 97])
600 ## Special handling for logical arrays
601 %!assert (min ([true false], false), [false false])
602 %!assert (min (true, [true false]), [true false])
603 ## Single values
604 %!test
605 %! x = single ([1, 2, 3, 4]); y = fliplr (x);
606 %! assert (min (x, y), single ([1 2 2 1]));
607 %! assert (min (x, 3), single ([1 2 3 3]));
608 %! assert (min (2, x), single ([1 2 2 2]));
609 %! assert (min (x, 2.1i), single ([1 2 2.1i 2.1i]));
610 ## Integer values
611 %!test
612 %! x = uint8 ([1, 2, 3, 4]); y = fliplr (x);
613 %! assert (min (x, y), uint8 ([1 2 2 1]));
614 %! assert (min (x, 3), uint8 ([1 2 3 3]));
615 %! assert (min (2, x), uint8 ([1 2 2 2]));
616 %! x = int8 ([1, 2, 3, 4]); y = fliplr (x);
617 %! assert (min (x, y), int8 ([1 2 2 1]));
618 %! assert (min (x, 3), int8 ([1 2 3 3]));
619 %! assert (min (2, x), int8 ([1 2 2 2]));
620 %! x = uint16 ([1, 2, 3, 4]); y = fliplr (x);
621 %! assert (min (x, y), uint16 ([1 2 2 1]));
622 %! assert (min (x, 3), uint16 ([1 2 3 3]));
623 %! assert (min (2, x), uint16 ([1 2 2 2]));
624 %! x = int16 ([1, 2, 3, 4]); y = fliplr (x);
625 %! assert (min (x, y), int16 ([1 2 2 1]));
626 %! assert (min (x, 3), int16 ([1 2 3 3]));
627 %! assert (min (2, x), int16 ([1 2 2 2]));
628 %! x = uint32 ([1, 2, 3, 4]); y = fliplr (x);
629 %! assert (min (x, y), uint32 ([1 2 2 1]));
630 %! assert (min (x, 3), uint32 ([1 2 3 3]));
631 %! assert (min (2, x), uint32 ([1 2 2 2]));
632 %! x = int32 ([1, 2, 3, 4]); y = fliplr (x);
633 %! assert (min (x, y), int32 ([1 2 2 1]));
634 %! assert (min (x, 3), int32 ([1 2 3 3]));
635 %! assert (min (2, x), int32 ([1 2 2 2]));
636 %! x = uint64 ([1, 2, 3, 4]); y = fliplr (x);
637 %! assert (min (x, y), uint64 ([1 2 2 1]));
638 %! assert (min (x, 3), uint64 ([1 2 3 3]));
639 %! assert (min (2, x), uint64 ([1 2 2 2]));
640 %! x = int64 ([1, 2, 3, 4]); y = fliplr (x);
641 %! assert (min (x, y), int64 ([1 2 2 1]));
642 %! assert (min (x, 3), int64 ([1 2 3 3]));
643 %! assert (min (2, x), int64 ([1 2 2 2]));
644 ## Sparse double values
645 %!test
646 %! x = sparse ([1, 2, 3, 4]); y = fliplr (x);
647 %! assert (min (x, y), sparse ([1 2 2 1]));
648 %! assert (min (x, 3), sparse ([1 2 3 3]));
649 %! assert (min (2, x), sparse ([1 2 2 2]));
650 %! assert (min (x, 2.1i), sparse ([1 2 2.1i 2.1i]));
651 
652 %!error min ()
653 %!error min (1, 2, 3, 4)
654 %!error <DIM must be a valid dimension> min ([1 2; 3 4], [], -3)
655 %!warning <second argument is ignored> min ([1 2 3 4], 2, 2);
656 %!error <wrong type argument 'cell'> min ({1 2 3 4})
657 %!error <cannot compute min \‍(cell, scalar\‍)> min ({1, 2, 3}, 2)
658 */
659 
660 DEFUN (max, args, nargout,
661  doc: /* -*- texinfo -*-
662 @deftypefn {} {} max (@var{x})
663 @deftypefnx {} {} max (@var{x}, [], @var{dim})
664 @deftypefnx {} {[@var{w}, @var{iw}] =} max (@var{x})
665 @deftypefnx {} {} max (@var{x}, @var{y})
666 Find maximum values in the array @var{x}.
667 
668 For a vector argument, return the maximum value. For a matrix argument,
669 return a row vector with the maximum value of each column. For a
670 multi-dimensional array, @code{max} operates along the first non-singleton
671 dimension.
672 
673 If the optional third argument @var{dim} is present then operate along
674 this dimension. In this case the second argument is ignored and should be
675 set to the empty matrix.
676 
677 For two inputs (@var{x} and @var{y}), return the pairwise maximum according to
678 the rules for @ref{Broadcasting}.
679 
680 Thus,
681 
682 @example
683 max (max (@var{x}))
684 @end example
685 
686 @noindent
687 returns the largest element of the 2-D matrix @var{x}, and
688 
689 @example
690 @group
691 max (2:5, pi)
692  @result{} 3.1416 3.1416 4.0000 5.0000
693 @end group
694 @end example
695 
696 @noindent
697 compares each element of the range @code{2:5} with @code{pi}, and returns a
698 row vector of the maximum values.
699 
700 For complex arguments, the magnitude of the elements are used for
701 comparison. If the magnitudes are identical, then the results are ordered
702 by phase angle in the range (-pi, pi]. Hence,
703 
704 @example
705 @group
706 max ([-1 i 1 -i])
707  @result{} -1
708 @end group
709 @end example
710 
711 @noindent
712 because all entries have magnitude 1, but -1 has the largest phase angle
713 with value pi.
714 
715 If called with one input and two output arguments, @code{max} also returns
716 the first index of the maximum value(s). Thus,
717 
718 @example
719 @group
720 [x, ix] = max ([1, 3, 5, 2, 5])
721  @result{} x = 5
722  ix = 3
723 @end group
724 @end example
725 @seealso{min, cummax, cummin}
726 @end deftypefn */)
727 {
728  return do_minmax_body (args, nargout, false);
729 }
730 
731 /*
732 ## Test generic double class
733 %!assert (max ([1, 4, 2, 3]), 4)
734 %!assert (max ([1; -10; 5; -2]), 5)
735 %!assert (max ([4, 2i 4.999; -2, 2, 3+4i]), [4, 2i, 3+4i])
736 ## Special routines for char arrays
737 %!assert (max (["abc", "ABC"]), 99)
738 %!assert (max (["abc"; "CBA"]), [97 98 99])
739 ## Special routines for logical arrays
740 %!assert (max (logical ([])), logical ([]))
741 %!assert (max (logical ([0 0 1 0])), true)
742 %!assert (max (logical ([0 0 1 0; 0 1 0 0])), logical ([0 1 1 0]))
743 ## Single values
744 %!assert (max (single ([1, 4, 2, 3])), single (4))
745 %!assert (max (single ([1; -10; 5; -2])), single (5))
746 %!assert (max (single ([4, 2i 4.999; -2, 2, 3+4i])), single ([4, 2i, 3+4i]))
747 ## Integer values
748 %!assert (max (uint8 ([1, 4, 2, 3])), uint8 (4))
749 %!assert (max (uint8 ([1; -10; 5; -2])), uint8 (5))
750 %!assert (max (int8 ([1, 4, 2, 3])), int8 (4))
751 %!assert (max (int8 ([1; -10; 5; -2])), int8 (5))
752 %!assert (max (uint16 ([1, 4, 2, 3])), uint16 (4))
753 %!assert (max (uint16 ([1; -10; 5; -2])), uint16 (5))
754 %!assert (max (int16 ([1, 4, 2, 3])), int16 (4))
755 %!assert (max (int16 ([1; -10; 5; -2])), int16 (5))
756 %!assert (max (uint32 ([1, 4, 2, 3])), uint32 (4))
757 %!assert (max (uint32 ([1; -10; 5; -2])), uint32 (5))
758 %!assert (max (int32 ([1, 4, 2, 3])), int32 (4))
759 %!assert (max (int32 ([1; -10; 5; -2])), int32 (5))
760 %!assert (max (uint64 ([1, 4, 2, 3])), uint64 (4))
761 %!assert (max (uint64 ([1; -10; 5; -2])), uint64 (5))
762 %!assert (max (int64 ([1, 4, 2, 3])), int64 (4))
763 %!assert (max (int64 ([1; -10; 5; -2])), int64 (5))
764 ## Sparse double values
765 %!assert (max (sparse ([1, 4, 2, 3])), sparse (4))
766 %!assert (max (sparse ([1; -10; 5; -2])), sparse(5))
767 %!assert (max (sparse ([4, 2i 4.999; -2, 2, 3+4i])), sparse ([4, 2i, 3+4i]))
768 
769 ## Test dimension argument
770 %!test
771 %! x = reshape (1:8, [2,2,2]);
772 %! assert (min (x, [], 1), reshape ([1, 3, 5, 7], [1,2,2]));
773 %! assert (min (x, [], 2), reshape ([1, 2, 5, 6], [2,1,2]));
774 %! [y, i] = min (x, [], 3);
775 %! assert (ndims (y), 2);
776 %! assert (y, [1, 3; 2, 4]);
777 %! assert (ndims (i), 2);
778 %! assert (i, [1, 1; 1, 1]);
779 
780 ## Test 2-output forms for various arg types
781 ## Special routines for char arrays
782 %!test
783 %! [y, i] = max (["abc", "ABC"]);
784 %! assert (y, 99);
785 %! assert (i, 3);
786 ## Special routines for logical arrays
787 %!test
788 %! x = logical ([0 0 1 0]);
789 %! [y, i] = max (x);
790 %! assert (y, true);
791 %! assert (i, 3);
792 ## Special handling of ranges
793 %!test
794 %! rng = 1:2:10;
795 %! [y, i] = max (rng);
796 %! assert (y, 9);
797 %! assert (i, 5);
798 %! rng = 10:-2:1;
799 %! [y, i] = max (rng);
800 %! assert (y, 10);
801 %! assert (i, 1);
802 
803 ## Test 2-input calling form for various arg types
804 ## Test generic double class
805 %!test
806 %! x = [1, 2, 3, 4]; y = fliplr (x);
807 %! assert (max (x, y), [4 3 3 4]);
808 %! assert (max (x, 3), [3 3 3 4]);
809 %! assert (max (2, x), [2 2 3 4]);
810 %! assert (max (x, 2.1i), [2.1i 2.1i 3 4]);
811 ## FIXME: Ordering of complex results with equal magnitude is not by phase
812 ## angle in the 2-input form. Instead, it is in the order in which it
813 ## appears in the argument list.
814 %!test <51307>
815 %! x = [1, 2, 3, 4]; y = fliplr (x);
816 %! assert (max (x, 2i), [2i 2i 3 4]);
817 ## Special routines for char arrays
818 %!assert (max ("abc", "b"), [98 98 99])
819 %!assert (max ("b", "cba"), [99 98 98])
820 ## Special handling for logical arrays
821 %!assert (max ([true false], false), [true false])
822 %!assert (max (true, [false false]), [true true])
823 ## Single values
824 %!test
825 %! x = single ([1, 2, 3, 4]); y = fliplr (x);
826 %! assert (max (x, y), single ([4 3 3 4]));
827 %! assert (max (x, 3), single ([3 3 3 4]));
828 %! assert (max (2, x), single ([2 2 3 4]));
829 %! assert (max (x, 2.1i), single ([2.1i 2.1i 3 4]));
830 ## Integer values
831 %!test
832 %! x = uint8 ([1, 2, 3, 4]); y = fliplr (x);
833 %! assert (max (x, y), uint8 ([4 3 3 4]));
834 %! assert (max (x, 3), uint8 ([3 3 3 4]));
835 %! assert (max (2, x), uint8 ([2 2 3 4]));
836 %! x = int8 ([1, 2, 3, 4]); y = fliplr (x);
837 %! assert (max (x, y), int8 ([4 3 3 4]));
838 %! assert (max (x, 3), int8 ([3 3 3 4]));
839 %! assert (max (2, x), int8 ([2 2 3 4]));
840 %! x = uint16 ([1, 2, 3, 4]); y = fliplr (x);
841 %! assert (max (x, y), uint16 ([4 3 3 4]));
842 %! assert (max (x, 3), uint16 ([3 3 3 4]));
843 %! assert (max (2, x), uint16 ([2 2 3 4]));
844 %! x = int16 ([1, 2, 3, 4]); y = fliplr (x);
845 %! assert (max (x, y), int16 ([4 3 3 4]));
846 %! assert (max (x, 3), int16 ([3 3 3 4]));
847 %! assert (max (2, x), int16 ([2 2 3 4]));
848 %! x = uint32 ([1, 2, 3, 4]); y = fliplr (x);
849 %! assert (max (x, y), uint32 ([4 3 3 4]));
850 %! assert (max (x, 3), uint32 ([3 3 3 4]));
851 %! assert (max (2, x), uint32 ([2 2 3 4]));
852 %! x = int32 ([1, 2, 3, 4]); y = fliplr (x);
853 %! assert (max (x, y), int32 ([4 3 3 4]));
854 %! assert (max (x, 3), int32 ([3 3 3 4]));
855 %! assert (max (2, x), int32 ([2 2 3 4]));
856 %! x = uint64 ([1, 2, 3, 4]); y = fliplr (x);
857 %! assert (max (x, y), uint64 ([4 3 3 4]));
858 %! assert (max (x, 3), uint64 ([3 3 3 4]));
859 %! assert (max (2, x), uint64 ([2 2 3 4]));
860 %! x = int64 ([1, 2, 3, 4]); y = fliplr (x);
861 %! assert (max (x, y), int64 ([4 3 3 4]));
862 %! assert (max (x, 3), int64 ([3 3 3 4]));
863 %! assert (max (2, x), int64 ([2 2 3 4]));
864 ## Sparse double values
865 %!test
866 %! x = sparse ([1, 2, 3, 4]); y = fliplr (x);
867 %! assert (max (x, y), sparse ([4 3 3 4]));
868 %! assert (max (x, 3), sparse ([3 3 3 4]));
869 %! assert (max (2, x), sparse ([2 2 3 4]));
870 %! assert (max (x, 2.1i), sparse ([2.1i 2.1i 3 4]));
871 
872 ## Test for bug #40743
873 %!assert <*40743> (max (zeros (1,0), ones (1,1)), zeros (1,0))
874 %!assert <*40743> (max (sparse (zeros (1,0)), sparse (ones (1,1))),
875 %! sparse (zeros (1,0)))
876 
877 %!error max ()
878 %!error max (1, 2, 3, 4)
879 %!error <DIM must be a valid dimension> max ([1 2; 3 4], [], -3)
880 %!warning <second argument is ignored> max ([1 2 3 4], 2, 2);
881 %!error <wrong type argument 'cell'> max ({1 2 3 4})
882 %!error <cannot compute max \‍(cell, scalar\‍)> max ({1, 2, 3}, 2)
883 
884 */
885 
886 template <typename ArrayType>
887 static octave_value_list
889  int nargout, int dim, bool ismin)
890 {
891  octave_value_list retval (nargout > 1 ? 2 : 1);
892  ArrayType array = octave_value_extract<ArrayType> (arg);
893 
894  if (nargout <= 1)
895  {
896  if (ismin)
897  retval(0) = array.cummin (dim);
898  else
899  retval(0) = array.cummax (dim);
900  }
901  else
902  {
903  retval.resize (2);
905  if (ismin)
906  retval(0) = array.cummin (idx, dim);
907  else
908  retval(0) = array.cummax (idx, dim);
909 
910  retval(1) = octave_value (idx, true, true);
911  }
912 
913  return retval;
914 }
915 
916 static octave_value_list
918  int nargout, bool ismin)
919 {
920  int nargin = args.length ();
921 
922  if (nargin < 1 || nargin > 2)
923  print_usage ();
924 
925  const char *func = (ismin ? "cummin" : "cummax");
926 
927  octave_value arg = args(0);
928  int dim = -1;
929  if (nargin == 2)
930  {
931  dim = args(1).int_value (true) - 1;
932 
933  if (dim < 0)
934  error ("%s: DIM must be a valid dimension", func);
935  }
936 
938 
939  switch (arg.builtin_type ())
940  {
941  case btyp_double:
942  retval = do_cumminmax_red_op<NDArray> (arg, nargout, dim, ismin);
943  break;
944 
945  case btyp_complex:
946  retval = do_cumminmax_red_op<ComplexNDArray> (arg, nargout, dim,
947  ismin);
948  break;
949 
950  case btyp_float:
951  retval = do_cumminmax_red_op<FloatNDArray> (arg, nargout, dim, ismin);
952  break;
953 
954  case btyp_float_complex:
955  retval = do_cumminmax_red_op<FloatComplexNDArray> (arg, nargout, dim,
956  ismin);
957  break;
958 
959 #define MAKE_INT_BRANCH(X) \
960  case btyp_ ## X: \
961  retval = do_cumminmax_red_op<X ## NDArray> (arg, nargout, dim, ismin); \
962  break;
963 
964  MAKE_INT_BRANCH (int8);
965  MAKE_INT_BRANCH (int16);
966  MAKE_INT_BRANCH (int32);
967  MAKE_INT_BRANCH (int64);
968  MAKE_INT_BRANCH (uint8);
969  MAKE_INT_BRANCH (uint16);
970  MAKE_INT_BRANCH (uint32);
971  MAKE_INT_BRANCH (uint64);
972 
973 #undef MAKE_INT_BRANCH
974 
975  case btyp_bool:
976  {
977  retval = do_cumminmax_red_op<int8NDArray> (arg, nargout, dim,
978  ismin);
979  if (retval.length () > 0)
980  retval(0) = retval(0).bool_array_value ();
981  }
982  break;
983 
984  default:
985  err_wrong_type_arg (func, arg);
986  }
987 
988  return retval;
989 }
990 
991 DEFUN (cummin, args, nargout,
992  doc: /* -*- texinfo -*-
993 @deftypefn {} {} cummin (@var{x})
994 @deftypefnx {} {} cummin (@var{x}, @var{dim})
995 @deftypefnx {} {[@var{w}, @var{iw}] =} cummin (@var{x})
996 Return the cumulative minimum values along dimension @var{dim}.
997 
998 If @var{dim} is unspecified it defaults to column-wise operation. For
999 example:
1000 
1001 @example
1002 @group
1003 cummin ([5 4 6 2 3 1])
1004  @result{} 5 4 4 2 2 1
1005 @end group
1006 @end example
1007 
1008 If called with two output arguments the index of the minimum value is also
1009 returned.
1010 
1011 @example
1012 @group
1013 [w, iw] = cummin ([5 4 6 2 3 1])
1014 @result{}
1015 w = 5 4 4 2 2 1
1016 iw = 1 2 2 4 4 6
1017 @end group
1018 @end example
1019 
1020 @seealso{cummax, min, max}
1021 @end deftypefn */)
1022 {
1023  return do_cumminmax_body (args, nargout, true);
1024 }
1025 
1026 /*
1027 %!assert (cummin ([1, 4, 2, 3]), [1 1 1 1])
1028 %!assert (cummin ([1; -10; 5; -2]), [1; -10; -10; -10])
1029 %!assert (cummin ([4, i; -2, 2]), [4, i; -2, i])
1030 %!assert (cummin ([1 2; NaN 1], 2), [1 1; NaN 1])
1031 
1032 %!test
1033 %! x = reshape (1:8, [2,2,2]);
1034 %! assert (cummin (x, 1), reshape ([1 1 3 3 5 5 7 7], [2,2,2]));
1035 %! assert (cummin (x, 2), reshape ([1 2 1 2 5 6 5 6], [2,2,2]));
1036 %! [w, iw] = cummin (x, 3);
1037 %! assert (ndims (w), 3);
1038 %! assert (w, repmat ([1 3; 2 4], [1 1 2]));
1039 %! assert (ndims (iw), 3);
1040 %! assert (iw, ones (2,2,2));
1041 
1042 %!error cummin ()
1043 %!error cummin (1, 2, 3)
1044 */
1045 
1046 DEFUN (cummax, args, nargout,
1047  doc: /* -*- texinfo -*-
1048 @deftypefn {} {} cummax (@var{x})
1049 @deftypefnx {} {} cummax (@var{x}, @var{dim})
1050 @deftypefnx {} {[@var{w}, @var{iw}] =} cummax (@dots{})
1051 Return the cumulative maximum values along dimension @var{dim}.
1052 
1053 If @var{dim} is unspecified it defaults to column-wise operation. For
1054 example:
1055 
1056 @example
1057 @group
1058 cummax ([1 3 2 6 4 5])
1059  @result{} 1 3 3 6 6 6
1060 @end group
1061 @end example
1062 
1063 If called with two output arguments the index of the maximum value is also
1064 returned.
1065 
1066 @example
1067 @group
1068 [w, iw] = cummax ([1 3 2 6 4 5])
1069 @result{}
1070 w = 1 3 3 6 6 6
1071 iw = 1 2 2 4 4 4
1072 @end group
1073 @end example
1074 
1075 @seealso{cummin, max, min}
1076 @end deftypefn */)
1077 {
1078  return do_cumminmax_body (args, nargout, false);
1079 }
1080 
1081 /*
1082 %!assert (cummax ([1, 4, 2, 3]), [1 4 4 4])
1083 %!assert (cummax ([1; -10; 5; -2]), [1; 1; 5; 5])
1084 %!assert (cummax ([4, i 4.9, -2, 2, 3+4i]), [4, 4, 4.9, 4.9, 4.9, 3+4i])
1085 %!assert (cummax ([1 NaN 0; NaN NaN 1], 2), [1 1 1; NaN NaN 1])
1086 
1087 %!test
1088 %! x = reshape (8:-1:1, [2,2,2]);
1089 %! assert (cummax (x, 1), reshape ([8 8 6 6 4 4 2 2], [2,2,2]));
1090 %! assert (cummax (x, 2), reshape ([8 7 8 7 4 3 4 3], [2,2,2]));
1091 %! [w, iw] = cummax (x, 3);
1092 %! assert (ndims (w), 3);
1093 %! assert (w, repmat ([8 6; 7 5], [1 1 2]));
1094 %! assert (ndims (iw), 3);
1095 %! assert (iw, ones (2,2,2));
1096 
1097 %!error cummax ()
1098 %!error cummax (1, 2, 3)
1099 */
charNDArray max(char d, const charNDArray &m)
Definition: chNDArray.cc:230
charNDArray min(char d, const charNDArray &m)
Definition: chNDArray.cc:207
void resize(const dim_vector &dv, const T &rfv)
Size of the specified dimension.
Definition: Array.cc:1011
octave_idx_type numel(void) const
Number of elements in the array.
Definition: Array.h:377
bool isempty(void) const
Size of the specified dimension.
Definition: Array.h:572
Definition: Range.h:40
double inc(void) const
Definition: Range.h:85
octave_idx_type numel(void) const
Definition: Range.h:87
double max(void) const
Definition: Range.cc:230
double min(void) const
Definition: Range.cc:209
boolNDArray any(int dim=-1) const
Definition: boolNDArray.cc:67
boolNDArray all(int dim=-1) const
Definition: boolNDArray.cc:61
charNDArray max(int dim=-1) const
Definition: chNDArray.cc:142
charNDArray min(int dim=-1) const
Definition: chNDArray.cc:154
octave_idx_type length(void) const
Definition: ovl.h:113
boolNDArray bool_array_value(bool warn=false) const
Definition: ov.h:844
bool issparse(void) const
Definition: ov.h:706
builtin_type_t builtin_type(void) const
Definition: ov.h:643
bool is_scalar_type(void) const
Definition: ov.h:697
Range range_value(void) const
Definition: ov.h:938
bool is_range(void) const
Definition: ov.h:602
std::string type_name(void) const
Definition: ov.h:1254
OCTINTERP_API void print_usage(void)
Definition: defun.cc:53
#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:1050
void error(const char *fmt,...)
Definition: error.cc:968
void err_wrong_type_arg(const char *name, const char *s)
Definition: errwarn.cc:166
F77_RET_T const F77_DBLE * x
static octave_value_list do_minmax_body(const octave_value_list &args, int nargout, bool ismin)
Definition: max.cc:239
octave_value do_minmax_bin_op< charNDArray >(const octave_value &argx, const octave_value &argy, bool ismin)
Definition: max.cc:208
#define MAKE_INT_BRANCH(X)
static octave_value_list do_cumminmax_red_op(const octave_value &arg, int nargout, int dim, bool ismin)
Definition: max.cc:888
octave_value_list do_minmax_red_op< boolNDArray >(const octave_value &arg, int nargout, int dim, bool ismin)
Definition: max.cc:116
octave_value_list do_minmax_red_op< charNDArray >(const octave_value &arg, int nargout, int dim, bool ismin)
Definition: max.cc:86
static octave_value_list do_cumminmax_body(const octave_value_list &args, int nargout, bool ismin)
Definition: max.cc:917
static octave_value_list do_minmax_red_op(const octave_value &arg, int nargout, int dim, bool ismin)
Definition: max.cc:49
static octave_value do_minmax_bin_op(const octave_value &argx, const octave_value &argy, bool ismin)
Definition: max.cc:157
return octave_value(v1.char_array_value() . concat(v2.char_array_value(), ra_idx),((a1.is_sq_string()||a2.is_sq_string()) ? '\'' :'"'))
builtin_type_t btyp_mixed_numeric(builtin_type_t x, builtin_type_t y)
Determine the resulting type for a possible mixed-type operation.
Definition: ov-base.cc:65
builtin_type_t
Definition: ov-base.h:72
@ btyp_float_complex
Definition: ov-base.h:76
@ btyp_double
Definition: ov-base.h:73
@ btyp_float
Definition: ov-base.h:74
@ btyp_bool
Definition: ov-base.h:85
@ btyp_char
Definition: ov-base.h:86
@ btyp_complex
Definition: ov-base.h:75
octave_value::octave_value(const Array< char > &chm, char type) return retval
Definition: ov.cc:811
charNDArray octave_value_extract< charNDArray >(const octave_value &v)
Definition: ov.h:1637