view scripts/statistics/discrete_rnd.m @ 30875:5d3faba0342e

doc: Ensure documentation lists output argument when it exists for all m-files. For new users of Octave it is best to show explicit calling forms in the documentation and to show a return argument when it exists. * bp-table.cc, shift.m, accumarray.m, accumdim.m, bincoeff.m, bitcmp.m, bitget.m, bitset.m, blkdiag.m, celldisp.m, cplxpair.m, dblquad.m, flip.m, fliplr.m, flipud.m, idivide.m, int2str.m, interpft.m, logspace.m, num2str.m, polyarea.m, postpad.m, prepad.m, randi.m, repmat.m, rng.m, rot90.m, rotdim.m, structfun.m, triplequad.m, uibuttongroup.m, uicontrol.m, uipanel.m, uipushtool.m, uitoggletool.m, uitoolbar.m, waitforbuttonpress.m, help.m, __additional_help_message__.m, hsv.m, im2double.m, im2frame.m, javachk.m, usejava.m, argnames.m, char.m, formula.m, inline.m, __vectorize__.m, findstr.m, flipdim.m, strmatch.m, vectorize.m, commutation_matrix.m, cond.m, cross.m, duplication_matrix.m, expm.m, orth.m, rank.m, rref.m, trace.m, vech.m, cast.m, compare_versions.m, delete.m, dir.m, fileattrib.m, grabcode.m, gunzip.m, inputname.m, license.m, list_primes.m, ls.m, mexext.m, movefile.m, namelengthmax.m, nargoutchk.m, nthargout.m, substruct.m, swapbytes.m, ver.m, verLessThan.m, what.m, fminunc.m, fsolve.m, fzero.m, optimget.m, __fdjac__.m, matlabroot.m, savepath.m, campos.m, camroll.m, camtarget.m, camup.m, camva.m, camzoom.m, clabel.m, diffuse.m, legend.m, orient.m, rticks.m, specular.m, thetaticks.m, xlim.m, xtickangle.m, xticklabels.m, xticks.m, ylim.m, ytickangle.m, yticklabels.m, yticks.m, zlim.m, ztickangle.m, zticklabels.m, zticks.m, ellipsoid.m, isocolors.m, isonormals.m, stairs.m, surfnorm.m, __actual_axis_position__.m, __pltopt__.m, close.m, graphics_toolkit.m, pan.m, print.m, printd.m, __ghostscript__.m, __gnuplot_print__.m, __opengl_print__.m, rotate3d.m, subplot.m, zoom.m, compan.m, conv.m, poly.m, polyaffine.m, polyder.m, polyint.m, polyout.m, polyreduce.m, polyvalm.m, roots.m, prefdir.m, prefsfile.m, profexplore.m, profexport.m, profshow.m, powerset.m, unique.m, arch_rnd.m, arma_rnd.m, autoreg_matrix.m, bartlett.m, blackman.m, detrend.m, durbinlevinson.m, fftconv.m, fftfilt.m, fftshift.m, fractdiff.m, hamming.m, hanning.m, hurst.m, ifftshift.m, rectangle_lw.m, rectangle_sw.m, triangle_lw.m, sinc.m, sinetone.m, sinewave.m, spectral_adf.m, spectral_xdf.m, spencer.m, ilu.m, __sprand__.m, sprand.m, sprandn.m, sprandsym.m, treelayout.m, beta.m, betainc.m, betaincinv.m, betaln.m, cosint.m, expint.m, factorial.m, gammainc.m, gammaincinv.m, lcm.m, nthroot.m, perms.m, reallog.m, realpow.m, realsqrt.m, sinint.m, hadamard.m, hankel.m, hilb.m, invhilb.m, magic.m, pascal.m, rosser.m, toeplitz.m, vander.m, wilkinson.m, center.m, corr.m, cov.m, discrete_cdf.m, discrete_inv.m, discrete_pdf.m, discrete_rnd.m, empirical_cdf.m, empirical_inv.m, empirical_pdf.m, empirical_rnd.m, kendall.m, kurtosis.m, mad.m, mean.m, meansq.m, median.m, mode.m, moment.m, range.m, ranks.m, run_count.m, skewness.m, spearman.m, statistics.m, std.m, base2dec.m, bin2dec.m, blanks.m, cstrcat.m, deblank.m, dec2base.m, dec2bin.m, dec2hex.m, hex2dec.m, index.m, regexptranslate.m, rindex.m, strcat.m, strjust.m, strtrim.m, strtrunc.m, substr.m, untabify.m, __have_feature__.m, __prog_output_assert__.m, __run_test_suite__.m, example.m, fail.m, asctime.m, calendar.m, ctime.m, date.m, etime.m: Add return arguments to @deftypefn macros where they were missing. Rename variables in functions (particularly generic "retval") to match documentation. Rename some return variables for (hopefully) better clarity (e.g., 'ax' to 'hax' to indicate it is a graphics handle to an axes object).
author Rik <rik@octave.org>
date Wed, 30 Mar 2022 20:40:27 -0700
parents 796f54d4ddbf
children 5d6b058a22dc
line wrap: on
line source

########################################################################
##
## Copyright (C) 1996-2022 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{rnd} =} discrete_rnd (@var{v}, @var{p})
## @deftypefnx {} {@var{rnd} =} discrete_rnd (@var{v}, @var{p}, @var{r})
## @deftypefnx {} {@var{rnd} =} discrete_rnd (@var{v}, @var{p}, @var{r}, @var{c}, @dots{})
## @deftypefnx {} {@var{rnd} =} discrete_rnd (@var{v}, @var{p}, [@var{sz}])
## Return a matrix of random samples from the univariate distribution which
## assumes the values in @var{v} with probabilities @var{p}.
##
## When called with a single size argument, return a square matrix with
## the dimension specified.  When called with more than one scalar argument the
## first two arguments are taken as the number of rows and columns and any
## further arguments specify additional matrix dimensions.  The size may also
## be specified with a vector of dimensions @var{sz}.
##
## If no size arguments are given then the result matrix is the common size of
## @var{v} and @var{p}.
## @end deftypefn

function rnd = discrete_rnd (v, p, varargin)

  if (nargin < 2)
    print_usage ();
  endif

  if (! isvector (v))
    error ("discrete_rnd: V must be a vector");
  elseif (! isvector (p) || (length (p) != length (v)))
    error ("discrete_rnd: P must be a vector with length (V) elements");
  elseif (any (isnan (p)))
    error ("discrete_rnd: P must not have any NaN elements");
  elseif (! (all (p >= 0) && any (p)))
    error ("discrete_rnd: P must be a nonzero, non-negative vector");
  endif

  if (nargin == 2)
    sz = size (v);
  elseif (nargin == 3)
    if (isscalar (varargin{1}) && varargin{1} >= 0)
      sz = [varargin{1}, varargin{1}];
    elseif (isrow (varargin{1}) && all (varargin{1} >= 0))
      sz = varargin{1};
    else
      error ("discrete_rnd: dimension vector must be row vector of non-negative integers");
    endif
  elseif (nargin > 3)
    if (any (cellfun (@(x) (! isscalar (x) || x < 0), varargin)))
      error ("discrete_rnd: dimensions must be non-negative integers");
    endif
    sz = [varargin{:}];
  endif

  rnd = v(lookup (cumsum (p(1:end-1)) / sum (p), rand (sz)) + 1);
  rnd = reshape (rnd, sz);

endfunction


%!assert (size (discrete_rnd (1:2, 1:2, 3)), [3, 3])
%!assert (size (discrete_rnd (1:2, 1:2, [4 1])), [4, 1])
%!assert (size (discrete_rnd (1:2, 1:2, 4, 1)), [4, 1])

## Test class of input preserved
%!assert (class (discrete_rnd (1:2, 1:2)), "double")
%!assert (class (discrete_rnd (single (1:2), 1:2)), "single")
## FIXME: Maybe this should work, maybe it shouldn't.
#%!assert (class (discrete_rnd (1:2, single(1:2))), "single")

## Test input validation
%!error <Invalid call> discrete_rnd ()
%!error <Invalid call> discrete_rnd (1)
%!error discrete_rnd (1:2,1:2, -1)
%!error discrete_rnd (1:2,1:2, ones (2))
%!error discrete_rnd (1:2,1:2, [2 -1 2])
%!error discrete_rnd (1:2,1:2, 1, ones (2))
%!error discrete_rnd (1:2,1:2, 1, -1)
## test v,p verification
%!error discrete_rnd (1, ones (2), ones (2,1))
%!error discrete_rnd (1, ones (2,1), ones (1,1))
%!error discrete_rnd (1, ones (2,1), [1 -1])
%!error discrete_rnd (1, ones (2,1), [1 NaN])
%!error discrete_rnd (1, ones (2,1), [0  0])