Mercurial > octave
annotate libinterp/dldfcn/convhulln.cc @ 30564:796f54d4ddbf stable
update Octave Project Developers copyright for the new year
In files that have the "Octave Project Developers" copyright notice,
update for 2021.
In all .txi and .texi files except gpl.txi and gpl.texi in the
doc/liboctave and doc/interpreter directories, change the copyright
to "Octave Project Developers", the same as used for other source
files. Update copyright notices for 2022 (not done since 2019). For
gpl.txi and gpl.texi, change the copyright notice to be "Free Software
Foundation, Inc." and leave the date at 2007 only because this file
only contains the text of the GPL, not anything created by the Octave
Project Developers.
Add Paul Thomas to contributors.in.
author | John W. Eaton <jwe@octave.org> |
---|---|
date | Tue, 28 Dec 2021 18:22:40 -0500 |
parents | 7d6709900da7 |
children | 83f9f8bda883 |
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 // |
30564
796f54d4ddbf
update Octave Project Developers copyright for the new year
John W. Eaton <jwe@octave.org>
parents:
29961
diff
changeset
|
3 // Copyright (C) 2000-2022 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 | 25 |
26 /* | |
27 29. July 2000 - Kai Habel: first release | |
28 2002-04-22 Paul Kienzle | |
29 * Use warning(...) function rather than writing to cerr | |
30 2006-05-01 Tom Holroyd | |
31 * add support for consistent winding in all dimensions; output is | |
32 * guaranteed to be simplicial. | |
33 */ | |
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 | 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 | 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 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
56 # if defined (NEED_QHULL_R_VERSION) |
6823 | 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 | 59 |
14339
3e4350f09a55
close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents:
14310
diff
changeset
|
60 static void |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
61 free_qhull_memory (qhT *qh) |
24260
0f2dc8d6c34d
Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents:
23219
diff
changeset
|
62 { |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
63 qh_freeqhull (qh, ! qh_ALL); |
24260
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; |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
66 qh_memfreeshort (qh, &curlong, &totlong); |
24260
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 |
29958
32c3a5805893
move DEFUN and DEFMETHOD functions inside octave namespace
John W. Eaton <jwe@octave.org>
parents:
29629
diff
changeset
|
89 OCTAVE_NAMESPACE_BEGIN |
32c3a5805893
move DEFUN and DEFMETHOD functions inside octave namespace
John W. Eaton <jwe@octave.org>
parents:
29629
diff
changeset
|
90 |
8786
dbd428efbf56
Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents:
7520
diff
changeset
|
91 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
|
92 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
|
93 @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
|
94 @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
|
95 @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
|
96 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
|
97 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21942
diff
changeset
|
98 @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
|
99 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
|
100 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21942
diff
changeset
|
101 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
|
102 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
|
103 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21942
diff
changeset
|
104 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
|
105 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
|
106 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
|
107 @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
|
108 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
|
109 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21942
diff
changeset
|
110 @itemize |
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 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
|
112 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21942
diff
changeset
|
113 @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
|
114 @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
|
115 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21942
diff
changeset
|
116 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
|
117 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
|
118 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
|
119 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
|
120 |
112b20240c87
move docstrings in C++ files out of C strings and into comments
John W. Eaton <jwe@octave.org>
parents:
21942
diff
changeset
|
121 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
|
122 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
|
123 @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
|
124 @end deftypefn */) |
6823 | 125 { |
14043
f913363318e0
handle new names and locations of qhull include files (bug #33712)
John W. Eaton <jwe@octave.org>
parents:
13879
diff
changeset
|
126 #if defined (HAVE_QHULL) |
6823 | 127 |
6959 | 128 int nargin = args.length (); |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
129 |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11553
diff
changeset
|
130 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
|
131 print_usage (); |
6823 | 132 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
133 octave_value_list retval; |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
134 |
13877
0d32a681d943
* convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
135 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
|
136 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
|
137 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
|
138 |
18077
ac74b0c4c564
avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents:
17787
diff
changeset
|
139 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
|
140 return retval; |
ac74b0c4c564
avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents:
17787
diff
changeset
|
141 |
13877
0d32a681d943
* convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
142 points = points.transpose (); |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
143 |
13879
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13877
diff
changeset
|
144 std::string options; |
13877
0d32a681d943
* convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
145 |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
146 if (dim <= 4) |
13879
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13877
diff
changeset
|
147 options = " Qt"; |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
148 else |
13879
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13877
diff
changeset
|
149 options = " Qt Qx"; |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
150 |
13879
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13877
diff
changeset
|
151 if (nargin == 2) |
6823 | 152 { |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
153 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
|
154 options = ' ' + args(1).string_value (); |
23577
80c42f4cca13
maint: Deprecate is_empty and replace with isempty.
Rik <rik@octave.org>
parents:
23575
diff
changeset
|
155 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
|
156 ; // Use default options. |
23575
e95738a119da
maint: Deprecate is_cellstr and replace with iscellstr.
Rik <rik@octave.org>
parents:
23220
diff
changeset
|
157 else if (args(1).iscellstr ()) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
158 { |
13879
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13877
diff
changeset
|
159 options = ""; |
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13877
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 Array<std::string> tmp = args(1).cellstr_value (); |
6823 | 162 |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
163 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
|
164 options += ' ' + tmp(i); |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
165 } |
6823 | 166 else |
20825
66cd994d1b79
eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents:
20790
diff
changeset
|
167 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
|
168 } |
6823 | 169 |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
170 boolT ismalloc = false; |
10257
cd550069240e
assume vsnprintf from gnulib; use sstream instead of snprintf
John W. Eaton <jwe@octave.org>
parents:
10154
diff
changeset
|
171 |
29517
78ccd8bf439c
Suppress extraneous error messages from Qhull (bug #57727).
Rik <rik@octave.org>
parents:
29359
diff
changeset
|
172 // Set the outfile pointer to stdout for status information. |
78ccd8bf439c
Suppress extraneous error messages from Qhull (bug #57727).
Rik <rik@octave.org>
parents:
29359
diff
changeset
|
173 FILE *outfile = nullptr; |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
174 FILE *errfile = stderr; |
14339
3e4350f09a55
close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents:
14310
diff
changeset
|
175 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
176 qhT context = { }; |
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
177 qhT *qh = &context; |
13879
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13877
diff
changeset
|
178 |
440d7914cf01
fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents:
13877
diff
changeset
|
179 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
|
180 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
181 int exitcode = qh_new_qhull (qh, dim, num_points, points.fortran_vec (), |
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
182 ismalloc, &cmd[0], outfile, errfile); |
29961
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
183 |
7d6709900da7
eliminate octave:: namespace tags in DEFUN and DEFMETHOD and more
John W. Eaton <jwe@octave.org>
parents:
29958
diff
changeset
|
184 unwind_action free_memory ([qh] () { free_qhull_memory (qh); }); |
24260
0f2dc8d6c34d
Eliminate possible segfaults related to not cleaning up Qhull workspace.
Rik <rik@octave.org>
parents:
23219
diff
changeset
|
185 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
186 if (exitcode) |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
187 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
|
188 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
189 bool nonsimp_seen = false; |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
190 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
191 octave_idx_type nf = qh->num_facets; |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
192 |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
193 Matrix idx (nf, dim + 1); |
6823 | 194 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
195 facetT *facet; |
6823 | 196 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
197 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
|
198 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
199 FORALLfacets |
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 j = 0; |
13877
0d32a681d943
* convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
202 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
203 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
|
204 { |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
205 nonsimp_seen = true; |
13746
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
206 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
207 if (cmd.find ("QJ") != std::string::npos) |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
208 // Should never happen with QJ. |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
209 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
|
210 } |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
211 |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
212 if (dim == 3) |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
213 { |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
214 setT *vertices = qh_facet3vertex (qh, facet); |
6823 | 215 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
216 vertexT *vertex, **vertexp; |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
217 |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
218 FOREACHvertex_ (vertices) |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
219 idx(i, j++) = 1 + qh_pointid(qh, vertex->point); |
13877
0d32a681d943
* convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
220 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
221 qh_settempfree (qh, &vertices); |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
222 } |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
223 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
224 { |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
225 if (facet->toporient ^ qh_ORIENTclock) |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
226 { |
13877
0d32a681d943
* convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
227 vertexT *vertex, **vertexp; |
0d32a681d943
* convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
228 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
229 FOREACHvertex_ (facet->vertices) |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
230 idx(i, j++) = 1 + qh_pointid(qh, vertex->point); |
11586
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11553
diff
changeset
|
231 } |
12df7854fa7c
strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents:
11553
diff
changeset
|
232 else |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
233 { |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
234 vertexT *vertex, **vertexp; |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
235 |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
236 FOREACHvertexreverse12_ (facet->vertices) |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
237 idx(i, j++) = 1 + qh_pointid(qh, vertex->point); |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
238 } |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
239 } |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
240 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
|
241 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
|
242 " 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
|
243 " vertices", i, j); |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
244 |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
245 i++; |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
246 } |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
247 |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
248 // 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
|
249 |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
250 if (! nonsimp_seen) |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
251 idx.resize (nf, dim, 0.0); |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
252 |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
253 if (nargout == 2) |
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 // 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
|
256 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
257 realT area; |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
258 realT dist; |
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 FORALLfacets |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
261 { |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
262 if (! facet->normal) |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
263 continue; |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
264 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
265 if (facet->upperdelaunay && qh->ATinfinity) |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
266 continue; |
13877
0d32a681d943
* convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents:
13746
diff
changeset
|
267 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
268 facet->f.area = area = qh_facetarea (qh, facet); |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
269 facet->isarea = True; |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
270 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
271 if (qh->DELAUNAY) |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
272 { |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
273 if (facet->upperdelaunay == qh->UPPERdelaunay) |
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
274 qh->totarea += area; |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
275 } |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
276 else |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
277 { |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
278 qh->totarea += area; |
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
279 qh_distplane (qh, qh->interior_point, facet, &dist); |
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
280 qh->totvol += -dist * area / qh->hull_dim; |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
281 } |
10154
40dfc0c99116
DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents:
9786
diff
changeset
|
282 } |
8786
dbd428efbf56
Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents:
7520
diff
changeset
|
283 |
29629
93c8df989ea0
qhull: Use reentrant libqhull_r (bug #60016).
Stefan Brüns <stefan.bruens@rwth-aachen.de>
parents:
29517
diff
changeset
|
284 retval(1) = octave_value (qh->totvol); |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
285 } |
8786
dbd428efbf56
Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents:
7520
diff
changeset
|
286 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
287 retval(0) = idx; |
6959 | 288 |
20898
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
289 return retval; |
8da80da1ac37
maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents:
20853
diff
changeset
|
290 |
6823 | 291 #else |
21691
263d18409fdf
Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents:
21547
diff
changeset
|
292 |
263d18409fdf
Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents:
21547
diff
changeset
|
293 octave_unused_parameter (args); |
263d18409fdf
Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents:
21547
diff
changeset
|
294 octave_unused_parameter (nargout); |
263d18409fdf
Eliminate unused variable warnings for conditionally compiled code.
John W. Eaton <jwe@octave.org>
parents:
21547
diff
changeset
|
295 |
21109
bd1752782e56
Use err_disabled_feature, warn_disabled_feature throughout code base.
Rik <rik@octave.org>
parents:
20940
diff
changeset
|
296 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
|
297 |
6823 | 298 #endif |
299 } | |
8786
dbd428efbf56
Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents:
7520
diff
changeset
|
300 |
dbd428efbf56
Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents:
7520
diff
changeset
|
301 /* |
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
|
302 %!testif HAVE_QHULL |
8786
dbd428efbf56
Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents:
7520
diff
changeset
|
303 %! 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 | 304 %! [h, v] = convhulln (cube, "Qt"); |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
305 %! assert (size (h), [12 3]); |
14310 | 306 %! h = sortrows (sort (h, 2), [1:3]); |
307 %! 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
|
308 %! 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
|
309 %! [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
|
310 %! assert (size (h2), size (h)); |
14310 | 311 %! h2 = sortrows (sort (h2, 2), [1:3]); |
312 %! assert (h2, h); | |
313 %! 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
|
314 |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
315 %!testif HAVE_QHULL |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
316 %! 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
|
317 %! [h, v] = convhulln (cube, "QJ"); |
17787
175b392e91fe
Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents:
17744
diff
changeset
|
318 %! 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
|
319 %! 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
|
320 %! 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
|
321 |
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
|
322 %!testif HAVE_QHULL |
8786
dbd428efbf56
Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents:
7520
diff
changeset
|
323 %! 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
|
324 %! [h, v] = convhulln (tetrahedron); |
7ff0bdc3dc4c
Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents:
11586
diff
changeset
|
325 %! 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
|
326 %! 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
|
327 %! 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
|
328 |
065bc7944335
fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents:
14339
diff
changeset
|
329 %!testif HAVE_QHULL |
22737
7abc25e6206a
maint: Clean up code base to follow Octave coding conventions.
Rik <rik@octave.org>
parents:
22407
diff
changeset
|
330 %! 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
|
331 %! h = convhulln (triangle); |
065bc7944335
fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents:
14339
diff
changeset
|
332 %! assert (size (h), [3 2]); |
8786
dbd428efbf56
Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents:
7520
diff
changeset
|
333 */ |
29958
32c3a5805893
move DEFUN and DEFMETHOD functions inside octave namespace
John W. Eaton <jwe@octave.org>
parents:
29629
diff
changeset
|
334 |
32c3a5805893
move DEFUN and DEFMETHOD functions inside octave namespace
John W. Eaton <jwe@octave.org>
parents:
29629
diff
changeset
|
335 OCTAVE_NAMESPACE_END |