]> Creatis software - CreaPhase.git/blob - octave_packages/m/geometry/delaunay.m
update packages
[CreaPhase.git] / octave_packages / m / geometry / delaunay.m
1 ## Copyright (C) 1999-2012 Kai Habel
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
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 ## -*- texinfo -*-
20 ## @deftypefn  {Function File} {} delaunay (@var{x}, @var{y})
21 ## @deftypefnx {Function File} {@var{tri} =} delaunay (@var{x}, @var{y})
22 ## @deftypefnx {Function File} {@var{tri} =} delaunay (@var{x}, @var{y}, @var{options})
23 ## Compute the Delaunay triangulation for a 2-D set of points.
24 ## The return value @var{tri} is a set of triangles which satisfies the
25 ## Delaunay circum-circle criterion, i.e., only a single data point from
26 ## [@var{x}, @var{y}] is within the circum-circle of the defining triangle.
27 ##
28 ## The set of triangles @var{tri} is a matrix of size [n, 3].  Each
29 ## row defines a triangle and the three columns are the three vertices
30 ## of the triangle.  The value of @code{@var{tri}(i,j)} is an index into
31 ## @var{x} and @var{y} for the location of the j-th vertex of the i-th
32 ## triangle.
33 ##
34 ## An optional third argument, which must be a string or cell array of strings,
35 ## contains options passed to the underlying qhull command.
36 ## See the documentation for the Qhull library for details
37 ## @url{http://www.qhull.org/html/qh-quick.htm#options}.
38 ## The default options are @code{@{"Qt", "Qbb", "Qc", "Qz"@}}.
39 ##
40 ## If @var{options} is not present or @code{[]} then the default arguments are
41 ## used.  Otherwise, @var{options} replaces the default argument list. 
42 ## To append user options to the defaults it is necessary to repeat the 
43 ## default arguments in @var{options}.  Use a null string to pass no arguments.
44 ##
45 ## If no output argument is specified the resulting Delaunay triangulation 
46 ## is plotted along with the original set of points.
47 ##
48 ## @example
49 ## @group
50 ## x = rand (1, 10);
51 ## y = rand (1, 10);
52 ## T = delaunay (x, y);
53 ## VX = [ x(T(:,1)); x(T(:,2)); x(T(:,3)); x(T(:,1)) ];
54 ## VY = [ y(T(:,1)); y(T(:,2)); y(T(:,3)); y(T(:,1)) ];
55 ## axis ([0,1,0,1]);
56 ## plot (VX, VY, "b", x, y, "r*");
57 ## @end group
58 ## @end example
59 ## @seealso{delaunay3, delaunayn, convhull, voronoi}
60 ## @end deftypefn
61
62 ## Author: Kai Habel <kai.habel@gmx.de>
63
64 function tri = delaunay (x, y, options)
65
66   if (nargin != 2 && nargin != 3)
67     print_usage ();
68   endif
69
70   if (! (isvector (x) && isvector (y) && length (x) == length (y))
71       && ! size_equal (x, y))
72     error ("delaunay: X and Y must be the same size");
73   elseif (nargin == 3 && ! (ischar (options) || iscellstr (options)))
74     error ("delaunay: OPTIONS must be a string or cell array of strings");
75   endif
76
77   if (nargin == 2)
78     T = delaunayn ([x(:), y(:)]);
79   else
80     T = delaunayn ([x(:), y(:)], options);
81   endif
82
83   if (nargout == 0)
84     x = x(:).';
85     y = y(:).';
86     VX = [ x(T(:,1)); x(T(:,2)); x(T(:,3)); x(T(:,1)) ];
87     VY = [ y(T(:,1)); y(T(:,2)); y(T(:,3)); y(T(:,1)) ];
88     plot (VX, VY, "b", x, y, "r*");
89   else
90     tri = T;
91   endif
92
93 endfunction
94
95
96 %!demo
97 %! old_state = rand ("state");
98 %! restore_state = onCleanup (@() rand ("state", old_state));
99 %! rand ("state", 1);
100 %! x = rand (1,10);
101 %! y = rand (1,10);
102 %! T = delaunay (x,y);
103 %! VX = [ x(T(:,1)); x(T(:,2)); x(T(:,3)); x(T(:,1)) ];
104 %! VY = [ y(T(:,1)); y(T(:,2)); y(T(:,3)); y(T(:,1)) ];
105 %! axis ([0,1,0,1]);
106 %! plot (VX,VY,"b", x,y,"r*");
107
108 %!testif HAVE_QHULL
109 %! x = [-1, 0, 1, 0];
110 %! y = [0, 1, 0, -1];
111 %! assert (sortrows (sort (delaunay (x, y), 2)), [1,2,4;2,3,4]);
112
113 %!testif HAVE_QHULL
114 %! x = [-1, 0, 1, 0, 0];
115 %! y = [0, 1, 0, -1, 0];
116 %! assert (sortrows (sort (delaunay (x, y), 2)), [1,2,5;1,4,5;2,3,5;3,4,5]);
117
118 %% FIXME: Need input validation tests
119