00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032 #ifdef HAVE_CONFIG_H
00033 #include <config.h>
00034 #endif
00035
00036 #include <sstream>
00037
00038 #include "Cell.h"
00039 #include "defun-dld.h"
00040 #include "error.h"
00041 #include "oct-obj.h"
00042 #include "parse.h"
00043 #include "unwind-prot.h"
00044
00045 #if defined (HAVE_QHULL)
00046 # include "oct-qhull.h"
00047 # if defined (NEED_QHULL_VERSION)
00048 char qh_version[] = "convhulln.oct 2007-07-24";
00049 # endif
00050 #endif
00051
00052 static void
00053 close_fcn (FILE *f)
00054 {
00055 gnulib::fclose (f);
00056 }
00057
00058 DEFUN_DLD (convhulln, args, nargout,
00059 "-*- texinfo -*-\n\
00060 @deftypefn {Loadable Function} {@var{h} =} convhulln (@var{pts})\n\
00061 @deftypefnx {Loadable Function} {@var{h} =} convhulln (@var{pts}, @var{options})\n\
00062 @deftypefnx {Loadable Function} {[@var{h}, @var{v}] =} convhulln (@dots{})\n\
00063 Compute the convex hull of the set of points @var{pts} which is a matrix\n\
00064 of size [n, dim] containing n points in a space of dimension dim.\n\
00065 The hull @var{h} is an index vector into the set of points and specifies\n\
00066 which points form the enclosing hull.\n\
00067 \n\
00068 An optional second argument, which must be a string or cell array of strings,\n\
00069 contains options passed to the underlying qhull command.\n\
00070 See the documentation for the Qhull library for details\n\
00071 @url{http://www.qhull.org/html/qh-quick.htm#options}.\n\
00072 The default options depend on the dimension of the input:\n\
00073 \n\
00074 @itemize\n\
00075 @item 2D, 3D, 4D: @var{options} = @code{@{\"Qt\"@}}\n\
00076 \n\
00077 @item 5D and higher: @var{options} = @code{@{\"Qt\", \"Qx\"@}}\n\
00078 @end itemize\n\
00079 \n\
00080 If @var{options} is not present or @code{[]} then the default arguments are\n\
00081 used. Otherwise, @var{options} replaces the default argument list.\n\
00082 To append user options to the defaults it is necessary to repeat the\n\
00083 default arguments in @var{options}. Use a null string to pass no arguments.\n\
00084 \n\
00085 If the second output @var{v} is requested the volume of the enclosing\n\
00086 convex hull is calculated.\n\n\
00087 @seealso{convhull, delaunayn, voronoin}\n\
00088 @end deftypefn")
00089 {
00090 octave_value_list retval;
00091
00092 #if defined (HAVE_QHULL)
00093
00094 int nargin = args.length ();
00095 if (nargin < 1 || nargin > 2)
00096 {
00097 print_usage ();
00098 return retval;
00099 }
00100
00101 Matrix points (args(0).matrix_value ());
00102 const octave_idx_type dim = points.columns ();
00103 const octave_idx_type num_points = points.rows ();
00104
00105 points = points.transpose ();
00106
00107 std::string options;
00108
00109 if (dim <= 4)
00110 options = " Qt";
00111 else
00112 options = " Qt Qx";
00113
00114 if (nargin == 2)
00115 {
00116 if (args(1).is_string ())
00117 options = " " + args(1).string_value ();
00118 else if (args(1).is_empty ())
00119 ;
00120 else if (args(1).is_cellstr ())
00121 {
00122 options = "";
00123
00124 Array<std::string> tmp = args(1).cellstr_value ();
00125
00126 for (octave_idx_type i = 0; i < tmp.numel (); i++)
00127 options += " " + tmp(i);
00128 }
00129 else
00130 {
00131 error ("convhulln: OPTIONS must be a string, cell array of strings, or empty");
00132 return retval;
00133 }
00134 }
00135
00136 boolT ismalloc = false;
00137
00138 unwind_protect frame;
00139
00140
00141 #if defined (OCTAVE_HAVE_WINDOWS_FILESYSTEM) && ! defined (OCTAVE_HAVE_POSIX_FILESYSTEM)
00142 FILE *outfile = gnulib::fopen ("NUL", "w");
00143 #else
00144 FILE *outfile = gnulib::fopen ("/dev/null", "w");
00145 #endif
00146 FILE *errfile = stderr;
00147
00148 if (outfile)
00149 frame.add_fcn (close_fcn, outfile);
00150 else
00151 {
00152 error ("convhulln: unable to create temporary file for output");
00153 return retval;
00154 }
00155
00156
00157
00158 std::string cmd = "qhull" + options;
00159
00160 OCTAVE_LOCAL_BUFFER (char, cmd_str, cmd.length () + 1);
00161
00162 strcpy (cmd_str, cmd.c_str ());
00163
00164 int exitcode = qh_new_qhull (dim, num_points, points.fortran_vec (),
00165 ismalloc, cmd_str, outfile, errfile);
00166 if (! exitcode)
00167 {
00168 bool nonsimp_seen = false;
00169
00170 octave_idx_type nf = qh num_facets;
00171
00172 Matrix idx (nf, dim + 1);
00173
00174 facetT *facet;
00175
00176 octave_idx_type i = 0;
00177
00178 FORALLfacets
00179 {
00180 octave_idx_type j = 0;
00181
00182 if (! (nonsimp_seen || facet->simplicial || qh hull_dim == 2))
00183 {
00184 nonsimp_seen = true;
00185
00186 if (cmd.find ("QJ") != std::string::npos)
00187 {
00188
00189 error ("convhulln: qhull failed: option 'QJ' returned non-simplicial facet");
00190 return retval;
00191 }
00192 }
00193
00194 if (dim == 3)
00195 {
00196 setT *vertices = qh_facet3vertex (facet);
00197
00198 vertexT *vertex, **vertexp;
00199
00200 FOREACHvertex_ (vertices)
00201 idx(i, j++) = 1 + qh_pointid(vertex->point);
00202
00203 qh_settempfree (&vertices);
00204 }
00205 else
00206 {
00207 if (facet->toporient ^ qh_ORIENTclock)
00208 {
00209 vertexT *vertex, **vertexp;
00210
00211 FOREACHvertex_ (facet->vertices)
00212 idx(i, j++) = 1 + qh_pointid(vertex->point);
00213 }
00214 else
00215 {
00216 vertexT *vertex, **vertexp;
00217
00218 FOREACHvertexreverse12_ (facet->vertices)
00219 idx(i, j++) = 1 + qh_pointid(vertex->point);
00220 }
00221 }
00222 if (j < dim)
00223 warning ("convhulln: facet %d only has %d vertices", i, j);
00224
00225 i++;
00226 }
00227
00228
00229
00230 if (! nonsimp_seen)
00231 idx.resize (nf, dim, 0.0);
00232
00233 if (nargout == 2)
00234 {
00235
00236
00237 realT area;
00238 realT dist;
00239
00240 FORALLfacets
00241 {
00242 if (! facet->normal)
00243 continue;
00244
00245 if (facet->upperdelaunay && qh ATinfinity)
00246 continue;
00247
00248 facet->f.area = area = qh_facetarea (facet);
00249 facet->isarea = True;
00250
00251 if (qh DELAUNAY)
00252 {
00253 if (facet->upperdelaunay == qh UPPERdelaunay)
00254 qh totarea += area;
00255 }
00256 else
00257 {
00258 qh totarea += area;
00259 qh_distplane (qh interior_point, facet, &dist);
00260 qh totvol += -dist * area/ qh hull_dim;
00261 }
00262 }
00263
00264 retval(1) = octave_value (qh totvol);
00265 }
00266
00267 retval(0) = idx;
00268 }
00269 else
00270 error ("convhulln: qhull failed");
00271
00272
00273 qh_freeqhull (! qh_ALL);
00274
00275 int curlong, totlong;
00276 qh_memfreeshort (&curlong, &totlong);
00277
00278 if (curlong || totlong)
00279 warning ("convhulln: did not free %d bytes of long memory (%d pieces)",
00280 totlong, curlong);
00281
00282 #else
00283 error ("convhulln: not available in this version of Octave");
00284 #endif
00285
00286 return retval;
00287 }
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321