Mercurial > octave
view scripts/sparse/private/__sprand__.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 | e1788b1a315f |
line wrap: on
line source
######################################################################## ## ## Copyright (C) 2004-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/>. ## ######################################################################## ## ## Original version by Paul Kienzle distributed as free software in the ## public domain. ## -*- texinfo -*- ## @deftypefn {} {@var{S} =} __sprand__ (@var{s}, @var{randfun}) ## @deftypefnx {} {@var{S} =} __sprand__ (@var{m}, @var{n}, @var{d}, @var{fcnname}, @var{randfun}) ## @deftypefnx {} {@var{S} =} __sprand__ (@var{m}, @var{n}, @var{d}, @var{rc}, @var{fcnname}, @var{randfun}) ## Undocumented internal function. ## @end deftypefn ## Actual implementation of sprand and sprandn happens here. function S = __sprand__ (varargin) if (nargin == 2) [m, randfun] = deal (varargin{1:2}); [i, j] = find (m); [nr, nc] = size (m); S = sparse (i, j, randfun (size (i)), nr, nc); else if (nargin == 5) [m, n, d, fcnname, randfun] = deal (varargin{:}); else [m, n, d, rc, fcnname, randfun] = deal (varargin{:}); endif if (! (isscalar (m) && m == fix (m) && m >= 0)) error ("%s: M must be a non-negative integer", fcnname); endif if (! (isscalar (n) && n == fix (n) && n >= 0)) error ("%s: N must be a non-negative integer", fcnname); endif if (d < 0 || d > 1) error ("%s: density D must be between 0 and 1", fcnname); endif if (m == 0 || n == 0) S = sparse (m, n); return; endif if (nargin == 5) mn = m*n; k = round (d*mn); if (mn > sizemax ()) ## randperm will overflow, so use alternative methods idx = unique (fix (rand (1.01*k, 1) * mn)) + 1; ## idx contains random numbers in [1,mn] ## Generate 1% more random values than necessary in order to reduce the ## probability that there are less than k distinct values; maybe a ## better strategy could be used but I don't think it's worth the price. ## actual number of entries in S k = min (length (idx), k); j = floor ((idx(1:k) - 1) / m); i = idx(1:k) - j * m; j += 1; else idx = randperm (mn, k); [i, j] = ind2sub ([m, n], idx); endif S = sparse (i, j, randfun (k, 1), m, n); elseif (nargin == 6) ## Create a matrix with specified reciprocal condition number. if (! isscalar (rc) && ! isvector (rc)) error ("%s: RC must be a scalar or vector", fcnname); endif ## We want to reverse singular valued decomposition A=U*S*V'. ## First, first S is constructed and then U = U1*U2*..Un and ## V' = V1*V2*..Vn are seen as Jacobi rotation matrices with angles and ## planes of rotation randomized. Repeatedly apply rotations until the ## required density for A is achieved. if (isscalar (rc)) if (rc < 0 || rc > 1) error ("%s: reciprocal condition number RC must be between 0 and 1", fcnname); endif ## Reciprocal condition number is ratio of smallest SV to largest SV ## Generate singular values randomly and sort them to build S ## Random singular values in range [rc, 1]. v = rand (1, min (m,n)) * (1 - rc) + rc; v(1) = 1; v(end) = rc; v = sort (v, "descend"); S = sparse (diag (v, m, n)); else ## Only the min (m, n) greater singular values from rc vector are used. if (length (rc) > min (m,n)) rc = rc(1:min (m, n)); endif S = sparse (diag (sort (rc, "descend"), m, n)); endif Uinit = speye (m); Vinit = speye (n); k = round (d*m*n); while (nnz (S) < k) if (m > 1) ## Construct U randomized rotation matrix rot_angleu = 2 * pi * rand (); cu = cos (rot_angleu); su = sin (rot_angleu); rndtmp = randperm (m, 2); i = rndtmp(1); j = rndtmp(2); U = Uinit; U(i, i) = cu; U(i, j) = -su; U(j, i) = su; U(j, j) = cu; S = U * S; endif if (n > 1) ## Construct V' randomized rotation matrix rot_anglev = 2 * pi * rand (); cv = cos (rot_anglev); sv = sin (rot_anglev); rndtmp = randperm (n, 2); i = rndtmp(1); j = rndtmp(2); V = Vinit; V(i, i) = cv; V(i, j) = sv; V(j, i) = -sv; V(j, j) = cv; S *= V; endif endwhile endif endif endfunction