comparison scripts/sparse/sprand.m @ 18725: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 c53e11fab75f
comparison
equal deleted inserted replaced
18724:35a5e7740a6d 18725:54a1e95365e1
19 ## Original version by Paul Kienzle distributed as free software in the 19 ## Original version by Paul Kienzle distributed as free software in the
20 ## public domain. 20 ## public domain.
21 21
22 ## -*- texinfo -*- 22 ## -*- texinfo -*-
23 ## @deftypefn {Function File} {} sprand (@var{m}, @var{n}, @var{d}) 23 ## @deftypefn {Function File} {} sprand (@var{m}, @var{n}, @var{d})
24 ## @deftypefnx {Function File} {} sprand (@var{m}, @var{n}, @var{d}, @var{rc}) 24 ## @deftypefnx {Function File} {} sprand (@var{m}, @var{n}, @var{d}, @var{rc})
25 ## @deftypefnx {Function File} {} sprand (@var{s}) 25 ## @deftypefnx {Function File} {} sprand (@var{s})
26 ## Generate a random sparse matrix. The size of the matrix will be 26 ## Generate a sparse matrix with uniformly distributed random values.
27 ## @var{m} by @var{n}, with a density of values given by @var{d}.
28 ## @var{d} should be between 0 and 1. Values will be uniformly
29 ## distributed between 0 and 1.
30 ## 27 ##
31 ## If called with a single matrix argument, a random sparse matrix is 28 ## The size of the matrix is @var{m}x@var{n} with a density of values @var{d}.
32 ## generated wherever the matrix @var{S} is non-zero. 29 ## @var{d} must be between 0 and 1. Values will be uniformly distributed on
30 ## the interval (0, 1).
33 ## 31 ##
34 ## If called with the rc parameter as a scalar, a random sparse matrix with a 32 ## If called with a single matrix argument, a sparse matrix is generated with
35 ## reciprocal condition number rc is generated. If rc is a vector, then it 33 ## random values wherever the matrix @var{s} is non-zero.
36 ## contains the first singular values of the matrix generated (length(rc) <= min(m, n)).
37 ## 34 ##
38 ## @seealso{sprandn, sprandsym} 35 ## If called with a scalar fourth argument @var{rc}, a random sparse matrix
36 ## with reciprocal condition number @var{rc} is generated. If @var{rc} is
37 ## a vector, then it specifies the first singular values of the generated
38 ## matrix (@code{length (@var{rc}) <= min (@var{m}, @var{n})}).
39 ##
40 ## @seealso{sprandn, sprandsym, rand}
39 ## @end deftypefn 41 ## @end deftypefn
40 42
41 ## Author: Paul Kienzle <pkienzle@users.sf.net> 43 ## Author: Paul Kienzle <pkienzle@users.sf.net>
42 ## 44 ##
43 ## Changelog: 45 ## Changelog:
46 ## 2004-09-27 use Paul's hint to allow larger random matrices 48 ## 2004-09-27 use Paul's hint to allow larger random matrices
47 ## at the price of sometimes lower density than desired 49 ## at the price of sometimes lower density than desired
48 ## David Bateman 50 ## David Bateman
49 ## 2004-10-20 Texinfo help and copyright message 51 ## 2004-10-20 Texinfo help and copyright message
50 52
51 function S = sprand (m, n, d, rc) 53 function s = sprand (m, n, d, rc)
52 54
53 if (nargin == 1 ) 55 if (nargin == 1 )
54 S = __sprand_impl__ (m, @rand); 56 s = __sprand_impl__ (m, @rand);
55 elseif ( nargin == 3) 57 elseif ( nargin == 3)
56 S = __sprand_impl__ (m, n, d, "sprand", @rand); 58 s = __sprand_impl__ (m, n, d, "sprand", @rand);
57 elseif (nargin == 4) 59 elseif (nargin == 4)
58 S = __sprand_impl__ (m, n, d, rc, "sprand", @rand); 60 s = __sprand_impl__ (m, n, d, rc, "sprand", @rand);
59 else 61 else
60 print_usage (); 62 print_usage ();
61 endif 63 endif
62 64
63 endfunction 65 endfunction
64 66
65 67
68 %% Test 3-input calling form
66 %!test 69 %!test
67 %! s = sprand (4, 10, 0.1); 70 %! s = sprand (4, 10, 0.1);
68 %! assert (size (s), [4, 10]); 71 %! assert (size (s), [4, 10]);
69 %! assert (nnz (s) / numel (s), 0.1); 72 %! assert (nnz (s) / numel (s), 0.1);
70 73
71 %% Test 1-input calling form
72 %!test
73 %! s = sprand (sparse ([1 2 3], [3 2 3], [2 2 2]));
74 %! [i, j, v] = find (s);
75 %! assert (sort (i), [1 2 3]');
76 %! assert (sort (j), [2 3 3]');
77 %! assert (all (v > 0 & v < 1));
78
79 %% Test 4-input calling form 74 %% Test 4-input calling form
80 %!test 75 %!test
81 %! d = rand(); 76 %! d = rand ();
82 %! s1 = sprand (100, 100, d, 0.4); 77 %! s1 = sprand (100, 100, d, 0.4);
83 %! rc = [5, 4, 3, 2, 1, 0.1]; 78 %! rc = [5, 4, 3, 2, 1, 0.1];
84 %! s2 = sprand (100, 100, d, rc); 79 %! s2 = sprand (100, 100, d, rc);
85 %! s3 = sprand (6, 4, d, rc); 80 %! s3 = sprand (6, 4, d, rc);
86 %! assert (svd (s2)'(1:length (rc)), rc, sqrt (eps)); 81 %! assert (svd (s2)'(1:length (rc)), rc, sqrt (eps));
87 %! assert (1/cond (s1), 0.4, sqrt (eps)); 82 %! assert (1/cond (s1), 0.4, sqrt (eps));
88 %! assert (nnz (s1) / (100*100), d, 0.02); 83 %! assert (nnz (s1) / (100*100), d, 0.02);
89 %! assert (nnz (s2) / (100*100), d, 0.02); 84 %! assert (nnz (s2) / (100*100), d, 0.02);
90 %! assert (svd (s3)', [5 4 3 2], sqrt (eps)); 85 %! assert (svd (s3)', [5 4 3 2], sqrt (eps));
91 86
87 %% Test 1-input calling form
88 %!test
89 %! s = sprand (sparse ([1 2 3], [3 2 3], [2 2 2]));
90 %! [i, j, v] = find (s);
91 %! assert (sort (i), [1 2 3]');
92 %! assert (sort (j), [2 3 3]');
93 %! assert (all (v > 0 & v < 1));
94
92 %% Test input validation 95 %% Test input validation
93 %!error sprand () 96 %!error sprand ()
94 %!error sprand (1, 2) 97 %!error sprand (1, 2)
95 %!error sprand (1, 2, 3, 4) 98 %!error sprand (1, 2, 3, 4)
96 %!error sprand (ones (3), 3, 0.5) 99 %!error <M must be an integer greater than 0> sprand (ones (3), 3, 0.5)
97 %!error sprand (3.5, 3, 0.5) 100 %!error <M must be an integer greater than 0> sprand (3.5, 3, 0.5)
98 %!error sprand (0, 3, 0.5) 101 %!error <M must be an integer greater than 0> sprand (0, 3, 0.5)
99 %!error sprand (3, ones (3), 0.5) 102 %!error <N must be an integer greater than 0> sprand (3, ones (3), 0.5)
100 %!error sprand (3, 3.5, 0.5) 103 %!error <N must be an integer greater than 0> sprand (3, 3.5, 0.5)
101 %!error sprand (3, 0, 0.5) 104 %!error <N must be an integer greater than 0> sprand (3, 0, 0.5)
102 %!error sprand (3, 3, -1) 105 %!error <D must be between 0 and 1> sprand (3, 3, -1)
103 %!error sprand (3, 3, 2) 106 %!error <D must be between 0 and 1> sprand (3, 3, 2)
104 %!error sprand (2, 2, 0.2, -1) 107 %!error <RC must be a scalar or vector> sprand (2, 2, 0.2, ones (3,3))
105 %!error sprand (2, 2, 0.2, 2) 108 %!error <RC must be between 0 and 1> sprand (2, 2, 0.2, -1)
109 %!error <RC must be between 0 and 1> sprand (2, 2, 0.2, 2)
106 110
107 %% Test very large, very low density matrix doesn't fail 111 %% Test very large, very low density matrix doesn't fail
108 %!test 112 %!test
109 %! s = sprand(1e6,1e6,1e-7); 113 %! s = sprand (1e6, 1e6, 1e-7);
110 114