Mercurial > hg > octave-avbm
annotate scripts/geometry/voronoi.m @ 14868:5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
* lin2mu.m, loadaudio.m, wavread.m, accumarray.m, bicubic.m, celldisp.m,
colon.m, cplxpair.m, dblquad.m, divergence.m, genvarname.m, gradient.m,
int2str.m, interp1.m, interp1q.m, interp2.m, interpn.m, loadobj.m, nthargout.m,
__isequal__.m, __splinen__.m, quadgk.m, quadl.m, quadv.m, rat.m, rot90.m,
rotdim.m, saveobj.m, subsindex.m, triplequad.m, delaunay3.m, griddata.m,
inpolygon.m, tsearchn.m, voronoi.m, get_first_help_sentence.m, which.m,
gray2ind.m, pink.m, dlmwrite.m, strread.m, textread.m, textscan.m, housh.m,
ishermitian.m, issymmetric.m, krylov.m, logm.m, null.m, rref.m,
compare_versions.m, copyfile.m, dump_prefs.m, edit.m, fileparts.m,
getappdata.m, isappdata.m, movefile.m, orderfields.m, parseparams.m,
__xzip__.m, rmappdata.m, setappdata.m, swapbytes.m, unpack.m, ver.m, fminbnd.m,
fminunc.m, fsolve.m, glpk.m, lsqnonneg.m, qp.m, sqp.m, configure_make.m,
copy_files.m, describe.m, get_description.m, get_forge_pkg.m, install.m,
installed_packages.m, is_architecture_dependent.m, load_package_dirs.m,
print_package_description.m, rebuild.m, repackage.m, save_order.m, shell.m,
allchild.m, ancestor.m, area.m, axes.m, axis.m, clabel.m, close.m, colorbar.m,
comet.m, comet3.m, contour.m, cylinder.m, ezmesh.m, ezsurf.m, findobj.m,
fplot.m, hist.m, isocolors.m, isonormals.m, isosurface.m, isprop.m, legend.m,
mesh.m, meshz.m, pareto.m, pcolor.m, peaks.m, plot3.m, plotmatrix.m, plotyy.m,
polar.m, print.m, __add_datasource__.m, __add_default_menu__.m,
__axes_limits__.m, __bar__.m, __clabel__.m, __contour__.m, __errcomm__.m,
__errplot__.m, __ezplot__.m, __file_filter__.m, __fltk_print__.m,
__ghostscript__.m, __gnuplot_print__.m, __go_draw_axes__.m,
__go_draw_figure__.m, __interp_cube__.m, __marching_cube__.m, __patch__.m,
__pie__.m, __plt__.m, __print_parse_opts__.m, __quiver__.m, __scatter__.m,
__stem__.m, __tight_eps_bbox__.m, __uigetdir_fltk__.m, __uigetfile_fltk__.m,
__uiputfile_fltk__.m, quiver.m, quiver3.m, rectangle.m, refreshdata.m,
ribbon.m, scatter.m, semilogy.m, shading.m, slice.m, subplot.m, surface.m,
surfl.m, surfnorm.m, text.m, uigetfile.m, uiputfile.m, whitebg.m, deconv.m,
mkpp.m, pchip.m, polyaffine.m, polyder.m, polygcd.m, polyout.m, polyval.m,
ppint.m, ppjumps.m, ppval.m, residue.m, roots.m, spline.m, splinefit.m,
addpref.m, getpref.m, setpref.m, ismember.m, setxor.m, arch_fit.m, arch_rnd.m,
arch_test.m, autoreg_matrix.m, diffpara.m, fftconv.m, filter2.m, hanning.m,
hurst.m, periodogram.m, triangle_sw.m, sinc.m, spectral_xdf.m, spencer.m,
stft.m, synthesis.m, unwrap.m, yulewalker.m, bicgstab.m, gmres.m, pcg.m, pcr.m,
__sprand_impl__.m, speye.m, spfun.m, sprandn.m, spstats.m, svds.m,
treelayout.m, treeplot.m, bessel.m, factor.m, legendre.m, perms.m, primes.m,
magic.m, toeplitz.m, corr.m, cov.m, mean.m, median.m, mode.m, qqplot.m,
quantile.m, ranks.m, zscore.m, logistic_regression_likelihood.m,
bartlett_test.m, chisquare_test_homogeneity.m, chisquare_test_independence.m,
kolmogorov_smirnov_test.m, run_test.m, u_test.m, wilcoxon_test.m, z_test.m,
z_test_2.m, bin2dec.m, dec2base.m, mat2str.m, strcat.m, strchr.m, strjust.m,
strtok.m, substr.m, untabify.m, assert.m, demo.m, example.m, fail.m, speed.m,
test.m, now.m: Use Octave coding conventions for cuddling parentheses in
scripts directory.
author | Rik <octave@nomad.inbox5.com> |
---|---|
date | Tue, 17 Jul 2012 07:08:39 -0700 |
parents | ce2b59a6d0e5 |
children | 02952657182e |
rev | line source |
---|---|
14138
72c96de7a403
maint: update copyright notices for 2012
John W. Eaton <jwe@octave.org>
parents:
14001
diff
changeset
|
1 ## Copyright (C) 2000-2012 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 ## |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
29 ## If "linespec" is given it is used to set the color and line style of the |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
30 ## plot. If an axis graphics handle @var{hax} is supplied then the Voronoi |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
31 ## diagram is drawn on the specified axis rather than in a new figure. |
6823 | 32 ## |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
33 ## 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
|
34 ## 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
|
35 ## 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
|
36 ## @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
|
37 ## |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
38 ## 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
|
39 ## 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
|
40 ## [@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
|
41 ## instead of plotting the diagram. |
6823 | 42 ## |
43 ## @example | |
44 ## @group | |
14327
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
45 ## x = rand (10, 1); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
46 ## y = rand (size (x)); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
47 ## h = convhull (x, y); |
4d917a6a858b
doc: Use Octave coding conventions in @example blocks of docstrings.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
48 ## [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
|
49 ## 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
|
50 ## legend ("", "points", "hull"); |
6823 | 51 ## @end group |
52 ## @end example | |
53 ## | |
54 ## @seealso{voronoin, delaunay, convhull} | |
55 ## @end deftypefn | |
56 | |
57 ## Author: Kai Habel <kai.habel@gmx.de> | |
58 ## First Release: 20/08/2000 | |
59 | |
60 ## 2002-01-04 Paul Kienzle <pkienzle@users.sf.net> | |
61 ## * limit the default graph to the input points rather than the whole diagram | |
62 ## * provide example | |
63 ## * use unique(x,"rows") rather than __unique_rows__ | |
64 | |
65 ## 2003-12-14 Rafael Laboissiere <rafael@laboissiere.net> | |
66 ## Added optional fourth argument to pass options to the underlying | |
67 ## qhull command | |
68 | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
69 function [vx, vy] = voronoi (varargin) |
6823 | 70 |
71 if (nargin < 1) | |
72 print_usage (); | |
73 endif | |
74 | |
75 narg = 1; | |
76 if (isscalar (varargin{1}) && ishandle (varargin{1})) | |
77 handl = varargin{1}; | |
6826 | 78 if (! strcmp (get (handl, "type"), "axes")) |
6823 | 79 error ("voronoi: expecting first argument to be an axes object"); |
80 endif | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
81 narg++; |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
82 elseif (nargout < 2) |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
83 handl = gca (); |
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 | |
108 lx = length (x); | |
109 ly = length (y); | |
110 | |
111 if (lx != ly) | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
112 error ("voronoi: X and Y must be vectors of the same length"); |
6823 | 113 endif |
114 | |
8440
e792c736b1ac
Fix for dense grids and speed up for voronoi function
David Bateman <dbateman@free.fr>
parents:
8347
diff
changeset
|
115 ## 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
|
116 ## 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
|
117 ## 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
|
118 ## delta below however |
6852 | 119 xmax = max (x(:)); |
120 xmin = min (x(:)); | |
121 ymax = max (y(:)); | |
122 ymin = min (y(:)); | |
123 xdelta = xmax - xmin; | |
124 ydelta = ymax - ymin; | |
8440
e792c736b1ac
Fix for dense grids and speed up for voronoi function
David Bateman <dbateman@free.fr>
parents:
8347
diff
changeset
|
125 scale = 2; |
6852 | 126 |
127 xbox = [xmin - scale * xdelta; xmin - scale * xdelta; ... | |
10549 | 128 xmax + scale * xdelta; xmax + scale * xdelta]; |
6852 | 129 ybox = [xmin - scale * xdelta; xmax + scale * xdelta; ... |
10549 | 130 xmax + scale * xdelta; xmin - scale * xdelta]; |
6852 | 131 |
13879
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
132 [p, c, infi] = __voronoi__ ("voronoi", |
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
133 [[x(:) ; xbox(:)], [y(:); ybox(:)]], |
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
134 opts{:}); |
6823 | 135 |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
136 idx = find (! infi); |
6823 | 137 ll = length (idx); |
8440
e792c736b1ac
Fix for dense grids and speed up for voronoi function
David Bateman <dbateman@free.fr>
parents:
8347
diff
changeset
|
138 c = c(idx).'; |
12931
cefd568ea073
Replace function handles with function names in cellfun calls for 15% speedup.
Rik <octave@nomad.inbox5.com>
parents:
12824
diff
changeset
|
139 k = sum (cellfun ("length", c)); |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
140 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
|
141 "uniformoutput", false)); |
6823 | 142 |
143 ## Identify the unique edges of the Voronoi diagram | |
144 edges = sortrows (sort (edges).').'; | |
6852 | 145 edges = edges (:, [(edges(1, 1: end - 1) != edges(1, 2 : end) | ... |
10549 | 146 edges(2, 1 :end - 1) != edges(2, 2 : end)), true]); |
6852 | 147 |
148 ## Eliminate the edges of the diagram representing the box | |
14868
5d3a684236b0
maint: Use Octave coding conventions for cuddling parentheses in scripts directory
Rik <octave@nomad.inbox5.com>
parents:
14335
diff
changeset
|
149 poutside = (1 : rows (p)) ... |
6852 | 150 (p (:, 1) < xmin - xdelta | p (:, 1) > xmax + xdelta | ... |
151 p (:, 2) < ymin - ydelta | p (:, 2) > ymax + ydelta); | |
152 edgeoutside = ismember (edges (1, :), poutside) & ... | |
153 ismember (edges (2, :), poutside); | |
154 edges (:, edgeoutside) = []; | |
6823 | 155 |
156 ## 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
|
157 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
|
158 Vvy = reshape (p(edges, 2), size (edges)); |
6823 | 159 |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
160 if (nargout < 2) |
6852 | 161 lim = [xmin, xmax, ymin, ymax]; |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
162 h = plot (handl, Vvx, Vvy, linespec{:}, x, y, '+'); |
6852 | 163 axis (lim + 0.1 * [[-1, 1] * (lim (2) - lim (1)), ... |
10549 | 164 [-1, 1] * (lim (4) - lim (3))]); |
6823 | 165 if (nargout == 1) |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
166 vx = h; |
6823 | 167 endif |
168 else | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
169 vx = Vvx; |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
170 vy = Vvy; |
6823 | 171 endif |
172 | |
173 endfunction | |
12824
819a60a05a65
codesprint: add test and demo for voronoi.m
Kai Habel <kai.habel@gmx.de>
parents:
12575
diff
changeset
|
174 |
819a60a05a65
codesprint: add test and demo for voronoi.m
Kai Habel <kai.habel@gmx.de>
parents:
12575
diff
changeset
|
175 |
819a60a05a65
codesprint: add test and demo for voronoi.m
Kai Habel <kai.habel@gmx.de>
parents:
12575
diff
changeset
|
176 %!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
|
177 %! 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
|
178 |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
179 %!testif HAVE_QHULL |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
180 %! 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
|
181 %! [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
|
182 %! [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
|
183 %! 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
|
184 %! 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
|
185 |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
186 %% FIXME: Need input validation tests |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
12931
diff
changeset
|
187 |