3Copyright (C) 2000, 2007, 2008, 2009 Kai Habel
5This file is part of Octave.
7Octave is free software; you can redistribute it and/or modify it
8under the terms of the GNU General Public License as published by the
9Free Software Foundation; either version 3 of the License, or (at your
10option) any later version.
12Octave is distributed in the hope that it will be useful, but WITHOUT
13ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
14FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
17You should have received a copy of the GNU General Public License
18along with Octave; see the file COPYING. If not, see
19<http://www.gnu.org/licenses/>.
2429. July 2000 - Kai Habel: first release
252002-04-22 Paul Kienzle
26* Use warning(...) function rather than writing to cerr
28* add support for consistent winding in all dimensions; output is
29* guaranteed to be simplicial.
45#include <qhull/qhull_a.h>
48#if defined (HAVE_QHULL) && defined (NEED_QHULL_VERSION)
49char qh_version[] = "convhulln.oct 2007-07-24";
52DEFUN_DLD (convhulln, args, nargout,
54@deftypefn {Loadable Function} {@var{h} =} convhulln (@var{p})\n\
55@deftypefnx {Loadable Function} {@var{h} =} convhulln (@var{p}, @var{opt})\n\
56@deftypefnx {Loadable Function} {[@var{h}, @var{v}] =} convhulln (@dots{})\n\
57Return an index vector to the points of the enclosing convex hull.\n\
58The input matrix of size [n, dim] contains n points of dimension dim.\n\n\
59If a second optional argument is given, it must be a string or cell array\n\
60of strings containing options for the underlying qhull command. (See\n\
61the Qhull documentation for the available options.) The default options\n\
63If the second output @var{V} is requested the volume of the convex hull is\n\
65@seealso{convhull, delaunayn}\n\
68 octave_value_list retval;
73 int nargin = args.length ();
74 if (nargin < 1 || nargin > 2)
82 if (args (1).is_string ())
83 options = args(1).string_value ();
84 else if (args(1).is_cell ())
86 Cell c = args(1).cell_value ();
88 for (octave_idx_type i = 0; i < c.numel (); i++)
90 if (! c.elem(i).is_string ())
92 error ("convhulln: second argument must be a string or cell array of strings");
96 options = options + c.elem(i).string_value() + " ";
101 error ("convhulln: second argument must be a string or cell array of strings");
106 // turn on some consistency checks
107 options = "s Qci Tcv";
109 Matrix p (args(0).matrix_value ());
111 const octave_idx_type dim = p.columns ();
112 const octave_idx_type n = p.rows ();
115 double *pt_array = p.fortran_vec ();
117 boolT ismalloc = False;
119 std::ostringstream buf;
121 buf << "qhull QJ " << options;
123 std::string buf_string = buf.str ();
125 // FIXME -- we can't just pass buf_string.c_str () to qh_new_qhull
126 // because the argument is not declared const. Ugh. Unless
127 // qh_new_qhull really needs to modify this argument, someone should
130 OCTAVE_LOCAL_BUFFER (char, flags, buf_string.length () + 1);
132 strcpy (flags, buf_string.c_str ());
134 if (! qh_new_qhull (dim, n, pt_array, ismalloc, flags, 0, stderr))
136 // If you want some debugging information replace the NULL
137 // pointer with stdout
139 vertexT *vertex, **vertexp;
142 unsigned int nf = qh num_facets;
144 Matrix idx (nf, dim);
146 octave_idx_type j, i = 0;
150 if (! facet->simplicial)
151 // should never happen with QJ
152 error ("convhulln: non-simplicial facet");
156 vertices = qh_facet3vertex (facet);
157 FOREACHvertex_ (vertices)
158 idx(i, j++) = 1 + qh_pointid(vertex->point);
159 qh_settempfree (&vertices);
163 if (facet->toporient ^ qh_ORIENTclock)
165 FOREACHvertex_ (facet->vertices)
166 idx(i, j++) = 1 + qh_pointid(vertex->point);
170 FOREACHvertexreverse12_ (facet->vertices)
171 idx(i, j++) = 1 + qh_pointid(vertex->point);
175 // likewise but less fatal
176 warning ("facet %d only has %d vertices", i, j);
181 // calculate volume of convex hull
182 // taken from qhull src/geom2.c
192 if (facet->upperdelaunay && qh ATinfinity)
195 facet->f.area = area = qh_facetarea (facet);
196 facet->isarea = True;
200 if (facet->upperdelaunay == qh UPPERdelaunay)
206 qh_distplane (qh interior_point, facet, &dist);
207 qh totvol += -dist * area/ qh hull_dim;
211 retval(1) = octave_value (qh totvol);
214 retval(0) = octave_value (idx);
218 qh_freeqhull (! qh_ALL);
220 // free short memory and memory allocator
221 int curlong, totlong;
222 qh_memfreeshort (&curlong, &totlong);
224 if (curlong || totlong)
225 warning ("convhulln: did not free %d bytes of long memory (%d pieces)",
228 error ("convhulln: not available in this version of Octave");
236%! 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];
237%! [h, v] = convhulln(cube,'Pp');
238%! assert (v, 1.0, 1e6*eps);
240%! tetrahedron = [1 1 1;-1 -1 1;-1 1 -1;1 -1 -1];
241%! [h, v] = convhulln(tetrahedron,'Pp');
242%! assert (v, 8/3, 1e6*eps);