changeset 18624:54a1e95365e1

Overhaul sprand, sprandn functions. * __sprand_impl__: Rename variable "funname" to "fcnname". Add comments to Reciprocal Condition number calculation. Rename "mynnz" to "k" to match rest of code. Add input validation test that RC is scalar or vector. Use double quotes instead of single quotes per Octave guidelines. Check for special case of output vector to avoid problems. Use randperm to replace do/until loop for speed. Pre-calculate speye() value instead of doing per loop iteration. * sprand.m: Improve docstring. Match function output variable name to documentation. Add check string to %!error tests. * sprandn.m: Improve docstring. Match function output variable name to documentation. Add check string to %!error tests.
author Rik <rik@octave.org>
date Sun, 23 Mar 2014 20:35:22 -0700
parents 35a5e7740a6d
children d0d9f6daa4b6
files scripts/sparse/private/__sprand_impl__.m scripts/sparse/sprand.m scripts/sparse/sprandn.m
diffstat 3 files changed, 148 insertions(+), 137 deletions(-) [+]
line wrap: on
line diff
--- a/scripts/sparse/private/__sprand_impl__.m	Sat Mar 22 13:23:41 2014 +0100
+++ b/scripts/sparse/private/__sprand_impl__.m	Sun Mar 23 20:35:22 2014 -0700
@@ -22,8 +22,8 @@
 
 ## -*- texinfo -*-
 ## @deftypefn  {Function File} {} __sprand_impl__ (@var{s}, @var{randfun})
-## @deftypefnx {Function File} {} __sprand_impl__ (@var{m}, @var{n}, @var{d}, @var{funname}, @var{randfun})
-## @deftypefnx {Function File} {} __sprand_impl__ (@var{m}, @var{n}, @var{d}, @var{rc}, @var{funname}, @var{randfun})
+## @deftypefnx {Function File} {} __sprand_impl__ (@var{m}, @var{n}, @var{d}, @var{fcnname}, @var{randfun})
+## @deftypefnx {Function File} {} __sprand_impl__ (@var{m}, @var{n}, @var{d}, @var{rc}, @var{fcnname}, @var{randfun})
 ## Undocumented internal function.
 ## @end deftypefn
 
@@ -32,46 +32,40 @@
 function S = __sprand_impl__ (varargin)
 
   if (nargin == 2)
-    m = varargin{1};
-    randfun = varargin{2};
+    [m, randfun] = deal (varargin{1:2});
     [i, j] = find (m);
     [nr, nc] = size (m);
     S = sparse (i, j, randfun (size (i)), nr, nc);
-    return;
   else
     if (nargin == 5)
-      [m, n, d, funname, randfun] = deal(varargin{:});
-    else 
-      [m, n, d, rc, funname, randfun] = deal(varargin{:});
+      [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 an integer greater than 0", funname);
+    if (! (isscalar (m) && m == fix (m) && m > 0))
+      error ("%s: M must be an integer greater than 0", fcnname);
     endif
-
-    if (!(isscalar (n) && n == fix (n) && n > 0))
-      error ("%s: N must be an integer greater than 0", funname);
+    if (! (isscalar (n) && n == fix (n) && n > 0))
+      error ("%s: N must be an integer greater than 0", fcnname);
     endif
-
     if (d < 0 || d > 1)
-      error ("%s: density D must be between 0 and 1", funname);
+      error ("%s: density D must be between 0 and 1", fcnname);
     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 (min (k*1.01, k+10), 1) * mn)) + 1;
+        idx = unique (fix (rand (1.01*k, 1) * mn)) + 1;
 
         ## idx contains random numbers in [1,mn]
-        ## generate 1% or 10 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
-        
+        ## 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);
@@ -83,59 +77,69 @@
       endif
 
       S = sparse (i, j, randfun (k, 1), m, n);
+
     elseif (nargin == 6)
-      ## We assume that we want to reverse A=U*S*V' so firstly 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. In the nth step the density required for A is achieved.
+      ## 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.
 
-      mynnz = round (m * n * d);
-      if (!isscalar(rc)) 
-        ## Only the min(m, n) greater singular values from rc vector are used. Needed to be compliant.
+      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));
-      else
-        if(rc < 0 || rc > 1)
-          error ("%s: reciprocal condition number rc must be between 0 and 1", funname);
-        endif
-        ## Generate the singular values randomly and sort them to build S
-        for (i = 1:min(m, n))
-          ## Randon singular values between 1 and rc.
-          v(i) = rand () * (1 - rc) + rc; 
-        endfor
-        v(1) = 1;
-        v(end) = rc;
-        v = sort (v, 'descend');
-        S = sparse (diag (v, m, n));
+        S = sparse (diag (sort (rc, "descend"), m, n));
       endif
-      while (nnz(S) < mynnz) 
-        [mm, nn] = size(S);
-        rot_angleu = 2 * randfun () * pi;
-        rot_anglev = 2 * randfun () * pi;
-        cu = cos (rot_angleu); cv = cos (rot_anglev); 
-        su = sin (rot_angleu); sv = sin (rot_anglev);
-        ## Rotation related with U
-        i = fix (rand () * m) + 1;
-        do
-          ## If j==i rotation matrix would be no longer that kind
-          j = fix (rand () * m) + 1;   
-        until (j != i)
-        U = sparse (eye (m,m));
-        U(i, i) = cu; U(i, j) = -su;
-        U(j, i) = su; U(j, j) = cu;
-        S = U * S;
-        ## Rotation related with V'
-        i = fix (rand () * nn) + 1;
-        do
-          j = fix (rand () * nn) + 1; 
-        until(j != i)
-        V = sparse (eye (n, n));
-        V(i, i) = cv; V(i, j) = sv;
-        V(j, i) = -sv; V(j, j) = cv;
-        S = S * V;
+
+      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 = S * V;
+        endif
       endwhile
     endif
   endif
+
 endfunction
 
--- a/scripts/sparse/sprand.m	Sat Mar 22 13:23:41 2014 +0100
+++ b/scripts/sparse/sprand.m	Sun Mar 23 20:35:22 2014 -0700
@@ -21,21 +21,23 @@
 
 ## -*- texinfo -*-
 ## @deftypefn  {Function File} {} sprand (@var{m}, @var{n}, @var{d})
-## @deftypefnx  {Function File} {} sprand (@var{m}, @var{n}, @var{d}, @var{rc})
+## @deftypefnx {Function File} {} sprand (@var{m}, @var{n}, @var{d}, @var{rc})
 ## @deftypefnx {Function File} {} sprand (@var{s})
-## Generate a random sparse matrix.  The size of the matrix will be
-## @var{m} by @var{n}, with a density of values given by @var{d}.
-## @var{d} should be between 0 and 1.  Values will be uniformly
-## distributed between 0 and 1.
+## Generate a sparse matrix with uniformly 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 uniformly distributed on
+## the interval (0, 1).
 ##
-## If called with a single matrix argument, a random sparse matrix is
-## generated wherever the matrix @var{S} is non-zero.
+## If called with a single matrix argument, a sparse matrix is generated with
+## random values wherever the matrix @var{s} is non-zero.
 ##
-## If called with the rc parameter as a scalar, a random sparse matrix with a
-## reciprocal condition number rc is generated. If rc is a vector, then it 
-## contains the first singular values of the matrix generated (length(rc) <= min(m, n)).
+## 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{sprandn, sprandsym}
+## @seealso{sprandn, sprandsym, rand}
 ## @end deftypefn
 
 ## Author: Paul Kienzle <pkienzle@users.sf.net>
@@ -48,14 +50,14 @@
 ## David Bateman
 ##      2004-10-20      Texinfo help and copyright message
 
-function S = sprand (m, n, d, rc)
+function s = sprand (m, n, d, rc)
 
   if (nargin == 1 )
-    S = __sprand_impl__ (m, @rand);
+    s = __sprand_impl__ (m, @rand);
   elseif ( nargin == 3)
-    S = __sprand_impl__ (m, n, d, "sprand", @rand);
+    s = __sprand_impl__ (m, n, d, "sprand", @rand);
   elseif (nargin == 4)
-    S = __sprand_impl__ (m, n, d, rc, "sprand", @rand);
+    s = __sprand_impl__ (m, n, d, rc, "sprand", @rand);
   else
     print_usage ();
   endif
@@ -63,22 +65,15 @@
 endfunction
 
 
+%% Test 3-input calling form
 %!test
 %! s = sprand (4, 10, 0.1);
 %! assert (size (s), [4, 10]);
 %! assert (nnz (s) / numel (s), 0.1);
 
-%% Test 1-input calling form
-%!test
-%! s = sprand (sparse ([1 2 3], [3 2 3], [2 2 2]));
-%! [i, j, v] = find (s);
-%! assert (sort (i), [1 2 3]');
-%! assert (sort (j), [2 3 3]');
-%! assert (all (v > 0 & v < 1));
-
 %% Test 4-input calling form
 %!test
-%! d = rand();
+%! d = rand ();
 %! s1 = sprand (100, 100, d, 0.4);
 %! rc = [5, 4, 3, 2, 1, 0.1];
 %! s2 = sprand (100, 100, d, rc);
@@ -89,22 +84,31 @@
 %! assert (nnz (s2) / (100*100), d, 0.02); 
 %! assert (svd (s3)', [5 4 3 2], sqrt (eps));
 
+%% Test 1-input calling form
+%!test
+%! s = sprand (sparse ([1 2 3], [3 2 3], [2 2 2]));
+%! [i, j, v] = find (s);
+%! assert (sort (i), [1 2 3]');
+%! assert (sort (j), [2 3 3]');
+%! assert (all (v > 0 & v < 1));
+
 %% Test input validation
 %!error sprand ()
 %!error sprand (1, 2)
 %!error sprand (1, 2, 3, 4)
-%!error sprand (ones (3), 3, 0.5)
-%!error sprand (3.5, 3, 0.5)
-%!error sprand (0, 3, 0.5)
-%!error sprand (3, ones (3), 0.5)
-%!error sprand (3, 3.5, 0.5)
-%!error sprand (3, 0, 0.5)
-%!error sprand (3, 3, -1)
-%!error sprand (3, 3, 2)
-%!error sprand (2, 2, 0.2, -1)
-%!error sprand (2, 2, 0.2, 2)
+%!error <M must be an integer greater than 0> sprand (ones (3), 3, 0.5)
+%!error <M must be an integer greater than 0> sprand (3.5, 3, 0.5)
+%!error <M must be an integer greater than 0> sprand (0, 3, 0.5)
+%!error <N must be an integer greater than 0> sprand (3, ones (3), 0.5)
+%!error <N must be an integer greater than 0> sprand (3, 3.5, 0.5)
+%!error <N must be an integer greater than 0> sprand (3, 0, 0.5)
+%!error <D must be between 0 and 1> sprand (3, 3, -1)
+%!error <D must be between 0 and 1> sprand (3, 3, 2)
+%!error <RC must be a scalar or vector> sprand (2, 2, 0.2, ones (3,3))
+%!error <RC must be between 0 and 1> sprand (2, 2, 0.2, -1)
+%!error <RC must be between 0 and 1> sprand (2, 2, 0.2, 2)
 
 %% Test very large, very low density matrix doesn't fail 
 %!test
-%! s = sprand(1e6,1e6,1e-7);
+%! s = sprand (1e6, 1e6, 1e-7);
 
--- a/scripts/sparse/sprandn.m	Sat Mar 22 13:23:41 2014 +0100
+++ b/scripts/sparse/sprandn.m	Sun Mar 23 20:35:22 2014 -0700
@@ -23,32 +23,33 @@
 ## @deftypefn  {Function File} {} sprandn (@var{m}, @var{n}, @var{d})
 ## @deftypefnx {Function File} {} sprandn (@var{m}, @var{n}, @var{d}, @var{rc})
 ## @deftypefnx {Function File} {} sprandn (@var{s})
-## Generate a random sparse matrix.  The size of the matrix will be
-## @var{m} by @var{n}, with a density of values given by @var{d}.
-## @var{d} should be between 0 and 1.  Values will be normally
-## distributed with mean of zero and variance 1.
+## 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 random sparse matrix is
-## generated wherever the matrix @var{S} is non-zero.
+## If called with a single matrix argument, a sparse matrix is generated with
+## random values wherever the matrix @var{s} is non-zero.
+##
+## 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})}).
 ## 
-## If called with the rc parameter, then a random sparse matrix with a
-## reciprocal condition number rc is generated if rc is scalar. If rc 
-## is a vector, then it contains the first singular values of the matrix
-## generated (length(rc) <= min(m, n)).
-## 
-## @seealso{sprand, sprandsym}
+## @seealso{sprand, sprandsym, randn}
 ## @end deftypefn
 
 ## Author: Paul Kienzle <pkienzle@users.sf.net>
 
-function S = sprandn (m, n, d, rc)
+function s = sprandn (m, n, d, rc)
 
   if (nargin == 1 )
-    S = __sprand_impl__ (m, @randn);
+    s = __sprand_impl__ (m, @randn);
   elseif ( nargin == 3)
-    S = __sprand_impl__ (m, n, d, "sprandn", @randn);
+    s = __sprand_impl__ (m, n, d, "sprandn", @randn);
   elseif (nargin == 4)
-    S = __sprand_impl__ (m, n, d, rc, "sprandn", @randn);
+    s = __sprand_impl__ (m, n, d, rc, "sprandn", @randn);
   else
     print_usage ();
   endif
@@ -56,21 +57,15 @@
 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 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 4-input calling form
 %!test
-%! d = rand();
+%! d = rand ();
 %! s1 = sprandn (100, 100, d, 0.4);
 %! rc = [5, 4, 3, 2, 1, 0.1];
 %! s2 = sprandn (100, 100, d, rc);
@@ -81,22 +76,30 @@
 %! 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 input validation
 %!error sprandn ()
 %!error sprandn (1, 2)
 %!error sprandn (1, 2, 3, 4)
-%!error sprandn (ones (3), 3, 0.5)
-%!error sprandn (3.5, 3, 0.5)
-%!error sprandn (0, 3, 0.5)
-%!error sprandn (3, ones (3), 0.5)
-%!error sprandn (3, 3.5, 0.5)
-%!error sprandn (3, 0, 0.5)
-%!error sprandn (3, 3, -1)
-%!error sprandn (3, 3, 2)
-%!error sprandn (2, 2, 0.2, -1)
-%!error sprandn (2, 2, 0.2, 2)
+%!error <M must be an integer greater than 0> sprandn (ones (3), 3, 0.5)
+%!error <M must be an integer greater than 0> sprandn (3.5, 3, 0.5)
+%!error <M must be an integer greater than 0> sprandn (0, 3, 0.5)
+%!error <N must be an integer greater than 0> sprandn (3, ones (3), 0.5)
+%!error <N must be an integer greater than 0> sprandn (3, 3.5, 0.5)
+%!error <N must be an integer greater than 0> sprandn (3, 0, 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)
 
 %% Test very large, very low density matrix doesn't fail 
 %!test
-%! s = sprandn(1e6,1e6,1e-7);
+%! s = sprandn (1e6,1e6,1e-7);