annotate libinterp/dldfcn/convhulln.cc @ 29359:7854d5752dd2

maint: merge stable to default.
author John W. Eaton <jwe@octave.org>
date Wed, 10 Feb 2021 10:10:40 -0500
parents 1ac5a76ae91d 0a5b15007766
children 78ccd8bf439c
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
1 ////////////////////////////////////////////////////////////////////////
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
2 //
29358
0a5b15007766 update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents: 27923
diff changeset
3 // Copyright (C) 2000-2021 The Octave Project Developers
27923
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
4 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
5 // See the file COPYRIGHT.md in the top-level directory of this
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
6 // distribution or <https://octave.org/copyright/>.
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
7 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
8 // This file is part of Octave.
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
9 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
10 // Octave is free software: you can redistribute it and/or modify it
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
11 // under the terms of the GNU General Public License as published by
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
12 // the Free Software Foundation, either version 3 of the License, or
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
13 // (at your option) any later version.
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
14 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
15 // Octave is distributed in the hope that it will be useful, but
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
16 // WITHOUT ANY WARRANTY; without even the implied warranty of
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
17 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
18 // GNU General Public License for more details.
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
19 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
20 // You should have received a copy of the GNU General Public License
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
21 // along with Octave; see the file COPYING. If not, see
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
22 // <https://www.gnu.org/licenses/>.
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
23 //
bd51beb6205e update formatting of copyright notices
John W. Eaton <jwe@octave.org>
parents: 27919
diff changeset
24 ////////////////////////////////////////////////////////////////////////
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
25
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
26 /*
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
27 29. July 2000 - Kai Habel: first release
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
28 2002-04-22 Paul Kienzle
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
29 * Use warning(...) function rather than writing to cerr
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
30 2006-05-01 Tom Holroyd
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
31 * add support for consistent winding in all dimensions; output is
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
32 * guaranteed to be simplicial.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
33 */
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
34
21724
aba2e6293dd8 use "#if ..." consistently instead of "#ifdef" and "#ifndef"
John W. Eaton <jwe@octave.org>
parents: 21691
diff changeset
35 #if defined (HAVE_CONFIG_H)
21301
40de9f8f23a6 Use '#include "config.h"' rather than <config.h>.
Rik <rik@octave.org>
parents: 21200
diff changeset
36 # include "config.h"
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
37 #endif
9786
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
38
23024
a6a7b054e4ba Rationalize #includes in libinterp/dldfcn directory.
Rik <rik@octave.org>
parents: 22755
diff changeset
39 #include <limits>
a6a7b054e4ba Rationalize #includes in libinterp/dldfcn directory.
Rik <rik@octave.org>
parents: 22755
diff changeset
40 #include <string>
9786
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
41
23024
a6a7b054e4ba Rationalize #includes in libinterp/dldfcn directory.
Rik <rik@octave.org>
parents: 22755
diff changeset
42 #include "Array.h"
a6a7b054e4ba Rationalize #includes in libinterp/dldfcn directory.
Rik <rik@octave.org>
parents: 22755
diff changeset
43 #include "dMatrix.h"
a6a7b054e4ba Rationalize #includes in libinterp/dldfcn directory.
Rik <rik@octave.org>
parents: 22755
diff changeset
44 #include "oct-locbuf.h"
a6a7b054e4ba Rationalize #includes in libinterp/dldfcn directory.
Rik <rik@octave.org>
parents: 22755
diff changeset
45 #include "unwind-prot.h"
a6a7b054e4ba Rationalize #includes in libinterp/dldfcn directory.
Rik <rik@octave.org>
parents: 22755
diff changeset
46
9786
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
47 #include "defun-dld.h"
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
48 #include "error.h"
21126
ba0a4b4f021d include errwarn.h in files that use err_disabled_feature conditionally
John W. Eaton <jwe@octave.org>
parents: 21109
diff changeset
49 #include "errwarn.h"
23024
a6a7b054e4ba Rationalize #includes in libinterp/dldfcn directory.
Rik <rik@octave.org>
parents: 22755
diff changeset
50 #include "ov.h"
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
51
14043
f913363318e0 handle new names and locations of qhull include files (bug #33712)
John W. Eaton <jwe@octave.org>
parents: 13879
diff changeset
52 #if defined (HAVE_QHULL)
21691
263d18409fdf Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents: 21547
diff changeset
53
21200
fcac5dbbf9ed maint: Indent #ifdef blocks in libinterp.
Rik <rik@octave.org>
parents: 21126
diff changeset
54 # include "oct-qhull.h"
21691
263d18409fdf Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents: 21547
diff changeset
55
21200
fcac5dbbf9ed maint: Indent #ifdef blocks in libinterp.
Rik <rik@octave.org>
parents: 21126
diff changeset
56 # if defined (NEED_QHULL_VERSION)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
57 char qh_version[] = "convhulln.oct 2007-07-24";
21200
fcac5dbbf9ed maint: Indent #ifdef blocks in libinterp.
Rik <rik@octave.org>
parents: 21126
diff changeset
58 # endif
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
59
14339
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
60 static void
24260
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
61 free_qhull_memory ()
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
62 {
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
63 qh_freeqhull (! qh_ALL);
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
64
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
65 int curlong, totlong;
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
66 qh_memfreeshort (&curlong, &totlong);
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
67
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
68 if (curlong || totlong)
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
69 warning ("convhulln: did not free %d bytes of long memory (%d pieces)",
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
70 totlong, curlong);
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
71 }
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
72
18077
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
73 static bool
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
74 octave_qhull_dims_ok (octave_idx_type dim, octave_idx_type n, const char *who)
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
75 {
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
76 if (sizeof (octave_idx_type) > sizeof (int))
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
77 {
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
78 int maxval = std::numeric_limits<int>::max ();
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
79
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
80 if (dim > maxval || n > maxval)
20825
66cd994d1b79 eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents: 20790
diff changeset
81 error ("%s: dimension too large for Qhull", who);
18077
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
82 }
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
83
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
84 return true;
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
85 }
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
86
21691
263d18409fdf Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents: 21547
diff changeset
87 #endif
263d18409fdf Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents: 21547
diff changeset
88
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
89 DEFUN_DLD (convhulln, args, nargout,
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
90 doc: /* -*- texinfo -*-
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
91 @deftypefn {} {@var{h} =} convhulln (@var{pts})
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
92 @deftypefnx {} {@var{h} =} convhulln (@var{pts}, @var{options})
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
93 @deftypefnx {} {[@var{h}, @var{v}] =} convhulln (@dots{})
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
94 Compute the convex hull of the set of points @var{pts}.
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
95
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
96 @var{pts} is a matrix of size [n, dim] containing n points in a space of
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
97 dimension dim.
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
98
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
99 The hull @var{h} is an index vector into the set of points and specifies
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
100 which points form the enclosing hull.
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
101
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
102 An optional second argument, which must be a string or cell array of
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
103 strings, contains options passed to the underlying qhull command. See the
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
104 documentation for the Qhull library for details
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
105 @url{http://www.qhull.org/html/qh-quick.htm#options}.
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
106 The default options depend on the dimension of the input:
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
107
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
108 @itemize
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
109 @item 2D, 3D, 4D: @var{options} = @code{@{"Qt"@}}
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
110
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
111 @item 5D and higher: @var{options} = @code{@{"Qt", "Qx"@}}
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
112 @end itemize
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
113
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
114 If @var{options} is not present or @code{[]} then the default arguments are
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
115 used. Otherwise, @var{options} replaces the default argument list.
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
116 To append user options to the defaults it is necessary to repeat the
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
117 default arguments in @var{options}. Use a null string to pass no arguments.
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
118
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
119 If the second output @var{v} is requested the volume of the enclosing
21968
973db845cb43 doc: Remove stray '\n' from docstrings.
Rik <rik@octave.org>
parents: 21966
diff changeset
120 convex hull is calculated.
21966
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
121 @seealso{convhull, delaunayn, voronoin}
112b20240c87 move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents: 21942
diff changeset
122 @end deftypefn */)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
123 {
14043
f913363318e0 handle new names and locations of qhull include files (bug #33712)
John W. Eaton <jwe@octave.org>
parents: 13879
diff changeset
124 #if defined (HAVE_QHULL)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
125
6959
47f4f4e88166 [project @ 2007-10-04 20:43:32 by jwe]
jwe
parents: 6882
diff changeset
126 int nargin = args.length ();
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
127
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11553
diff changeset
128 if (nargin < 1 || nargin > 2)
20790
c2d9556d51d0 eliminate return statements after calls to print_usage
John W. Eaton <jwe@octave.org>
parents: 20163
diff changeset
129 print_usage ();
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
130
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
131 octave_value_list retval;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
132
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
133 Matrix points (args(0).matrix_value ());
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
134 const octave_idx_type dim = points.columns ();
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
135 const octave_idx_type num_points = points.rows ();
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
136
18077
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
137 if (! octave_qhull_dims_ok (dim, num_points, "convhulln"))
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
138 return retval;
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
139
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
140 points = points.transpose ();
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
141
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
142 std::string options;
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
143
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
144 if (dim <= 4)
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
145 options = " Qt";
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
146 else
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
147 options = " Qt Qx";
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
148
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
149 if (nargin == 2)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
150 {
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
151 if (args(1).is_string ())
23807
336f89b6208b Use character literals 'c' rather than string literals "c" when possible.
Rik <rik@octave.org>
parents: 23577
diff changeset
152 options = ' ' + args(1).string_value ();
23577
80c42f4cca13 maint: Deprecate is_empty and replace with isempty.
Rik <rik@octave.org>
parents: 23575
diff changeset
153 else if (args(1).isempty ())
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
154 ; // Use default options.
23575
e95738a119da maint: Deprecate is_cellstr and replace with iscellstr.
Rik <rik@octave.org>
parents: 23220
diff changeset
155 else if (args(1).iscellstr ())
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
156 {
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
157 options = "";
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
158
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
159 Array<std::string> tmp = args(1).cellstr_value ();
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
160
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
161 for (octave_idx_type i = 0; i < tmp.numel (); i++)
23807
336f89b6208b Use character literals 'c' rather than string literals "c" when possible.
Rik <rik@octave.org>
parents: 23577
diff changeset
162 options += ' ' + tmp(i);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
163 }
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
164 else
20825
66cd994d1b79 eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents: 20790
diff changeset
165 error ("convhulln: OPTIONS must be a string, cell array of strings, or empty");
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
166 }
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
167
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
168 boolT ismalloc = false;
10257
cd550069240e assume vsnprintf from gnulib; use sstream instead of snprintf
John W. Eaton <jwe@octave.org>
parents: 10154
diff changeset
169
14309
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
170 // Replace the outfile pointer with stdout for debugging information.
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
171 #if defined (OCTAVE_HAVE_WINDOWS_FILESYSTEM) && ! defined (OCTAVE_HAVE_POSIX_FILESYSTEM)
21942
aab79a1885cc limit gnulib headers to liboctave/wrappers directory
John W. Eaton <jwe@octave.org>
parents: 21743
diff changeset
172 FILE *outfile = std::fopen ("NUL", "w");
14309
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
173 #else
21942
aab79a1885cc limit gnulib headers to liboctave/wrappers directory
John W. Eaton <jwe@octave.org>
parents: 21743
diff changeset
174 FILE *outfile = std::fopen ("/dev/null", "w");
14309
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
175 #endif
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
176 FILE *errfile = stderr;
14339
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
177
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
178 if (! outfile)
20825
66cd994d1b79 eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents: 20790
diff changeset
179 error ("convhulln: unable to create temporary file for output");
14309
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
180
28851
1ac5a76ae91d use [=] capture default specification where possible
John W. Eaton <jwe@octave.org>
parents: 28850
diff changeset
181 octave::unwind_action close_outfile ([=] () { std::fclose (outfile); });
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
182
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
183 // qh_new_qhull command and points arguments are not const...
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
184
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
185 std::string cmd = "qhull" + options;
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
186
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
187 OCTAVE_LOCAL_BUFFER (char, cmd_str, cmd.length () + 1);
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
188
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
189 strcpy (cmd_str, cmd.c_str ());
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
190
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
191 int exitcode = qh_new_qhull (dim, num_points, points.fortran_vec (),
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
192 ismalloc, cmd_str, outfile, errfile);
24260
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
193
28823
26cfccfee9a0 Replace unwind_protect with more efficient constructs (bug #59192).
Rik <rik@octave.org>
parents: 27923
diff changeset
194 octave::unwind_action free_memory ([] () { free_qhull_memory (); });
24260
0f2dc8d6c34d Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents: 23219
diff changeset
195
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
196 if (exitcode)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
197 error ("convhulln: qhull failed");
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
198
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
199 bool nonsimp_seen = false;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
200
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
201 octave_idx_type nf = qh num_facets;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
202
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
203 Matrix idx (nf, dim + 1);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
204
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
205 facetT *facet;
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
206
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
207 octave_idx_type i = 0;
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
208
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
209 FORALLfacets
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
210 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
211 octave_idx_type j = 0;
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
212
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
213 if (! (nonsimp_seen || facet->simplicial || qh hull_dim == 2))
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
214 {
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
215 nonsimp_seen = true;
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
216
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
217 if (cmd.find ("QJ") != std::string::npos)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
218 // Should never happen with QJ.
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
219 error ("convhulln: qhull failed: option 'QJ' returned non-simplicial facet");
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
220 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
221
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
222 if (dim == 3)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
223 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
224 setT *vertices = qh_facet3vertex (facet);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
225
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
226 vertexT *vertex, **vertexp;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
227
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
228 FOREACHvertex_ (vertices)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
229 idx(i, j++) = 1 + qh_pointid(vertex->point);
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
230
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
231 qh_settempfree (&vertices);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
232 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
233 else
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
234 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
235 if (facet->toporient ^ qh_ORIENTclock)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
236 {
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
237 vertexT *vertex, **vertexp;
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
238
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
239 FOREACHvertex_ (facet->vertices)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
240 idx(i, j++) = 1 + qh_pointid(vertex->point);
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11553
diff changeset
241 }
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11553
diff changeset
242 else
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
243 {
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
244 vertexT *vertex, **vertexp;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
245
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
246 FOREACHvertexreverse12_ (facet->vertices)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
247 idx(i, j++) = 1 + qh_pointid(vertex->point);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
248 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
249 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
250 if (j < dim)
26164
7f6a50f73625 Silence compiler warnings about format identifier for octave_idx_type (bug #55046).
Markus Mützel <markus.muetzel@gmx.de>
parents: 25438
diff changeset
251 warning ("convhulln: facet %" OCTAVE_IDX_TYPE_FORMAT
7f6a50f73625 Silence compiler warnings about format identifier for octave_idx_type (bug #55046).
Markus Mützel <markus.muetzel@gmx.de>
parents: 25438
diff changeset
252 " only has %" OCTAVE_IDX_TYPE_FORMAT
7f6a50f73625 Silence compiler warnings about format identifier for octave_idx_type (bug #55046).
Markus Mützel <markus.muetzel@gmx.de>
parents: 25438
diff changeset
253 " vertices", i, j);
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
254
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
255 i++;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
256 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
257
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
258 // Remove extra dimension if all facets were simplicial.
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
259
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
260 if (! nonsimp_seen)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
261 idx.resize (nf, dim, 0.0);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
262
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
263 if (nargout == 2)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
264 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
265 // Calculate volume of convex hull, taken from qhull src/geom2.c.
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
266
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
267 realT area;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
268 realT dist;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
269
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
270 FORALLfacets
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
271 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
272 if (! facet->normal)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
273 continue;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
274
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
275 if (facet->upperdelaunay && qh ATinfinity)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
276 continue;
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
277
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
278 facet->f.area = area = qh_facetarea (facet);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
279 facet->isarea = True;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
280
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
281 if (qh DELAUNAY)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
282 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
283 if (facet->upperdelaunay == qh UPPERdelaunay)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
284 qh totarea += area;
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
285 }
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
286 else
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
287 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
288 qh totarea += area;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
289 qh_distplane (qh interior_point, facet, &dist);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
290 qh totvol += -dist * area/ qh hull_dim;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
291 }
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
292 }
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
293
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
294 retval(1) = octave_value (qh totvol);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
295 }
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
296
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
297 retval(0) = idx;
6959
47f4f4e88166 [project @ 2007-10-04 20:43:32 by jwe]
jwe
parents: 6882
diff changeset
298
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
299 return retval;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
300
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
301 #else
21691
263d18409fdf Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents: 21547
diff changeset
302
263d18409fdf Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents: 21547
diff changeset
303 octave_unused_parameter (args);
263d18409fdf Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents: 21547
diff changeset
304 octave_unused_parameter (nargout);
263d18409fdf Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents: 21547
diff changeset
305
21109
bd1752782e56 Use err_disabled_feature, warn_disabled_feature throughout code base.
Rik <rik@octave.org>
parents: 20940
diff changeset
306 err_disabled_feature ("convhulln", "Qhull");
21691
263d18409fdf Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents: 21547
diff changeset
307
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
308 #endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
309 }
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
310
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
311 /*
10474
b2143d97c002 Neither include qhull/qhull_a.h nor run tests unless we HAVE_QHULL.
David Grundberg <davidg@cs.umu.se>
parents: 10257
diff changeset
312 %!testif HAVE_QHULL
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
313 %! cube = [0 0 0;1 0 0;1 1 0;0 1 0;0 0 1;1 0 1;1 1 1;0 1 1];
14310
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
314 %! [h, v] = convhulln (cube, "Qt");
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
315 %! assert (size (h), [12 3]);
14310
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
316 %! h = sortrows (sort (h, 2), [1:3]);
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
317 %! assert (h, [1 2 4; 1 2 6; 1 4 8; 1 5 6; 1 5 8; 2 3 4; 2 3 7; 2 6 7; 3 4 7; 4 7 8; 5 6 7; 5 7 8]);
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
318 %! assert (v, 1, 10*eps);
22737
7abc25e6206a maint: Clean up code base to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22407
diff changeset
319 %! [h2, v2] = convhulln (cube); # Test default option = "Qt"
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14346
diff changeset
320 %! assert (size (h2), size (h));
14310
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
321 %! h2 = sortrows (sort (h2, 2), [1:3]);
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
322 %! assert (h2, h);
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
323 %! assert (v2, v, 10*eps);
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
324
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
325 %!testif HAVE_QHULL
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
326 %! cube = [0 0 0;1 0 0;1 1 0;0 1 0;0 0 1;1 0 1;1 1 1;0 1 1];
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
327 %! [h, v] = convhulln (cube, "QJ");
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
328 %! assert (size (h), [12 3]);
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
329 %! assert (sortrows (sort (h, 2), [1:3]), [1 2 4; 1 2 5; 1 4 5; 2 3 4; 2 3 6; 2 5 6; 3 4 8; 3 6 7; 3 7 8; 4 5 8; 5 6 8; 6 7 8]);
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
330 %! assert (v, 1.0, 1e6*eps);
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
331
10474
b2143d97c002 Neither include qhull/qhull_a.h nor run tests unless we HAVE_QHULL.
David Grundberg <davidg@cs.umu.se>
parents: 10257
diff changeset
332 %!testif HAVE_QHULL
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
333 %! tetrahedron = [1 1 1;-1 -1 1;-1 1 -1;1 -1 -1];
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
334 %! [h, v] = convhulln (tetrahedron);
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
335 %! h = sortrows (sort (h, 2), [1 2 3]);
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
336 %! assert (h, [1 2 3;1 2 4; 1 3 4; 2 3 4]);
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
337 %! assert (v, 8/3, 10*eps);
15885
065bc7944335 fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents: 14339
diff changeset
338
065bc7944335 fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents: 14339
diff changeset
339 %!testif HAVE_QHULL
22737
7abc25e6206a maint: Clean up code base to follow Octave coding conventions.
Rik <rik@octave.org>
parents: 22407
diff changeset
340 %! triangle = [0 0; 1 1; 1 0; 1 2];
15885
065bc7944335 fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents: 14339
diff changeset
341 %! h = convhulln (triangle);
065bc7944335 fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents: 14339
diff changeset
342 %! assert (size (h), [3 2]);
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
343 */