Mercurial > hg > octave-nkf
annotate scripts/plot/hist.m @ 17122:eaab03308c0b
doc: Rewrite docstrings for most plot functions.
Emphasize clarity, use common "voice", and increase density of seealso links.
* doc/interpreter/plot.txi: Add @findex entries that were in xlim.m
* scripts/miscellaneous/getappdata.m scripts/miscellaneous/setappdata.m,
scripts/plot/allchild.m, scripts/plot/ancestor.m, scripts/plot/area.m,
scripts/plot/axes.m, scripts/plot/axis.m, scripts/plot/bar.m,
scripts/plot/barh.m, scripts/plot/box.m, scripts/plot/caxis.m,
scripts/plot/cla.m, scripts/plot/clabel.m, scripts/plot/clf.m,
scripts/plot/close.m, scripts/plot/closereq.m, scripts/plot/colorbar.m,
scripts/plot/comet.m, scripts/plot/comet3.m, scripts/plot/compass.m,
scripts/plot/contour.m, scripts/plot/contour3.m, scripts/plot/contourc.m,
scripts/plot/contourf.m, scripts/plot/copyobj.m, scripts/plot/cylinder.m,
scripts/plot/daspect.m, scripts/plot/diffuse.m, scripts/plot/ellipsoid.m,
scripts/plot/errorbar.m, scripts/plot/ezcontour.m, scripts/plot/ezcontourf.m,
scripts/plot/ezmesh.m, scripts/plot/ezmeshc.m, scripts/plot/ezplot.m,
scripts/plot/ezplot3.m, scripts/plot/ezpolar.m, scripts/plot/ezsurf.m,
scripts/plot/ezsurfc.m, scripts/plot/feather.m, scripts/plot/figure.m,
scripts/plot/fill.m, scripts/plot/findall.m, scripts/plot/findobj.m,
scripts/plot/fplot.m, scripts/plot/gca.m, scripts/plot/gcbf.m,
scripts/plot/gcbo.m, scripts/plot/gcf.m, scripts/plot/gco.m,
scripts/plot/ginput.m, scripts/plot/graphics_toolkit.m, scripts/plot/grid.m,
scripts/plot/gtext.m, scripts/plot/guidata.m, scripts/plot/guihandles.m,
scripts/plot/hdl2struct.m, scripts/plot/hggroup.m, scripts/plot/hidden.m,
scripts/plot/hist.m, scripts/plot/hold.m, scripts/plot/ishghandle.m,
scripts/plot/ishold.m, scripts/plot/isocolors.m, scripts/plot/isprop.m,
scripts/plot/legend.m, scripts/plot/line.m, scripts/plot/linkprop.m,
scripts/plot/loglog.m, scripts/plot/loglogerr.m, scripts/plot/mesh.m,
scripts/plot/meshc.m, scripts/plot/meshgrid.m, scripts/plot/meshz.m,
scripts/plot/newplot.m, scripts/plot/orient.m, scripts/plot/pareto.m,
scripts/plot/patch.m, scripts/plot/pcolor.m, scripts/plot/peaks.m,
scripts/plot/pie.m, scripts/plot/pie3.m, scripts/plot/plot.m,
scripts/plot/plot3.m, scripts/plot/plotmatrix.m, scripts/plot/plotyy.m,
scripts/plot/polar.m, scripts/plot/print.m, scripts/plot/quiver.m,
scripts/plot/quiver3.m, scripts/plot/rectangle.m, scripts/plot/refresh.m,
scripts/plot/refreshdata.m, scripts/plot/ribbon.m, scripts/plot/rose.m,
scripts/plot/saveas.m, scripts/plot/scatter.m, scripts/plot/scatter3.m,
scripts/plot/semilogx.m, scripts/plot/semilogxerr.m, scripts/plot/semilogy.m,
scripts/plot/semilogyerr.m, scripts/plot/shading.m, scripts/plot/shg.m,
scripts/plot/shrinkfaces.m, scripts/plot/slice.m, scripts/plot/specular.m,
scripts/plot/sphere.m, scripts/plot/stairs.m, scripts/plot/stem.m,
scripts/plot/stem3.m, scripts/plot/struct2hdl.m, scripts/plot/subplot.m,
scripts/plot/surf.m, scripts/plot/surface.m, scripts/plot/surfc.m,
scripts/plot/surfl.m, scripts/plot/tetramesh.m, scripts/plot/text.m,
scripts/plot/title.m, scripts/plot/trimesh.m, scripts/plot/triplot.m,
scripts/plot/trisurf.m, scripts/plot/view.m, scripts/plot/waitbar.m,
scripts/plot/waitforbuttonpress.m, scripts/plot/waterfall.m,
scripts/plot/whitebg.m, scripts/plot/xlabel.m, scripts/plot/xlim.m,
scripts/plot/ylabel.m, scripts/plot/ylim.m, scripts/plot/zlabel.m,
scripts/plot/zlim.m: Rewrite docstrings for most plot functions.
Emphasize clarity, use common "voice", and increase density of seealso links.
author | Rik <rik@octave.org> |
---|---|
date | Wed, 31 Jul 2013 13:53:30 -0700 |
parents | c2dbdeaa25df |
children | b491ef539071 |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
11587
diff
changeset
|
1 ## Copyright (C) 1994-2012 John W. Eaton |
2313 | 2 ## |
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 | |
7016 | 7 ## the Free Software Foundation; either version 3 of the License, or (at |
8 ## your option) any later version. | |
2313 | 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 | |
7016 | 16 ## along with Octave; see the file COPYING. If not, see |
17 ## <http://www.gnu.org/licenses/>. | |
724 | 18 |
3368 | 19 ## -*- texinfo -*- |
11563
3c6e8aaa9555
Grammarcheck m-files before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11523
diff
changeset
|
20 ## @deftypefn {Function File} {} hist (@var{y}) |
11351
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
21 ## @deftypefnx {Function File} {} hist (@var{y}, @var{x}) |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
22 ## @deftypefnx {Function File} {} hist (@var{y}, @var{nbins}) |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
23 ## @deftypefnx {Function File} {} hist (@var{y}, @var{x}, @var{norm}) |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
24 ## @deftypefnx {Function File} {[@var{nn}, @var{xx}] =} hist (@dots{}) |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
25 ## @deftypefnx {Function File} {[@dots{}] =} hist (@dots{}, @var{prop}, @var{val}) |
2311 | 26 ## Produce histogram counts or plots. |
3426 | 27 ## |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11575
diff
changeset
|
28 ## With one vector input argument, @var{y}, plot a histogram of the values |
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11575
diff
changeset
|
29 ## with 10 bins. The range of the histogram bins is determined by the |
11351
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
30 ## range of the data. With one matrix input argument, @var{y}, plot a |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
31 ## histogram where each bin contains a bar per input column. |
3426 | 32 ## |
11351
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
33 ## Given a second vector argument, @var{x}, use that as the centers of |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11575
diff
changeset
|
34 ## the bins, with the width of the bins determined from the adjacent |
11351
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
35 ## values in the vector. |
3426 | 36 ## |
11575
d6619410e79c
Spellcheck documentation before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11563
diff
changeset
|
37 ## If scalar, the second argument, @var{nbins}, defines the number of bins. |
3426 | 38 ## |
11575
d6619410e79c
Spellcheck documentation before 3.4 release.
Rik <octave@nomad.inbox5.com>
parents:
11563
diff
changeset
|
39 ## If a third argument is provided, the histogram is normalized such that |
3597 | 40 ## the sum of the bars is equal to @var{norm}. |
41 ## | |
2311 | 42 ## Extreme values are lumped in the first and last bins. |
3426 | 43 ## |
3368 | 44 ## With two output arguments, produce the values @var{nn} and @var{xx} such |
45 ## that @code{bar (@var{xx}, @var{nn})} will plot the histogram. | |
11351
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
46 ## |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
47 ## The histogram's appearance may be modified by specifying property/value |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11575
diff
changeset
|
48 ## pairs, @var{prop} and @var{val} pairs. For example the face and edge |
11351
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
49 ## color may be modified. |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
50 ## |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
51 ## @example |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
52 ## @group |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
53 ## hist (randn (1, 100), 25, "facecolor", "r", "edgecolor", "b"); |
11351
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
54 ## @end group |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
55 ## @end example |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
56 ## |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
57 ## @noindent |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
58 ## The histograms colors also depend upon the colormap. |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
59 ## |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
60 ## @example |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
61 ## @group |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
62 ## hist (rand (10, 3)); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
63 ## colormap (summer ()); |
11351
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
64 ## @end group |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
65 ## @end example |
bff585d759cf
improved the help text in hist.m
Doug Stewart <doug.dastew@gmail.com>
parents:
10549
diff
changeset
|
66 ## |
17122
eaab03308c0b
doc: Rewrite docstrings for most plot functions.
Rik <rik@octave.org>
parents:
14872
diff
changeset
|
67 ## @seealso{histc, bar, pie, rose} |
3368 | 68 ## @end deftypefn |
724 | 69 |
2314 | 70 ## Author: jwe |
71 | |
7112 | 72 function [nn, xx] = hist (y, varargin) |
724 | 73 |
7112 | 74 if (nargin < 1) |
6046 | 75 print_usage (); |
724 | 76 endif |
2325 | 77 |
5443 | 78 arg_is_vector = isvector (y); |
5065 | 79 |
80 if (rows (y) == 1) | |
4880 | 81 y = y(:); |
82 endif | |
83 | |
84 if (isreal (y)) | |
7112 | 85 max_val = max (y(:)); |
86 min_val = min (y(:)); | |
724 | 87 else |
8427
26b899d309f6
help and error string corrected in hist.m.
Francesco Potortì <pot@gnu.org>
parents:
7566
diff
changeset
|
88 error ("hist: first argument must be real valued"); |
724 | 89 endif |
90 | |
7112 | 91 iarg = 1; |
92 if (nargin == 1 || ischar (varargin{iarg})) | |
724 | 93 n = 10; |
4880 | 94 x = [0.5:n]'/n; |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14613
diff
changeset
|
95 x = x * (max_val - min_val) + ones (size (x)) * min_val; |
3597 | 96 else |
97 ## nargin is either 2 or 3 | |
7208 | 98 x = varargin{iarg++}; |
4030 | 99 if (isscalar (x)) |
724 | 100 n = x; |
101 if (n <= 0) | |
102 error ("hist: number of bins must be positive"); | |
103 endif | |
4880 | 104 x = [0.5:n]'/n; |
7191 | 105 x = x * (max_val - min_val) + ones (size (x)) * min_val; |
4880 | 106 elseif (isreal (x)) |
107 if (isvector (x)) | |
10549 | 108 x = x(:); |
4880 | 109 endif |
724 | 110 tmp = sort (x); |
111 if (any (tmp != x)) | |
904 | 112 warning ("hist: bin values not sorted on input"); |
724 | 113 x = tmp; |
114 endif | |
115 else | |
116 error ("hist: second argument must be a scalar or a vector"); | |
117 endif | |
118 endif | |
119 | |
7189 | 120 ## Avoid issues with integer types for x and y |
121 x = double (x); | |
122 y = double (y); | |
123 | |
4880 | 124 cutoff = (x(1:end-1,:) + x(2:end,:)) / 2; |
125 n = rows (x); | |
7566
b3acdf1c41a5
hist: avoid temps; allow matrix args when number of bins > 30
John W. Eaton <jwe@octave.org>
parents:
7208
diff
changeset
|
126 y_nc = columns (y); |
4880 | 127 if (n < 30 && columns (x) == 1) |
4407 | 128 ## The following algorithm works fastest for n less than about 30. |
7566
b3acdf1c41a5
hist: avoid temps; allow matrix args when number of bins > 30
John W. Eaton <jwe@octave.org>
parents:
7208
diff
changeset
|
129 chist = zeros (n+1, y_nc); |
4407 | 130 for i = 1:n-1 |
4880 | 131 chist(i+1,:) = sum (y <= cutoff(i)); |
4407 | 132 endfor |
5746 | 133 chist(n+1,:) = sum (! isnan (y)); |
4407 | 134 else |
135 ## The following algorithm works fastest for n greater than about 30. | |
136 ## Put cutoff elements between boundaries, integrate over all | |
137 ## elements, keep totals at boundaries. | |
7566
b3acdf1c41a5
hist: avoid temps; allow matrix args when number of bins > 30
John W. Eaton <jwe@octave.org>
parents:
7208
diff
changeset
|
138 [s, idx] = sort ([y; repmat(cutoff, 1, y_nc)]); |
4880 | 139 len = rows (y); |
140 chist = cumsum (idx <= len); | |
7566
b3acdf1c41a5
hist: avoid temps; allow matrix args when number of bins > 30
John W. Eaton <jwe@octave.org>
parents:
7208
diff
changeset
|
141 chist = [(zeros (1, y_nc)); |
10549 | 142 (reshape (chist(idx > len), rows (cutoff), y_nc)); |
143 (chist(end,:) - sum (isnan (y)))]; | |
4407 | 144 endif |
145 | |
4880 | 146 freq = diff (chist); |
724 | 147 |
7191 | 148 if (nargin > 2 && ! ischar (varargin{iarg})) |
3597 | 149 ## Normalise the histogram. |
7112 | 150 norm = varargin{iarg++}; |
14613
e7c8e31f8e5d
hist.m: Add test to check for correct NaN normalising
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
14612
diff
changeset
|
151 freq = freq / sum(! isnan (y)) * norm; |
3597 | 152 endif |
153 | |
6586 | 154 if (nargout > 0) |
5065 | 155 if (arg_is_vector) |
4880 | 156 nn = freq'; |
157 xx = x'; | |
158 else | |
159 nn = freq; | |
160 xx = x; | |
161 endif | |
14872
c2dbdeaa25df
maint: use rows() and columns() to clarify m-files.
Rik <octave@nomad.inbox5.com>
parents:
14868
diff
changeset
|
162 elseif (columns (freq) != 1) |
7112 | 163 bar (x, freq, 0.8, varargin{iarg:end}); |
736 | 164 else |
7112 | 165 bar (x, freq, 1.0, varargin{iarg:end}); |
724 | 166 endif |
167 | |
168 endfunction | |
4811 | 169 |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
170 |
4811 | 171 %!test |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
172 %! [nn,xx] = hist ([1:4], 3); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
173 %! assert (xx, [1.5,2.5,3.5]); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
174 %! assert (nn, [2,1,1]); |
4880 | 175 %!test |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
176 %! [nn,xx] = hist ([1:4]', 3); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
177 %! assert (xx, [1.5,2.5,3.5]); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
178 %! assert (nn, [2,1,1]); |
4880 | 179 %!test |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
180 %! [nn,xx] = hist ([1 1 1 NaN NaN NaN 2 2 3],[1 2 3]); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
181 %! assert (xx, [1,2,3]); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
182 %! assert (nn, [3,2,1]); |
5746 | 183 %!test |
14613
e7c8e31f8e5d
hist.m: Add test to check for correct NaN normalising
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
14612
diff
changeset
|
184 %! [nn,xx] = hist ([1 1 1 NaN NaN NaN 2 2 3],[1 2 3], 6); |
e7c8e31f8e5d
hist.m: Add test to check for correct NaN normalising
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
14612
diff
changeset
|
185 %! assert (xx, [1,2,3]); |
e7c8e31f8e5d
hist.m: Add test to check for correct NaN normalising
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
14612
diff
changeset
|
186 %! assert (nn, [3,2,1]); |
e7c8e31f8e5d
hist.m: Add test to check for correct NaN normalising
Jordi Gutiérrez Hermoso <jordigh@octave.org>
parents:
14612
diff
changeset
|
187 %!test |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
188 %! [nn,xx] = hist ([[1:4]', [1:4]'], 3); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
189 %! assert (xx, [1.5;2.5;3.5]); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
190 %! assert (nn, [[2,1,1]',[2,1,1]']); |
4880 | 191 %!test |
14363
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
192 %! for n = [10, 30, 100, 1000] |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
193 %! assert (sum (hist ([1:n], n)), n); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
194 %! assert (sum (hist ([1:n], [2:n-1])), n); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
195 %! assert (sum (hist ([1:n], [1:n])), n); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
196 %! assert (sum (hist ([1:n], 29)), n); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
197 %! assert (sum (hist ([1:n], 30)), n); |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
198 %! endfor |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
199 %!assert (hist (1,1), 1) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
200 %!assert (size (hist (randn (750,240), 200)), [200,240]) |
f3d52523cde1
Use Octave coding conventions in all m-file %!test blocks
Rik <octave@nomad.inbox5.com>
parents:
14327
diff
changeset
|
201 |