changelog shortlog tags changeset files revisions annotate raw

src/DLD-FUNCTIONS/convhulln.cc

changeset 10289: 4b124317dc38
parent:40dfc0c99116
author: John W. Eaton <jwe@octave.org>
date: Tue Feb 09 20:58:55 2010 -0500 (37 minutes ago)
permissions: -rw-r--r--
description: base_properties::set_children: account for hidden children
1/*
2
3Copyright (C) 2000, 2007, 2008, 2009 Kai Habel
4
5This file is part of Octave.
6
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.
11
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
15for more details.
16
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/>.
20
21*/
22
23/*
2429. July 2000 - Kai Habel: first release
252002-04-22 Paul Kienzle
26* Use warning(...) function rather than writing to cerr
272006-05-01 Tom Holroyd
28* add support for consistent winding in all dimensions; output is
29* guaranteed to be simplicial.
30*/
31
32#ifdef HAVE_CONFIG_H
33#include <config.h>
34#endif
35
36#include <sstream>
37
38#include "Cell.h"
39#include "defun-dld.h"
40#include "error.h"
41#include "oct-obj.h"
42#include "parse.h"
43
44extern "C" {
45#include <qhull/qhull_a.h>
46}
47
48#if defined (HAVE_QHULL) && defined (NEED_QHULL_VERSION)
49char qh_version[] = "convhulln.oct 2007-07-24";
50#endif
51
52DEFUN_DLD (convhulln, args, nargout,
53 "-*- texinfo -*-\n\
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\
62are \"s Qci Tcv\".\n\
63If the second output @var{V} is requested the volume of the convex hull is\n\
64calculated.\n\n\
65@seealso{convhull, delaunayn}\n\
66@end deftypefn")
67{
68 octave_value_list retval;
69
70#ifdef HAVE_QHULL
71 std::string options;
72
73 int nargin = args.length ();
74 if (nargin < 1 || nargin > 2)
75 {
76 print_usage ();
77 return retval;
78 }
79
80 if (nargin == 2)
81 {
82 if (args (1).is_string ())
83 options = args(1).string_value ();
84 else if (args(1).is_cell ())
85 {
86 Cell c = args(1).cell_value ();
87 options = "";
88 for (octave_idx_type i = 0; i < c.numel (); i++)
89 {
90 if (! c.elem(i).is_string ())
91 {
92 error ("convhulln: second argument must be a string or cell array of strings");
93 return retval;
94 }
95
96 options = options + c.elem(i).string_value() + " ";
97 }
98 }
99 else
100 {
101 error ("convhulln: second argument must be a string or cell array of strings");
102 return retval;
103 }
104 }
105 else
106 // turn on some consistency checks
107 options = "s Qci Tcv";
108
109 Matrix p (args(0).matrix_value ());
110
111 const octave_idx_type dim = p.columns ();
112 const octave_idx_type n = p.rows ();
113 p = p.transpose ();
114
115 double *pt_array = p.fortran_vec ();
116
117 boolT ismalloc = False;
118
119 std::ostringstream buf;
120
121 buf << "qhull QJ " << options;
122
123 std::string buf_string = buf.str ();
124
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
128 // fix QHULL.
129
130 OCTAVE_LOCAL_BUFFER (char, flags, buf_string.length () + 1);
131
132 strcpy (flags, buf_string.c_str ());
133
134 if (! qh_new_qhull (dim, n, pt_array, ismalloc, flags, 0, stderr))
135 {
136 // If you want some debugging information replace the NULL
137 // pointer with stdout
138
139 vertexT *vertex, **vertexp;
140 facetT *facet;
141 setT *vertices;
142 unsigned int nf = qh num_facets;
143
144 Matrix idx (nf, dim);
145
146 octave_idx_type j, i = 0;
147 FORALLfacets
148 {
149 j = 0;
150 if (! facet->simplicial)
151 // should never happen with QJ
152 error ("convhulln: non-simplicial facet");
153
154 if (dim == 3)
155 {
156 vertices = qh_facet3vertex (facet);
157 FOREACHvertex_ (vertices)
158 idx(i, j++) = 1 + qh_pointid(vertex->point);
159 qh_settempfree (&vertices);
160 }
161 else
162 {
163 if (facet->toporient ^ qh_ORIENTclock)
164 {
165 FOREACHvertex_ (facet->vertices)
166 idx(i, j++) = 1 + qh_pointid(vertex->point);
167 }
168 else
169 {
170 FOREACHvertexreverse12_ (facet->vertices)
171 idx(i, j++) = 1 + qh_pointid(vertex->point);
172 }
173 }
174 if (j < dim)
175 // likewise but less fatal
176 warning ("facet %d only has %d vertices", i, j);
177 i++;
178 }
179
180 if (nargout == 2)
181 // calculate volume of convex hull
182 // taken from qhull src/geom2.c
183 {
184 realT area;
185 realT dist;
186
187 FORALLfacets
188 {
189 if (! facet->normal)
190 continue;
191
192 if (facet->upperdelaunay && qh ATinfinity)
193 continue;
194
195 facet->f.area = area = qh_facetarea (facet);
196 facet->isarea = True;
197
198 if (qh DELAUNAY)
199 {
200 if (facet->upperdelaunay == qh UPPERdelaunay)
201 qh totarea += area;
202 }
203 else
204 {
205 qh totarea += area;
206 qh_distplane (qh interior_point, facet, &dist);
207 qh totvol += -dist * area/ qh hull_dim;
208 }
209 }
210
211 retval(1) = octave_value (qh totvol);
212 }
213
214 retval(0) = octave_value (idx);
215 }
216
217 // free long memory
218 qh_freeqhull (! qh_ALL);
219
220 // free short memory and memory allocator
221 int curlong, totlong;
222 qh_memfreeshort (&curlong, &totlong);
223
224 if (curlong || totlong)
225 warning ("convhulln: did not free %d bytes of long memory (%d pieces)",
226 totlong, curlong);
227#else
228 error ("convhulln: not available in this version of Octave");
229#endif
230
231 return retval;
232}
233
234/*
235%!test
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);
239%!test
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);
243*/