File: quantile.m

package info (click to toggle)
octave 11.0.92-1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 148,624 kB
  • sloc: cpp: 347,499; ansic: 85,112; fortran: 20,693; objc: 10,276; sh: 8,747; lex: 4,496; yacc: 4,406; perl: 1,544; java: 1,365; awk: 1,282; makefile: 666; xml: 192
file content (766 lines) | stat: -rw-r--r-- 26,442 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
########################################################################
##
## Copyright (C) 2008-2026 The Octave Project Developers
##
## See the file COPYRIGHT.md in the top-level directory of this
## distribution or <https://octave.org/copyright/>.
##
## This file is part of Octave.
##
## Octave is free software: you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation, either version 3 of the License, or
## (at your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
## GNU General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <https://www.gnu.org/licenses/>.
##
########################################################################

## -*- texinfo -*-
## @deftypefn  {} {@var{q} =} quantile (@var{x})
## @deftypefnx {} {@var{q} =} quantile (@var{x}, @var{p})
## @deftypefnx {} {@var{q} =} quantile (@var{x}, @var{n})
## @deftypefnx {} {@var{q} =} quantile (@var{x}, @dots{}, @var{dim})
## @deftypefnx {} {@var{q} =} quantile (@var{x}, @dots{}, @var{vecdim})
## @deftypefnx {} {@var{q} =} quantile (@var{x}, @dots{}, "all")
## @deftypefnx {} {@var{q} =} quantile (@var{x}, @var{p}, @dots{}, @var{method})
## @deftypefnx {} {@var{q} =} quantile (@var{x}, @var{n}, @dots{}, @var{method})
## Compute the quantiles of the input data @var{x}.
##
## If @var{x} is a vector, then @code{quantile (@var{x})} computes the quantiles
## specified by @var{p} of the data in @var{x}.
##
## If @var{x} is a matrix, then @code{quantile (@var{x})} returns a matrix such
## that the i-th row of @var{q} contains the @var{p}(i)th quantiles of each
## column of @var{x}.
##
## If @var{x} is an array, then @code{quantile (@var{x})} computes the quantiles
## specified by @var{p} along the first non-singleton dimension of @var{x}.
##
## The data in @var{x} must be numeric and any NaN values are ignored.  The
## size of @var{q} is equal to the size of @var{x} except for the operating
## dimension, which equals to the number of quantiles specified by @var{p}
## or @var{n}.
##
## @var{p} is a numeric vector specifying the percentiles to be computed, which
## correspond to the cumulative probabilities of the data .  All elements of
## @var{p} must be in the range from 0 to 1.  If @var{p} is unspecified, return
## the percentiles for @code{[0.00 0.25 0.50 0.75 1.00]}.  Alternatively, the
## second input argument may be specified as a positive integer value @var{n},
## in which case @code{quantile} returns the quantiles for @var{n} evenly
## spaced cumulative probabilities computed as (1/(@var{n} + 1), 2/(@var{n}
## + 1), @dots{}, @var{n}/(@var{n} + 1)) for @code{@var{n} > 1}.
##
## The optional input @var{dim} specifies the dimension to operate on and must
## be a positive integer.  Specifying any singleton dimension of @var{x},
## including any dimension exceeding @code{ndims (@var{x})}, will return N
## copies of @var{x} along the operating dimension, where N is the number of
## specified quantiles.
##
## Specifying multiple dimensions with input @var{vecdim}, a vector of
## non-repeating dimensions, will operate along the array slice defined by
## @var{vecdim}.  If @var{vecdim} indexes all dimensions of @var{x}, then it is
## equivalent to the option @qcode{"all"}.  Any dimension in @var{vecdim}
## greater than @code{ndims (@var{x})} is ignored.  If all dimensions in
## @var{vecdim} are greater than @code{ndims (@var{x})}, then @code{quantile}
## will return N copies of @var{x} along the smallest dimension in
## @var{vecdim}.
##
## Specifying the dimension as @qcode{"all"} will cause @code{iqr} to operate
## on all elements of @var{x}, and is equivalent to @code{iqr (@var{x}(:))}.
##
## The fourth input argument, @var{methods}, determines the method to calculate
## the quantiles specified by @var{p} or @var{n}.  The methods available to
## calculate sample quantiles are the nine methods used by R
## (@url{https://www.r-project.org/}) and can be specified by the corresponding
## integer value.  The default value is @w{@var{method} = 5}.
##
## Discontinuous sample quantile methods 1, 2, and 3
##
## @enumerate 1
## @item Method 1: Inverse of empirical distribution function.
##
## @item Method 2: Similar to method 1 but with averaging at discontinuities.
##
## @item Method 3: SAS definition: nearest even order statistic.
## @end enumerate
##
## Continuous sample quantile methods 4 through 9, where
## @tex
## $p(k)$
## @end tex
## @ifnottex
## @var{p}(k)
## @end ifnottex
## is the linear
## interpolation function respecting each method's representative cdf.
##
## @enumerate 4
## @item Method 4:
## @tex
## $p(k) = k / N$.
## @end tex
## @ifnottex
## @var{p}(k) = k / N.
## @end ifnottex
## That is, linear interpolation of the empirical cdf, where @math{N} is the
## length of @var{P}.
##
## @item Method 5:
## @tex
## $p(k) = (k - 0.5) / N$.
## @end tex
## @ifnottex
## @var{p}(k) = (k - 0.5) / N.
## @end ifnottex
## That is, a piecewise linear function where the knots are the values midway
## through the steps of the empirical cdf.
##
## @item Method 6:
## @tex
## $p(k) = k / (N + 1)$.
## @end tex
## @ifnottex
## @var{p}(k) = k / (N + 1).
## @end ifnottex
##
## @item Method 7:
## @tex
## $p(k) = (k - 1) / (N - 1)$.
## @end tex
## @ifnottex
## @var{p}(k) = (k - 1) / (N - 1).
## @end ifnottex
##
## @item Method 8:
## @tex
## $p(k) = (k - 1/3) / (N + 1/3)$.
## @end tex
## @ifnottex
## @var{p}(k) = (k - 1/3) / (N + 1/3).
## @end ifnottex
## The resulting quantile estimates are approximately median-unbiased
## regardless of the distribution of @var{x}.
##
## @item Method 9:
## @tex
## $p(k) = (k - 3/8) / (N + 1/4)$.
## @end tex
## @ifnottex
## @var{p}(k) = (k - 3/8) / (N + 1/4).
## @end ifnottex
## The resulting quantile estimates are approximately unbiased for the
## expected order statistics if @var{x} is normally distributed.
## @end enumerate
##
## @nospell{Hyndman and Fan} (1996) recommend method 8.  Maxima, S, and R
## (versions prior to 2.0.0) use 7 as their default.  Minitab and SPSS
## use method 6.  @sc{matlab} uses method 5.
##
## References:
##
## @itemize @bullet
## @item @nospell{R. A. Becker, J. M. Chambers, and A. R. Wilks},
## @cite{The New S Language}, @nospell{Wadsworth & Brooks/Cole}, 1988.
##
## @item @nospell{R. J. Hyndman, and Y. Fan}, "Sample quantiles in statistical
## packages", @cite{American Statistician}, 50, @w{pp.@: 361}--365, 1996.
##
## @item @cite{R: A Language and Environment for Statistical Computing},
## @url{https://cran.r-project.org/doc/manuals/fullrefman.pdf}.
## @end itemize
##
## Examples:
## @c Set example in small font to prevent overfull line
##
## @smallexample
## @group
## x = randi (1000, [10, 1]);  # Create empirical data in range 1-1000
## q = quantile (x, [0, 1]);   # Return minimum, maximum of distribution
## q = quantile (x, [0.25 0.5 0.75]); # Return quartiles of distribution
## @end group
## @end smallexample
## @seealso{prctile}
## @end deftypefn

function q = quantile (x, p = [], dim, method = 5)

  if (nargin < 1)
    print_usage ();
  endif

  if (! (isnumeric (x)))
    error ("quantile: X must be a numeric array");
  endif

  if (isempty (p))
    p = [0.00, 0.25, 0.50, 0.75, 1.00];
  endif

  if (! (isnumeric (p) && isvector (p)))
    error ("quantile: P must be a numeric vector");
  endif

  if (isscalar (p) && fix (p) == p && p > 1)
    p = [1:p] ./ (p + 1);
  elseif (any (p < 0 | p > 1))
    error (strcat ("quantile: P values must range from 0 to 1, unless", ...
                   " specifying N evenly spaced cumulative probabilities"));
  endif

  do_perm = false;
  nd = ndims (x);
  sz = size (x);
  empty_x = isempty (x);

  if (nargin < 3)
    ## Find the first non-singleton dimension.
    (dim = find (size (x) > 1, 1)) || (dim = 1);

    ## Return immediately for an empty matrix.
    if (empty_x)
      if (nd == 2 && max (sz) <= 1)
        ## Return the size of vector P
        sz = size (p);
      else
        ## Reduce operating DIM to length of P
        sz(dim) = numel (p);
      endif
      q = NaN (sz);
      return;
    endif
  endif

  if (isnumeric (dim))

    ## Check for DIM argument
    if (isscalar (dim))
      if (! (dim == fix (dim) && dim > 0))
        error ("quantile: DIM must be a positive integer");
      endif

      ## Return immediately for an empty matrix.
      if (empty_x)
        ## Reduce operating DIM to length of P
        sz(dim) = numel (p);
        ## Mask existing dims, new zero-dims must be 1
        mask = ones (1, dim);
        mask(1:nd) = 0;
        sz(mask & sz == 0) = 1;
        q = NaN (sz);
        return;
      endif

      ## Return numel (p) copies of X if DIM > nd
      if (dim > nd)
        sz = ones (1, dim);
        sz(dim) = numel (p);
        q = repmat (x, sz);
        return;
      endif

      ## Set the permutation vector.
      perm = 1:(max (ndims (x), dim));
      perm(1) = dim;
      perm(dim) = 1;
      do_perm = true;

      ## Permute dim to the 1st index.
      x = permute (x, perm);

      ## Save the size of the permuted x N-D array.
      sx = size (x);

      ## Reshape to a 2-D array.
      x = reshape (x, sx(1), []);

    ## Check for proper VECDIM (more than 1 dim, no repeats)
    elseif (isvector (dim) && isindex (dim) && all (diff (sort (dim))))

      ## Discard exceeding dims, unless all dims > nd so keep smallest
      vecdim = dim(dim <= nd);
      if (isempty (vecdim))
        dim = min (dim);
      else
        dim = vecdim;
      endif

      ## Return immediately for an empty matrix.
      if (empty_x)
        sz(dim(1)) = numel (p);   # reduce first operating VECDIM to P
        sz(dim(2:end)) = 1;       # reduce other operating VECDIM to 1
        ## Mask existing dims, new zero-dims must be 1
        mask = ones (1, max (dim));
        mask(1:nd) = 0;
        sz(mask & sz == 0) = 1;
        q = NaN (sz);
        return;
      endif

      ## Return numel (p) copies of X if remaining DIM > nd
      if (dim > nd)
        sz = ones (1, dim);
        sz(dim) = numel (p);
        q = repmat (x, sz);
        return;
      endif

      ## Return X with X(dim(1)) expanded to P, if all DIM == 1
      if (all (sz(dim) == 1))
        sz = ones (1, nd);
        sz(dim(1)) = numel (p);   # reduce first operating VECDIM to P
        sz(dim(2:end)) = 1;       # reduce other operating VECDIM to 1
        q = repmat (x, sz);
        return;
      endif

      ## Detect trivial case of DIM being all dimensions (same as "all").
      vecdims = numel (dim);
      max_dim = max (nd, max (dim));
      if (vecdims == nd && max_dim == nd)
        x = x(:);
        sx = size (x);
        dim = 1;
      else
        ## Algorithm: Move dimensions for operation to the front, keeping the
        ## order of the remaining dimensions.  Reshape the moved dims into a
        ## single dimension (row).  Calculate with __quantile__ along dim1 of
        ## X, then reshape to correct dimensions.

        dim = dim(:).';  # Force row vector

        ## Permutation vector with DIM at front
        perm = [1:max_dim];
        perm(dim) = [];
        perm = [dim, perm];
        do_perm = true;

        ## Reshape X to put dims to process at front.
        x = permute (x, perm);
        sx = size (x);

        ## Preserve trailing singletons when dim > ndims (x).
        sx = [sx, ones(1, max_dim - numel (sx))];
        sx = [prod(sx(1:vecdims)), ones(1, (vecdims-1)), sx((vecdims+1):end)];

        ## Size must always have 2 dimensions.
        if (isscalar (sx))
          sx = [sx, 1];
        endif

        ## Collapse dimensions to be processsed into single column.
        x = reshape (x, sx);
      endif

    else
      error ("quantile: VECDIM must be a vector of non-repeating positive integers");
    endif

  elseif (strcmpi (dim, "all"))

    ## Return immediately for an empty matrix
    if (empty_x)
      ## Always return a column vector.
      q = NaN (numel (p), 1);
      return;
    endif

    ## "all" simplifies to collapsing all elements to single vector.
    x = x(:);
    sx = size (x);
    dim = 1;

  else
    error ("quantile: DIM must be a positive integer scalar, vector, or 'all'");
  endif

  ## Calculate the quantiles.
  q = __quantile__ (x, p, method);

  ## Return the shape to the original N-D array.
  q = reshape (q, [numel(p), sx(2:end)]);

  ## Permute the 1st index back to dim.
  if (do_perm)
    q = ipermute (q, perm);
  endif

  ## For Matlab compatibility, return vectors with the same orientation as p
  if (isvector (q) && ! isscalar (q) && ! isscalar (p))
    if (isrow (p))
      q = reshape (q, 1, []);
    else
      q = reshape (q, [], 1);
    endif
  endif

endfunction


%!test
%! p = 0.50;
%! q = quantile (1:4, p);
%! qa = 2.5;
%! assert (q, qa);
%! q = quantile (1:4, p, 1);
%! qa = [1, 2, 3, 4];
%! assert (q, qa);
%! q = quantile (1:4, p, 2);
%! qa = 2.5;
%! assert (q, qa);

%!test
%! p = [0.50 0.75];
%! q = quantile (1:4, p);
%! qa = [2.5 3.5];
%! assert (q, qa);
%! q = quantile (1:4, p, 1);
%! qa = [1, 2, 3, 4; 1, 2, 3, 4];
%! assert (q, qa);
%! q = quantile (1:4, p, 2);
%! qa = [2.5 3.5];
%! assert (q, qa);

%!test
%! p = 0.5;
%! x = sort (rand (11));
%! q = quantile (x, p);
%! assert (q, x(6,:));
%! x = x.';
%! q = quantile (x, p, 2);
%! assert (q, x(:,6));

%!test
%! p = [0.00, 0.25, 0.50, 0.75, 1.00];
%! x = [1; 2; 3; 4];
%! a = [1.0000   1.0000   2.0000   3.0000   4.0000
%!      1.0000   1.5000   2.5000   3.5000   4.0000
%!      1.0000   1.0000   2.0000   3.0000   4.0000
%!      1.0000   1.0000   2.0000   3.0000   4.0000
%!      1.0000   1.5000   2.5000   3.5000   4.0000
%!      1.0000   1.2500   2.5000   3.7500   4.0000
%!      1.0000   1.7500   2.5000   3.2500   4.0000
%!      1.0000   1.4167   2.5000   3.5833   4.0000
%!      1.0000   1.4375   2.5000   3.5625   4.0000];
%! for m = 1:9
%!   q = quantile (x, p, 1, m);
%!   assert (q, a(m,:), 0.0001);
%! endfor

%!test
%! p = [0.00, 0.25, 0.50, 0.75, 1.00];
%! x = [1; 2; 3; 4; 5];
%! a = [1.0000   2.0000   3.0000   4.0000   5.0000
%!      1.0000   2.0000   3.0000   4.0000   5.0000
%!      1.0000   1.0000   2.0000   4.0000   5.0000
%!      1.0000   1.2500   2.5000   3.7500   5.0000
%!      1.0000   1.7500   3.0000   4.2500   5.0000
%!      1.0000   1.5000   3.0000   4.5000   5.0000
%!      1.0000   2.0000   3.0000   4.0000   5.0000
%!      1.0000   1.6667   3.0000   4.3333   5.0000
%!      1.0000   1.6875   3.0000   4.3125   5.0000];
%! for m = 1:9
%!   q = quantile (x, p, 1, m);
%!   assert (q, a(m,:), 0.0001);
%! endfor

%!test
%! p = [0.00, 0.25, 0.50, 0.75, 1.00];
%! x = [1; 2; 5; 9];
%! a = [1.0000   1.0000   2.0000   5.0000   9.0000
%!      1.0000   1.5000   3.5000   7.0000   9.0000
%!      1.0000   1.0000   2.0000   5.0000   9.0000
%!      1.0000   1.0000   2.0000   5.0000   9.0000
%!      1.0000   1.5000   3.5000   7.0000   9.0000
%!      1.0000   1.2500   3.5000   8.0000   9.0000
%!      1.0000   1.7500   3.5000   6.0000   9.0000
%!      1.0000   1.4167   3.5000   7.3333   9.0000
%!      1.0000   1.4375   3.5000   7.2500   9.0000];
%! for m = 1:9
%!   q = quantile (x, p, 1, m);
%!   assert (q, a(m,:), 0.0001);
%! endfor

%!test
%! p = [0.00, 0.25, 0.50, 0.75, 1.00];
%! x = [1; 2; 5; 9; 11];
%! a = [1.0000    2.0000    5.0000    9.0000   11.0000
%!      1.0000    2.0000    5.0000    9.0000   11.0000
%!      1.0000    1.0000    2.0000    9.0000   11.0000
%!      1.0000    1.2500    3.5000    8.0000   11.0000
%!      1.0000    1.7500    5.0000    9.5000   11.0000
%!      1.0000    1.5000    5.0000   10.0000   11.0000
%!      1.0000    2.0000    5.0000    9.0000   11.0000
%!      1.0000    1.6667    5.0000    9.6667   11.0000
%!      1.0000    1.6875    5.0000    9.6250   11.0000];
%! for m = 1:9
%!   q = quantile (x, p, 1, m);
%!   assert (q, a(m,:), 0.0001);
%! endfor

%!test
%! p = [0.00, 0.25, 0.50, 0.75, 1.00];
%! x = [16; 11; 15; 12; 15;  8; 11; 12;  6; 10];
%! a = [6.0000   10.0000   11.0000   15.0000   16.0000
%!      6.0000   10.0000   11.5000   15.0000   16.0000
%!      6.0000    8.0000   11.0000   15.0000   16.0000
%!      6.0000    9.0000   11.0000   13.5000   16.0000
%!      6.0000   10.0000   11.5000   15.0000   16.0000
%!      6.0000    9.5000   11.5000   15.0000   16.0000
%!      6.0000   10.2500   11.5000   14.2500   16.0000
%!      6.0000    9.8333   11.5000   15.0000   16.0000
%!      6.0000    9.8750   11.5000   15.0000   16.0000];
%! for m = 1:9
%!   q = quantile (x, p, 1, m);
%!   assert (q, a(m,:), 0.0001);
%! endfor

%!test
%! p = [0.00, 0.25, 0.50, 0.75, 1.00];
%! x = [-0.58851;  0.40048;  0.49527; -2.551500; -0.52057; ...
%!      -0.17841; 0.057322; -0.62523;  0.042906;  0.12337];
%! a = [-2.551474  -0.588505  -0.178409   0.123366   0.495271
%!      -2.551474  -0.588505  -0.067751   0.123366   0.495271
%!      -2.551474  -0.625231  -0.178409   0.123366   0.495271
%!      -2.551474  -0.606868  -0.178409   0.090344   0.495271
%!      -2.551474  -0.588505  -0.067751   0.123366   0.495271
%!      -2.551474  -0.597687  -0.067751   0.192645   0.495271
%!      -2.551474  -0.571522  -0.067751   0.106855   0.495271
%!      -2.551474  -0.591566  -0.067751   0.146459   0.495271
%!      -2.551474  -0.590801  -0.067751   0.140686   0.495271];
%! for m = 1:9
%!   q = quantile (x, p, 1, m);
%!   assert (q, a(m,:), 0.0001);
%! endfor

%!test
%! p = 0.5;
%! x = [0.112600, 0.114800, 0.052100, 0.236400, 0.139300
%!      0.171800, 0.727300, 0.204100, 0.453100, 0.158500
%!      0.279500, 0.797800, 0.329600, 0.556700, 0.730700
%!      0.428800, 0.875300, 0.647700, 0.628700, 0.816500
%!      0.933100, 0.931200, 0.963500, 0.779600, 0.846100];
%! tol = 0.00001;
%! x(5,5) = NaN;
%! assert (quantile (x, p, 1),
%!         [0.27950, 0.79780, 0.32960, 0.55670, 0.44460], tol);
%! x(1,1) = NaN;
%! assert (quantile (x, p, 1),
%!         [0.35415, 0.79780, 0.32960, 0.55670, 0.44460], tol);
%! x(3,3) = NaN;
%! assert (quantile (x, p, 1),
%!         [0.35415, 0.79780, 0.42590, 0.55670, 0.44460], tol);

%!test
%! sx = [2, 3, 4];
%! x = rand (sx);
%! dim = 2;
%! p = 0.5;
%! yobs = quantile (x, p, dim);
%! yexp = median (x, dim);
%! assert (yobs, yexp);

%!assert <*45455> (quantile ([1 3 2], 0.5, 1), [1 3 2])
%!assert <*54421> (quantile ([1:10], 0.5, 1), 1:10)
%!assert <*54421> (quantile ([1:10]', 0.5, 2), [1:10]')
%!assert <*54421> (quantile ([1:10], [0.25, 0.75]), [3, 8])
%!assert <*54421> (quantile ([1:10], [0.25, 0.75]'), [3; 8])
%!assert (quantile ([1:10], 1, 3), [1:10])

## Test empty input arrays
%!assert (quantile ([], [0.2, 0.5, 0.7]), NaN (1, 3))
%!assert (quantile (ones (1, 0), [0.2, 0.5, 0.7]), NaN (1, 3))
%!assert (quantile ([], [0.2, 0.5, 0.7, 0.9]), NaN (1, 4))
%!assert (quantile (ones (1, 0), [0.2, 0.5, 0.7, 0.9]), NaN (1, 4))
%!assert (quantile ([], [0.2, 0.5, 0.7], 2), NaN (0, 3))
%!assert (quantile (ones (1, 0), [0.2, 0.5, 0.7], 2), NaN (1, 3))
%!assert (quantile (ones (0, 1), [0.2, 0.5, 0.7], 2), NaN (0, 3))
%!assert (quantile (ones (1, 0), [0.2, 0.5, 0.7], 1), NaN (3, 0))
%!assert (quantile (ones (1, 0), [0.2, 0.5, 0.7], 3), NaN (1, 0, 3))
%!assert (quantile (ones (1, 0, 1), [0.2, 0.5, 0.7], 3), NaN (1, 0, 3))
%!assert (quantile (ones (3, 0, 1, 2), [0.2, 0.5, 0.7]), NaN (3, 0, 1, 2))
%!assert (quantile (ones (3, 0, 1, 2), [0.2, 0.5, 0.7], 2), NaN (3, 3, 1, 2))
%!assert (quantile (ones (3, 0, 1, 2), [0.2; 0.5; 0.7]), NaN (3, 0, 1, 2))
%!assert (quantile (ones (3, 0, 1, 2), [0.2; 0.5; 0.7], 2), NaN (3, 3, 1, 2))
%!assert (quantile (ones (1, 0), [0.2; 0.5; 0.7]), NaN (3, 1))
%!assert (quantile (ones (0, 1), [0.2, 0.5, 0.7]), NaN (1, 3))
%!assert (quantile (ones (5, 0, 1, 2), [0.2, 0.5, 0.7]), NaN (3, 0, 1, 2))
%!assert (quantile (ones (5, 0), [0.2, 0.5, 0.7]), NaN (3, 0))
%!assert (quantile (ones (1, 0), [0.2, 0.5, 0.7], 4), NaN (1, 0, 1, 3))
%!assert (quantile (ones (1, 0), [0.2, 0.5, 0.7], 5), NaN (1, 0, 1, 1, 3))
%!assert (quantile (ones (5, 0, 1), [0.2, 0.5, 0.7], [4:6]), NaN (5, 0, 1, 3))
%!assert (quantile (ones (5, 0, 1), [0.2, 0.5, 0.7], [3, 4]), NaN (5, 0, 3))
%!assert (quantile (ones (5, 0, 2, 2), [0.2, 0.5, 0.7], [2, 3]), NaN (5, 3, 1, 2))
%!assert (quantile (ones (5, 0, 2, 2), [0.2, 0.5, 0.7], [1, 3]), NaN (3, 0, 1, 2))
%!assert (quantile (ones (5, 0, 2, 2), [0.2, 0.5, 0.7], 'all'), NaN (3, 1))
%!assert (quantile (ones (5, 0, 2, 2), [0.2; 0.5; 0.7], 'all'), NaN (3, 1))
%!assert (quantile (ones (0, 1), [0.2, 0.5, 0.7], 'all'), NaN (3, 1))
%!assert (quantile (ones (0, 1), [0.2; 0.5; 0.7], 'all'), NaN (3, 1))
%!assert (quantile (ones (1, 0), [0.2, 0.5, 0.7], 'all'), NaN (3, 1))
%!assert (quantile (ones (1, 0), [0.2; 0.5; 0.7], 'all'), NaN (3, 1))

## Test DIM and VECDIM with exceeding dimensions
%!assert (quantile (repmat ([1:10], 5, 1), [0.2, 0.5, 0.7], 2), ...
%!        repmat ([2.5, 5.5, 7.5], 5, 1))
%!assert (quantile (repmat ([1:10], 5, 1), [0.2, 0.5, 0.7], 3), ...
%!        repmat ([1:10], 5, 1, 3))
%!assert (quantile (repmat ([1:10], 5, 1), [0.2, 0.5, 0.7], 4), ...
%!        repmat ([1:10], 5, 1, 1, 3))
%!assert (quantile (repmat ([1:10], 5, 1), [0.2, 0.5, 0.7], [2, 4]), ...
%!        repmat ([2.5, 5.5, 7.5], 5, 1))
%!assert (quantile (repmat ([1:10], 5, 1), [0.2, 0.5, 0.7], [3, 5]), ...
%!        repmat ([1:10], 5, 1, 3))
%!assert (quantile (repmat ([1:10], 5, 1), [0.2, 0.5, 0.7], [4, 6]), ...
%!        repmat ([1:10], 5, 1, 1, 3))

## Test DIM and VECCDIM with dimensions of length 1
%!assert (quantile (ones (5, 1, 2, 2), [0.2, 0.5, 0.7], 2), ones (5, 3, 2, 2))
%!assert (quantile (ones (5, 1, 1, 2), [0.2, 0.5, 0.7], [2, 3]), ones (5, 3, 1, 2))

## Test direction of P vector
%!assert (quantile ([1:10], [0.2; 0.5; 0.7]), [2.5; 5.5; 7.5])
%!assert (quantile ([1:10]', [0.2; 0.5; 0.7]), [2.5; 5.5; 7.5])
%!assert (quantile ([1:10], [0.2, 0.5, 0.7]), [2.5, 5.5, 7.5])
%!assert (quantile ([1:10]', [0.2, 0.5, 0.7]), [2.5, 5.5, 7.5])

## Test N evenly spaced cummulative probabilities
%!assert (quantile ([1:10], 3), [3, 5.5, 8])
%!assert (quantile ([1:10]', 3), [3, 5.5, 8])
%!assert (quantile ([1:10], 9), [1.5, 2.5, 3.5, 4.5, 5.5, 6.5, 7.5, 8.5, 9.5])
%!assert (quantile ([1:10]', 9), [1.5, 2.5, 3.5, 4.5, 5.5, 6.5, 7.5, 8.5, 9.5])

## Test input validation
%!error <Invalid call> quantile ()
%!error <quantile: X must be a numeric array> quantile (['A'; 'B'], 10)
%!error <quantile: X must be a numeric array> quantile ([true; false], 10)
%!error <quantile: X must be a numeric array> quantile ({1, 2}, [0.2, 0.5, 0.8])
%!error <P must be a numeric vector> quantile (1:10, [true, false])
%!error <P must be a numeric vector> quantile (1:10, ones (2,2))
%!error <quantile: P values must range from 0 to 1, unless specifying N evenly spaced cumulative probabilities> ...
%!       quantile (1:10, -1)
%!error <quantile: P values must range from 0 to 1, unless specifying N evenly spaced cumulative probabilities> ...
%!       quantile (1:10, [0.2, 0.5, -0.8])
%!error <quantile: DIM must be a positive integer> quantile (1, 1, 1.5)
%!error <quantile: DIM must be a positive integer> quantile (1, 1, 0)
%!error <quantile: VECDIM must be a vector of non-repeating positive integers> ...
%!       quantile (1, 1, [1, 2, 2])
%!error <quantile: VECDIM must be a vector of non-repeating positive integers> ...
%!       quantile (1, 1, [1, 2, 0])
%!error <quantile: DIM must be a positive integer scalar, vector, or 'all'> ...
%!       quantile (1, 1, "some")
%!error quantile ((1:5)', 0.5, 1, 0)
%!error quantile ((1:5)', 0.5, 1, 10)

## For the cumulative probability values in @var{p}, compute the
## quantiles, @var{q} (the inverse of the cdf), for the sample, @var{x}.
##
## The optional input, @var{method}, refers to nine methods available in R
## (https://www.r-project.org/).  The default is @var{method} = 5.
## @seealso{prctile, quantile, statistics}

## Description: Quantile function of empirical samples
function inv = __quantile__ (x, p, method = 5)

  if (nargin < 2)
    print_usage ("quantile");
  endif

  if (isinteger (x) || islogical (x))
    x = double (x);
  endif

  ## set shape of quantiles to column vector.
  p = p(:);

  ## Save length and set shape of samples.
  x = sort (x, 1);
  m = sum (! isnan (x));
  [xr, xc] = size (x);

  ## Initialize output values.
  inv = Inf (class (x)) * (-(p < 0) + (p > 1));
  inv = repmat (inv, 1, xc);

  ## Do the work.
  if (any (k = find ((p >= 0) & (p <= 1))))
    n = length (k);
    p = p(k);
    ## Special case of 1 row.
    if (xr == 1)
      inv(k,:) = repmat (x, n, 1);
      return;
    endif

    ## The column-distribution indices.
    pcd = kron (ones (n, 1), xr*(0:xc-1));
    mm = kron (ones (n, 1), m);
    switch (method)
      case {1, 2, 3}
        switch (method)
          case 1
            p = max (ceil (kron (p, m)), 1);
            inv(k,:) = x(p + pcd);

          case 2
            p = kron (p, m);
            p_lr = max (ceil (p), 1);
            p_rl = min (floor (p + 1), mm);
            inv(k,:) = (x(p_lr + pcd) + x(p_rl + pcd))/2;

          case 3
           ## Used by SAS, method PCTLDEF=2.
           ## http://support.sas.com/onlinedoc/913/getDoc/en/statug.hlp/stdize_sect14.htm
            t = max (kron (p, m), 1);
            t = roundb (t);
            inv(k,:) = x(t + pcd);
        endswitch

      otherwise
        switch (method)
          case 4
            p = kron (p, m);

          case 5
            ## Used by Matlab.
            p = kron (p, m) + 0.5;

          case 6
            ## Used by Minitab and SPSS.
            p = kron (p, m+1);

          case 7
            ## Used by S and R.
            p = kron (p, m-1) + 1;

          case 8
            ## Median unbiased.
            p = kron (p, m+1/3) + 1/3;

          case 9
            ## Approximately unbiased respecting order statistics.
            p = kron (p, m+0.25) + 0.375;

          otherwise
            error ("quantile: Unknown METHOD, '%d'", method);
        endswitch

        ## Duplicate single values.
        imm1 = (mm(1,:) == 1);
        x(2,imm1) = x(1,imm1);

        ## Interval indices.
        pi = max (min (floor (p), mm-1), 1);
        pr = max (min (p - pi, 1), 0);
        pi += pcd;
        inv(k,:) = (1-pr) .* x(pi) + pr .* x(pi+1);
    endswitch
  endif

endfunction