GNU Octave 11.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
max.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 1996-2026 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 "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
48
49static void
50get_dim_vecdim_all (const octave_value& dimarg, octave_value& arg,
51 int& dim, Array<int>& perm_vec, bool& do_perm,
52 bool& allflag, const char *fcn)
53{
54 if (dimarg.is_scalar_type ())
55 {
56 dim = dimarg.int_value () - 1;
57 if (dim < 0)
58 error ("%s: invalid dimension DIM = %d", fcn, dim + 1);
59 }
60 else
61 {
62 Array<int> vec = dimarg.int_vector_value ();
63 std::vector<int> vecdim;
64 dim_vector sz = arg.dims ();
65 int ndims = arg.ndims ();
66 // Check for invalid dims and ignore any dims larger than actual ndims
67 for (int i = 0; i < vec.numel (); i++)
68 {
69 vec(i)--;
70 if (vec(i) < 0)
71 error ("%s: invalid dimension in VECDIM = %d", fcn, vec(i)+1);
72 if (vec(i) < ndims)
73 vecdim.push_back (vec(i));
74 }
75 int n = vecdim.size ();
76 // If no dimensions left, set DIM = ndims to return input as is
77 if (n == 0)
78 {
79 dim = ndims;
80 return;
81 }
82 octave_idx_type szvecdim = 1;
83 // Check for duplicate dims and add VECDIM to permutation vector
84 perm_vec.resize (dim_vector (1, ndims));
85 std::sort (vecdim.begin (), vecdim.end ());
86 // Check for duplicates FIRST before any array writes
87 auto dup = std::adjacent_find (vecdim.begin (), vecdim.end ());
88 if (dup != vecdim.end ())
89 error ("%s: duplicate dimension in VECDIM = %d", fcn, *dup + 1);
90
91 // Verified vecdim has unique entries in [0, ndims-1], hence n <= ndims
92 int out_pos = ndims - n;
93 for (int d : vecdim)
94 {
95 szvecdim *= sz (d);
96 perm_vec(out_pos++) = d;
97 }
98
99 // Parse vecdim
100 if (n == 1)
101 // Only one dimension given, treat as if dim were specified instead
102 dim = vecdim[0];
103 else if (ndims == n)
104 // vecdim contains all dimensions, treat as if "all" flag given.
105 allflag = true;
106 else
107 {
108 dim_vector new_sz;
109 new_sz.resize (ndims - n + 1);
110 int idx = 0;
111 // Add remaining dims to permutation vector
112 for (int i = 0; i < ndims; i++)
113 {
114 if (std::find (vecdim.begin (), vecdim.end (), i)
115 == vecdim.end ())
116 {
117 perm_vec(idx) = i;
118 new_sz(idx) = sz(i);
119 idx++;
120 }
121 }
122 new_sz(idx) = szvecdim;
123 arg = arg.permute (perm_vec, false);
124 arg = arg.reshape (new_sz);
125 do_perm = true;
126 dim = idx;
127 }
128 }
129}
130
131static dim_vector
132custom_ind2sub (const dim_vector szx, octave_idx_type idx, int ndims)
133{
135 szc.resize (dim_vector (1, ndims));
136 for (int i = 0; i < ndims; i++)
137 szc(i) = szx(i);
138 if (ndims > 1)
139 for (int i = 1; i < ndims; i++)
140 szc(i) *= szc(i-1);
141
142 dim_vector sub;
143 sub.resize (ndims);
144 for (int i = 0; i < ndims; i++)
145 sub(i) = 1;
146
147 for (octave_idx_type i = ndims - 1; i >= 0; i--)
148 {
149 int dsz = szc(i);
150 if (idx > dsz)
151 {
152 int tmp = idx / dsz; // dzc is never 0 at this point
153 if (idx % dsz > 0)
154 {
155 sub(i+1) = tmp + 1;
156 idx -= tmp * dsz;
157 }
158 else
159 {
160 sub(i+1) = tmp;
161 idx -= (tmp - 1) * dsz;
162 }
163 }
164 if (i == 0)
165 sub(i) = idx;
166 }
167 return sub;
168}
169
170static octave_idx_type
171custom_sub2ind (const dim_vector& szx, const dim_vector& dims, int ndims)
172{
173 octave_idx_type ind = 1;
175 szc.resize (dim_vector (1, ndims));
176 for (int i = 0; i < ndims; i++)
177 szc(i) = szx(i);
178 for (octave_idx_type i = 1; i < ndims; i++)
179 szc(i) *= szc(i-1);
180
182 szce.resize (dim_vector (1, ndims + 1));
183 szce(0) = 1;
184 for (octave_idx_type i = 0; i < ndims; i++)
185 szce(i+1) = szc(i);
186
187 for (octave_idx_type i = 0; i < ndims; i++)
188 if (szx(i) > 1)
189 ind += szce(i) * (dims(i) - 1);
190
191 return ind;
192}
193
194static octave_value
195get_linear_index (const octave_value& index, const octave_value& arg,
196 const Array<int>& vecdim, bool cumop)
197{
198 // Get index values and original index size
200 dim_vector sz_idx = index.dims ();
201
202 // Get size and number of dimensions of input argument
203 dim_vector sz_arg = arg.dims ();
204 int nd_arg = arg.ndims ();
205
206 // Invalid input in VECDIM has already been ruled out
207 // Ignore any exceeding dimensions and get usable dimensions in VECDIM
208 std::vector<int> dim_page;
209 for (octave_idx_type i = 0; i < vecdim.numel (); i++)
210 if (vecdim(i) <= nd_arg)
211 dim_page.push_back (vecdim(i));
212 int nd_page = dim_page.size ();
213
214 // If no dimensions left, return linear index as [1:numel(x)]
215 if (nd_page == 0)
216 for (octave_idx_type i = 0; i < idx.numel (); i++)
217 idx(i) += i;
218 else
219 {
220 // Sort VECDIM
221 std::sort (dim_page.begin (), dim_page.end ());
222
223 // Calculate page size
224 dim_vector sz_page;
225 sz_page.resize (nd_page);
226 for (int i = 0; i < nd_page; i++)
227 sz_page(i) = sz_arg(dim_page[i]-1);
228
229 // Calculate operating index vector and reduced index dimensions
230 std::vector<int> dim_index;
231 dim_vector dim_vec;
232 dim_vec.resize (nd_arg);
233 for (int i = 0; i < nd_arg; i++)
234 {
235 dim_vec(i) = 1;
236 bool addindex = true;
237 for (int j = 0; j < nd_page; j++)
238 if (i+1 == dim_page[j])
239 addindex = false;
240 // Ignore singleton indexing dimensions
241 if (addindex && sz_arg(i) > 1)
242 dim_index.push_back (i + 1);
243 }
244 // Calculate reduced indexing size
245 int nd_index = dim_index.size ();
246 dim_vector sz_index;
247 sz_index.resize (nd_index);
248 for (int i = 0; i < nd_index; i++)
249 sz_index(i) = sz_arg(dim_index[i]-1);
250
251 // In cumulative operations, the first dimension of the operating page is
252 // expanded to 'prod (vecdim)'. If the first dimension of the operating
253 // page is less than the first dimension in dim_index, every prod_vecdim
254 // number of index elements reuse the corresponding index of the first
255 // element. Otherwise, they get reused after indexing the elements whose
256 // dim_index is less than first_vecdim.
257 octave_idx_type prod_vecdim = 1;
258 octave_idx_type group_index = 1;
259 if (cumop)
260 {
261 for (octave_idx_type i = 0; i < nd_page; i++)
262 prod_vecdim *= sz_page(i);
263 if (dim_page[0] > dim_index[0])
264 for (int i = 0; i < nd_index; i++)
265 if (dim_page[0] > dim_index[i])
266 group_index *= sz_index(i);
267 }
268
269 // Process each index value
270 dim_vector tmp;
272 for (octave_idx_type i = 0; i < idx.numel (); i++)
273 {
274 if (cumop)
275 {
276 if (group_index == 1)
277 ii = i / prod_vecdim;
278 else
279 {
280 octave_idx_type igm = i % group_index;
281 octave_idx_type igf = i / group_index;
282 ii = igm + ((igf / prod_vecdim) * group_index);
283 }
284 }
285 else
286 ii = i;
287 // Get index position and add it to corresponding dimensions
288 tmp = custom_ind2sub (sz_index, ii+1, nd_index);
289 for (octave_idx_type j = 0; j < nd_index; j++)
290 dim_vec(dim_index[j]-1) = tmp(j);
291
292 // Get page position and add it to corresponding dimensions
293 tmp = custom_ind2sub (sz_page, idx(i), nd_page);
294 for (octave_idx_type j = 0; j < nd_page; j++)
295 dim_vec(dim_page[j]-1) = tmp(j);
296
297 // Calculate linear index
298 idx(i) = custom_sub2ind (sz_arg, dim_vec, nd_arg);
299 }
300 }
301 return octave_value (idx).reshape (sz_idx);
302}
303
304template <typename ArrayType>
306do_minmax_red_op (const octave_value& arg,
307 int nargout, int dim, bool ismin)
308{
309 octave_value_list retval (nargout > 1 ? 2 : 1);
310 ArrayType array = octave_value_extract<ArrayType> (arg);
311
312 if (nargout <= 1)
313 {
314 if (ismin)
315 retval(0) = array.min (dim);
316 else
317 retval(0) = array.max (dim);
318 }
319 else
320 {
322 if (ismin)
323 retval(0) = array.min (idx, dim);
324 else
325 retval(0) = array.max (idx, dim);
326
327 retval(1) = octave_value (idx, true, true);
328 }
329
330 return retval;
331}
332
333template <typename ArrayType>
335do_minmax_red_op (const octave_value& arg, int nargout,
336 int dim, bool nanflag, bool ismin)
337{
338 octave_value_list retval (nargout > 1 ? 2 : 1);
339 ArrayType array = octave_value_extract<ArrayType> (arg);
340
341 if (nargout <= 1)
342 {
343 if (ismin)
344 retval(0) = array.min (dim, nanflag);
345 else
346 retval(0) = array.max (dim, nanflag);
347 }
348 else
349 {
351 if (ismin)
352 retval(0) = array.min (idx, dim, nanflag);
353 else
354 retval(0) = array.max (idx, dim, nanflag);
355
356 retval(1) = octave_value (idx, true, true);
357 }
358
359 return retval;
360}
361
362template <typename ArrayType>
364do_minmax_red_op (const octave_value& arg, int nargout,
365 int dim, bool nanflag, bool realabs, bool ismin)
366{
367 octave_value_list retval (nargout > 1 ? 2 : 1);
368 ArrayType array = octave_value_extract<ArrayType> (arg);
369
370 if (nargout <= 1)
371 {
372 if (ismin)
373 retval(0) = array.min (dim, nanflag, realabs);
374 else
375 retval(0) = array.max (dim, nanflag, realabs);
376 }
377 else
378 {
380 if (ismin)
381 retval(0) = array.min (idx, dim, nanflag, realabs);
382 else
383 retval(0) = array.max (idx, dim, nanflag, realabs);
384
385 retval(1) = octave_value (idx, true, true);
386 }
387
388 return retval;
389}
390
391// Matlab returns double arrays for min/max operations on character
392// arrays, so we specialize here to get that behavior. Other possible
393// solutions are to convert the argument to double here and call the
394// code for double, but that could waste memory, or to have the
395// underlying charNDArray::min/max functions return NDArray instead of
396// charNDArray, but that is inconsistent with the way other min/max
397// functions work.
398
399template <>
402 int nargout, int dim, bool ismin)
403{
404 octave_value_list retval (nargout > 1 ? 2 : 1);
406
407 if (nargout <= 1)
408 {
409 if (ismin)
410 retval(0) = NDArray (array.min (dim));
411 else
412 retval(0) = NDArray (array.max (dim));
413 }
414 else
415 {
417 if (ismin)
418 retval(0) = NDArray (array.min (idx, dim));
419 else
420 retval(0) = NDArray (array.max (idx, dim));
421
422 retval(1) = octave_value (idx, true, true);
423 }
424
425 return retval;
426}
427
428// Specialization for bool arrays (dense or sparse).
429template <>
432 int nargout, int dim, bool ismin)
433{
434 octave_value_list retval;
435
436 if (! arg.issparse ())
437 {
438 if (nargout <= 1)
439 {
440 // This case can be handled using any/all.
441 boolNDArray array = arg.bool_array_value ();
442
443 if (array.isempty ())
444 retval(0) = array;
445 else if (ismin)
446 retval(0) = array.all (dim);
447 else
448 retval(0) = array.any (dim);
449 }
450 else
451 {
452 // any/all don't have indexed versions, so do it via a conversion.
453 retval = do_minmax_red_op<int8NDArray> (arg, nargout, dim, ismin);
454
455 retval(0) = retval(0).bool_array_value ();
456 }
457 }
458 else
459 {
460 // Sparse: Don't use any/all trick, as full matrix could exceed memory.
461 // Instead, convert to double.
462 retval = do_minmax_red_op<SparseMatrix> (arg, nargout, dim, ismin);
463
464 retval(0) = retval(0).sparse_bool_matrix_value ();
465 }
466
467 return retval;
468}
469
470template <typename ArrayType>
471static octave_value
472do_minmax_bin_op (const octave_value& argx, const octave_value& argy,
473 bool ismin)
474{
475 typedef typename ArrayType::element_type ScalarType;
476
477 octave_value retval;
478
479 if (argx.is_scalar_type ())
480 {
481 ScalarType x = octave_value_extract<ScalarType> (argx);
482 ArrayType y = octave_value_extract<ArrayType> (argy);
483
484 if (ismin)
485 retval = min (x, y);
486 else
487 retval = max (x, y);
488 }
489 else if (argy.is_scalar_type ())
490 {
491 ArrayType x = octave_value_extract<ArrayType> (argx);
492 ScalarType y = octave_value_extract<ScalarType> (argy);
493
494 if (ismin)
495 retval = min (x, y);
496 else
497 retval = max (x, y);
498 }
499 else
500 {
501 ArrayType x = octave_value_extract<ArrayType> (argx);
502 ArrayType y = octave_value_extract<ArrayType> (argy);
503
504 if (ismin)
505 retval = min (x, y);
506 else
507 retval = max (x, y);
508 }
509
510 return retval;
511}
512
513template <typename ArrayType>
514static octave_value
515do_minmax_bin_op (const octave_value& argx, const octave_value& argy,
516 bool nanflag, bool ismin)
517{
518 typedef typename ArrayType::element_type ScalarType;
519
520 octave_value retval;
521
522 if (argx.is_scalar_type ())
523 {
524 ScalarType x = octave_value_extract<ScalarType> (argx);
525 ArrayType y = octave_value_extract<ArrayType> (argy);
526
527 if (ismin)
528 retval = min (x, y, nanflag);
529 else
530 retval = max (x, y, nanflag);
531 }
532 else if (argy.is_scalar_type ())
533 {
534 ArrayType x = octave_value_extract<ArrayType> (argx);
535 ScalarType y = octave_value_extract<ScalarType> (argy);
536
537 if (ismin)
538 retval = min (x, y, nanflag);
539 else
540 retval = max (x, y, nanflag);
541 }
542 else
543 {
544 ArrayType x = octave_value_extract<ArrayType> (argx);
545 ArrayType y = octave_value_extract<ArrayType> (argy);
546
547 if (ismin)
548 retval = min (x, y, nanflag);
549 else
550 retval = max (x, y, nanflag);
551 }
552
553 return retval;
554}
555
556template <typename ArrayType>
557static octave_value
558do_minmax_bin_op (const octave_value& argx, const octave_value& argy,
559 bool nanflag, bool realabs, bool ismin)
560{
561 typedef typename ArrayType::element_type ScalarType;
562
563 octave_value retval;
564
565 if (argx.is_scalar_type ())
566 {
567 ScalarType x = octave_value_extract<ScalarType> (argx);
568 ArrayType y = octave_value_extract<ArrayType> (argy);
569
570 if (ismin)
571 retval = min (x, y, nanflag, realabs);
572 else
573 retval = max (x, y, nanflag, realabs);
574 }
575 else if (argy.is_scalar_type ())
576 {
577 ArrayType x = octave_value_extract<ArrayType> (argx);
578 ScalarType y = octave_value_extract<ScalarType> (argy);
579
580 if (ismin)
581 retval = min (x, y, nanflag, realabs);
582 else
583 retval = max (x, y, nanflag, realabs);
584 }
585 else
586 {
587 ArrayType x = octave_value_extract<ArrayType> (argx);
588 ArrayType y = octave_value_extract<ArrayType> (argy);
589
590 if (ismin)
591 retval = min (x, y, nanflag, realabs);
592 else
593 retval = max (x, y, nanflag, realabs);
594 }
595
596 return retval;
597}
598
599// Matlab returns double arrays for min/max operations on character
600// arrays, so we specialize here to get that behavior. Other possible
601// solutions are to convert the arguments to double here and call the
602// code for double, but that could waste a lot of memory, or to have the
603// underlying charNDArray::min/max functions return NDArray instead of
604// charNDArray, but that is inconsistent with the way other min/max
605// functions work.
606
607template <>
610 const octave_value& argy, bool ismin)
611{
612 octave_value retval;
613
616
617 if (ismin)
618 {
619 if (x.numel () == 1)
620 retval = NDArray (min (x(0), y));
621 else if (y.numel () == 1)
622 retval = NDArray (min (x, y(0)));
623 else
624 retval = NDArray (min (x, y));
625 }
626 else
627 {
628 if (x.numel () == 1)
629 retval = NDArray (max (x(0), y));
630 else if (y.numel () == 1)
631 retval = NDArray (max (x, y(0)));
632 else
633 retval = NDArray (max (x, y));
634 }
635
636 return retval;
637}
638
640do_minmax_body (const octave_value_list& args, int nargout, bool ismin)
641{
642 int nargin = args.length ();
643 int orig_nargin = nargin;
644
645 const char *fcn = (ismin ? "min" : "max");
646
647 bool do_perm = false;
648 bool allflag = false;
649 bool nanflag = true; // NaNs are ignored by default
650 bool cmethod = true; // resolves to "auto"
651 bool realabs = true;
652 bool red_op = false;
653 bool linear = false;
654
655 // Get "ComparisonMethod" optional argument first
656 if (nargin > 2
657 && args(nargin - 1).is_string () && args(nargin - 2).is_string ())
658 {
659 std::string name = args(nargin - 2).string_value ();
660 std::string sval = args(nargin - 1).string_value ();
661 if (name == "ComparisonMethod" || name == "comparisonmethod")
662 {
663 if (sval == "auto")
664 cmethod = true;
665 else if (sval == "real")
666 cmethod = false;
667 else if (sval == "abs")
668 {
669 cmethod = false;
670 realabs = false;
671 }
672 else
673 error ("%s: invalid comparison method '%s'", fcn, sval.c_str ());
674
675 nargin -= 2;
676 }
677 }
678
679 while (nargin > 2 && args(nargin - 1).is_string ())
680 {
681 std::string str = args(nargin - 1).string_value ();
682
683 if (str == "all")
684 {
685 allflag = true;
686 red_op = true;
687 }
688 else if (str == "omitnan" || str == "omitmissing")
689 {
690 if (args(0).is_double_type () || args(0).is_single_type ())
691 nanflag = true;
692 }
693 else if (str == "includenan" || str == "includemissing")
694 nanflag = false;
695 else if (str == "linear")
696 linear = true;
697 else
698 error ("%s: unrecognized optional argument '%s'", fcn, str.c_str ());
699
700 nargin--;
701 }
702
703 if (nargin < 1 || nargin > 3)
704 print_usage ();
705 if (allflag && nargin > 2)
706 error ("%s: cannot set DIM or VECDIM with 'all' flag", fcn);
707 if (nargin > 2)
708 red_op = true;
709 if (nargin > 1)
710 if (orig_nargin > 2 && args(1).rows () == 0 && args(1). columns () == 0)
711 red_op = true;
712
713 octave_value_list retval (nargout > 1 ? 2 : 1);
714 Array<int> vecdim;
715
716 if (nargin == 1 || red_op)
717 {
718 octave_value arg = args(0);
719 // Handle DIM, VECDIM
720 int dim = -1;
721 Array<int> perm_vec;
722 if (nargin == 3)
723 {
724 octave_value dimarg = args(2);
725 get_dim_vecdim_all (dimarg, arg, dim, perm_vec, do_perm, allflag, fcn);
726 vecdim = dimarg.int_vector_value ();
727
728 if (! args(1).isempty ())
729 warning ("%s: second argument is ignored", fcn);
730 }
731 else
732 {
733 dim_vector dims = arg.dims ();
734 vecdim.resize (dim_vector (1, 1));
735 vecdim(0) = dims.first_non_singleton () + 1;
736 }
737 // Handle allflag
738 if (allflag)
739 arg = arg.reshape (dim_vector (arg.numel (), 1));
740
741 switch (arg.builtin_type ())
742 {
743 case btyp_double:
744 {
745 if (arg.is_range () && (dim == -1 || dim == 1))
746 {
748 if (range.numel () < 1)
749 {
750 retval(0) = arg;
751 if (nargout > 1)
752 retval(1) = arg;
753 }
754 else if (ismin)
755 {
756 retval(0) = range.min ();
757 if (nargout > 1)
758 retval(1) = static_cast<double>
759 (range.increment () < 0 ? range.numel () : 1);
760 }
761 else
762 {
763 retval(0) = range.max ();
764 if (nargout > 1)
765 retval(1) = static_cast<double>
766 (range.increment () >= 0 ? range.numel ()
767 : 1);
768 }
769 }
770 else if (arg.issparse ())
771 {
772 if (cmethod)
773 retval = do_minmax_red_op<SparseMatrix>
774 (arg, nargout, dim, nanflag, ismin);
775 else
776 retval = do_minmax_red_op<SparseMatrix>
777 (arg, nargout, dim, nanflag, realabs, ismin);
778 }
779 else
780 {
781 if (cmethod)
782 retval = do_minmax_red_op<NDArray>
783 (arg, nargout, dim, nanflag, ismin);
784 else
785 retval = do_minmax_red_op<NDArray>
786 (arg, nargout, dim, nanflag, realabs, ismin);
787 }
788 }
789 break;
790
791 case btyp_complex:
792 {
793 if (arg.issparse ())
794 {
795 if (cmethod)
796 retval = do_minmax_red_op<SparseComplexMatrix>
797 (arg, nargout, dim, nanflag, ismin);
798 else
799 retval = do_minmax_red_op<SparseComplexMatrix>
800 (arg, nargout, dim, nanflag, realabs, ismin);
801 }
802 else
803 {
804 if (cmethod)
805 retval = do_minmax_red_op<ComplexNDArray>
806 (arg, nargout, dim, nanflag, ismin);
807 else
808 retval = do_minmax_red_op<ComplexNDArray>
809 (arg, nargout, dim, nanflag, realabs, ismin);
810 }
811 }
812 break;
813
814 case btyp_float:
815 {
816 if (cmethod)
817 retval = do_minmax_red_op<FloatNDArray>
818 (arg, nargout, dim, nanflag, ismin);
819 else
820 retval = do_minmax_red_op<FloatNDArray>
821 (arg, nargout, dim, nanflag, realabs, ismin);
822 }
823 break;
824
826 {
827 if (cmethod)
828 retval = do_minmax_red_op<FloatComplexNDArray>
829 (arg, nargout, dim, nanflag, ismin);
830 else
831 retval = do_minmax_red_op<FloatComplexNDArray>
832 (arg, nargout, dim, nanflag, realabs, ismin);
833 }
834 break;
835
836 case btyp_char:
837 retval = do_minmax_red_op<charNDArray> (arg, nargout, dim, ismin);
838 break;
839
840#define MAKE_INT_BRANCH(X) \
841 case btyp_ ## X: \
842 { \
843 if (cmethod) \
844 retval = do_minmax_red_op<X ## NDArray> (arg, \
845 nargout, dim, ismin); \
846 else \
847 retval = do_minmax_red_op<X ## NDArray> (arg, \
848 nargout, dim, nanflag, realabs, ismin); \
849 } \
850 break;
851
852 MAKE_INT_BRANCH (int8);
853 MAKE_INT_BRANCH (int16);
854 MAKE_INT_BRANCH (int32);
855 MAKE_INT_BRANCH (int64);
856 MAKE_INT_BRANCH (uint8);
857 MAKE_INT_BRANCH (uint16);
858 MAKE_INT_BRANCH (uint32);
859 MAKE_INT_BRANCH (uint64);
860
861#undef MAKE_INT_BRANCH
862
863 case btyp_bool:
864 retval = do_minmax_red_op<boolNDArray> (arg, nargout, dim, ismin);
865 break;
866
867 default:
868 err_wrong_type_arg (fcn, arg);
869 }
870
871 if (do_perm)
872 {
873 retval(0) = retval(0).permute (perm_vec, true);
874 if (nargout > 1)
875 retval(1) = retval(1).permute (perm_vec, true);
876 }
877 }
878 else
879 {
880 if (nargout > 1)
881 error ("%s: two output arguments are not supported for two input arrays", fcn);
882 if (linear)
883 error ("%s: 'linear' is not supported for two input arrays", fcn);
884 octave_value argx = args(0);
885 octave_value argy = args(1);
886 builtin_type_t xtyp = argx.builtin_type ();
887 builtin_type_t ytyp = argy.builtin_type ();
888 builtin_type_t rtyp;
889 if (xtyp == btyp_char && ytyp == btyp_char)
890 rtyp = btyp_char;
891 // FIXME: This is what should happen when boolNDArray has max()
892 // else if (xtyp == btyp_bool && ytyp == btyp_bool)
893 // rtyp = btyp_bool;
894 else
895 rtyp = btyp_mixed_numeric (xtyp, ytyp);
896
897 switch (rtyp)
898 {
899 case btyp_double:
900 {
901 if ((argx.issparse ()
902 && (argy.issparse () || argy.is_scalar_type ()))
903 || (argy.issparse () && argx.is_scalar_type ()))
904 {
905 if (cmethod)
906 retval = do_minmax_bin_op<SparseMatrix>
907 (argx, argy, nanflag, ismin);
908 else
909 retval = do_minmax_bin_op<SparseMatrix>
910 (argx, argy, nanflag, realabs, ismin);
911 }
912 else
913 {
914 if (cmethod)
915 retval = do_minmax_bin_op<NDArray>
916 (argx, argy, nanflag, ismin);
917 else
918 retval = do_minmax_bin_op<NDArray>
919 (argx, argy, nanflag, realabs, ismin);
920 }
921 }
922 break;
923
924 case btyp_complex:
925 {
926 if ((argx.issparse ()
927 && (argy.issparse () || argy.is_scalar_type ()))
928 || (argy.issparse () && argx.is_scalar_type ()))
929 {
930 if (cmethod)
931 retval = do_minmax_bin_op<SparseComplexMatrix>
932 (argx, argy, nanflag, ismin);
933 else
934 retval = do_minmax_bin_op<SparseComplexMatrix>
935 (argx, argy, nanflag, realabs, ismin);
936 }
937 else
938 {
939 if (cmethod)
940 retval = do_minmax_bin_op<ComplexNDArray>
941 (argx, argy, nanflag, ismin);
942 else
943 retval = do_minmax_bin_op<ComplexNDArray>
944 (argx, argy, nanflag, realabs, ismin);
945 }
946 }
947 break;
948
949 case btyp_float:
950 {
951 if (cmethod)
952 retval = do_minmax_bin_op<FloatNDArray>
953 (argx, argy, nanflag, ismin);
954 else
955 retval = do_minmax_bin_op<FloatNDArray>
956 (argx, argy, nanflag, realabs, ismin);
957 }
958 break;
959
961 {
962 if (cmethod)
963 retval = do_minmax_bin_op<FloatComplexNDArray>
964 (argx, argy, nanflag, ismin);
965 else
966 retval = do_minmax_bin_op<FloatComplexNDArray>
967 (argx, argy, nanflag, realabs, ismin);
968 }
969 break;
970
971 case btyp_char:
972 retval = do_minmax_bin_op<charNDArray> (argx, argy, ismin);
973 break;
974
975#define MAKE_INT_BRANCH(X) \
976 case btyp_ ## X: \
977 { \
978 if (cmethod) \
979 retval = do_minmax_bin_op<X ## NDArray> \
980 (argx, argy, ismin); \
981 else \
982 retval = do_minmax_bin_op<X ## NDArray> \
983 (argx, argy, nanflag, realabs, ismin); \
984 } \
985 break;
986
987 MAKE_INT_BRANCH (int8);
988 MAKE_INT_BRANCH (int16);
989 MAKE_INT_BRANCH (int32);
990 MAKE_INT_BRANCH (int64);
991 MAKE_INT_BRANCH (uint8);
992 MAKE_INT_BRANCH (uint16);
993 MAKE_INT_BRANCH (uint32);
994 MAKE_INT_BRANCH (uint64);
995
996#undef MAKE_INT_BRANCH
997
998 // FIXME: This is what should happen when boolNDArray has max()
999 // case btyp_bool:
1000 // retval = do_minmax_bin_op<boolNDArray> (argx, argy, ismin);
1001 // break;
1002
1003 default:
1004 error ("%s: cannot compute %s (%s, %s)", fcn, fcn,
1005 argx.type_name ().c_str (), argy.type_name ().c_str ());
1006 }
1007
1008 // FIXME: Delete when boolNDArray has max()
1009 if (xtyp == btyp_bool && ytyp == btyp_bool)
1010 retval(0) = retval(0).bool_array_value ();
1011
1012 }
1013
1014 // Process "linear" option
1015 if (linear && nargout > 1 && ! allflag && ! args(0).isempty ())
1016 retval(1) = get_linear_index (retval(1), args(0), vecdim, false);
1017
1018 return retval;
1019}
1020
1021DEFUN (min, args, nargout,
1022 doc: /* -*- texinfo -*-
1023@deftypefn {} {@var{m} =} min (@var{x})
1024@deftypefnx {} {@var{m} =} min (@var{x}, [], @var{dim})
1025@deftypefnx {} {@var{m} =} min (@var{x}, [], @var{vecdim})
1026@deftypefnx {} {@var{m} =} min (@var{x}, [], "all")
1027@deftypefnx {} {@var{m} =} min (@var{x}, [], @var{nanflag})
1028@deftypefnx {} {[@var{m}, @var{im}] =} min (@dots{})
1029@deftypefnx {} {[@var{m}, @var{im}] =} min (@dots{}, "linear")
1030@deftypefnx {} {@var{m} =} min (@var{x}, @var{y})
1031@deftypefnx {} {@var{m} =} min (@var{x}, @var{y}, @var{nanflag})
1032@deftypefnx {} {@dots{} =} min (@dots{}, "ComparisonMethod", @var{method})
1033Find minimum values in the array @var{x}.
1034
1035If @var{x} is a vector, then @code{min (@var{x})} returns the minimum value of
1036the elements in @var{x}.
1037
1038If @var{x} is a matrix, then @code{min (@var{x})} returns a row vector with
1039each element containing the minimum value of the corresponding column in
1040@var{x}.
1041
1042If @var{x} is an array, then @code{min (@var{x})} computes the minimum value
1043along the first non-singleton dimension of @var{x}.
1044
1045The optional input @var{dim} specifies the dimension to operate on and must be
1046a positive integer. Specifying any singleton dimension of @var{x}, including
1047any dimension exceeding @code{ndims (@var{x})}, will return @var{x}.
1048
1049Specifying multiple dimensions with input @var{vecdim}, a vector of
1050non-repeating dimensions, will operate along the array slice defined by
1051@var{vecdim}. If @var{vecdim} indexes all dimensions of @var{x}, then it is
1052equivalent to the option @qcode{"all"}. Any dimension in @var{vecdim} greater
1053than @code{ndims (@var{x})} is ignored.
1054
1055Specifying the dimension as @qcode{"all"} will cause @code{min} to operate on
1056on all elements of @var{x}, and is equivalent to @code{min (@var{x}(:))}.
1057
1058If called with two output arguments, @code{min} also returns the first index of
1059the minimum value(s) in @var{x}. The second output argument is only valid when
1060@code{min} operates on a single input array. Setting the @qcode{"linear"} flag
1061returns the linear index to the corresponding minimum values in @var{x}.
1062
1063If called with two input arrays (@var{x} and @var{y}), @code{min} return the
1064pairwise minimum according to the rules for @ref{Broadcasting}.
1065
1066The optional variable @var{nanflag} specifies whether to include or exclude
1067@code{NaN} values from the calculation using any of the previously specified
1068input argument combinations. The default value for @var{nanflag} is
1069@qcode{"omitnan"} which ignores @code{NaN} values in the calculation. To
1070include @code{NaN} values, set the value of @var{nanflag} to
1071@qcode{"includenan"}, in which case @code{min} will return @code{NaN}, if any
1072element along the operating dimension is @code{NaN}.
1073
1074The optional @qcode{"ComparisonMethod"} paired argument specifies the
1075comparison method for numeric input and it applies to both one input and two
1076input arrays. @var{method} can take any of the following values:
1077
1078@table @asis
1079@item @qcode{'auto'}
1080This is the default method, which compares elements by @code{real (@var{x})}
1081when @var{x} is real, and by @code{abs (@var{x})} when @var{x} is complex.
1082
1083@item @qcode{'real'}
1084Compares elements by @code{real (@var{x})} when @var{x} is real or complex.
1085For elements with equal real parts, a second comparison by
1086@code{imag (@var{x})} is performed.
1087
1088@item @qcode{'abs'}
1089Compares elements by @code{abs (@var{x})} when @var{x} is real or complex. For
1090elements with equal magnitude, a second comparison by @code{angle (@var{x})} in
1091the interval from @math{-\pi} to @math{\pi} is performed.
1092@end table
1093@seealso{max, cummin, cummax}
1094@end deftypefn */)
1095{
1096 return do_minmax_body (args, nargout, true);
1097}
1098
1099/*
1100## Test generic double class
1101%!assert (min ([1, 4, 2, 3]), 1)
1102%!assert (min ([1; -10; 5; -2]), -10)
1103%!assert (min ([4, 2i 4.999; -2, 2, 3+4i]), [-2, 2, 4.999])
1104## Special routines for char arrays
1105%!assert (min (["abc", "ABC"]), 65)
1106%!assert (min (["abc"; "CBA"]), [67 66 65])
1107## Special routines for logical arrays
1108%!assert (min (logical ([])), logical ([]))
1109%!assert (min (logical ([0 0 1 0])), false)
1110%!assert (min (logical ([0 0 1 0; 0 1 1 0])), logical ([0 0 1 0]))
1111## Single values
1112%!assert (min (single ([1, 4, 2, 3])), single (1))
1113%!assert (min (single ([1; -10; 5; -2])), single (-10))
1114%!assert (min (single ([4, 2i 4.999; -2, 2, 3+4i])), single ([-2, 2, 4.999]))
1115## Integer values
1116%!assert (min (uint8 ([1, 4, 2, 3])), uint8 (1))
1117%!assert (min (uint8 ([1; -10; 5; -2])), uint8 (-10))
1118%!assert (min (int8 ([1, 4, 2, 3])), int8 (1))
1119%!assert (min (int8 ([1; -10; 5; -2])), int8 (-10))
1120%!assert (min (uint16 ([1, 4, 2, 3])), uint16 (1))
1121%!assert (min (uint16 ([1; -10; 5; -2])), uint16 (-10))
1122%!assert (min (int16 ([1, 4, 2, 3])), int16 (1))
1123%!assert (min (int16 ([1; -10; 5; -2])), int16 (-10))
1124%!assert (min (uint32 ([1, 4, 2, 3])), uint32 (1))
1125%!assert (min (uint32 ([1; -10; 5; -2])), uint32 (-10))
1126%!assert (min (int32 ([1, 4, 2, 3])), int32 (1))
1127%!assert (min (int32 ([1; -10; 5; -2])), int32 (-10))
1128%!assert (min (uint64 ([1, 4, 2, 3])), uint64 (1))
1129%!assert (min (uint64 ([1; -10; 5; -2])), uint64 (-10))
1130%!assert (min (int64 ([1, 4, 2, 3])), int64 (1))
1131%!assert (min (int64 ([1; -10; 5; -2])), int64 (-10))
1132## Sparse double values
1133%!assert (min (sparse ([1, 4, 2, 3])), sparse (1))
1134%!assert (min (sparse ([1; -10; 5; -2])), sparse(-10))
1135## Order sparse complex values by phase angle
1136%!test <*51307>
1137%! assert (min (sparse ([4, 2i 4.999; -2, 2, 3+4i])), sparse ([-2, 2, 4.999]));
1138%! assert (min (sparse ([4, -2i 4.99; -2, 2, 3+4i])), sparse ([-2, -2i, 4.99]));
1139%!assert (min (sparse ([4, 2i 4.999]), sparse ([-2, 2, 3+4i])),
1140%! sparse ([-2, 2, 4.999]))
1141%!assert (min (sparse ([4, -2i 4.999]), sparse ([-2, 2, 3+4i])),
1142%! sparse ([-2, -2i, 4.999]))
1143
1144## Test dimension argument
1145%!test
1146%! x = reshape (1:8, [2,2,2]);
1147%! assert (min (x, [], 1), reshape ([1, 3, 5, 7], [1,2,2]));
1148%! assert (min (x, [], 2), reshape ([1, 2, 5, 6], [2,1,2]));
1149%! [y, i] = min (x, [], 3);
1150%! assert (ndims (y), 2);
1151%! assert (y, [1, 3; 2, 4]);
1152%! assert (ndims (i), 2);
1153%! assert (i, [1, 1; 1, 1]);
1154
1155## Test vecdim argument
1156%!test
1157%! x = reshape (1:8, [2,2,2]);
1158%! assert (min (x, [], [1, 3]), [1, 3]);
1159%! assert (min (x, [], [2, 3]), [1; 2]);
1160%! assert (min (x, [], [1, 2]), reshape ([1, 5], [1,1,2]));
1161%! assert (min (x, [], [1, 2, 3]), min (x, [], "all"));
1162
1163## Test 2-output forms for various arg types
1164## Special routines for char arrays
1165%!test
1166%! [y, i] = min (["abc", "ABC"]);
1167%! assert (y, 65);
1168%! assert (i, 4);
1169## Special routines for logical arrays
1170%!test
1171%! x = logical ([0 0 1 0]);
1172%! [y, i] = min (x);
1173%! assert (y, false);
1174%! assert (i, 1);
1175## Special handling of ranges
1176%!test
1177%! rng = 1:2:10;
1178%! [y, i] = min (rng);
1179%! assert (y, 1);
1180%! assert (i, 1);
1181%! rng = 10:-2:1;
1182%! [y, i] = min (rng);
1183%! assert (y, 2);
1184%! assert (i, 5);
1185
1186## Test 2-input calling form for various arg types
1187## Test generic double class
1188%!test
1189%! x = [1, 2, 3, 4]; y = fliplr (x);
1190%! assert (min (x, y), [1 2 2 1]);
1191%! assert (min (x, 3), [1 2 3 3]);
1192%! assert (min (2, x), [1 2 2 2]);
1193%! assert (min (x, 2.1i), [1 2 2.1i 2.1i]);
1194## Test ordering of complex results with equal magnitude
1195%!test <*51307>
1196%! x = [1, 2, 3, 4];
1197%! assert (min (x, 2i), [1 2 2i 2i]);
1198%! assert (min (x, -2i), [1 -2i -2i -2i]);
1199## Special routines for char arrays
1200%!assert (min ("abc", "b"), [97 98 98])
1201%!assert (min ("b", "cba"), [98 98 97])
1202## Special handling for logical arrays
1203%!assert (min ([true false], false), [false false])
1204%!assert (min (true, [true false]), [true false])
1205## Single values
1206%!test
1207%! x = single ([1, 2, 3, 4]); y = fliplr (x);
1208%! assert (min (x, y), single ([1 2 2 1]));
1209%! assert (min (x, 3), single ([1 2 3 3]));
1210%! assert (min (2, x), single ([1 2 2 2]));
1211%! assert (min (x, 2.1i), single ([1 2 2.1i 2.1i]));
1212## Integer values
1213%!test
1214%! x = uint8 ([1, 2, 3, 4]); y = fliplr (x);
1215%! assert (min (x, y), uint8 ([1 2 2 1]));
1216%! assert (min (x, 3), uint8 ([1 2 3 3]));
1217%! assert (min (2, x), uint8 ([1 2 2 2]));
1218%! x = int8 ([1, 2, 3, 4]); y = fliplr (x);
1219%! assert (min (x, y), int8 ([1 2 2 1]));
1220%! assert (min (x, 3), int8 ([1 2 3 3]));
1221%! assert (min (2, x), int8 ([1 2 2 2]));
1222%! x = uint16 ([1, 2, 3, 4]); y = fliplr (x);
1223%! assert (min (x, y), uint16 ([1 2 2 1]));
1224%! assert (min (x, 3), uint16 ([1 2 3 3]));
1225%! assert (min (2, x), uint16 ([1 2 2 2]));
1226%! x = int16 ([1, 2, 3, 4]); y = fliplr (x);
1227%! assert (min (x, y), int16 ([1 2 2 1]));
1228%! assert (min (x, 3), int16 ([1 2 3 3]));
1229%! assert (min (2, x), int16 ([1 2 2 2]));
1230%! x = uint32 ([1, 2, 3, 4]); y = fliplr (x);
1231%! assert (min (x, y), uint32 ([1 2 2 1]));
1232%! assert (min (x, 3), uint32 ([1 2 3 3]));
1233%! assert (min (2, x), uint32 ([1 2 2 2]));
1234%! x = int32 ([1, 2, 3, 4]); y = fliplr (x);
1235%! assert (min (x, y), int32 ([1 2 2 1]));
1236%! assert (min (x, 3), int32 ([1 2 3 3]));
1237%! assert (min (2, x), int32 ([1 2 2 2]));
1238%! x = uint64 ([1, 2, 3, 4]); y = fliplr (x);
1239%! assert (min (x, y), uint64 ([1 2 2 1]));
1240%! assert (min (x, 3), uint64 ([1 2 3 3]));
1241%! assert (min (2, x), uint64 ([1 2 2 2]));
1242%! x = int64 ([1, 2, 3, 4]); y = fliplr (x);
1243%! assert (min (x, y), int64 ([1 2 2 1]));
1244%! assert (min (x, 3), int64 ([1 2 3 3]));
1245%! assert (min (2, x), int64 ([1 2 2 2]));
1246## Sparse double values
1247%!test
1248%! x = sparse ([1, 2, 3, 4]); y = fliplr (x);
1249%! assert (min (x, y), sparse ([1 2 2 1]));
1250%! assert (min (x, 3), sparse ([1 2 3 3]));
1251%! assert (min (2, x), sparse ([1 2 2 2]));
1252%! assert (min (x, 2.1i), sparse ([1 2 2.1i 2.1i]));
1253%!assert (min (sparse ([4, 2i, 4.999; -2, 2, 3+4i])), sparse ([-2, 2, 4.999]))
1254%!assert (min (sparse ([4, 2+i, 4.999; -2, 2, 3+4i])), sparse ([-2, 2, 4.999]))
1255%!assert (min (sparse ([4, 2i, 4.999; -2, 2-i, 3+4i])),
1256%! sparse ([-2, 2i, 4.999]))
1257
1258## Sparse with NaNs
1259%!assert (min (sparse ([NaN, 1, 2, 1])), sparse (1))
1260%!assert (min (sparse ([NaN; 1; 2; 1])), sparse (1))
1261%!assert (min (sparse ([0, NaN, 1, 2, 1])), sparse (0))
1262%!assert (min (sparse ([0; NaN; 1; 2; 1])), sparse (0))
1263%!assert (min (sparse ([0; NaN; 1; 2; 1]), [], "includenan"), sparse (NaN))
1264%!assert (min (sparse ([0; NaN; 1; 2; 1]), [], 1, "includenan"), sparse (NaN))
1265%!assert (min (sparse ([0; NaN; 1; 2; 1]), [], 2, "includenan"),
1266%! sparse ([0; NaN; 1; 2; 1]))
1267%!assert (min (sparse ([NaN, 1i, 2, 1])), sparse (1))
1268%!assert (min (sparse ([NaN, 1i, 2, 1]), "ComparisonMethod", "real"),
1269%! sparse (1i))
1270%!assert (min (sparse (single ([NaN, 1, 2, 1]))), sparse (1))
1271%!assert (min (sparse (single ([NaN; 1; 2; 1]))), sparse (1))
1272%!assert (min (sparse (single ([0, NaN, 1, 2, 1]))), sparse (0))
1273%!assert (min (sparse (single ([0; NaN; 1; 2; 1]))), sparse (0))
1274%!assert (min (sparse (single ([0; NaN; 1; 2; 1])), [], "includenan"),
1275%! sparse (NaN))
1276%!assert (min (sparse (single ([0; NaN; 1; 2; 1])), [], 1, "includenan"),
1277%! sparse (NaN))
1278%!assert (min (sparse (single ([NaN, 1i, 2, 1]))), sparse (1))
1279%!assert (min (sparse (single ([NaN, 1i, 2, 1])), "ComparisonMethod", "real"),
1280%! sparse (1i))
1281%!assert (min (sparse (single ([NaN, 1i, 2, 1]))), sparse (1))
1282
1283## Test broadcasting of empty matrices with min/max functions
1284%!assert (min (sparse (zeros (0,1)), sparse ([1, 2, 3, 4])),
1285%! sparse (zeros (0,4)))
1286%!error min (sparse (zeros (0,2)), sparse ([1, 2, 3, 4]))
1287%!assert (min (sparse (zeros (1,0)), sparse ([1; 2; 3; 4])),
1288%! sparse (zeros (4,0)))
1289%!error min (sparse (zeros (2,0)), sparse ([1; 2; 3; 4]))
1290%!assert (min (sparse (zeros (0,1)), sparse ([1, 2, 3, 4i])),
1291%! sparse (zeros (0,4)))
1292%!error min (sparse (zeros (0,2)), sparse ([1, 2, 3, 4i]))
1293%!assert (min (sparse (zeros (1,0)), sparse ([1; 2; 3; 4i])),
1294%! sparse (zeros (4,0)))
1295%!error min (sparse (zeros (2,0)), sparse ([1; 2; 3; 4i]))
1296
1297## NaNs and complex numbers
1298%!assert (min (NaN, 0), 0)
1299%!assert (min (NaN+1i, 0), 0)
1300%!assert (min (2+2i, 2-2i), 2-2i)
1301%!assert (min (2-2i, 2+2i), 2-2i)
1302%!test
1303%! [M1, I1] = min ([2+2i, 2-2i]);
1304%! [M2, I2] = min ([2-2i, 2+2i]);
1305%! assert (M1, M2);
1306%! assert (I1, 2);
1307%! assert (I2, 1);
1308%!assert (min (NaN, 0, "ComparisonMethod", "real"),
1309%! min (NaN+1i, 0, "ComparisonMethod", "real"))
1310%!assert (min (2+i, 1+10i), 2+1i)
1311%!assert (min (2+i, 1+10i, "ComparisonMethod", "real"), 1+10i)
1312%!assert (min (2+i, 1+10i, "ComparisonMethod", "abs"), 2+1i)
1313%!assert (min (2+i, -1+10i, "ComparisonMethod", "abs"), 2+1i)
1314%!assert (min (2+i, -2+i, "ComparisonMethod", "abs"), 2+1i)
1315%!assert (min (2+i, 2-i, "ComparisonMethod", "abs"), 2-1i)
1316%!assert (min (2+i, 2-i, "ComparisonMethod", "real"), 2-1i)
1317%!assert (min ([1i 2 -3 4]), 1i)
1318%!assert (min ([-2+i, 2-i]), 2-1i)
1319
1320## Test nanflag with dense arrays
1321%!test
1322%! [m,i] = min ([1,2,3;4,3,NaN;4,5,6], [], 2, "includenan");
1323%! assert (m, [1; NaN; 4]);
1324%! assert (i, [1; 3; 1]);
1325%! [m,i] = min ([1,2,3;4,NaN,NaN;4,5,6], [], 2, "includenan");
1326%! assert (m, [1; NaN; 4]);
1327%! assert (i, [1; 2; 1]);
1328%!test
1329%! x = magic (3);
1330%! x(2, 3) = NaN;
1331%! assert (min (x, [], 2, "includenan"), [1; NaN; 2]);
1332
1333## Test empty matrices and those with 0 dimensions (including catching errors)
1334%!assert (size (min (0, zeros (0, 1))), [0, 1])
1335%!assert (size (min ([], zeros (0, 1))), [0, 0])
1336%!assert (size (min (zeros (0, 1), [])), [0, 0])
1337%!assert (size (min ([], zeros (1, 0))), [0, 0])
1338%!assert (size (min (zeros (1, 0), [])), [0, 0])
1339%!error min ([1, 2, 3, 4], zeros (1, 0))
1340%!assert (size (min ([1, 2, 3, 4], zeros (0, 1))), [0, 4])
1341%!assert (min (0, sparse (zeros (1,0))), sparse (zeros (1, 0)))
1342%!assert (min (sparse (zeros (1,0)), 0), sparse (zeros (1, 0)))
1343%!error min (sparse (zeros (1,0)), sparse ([1, 2, 3, 4]))
1344%!error min ([1, 2, NaN, 4], zeros (1, 0), 'includenan')
1345%!assert (min ([1, 2, NaN, 4], zeros (0, 1), 'includenan'), zeros (0, 4))
1346%!error min (sparse (zeros (1,0)), sparse ([1, 2, NaN 4]), 'includenan')
1347%!assert (min (sparse (zeros (1,0)), sparse ([1; 2; 3; 4])),
1348%! sparse (zeros (4, 0)))
1349%!assert (min (zeros(0,3), zeros (0, 1)), zeros (0, 3))
1350%!error min (zeros(3, 0), zeros (0, 1))
1351%!assert (min (zeros(3, 0), zeros (1, 0)), zeros (3, 0))
1352%!error min (zeros(3, 0), zeros (0, 0))
1353%!error min (zeros(3, 0), [])
1354%!error min ([], zeros(3, 0))
1355%!error min (zeros(0,1), zeros(3, 0))
1356%!assert (min (zeros(1,0), zeros(3, 0)), zeros (3, 0))
1357
1358## Test "linear" option
1359%!shared x
1360%! x = repmat ([4,3,2,4,1,5,3,2], 3, 2, 5, 2);
1361%!test
1362%! [m, i] = min (x, [], [2, 3], 'linear');
1363%! assert (m, ones (3, 1, 1, 2));
1364%! assert (i(:,:,1,1), [13; 14; 15]);
1365%! assert (i(:,:,1,2), [253; 254; 255]);
1366%!test
1367%! [m, i] = min (x, [], [1, 3], 'linear');
1368%! assert (m, x(1,:,1,:));
1369%! assert (i(:,:,1,1), [1:3:46]);
1370%! assert (i(:,:,1,2), [241:3:286]);
1371%!test
1372%! [m, i] = min (x, [], [2, 4], 'linear');
1373%! assert (m, ones (3, 1, 5));
1374%! assert (i(:,:,1), [13; 14; 15]);
1375%! assert (i(:,:,2), [61; 62; 63]);
1376%! assert (i(:,:,3), [109; 110; 111]);
1377%! assert (i(:,:,4), [157; 158; 159]);
1378%! assert (i(:,:,5), [205; 206; 207]);
1379%!test
1380%! [m, i] = min (x, [], [1, 4], 'linear');
1381%! assert (m, x(1, :,:,1));
1382%! assert (i(:,:,1), [1:3:46]);
1383%! assert (i(:,:,2), [49:3:94]);
1384%! assert (i(:,:,3), [97:3:142]);
1385%! assert (i(:,:,4), [145:3:190]);
1386%! assert (i(:,:,5), [193:3:238]);
1387%!test
1388%! [m, i] = min (x, [], [1, 2, 3], 'linear');
1389%! assert (m, ones (1, 1, 1, 2));
1390%! assert (i(:,:,1,1), 13);
1391%! assert (i(:,:,1,2), 253);
1392%!test
1393%! [m, i] = min (x, [], [2, 3, 4], 'linear');
1394%! assert (m, [1; 1; 1]);
1395%! assert (i, [13; 14; 15]);
1396%!test
1397%! [m, i] = min ([1, 2, 3; 4, 5, 6], [], 1, 'linear');
1398%! assert (m, [1, 2, 3]);
1399%! assert (i, [1, 3, 5]);
1400%!test
1401%! [m, i] = min ([1, 2, 3; 4, 5, 6], [], 2, 'linear');
1402%! assert (m, [1; 4]);
1403%! assert (i, [1; 2]);
1404
1405## Test input validation
1406%!error <Invalid call> min ()
1407%!error <Invalid call> min (1, 2, 3, 4)
1408%!error <unrecognized optional argument 'foobar'> min (1, [], "foobar")
1409%!error <cannot set DIM or VECDIM with 'all' flag>
1410%! min (ones (3,3), [], 1, "all")
1411%!error <cannot set DIM or VECDIM with 'all' flag>
1412%! min (ones (3,3), [], [1, 2], "all")
1413%!error <invalid dimension DIM = 0> min (ones (3,3), [], 0)
1414%!error <invalid dimension DIM = -1> min (ones (3,3), [], -1)
1415%!error <invalid dimension in VECDIM = -2> min (ones (3,3), [], [1 -2])
1416%!error <duplicate dimension in VECDIM = 2> min (ones (3,3), [], [1 2 2])
1417%!error <duplicate dimension in VECDIM = 1> min (ones (3,3), [], [1 1 2])
1418%!warning <second argument is ignored> min ([1 2 3 4], 2, 2);
1419%!error <wrong type argument 'cell'> min ({1 2 3 4})
1420%!error <cannot compute min \‍(cell, scalar\‍)> min ({1, 2, 3}, 2)
1421%!error <two output arguments are not supported for two input arrays>
1422%! [m, i] = min ([], [])
1423%!error <'linear' is not supported for two input arrays>
1424%! min ([1 2 3 4], 1, "linear")
1425%! [m, i] = min ([1 2 3 4], [], "linear");
1426*/
1427
1428DEFUN (max, args, nargout,
1429 doc: /* -*- texinfo -*-
1430@deftypefn {} {@var{m} =} max (@var{x})
1431@deftypefnx {} {@var{m} =} max (@var{x}, [], @var{dim})
1432@deftypefnx {} {@var{m} =} max (@var{x}, [], @var{vecdim})
1433@deftypefnx {} {@var{m} =} max (@var{x}, [], "all")
1434@deftypefnx {} {@var{m} =} max (@var{x}, [], @var{nanflag})
1435@deftypefnx {} {[@var{m}, @var{im}] =} max (@dots{})
1436@deftypefnx {} {[@var{m}, @var{im}] =} max (@dots{}, "linear")
1437@deftypefnx {} {@var{m} =} max (@var{x}, @var{y})
1438@deftypefnx {} {@var{m} =} max (@var{x}, @var{y}, @var{nanflag})
1439@deftypefnx {} {@dots{} =} max (@dots{}, "ComparisonMethod", @var{method})
1440Find maximum values in the array @var{x}.
1441
1442If @var{x} is a vector, then @code{max (@var{x})} returns the maximum value of
1443the elements in @var{x}.
1444
1445If @var{x} is a matrix, then @code{max (@var{x})} returns a row vector with
1446each element containing the maximum value of the corresponding column in
1447@var{x}.
1448
1449If @var{x} is an array, then @code{max (@var{x})} computes the maximum value
1450along the first non-singleton dimension of @var{x}.
1451
1452The optional input @var{dim} specifies the dimension to operate on and must be
1453a positive integer. Specifying any singleton dimension of @var{x}, including
1454any dimension exceeding @code{ndims (@var{x})}, will return @var{x}.
1455
1456Specifying multiple dimensions with input @var{vecdim}, a vector of
1457non-repeating dimensions, will operate along the array slice defined by
1458@var{vecdim}. If @var{vecdim} indexes all dimensions of @var{x}, then it is
1459equivalent to the option @qcode{"all"}. Any dimension in @var{vecdim} greater
1460than @code{ndims (@var{x})} is ignored.
1461
1462Specifying the dimension as @qcode{"all"} will cause @code{max} to operate on
1463on all elements of @var{x}, and is equivalent to @code{max (@var{x}(:))}.
1464
1465If called with two output arguments, @code{max} also returns the first index of
1466the maximum value(s) in @var{x}. The second output argument is only valid when
1467@code{max} operates on a single input array. Setting the @qcode{"linear"} flag
1468returns the linear index to the corresponding minimum values in @var{x}.
1469
1470If called with two input arrays (@var{x} and @var{y}), @code{max} return the
1471pairwise maximum according to the rules for @ref{Broadcasting}.
1472
1473The optional variable @var{nanflag} specifies whether to include or exclude
1474@code{NaN} values from the calculation using any of the previously specified
1475input argument combinations. The default value for @var{nanflag} is
1476@qcode{"omitnan"} which ignores @code{NaN} values in the calculation. To
1477include @code{NaN} values, set the value of @var{nanflag} to
1478@qcode{"includenan"}, in which case @code{max} will return @code{NaN}, if any
1479element along the operating dimension is @code{NaN}.
1480
1481The optional @qcode{"ComparisonMethod"} paired argument specifies the
1482comparison method for numeric input and it applies to both one input and two
1483input arrays. @var{method} can take any of the following values:
1484
1485@table @asis
1486@item @qcode{'auto'}
1487This is the default method, which compares elements by @code{real (@var{x})}
1488when @var{x} is real, and by @code{abs (@var{x})} when @var{x} is complex.
1489
1490@item @qcode{'real'}
1491Compares elements by @code{real (@var{x})} when @var{x} is real or complex.
1492For elements with equal real parts, a second comparison by
1493@code{imag (@var{x})} is performed.
1494
1495@item @qcode{'abs'}
1496Compares elements by @code{abs (@var{x})} when @var{x} is real or complex. For
1497elements with equal magnitude, a second comparison by @code{angle (@var{x})} in
1498the interval from @math{-\pi} to @math{\pi} is performed.
1499@end table
1500@seealso{min, cummax, cummin}
1501@end deftypefn */)
1502{
1503 return do_minmax_body (args, nargout, false);
1504}
1505
1506/*
1507## Test generic double class
1508%!assert (max ([1, 4, 2, 3]), 4)
1509%!assert (max ([1; -10; 5; -2]), 5)
1510%!assert (max ([4, 2i 4.999; -2, 2, 3+4i]), [4, 2i, 3+4i])
1511## Special routines for char arrays
1512%!assert (max (["abc", "ABC"]), 99)
1513%!assert (max (["abc"; "CBA"]), [97 98 99])
1514## Special routines for logical arrays
1515%!assert (max (logical ([])), logical ([]))
1516%!assert (max (logical ([0 0 1 0])), true)
1517%!assert (max (logical ([0 0 1 0; 0 1 0 0])), logical ([0 1 1 0]))
1518## Single values
1519%!assert (max (single ([1, 4, 2, 3])), single (4))
1520%!assert (max (single ([1; -10; 5; -2])), single (5))
1521%!assert (max (single ([4, 2i 4.999; -2, 2, 3+4i])), single ([4, 2i, 3+4i]))
1522## Integer values
1523%!assert (max (uint8 ([1, 4, 2, 3])), uint8 (4))
1524%!assert (max (uint8 ([1; -10; 5; -2])), uint8 (5))
1525%!assert (max (int8 ([1, 4, 2, 3])), int8 (4))
1526%!assert (max (int8 ([1; -10; 5; -2])), int8 (5))
1527%!assert (max (uint16 ([1, 4, 2, 3])), uint16 (4))
1528%!assert (max (uint16 ([1; -10; 5; -2])), uint16 (5))
1529%!assert (max (int16 ([1, 4, 2, 3])), int16 (4))
1530%!assert (max (int16 ([1; -10; 5; -2])), int16 (5))
1531%!assert (max (uint32 ([1, 4, 2, 3])), uint32 (4))
1532%!assert (max (uint32 ([1; -10; 5; -2])), uint32 (5))
1533%!assert (max (int32 ([1, 4, 2, 3])), int32 (4))
1534%!assert (max (int32 ([1; -10; 5; -2])), int32 (5))
1535%!assert (max (uint64 ([1, 4, 2, 3])), uint64 (4))
1536%!assert (max (uint64 ([1; -10; 5; -2])), uint64 (5))
1537%!assert (max (int64 ([1, 4, 2, 3])), int64 (4))
1538%!assert (max (int64 ([1; -10; 5; -2])), int64 (5))
1539## Sparse double values
1540%!assert (max (sparse ([1, 4, 2, 3])), sparse (4))
1541%!assert (max (sparse ([1; -10; 5; -2])), sparse(5))
1542## Order sparse complex values by phase angle
1543%!test <*51307>
1544%! assert (max (sparse ([4, 2i 4.999; -2, 2, 3+4i])), sparse ([4, 2i, 3+4i]));
1545%! assert (max (sparse ([4, -2i 4.999; -2, 2, 3+4i])), sparse ([4, 2, 3+4i]));
1546%!assert (max (sparse ([4, 2i 4.999]), sparse ([-2, 2, 3+4i])),
1547%! sparse ([4, 2i, 3+4i]))
1548%!assert (max (sparse ([4, -2i 4.999]), sparse ([-2, 2, 3+4i])),
1549%! sparse ([4, 2, 3+4i]))
1550
1551## Test dimension argument
1552%!test
1553%! x = reshape (1:8, [2,2,2]);
1554%! assert (max (x, [], 1), reshape ([2, 4, 6, 8], [1,2,2]));
1555%! assert (max (x, [], 2), reshape ([3, 4, 7, 8], [2,1,2]));
1556%! [y, i] = max (x, [], 3);
1557%! assert (ndims (y), 2);
1558%! assert (y, [5, 7; 6, 8]);
1559%! assert (ndims (i), 2);
1560%! assert (i, [2, 2; 2, 2]);
1561
1562## Test vecdim argument
1563%!test
1564%! x = reshape (1:8, [2,2,2]);
1565%! assert (max (x, [], [1, 3]), [6, 8]);
1566%! assert (max (x, [], [2, 3]), [7; 8]);
1567%! assert (max (x, [], [1, 2]), reshape ([4, 8], [1,1,2]));
1568%! assert (max (x, [], [1, 2, 3]), max (x, [], "all"));
1569
1570## Test 2-output forms for various arg types
1571## Special routines for char arrays
1572%!test
1573%! [y, i] = max (["abc", "ABC"]);
1574%! assert (y, 99);
1575%! assert (i, 3);
1576## Special routines for logical arrays
1577%!test
1578%! x = logical ([0 0 1 0]);
1579%! [y, i] = max (x);
1580%! assert (y, true);
1581%! assert (i, 3);
1582## Special handling of ranges
1583%!test
1584%! rng = 1:2:10;
1585%! [y, i] = max (rng);
1586%! assert (y, 9);
1587%! assert (i, 5);
1588%! rng = 10:-2:1;
1589%! [y, i] = max (rng);
1590%! assert (y, 10);
1591%! assert (i, 1);
1592
1593## Test 2-input calling form for various arg types
1594## Test generic double class
1595%!test
1596%! x = [1, 2, 3, 4]; y = fliplr (x);
1597%! assert (max (x, y), [4 3 3 4]);
1598%! assert (max (x, 3), [3 3 3 4]);
1599%! assert (max (2, x), [2 2 3 4]);
1600%! assert (max (x, 2.1i), [2.1i 2.1i 3 4]);
1601## Test ordering of complex results with equal magnitude
1602%!test <*51307>
1603%! x = [1, 2, 3, 4];
1604%! assert (max (x, 2i), [2i 2i 3 4]);
1605%! assert (max (x, -2i), [-2i 2 3 4]);
1606## Special routines for char arrays
1607%!assert (max ("abc", "b"), [98 98 99])
1608%!assert (max ("b", "cba"), [99 98 98])
1609## Special handling for logical arrays
1610%!assert (max ([true false], false), [true false])
1611%!assert (max (true, [false false]), [true true])
1612## Single values
1613%!test
1614%! x = single ([1, 2, 3, 4]); y = fliplr (x);
1615%! assert (max (x, y), single ([4 3 3 4]));
1616%! assert (max (x, 3), single ([3 3 3 4]));
1617%! assert (max (2, x), single ([2 2 3 4]));
1618%! assert (max (x, 2.1i), single ([2.1i 2.1i 3 4]));
1619## Integer values
1620%!test
1621%! x = uint8 ([1, 2, 3, 4]); y = fliplr (x);
1622%! assert (max (x, y), uint8 ([4 3 3 4]));
1623%! assert (max (x, 3), uint8 ([3 3 3 4]));
1624%! assert (max (2, x), uint8 ([2 2 3 4]));
1625%! x = int8 ([1, 2, 3, 4]); y = fliplr (x);
1626%! assert (max (x, y), int8 ([4 3 3 4]));
1627%! assert (max (x, 3), int8 ([3 3 3 4]));
1628%! assert (max (2, x), int8 ([2 2 3 4]));
1629%! x = uint16 ([1, 2, 3, 4]); y = fliplr (x);
1630%! assert (max (x, y), uint16 ([4 3 3 4]));
1631%! assert (max (x, 3), uint16 ([3 3 3 4]));
1632%! assert (max (2, x), uint16 ([2 2 3 4]));
1633%! x = int16 ([1, 2, 3, 4]); y = fliplr (x);
1634%! assert (max (x, y), int16 ([4 3 3 4]));
1635%! assert (max (x, 3), int16 ([3 3 3 4]));
1636%! assert (max (2, x), int16 ([2 2 3 4]));
1637%! x = uint32 ([1, 2, 3, 4]); y = fliplr (x);
1638%! assert (max (x, y), uint32 ([4 3 3 4]));
1639%! assert (max (x, 3), uint32 ([3 3 3 4]));
1640%! assert (max (2, x), uint32 ([2 2 3 4]));
1641%! x = int32 ([1, 2, 3, 4]); y = fliplr (x);
1642%! assert (max (x, y), int32 ([4 3 3 4]));
1643%! assert (max (x, 3), int32 ([3 3 3 4]));
1644%! assert (max (2, x), int32 ([2 2 3 4]));
1645%! x = uint64 ([1, 2, 3, 4]); y = fliplr (x);
1646%! assert (max (x, y), uint64 ([4 3 3 4]));
1647%! assert (max (x, 3), uint64 ([3 3 3 4]));
1648%! assert (max (2, x), uint64 ([2 2 3 4]));
1649%! x = int64 ([1, 2, 3, 4]); y = fliplr (x);
1650%! assert (max (x, y), int64 ([4 3 3 4]));
1651%! assert (max (x, 3), int64 ([3 3 3 4]));
1652%! assert (max (2, x), int64 ([2 2 3 4]));
1653## Sparse double values
1654%!test
1655%! x = sparse ([1, 2, 3, 4]); y = fliplr (x);
1656%! assert (max (x, y), sparse ([4 3 3 4]));
1657%! assert (max (x, 3), sparse ([3 3 3 4]));
1658%! assert (max (2, x), sparse ([2 2 3 4]));
1659%! assert (max (x, 2.1i), sparse ([2.1i 2.1i 3 4]));
1660%!assert (max (sparse ([4, 2i, 4.999; -2, 2, 3+4i])), sparse ([4, 2i, 3+4i]))
1661%!assert (max (sparse ([4, 2+i, 4.999; -2, 2, 3+4i])), sparse ([4, 2+i, 3+4i]))
1662%!assert (max (sparse ([4, 2i, 4.999; -2, 2-i, 3+4i])), sparse ([4, 2-i, 3+4i]))
1663
1664## Test for bug #40743
1665%!assert <*40743> (max (zeros (1,0), ones (1,1)), zeros (1,0))
1666%!assert <*40743> (max (sparse (zeros (1,0)), sparse (ones (1,1))),
1667%! sparse (zeros (1,0)))
1668
1669## Sparse with NaNs
1670%!assert (max (sparse ([NaN, 1, 2, 1])), sparse (2))
1671%!assert (max (sparse ([NaN; 1; 2; 1])), sparse (2))
1672%!assert (max (sparse ([0, NaN, 1, 2, 1])), sparse (2))
1673%!assert (max (sparse ([0; NaN; 1; 2; 1])), sparse (2))
1674%!assert (max (sparse ([0; NaN; 1; 2; 1]), [], "includenan"), sparse (NaN))
1675%!assert (max (sparse ([0; NaN; 1; 2; 1]), [], 1, "includenan"), sparse (NaN))
1676%!assert (max (sparse ([0; NaN; 1; 2; 1]), [], 2, "includenan"),
1677%! sparse ([0; NaN; 1; 2; 1]))
1678%!assert (max (sparse ([NaN, 1i, 2, 1])), sparse (2))
1679%!assert (max (sparse ([NaN, 1i, 2, 1]), "ComparisonMethod", "real"),
1680%! sparse (2))
1681%!assert (max (sparse (single ([NaN, 1, 2, 1]))), sparse (2))
1682%!assert (max (sparse (single ([NaN; 1; 2; 1]))), sparse (2))
1683%!assert (max (sparse (single ([0, NaN, 1, 2, 1]))), sparse (2))
1684%!assert (max (sparse (single ([0; NaN; 1; 2; 1]))), sparse (2))
1685%!assert (max (sparse (single ([0; NaN; 1; 2; 1])), [], "includenan"),
1686%! sparse (NaN))
1687%!assert (max (sparse (single ([0; NaN; 1; 2; 1])), [], 1, "includenan"),
1688%! sparse (NaN))
1689%!assert (max (sparse (single ([NaN, 1i, 2, 1]))), sparse (2))
1690%!assert (max (sparse (single ([NaN, 1i, 2, 1])), "ComparisonMethod", "real"),
1691%! sparse (2))
1692%!assert (max (sparse (single ([NaN, 1i, 1]))), sparse (1i))
1693
1694## Test broadcasting of empty matrices with min/max functions
1695%!assert (max (sparse (zeros (0,1)), sparse ([1, 2, 3, 4])),
1696%! sparse (zeros (0,4)))
1697%!error max (sparse (zeros (0,2)), sparse ([1, 2, 3, 4]))
1698%!assert (max (sparse (zeros (1,0)), sparse ([1; 2; 3; 4])),
1699%! sparse (zeros (4,0)))
1700%!error max (sparse (zeros (2,0)), sparse ([1; 2; 3; 4]))
1701%!assert (max (sparse (zeros (0,1)), sparse ([1, 2, 3, 4i])),
1702%! sparse (zeros (0,4)))
1703%!error max (sparse (zeros (0,2)), sparse ([1, 2, 3, 4i]))
1704%!assert (max (sparse (zeros (1,0)), sparse ([1; 2; 3; 4i])),
1705%! sparse (zeros (4,0)))
1706%!error max (sparse (zeros (2,0)), sparse ([1; 2; 3; 4i]))
1707
1708## NaNs and complex numbers
1709%!assert (max (NaN, 0), 0)
1710%!assert (max (NaN+1i, 0), 0)
1711%!assert (max (2+2i, 2-2i), 2+2i)
1712%!assert (max (2-2i, 2+2i), 2+2i)
1713%!test
1714%! [M1, I1] = max ([2+2i, 2-2i]);
1715%! [M2, I2] = max ([2-2i, 2+2i]);
1716%! assert (M1, M2);
1717%! assert (I1, 1);
1718%! assert (I2, 2);
1719%!assert (max (NaN, 0, "ComparisonMethod", "real"),
1720%! max (NaN+1i, 0, "ComparisonMethod", "real"))
1721%!assert (max (2+i, 1+10i), 1+10i)
1722%!assert (max (2+i, 1+10i, "ComparisonMethod", "real"), 2+i)
1723%!assert (max (2+i, 1+10i, "ComparisonMethod", "abs"), 1+10i)
1724%!assert (max (2+i, -1+10i, "ComparisonMethod", "abs"), -1+10i)
1725%!assert (max (2+i, -2+i, "ComparisonMethod", "abs"), -2+1i)
1726%!assert (max (2+i, 2-i, "ComparisonMethod", "abs"), 2+1i)
1727%!assert (max (2+i, 2-i, "ComparisonMethod", "real"), 2+1i)
1728%!assert (max ([1i 2 -3 4]), 4)
1729%!assert (max ([-2+i, 2-i]), -2+1i)
1730
1731## Test nanflag with dense arrays
1732%!test
1733%! [m,i] = max ([1,2,3;4,3,NaN;4,5,6], [], 2, "includenan");
1734%! assert (m, [3; NaN; 6]);
1735%! assert (i, [3; 3; 3]);
1736%! [m,i] = max ([1,2,3;4,NaN,NaN;4,5,6], [], 2, "includenan");
1737%! assert (m, [3; NaN; 6]);
1738%! assert (i, [3; 2; 3]);
1739%!test
1740%! x = magic (3);
1741%! x(2, 3) = NaN;
1742%! assert (max (x, [], 2, "includenan"), [8; NaN; 9]);
1743
1744## Test empty matrices and those with 0 dimensions (including catching errors)
1745%!assert (size (max (0, zeros (0, 1))), [0, 1])
1746%!assert (size (max ([], zeros (0, 1))), [0, 0])
1747%!assert (size (max (zeros (0, 1), [])), [0, 0])
1748%!assert (size (max ([], zeros (1, 0))), [0, 0])
1749%!assert (size (max (zeros (1, 0), [])), [0, 0])
1750%!error max ([1, 2, 3, 4], zeros (1, 0))
1751%!assert (size (max ([1, 2, 3, 4], zeros (0, 1))), [0, 4])
1752%!assert (max (0, sparse (zeros (1,0))), sparse (zeros (1, 0)))
1753%!assert (max (sparse (zeros (1,0)), 0), sparse (zeros (1, 0)))
1754%!error max (sparse (zeros (1,0)), sparse ([1, 2, 3, 4]))
1755%!error max ([1, 2, NaN, 4], zeros (1, 0), 'includenan')
1756%!assert (max ([1, 2, NaN, 4], zeros (0, 1), 'includenan'), zeros (0, 4))
1757%!error max (sparse (zeros (1,0)), sparse ([1, 2, NaN 4]), 'includenan')
1758%!assert (max (sparse (zeros (1,0)), sparse ([1; 2; 3; 4])),
1759%! sparse (zeros (4, 0)))
1760%!assert (max (zeros(0,3), zeros (0, 1)), zeros (0, 3))
1761%!error max (zeros(3, 0), zeros (0, 1))
1762%!assert (max (zeros(3, 0), zeros (1, 0)), zeros (3, 0))
1763%!error max (zeros(3, 0), zeros (0, 0))
1764%!error max (zeros(3, 0), [])
1765%!error max ([], zeros(3, 0))
1766%!error max (zeros(0,1), zeros(3, 0))
1767%!assert (max (zeros(1,0), zeros(3, 0)), zeros (3, 0))
1768
1769## Test "linear" option
1770%!shared x
1771%! x = repmat ([4,3,2,4,1,5,3,2], 3, 2, 5, 2);
1772%!test
1773%! [m, i] = max (x, [], [2, 3], 'linear');
1774%! assert (m, 5 * ones (3, 1, 1, 2));
1775%! assert (i(:,:,1,1), [16; 17; 18]);
1776%! assert (i(:,:,1,2), [256; 257; 258]);
1777%!test
1778%! [m, i] = max (x, [], [1, 3], 'linear');
1779%! assert (m, x(1,:,1,:));
1780%! assert (i(:,:,1,1), [1:3:46]);
1781%! assert (i(:,:,1,2), [241:3:286]);
1782%!test
1783%! [m, i] = max (x, [], [2, 4], 'linear');
1784%! assert (m, 5 * ones (3, 1, 5));
1785%! assert (i(:,:,1), [16; 17; 18]);
1786%! assert (i(:,:,2), [64; 65; 66]);
1787%! assert (i(:,:,3), [112; 113; 114]);
1788%! assert (i(:,:,4), [160; 161; 162]);
1789%! assert (i(:,:,5), [208; 209; 210]);
1790%!test
1791%! [m, i] = max (x, [], [1, 4], 'linear');
1792%! assert (m, x(1, :,:,1));
1793%! assert (i(:,:,1), [1:3:46]);
1794%! assert (i(:,:,2), [49:3:94]);
1795%! assert (i(:,:,3), [97:3:142]);
1796%! assert (i(:,:,4), [145:3:190]);
1797%! assert (i(:,:,5), [193:3:238]);
1798%!test
1799%! [m, i] = max (x, [], [1, 2, 3], 'linear');
1800%! assert (m, 5 * ones (1, 1, 1, 2));
1801%! assert (i(:,:,1,1), 16);
1802%! assert (i(:,:,1,2), 256);
1803%!test
1804%! [m, i] = max (x, [], [2, 3, 4], 'linear');
1805%! assert (m, [5; 5; 5]);
1806%! assert (i, [16; 17; 18]);
1807%!test
1808%! [m, i] = max ([1, 2, 3; 4, 5, 6], [], 1, 'linear');
1809%! assert (m, [4, 5, 6]);
1810%! assert (i, [2, 4, 6]);
1811%!test
1812%! [m, i] = max ([1, 2, 3; 4, 5, 6], [], 2, 'linear');
1813%! assert (m, [3; 6]);
1814%! assert (i, [5; 6]);
1815
1816## Test input validation
1817%!error <Invalid call> max ()
1818%!error <Invalid call> max (1, 2, 3, 4)
1819%!error <unrecognized optional argument 'foobar'> max (1, [], "foobar")
1820%!error <cannot set DIM or VECDIM with 'all' flag>
1821%! max (ones (3,3), [], 1, "all")
1822%!error <cannot set DIM or VECDIM with 'all' flag>
1823%! max (ones (3,3), [], [1, 2], "all")
1824%!error <invalid dimension DIM = 0> max (ones (3,3), [], 0)
1825%!error <invalid dimension DIM = -1> max (ones (3,3), [], -1)
1826%!error <invalid dimension in VECDIM = -2> max (ones (3,3), [], [1 -2])
1827%!error <duplicate dimension in VECDIM = 2> max (ones (3,3), [], [1 2 2])
1828%!error <duplicate dimension in VECDIM = 1> max (ones (3,3), [], [1 1 2])
1829%!warning <second argument is ignored> max ([1 2 3 4], 2, 2);
1830%!error <wrong type argument 'cell'> max ({1 2 3 4})
1831%!error <cannot compute max \‍(cell, scalar\‍)> max ({1, 2, 3}, 2)
1832%!error <two output arguments are not supported for two input arrays>
1833%! [m, i] = max ([], [])
1834%!error <'linear' is not supported for two input arrays>
1835%! max ([1 2 3 4], 1, "linear")
1836%! [m ,i] = max ([1 2 3 4], [], "linear");
1837*/
1838
1840reverse_index (Array<octave_idx_type> idx, const dim_vector array_size, int dim)
1841{
1842 octave_idx_type n = idx.numel ();
1843 if (dim == -1)
1844 dim = array_size.first_non_singleton ();
1845 octave_idx_type dim_size = array_size(dim);
1846 for (octave_idx_type i = 0; i < n; i++)
1847 idx(i) = dim_size - idx(i) - 1;
1848 return idx;
1849}
1850
1851template <typename ArrayType>
1852static octave_value_list
1853do_cumminmax_red_op (const octave_value& arg, int nargout, int dim,
1854 bool nanflag, bool ismin, bool direction)
1855{
1856 octave_value_list retval (nargout > 1 ? 2 : 1);
1857 ArrayType array = octave_value_extract<ArrayType> (arg);
1858
1859 if (nargout <= 1)
1860 {
1861 if (ismin)
1862 {
1863 if (direction)
1864 retval(0) = array.flip (dim).cummin (dim, nanflag).flip (dim);
1865 else
1866 retval(0) = array.cummin (dim, nanflag);
1867 }
1868 else
1869 {
1870 if (direction)
1871 retval(0) = array.flip (dim).cummax (dim, nanflag).flip (dim);
1872 else
1873 retval(0) = array.cummax (dim, nanflag);
1874 }
1875 }
1876 else
1877 {
1878 retval.resize (2);
1880 if (ismin)
1881 {
1882 if (direction)
1883 {
1884 retval(0) = array.flip (dim).cummin (idx, dim, nanflag).flip (dim);
1885 idx = reverse_index (idx, array.dims (), dim);
1886 }
1887 else
1888 retval(0) = array.cummin (idx, dim, nanflag);
1889 }
1890 else
1891 {
1892 if (direction)
1893 {
1894 retval(0) = array.flip (dim).cummax (idx, dim, nanflag).flip (dim);
1895 idx = reverse_index (idx, array.dims (), dim);
1896 }
1897 else
1898 retval(0) = array.cummax (idx, dim, nanflag);
1899 }
1900
1901 retval(1) = octave_value (idx, true, true);
1902 if (direction)
1903 retval(1) = retval(1).array_value ().flip (dim);
1904 }
1905
1906 return retval;
1907}
1908
1909template <typename ArrayType>
1910static octave_value_list
1911do_cumminmax_red_op (const octave_value& arg, int nargout, int dim,
1912 bool nanflag, bool realabs, bool ismin, bool direction)
1913{
1914 octave_value_list retval (nargout > 1 ? 2 : 1);
1915 ArrayType array = octave_value_extract<ArrayType> (arg);
1916
1917 if (nargout <= 1)
1918 {
1919 if (ismin)
1920 {
1921 if (direction)
1922 retval(0) = array.flip (dim).cummin (dim, nanflag, realabs).flip (dim);
1923 else
1924 retval(0) = array.cummin (dim, nanflag, realabs);
1925 }
1926 else
1927 {
1928 if (direction)
1929 retval(0) = array.flip (dim).cummax (dim, nanflag, realabs).flip (dim);
1930 else
1931 retval(0) = array.cummax (dim, nanflag, realabs);
1932 }
1933 }
1934 else
1935 {
1936 retval.resize (2);
1938 if (ismin)
1939 {
1940 if (direction)
1941 {
1942 retval(0) = array.flip (dim).cummin (idx, dim, nanflag, realabs).flip (dim);
1943 idx = reverse_index (idx, array.dims (), dim);
1944 }
1945 else
1946 retval(0) = array.cummin (idx, dim, nanflag, realabs);
1947 }
1948 else
1949 {
1950 if (direction)
1951 {
1952 retval(0) = array.flip (dim).cummax (idx, dim, nanflag, realabs).flip (dim);
1953 idx = reverse_index (idx, array.dims (), dim);
1954 }
1955 else
1956 retval(0) = array.cummax (idx, dim, nanflag, realabs);
1957 }
1958
1959 retval(1) = octave_value (idx, true, true);
1960 if (direction)
1961 retval(1) = retval(1).array_value ().flip (dim);
1962 }
1963
1964 return retval;
1965}
1966
1967static octave_value_list
1968do_cumminmax_body (const octave_value_list& args,
1969 int nargout, bool ismin)
1970{
1971 int nargin = args.length ();
1972
1973 const char *fcn = (ismin ? "cummin" : "cummax");
1974
1975 bool direction = false;
1976 bool do_perm = false;
1977 bool allflag = false;
1978 bool nanflag = true; // NaNs are ignored by default
1979 bool cmethod = true; // resolves to "auto"
1980 bool realabs = true;
1981 bool linear = false;
1982
1983 // Get "ComparisonMethod" optional argument first
1984 if (nargin > 2
1985 && args(nargin - 1).is_string () && args(nargin - 2).is_string ())
1986 {
1987 std::string name = args(nargin - 2).string_value ();
1988 std::string sval = args(nargin - 1).string_value ();
1989 if (name == "ComparisonMethod" || name == "comparisonmethod")
1990 {
1991 if (sval == "auto")
1992 cmethod = true;
1993 else if (sval == "real")
1994 cmethod = false;
1995 else if (sval == "abs")
1996 {
1997 cmethod = false;
1998 realabs = false;
1999 }
2000 else
2001 error ("%s: invalid comparison method '%s'", fcn, sval.c_str ());
2002
2003 nargin -= 2;
2004 }
2005 }
2006
2007 while (nargin > 1 && args(nargin - 1).is_string ())
2008 {
2009 std::string str = args(nargin - 1).string_value ();
2010
2011 if (str == "forward")
2012 direction = false;
2013 else if (str == "reverse")
2014 direction = true;
2015 else if (str == "all")
2016 allflag = true;
2017 else if (str == "omitnan" || str == "omitmissing")
2018 {
2019 if (args(0).is_double_type () || args(0).is_single_type ())
2020 nanflag = true;
2021 }
2022 else if (str == "includenan" || str == "includemissing")
2023 nanflag = false;
2024 else if (str == "linear")
2025 linear = true;
2026 else
2027 error ("%s: unrecognized optional argument '%s'", fcn, str.c_str ());
2028
2029 nargin--;
2030 }
2031
2032 if (nargin < 1 || nargin > 2)
2033 print_usage ();
2034 if (allflag && nargin > 1)
2035 error ("%s: cannot set DIM or VECDIM with 'all' flag", fcn);
2036
2037 octave_value arg = args(0);
2038 Array<int> vecdim;
2039
2040 // Handle DIM, VECDIM
2041 int dim = -1;
2042 Array<int> perm_vec;
2043 if (nargin == 2)
2044 {
2045 octave_value dimarg = args(1);
2046 get_dim_vecdim_all (dimarg, arg, dim, perm_vec, do_perm, allflag, fcn);
2047 vecdim = dimarg.int_vector_value ();
2048 }
2049 else
2050 {
2051 dim_vector dims = arg.dims ();
2052 vecdim.resize (dim_vector (1, 1));
2053 vecdim(0) = dims.first_non_singleton () + 1;
2054 }
2055
2056 // Handle allflag
2057 if (allflag)
2058 arg = arg.reshape (dim_vector (arg.numel (), 1));
2059
2060 octave_value_list retval;
2061
2062 switch (arg.builtin_type ())
2063 {
2064 case btyp_double:
2065 {
2066 if (cmethod)
2067 retval = do_cumminmax_red_op<NDArray>
2068 (arg, nargout, dim, nanflag, ismin, direction);
2069 else
2070 retval = do_cumminmax_red_op<NDArray>
2071 (arg, nargout, dim, nanflag, realabs, ismin, direction);
2072 }
2073 break;
2074
2075 case btyp_complex:
2076 {
2077 if (cmethod)
2078 retval = do_cumminmax_red_op<ComplexNDArray>
2079 (arg, nargout, dim, nanflag, ismin, direction);
2080 else
2081 retval = do_cumminmax_red_op<ComplexNDArray>
2082 (arg, nargout, dim, nanflag, realabs, ismin, direction);
2083 }
2084 break;
2085
2086 case btyp_float:
2087 {
2088 if (cmethod)
2089 retval = do_cumminmax_red_op<FloatNDArray>
2090 (arg, nargout, dim, nanflag, ismin, direction);
2091 else
2092 retval = do_cumminmax_red_op<FloatNDArray>
2093 (arg, nargout, dim, nanflag, realabs, ismin, direction);
2094 }
2095 break;
2096
2097 case btyp_float_complex:
2098 {
2099 if (cmethod)
2100 retval = do_cumminmax_red_op<FloatComplexNDArray>
2101 (arg, nargout, dim, nanflag, ismin, direction);
2102 else
2103 retval = do_cumminmax_red_op<FloatComplexNDArray>
2104 (arg, nargout, dim, nanflag, realabs, ismin, direction);
2105 }
2106 break;
2107
2108#define MAKE_INT_BRANCH(X) \
2109 case btyp_ ## X: \
2110 { \
2111 if (cmethod) \
2112 retval = do_cumminmax_red_op<X ## NDArray> (arg, nargout, \
2113 dim, nanflag, ismin, direction); \
2114 else \
2115 retval = do_cumminmax_red_op<X ## NDArray> (arg, nargout, \
2116 dim, nanflag, realabs, ismin, direction); \
2117 } \
2118 break;
2119
2120 MAKE_INT_BRANCH (int8);
2121 MAKE_INT_BRANCH (int16);
2122 MAKE_INT_BRANCH (int32);
2123 MAKE_INT_BRANCH (int64);
2124 MAKE_INT_BRANCH (uint8);
2125 MAKE_INT_BRANCH (uint16);
2126 MAKE_INT_BRANCH (uint32);
2127 MAKE_INT_BRANCH (uint64);
2128
2129#undef MAKE_INT_BRANCH
2130
2131 case btyp_bool:
2132 {
2133 if (cmethod)
2134 retval = do_cumminmax_red_op<int8NDArray>
2135 (arg, nargout, dim, nanflag, ismin, direction);
2136 else
2137 retval = do_cumminmax_red_op<int8NDArray>
2138 (arg, nargout, dim, nanflag, realabs, ismin, direction);
2139 if (retval.length () > 0)
2140 retval(0) = retval(0).bool_array_value ();
2141 }
2142 break;
2143
2144 default:
2145 err_wrong_type_arg (fcn, arg);
2146 }
2147
2148 if (do_perm)
2149 {
2150 retval(0) = retval(0).permute (perm_vec, true);
2151 if (nargout > 1)
2152 retval(1) = retval(1).permute (perm_vec, true);
2153 }
2154
2155 // Process "linear" option
2156 if (linear && nargout > 1 && ! allflag && ! args(0).isempty ())
2157 retval(1) = get_linear_index (retval(1), args(0), vecdim, true);
2158
2159 return retval;
2160}
2161
2162DEFUN (cummin, args, nargout,
2163 doc: /* -*- texinfo -*-
2164@deftypefn {} {@var{M} =} cummin (@var{x})
2165@deftypefnx {} {@var{m} =} cummin (@var{x}, @var{dim})
2166@deftypefnx {} {@var{m} =} cummin (@var{x}, @var{vecdim})
2167@deftypefnx {} {@var{m} =} cummin (@var{x}, "all")
2168@deftypefnx {} {@var{m} =} cummin (@dots{}, @var{direction})
2169@deftypefnx {} {@var{m} =} cummin (@var{x}, @var{nanflag})
2170@deftypefnx {} {[@var{m}, @var{im}] =} cummin (@dots{})
2171@deftypefnx {} {[@var{m}, @var{im}] =} cummin (@dots{}, "linear")
2172@deftypefnx {} {@dots{} =} cummin (@dots{}, "ComparisonMethod", @var{method})
2173Return the cumulative minimum values from the array @var{x}.
2174
2175If @var{x} is a vector, then @code{cummin (@var{x})} returns a vector of the
2176same size with the cumulative minimum values of @var{x}.
2177
2178If @var{x} is a matrix, then @code{cummin (@var{x})} returns a matrix of the
2179same size with the cumulative minimum values along each column of @var{x}.
2180
2181If @var{x} is an array, then @code{cummin(@var{x})} returns an array of the
2182same size with the cumulative product along the first non-singleton dimension
2183of @var{x}.
2184
2185The class of output @var{y} is the same as the class of input @var{x}, unless
2186@var{x} is logical, in which case @var{y} is double.
2187
2188The optional input @var{dim} specifies the dimension to operate on and must be
2189a positive integer. Specifying any singleton dimension in @var{x}, including
2190any dimension exceeding @code{ndims (@var{x})}, will return @var{x}.
2191
2192Specifying multiple dimensions with input @var{vecdim}, a vector of
2193non-repeating dimensions, will operate along the array slice defined by
2194@var{vecdim}. If @var{vecdim} indexes all dimensions of @var{x}, then it is
2195equivalent to the option @qcode{"all"}. Any dimension in @var{vecdim} greater
2196than @code{ndims (@var{x})} is ignored.
2197
2198Specifying the dimension as @qcode{"all"} will cause @code{cummin} to operate
2199on all elements of @var{x}, and is equivalent to @code{cummin (@var{x}(:))}.
2200
2201The optional input @var{direction} specifies how the operating dimension is
2202traversed and can take the following values:
2203
2204@table @asis
2205@item @qcode{"forward"} (default)
2206
2207The cumulative minimum values are computed from beginning (index 1) to end
2208along the operating dimension.
2209
2210@item @qcode{"reverse"}
2211
2212The cumulative minimum values are computed from end to beginning along the
2213operating dimension.
2214@end table
2215
2216The optional variable @var{nanflag} specifies whether to include or exclude
2217@code{NaN} values from the calculation using any of the previously specified
2218input argument combinations. The default value for @var{nanflag} is
2219@qcode{"omitnan"} which ignores @code{NaN} values in the calculation. To
2220include @code{NaN} values, set the value of @var{nanflag} to
2221@qcode{"includenan"}, in which case @code{max} will return @code{NaN}, if any
2222element along the operating dimension is @code{NaN}.
2223
2224If called with two output arguments, @code{cummin} also returns the first index
2225of the minimum value(s) in @var{x}. Setting the @qcode{"linear"} flag returns
2226the linear index to the corresponding minimum values in @var{x}.
2227
2228The optional @qcode{"ComparisonMethod"} paired argument specifies the
2229comparison method for numeric input and it applies to both one input and two
2230input arrays. @var{method} can take any of the following values:
2231
2232@table @asis
2233@item @qcode{'auto'} (default)
2234Compare elements by @code{real (@var{x})} when @var{x} is real, and by
2235@code{abs (@var{x})} when @var{x} is complex.
2236
2237@item @qcode{'real'}
2238Compare elements by @code{real (@var{x})} when @var{x} is real or complex. For
2239elements with equal real parts, a second comparison by @code{imag (@var{x})} is
2240performed.
2241
2242@item @qcode{'abs'}
2243Compare elements by @code{abs (@var{x})} when @var{x} is real or complex. For
2244elements with equal magnitude, a second comparison by @code{angle (@var{x})} in
2245the interval from @math{-\pi} to @math{\pi} is performed.
2246@end table
2247@seealso{cummax, min, max}
2248@end deftypefn */)
2249{
2250 return do_cumminmax_body (args, nargout, true);
2251}
2252
2253/*
2254%!assert (cummin ([1, 4, 2, 3]), [1 1 1 1])
2255%!assert (cummin ([1; -10; 5; -2]), [1; -10; -10; -10])
2256%!assert (cummin ([4, i; -2, 2]), [4, i; -2, i])
2257%!assert (cummin ([1, 2; NaN, 1], 2), [1, 1; NaN, 1])
2258
2259%!test
2260%! x = reshape (1:8, [2,2,2]);
2261%! assert (cummin (x, 1), reshape ([1 1 3 3 5 5 7 7], [2,2,2]));
2262%! assert (cummin (x, 2), reshape ([1 2 1 2 5 6 5 6], [2,2,2]));
2263%! [w, iw] = cummin (x, 3);
2264%! assert (ndims (w), 3);
2265%! assert (w, repmat ([1 3; 2 4], [1 1 2]));
2266%! assert (ndims (iw), 3);
2267%! assert (iw, ones (2,2,2));
2268
2269## Test "includenan" and "reverse" options
2270%!assert (cummin ([1; -10; NaN; -2], 'includenan'), [1; -10; NaN; NaN])
2271%!assert (cummin ([1; -10; NaN; -2], 'reverse', 'includenan'), [NaN(3,1); -2])
2272%!assert (cummin ([1, -10, NaN, -2], 'includenan'), [1, -10, NaN, NaN])
2273%!assert (cummin ([1, -10, NaN, -2], 'reverse', 'includenan'), [NaN(1,3), -2])
2274%!shared A, C
2275%! A = [3, 5, NaN, 4, 2; 2, 6, 2, 9, 4; 1, 3, 0, NaN, 1; 5, 3, 4, 2, 0];
2276%! C = [3+2i, 5i, NaN, 4+i, 2i; 2, 6i, 2, 9+2i, 4i; ...
2277%! 1, 2+3i, 0, NaN, 1i; 5, 3i, 4, 2i, 0+i];
2278%!test
2279%! [m, i] = cummin (A, 'includenan');
2280%! m_exp = [3, 5, NaN, 4, 2; 2, 5, NaN, 4, 2; ...
2281%! 1, 3, NaN, NaN, 1; 1, 3, NaN, NaN, 0];
2282%! i_exp = [1, 1, 1, 1, 1; 2, 1, 1, 1, 1; 3, 3, 1, 3, 3; 3, 3, 1, 3, 4];
2283%! assert (m, m_exp);
2284%! assert (i, i_exp);
2285%!test
2286%! [m, i] = cummin (A, 'includenan', 'reverse');
2287%! m_exp = [1, 3, NaN, NaN, 0; 1, 3, 0, NaN, 0; ...
2288%! 1, 3, 0, NaN, 0; 5, 3, 4, 2, 0];
2289%! i_exp = [3, 4, 1, 3, 4; 3, 4, 3, 3, 4; 3, 4, 3, 3, 4; 4, 4, 4, 4, 4];
2290%! assert (m, m_exp);
2291%! assert (i, i_exp);
2292%!test
2293%! [m, i] = cummin (A, 'reverse');
2294%! m_exp = [1, 3, 0, 2, 0; 1, 3, 0, 2, 0; 1, 3, 0, 2, 0; 5, 3, 4, 2, 0];
2295%! i_exp = [3, 4, 3, 4, 4; 3, 4, 3, 4, 4; 3, 4, 3, 4, 4; 4, 4, 4, 4, 4];
2296%! assert (m, m_exp);
2297%! assert (i, i_exp);
2298%!test
2299%! [m, i] = cummin (A);
2300%! m_exp = [3, 5, NaN, 4, 2; 2, 5, 2, 4, 2; 1, 3, 0, 4, 1; 1, 3, 0, 2, 0];
2301%! i_exp = [1, 1, 1, 1, 1; 2, 1, 2, 1, 1; 3, 3, 3, 1, 3; 3, 3, 3, 4, 4];
2302%! assert (m, m_exp);
2303%! assert (i, i_exp);
2304%!test
2305%! [m, i] = cummin (A, 2, 'includenan');
2306%! m_exp = [3, 3, NaN, NaN, NaN; 2, 2, 2, 2, 2; 1, 1, 0, NaN, NaN; 5, 3, 3, 2, 0];
2307%! i_exp = [1, 1, 3, 3, 3; 1, 1, 1, 1, 1; 1, 1, 3, 4, 4; 1, 2, 2, 4, 5];
2308%! assert (m, m_exp);
2309%! assert (i, i_exp);
2310%!test
2311%! [m, i] = cummin (A, 2, 'includenan', 'reverse');
2312%! m_exp = [NaN, NaN, NaN, 2, 2; 2, 2, 2, 4, 4; ...
2313%! NaN, NaN, NaN, NaN, 1; 0, 0, 0, 0, 0];
2314%! i_exp = [3, 3, 3, 5, 5; 3, 3, 3, 5, 5; 4, 4, 4, 4, 5; 5, 5, 5, 5, 5];
2315%! assert (m, m_exp);
2316%! assert (i, i_exp);
2317%!test
2318%! [m, i] = cummin (A, 2, 'reverse');
2319%! m_exp = [2, 2, 2, 2, 2; 2, 2, 2, 4, 4; 0, 0, 0, 1, 1; 0, 0, 0, 0, 0];
2320%! i_exp = [5, 5, 5, 5, 5; 3, 3, 3, 5, 5; 3, 3, 3, 5, 5; 5, 5, 5, 5, 5];
2321%! assert (m, m_exp);
2322%! assert (i, i_exp);
2323%!test
2324%! [m, i] = cummin (A, 2);
2325%! m_exp = [3, 3, 3, 3, 2; 2, 2, 2, 2, 2; 1, 1, 0, 0, 0; 5, 3, 3, 2, 0];
2326%! i_exp = [1, 1, 1, 1, 5; 1, 1, 1, 1, 1; 1, 1, 3, 3, 3; 1, 2, 2, 4, 5];
2327%! assert (m, m_exp);
2328%! assert (i, i_exp);
2329
2330## Test single output against linear option with previous examples
2331%!test
2332%! [~, i_lin] = cummin (A, 'includenan', 'linear');
2333%! assert (cummin (A, 'includenan'), A(i_lin));
2334%!test
2335%! [~, i_lin] = cummin (A, 'includenan', 'reverse', 'linear');
2336%! assert (cummin (A, 'includenan', 'reverse'), A(i_lin));
2337%!test
2338%! [~, i_lin] = cummin (A, 'reverse', 'linear');
2339%! assert (cummin (A, 'reverse'), A(i_lin));
2340%!test
2341%! [~, i_lin] = cummin (A, 'linear');
2342%! assert (cummin (A), A(i_lin));
2343%!test
2344%! [~, i_lin] = cummin (A, 2, 'includenan', 'linear');
2345%! assert (cummin (A, 2, 'includenan'), A(i_lin));
2346%!test
2347%! [~, i_lin] = cummin (A, 2, 'includenan', 'reverse', 'linear');
2348%! assert (cummin (A, 2, 'includenan', 'reverse'), A(i_lin));
2349%!test
2350%! [~, i_lin] = cummin (A, 2, 'reverse', 'linear');
2351%! assert (cummin (A, 2, 'reverse'), A(i_lin));
2352%!test
2353%! [~, i_lin] = cummin (A, 2, 'linear');
2354%! assert (cummin (A, 2), A(i_lin));
2355
2356## Test "includenan" and "reverse" options with complex input
2357%!test
2358%! [m, idx] = cummin (C, 'includenan');
2359%! m_exp = [3+2i, 5i, NaN, 4+i, 2i; 2, 5i, NaN, 4+i, 2i; ...
2360%! 1, 2+3i, NaN, NaN, 1i; 1, 3i, NaN, NaN, 1i];
2361%! i_exp = [1, 1, 1, 1, 1; 2, 1, 1, 1, 1; 3, 3, 1, 3, 3; 3, 4, 1, 3, 3];
2362%! assert (m, m_exp);
2363%! assert (idx, i_exp);
2364%!test
2365%! [m, idx] = cummin (C, 'includenan', 'reverse');
2366%! m_exp = [1, 3i, NaN, NaN, 1i; 1, 3i, 0, NaN, 1i; ...
2367%! 1, 3i, 0, NaN, 1i; 5, 3i, 4, 2i, 1i];
2368%! i_exp = [3, 4, 1, 3, 4; 3, 4, 3, 3, 4; 3, 4, 3, 3, 4; 4, 4, 4, 4, 4];
2369%! assert (m, m_exp);
2370%! assert (idx, i_exp);
2371%!test
2372%! [m, idx] = cummin (C, 'reverse');
2373%! m_exp = [1, 3i, 0, 2i, 1i; 1, 3i, 0, 2i, 1i; ...
2374%! 1, 3i, 0, 2i, 1i; 5, 3i, 4, 2i, 1i];
2375%! i_exp = [3, 4, 3, 4, 4; 3, 4, 3, 4, 4; 3, 4, 3, 4, 4; 4, 4, 4, 4, 4];
2376%! assert (m, m_exp);
2377%! assert (idx, i_exp);
2378%!test
2379%! [m, idx] = cummin (C);
2380%! m_exp = [3+2i, 5i, NaN, 4+i, 2i; 2, 5i, 2, 4+i, 2i; ...
2381%! 1, 2+3i, 0, 4+i, 1i; 1, 3i, 0, 2i, 1i];
2382%! i_exp = [1, 1, 1, 1, 1; 2, 1, 2, 1, 1; 3, 3, 3, 1, 3; 3, 4, 3, 4, 3];
2383%! assert (m, m_exp);
2384%! assert (idx, i_exp);
2385%!test
2386%! [m, idx] = cummin (C, 2, 'includenan');
2387%! m_exp = [3+2i, 3+2i, NaN, NaN, NaN; 2, 2, 2, 2, 2; ...
2388%! 1, 1, 0, NaN, NaN; 5, 3i, 3i, 2i, 1i];
2389%! i_exp = [1, 1, 3, 3, 3; 1, 1, 1, 1, 1; 1, 1, 3, 4, 4; 1, 2, 2, 4, 5];
2390%! assert (m, m_exp);
2391%! assert (idx, i_exp);
2392%!test
2393%! [m, idx] = cummin (C, 2, 'includenan', 'reverse');
2394%! m_exp = [NaN, NaN, NaN, 2i, 2i; 2, 2, 2, 4i, 4i; ...
2395%! NaN, NaN, NaN, NaN, 1i; 1i, 1i, 1i, 1i, 1i];
2396%! i_exp = [3, 3, 3, 5, 5; 3, 3, 3, 5, 5; 4, 4, 4, 4, 5; 5, 5, 5, 5, 5];
2397%! assert (m, m_exp);
2398%! assert (idx, i_exp);
2399%!test
2400%! [m, idx] = cummin (C, 2, 'reverse');
2401%! m_exp = [2i, 2i, 2i, 2i, 2i; 2, 2, 2, 4i, 4i; ...
2402%! 0, 0, 0, 1i, 1i; 1i, 1i, 1i, 1i, 1i];
2403%! i_exp = [5, 5, 5, 5, 5; 3, 3, 3, 5, 5; 3, 3, 3, 5, 5; 5, 5, 5, 5, 5];
2404%! assert (m, m_exp);
2405%! assert (idx, i_exp);
2406%!test
2407%! [m, idx] = cummin (C, 2);
2408%! m_exp = [3+2i, 3+2i, 3+2i, 3+2i, 2i; 2, 2, 2, 2, 2; ...
2409%! 1, 1, 0, 0, 0; 5, 3i, 3i, 2i, 1i];
2410%! i_exp = [1, 1, 1, 1, 5; 1, 1, 1, 1, 1; 1, 1, 3, 3, 3; 1, 2, 2, 4, 5];
2411%! assert (m, m_exp);
2412%! assert (idx, i_exp);
2413
2414## Test 'ComparisonMethod' with complex input
2415%!test
2416%! [m, idx] = cummin (C, 2, 'ComparisonMethod', 'real');
2417%! m_exp = [3+2i, 5i, 5i, 5i, 2i; 2, 6i, 6i, 6i, 4i; ...
2418%! 1, 1, 0, 0, 0; 5, 3i, 3i, 2i, 1i];
2419%! i_exp = [1, 2, 2, 2, 5; 1, 2, 2, 2, 5; 1, 1, 3, 3, 3; 1, 2, 2, 4, 5];
2420%! assert (m, m_exp);
2421%! assert (idx, i_exp);
2422%!test
2423%! [m, idx] = cummin (C, 'includenan', 'ComparisonMethod', 'real');
2424%! m_exp = [3+2i, 5i, NaN, 4+i, 2i; 2, 5i, NaN, 4+i, 2i; ...
2425%! 1, 5i, NaN, NaN, 1i; 1, 3i, NaN, NaN, 1i];
2426%! i_exp = [1, 1, 1, 1, 1; 2, 1, 1, 1, 1; 3, 1, 1, 3, 3; 3, 4, 1, 3, 3];
2427%! assert (m, m_exp);
2428%! assert (idx, i_exp);
2429
2430## Test "linear" option with ND array
2431%!shared x
2432%! x = randi ([-10, 10], 3, 4, 5, 2);
2433%!test
2434%! [m, i] = cummin (x, [2, 3], 'linear');
2435%! assert (m, x(i));
2436%!test
2437%! [m, i] = cummin (x, [1, 3], 'linear');
2438%! assert (m, x(i));
2439%!test
2440%! [m, i] = cummin (x, [2, 4], 'linear');
2441%! assert (m, x(i));
2442%!test
2443%! [m, i] = cummin (x, [1, 4], 'linear');
2444%! assert (m, x(i));
2445%!test
2446%! [m, i] = cummin (x, [1, 2, 3], 'linear');
2447%! assert (m, x(i));
2448%!test
2449%! [m, i] = cummin (x, [2, 3, 4], 'linear');
2450%! assert (m, x(i));
2451%!test
2452%! [m, i] = cummin ([1, 2, 3; 4, 5, 6], 1, 'linear');
2453%! assert (m, [1, 2, 3; 1, 2, 3]);
2454%! assert (i, [1, 3, 5; 1, 3, 5]);
2455%!test
2456%! [m, i] = cummin ([1, 2, 3; 4, 5, 6], 2, 'linear');
2457%! assert (m, [1, 1, 1; 4, 4, 4]);
2458%! assert (i, [1, 1, 1; 2, 2, 2]);
2459%!test
2460%! [m, i] = cummin ([1, 2, 3; 4, 5, 6], 1, 'linear', 'reverse');
2461%! assert (m, [1, 2, 3; 4, 5, 6]);
2462%! assert (i, [1, 3, 5; 2, 4, 6]);
2463%!test
2464%! [m, i] = cummin ([1, 2, 3; 4, 5, 6], 2, 'linear', 'reverse');
2465%! assert (m, [1, 2, 3; 4, 5, 6]);
2466%! assert (i, [1, 3, 5; 2, 4, 6]);
2467
2468%!error <Invalid call> cummin ()
2469%!error <Invalid call> cummin (1, 2, 3)
2470%!error <unrecognized optional argument 'foobar'> cummin (1, "foobar")
2471%!error <cannot set DIM or VECDIM with 'all' flag>
2472%! cummin (ones (3,3), 1, "all")
2473%!error <cannot set DIM or VECDIM with 'all' flag>
2474%! cummin (ones (3,3), [1, 2], "all")
2475%!error <invalid dimension DIM = 0> cummin (ones (3,3), 0)
2476%!error <invalid dimension DIM = -1> cummin (ones (3,3), -1)
2477%!error <invalid dimension in VECDIM = -2> cummin (ones (3,3), [1 -2])
2478%!error <duplicate dimension in VECDIM = 2> cummin (ones (3,3), [1 2 2])
2479%!error <duplicate dimension in VECDIM = 1> cummin (ones (3,3), [1 1 2])
2480%!error <wrong type argument 'cell'> cummin ({1 2 3 4})
2481*/
2482
2483DEFUN (cummax, args, nargout,
2484 doc: /* -*- texinfo -*-
2485@deftypefn {} {@var{M} =} cummax (@var{x})
2486@deftypefnx {} {@var{m} =} cummax (@var{x}, @var{dim})
2487@deftypefnx {} {@var{m} =} cummax (@var{x}, @var{vecdim})
2488@deftypefnx {} {@var{m} =} cummax (@var{x}, "all")
2489@deftypefnx {} {@var{m} =} cummax (@dots{}, @var{direction})
2490@deftypefnx {} {@var{m} =} cummax (@var{x}, @var{nanflag})
2491@deftypefnx {} {[@var{m}, @var{im}] =} cummax (@dots{})
2492@deftypefnx {} {[@var{m}, @var{im}] =} cummax (@dots{}, "linear")
2493@deftypefnx {} {@dots{} =} cummax (@dots{}, "ComparisonMethod", @var{method})
2494Return the cumulative maximum values from the array @var{x}.
2495
2496If @var{x} is a vector, then @code{cummax (@var{x})} returns a vector of the
2497same size with the cumulative maximum values of @var{x}.
2498
2499If @var{x} is a matrix, then @code{cummax (@var{x})} returns a matrix of the
2500same size with the cumulative maximum values along each column of @var{x}.
2501
2502If @var{x} is an array, then @code{cummax(@var{x})} returns an array of the
2503same size with the cumulative product along the first non-singleton dimension
2504of @var{x}.
2505
2506The class of output @var{y} is the same as the class of input @var{x}, unless
2507@var{x} is logical, in which case @var{y} is double.
2508
2509The optional input @var{dim} specifies the dimension to operate on and must be
2510a positive integer. Specifying any singleton dimension in @var{x}, including
2511any dimension exceeding @code{ndims (@var{x})}, will return @var{x}.
2512
2513Specifying multiple dimensions with input @var{vecdim}, a vector of
2514non-repeating dimensions, will operate along the array slice defined by
2515@var{vecdim}. If @var{vecdim} indexes all dimensions of @var{x}, then it is
2516equivalent to the option @qcode{"all"}. Any dimension in @var{vecdim} greater
2517than @code{ndims (@var{x})} is ignored.
2518
2519Specifying the dimension as @qcode{"all"} will cause @code{cummax} to operate
2520on all elements of @var{x}, and is equivalent to @code{cummax (@var{x}(:))}.
2521
2522The optional input @var{direction} specifies how the operating dimension is
2523traversed and can take the following values:
2524
2525@table @asis
2526@item @qcode{"forward"} (default)
2527
2528The cumulative maximum values are computed from beginning (index 1) to end
2529along the operating dimension.
2530
2531@item @qcode{"reverse"}
2532
2533The cumulative maximum values are computed from end to beginning along the
2534operating dimension.
2535@end table
2536
2537The optional variable @var{nanflag} specifies whether to include or exclude
2538@code{NaN} values from the calculation using any of the previously specified
2539input argument combinations. The default value for @var{nanflag} is
2540@qcode{"omitnan"} which ignores @code{NaN} values in the calculation. To
2541include @code{NaN} values, set the value of @var{nanflag} to
2542@qcode{"includenan"}, in which case @code{max} will return @code{NaN}, if any
2543element along the operating dimension is @code{NaN}.
2544
2545If called with two output arguments, @code{cummax} also returns the first index
2546of the maximum value(s) in @var{x}. Setting the @qcode{"linear"} flag returns
2547the linear index to the corresponding minimum values in @var{x}.
2548
2549The optional @qcode{"ComparisonMethod"} paired argument specifies the
2550comparison method for numeric input and it applies to both one input and two
2551input arrays. @var{method} can take any of the following values:
2552
2553@table @asis
2554@item @qcode{'auto'} (default)
2555Compare elements by @code{real (@var{x})} when @var{x} is real, and by
2556@code{abs (@var{x})} when @var{x} is complex.
2557
2558@item @qcode{'real'}
2559Compare elements by @code{real (@var{x})} when @var{x} is real or complex. For
2560elements with equal real parts, a second comparison by @code{imag (@var{x})} is
2561performed.
2562
2563@item @qcode{'abs'}
2564Compare elements by @code{abs (@var{x})} when @var{x} is real or complex. For
2565elements with equal magnitude, a second comparison by @code{angle (@var{x})}
2566in the interval from @math{-\pi} to @math{\pi} is performed.
2567@end table
2568@seealso{cummin, max, min}
2569@end deftypefn */)
2570{
2571 return do_cumminmax_body (args, nargout, false);
2572}
2573
2574/*
2575%!assert (cummax ([1, 4, 2, 3]), [1, 4, 4, 4])
2576%!assert (cummax ([1; -10; 5; -2]), [1; 1; 5; 5])
2577%!assert (cummax ([4, i, 4.9, -2, 2, 3+4i]), [4, 4, 4.9, 4.9, 4.9, 3+4i])
2578%!assert (cummax ([1, NaN, 0; NaN, NaN, 1], 2), [1, 1, 1; NaN, NaN, 1])
2579
2580%!test
2581%! x = reshape (8:-1:1, [2,2,2]);
2582%! assert (cummax (x, 1), reshape ([8 8 6 6 4 4 2 2], [2,2,2]));
2583%! assert (cummax (x, 2), reshape ([8 7 8 7 4 3 4 3], [2,2,2]));
2584%! [w, iw] = cummax (x, 3);
2585%! assert (ndims (w), 3);
2586%! assert (w, repmat ([8 6; 7 5], [1 1 2]));
2587%! assert (ndims (iw), 3);
2588%! assert (iw, ones (2,2,2));
2589
2590## Test "includenan" and "reverse" options
2591%!assert (cummax ([1; -10; NaN; -2], 'includenan'), [1; 1; NaN; NaN])
2592%!assert (cummax ([1; -10; NaN; -2], 'reverse', 'includenan'), [NaN(3,1); -2])
2593%!assert (cummax ([1, -10, NaN, -2], 'includenan'), [1, 1, NaN, NaN])
2594%!assert (cummax ([1, -10, NaN, -2], 'reverse', 'includenan'), [NaN(1,3), -2])
2595%!shared A, C
2596%! A = [3, 5, NaN, 4, 2; 2, 6, 2, 9, 4; 1, 3, 0, NaN, 1; 5, 3, 4, 2, 0];
2597%! C = [3+2i, 5i, NaN, 4+i, 2i; 2, 6i, 2, 9+2i, 4i; ...
2598%! 1, 2+3i, 0, NaN, 1i; 5, 3i, 4, 2i, 0+i];
2599%!test
2600%! [m, i] = cummax (A, 'includenan');
2601%! m_exp = [3, 5, NaN, 4, 2; 3, 6, NaN, 9, 4; ...
2602%! 3, 6, NaN, NaN, 4; 5, 6, NaN, NaN, 4];
2603%! i_exp = [1, 1, 1, 1, 1; 1, 2, 1, 2, 2; 1, 2, 1, 3, 2; 4, 2, 1, 3, 2];
2604%! assert (m, m_exp);
2605%! assert (i, i_exp);
2606%!test
2607%! [m, i] = cummax (A, 'includenan', 'reverse');
2608%! m_exp = [5, 6, NaN, NaN, 4; 5, 6, 4, NaN, 4; ...
2609%! 5, 3, 4, NaN, 1; 5, 3, 4, 2, 0];
2610%! i_exp = [4, 2, 1, 3, 2; 4, 2, 4, 3, 2; 4, 4, 4, 3, 3; 4, 4, 4, 4, 4];
2611%! assert (m, m_exp);
2612%! assert (i, i_exp);
2613%!test
2614%! [m, i] = cummax (A, 'reverse');
2615%! m_exp = [5, 6, 4, 9, 4; 5, 6, 4, 9, 4; 5, 3, 4, 2, 1; 5, 3, 4, 2, 0];
2616%! i_exp = [4, 2, 4, 2, 2; 4, 2, 4, 2, 2; 4, 4, 4, 4, 3; 4, 4, 4, 4, 4];
2617%! assert (m, m_exp);
2618%! assert (i, i_exp);
2619%!test
2620%! [m, i] = cummax (A);
2621%! m_exp = [3, 5, NaN, 4, 2; 3, 6, 2, 9, 4; 3, 6, 2, 9, 4; 5, 6, 4, 9, 4];
2622%! i_exp = [1, 1, 1, 1, 1; 1, 2, 2, 2, 2; 1, 2, 2, 2, 2; 4, 2, 4, 2, 2];
2623%! assert (m, m_exp);
2624%! assert (i, i_exp);
2625%!test
2626%! [m, i] = cummax (A, 2, 'includenan');
2627%! m_exp = [3, 5, NaN, NaN, NaN; 2, 6, 6, 9, 9; 1, 3, 3, NaN, NaN; 5, 5, 5, 5, 5];
2628%! i_exp = [1, 2, 3, 3, 3; 1, 2, 2, 4, 4; 1, 2, 2, 4, 4; 1, 1, 1, 1, 1];
2629%! assert (m, m_exp);
2630%! assert (i, i_exp);
2631%!test
2632%! [m, i] = cummax (A, 2, 'includenan', 'reverse');
2633%! m_exp = [NaN, NaN, NaN, 4, 2; 9, 9, 9, 9, 4; ...
2634%! NaN, NaN, NaN, NaN, 1; 5, 4, 4, 2, 0];
2635%! i_exp = [3, 3, 3, 4, 5; 4, 4, 4, 4, 5; 4, 4, 4, 4, 5; 1, 3, 3, 4, 5];
2636%! assert (m, m_exp);
2637%! assert (i, i_exp);
2638%!test
2639%! [m, i] = cummax (A, 2, 'reverse');
2640%! m_exp = [5, 5, 4, 4, 2; 9, 9, 9, 9, 4; 3, 3, 1, 1, 1; 5, 4, 4, 2, 0];
2641%! i_exp = [2, 2, 4, 4, 5; 4, 4, 4, 4, 5; 2, 2, 5, 5, 5; 1, 3, 3, 4, 5];
2642%! assert (m, m_exp);
2643%! assert (i, i_exp);
2644%!test
2645%! [m, i] = cummax (A, 2);
2646%! m_exp = [3, 5, 5, 5, 5; 2, 6, 6, 9, 9; 1, 3, 3, 3, 3; 5, 5, 5, 5, 5];
2647%! i_exp = [1, 2, 2, 2, 2; 1, 2, 2, 4, 4; 1, 2, 2, 2, 2; 1, 1, 1, 1, 1];
2648%! assert (m, m_exp);
2649%! assert (i, i_exp);
2650
2651## Test single output against linear option with previous examples
2652%!test
2653%! [~, i_lin] = cummax (A, 'includenan', 'linear');
2654%! assert (cummax (A, 'includenan'), A(i_lin));
2655%!test
2656%! [~, i_lin] = cummax (A, 'includenan', 'reverse', 'linear');
2657%! assert (cummax (A, 'includenan', 'reverse'), A(i_lin));
2658%!test
2659%! [~, i_lin] = cummax (A, 'reverse', 'linear');
2660%! assert (cummax (A, 'reverse'), A(i_lin));
2661%!test
2662%! [~, i_lin] = cummax (A, 'linear');
2663%! assert (cummax (A), A(i_lin));
2664%!test
2665%! [~, i_lin] = cummax (A, 2, 'includenan', 'linear');
2666%! assert (cummax (A, 2, 'includenan'), A(i_lin));
2667%!test
2668%! [~, i_lin] = cummax (A, 2, 'includenan', 'reverse', 'linear');
2669%! assert (cummax (A, 2, 'includenan', 'reverse'), A(i_lin));
2670%!test
2671%! [~, i_lin] = cummax (A, 2, 'reverse', 'linear');
2672%! assert (cummax (A, 2, 'reverse'), A(i_lin));
2673%!test
2674%! [~, i_lin] = cummax (A, 2, 'linear');
2675%! assert (cummax (A, 2), A(i_lin));
2676
2677## Test "includenan" and "reverse" options with complex input
2678%!test
2679%! [m, idx] = cummax (C, 'includenan');
2680%! m_exp = [3+2i, 5i, NaN, 4+i, 2i; 3+2i, 6i, NaN, 9+2i, 4i; ...
2681%! 3+2i, 6i, NaN, NaN, 4i; 5, 6i, NaN, NaN, 4i];
2682%! i_exp = [1, 1, 1, 1, 1; 1, 2, 1, 2, 2; 1, 2, 1, 3, 2; 4, 2, 1, 3, 2];
2683%! assert (m, m_exp);
2684%! assert (idx, i_exp);
2685%!test
2686%! [m, idx] = cummax (C, 'includenan', 'reverse');
2687%! m_exp = [5, 6i, NaN, NaN, 4i; 5, 6i, 4, NaN, 4i; ...
2688%! 5, 2+3i, 4, NaN, 1i; 5, 3i, 4, 2i, 1i];
2689%! i_exp = [4, 2, 1, 3, 2; 4, 2, 4, 3, 2; 4, 3, 4, 3, 4; 4, 4, 4, 4, 4];
2690%! assert (m, m_exp);
2691%! assert (idx, i_exp);
2692%!test
2693%! [m, idx] = cummax (C, 'reverse');
2694%! m_exp = [5, 6i, 4, 9+2i, 4i; 5, 6i, 4, 9+2i, 4i; ...
2695%! 5, 2+3i, 4, 2i, 1i; 5, 3i, 4, 2i, 1i];
2696%! i_exp = [4, 2, 4, 2, 2; 4, 2, 4, 2, 2; 4, 3, 4, 4, 4; 4, 4, 4, 4, 4];
2697%! assert (m, m_exp);
2698%! assert (idx, i_exp);
2699%!test
2700%! [m, idx] = cummax (C);
2701%! m_exp = [3+2i, 5i, NaN, 4+i, 2i; 3+2i, 6i, 2, 9+2i, 4i; ...
2702%! 3+2i, 6i, 2, 9+2i, 4i; 5, 6i, 4, 9+2i, 4i];
2703%! i_exp = [1, 1, 1, 1, 1; 1, 2, 2, 2, 2; 1, 2, 2, 2, 2; 4, 2, 4, 2, 2];
2704%! assert (m, m_exp);
2705%! assert (idx, i_exp);
2706%!test
2707%! [m, idx] = cummax (C, 2, 'includenan');
2708%! m_exp = [3+2i, 5i, NaN, NaN, NaN; 2, 6i, 6i, 9+2i, 9+2i; ...
2709%! 1, 2+3i, 2+3i, NaN, NaN; 5, 5, 5, 5, 5];
2710%! i_exp = [1, 2, 3, 3, 3; 1, 2, 2, 4, 4; 1, 2, 2, 4, 4; 1, 1, 1, 1, 1];
2711%! assert (m, m_exp);
2712%! assert (idx, i_exp);
2713%!test
2714%! [m, idx] = cummax (C, 2, 'includenan', 'reverse');
2715%! m_exp = [NaN, NaN, NaN, 4+i, 2i; 9+2i, 9+2i, 9+2i, 9+2i, 4i; ...
2716%! NaN, NaN, NaN, NaN, 1i; 5, 4, 4, 2i, 1i];
2717%! i_exp = [3, 3, 3, 4, 5; 4, 4, 4, 4, 5; 4, 4, 4, 4, 5; 1, 3, 3, 4, 5];
2718%! assert (m, m_exp);
2719%! assert (idx, i_exp);
2720%!test
2721%! [m, idx] = cummax (C, 2, 'reverse');
2722%! m_exp = [5i, 5i, 4+i, 4+i, 2i; 9+2i, 9+2i, 9+2i, 9+2i, 4i; ...
2723%! 2+3i, 2+3i, 1i, 1i, 1i; 5, 4, 4, 2i, 1i];
2724%! i_exp = [2, 2, 4, 4, 5; 4, 4, 4, 4, 5; 2, 2, 5, 5, 5; 1, 3, 3, 4, 5];
2725%! assert (m, m_exp);
2726%! assert (idx, i_exp);
2727%!test
2728%! [m, idx] = cummax (C, 2);
2729%! m_exp = [3+2i, 5i, 5i, 5i, 5i; 2, 6i, 6i, 9+2i, 9+2i; ...
2730%! 1, 2+3i, 2+3i, 2+3i, 2+3i; 5, 5, 5, 5, 5];
2731%! i_exp = [1, 2, 2, 2, 2; 1, 2, 2, 4, 4; 1, 2, 2, 2, 2; 1, 1, 1, 1, 1];
2732%! assert (m, m_exp);
2733%! assert (idx, i_exp);
2734
2735## Test 'ComparisonMethod' with complex input
2736%!test
2737%! [m, idx] = cummax (C, 2, 'ComparisonMethod', 'real');
2738%! m_exp = [3+2i, 3+2i, 3+2i, 4+i, 4+i; 2, 2, 2, 9+2i, 9+2i; ...
2739%! 1, 2+3i, 2+3i, 2+3i, 2+3i; 5, 5, 5, 5, 5];
2740%! i_exp = [1, 1, 1, 4, 4; 1, 1, 1, 4, 4; 1, 2, 2, 2, 2; 1, 1, 1, 1, 1];
2741%! assert (m, m_exp);
2742%! assert (idx, i_exp);
2743%!test
2744%! [m, idx] = cummax (C, 'includenan', 'reverse', 'ComparisonMethod', 'real');
2745%! m_exp = [5, 2+3i, NaN, NaN, 4i; 5, 2+3i, 4, NaN, 4i; ...
2746%! 5, 2+3i, 4, NaN, 1i; 5, 3i, 4, 2i, 1i];
2747%! i_exp = [4, 3, 1, 3, 2; 4, 3, 4, 3, 2; 4, 3, 4, 3, 4; 4, 4, 4, 4, 4];
2748%! assert (m, m_exp);
2749%! assert (idx, i_exp);
2750
2751## Test "linear" option with ND array
2752%!shared x
2753%! x = randi ([-10, 10], 3, 4, 5, 2);
2754%!test
2755%! [m, i] = cummax (x, [2, 3], 'linear');
2756%! assert (m, x(i));
2757%!test
2758%! [m, i] = cummax (x, [1, 3], 'linear');
2759%! assert (m, x(i));
2760%!test
2761%! [m, i] = cummax (x, [2, 4], 'linear');
2762%! assert (m, x(i));
2763%!test
2764%! [m, i] = cummax (x, [1, 4], 'linear');
2765%! assert (m, x(i));
2766%!test
2767%! [m, i] = cummax (x, [1, 2, 3], 'linear');
2768%! assert (m, x(i));
2769%!test
2770%! [m, i] = cummax (x, [2, 3, 4], 'linear');
2771%! assert (m, x(i));
2772%!test
2773%! [m, i] = cummax ([1, 2, 3; 4, 5, 6], 1, 'linear');
2774%! assert (m, [1, 2, 3; 4, 5, 6]);
2775%! assert (i, [1, 3, 5; 2, 4, 6]);
2776%!test
2777%! [m, i] = cummax ([1, 2, 3; 4, 5, 6], 2, 'linear');
2778%! assert (m, [1, 2, 3; 4, 5, 6]);
2779%! assert (i, [1, 3, 5; 2, 4, 6]);
2780%!test
2781%! [m, i] = cummax ([1, 2, 3; 4, 5, 6], 1, 'linear', 'reverse');
2782%! assert (m, [4, 5, 6; 4, 5, 6]);
2783%! assert (i, [2, 4, 6; 2, 4, 6]);
2784%!test
2785%! [m, i] = cummax ([1, 2, 3; 4, 5, 6], 2, 'linear', 'reverse');
2786%! assert (m, [3, 3, 3; 6, 6, 6]);
2787%! assert (i, [5, 5, 5; 6, 6, 6]);
2788
2789%!error <Invalid call> cummax ()
2790%!error <Invalid call> cummax (1, 2, 3)
2791%!error <unrecognized optional argument 'foobar'> cummax (1, "foobar")
2792%!error <cannot set DIM or VECDIM with 'all' flag>
2793%! cummax (ones (3,3), 1, "all")
2794%!error <cannot set DIM or VECDIM with 'all' flag>
2795%! cummax (ones (3,3), [1, 2], "all")
2796%!error <invalid dimension DIM = 0> cummax (ones (3,3), 0)
2797%!error <invalid dimension DIM = -1> cummax (ones (3,3), -1)
2798%!error <invalid dimension in VECDIM = -2> cummax (ones (3,3), [1 -2])
2799%!error <duplicate dimension in VECDIM = 2> cummax (ones (3,3), [1 2 2])
2800%!error <duplicate dimension in VECDIM = 1> cummax (ones (3,3), [1 1 2])
2801%!error <wrong type argument 'cell'> cummax ({1 2 3 4})
2802*/
2803
2804OCTAVE_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
N Dimensional Array with copy-on-write semantics.
Definition Array-base.h:130
void resize(const dim_vector &dv, const T &rfv)
Size of the specified dimension.
bool isempty() const
Size of the specified dimension.
Definition Array-base.h:674
octave_idx_type numel() const
Number of elements in the array.
Definition Array-base.h:440
boolNDArray any(int dim=-1) const
boolNDArray all(int dim=-1) const
Vector representing the dimensions (size) of an Array.
Definition dim-vector.h:92
void resize(int n, int fill_value=0)
Definition dim-vector.h:278
int first_non_singleton(int def=0) const
Definition dim-vector.h:481
void resize(octave_idx_type n, const octave_value &rfv=octave_value())
Definition ovl.h:115
Array< octave_value > array_value() const
Definition ovl.h:88
octave_idx_type length() const
Definition ovl.h:111
boolNDArray bool_array_value(bool warn=false) const
Definition ov.h:898
int int_value(bool req_int=false, bool frc_str_conv=false) const
Definition ov.h:810
int ndims() const
Definition ov.h:549
bool is_scalar_type() const
Definition ov.h:742
octave_value permute(const Array< int > &vec, bool inv=false) const
Definition ov.h:572
bool is_range() const
Definition ov.h:644
bool issparse() const
Definition ov.h:751
octave_value reshape(const dim_vector &dv) const
Definition ov.h:569
octave_idx_type numel() const
Definition ov.h:557
Array< int > int_vector_value(bool req_int=false, bool frc_str_conv=false, bool frc_vec_conv=false) const
std::string type_name() const
Definition ov.h:1358
builtin_type_t builtin_type() const
Definition ov.h:688
octave::range< double > range_value() const
Definition ov.h:992
dim_vector dims() const
Definition ov.h:539
OCTAVE_BEGIN_NAMESPACE(octave) static octave_value daspk_fcn
void print_usage()
Definition defun-int.h:72
#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:1083
void error(const char *fmt,...)
Definition error.cc:1008
void err_wrong_type_arg(const char *name, const char *s)
Definition errwarn.cc:166
octave_value do_minmax_bin_op< charNDArray >(const octave_value &argx, const octave_value &argy, bool ismin)
Definition max.cc:609
#define MAKE_INT_BRANCH(X)
octave_value_list do_minmax_red_op< boolNDArray >(const octave_value &arg, int nargout, int dim, bool ismin)
Definition max.cc:431
octave_value_list do_minmax_red_op< charNDArray >(const octave_value &arg, int nargout, int dim, bool ismin)
Definition max.cc:401
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:67
builtin_type_t
Definition ov-base.h:83
@ btyp_float_complex
Definition ov-base.h:87
@ btyp_double
Definition ov-base.h:84
@ btyp_float
Definition ov-base.h:85
@ btyp_bool
Definition ov-base.h:96
@ btyp_char
Definition ov-base.h:97
@ btyp_complex
Definition ov-base.h:86
charNDArray octave_value_extract< charNDArray >(const octave_value &v)
Definition ov.h:1753
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
F77_RET_T const F77_DBLE * x