view scripts/sparse/sprandsym.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 597f3ee61a48
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/>.
##
########################################################################

## -*- texinfo -*-
## @deftypefn  {} {@var{S} =} sprandsym (@var{n}, @var{d})
## @deftypefnx {} {@var{S} =} sprandsym (@var{s})
## Generate a symmetric random sparse matrix.
##
## The size of the matrix will be @var{n}x@var{n}, with a density of values
## given by @var{d}.  @var{d} must be between 0 and 1 inclusive.  Values will
## be normally distributed with a mean of zero and a variance of 1.
##
## If called with a single matrix argument, a random sparse matrix is generated
## wherever the matrix @var{s} is nonzero in its lower triangular part.
## @seealso{sprand, sprandn, spones, sparse}
## @end deftypefn

function S = sprandsym (n, d)

  if (nargin < 1)
    print_usage ();
  endif

  if (nargin == 1)
    [i, j] = find (tril (n));
    [nr, nc] = size (n);
    S = sparse (i, j, randn (size (i)), nr, nc);
    S += tril (S, -1)';
    return;
  endif

  if (!(isscalar (n) && n == fix (n) && n >= 0))
    error ("sprandsym: N must be a non-negative integer 0");
  endif

  if (n == 0)
    S = sparse (n, n);
    return;
  endif

  if (d < 0 || d > 1)
    error ("sprandsym: density D must be between 0 and 1");
  endif

  ## Actual number of nonzero entries
  k = round (n^2*d);

  ## Diagonal nonzero entries, same parity as k
  r = pick_rand_diag (n, k);

  ## Off diagonal nonzero entries
  m = (k - r)/2;

  ondiag = randperm (n, r);
  offdiag = randperm (n*(n - 1)/2, m);

  ## Row index
  i = lookup (cumsum (0:n), offdiag - 1) + 1;

  ## Column index
  j = offdiag - (i - 1).*(i - 2)/2;

  diagvals = randn (1, r);
  offdiagvals = randn (1, m);

  S = sparse ([ondiag, i, j], [ondiag, j, i],
              [diagvals, offdiagvals, offdiagvals], n, n);

endfunction

function r = pick_rand_diag (n, k)

  ## Pick a random number R of entries for the diagonal of a sparse NxN
  ## symmetric square matrix with exactly K nonzero entries, ensuring
  ## that this R is chosen uniformly over all such matrices.
  ##
  ## Let D be the number of diagonal entries and M the number of
  ## off-diagonal entries.  Then K = D + 2*M.  Let A = N*(N-1)/2 be the
  ## number of available entries in the upper triangle of the matrix.
  ## Then, by a simple counting argument, there is a total of
  ##
  ##     T = nchoosek (N, D) * nchoosek (A, M)
  ##
  ## symmetric NxN matrices with a total of K nonzero entries and D on
  ## the diagonal.  Letting D range from mod (K,2) through min (N,K), and
  ## dividing by this sum, we obtain the probability P for D to be each
  ## of those values.
  ##
  ## However, we cannot use this form for computation, as the binomial
  ## coefficients become unmanageably large.  Instead, we use the
  ## successive quotients Q(i) = T(i+1)/T(i), which we easily compute to
  ## be
  ##
  ##               (N - D)*(N - D - 1)*M
  ##     Q =  -------------------------------
  ##            (D + 2)*(D + 1)*(A - M + 1)
  ##
  ## Then, after prepending 1, the cumprod of these quotients is
  ##
  ##      C = [ T(1)/T(1), T(2)/T(1), T(3)/T(1), ..., T(N)/T(1) ]
  ##
  ## Their sum is thus S = sum (T)/T(1), and then C(i)/S is the desired
  ## probability P(i) for i=1:N.  The cumsum will finally give the
  ## distribution function for computing the random number of entries on
  ## the diagonal R.
  ##
  ## Thanks to Zsbán Ambrus <ambrus@math.bme.hu> for most of the ideas
  ## of the implementation here, especially how to do the computation
  ## numerically to avoid overflow.

  ## Degenerate case
  if (k == 1)
    r = 1;
    return;
  endif

  ## Compute the stuff described above
  a = n*(n - 1)/2;
  d = [mod(k,2):2:min(n,k)-2];
  m = (k - d)/2;
  q = (n - d).*(n - d - 1).*m ./ (d + 2)./(d + 1)./(a - m + 1);

  ## Slight modification from discussion above: pivot around the max in
  ## order to avoid overflow (underflow is fine, just means effectively
  ## zero probabilities).
  [~, midx] = max (cumsum (log (q)));
  midx += 1;
  lc = fliplr (cumprod (1./q(midx-1:-1:1)));
  rc = cumprod (q(midx:end));

  ## Now c = t(i)/t(midx), so c > 1 == [].
  c = [lc, 1, rc];
  s = sum (c);
  p = c/s;

  ## Add final d
  d(end+1) = d(end) + 2;

  ## Pick a random r using this distribution
  r = d(sum (cumsum (p) < rand) + 1);

endfunction


%!test
%! s = sprandsym (10, 0.1);
%! assert (issparse (s));
%! assert (issymmetric (s));
%! assert (size (s), [10, 10]);
%! assert (nnz (s) / numel (s), 0.1, .01);

## Test 1-input calling form
%!test
%! s = sprandsym (sparse ([1 2 3], [3 2 3], [2 2 2]));
%! [i, j] = find (s);
%! assert (sort (i), [2 3]');
%! assert (sort (j), [2 3]');

## Test empty array creation
%!assert (size (sprandsym (0, 0.5)), [0, 0])

## Test input validation
%!error <Invalid call> sprandsym ()
%!error sprandsym (ones (3), 0.5)
%!error sprandsym (3.5, 0.5)
%!error sprandsym (-1, 0.5)
%!error sprandsym (3, -1)
%!error sprandsym (3, 2)