annotate libinterp/dldfcn/convhulln.cc @ 21109:bd1752782e56

Use err_disabled_feature, warn_disabled_feature throughout code base. In liboctave, use the same error text as err_disabled_feature, but call error_handler directly because err_disabled_feature is in libinterp. * errwarn.cc (err_disabled_feature): Don't print leading "%s:" if calling function is "". * errwarn.cc (warn_disabled_feature): New function. Same msg as err_disabled_feature but uses warning rather than error. * errwarn.h (warn_disabled_feature): prototype for new function. * file-io.cc, gl-render.cc, gl2ps-renderer.cc, load-save.cc, ls-mat5.cc, oct-hdf5-types.cc, pt-jit.cc, syscalls.cc, sysdep.cc, toplev.cc, __delaunayn__.cc, __eigs__.cc, __fltk_uigetfile__.cc, __init_fltk__.cc, __osmesa_print__.cc, __voronoi__.cc, amd.cc, ccolamd.cc, colamd.cc, convhulln.cc, dmperm.cc, fftw.cc, symbfact.cc: Replace calls to error about missing feature with calls to err_disabled_feature. Replace calls to warning about missing feature with calls to warn_disabled_feature. * CSparse.cc, dSparse.cc, SparseCmplxLU.cc, SparseCmplxQR.cc, SparseQR.cc, SparsedbleLU.cc, sparse-base-chol.cc, sparse-dmsolve.cc, oct-shlib.cc: Use same text of message from err_disabled_feature but call (*current_liboctave_error_handler) or (*current_liboctave_warning_with_id_handler).
author Rik <rik@octave.org>
date Tue, 19 Jan 2016 12:44:54 -0800
parents 48b2ad5ee801
children ba0a4b4f021d
Ignore whitespace changes - Everywhere: Within whitespace: At end of lines:
rev   line source
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
1 /*
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
2
19697
4197fc428c7d maint: Update copyright notices for 2015.
John W. Eaton <jwe@octave.org>
parents: 19269
diff changeset
3 Copyright (C) 2000-2015 Kai Habel
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
4
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
5 This file is part of Octave.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
6
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
7 Octave is free software; you can redistribute it and/or modify it
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
8 under the terms of the GNU General Public License as published by the
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6959
diff changeset
9 Free Software Foundation; either version 3 of the License, or (at your
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6959
diff changeset
10 option) any later version.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
11
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
12 Octave is distributed in the hope that it will be useful, but WITHOUT
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
13 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
14 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
15 for more details.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
16
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
17 You should have received a copy of the GNU General Public License
7016
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6959
diff changeset
18 along with Octave; see the file COPYING. If not, see
93c65f2a5668 [project @ 2007-10-12 06:40:56 by jwe]
jwe
parents: 6959
diff changeset
19 <http://www.gnu.org/licenses/>.
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
20
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
21 */
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
22
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
23 /*
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
24 29. July 2000 - Kai Habel: first release
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
25 2002-04-22 Paul Kienzle
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
26 * Use warning(...) function rather than writing to cerr
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
27 2006-05-01 Tom Holroyd
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
28 * add support for consistent winding in all dimensions; output is
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
29 * guaranteed to be simplicial.
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
30 */
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
31
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
32 #ifdef HAVE_CONFIG_H
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
33 #include <config.h>
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
34 #endif
9786
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
35
19269
65554f5847ac don't include oct-locbuf.h in header files unnecessarily
John W. Eaton <jwe@octave.org>
parents: 18077
diff changeset
36 #include "oct-locbuf.h"
9786
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
37
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
38 #include "Cell.h"
9786
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
39 #include "defun-dld.h"
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
40 #include "error.h"
20940
48b2ad5ee801 maint: Rename oct-obj.[cc|h] to ovl.[cc|h] for clarity.
Rik <rik@octave.org>
parents: 20898
diff changeset
41 #include "ovl.h"
9786
2c279308f6ab fix includes in some src/DLD-FUNCTIONS files
John W. Eaton <jwe@octave.org>
parents: 8920
diff changeset
42 #include "parse.h"
14339
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
43 #include "unwind-prot.h"
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
44
14043
f913363318e0 handle new names and locations of qhull include files (bug #33712)
John W. Eaton <jwe@octave.org>
parents: 13879
diff changeset
45 #if defined (HAVE_QHULL)
f913363318e0 handle new names and locations of qhull include files (bug #33712)
John W. Eaton <jwe@octave.org>
parents: 13879
diff changeset
46 # include "oct-qhull.h"
f913363318e0 handle new names and locations of qhull include files (bug #33712)
John W. Eaton <jwe@octave.org>
parents: 13879
diff changeset
47 # if defined (NEED_QHULL_VERSION)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
48 char qh_version[] = "convhulln.oct 2007-07-24";
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
49 # endif
14043
f913363318e0 handle new names and locations of qhull include files (bug #33712)
John W. Eaton <jwe@octave.org>
parents: 13879
diff changeset
50 #endif
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
51
14339
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
52 static void
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
53 close_fcn (FILE *f)
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
54 {
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
55 gnulib::fclose (f);
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
56 }
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
57
18077
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
58 static bool
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
59 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
60 {
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
61 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
62 {
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
63 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
64
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
65 if (dim > maxval || n > maxval)
20825
66cd994d1b79 eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents: 20790
diff changeset
66 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
67 }
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
68
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
69 return true;
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
70 }
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
71
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
72 DEFUN_DLD (convhulln, args, nargout,
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
73 "-*- texinfo -*-\n\
20853
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20825
diff changeset
74 @deftypefn {} {@var{h} =} convhulln (@var{pts})\n\
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20825
diff changeset
75 @deftypefnx {} {@var{h} =} convhulln (@var{pts}, @var{options})\n\
1142cf6abc0d 2015 Code Sprint: remove class of function from docstring for all C++ files.
Rik <rik@octave.org>
parents: 20825
diff changeset
76 @deftypefnx {} {[@var{h}, @var{v}] =} convhulln (@dots{})\n\
20163
075a5e2e1ba5 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19697
diff changeset
77 Compute the convex hull of the set of points @var{pts}.\n\
075a5e2e1ba5 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19697
diff changeset
78 \n\
075a5e2e1ba5 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19697
diff changeset
79 @var{pts} is a matrix of size [n, dim] containing n points in a space of\n\
075a5e2e1ba5 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19697
diff changeset
80 dimension dim.\n\
075a5e2e1ba5 doc: Update more docstrings to have one sentence summary as first line.
Rik <rik@octave.org>
parents: 19697
diff changeset
81 \n\
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
82 The hull @var{h} is an index vector into the set of points and specifies\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
83 which points form the enclosing hull.\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
84 \n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
85 An optional second argument, which must be a string or cell array of strings,\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
86 contains options passed to the underlying qhull command.\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
87 See the documentation for the Qhull library for details\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
88 @url{http://www.qhull.org/html/qh-quick.htm#options}.\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
89 The default options depend on the dimension of the input:\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
90 \n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
91 @itemize\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
92 @item 2D, 3D, 4D: @var{options} = @code{@{\"Qt\"@}}\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
93 \n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
94 @item 5D and higher: @var{options} = @code{@{\"Qt\", \"Qx\"@}}\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
95 @end itemize\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
96 \n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
97 If @var{options} is not present or @code{[]} then the default arguments are\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
98 used. Otherwise, @var{options} replaces the default argument list.\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
99 To append user options to the defaults it is necessary to repeat the\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
100 default arguments in @var{options}. Use a null string to pass no arguments.\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
101 \n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
102 If the second output @var{v} is requested the volume of the enclosing\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
103 convex hull is calculated.\n\n\
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
104 @seealso{convhull, delaunayn, voronoin}\n\
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
105 @end deftypefn")
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
106 {
14043
f913363318e0 handle new names and locations of qhull include files (bug #33712)
John W. Eaton <jwe@octave.org>
parents: 13879
diff changeset
107 #if defined (HAVE_QHULL)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
108
6959
47f4f4e88166 [project @ 2007-10-04 20:43:32 by jwe]
jwe
parents: 6882
diff changeset
109 int nargin = args.length ();
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
110
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11553
diff changeset
111 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
112 print_usage ();
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
113
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
114 octave_value_list retval;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
115
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
116 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
117 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
118 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
119
18077
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
120 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
121 return retval;
ac74b0c4c564 avoid overflow when passing problem dimensions to qhull with --enable-64
John W. Eaton <jwe@octave.org>
parents: 17787
diff changeset
122
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
123 points = points.transpose ();
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
124
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
125 std::string options;
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
126
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
127 if (dim <= 4)
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
128 options = " Qt";
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
129 else
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
130 options = " Qt Qx";
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
131
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
132 if (nargin == 2)
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
133 {
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
134 if (args(1).is_string ())
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
135 options = " " + args(1).string_value ();
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
136 else if (args(1).is_empty ())
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
137 ; // Use default options.
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
138 else if (args(1).is_cellstr ())
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
139 {
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
140 options = "";
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
141
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
142 Array<std::string> tmp = args(1).cellstr_value ();
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
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 for (octave_idx_type i = 0; i < tmp.numel (); i++)
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
145 options += " " + tmp(i);
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
146 }
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
147 else
20825
66cd994d1b79 eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents: 20790
diff changeset
148 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
149 }
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 boolT ismalloc = false;
10257
cd550069240e assume vsnprintf from gnulib; use sstream instead of snprintf
John W. Eaton <jwe@octave.org>
parents: 10154
diff changeset
152
14339
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
153 unwind_protect frame;
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
154
14309
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
155 // 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
156 #if defined (OCTAVE_HAVE_WINDOWS_FILESYSTEM) && ! defined (OCTAVE_HAVE_POSIX_FILESYSTEM)
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
157 FILE *outfile = gnulib::fopen ("NUL", "w");
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
158 #else
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
159 FILE *outfile = gnulib::fopen ("/dev/null", "w");
824e5d362aba Fix Qhull calling convention by passing true file pointer to qh_new_qhull
Rik <octave@nomad.inbox5.com>
parents: 14138
diff changeset
160 #endif
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
161 FILE *errfile = stderr;
14339
3e4350f09a55 close temporary files opened for Qhull
John W. Eaton <jwe@octave.org>
parents: 14310
diff changeset
162
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
163 if (! outfile)
20825
66cd994d1b79 eliminate return statements after calls to error
John W. Eaton <jwe@octave.org>
parents: 20790
diff changeset
164 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
165
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
166 frame.add_fcn (close_fcn, outfile);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
167
13879
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
168 // 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
169
440d7914cf01 fix regression in __voronoi__ and convhulln option processing
John W. Eaton <jwe@octave.org>
parents: 13877
diff changeset
170 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
171
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
172 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
173
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
174 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
175
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
176 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
177 ismalloc, cmd_str, outfile, errfile);
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
178 if (exitcode)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
179 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
180
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
181 bool nonsimp_seen = false;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
182
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
183 octave_idx_type nf = qh num_facets;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
184
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
185 Matrix idx (nf, dim + 1);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
186
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
187 facetT *facet;
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
188
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
189 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
190
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
191 FORALLfacets
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 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
194
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
195 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
196 {
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
197 nonsimp_seen = true;
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
198
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
199 if (cmd.find ("QJ") != std::string::npos)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
200 // Should never happen with QJ.
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
201 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
202 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
203
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
204 if (dim == 3)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
205 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
206 setT *vertices = qh_facet3vertex (facet);
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
207
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
208 vertexT *vertex, **vertexp;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
209
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
210 FOREACHvertex_ (vertices)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
211 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
212
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
213 qh_settempfree (&vertices);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
214 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
215 else
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
216 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
217 if (facet->toporient ^ qh_ORIENTclock)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
218 {
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
219 vertexT *vertex, **vertexp;
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
220
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
221 FOREACHvertex_ (facet->vertices)
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
222 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
223 }
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11553
diff changeset
224 else
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
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 FOREACHvertexreverse12_ (facet->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);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
230 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
231 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
232 if (j < dim)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
233 warning ("convhulln: facet %d only has %d vertices", i, j);
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 i++;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
236 }
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
237
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
238 // 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
239
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
240 if (! nonsimp_seen)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
241 idx.resize (nf, dim, 0.0);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
242
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
243 if (nargout == 2)
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 // 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
246
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
247 realT area;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
248 realT dist;
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 FORALLfacets
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
251 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
252 if (! facet->normal)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
253 continue;
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 if (facet->upperdelaunay && qh ATinfinity)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
256 continue;
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
257
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
258 facet->f.area = area = qh_facetarea (facet);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
259 facet->isarea = True;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
260
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
261 if (qh DELAUNAY)
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 (facet->upperdelaunay == qh UPPERdelaunay)
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
264 qh totarea += area;
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
265 }
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
266 else
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
267 {
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
268 qh totarea += area;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
269 qh_distplane (qh interior_point, facet, &dist);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
270 qh totvol += -dist * area/ qh hull_dim;
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
271 }
10154
40dfc0c99116 DLD-FUNCTIONS/*.cc: untabify
John W. Eaton <jwe@octave.org>
parents: 9786
diff changeset
272 }
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
273
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
274 retval(1) = octave_value (qh totvol);
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
275 }
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
276
20898
8da80da1ac37 maint: Use ovl() more places in the code.
Rik <rik@octave.org>
parents: 20853
diff changeset
277 retval(0) = idx;
6959
47f4f4e88166 [project @ 2007-10-04 20:43:32 by jwe]
jwe
parents: 6882
diff changeset
278
13877
0d32a681d943 * convhulln.cc: Clean up argument parsing and variable decls.
John W. Eaton <jwe@octave.org>
parents: 13746
diff changeset
279 // Free memory from Qhull
6959
47f4f4e88166 [project @ 2007-10-04 20:43:32 by jwe]
jwe
parents: 6882
diff changeset
280 qh_freeqhull (! qh_ALL);
47f4f4e88166 [project @ 2007-10-04 20:43:32 by jwe]
jwe
parents: 6882
diff changeset
281
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
282 int curlong, totlong;
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
283 qh_memfreeshort (&curlong, &totlong);
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
284
11586
12df7854fa7c strip trailing whitespace from source files
John W. Eaton <jwe@octave.org>
parents: 11553
diff changeset
285 if (curlong || totlong)
6959
47f4f4e88166 [project @ 2007-10-04 20:43:32 by jwe]
jwe
parents: 6882
diff changeset
286 warning ("convhulln: did not free %d bytes of long memory (%d pieces)",
13746
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
287 totlong, curlong);
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
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
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
291 #else
21109
bd1752782e56 Use err_disabled_feature, warn_disabled_feature throughout code base.
Rik <rik@octave.org>
parents: 20940
diff changeset
292 err_disabled_feature ("convhulln", "Qhull");
6823
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
293 #endif
9fddcc586065 [project @ 2007-08-24 08:27:27 by dbateman]
dbateman
parents:
diff changeset
294 }
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
295
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
296 /*
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
297 %!testif HAVE_QHULL
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
298 %! 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
299 %! [h, v] = convhulln (cube, "Qt");
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
300 %! assert (size (h), [12 3]);
14310
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
301 %! h = sortrows (sort (h, 2), [1:3]);
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
302 %! 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
303 %! assert (v, 1, 10*eps);
14310
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
304 %! [h2, v2] = convhulln (cube); % Test defaut option = "Qt"
14501
60e5cf354d80 Update %!tests in DLD-FUNCTIONS/ directory with Octave coding conventions.
Rik <octave@nomad.inbox5.com>
parents: 14346
diff changeset
305 %! assert (size (h2), size (h));
14310
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
306 %! h2 = sortrows (sort (h2, 2), [1:3]);
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
307 %! assert (h2, h);
decea31ea010 Fix qhull tests.
Ben Abbott <bpabbott@mac.com>
parents: 14309
diff changeset
308 %! 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
309
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
310 %!testif HAVE_QHULL
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
311 %! 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
312 %! [h, v] = convhulln (cube, "QJ");
17787
175b392e91fe Use GNU style coding conventions for code in libinterp/
Rik <rik@octave.org>
parents: 17744
diff changeset
313 %! 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
314 %! 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
315 %! 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
316
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
317 %!testif HAVE_QHULL
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
318 %! 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
319 %! [h, v] = convhulln (tetrahedron);
7ff0bdc3dc4c Revamp geometry functions dependent on Qhull (Bug #34604, Bug #33346)
Rik <octave@nomad.inbox5.com>
parents: 11586
diff changeset
320 %! 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
321 %! 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
322 %! 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
323
065bc7944335 fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents: 14339
diff changeset
324 %!testif HAVE_QHULL
065bc7944335 fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents: 14339
diff changeset
325 %! triangle=[0 0; 1 1; 1 0; 1 2];
065bc7944335 fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents: 14339
diff changeset
326 %! h = convhulln (triangle);
065bc7944335 fix incorrect results for convhulln in some cases (bug #38013)
John W. Eaton <jwe@octave.org>
parents: 14339
diff changeset
327 %! assert (size (h), [3 2]);
8786
dbd428efbf56 Add calculation convex hull volume
Kai Habel <kai.habel@gmx.de>
parents: 7520
diff changeset
328 */