Mercurial > octave-antonio
comparison scripts/geometry/delaunay.m @ 19164:ba167badef9f
Deprecate delaunay3 and extend delaunay to 3-D inputs.
* NEWS: Announce deprecation of delaunay3. Announce extension of delaunay to
3-D inputs.
* scripts/deprecated/delaunay3.m: Moved from geometry/. Print warning about
deprecation when executing function for the first time.
* scripts/deprecated/module.mk: Add deprecated delaunay3.m to build system.
* scripts/geometry/module.mk: Remove from geometry directory build system
* delaunay.m: Redo docstring to mention 2-D and 3-D inputs.
Overhaul function to accept 3-D inputs. Add %!error input validation tests.
* delaunayn, tetramesh.m: Remove seealso reference to delaunay3.
* geometry.txi, install.txi: Remove delaunay3 from manual.
author | Rik <rik@octave.org> |
---|---|
date | Fri, 26 Sep 2014 09:02:53 -0700 |
parents | d63878346099 |
children | 0e1f5a750d00 |
comparison
equal
deleted
inserted
replaced
19163:71da5cce37d6 | 19164:ba167badef9f |
---|---|
15 ## You should have received a copy of the GNU General Public License | 15 ## You should have received a copy of the GNU General Public License |
16 ## along with Octave; see the file COPYING. If not, see | 16 ## along with Octave; see the file COPYING. If not, see |
17 ## <http://www.gnu.org/licenses/>. | 17 ## <http://www.gnu.org/licenses/>. |
18 | 18 |
19 ## -*- texinfo -*- | 19 ## -*- texinfo -*- |
20 ## @deftypefn {Function File} {} delaunay (@var{x}, @var{y}) | 20 ## @deftypefn {Function File} {@var{tri} =} delaunay (@var{x}, @var{y}) |
21 ## @deftypefnx {Function File} {} delaunay (@var{x}) | 21 ## @deftypefnx {Function File} {@var{tetr} =} delaunay (@var{x}, @var{y}, @var{z}) |
22 ## @deftypefnx {Function File} {} delaunay (@dots{}, @var{options}) | 22 ## @deftypefnx {Function File} {@var{tri} =} delaunay (@var{x}) |
23 ## @deftypefnx {Function File} {@var{tri} =} delaunay (@dots{}) | 23 ## @deftypefnx {Function File} {@var{tri} =} delaunay (@dots{}, @var{options}) |
24 ## Compute the Delaunay triangulation for a 2-D set of points. | 24 ## Compute the Delaunay triangulation for a 2-D or 3-D set of points. |
25 ## The return value @var{tri} is a set of triangles which satisfies the | 25 ## |
26 ## Delaunay circum-circle criterion, i.e., only a single data point from | 26 ## For 2-D sets, the return value @var{tri} is a set of triangles which |
27 ## [@var{x}, @var{y}] is within the circum-circle of the defining triangle. | 27 ## satisfies the Delaunay circum-circle criterion, i.e., only a single data |
28 ## The input @var{x} may also be a matrix with two columns where the first | 28 ## point from [@var{x}, @var{y}] is within the circum-circle of the defining |
29 ## column contains x-data and the second y-data. | 29 ## triangle. The set of triangles @var{tri} is a matrix of size [n, 3]. Each |
30 ## | 30 ## row defines a triangle and the three columns are the three vertices of the |
31 ## The set of triangles @var{tri} is a matrix of size [n, 3]. Each | 31 ## triangle. The value of @code{@var{tri}(i,j)} is an index into @var{x} and |
32 ## row defines a triangle and the three columns are the three vertices | 32 ## @var{y} for the location of the j-th vertex of the i-th triangle. |
33 ## of the triangle. The value of @code{@var{tri}(i,j)} is an index into | 33 ## |
34 ## @var{x} and @var{y} for the location of the j-th vertex of the i-th | 34 ## For 3-D sets, the return value @var{tetr} is a set of tetrahedrons which |
35 ## triangle. | 35 ## satisfies the Delaunay circum-circle criterion, i.e., only a single data |
36 ## point from [@var{x}, @var{y}, @var{z}] is within the circum-circle of the | |
37 ## defining tetrahedron. The set of tetrahedrons is a matrix of size [n, 4]. | |
38 ## Each row defines a tetrahedron and the four columns are the four vertices of | |
39 ## the tetrahedron. The value of @code{@var{tetr}(i,j)} is an index into | |
40 ## @var{x}, @var{y}, @var{z} for the location of the j-th vertex of the i-th | |
41 ## tetrahedron. | |
42 ## | |
43 ## The input @var{x} may also be a matrix with two or three columns where the | |
44 ## first column contains x-data, the second y-data, and the optional third | |
45 ## column contains z-data. | |
36 ## | 46 ## |
37 ## The optional last argument, which must be a string or cell array of strings, | 47 ## The optional last argument, which must be a string or cell array of strings, |
38 ## contains options passed to the underlying qhull command. | 48 ## contains options passed to the underlying qhull command. |
39 ## See the documentation for the Qhull library for details | 49 ## See the documentation for the Qhull library for details |
40 ## @url{http://www.qhull.org/html/qh-quick.htm#options}. | 50 ## @url{http://www.qhull.org/html/qh-quick.htm#options}. |
43 ## If @var{options} is not present or @code{[]} then the default arguments are | 53 ## If @var{options} is not present or @code{[]} then the default arguments are |
44 ## used. Otherwise, @var{options} replaces the default argument list. | 54 ## used. Otherwise, @var{options} replaces the default argument list. |
45 ## To append user options to the defaults it is necessary to repeat the | 55 ## To append user options to the defaults it is necessary to repeat the |
46 ## default arguments in @var{options}. Use a null string to pass no arguments. | 56 ## default arguments in @var{options}. Use a null string to pass no arguments. |
47 ## | 57 ## |
48 ## If no output argument is specified the resulting Delaunay triangulation | |
49 ## is plotted along with the original set of points. | |
50 ## | |
51 ## @example | 58 ## @example |
52 ## @group | 59 ## @group |
53 ## x = rand (1, 10); | 60 ## x = rand (1, 10); |
54 ## y = rand (1, 10); | 61 ## y = rand (1, 10); |
55 ## T = delaunay (x, y); | 62 ## tri = delaunay (x, y); |
56 ## VX = [ x(T(:,1)); x(T(:,2)); x(T(:,3)); x(T(:,1)) ]; | 63 ## triplot (tri, x, y); |
57 ## VY = [ y(T(:,1)); y(T(:,2)); y(T(:,3)); y(T(:,1)) ]; | 64 ## hold on; |
65 ## plot (x, y, "r*"); | |
58 ## axis ([0,1,0,1]); | 66 ## axis ([0,1,0,1]); |
59 ## plot (VX, VY, "b", x, y, "r*"); | |
60 ## @end group | 67 ## @end group |
61 ## @end example | 68 ## @end example |
62 ## @seealso{delaunay3, delaunayn, convhull, voronoi, triplot, trimesh, trisurf} | 69 ## @seealso{delaunayn, convhull, voronoi, triplot, trimesh, tetramesh, trisurf} |
63 ## @end deftypefn | 70 ## @end deftypefn |
64 | 71 |
65 ## Author: Kai Habel <kai.habel@gmx.de> | 72 ## Author: Kai Habel <kai.habel@gmx.de> |
66 | 73 |
67 function tri = delaunay (varargin) | 74 function tri = delaunay (varargin) |
68 | 75 |
69 if (nargin < 1 || nargin > 3) | 76 if (nargin < 1 || nargin > 4) |
70 print_usage (); | 77 print_usage (); |
71 endif | 78 endif |
72 | 79 |
80 z = []; | |
73 options = []; | 81 options = []; |
74 | 82 |
75 switch (nargin) | 83 switch (nargin) |
76 | 84 |
77 case 1 | 85 case 1 |
78 if (! ismatrix (varargin{1}) || columns (varargin{1}) != 2) | 86 if (! ismatrix (varargin{1}) |
79 error ("delaunay: X must be a matrix with 2 columns"); | 87 || (columns (varargin{1}) != 2 && columns (varargin{1}) != 3)) |
88 error ("delaunay: X must be a matrix with 2 or 3 columns"); | |
80 else | 89 else |
81 x = varargin{1}(:,1); | 90 x = varargin{1}(:,1); |
82 y = varargin{1}(:,2); | 91 y = varargin{1}(:,2); |
92 if (columns (varargin{1}) == 3) | |
93 z = varargin{1}(:,3); | |
94 endif | |
83 endif | 95 endif |
84 | 96 |
85 case 2 | 97 case 2 |
86 if (isnumeric (varargin{2})) | 98 if (isnumeric (varargin{2})) |
87 x = varargin{1}; | 99 x = varargin{1}; |
88 y = varargin{2}; | 100 y = varargin{2}; |
89 elseif (ischar (varargin{2}) || iscellstr (varargin{2})) | 101 elseif (! (ischar (varargin{2}) || iscellstr (varargin{2}))) |
102 error ("delaunay: OPTIONS must be a string or cell array of strings"); | |
103 else | |
90 options = varargin{2}; | 104 options = varargin{2}; |
91 if (! ismatrix (varargin{1}) && columns (varargin{1}) != 2) | 105 ncols = columns (varargin{1}); |
92 error ("delaunay: X must be a matrix with 2 columns"); | 106 |
107 if (! ismatrix (varargin{1}) || (ncols != 2 && ncols != 3)) | |
108 error ("delaunay: X must be a matrix with 2 or 3 columns"); | |
93 else | 109 else |
94 x = varargin{1}(:,1); | 110 x = varargin{1}(:,1); |
95 y = varargin{1}(:,2); | 111 y = varargin{1}(:,2); |
112 if (ncols == 3) | |
113 z = varargin{1}(:,3); | |
114 endif | |
96 endif | 115 endif |
116 endif | |
117 | |
118 case 3 | |
119 if (isnumeric (varargin{3})) | |
120 x = varargin{1}; | |
121 y = varargin{2}; | |
122 z = varargin{3}; | |
123 elseif (! (ischar (varargin{3}) || iscellstr (varargin{3}))) | |
124 error ("delaunay: OPTIONS must be a string or cell array of strings"); | |
97 else | 125 else |
98 error ("delaunay: OPTIONS must be a string or cell array of strings"); | 126 x = varargin{1}; |
99 endif | 127 y = varargin{2}; |
100 | 128 options = varargin{3}; |
101 case 3 | 129 endif |
130 | |
131 case 4 | |
102 x = varargin{1}; | 132 x = varargin{1}; |
103 y = varargin{2}; | 133 y = varargin{2}; |
104 options = varargin{3}; | 134 z = varargin{3}; |
135 options = varargin{4}; | |
105 | 136 |
106 if (! (ischar (options) || iscellstr (options))) | 137 if (! (ischar (options) || iscellstr (options))) |
107 error ("delaunay: OPTIONS must be a string or cell array of strings"); | 138 error ("delaunay: OPTIONS must be a string or cell array of strings"); |
108 endif | 139 endif |
109 | 140 |
110 endswitch | 141 endswitch |
111 | 142 |
112 if (! (isequal(size(x),size(y)))) | 143 if (isempty (z)) |
113 error ("delaunay: X and Y must be the same size"); | 144 if (! size_equal (x, y)) |
114 endif | 145 error ("delaunay: X and Y must be the same size"); |
115 | 146 endif |
116 T = delaunayn ([x(:), y(:)], options); | 147 tri = delaunayn ([x(:), y(:)], options); |
117 | |
118 if (nargout == 0) | |
119 x = x(:).'; | |
120 y = y(:).'; | |
121 VX = [ x(T(:,1)); x(T(:,2)); x(T(:,3)); x(T(:,1)) ]; | |
122 VY = [ y(T(:,1)); y(T(:,2)); y(T(:,3)); y(T(:,1)) ]; | |
123 plot (VX, VY, "b", x, y, "r*"); | |
124 else | 148 else |
125 tri = T; | 149 if (! size_equal (x, y, z)) |
150 error ("delaunay: X, Y, and Z must be the same size"); | |
151 endif | |
152 tri = delaunayn ([x(:), y(:), z(:)], options); | |
126 endif | 153 endif |
127 | 154 |
128 endfunction | 155 endfunction |
129 | 156 |
130 | 157 |
132 %! old_state = rand ("state"); | 159 %! old_state = rand ("state"); |
133 %! restore_state = onCleanup (@() rand ("state", old_state)); | 160 %! restore_state = onCleanup (@() rand ("state", old_state)); |
134 %! rand ("state", 1); | 161 %! rand ("state", 1); |
135 %! x = rand (1,10); | 162 %! x = rand (1,10); |
136 %! y = rand (1,10); | 163 %! y = rand (1,10); |
137 %! T = delaunay (x,y); | 164 %! tri = delaunay (x,y); |
138 %! VX = [ x(T(:,1)); x(T(:,2)); x(T(:,3)); x(T(:,1)) ]; | |
139 %! VY = [ y(T(:,1)); y(T(:,2)); y(T(:,3)); y(T(:,1)) ]; | |
140 %! clf; | 165 %! clf; |
141 %! plot (VX,VY,"b", x,y,"r*"); | 166 %! triplot (tri, x, y); |
167 %! hold on; | |
168 %! plot (x, y, "r*"); | |
142 %! axis ([0,1,0,1]); | 169 %! axis ([0,1,0,1]); |
143 | 170 |
144 %!testif HAVE_QHULL | 171 %!testif HAVE_QHULL |
145 %! x = [-1, 0, 1, 0]; | 172 %! x = [-1, 0, 1, 0]; |
146 %! y = [0, 1, 0, -1]; | 173 %! y = [0, 1, 0, -1]; |
163 %!testif HAVE_QHULL | 190 %!testif HAVE_QHULL |
164 %! x = [1 5 2; 5 6 7]; | 191 %! x = [1 5 2; 5 6 7]; |
165 %! y = [5 7 8; 1 2 3]; | 192 %! y = [5 7 8; 1 2 3]; |
166 %! assert (sortrows (sort (delaunay (x, y), 2)), [1,2,4;1,3,4;1,3,5;3,4,6]); | 193 %! assert (sortrows (sort (delaunay (x, y), 2)), [1,2,4;1,3,4;1,3,5;3,4,6]); |
167 | 194 |
168 %% FIXME: Need input validation tests | 195 ## Test 3-D input |
169 | 196 %!testif HAVE_QHULL |
197 %! x = [-1, -1, 1, 0, -1]; y = [-1, 1, 1, 0, -1]; z = [0, 0, 0, 1, 1]; | |
198 %! assert (sortrows (sort (delaunay (x, y, z), 2)), [1,2,3,4;1,2,4,5]) | |
199 | |
200 %% Input validation tests | |
201 %!error delaunay () | |
202 %!error delaunay (1,2,3,4,5) | |
203 %!error <X must be a matrix with 2 or 3 columns> delaunay (ones (2,4)) | |
204 %!error <OPTIONS must be a string or cell array> delaunay (ones (2,2), struct()) | |
205 %!error <X must be a matrix with 2 or 3 columns> delaunay (ones (2,4), "") | |
206 %!error <OPTIONS must be a string or cell array> delaunay (ones (2,2), ones (2,2), struct()) | |
207 %!error <OPTIONS must be a string or cell array> delaunay (ones (2,2), ones (2,2), ones (2,2), struct()) | |
208 %!error <X and Y must be the same size> delaunay (1, [1 2]) | |
209 %!error <X, Y, and Z must be the same size> delaunay (1, [1 2], [1 2]) | |
210 |