Mercurial > octave-nkf
annotate scripts/geometry/delaunay.m @ 18549:43cc202335dc stable release-3-8-1
Version 3.8.1 released.
* configure.ac (OCTAVE_VERSION): Now 3.8.1.
(OCTAVE_MINOR_VERSION): Now 1.
(OCTAVE_RELEASE_DATE): Set to 2014-03-06.
Update copyright date for startup message.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Thu, 06 Mar 2014 14:40:35 -0500 |
parents | d63878346099 |
children | ba167badef9f 446c46af4b42 |
rev | line source |
---|---|
17744
d63878346099
maint: Update copyright notices for release.
John W. Eaton <jwe@octave.org>
parents:
17338
diff
changeset
|
1 ## Copyright (C) 1999-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 -*- | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
20 ## @deftypefn {Function File} {} delaunay (@var{x}, @var{y}) |
16932
c55df4e5e216
doc: Fix missing @deftypefnx entry in delaunay.m.
Rik <rik@octave.org>
parents:
16928
diff
changeset
|
21 ## @deftypefnx {Function File} {} delaunay (@var{x}) |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
22 ## @deftypefnx {Function File} {} delaunay (@dots{}, @var{options}) |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
23 ## @deftypefnx {Function File} {@var{tri} =} delaunay (@dots{}) |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
24 ## Compute the Delaunay triangulation for a 2-D set of points. |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
25 ## The return value @var{tri} is a set of triangles which satisfies the |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
26 ## Delaunay circum-circle criterion, i.e., only a single data point from |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
27 ## [@var{x}, @var{y}] is within the circum-circle of the defining triangle. |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
28 ## The input @var{x} may also be a matrix with two columns where the first |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
29 ## column contains x-data and the second y-data. |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
30 ## |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
31 ## The set of triangles @var{tri} is a matrix of size [n, 3]. Each |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
32 ## row defines a triangle and the three columns are the three vertices |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
33 ## of the triangle. The value of @code{@var{tri}(i,j)} is an index into |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
34 ## @var{x} and @var{y} for the location of the j-th vertex of the i-th |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
35 ## triangle. |
6823 | 36 ## |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
37 ## The optional last argument, which must be a string or cell array of strings, |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
38 ## 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:
11587
diff
changeset
|
39 ## 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:
11587
diff
changeset
|
40 ## @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:
11587
diff
changeset
|
41 ## The default options are @code{@{"Qt", "Qbb", "Qc", "Qz"@}}. |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
42 ## |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
43 ## If @var{options} is not present or @code{[]} then the default arguments are |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
44 ## used. Otherwise, @var{options} replaces the default argument list. |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
45 ## To append user options to the defaults it is necessary to repeat the |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
46 ## default arguments in @var{options}. Use a null string to pass no arguments. |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
47 ## |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
48 ## If no output argument is specified the resulting Delaunay triangulation |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
49 ## is plotted along with the original set of points. |
6823 | 50 ## |
51 ## @example | |
52 ## @group | |
6826 | 53 ## x = rand (1, 10); |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
54 ## y = rand (1, 10); |
6826 | 55 ## T = delaunay (x, y); |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
56 ## VX = [ x(T(:,1)); x(T(:,2)); x(T(:,3)); x(T(:,1)) ]; |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
57 ## VY = [ y(T(:,1)); y(T(:,2)); y(T(:,3)); y(T(:,1)) ]; |
6826 | 58 ## axis ([0,1,0,1]); |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
59 ## plot (VX, VY, "b", x, y, "r*"); |
6823 | 60 ## @end group |
61 ## @end example | |
14541
759944521fd6
Improve tetramesh docstring and add function to manual.
Rik <octave@nomad.inbox5.com>
parents:
14237
diff
changeset
|
62 ## @seealso{delaunay3, delaunayn, convhull, voronoi, triplot, trimesh, trisurf} |
6823 | 63 ## @end deftypefn |
64 | |
6826 | 65 ## Author: Kai Habel <kai.habel@gmx.de> |
6823 | 66 |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
67 function tri = delaunay (varargin) |
6823 | 68 |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
69 if (nargin < 1 || nargin > 3) |
6823 | 70 print_usage (); |
71 endif | |
11587
c792872f8942
all script files: untabify and strip trailing whitespace
John W. Eaton <jwe@octave.org>
parents:
11523
diff
changeset
|
72 |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
73 options = []; |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
74 |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
75 switch (nargin) |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
76 |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
77 case 1 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
78 if (! ismatrix (varargin{1}) || columns (varargin{1}) != 2) |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
79 error ("delaunay: X must be a matrix with 2 columns"); |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
80 else |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
81 x = varargin{1}(:,1); |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
82 y = varargin{1}(:,2); |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
83 endif |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
84 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
85 case 2 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
86 if (isnumeric (varargin{2})) |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
87 x = varargin{1}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
88 y = varargin{2}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
89 elseif (ischar (varargin{2}) || iscellstr (varargin{2})) |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
90 options = varargin{2}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
91 if (! ismatrix (varargin{1}) && columns (varargin{1}) != 2) |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
92 error ("delaunay: X must be a matrix with 2 columns"); |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
93 else |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
94 x = varargin{1}(:,1); |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
95 y = varargin{1}(:,2); |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
96 endif |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
97 else |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
98 error ("delaunay: OPTIONS must be a string or cell array of strings"); |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
99 endif |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
100 |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
101 case 3 |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
102 x = varargin{1}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
103 y = varargin{2}; |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
104 options = varargin{3}; |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
105 |
17174
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
106 if (! (ischar (options) || iscellstr (options))) |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
107 error ("delaunay: OPTIONS must be a string or cell array of strings"); |
c3c1ebfaa7dc
maint: Use common indentation for switch statement.
Rik <rik@octave.org>
parents:
17132
diff
changeset
|
108 endif |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
109 |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
110 endswitch |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
111 |
17132
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
112 if (! (isequal(size(x),size(y)))) |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
113 error ("delaunay: X and Y must be the same size"); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
114 endif |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
115 |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
116 T = delaunayn ([x(:), y(:)], options); |
6823 | 117 |
6826 | 118 if (nargout == 0) |
119 x = x(:).'; | |
120 y = y(:).'; | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
121 VX = [ x(T(:,1)); x(T(:,2)); x(T(:,3)); x(T(:,1)) ]; |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
122 VY = [ y(T(:,1)); y(T(:,2)); y(T(:,3)); y(T(:,1)) ]; |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
123 plot (VX, VY, "b", x, y, "r*"); |
6823 | 124 else |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
125 tri = T; |
6823 | 126 endif |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
127 |
6823 | 128 endfunction |
129 | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
130 |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
131 %!demo |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
132 %! old_state = rand ("state"); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
133 %! restore_state = onCleanup (@() rand ("state", old_state)); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
134 %! rand ("state", 1); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
135 %! x = rand (1,10); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
136 %! y = rand (1,10); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
137 %! T = delaunay (x,y); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
138 %! VX = [ x(T(:,1)); x(T(:,2)); x(T(:,3)); x(T(:,1)) ]; |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
139 %! VY = [ y(T(:,1)); y(T(:,2)); y(T(:,3)); y(T(:,1)) ]; |
14237
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
140 %! clf; |
11949c9795a0
Revamp %!demos in m-files to use Octave coding conventions on spacing, etc.
Rik <octave@nomad.inbox5.com>
parents:
14138
diff
changeset
|
141 %! plot (VX,VY,"b", x,y,"r*"); |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
142 %! axis ([0,1,0,1]); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
143 |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
144 %!testif HAVE_QHULL |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
145 %! x = [-1, 0, 1, 0]; |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
146 %! y = [0, 1, 0, -1]; |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
147 %! assert (sortrows (sort (delaunay (x, y), 2)), [1,2,4;2,3,4]); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
148 |
8153
ec0a13863eb7
Only run tests that depend on HDF5 and QHull if Octave was actually
Soren Hauberg <hauberg@gmail.com>
parents:
7017
diff
changeset
|
149 %!testif HAVE_QHULL |
16928
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
150 %! x = [-1, 0, 1, 0]; |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
151 %! y = [0, 1, 0, -1]; |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
152 %! assert (sortrows (sort (delaunay ([x(:) y(:)]), 2)), [1,2,4;2,3,4]); |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
153 |
7c4a6197e020
delaunay.m: Accept a single matrix as input.
Rik <rik@octave.org>
parents:
14541
diff
changeset
|
154 %!testif HAVE_QHULL |
6823 | 155 %! x = [-1, 0, 1, 0, 0]; |
156 %! y = [0, 1, 0, -1, 0]; | |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
157 %! assert (sortrows (sort (delaunay (x, y), 2)), [1,2,5;1,4,5;2,3,5;3,4,5]); |
6823 | 158 |
17132
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
159 %!testif HAVE_QHULL |
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
160 %! x = [-1, 0; 0, 1; 1, 0; 0, -1; 0, 0]; |
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
161 %! assert (sortrows (sort (delaunay (x), 2)), [1,2,5;1,4,5;2,3,5;3,4,5]); |
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
162 |
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
163 %!testif HAVE_QHULL |
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
164 %! x = [1 5 2; 5 6 7]; |
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
165 %! y = [5 7 8; 1 2 3]; |
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
166 %! assert (sortrows (sort (delaunay (x, y), 2)), [1,2,4;1,3,4;1,3,5;3,4,6]); |
f72ffae1fcc3
delaunay.m: Fixed matlab compatibility and input check for single matrix (bug #39644)
Andreas Weber <andreas.weber@hs-offenburg.de>
parents:
16932
diff
changeset
|
167 |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11587
diff
changeset
|
168 %% FIXME: Need input validation tests |
17338
1c89599167a6
maint: End m-files with 1 blank line.
Rik <rik@octave.org>
parents:
17174
diff
changeset
|
169 |