Mercurial > octave
view scripts/sparse/sprandn.m @ 28789:28de41192f3c
Eliminate unneeded verification of nargin, nargout in m-files.
* FIRfilter.m, FIRfilter_aggregation.m, get.m, polynomial.m,
polynomial_superiorto.m, polynomial2.m, makeUniqueStrings.m, base64decode.m,
base64encode.m, cd.m, lin2mu.m, record.m, sound.m, soundsc.m, accumarray.m,
accumdim.m, bitcmp.m, bitset.m, cart2pol.m, celldisp.m, circshift.m,
cplxpair.m, cumtrapz.m, flip.m, idivide.m, interpft.m, logspace.m, pol2cart.m,
polyarea.m, postpad.m, prepad.m, rat.m, rot90.m, rotdim.m, shift.m, shiftdim.m,
sortrows.m, trapz.m, dsearch.m, dsearchn.m, getappdata.m, getpixelposition.m,
guidata.m, guihandles.m, isappdata.m, listfonts.m, uigetdir.m,
waitforbuttonpress.m, __makeinfo__.m, doc.m, get_first_help_sentence.m,
autumn.m, bone.m, brighten.m, cmpermute.m, cmunique.m, colorcube.m, contrast.m,
cool.m, copper.m, cubehelix.m, flag.m, gray.m, gray2ind.m, hot.m, hsv.m,
im2double.m, im2frame.m, imformats.m, jet.m, lines.m, ocean.m, pink.m, prism.m,
rainbow.m, rgbplot.m, spinmap.m, spring.m, summer.m, viridis.m, white.m,
winter.m, beep.m, importdata.m, is_valid_file_id.m, javachk.m, javaclasspath.m,
findstr.m, genvarname.m, strmatch.m, bandwidth.m, commutation_matrix.m, cond.m,
cross.m, isdefinite.m, ishermitian.m, issymmetric.m, krylov.m, linsolve.m,
logm.m, lscov.m, null.m, ordeig.m, orth.m, rank.m, rref.m, vecnorm.m,
bunzip2.m, citation.m, computer.m, copyfile.m, dir.m, dos.m, fileattrib.m,
gunzip.m, inputParser.m, inputname.m, ismac.m, ispc.m, isunix.m, license.m,
list_primes.m, methods.m, mkdir.m, movefile.m, nargchk.m, news.m,
orderfields.m, recycle.m, tar.m, unix.m, unpack.m, untar.m, unzip.m, ver.m,
version.m, what.m, zip.m, decic.m, fminbnd.m, fminunc.m, fsolve.m, fzero.m,
glpk.m, humps.m, lsqnonneg.m, optimget.m, pqpnonneg.m, sqp.m, pathdef.m,
camlookat.m, hidden.m, specular.m, plotmatrix.m, smooth3.m, sombrero.m,
stemleaf.m, __gnuplot_drawnow__.m, __opengl_info__.m, ancestor.m, cla.m,
close.m, closereq.m, copyobj.m, gca.m, gcf.m, ginput.m, graphics_toolkit.m,
groot.m, hgload.m, hgsave.m, isgraphics.m, ishold.m, linkaxes.m, meshgrid.m,
newplot.m, refresh.m, refreshdata.m, rotate.m, saveas.m, struct2hdl.m, conv.m,
mkpp.m, mpoles.m, padecoef.m, pchip.m, polyder.m, polyfit.m, polygcd.m,
polyint.m, polyout.m, polyval.m, ppder.m, ppint.m, getpref.m, ispref.m,
rmpref.m, profexport.m, profshow.m, powerset.m, arch_fit.m, arma_rnd.m,
blackman.m, detrend.m, diffpara.m, fftconv.m, fftfilt.m, filter2.m, freqz.m,
freqz_plot.m, hamming.m, hanning.m, sinetone.m, sinewave.m, spectral_adf.m,
spectral_xdf.m, stft.m, unwrap.m, gplot.m, ichol.m, ilu.m, spdiags.m, sprand.m,
sprandn.m, spstats.m, svds.m, treelayout.m, treeplot.m, betainc.m,
betaincinv.m, ellipke.m, gammainc.m, gammaincinv.m, legendre.m, pow2.m,
hankel.m, pascal.m, rosser.m, toeplitz.m, bounds.m, corr.m, cov.m, histc.m,
kendall.m, kurtosis.m, mad.m, mode.m, moment.m, prctile.m, quantile.m, range.m,
ranks.m, run_count.m, skewness.m, spearman.m, std.m, var.m, zscore.m,
dec2base.m, dec2bin.m, dec2hex.m, index.m, mat2str.m, native2unicode.m,
ostrsplit.m, strjoin.m, strjust.m, strtok.m, substr.m, unicode2native.m,
untabify.m, __debug_octave__.m, demo.m, example.m, fail.m, oruntests.m,
dump_demos.m, speed.m, test.m, date.m, datenum.m, datestr.m, datevec.m,
is_leap_year.m, now.m, weekday.m:
Eliminate unneeded verification of nargin, nargout in m-files now that
the interpreter checks these values.
author | Rik <rik@octave.org> |
---|---|
date | Thu, 24 Sep 2020 14:44:58 -0700 |
parents | e9d57f6d6353 |
children | d8318c12d903 |
line wrap: on
line source
######################################################################## ## ## Copyright (C) 2004-2020 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 {} {} sprandn (@var{m}, @var{n}, @var{d}) ## @deftypefnx {} {} sprandn (@var{m}, @var{n}, @var{d}, @var{rc}) ## @deftypefnx {} {} sprandn (@var{s}) ## Generate a sparse matrix with normally distributed random values. ## ## The size of the matrix is @var{m}x@var{n} with a density of values @var{d}. ## @var{d} must be between 0 and 1. Values will be normally distributed with a ## mean of 0 and a variance of 1. ## ## If called with a single matrix argument, a sparse matrix is generated with ## random values wherever the matrix @var{s} is nonzero. ## ## If called with a scalar fourth argument @var{rc}, a random sparse matrix ## with reciprocal condition number @var{rc} is generated. If @var{rc} is ## a vector, then it specifies the first singular values of the generated ## matrix (@code{length (@var{rc}) <= min (@var{m}, @var{n})}). ## ## @seealso{sprand, sprandsym, randn} ## @end deftypefn function s = sprandn (m, n, d, rc) if (nargin == 1) s = __sprand__ (m, @randn); elseif (nargin == 3) s = __sprand__ (m, n, d, "sprandn", @randn); elseif (nargin == 4) s = __sprand__ (m, n, d, rc, "sprandn", @randn); else print_usage (); endif endfunction ## Test 3-input calling form %!test %! s = sprandn (4, 10, 0.1); %! assert (size (s), [4, 10]); %! assert (nnz (s) / numel (s), 0.1); ## Test 4-input calling form %!test %! d = rand (); %! s1 = sprandn (100, 100, d, 0.4); %! rc = [5, 4, 3, 2, 1, 0.1]; %! s2 = sprandn (100, 100, d, rc); %! s3 = sprandn (6, 4, d, rc); %! assert (svd (s2)'(1:length (rc)), rc, sqrt (eps)); %! assert (1/cond (s1), 0.4, sqrt (eps)); %! assert (nnz (s1) / (100*100), d, 0.02); %! assert (nnz (s2) / (100*100), d, 0.02); %! assert (svd (s3)', [5 4 3 2], sqrt (eps)); ## Test 1-input calling form %!test %! s = sprandn (sparse ([1 2 3], [3 2 3], [2 2 2])); %! [i, j] = find (s); %! assert (sort (i), [1 2 3]'); %! assert (sort (j), [2 3 3]'); ## Test very large, very low density matrix doesn't fail %!test %! s = sprandn (1e6,1e6,1e-7); ## Test empty array creation %!assert (size (sprandn (0, 0, 0.5)), [0, 0]) %!assert (size (sprandn (0, 3, 0.5)), [0, 3]) %!assert (size (sprandn (3, 0, 0.5)), [3, 0]) ## Test input validation %!error sprandn () %!error sprandn (1, 2) %!error sprandn (1, 2, 3, 4) %!error <M must be a non-negative integer> sprand (-1, -1, 0.5) %!error <M must be a non-negative integer> sprandn (ones (3), 3, 0.5) %!error <M must be a non-negative integer> sprandn (3.5, 3, 0.5) %!error <M must be a non-negative integer> sprandn (-1, 3, 0.5) %!error <N must be a non-negative integer> sprandn (3, ones (3), 0.5) %!error <N must be a non-negative integer> sprandn (3, 3.5, 0.5) %!error <N must be a non-negative integer> sprandn (3, -1, 0.5) %!error <D must be between 0 and 1> sprandn (3, 3, -1) %!error <D must be between 0 and 1> sprandn (3, 3, 2) %!error <RC must be a scalar or vector> sprandn (2, 2, 0.2, ones (3,3)) %!error <RC must be between 0 and 1> sprandn (2, 2, 0.2, -1) %!error <RC must be between 0 and 1> sprandn (2, 2, 0.2, 2)