annotate scripts/geometry/convhull.m @ 7016:93c65f2a5668

[project @ 2007-10-12 06:40:56 by jwe]
author jwe
date Fri, 12 Oct 2007 06:41:26 +0000
parents 7911a62a300d
children a1dbe9d80eee
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
1 ## Copyright (C) 2000 Kai Habel
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
2 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
3 ## This file is part of Octave.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
4 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
5 ## Octave is free software; you can redistribute it and/or modify it
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
6 ## under the terms of the GNU General Public License as published by
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6846
diff changeset
7 ## the Free Software Foundation; either version 3 of the License, or (at
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6846
diff changeset
8 ## your option) any later version.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
9 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
10 ## Octave is distributed in the hope that it will be useful, but
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
11 ## WITHOUT ANY WARRANTY; without even the implied warranty of
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
12 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
13 ## General Public License for more details.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
14 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
15 ## You should have received a copy of the GNU General Public License
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6846
diff changeset
16 ## along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6846
diff changeset
17 ## <http://www.gnu.org/licenses/>.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
18
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
19 ## -*- texinfo -*-
6846
7911a62a300d [project @ 2007-08-30 07:25:57 by dbateman]
dbateman
parents: 6826
diff changeset
20 ## @deftypefn {Function File} {@var{H} =} convhull (@var{x}, @var{y})
7911a62a300d [project @ 2007-08-30 07:25:57 by dbateman]
dbateman
parents: 6826
diff changeset
21 ## @deftypefnx {Function File} {@var{H} =} convhull (@var{x}, @var{y}, @var{opt})
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
22 ## Returns the index vector to the points of the enclosing convex hull. The
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
23 ## data points are defined by the x and y vectors.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
24 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
25 ## A third optional argument, which must be a string, contains extra options
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
26 ## passed to the underlying qhull command. See the documentation for the
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
27 ## Qhull library for details.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
28 ##
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
29 ## @seealso{delaunay, convhulln}
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
30 ## @end deftypefn
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
31
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
32 ## Author: Kai Habel <kai.habel@gmx.de>
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
33
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
34 function H = convhull (x, y, opt)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
35
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
36 if (nargin != 2 && nargin != 3)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
37 print_usage ();
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
38 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
39
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
40 if (isvector (x) && isvector (y) && length (x) == length (y))
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
41 if (nargin == 2)
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
42 i = convhulln ([x(:), y(:)]);
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
43 elseif (ischar (opt) || iscell (opt))
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
44 i = convhulln ([x(:), y(:)], opt);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
45 else
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
46 error ("convhull: third argument must be a string or cell array of strings");
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
47 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
48 else
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
49 error ("convhull: first two input arguments must be vectors of same size");
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
50 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
51
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
52 n = rows (i);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
53 i = i'(:);
6826
8618f29520c6 [project @ 2007-08-24 16:02:07 by jwe]
jwe
parents: 6823
diff changeset
54 H = zeros (n + 1, 1);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
55
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
56 H(1) = i(1);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
57 next_i = i(2);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
58 i(2) = 0;
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
59 for k = 2:n
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
60 next_idx = find (i == next_i);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
61
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
62 if (rem (next_idx, 2) == 0)
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
63 H(k) = i(next_idx);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
64 next_i = i(next_idx - 1);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
65 i(next_idx - 1) = 0;
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
66 else
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
67 H(k) = i(next_idx);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
68 next_i = i(next_idx + 1);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
69 i(next_idx + 1) = 0;
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
70 endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
71 endfor
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
72
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
73 H(n + 1) = H(1);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
74 endfunction
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
75
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
76 %!test
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
77 %! x = -3:0.5:3;
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
78 %! y = abs (sin (x));
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
79 %! assert (convhull (x, y, {"s","Qci","Tcv","Pp"}), [1;7;13;12;11;10;4;3;2;1])
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
80
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
81 %!demo
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
82 %! x = -3:0.05:3;
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
83 %! y = abs (sin (x));
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
84 %! k = convhull (x, y);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
85 %! plot (x(k),y(k),'r-',x,y,'b+');
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
86 %! axis ([-3.05, 3.05, -0.05, 1.05]);