Mercurial > hg > octave-nkf
annotate scripts/geometry/voronoi.m @ 18633:3f2a95a4b98d draft lyh-review
jit compiler: use existing int functions
author | Stefan Mahr <dac922@gmx.de> |
---|---|
date | Sun, 24 Nov 2013 22:46:32 +0100 |
parents | 5646f999245d |
children | 0e1f5a750d00 |
rev | line source |
---|---|
17744
d63878346099
maint: Update copyright notices for release.
John W. Eaton <jwe@octave.org>
parents:
17281
diff
changeset
|
1 ## Copyright (C) 2000-2013 Kai Habel |
6823 | 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. | |
6823 | 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/>. | |
6823 | 18 |
19 ## -*- texinfo -*- | |
10793
be55736a0783
Grammarcheck the documentation from m-files.
Rik <octave@nomad.inbox5.com>
parents:
10791
diff
changeset
|
20 ## @deftypefn {Function File} {} voronoi (@var{x}, @var{y}) |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
21 ## @deftypefnx {Function File} {} voronoi (@var{x}, @var{y}, @var{options}) |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
22 ## @deftypefnx {Function File} {} voronoi (@dots{}, "linespec") |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
23 ## @deftypefnx {Function File} {} voronoi (@var{hax}, @dots{}) |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
24 ## @deftypefnx {Function File} {@var{h} =} voronoi (@dots{}) |
6823 | 25 ## @deftypefnx {Function File} {[@var{vx}, @var{vy}] =} voronoi (@dots{}) |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
26 ## Plot the Voronoi diagram of points @code{(@var{x}, @var{y})}. |
10791
3140cb7a05a1
Add spellchecker scripts for Octave and run spellcheck of documentation
Rik <octave@nomad.inbox5.com>
parents:
10549
diff
changeset
|
27 ## The Voronoi facets with points at infinity are not drawn. |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
28 ## |
17281
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17125
diff
changeset
|
29 ## If @qcode{"linespec"} is given it is used to set the color and line style |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17125
diff
changeset
|
30 ## of the plot. If an axis graphics handle @var{hax} is supplied then the |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17125
diff
changeset
|
31 ## Voronoi diagram is drawn on the specified axis rather than in a new |
bc924baa2c4e
doc: Add new @qcode macro for code samples which are quoted.
Rik <rik@octave.org>
parents:
17125
diff
changeset
|
32 ## figure. |
6823 | 33 ## |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
34 ## The @var{options} argument, which must be a string or cell array of strings, |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
35 ## contains options passed to the underlying qhull command. |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
36 ## See the documentation for the Qhull library for details |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
37 ## @url{http://www.qhull.org/html/qh-quick.htm#options}. |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
38 ## |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
39 ## If a single output argument is requested then the Voronoi diagram will be |
14001
5f0bb45e615c
doc: Update documentation for functions returning a graphics handle h (Bug #34761)
Rik <octave@nomad.inbox5.com>
parents:
13879
diff
changeset
|
40 ## plotted and a graphics handle @var{h} to the plot is returned. |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
41 ## [@var{vx}, @var{vy}] = voronoi (@dots{}) returns the Voronoi vertices |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
42 ## instead of plotting the diagram. |
6823 | 43 ## |
44 ## @example | |
45 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
46 ## x = rand (10, 1); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
47 ## y = rand (size (x)); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
48 ## h = convhull (x, y); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
49 ## [vx, vy] = voronoi (x, y); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
50 ## plot (vx, vy, "-b", x, y, "o", x(h), y(h), "-g"); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
51 ## legend ("", "points", "hull"); |
6823 | 52 ## @end group |
53 ## @end example | |
54 ## | |
55 ## @seealso{voronoin, delaunay, convhull} | |
56 ## @end deftypefn | |
57 | |
58 ## Author: Kai Habel <kai.habel@gmx.de> | |
59 ## First Release: 20/08/2000 | |
60 | |
61 ## 2002-01-04 Paul Kienzle <pkienzle@users.sf.net> | |
62 ## * limit the default graph to the input points rather than the whole diagram | |
63 ## * provide example | |
64 ## * use unique(x,"rows") rather than __unique_rows__ | |
65 | |
66 ## 2003-12-14 Rafael Laboissiere <rafael@laboissiere.net> | |
67 ## Added optional fourth argument to pass options to the underlying | |
68 ## qhull command | |
69 | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
70 function [vx, vy] = voronoi (varargin) |
6823 | 71 |
72 if (nargin < 1) | |
73 print_usage (); | |
74 endif | |
75 | |
76 narg = 1; | |
18200
0ecd4618b1fc
voronoi.m: Fix input validation so it doesn't open blank figure window.
Rik <rik@octave.org>
parents:
18192
diff
changeset
|
77 hax = NaN; |
6823 | 78 if (isscalar (varargin{1}) && ishandle (varargin{1})) |
17125
b5d6314314fc
Change various plot functions to take advantage of new isaxes() function.
Rik <rik@octave.org>
parents:
15393
diff
changeset
|
79 hax = varargin{1}; |
18192
31d8e19a745d
voronoi.m: Add input validation test for 2 points (bug #40996)
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
80 if (! isaxes (hax)) |
31d8e19a745d
voronoi.m: Add input validation test for 2 points (bug #40996)
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
81 error ("voronoi: HAX argument must be an axes object"); |
6823 | 82 endif |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
83 narg++; |
6823 | 84 endif |
85 | |
86 if (nargin < 1 + narg || nargin > 3 + narg) | |
87 print_usage (); | |
88 endif | |
89 | |
90 x = varargin{narg++}; | |
91 y = varargin{narg++}; | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
92 |
6823 | 93 opts = {}; |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
94 if (narg <= nargin) |
6823 | 95 if (iscell (varargin{narg})) |
96 opts = varargin(narg++); | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
97 elseif (isnumeric (varargin{narg})) |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
98 ## Accept, but ignore, the triangulation |
6823 | 99 narg++; |
100 endif | |
101 endif | |
102 | |
103 linespec = {"b"}; | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
104 if (narg <= nargin && ischar (varargin{narg})) |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
105 linespec = varargin(narg); |
6823 | 106 endif |
107 | |
15393
b99c52303d0b
voronoi.m: Fix bug when there are multiple identical input points (bug #37270)
Rik <rik@octave.org>
parents:
14881
diff
changeset
|
108 if (length (x) != length (y)) |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
109 error ("voronoi: X and Y must be vectors of the same length"); |
18202
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
110 elseif (length (x) < 2) |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
111 error ("voronoi: minimum of 2 points needed"); |
6823 | 112 endif |
113 | |
8440
e792c736b1ac
Fix for dense grids and speed up for voronoi function
David Bateman <dbateman@free.fr>
parents:
8347
diff
changeset
|
114 ## Add box to approximate rays to infinity. For Voronoi diagrams the |
e792c736b1ac
Fix for dense grids and speed up for voronoi function
David Bateman <dbateman@free.fr>
parents:
8347
diff
changeset
|
115 ## box can (and should) be close to the points themselves. To make the |
e792c736b1ac
Fix for dense grids and speed up for voronoi function
David Bateman <dbateman@free.fr>
parents:
8347
diff
changeset
|
116 ## job of finding the exterior edges it should be at least two times the |
e792c736b1ac
Fix for dense grids and speed up for voronoi function
David Bateman <dbateman@free.fr>
parents:
8347
diff
changeset
|
117 ## delta below however |
6852 | 118 xmax = max (x(:)); |
119 xmin = min (x(:)); | |
120 ymax = max (y(:)); | |
121 ymin = min (y(:)); | |
122 xdelta = xmax - xmin; | |
123 ydelta = ymax - ymin; | |
8440
e792c736b1ac
Fix for dense grids and speed up for voronoi function
David Bateman <dbateman@free.fr>
parents:
8347
diff
changeset
|
124 scale = 2; |
6852 | 125 |
126 xbox = [xmin - scale * xdelta; xmin - scale * xdelta; ... | |
10549 | 127 xmax + scale * xdelta; xmax + scale * xdelta]; |
14865
70f86a64c412
* voronoi.m: Fix cut and paste error.
Nicholas Musolino <musolino@mit.edu>
parents:
14327
diff
changeset
|
128 ybox = [ymin - scale * ydelta; ymax + scale * ydelta; ... |
70f86a64c412
* voronoi.m: Fix cut and paste error.
Nicholas Musolino <musolino@mit.edu>
parents:
14327
diff
changeset
|
129 ymax + scale * ydelta; ymin - scale * ydelta]; |
6852 | 130 |
13879
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
131 [p, c, infi] = __voronoi__ ("voronoi", |
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
132 [[x(:) ; xbox(:)], [y(:); ybox(:)]], |
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
133 opts{:}); |
6823 | 134 |
15393
b99c52303d0b
voronoi.m: Fix bug when there are multiple identical input points (bug #37270)
Rik <rik@octave.org>
parents:
14881
diff
changeset
|
135 c = c(! infi).'; |
b99c52303d0b
voronoi.m: Fix bug when there are multiple identical input points (bug #37270)
Rik <rik@octave.org>
parents:
14881
diff
changeset
|
136 ## Delete null entries which cause problems in next cellfun function |
b99c52303d0b
voronoi.m: Fix bug when there are multiple identical input points (bug #37270)
Rik <rik@octave.org>
parents:
14881
diff
changeset
|
137 c(cellfun ("isempty", c)) = []; |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
138 edges = cell2mat (cellfun (@(x) [x ; [x(end), x(1:end-1)]], c, |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
139 "uniformoutput", false)); |
6823 | 140 |
141 ## Identify the unique edges of the Voronoi diagram | |
142 edges = sortrows (sort (edges).').'; | |
15393
b99c52303d0b
voronoi.m: Fix bug when there are multiple identical input points (bug #37270)
Rik <rik@octave.org>
parents:
14881
diff
changeset
|
143 edges = edges(:, [(edges(1, 1 :end - 1) != edges(1, 2 : end) | ... |
b99c52303d0b
voronoi.m: Fix bug when there are multiple identical input points (bug #37270)
Rik <rik@octave.org>
parents:
14881
diff
changeset
|
144 edges(2, 1 :end - 1) != edges(2, 2 : end)), true]); |
6852 | 145 |
18202
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
146 if (numel (x) > 2) |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
147 ## Eliminate the edges of the diagram representing the box |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
148 poutside = (1:rows (p)) ... |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
149 (p(:, 1) < xmin - xdelta | p(:, 1) > xmax + xdelta | ... |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
150 p(:, 2) < ymin - ydelta | p(:, 2) > ymax + ydelta); |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
151 edgeoutside = ismember (edges(1, :), poutside) & ... |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
152 ismember (edges(2, :), poutside); |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
153 edges(:, edgeoutside) = []; |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
154 else |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
155 ## look for the edge between the two given points |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
156 for edge = edges(1:2,:) |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
157 if (det ([[[1;1],p(edge,1:2)];1,x(1),y(1)]) |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
158 * det ([[[1;1],p(edge,1:2)];1,x(2),y(2)]) < 0) |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
159 edges = edge; |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
160 break; |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
161 endif |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
162 endfor |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
163 ## Use larger plot limits to make it more likely single bisector is shown. |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
164 xdelta = ydelta = max (xdelta, ydelta); |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
165 endif |
6823 | 166 |
167 ## Get points of the diagram | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
168 Vvx = reshape (p(edges, 1), size (edges)); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
169 Vvy = reshape (p(edges, 2), size (edges)); |
6823 | 170 |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
171 if (nargout < 2) |
18200
0ecd4618b1fc
voronoi.m: Fix input validation so it doesn't open blank figure window.
Rik <rik@octave.org>
parents:
18192
diff
changeset
|
172 if (isnan (hax)) |
0ecd4618b1fc
voronoi.m: Fix input validation so it doesn't open blank figure window.
Rik <rik@octave.org>
parents:
18192
diff
changeset
|
173 hax = gca (); |
0ecd4618b1fc
voronoi.m: Fix input validation so it doesn't open blank figure window.
Rik <rik@octave.org>
parents:
18192
diff
changeset
|
174 endif |
17125
b5d6314314fc
Change various plot functions to take advantage of new isaxes() function.
Rik <rik@octave.org>
parents:
15393
diff
changeset
|
175 h = plot (hax, Vvx, Vvy, linespec{:}, x, y, '+'); |
6852 | 176 lim = [xmin, xmax, ymin, ymax]; |
15393
b99c52303d0b
voronoi.m: Fix bug when there are multiple identical input points (bug #37270)
Rik <rik@octave.org>
parents:
14881
diff
changeset
|
177 axis (lim + 0.1 * [[-1, 1] * xdelta, [-1, 1] * ydelta]); |
6823 | 178 if (nargout == 1) |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
179 vx = h; |
6823 | 180 endif |
181 else | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
182 vx = Vvx; |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
183 vy = Vvy; |
6823 | 184 endif |
185 | |
186 endfunction | |
12824
819a60a05a65
codesprint: add test and demo for voronoi.m
Kai Habel <kai.habel@gmx.de>
parents:
12575
diff
changeset
|
187 |
819a60a05a65
codesprint: add test and demo for voronoi.m
Kai Habel <kai.habel@gmx.de>
parents:
12575
diff
changeset
|
188 |
819a60a05a65
codesprint: add test and demo for voronoi.m
Kai Habel <kai.habel@gmx.de>
parents:
12575
diff
changeset
|
189 %!demo |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
190 %! voronoi (rand (10,1), rand (10,1)); |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
191 |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
192 %!testif HAVE_QHULL |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
193 %! phi = linspace (-pi, 3/4*pi, 8); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
194 %! [x,y] = pol2cart (phi, 1); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
195 %! [vx,vy] = voronoi (x,y); |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
196 %! assert (vx(2,:), zeros (1, columns (vx)), eps); |
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
197 %! assert (vy(2,:), zeros (1, columns (vy)), eps); |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
198 |
18202
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
199 %!testif HAVE_QHULL |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
200 %! ## Special case of just 2 points |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
201 %! x = [0 1]; y = [1 0]; |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
202 %! [vx, vy] = voronoi (x,y); |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
203 %! assert (vx, [-0.7; 1.7], eps); |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
204 %! assert (vy, [-0.7; 1.7], eps); |
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
205 |
18192
31d8e19a745d
voronoi.m: Add input validation test for 2 points (bug #40996)
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
206 %% Input validation tests |
31d8e19a745d
voronoi.m: Add input validation test for 2 points (bug #40996)
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
207 %!error voronoi () |
31d8e19a745d
voronoi.m: Add input validation test for 2 points (bug #40996)
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
208 %!error voronoi (ones (3,1)) |
31d8e19a745d
voronoi.m: Add input validation test for 2 points (bug #40996)
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
209 %!error voronoi (ones (3,1), ones (3,1), "bogus1", "bogus2", "bogus3") |
31d8e19a745d
voronoi.m: Add input validation test for 2 points (bug #40996)
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
210 %!error <HAX argument must be an axes object> voronoi (0, ones (3,1), ones (3,1)) |
31d8e19a745d
voronoi.m: Add input validation test for 2 points (bug #40996)
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
211 %!error <X and Y must be vectors of the same length> voronoi (ones (3,1), ones (4,1)) |
18202
5646f999245d
voronoi.m: Add special handling for 2-point input (bug #40996).
Rik <rik@octave.org>
parents:
18200
diff
changeset
|
212 %!error <minimum of 2 points needed> voronoi (2.5, 3.5) |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
213 |