Mercurial > hg > octave-nkf
annotate scripts/sparse/sprandn.m @ 20788:7374a3a6d594
use new string_value method to handle value extraction errors
* urlwrite.cc: Use new string_value method.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Thu, 08 Oct 2015 17:26:40 -0400 |
parents | 26bd6008fc9c |
children |
rev | line source |
---|---|
19898
4197fc428c7d
maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents:
19794
diff
changeset
|
1 ## Copyright (C) 2004-2015 Paul Kienzle |
5164 | 2 ## |
7016 | 3 ## This file is part of Octave. |
4 ## | |
5 ## Octave is free software; you can redistribute it and/or modify it | |
6 ## under the terms of the GNU General Public License as published by | |
7 ## the Free Software Foundation; either version 3 of the License, or (at | |
8 ## your option) any later version. | |
9 ## | |
10 ## Octave is distributed in the hope that it will be useful, but | |
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of | |
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU | |
13 ## General Public License for more details. | |
14 ## | |
15 ## You should have received a copy of the GNU General Public License | |
16 ## along with Octave; see the file COPYING. If not, see | |
17 ## <http://www.gnu.org/licenses/>. | |
18 ## | |
19 ## Original version by Paul Kienzle distributed as free software in the | |
20 ## public domain. | |
5164 | 21 |
22 ## -*- texinfo -*- | |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
9245
diff
changeset
|
23 ## @deftypefn {Function File} {} sprandn (@var{m}, @var{n}, @var{d}) |
18724
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
24 ## @deftypefnx {Function File} {} sprandn (@var{m}, @var{n}, @var{d}, @var{rc}) |
5610 | 25 ## @deftypefnx {Function File} {} sprandn (@var{s}) |
18725 | 26 ## Generate a sparse matrix with normally distributed random values. |
27 ## | |
28 ## The size of the matrix is @var{m}x@var{n} with a density of values @var{d}. | |
29 ## @var{d} must be between 0 and 1. Values will be normally distributed with a | |
30 ## mean of 0 and a variance of 1. | |
5164 | 31 ## |
18725 | 32 ## If called with a single matrix argument, a sparse matrix is generated with |
19007
9ac2357f19bc
doc: Replace "non-zero" with "nonzero" to match existing usage.
Rik <rik@octave.org>
parents:
19006
diff
changeset
|
33 ## random values wherever the matrix @var{s} is nonzero. |
18725 | 34 ## |
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})}). | |
19794
db92e7e28e1f
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
19793
diff
changeset
|
39 ## |
18725 | 40 ## @seealso{sprand, sprandsym, randn} |
5164 | 41 ## @end deftypefn |
42 | |
43 ## Author: Paul Kienzle <pkienzle@users.sf.net> | |
44 | |
18725 | 45 function s = sprandn (m, n, d, rc) |
6498 | 46 |
13197
6db186dfdeaa
Refactor sprandn/sprand code, move common code to common function (bug #34352)
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
13064
diff
changeset
|
47 if (nargin == 1 ) |
20609
26bd6008fc9c
maint: Rename __sprand_impl__.m to __sprand__.m
Rik <rik@octave.org>
parents:
20038
diff
changeset
|
48 s = __sprand__ (m, @randn); |
13197
6db186dfdeaa
Refactor sprandn/sprand code, move common code to common function (bug #34352)
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
13064
diff
changeset
|
49 elseif ( nargin == 3) |
20609
26bd6008fc9c
maint: Rename __sprand_impl__.m to __sprand__.m
Rik <rik@octave.org>
parents:
20038
diff
changeset
|
50 s = __sprand__ (m, n, d, "sprandn", @randn); |
18724
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
51 elseif (nargin == 4) |
20609
26bd6008fc9c
maint: Rename __sprand_impl__.m to __sprand__.m
Rik <rik@octave.org>
parents:
20038
diff
changeset
|
52 s = __sprand__ (m, n, d, rc, "sprandn", @randn); |
13197
6db186dfdeaa
Refactor sprandn/sprand code, move common code to common function (bug #34352)
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
13064
diff
changeset
|
53 else |
6046 | 54 print_usage (); |
5164 | 55 endif |
13058
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
56 |
5164 | 57 endfunction |
13058
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
58 |
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
59 |
20038
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
60 ## Test 3-input calling form |
13058
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
61 %!test |
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
62 %! s = sprandn (4, 10, 0.1); |
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
63 %! assert (size (s), [4, 10]); |
13197
6db186dfdeaa
Refactor sprandn/sprand code, move common code to common function (bug #34352)
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
13064
diff
changeset
|
64 %! assert (nnz (s) / numel (s), 0.1); |
13058
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
65 |
20038
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
66 ## Test 4-input calling form |
18724
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
67 %!test |
18725 | 68 %! d = rand (); |
18724
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
69 %! s1 = sprandn (100, 100, d, 0.4); |
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
70 %! rc = [5, 4, 3, 2, 1, 0.1]; |
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
71 %! s2 = sprandn (100, 100, d, rc); |
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
72 %! s3 = sprandn (6, 4, d, rc); |
19794
db92e7e28e1f
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
19793
diff
changeset
|
73 %! assert (svd (s2)'(1:length (rc)), rc, sqrt (eps)); |
18724
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
74 %! assert (1/cond (s1), 0.4, sqrt (eps)); |
19794
db92e7e28e1f
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
19793
diff
changeset
|
75 %! assert (nnz (s1) / (100*100), d, 0.02); |
db92e7e28e1f
strip trailing whitespace from most source files
John W. Eaton <jwe@octave.org>
parents:
19793
diff
changeset
|
76 %! assert (nnz (s2) / (100*100), d, 0.02); |
18724
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
77 %! assert (svd (s3)', [5 4 3 2], sqrt (eps)); |
35a5e7740a6d
Added implementation for 4th argument of sprand/sprandn (bug #41839).
Eduardo Ramos (edu159) <eduradical951@gmail.com>
parents:
17744
diff
changeset
|
78 |
20038
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
79 ## Test 1-input calling form |
13058
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
80 %!test |
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
81 %! s = sprandn (sparse ([1 2 3], [3 2 3], [2 2 2])); |
13064
bae887ebea48
codesprint: Add input validation and tests for sprandsym.m
Rik <octave@nomad.inbox5.com>
parents:
13058
diff
changeset
|
82 %! [i, j] = find (s); |
13058
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
83 %! assert (sort (i), [1 2 3]'); |
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
84 %! assert (sort (j), [2 3 3]'); |
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
85 |
20038
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
86 ## Test very large, very low density matrix doesn't fail |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
87 %!test |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
88 %! s = sprandn (1e6,1e6,1e-7); |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
89 |
9fc020886ae9
maint: Clean up m-files to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
19898
diff
changeset
|
90 ## Test input validation |
13058
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
91 %!error sprandn () |
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
92 %!error sprandn (1, 2) |
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
93 %!error sprandn (1, 2, 3, 4) |
18725 | 94 %!error <M must be an integer greater than 0> sprandn (ones (3), 3, 0.5) |
95 %!error <M must be an integer greater than 0> sprandn (3.5, 3, 0.5) | |
96 %!error <M must be an integer greater than 0> sprandn (0, 3, 0.5) | |
97 %!error <N must be an integer greater than 0> sprandn (3, ones (3), 0.5) | |
98 %!error <N must be an integer greater than 0> sprandn (3, 3.5, 0.5) | |
99 %!error <N must be an integer greater than 0> sprandn (3, 0, 0.5) | |
100 %!error <D must be between 0 and 1> sprandn (3, 3, -1) | |
101 %!error <D must be between 0 and 1> sprandn (3, 3, 2) | |
102 %!error <RC must be a scalar or vector> sprandn (2, 2, 0.2, ones (3,3)) | |
103 %!error <RC must be between 0 and 1> sprandn (2, 2, 0.2, -1) | |
104 %!error <RC must be between 0 and 1> sprandn (2, 2, 0.2, 2) | |
13058
14422cc782b2
codesprint: Write input validation and tests for sprandn.m
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
105 |