qhull.h 35 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724
  1. /* qhull.h -- user-level header file for using qhull.a library
  2. see README, qhull_a.h
  3. see user.h for all user defineable constants
  4. copyright (c) 1993, 1997 The Geometry Center
  5. defines qh_qh, global data structure for qhull.
  6. NOTE: access to qh_qh is via the 'qh' macro. This allows
  7. qh_qh to be either a pointer or a structure. An example
  8. of using qh is "qh DROPdim" which accesses the DROPdim
  9. field of qh_qh. Similarly, access to qh_qhstat is via
  10. the 'qhstat' macro.
  11. includes function prototypes for qhull.c, geom.c, global.c, io.c, user.c
  12. use mem.h for mem.c
  13. use set.h for set.c
  14. see unix.c for an example of using qhull.h
  15. recompile qhull if you change this file
  16. */
  17. #ifndef qhDEFqhull
  18. #define qhDEFqhull 1
  19. /* =========================== -included files ============== */
  20. #include <setjmp.h>
  21. #include <float.h>
  22. #include <time.h>
  23. #if __MWERKS__ && __POWERPC__
  24. #include <SIOUX.h>
  25. #include <Files.h>
  26. #include <Desk.h>
  27. #endif
  28. #ifndef __STDC__
  29. #ifndef __cplusplus
  30. #error neither __STDC__ nor __cplusplus are defined. Please use ANSI C or C++ to compile Qhull.
  31. #endif
  32. #endif
  33. #include "user.h" /* user defineable constants */
  34. /* ============ -types- ====================
  35. Note: could use 'float' for data and 'double' for calculations (realT vs. coordT)
  36. This requires many type casts, and adjusted error bounds.
  37. Also C compilers may do expressions in double anyway.
  38. */
  39. #define coordT realT /* for stored coordinates and coefficients */
  40. #define pointT coordT /* array of hull_dim coordinates */
  41. #define flagT unsigned int /* Boolean flag */
  42. #define boolT unsigned int /* Boolean value, nothing much is portable */
  43. #ifdef False
  44. #undef False
  45. #endif
  46. #ifdef True
  47. #undef True
  48. #endif
  49. #define False 0
  50. #define True 1
  51. /* -CENTERtype to distinguish facet->center */
  52. typedef enum {qh_ASnone= 0, qh_ASvoronoi, qh_AScentrum} qh_CENTER;
  53. /* -output formats for printing (qh PRINTout).
  54. notes:
  55. some of these names are similar to qh names. The similar names are only
  56. used in switch statements in qh_printbegin() etc.
  57. */
  58. typedef enum {qh_PRINTnone= 0, qh_PRINTarea, qh_PRINTaverage,
  59. qh_PRINTcoplanars, qh_PRINTcentrums,
  60. qh_PRINTfacets, qh_PRINTfacets_xridge,
  61. qh_PRINTgeom, qh_PRINTids, qh_PRINTinner, qh_PRINTneighbors,
  62. qh_PRINTnormals, qh_PRINTouter, qh_PRINTincidences,
  63. qh_PRINTmathematica, qh_PRINTmerges, qh_PRINToff,
  64. qh_PRINToptions, qh_PRINTpointintersect, qh_PRINTpointnearest,
  65. qh_PRINTpoints, qh_PRINTqhull, qh_PRINTsize, qh_PRINTsummary,
  66. qh_PRINTtriangles, qh_PRINTvertices, qh_PRINTvneighbors,
  67. qh_PRINTEND} qh_PRINT;
  68. /*---------------------------
  69. arguments
  70. */
  71. #define qh_ALL True /* argument for printall and checkall parameters*/
  72. #define qh_NOupper True /* argument for qh_findbest */
  73. /*-------------------------------------------
  74. -ERR - qhull exit codes, for indicating errors
  75. */
  76. #define qh_ERRnone 0 /* no error occurred during qhull */
  77. #define qh_ERRinput 1 /* input inconsistency */
  78. #define qh_ERRsingular 2 /* singular input data */
  79. #define qh_ERRprec 3 /* precision error */
  80. #define qh_ERRmem 4 /* insufficient memory, matches mem.h */
  81. #define qh_ERRqhull 5 /* internal error detected, matches mem.h */
  82. /* ============ -structures- ====================
  83. each of the following structures is defined by a typedef
  84. all realT and coordT fields occur at the beginning of a structure
  85. (otherwise space may be wasted due to alignment)
  86. define all flags together and pack into 32-bit number
  87. */
  88. typedef struct vertexT vertexT;
  89. typedef struct ridgeT ridgeT;
  90. typedef struct facetT facetT;
  91. #ifndef DEFsetT
  92. #define DEFsetT 1
  93. typedef struct setT setT; /* defined in set.h */
  94. #endif
  95. /* ----------------------------------------------
  96. -facetT- specifies a facet.
  97. qhull() generates the hull as a list of facets.
  98. */
  99. struct facetT {
  100. #if !qh_COMPUTEfurthest
  101. coordT furthestdist;/* distance to furthest point of outsideset */
  102. #endif
  103. #if qh_MAXoutside
  104. coordT maxoutside; /* max computed distance of point to facet
  105. Before QHULLfinished this is an approximation
  106. since maxdist not always set for mergefacet
  107. Actual outer plane is +DISTround and
  108. computed outer plane is +2*DISTround */
  109. #endif
  110. coordT offset; /* exact offset of hyperplane from origin */
  111. coordT *normal; /* normal of hyperplane, hull_dim coefficients */
  112. union { /* in order of testing */
  113. realT area; /* area of facet, only in io.c if ->isarea */
  114. facetT *replace; /* replacement facet if ->visible and NEWfacets
  115. is NULL only if qh_mergedegen_redundant or interior */
  116. facetT *samecycle; /* cycle of facets from the same visible/horizon intersection,
  117. if ->newfacet */
  118. facetT *newcycle; /* in horizon facet, current samecycle of new facets */
  119. }f;
  120. coordT *center; /* centrum for convexity, qh CENTERtype == qh_AScentrum */
  121. /* Voronoi center, qh CENTERtype == qh_ASvoronoi */
  122. facetT *previous; /* previous facet in the facet_list */
  123. facetT *next; /* next facet in the facet_list */
  124. setT *vertices; /* vertices for this facet, inverse sorted by id
  125. if simplicial, 1st vertex was apex/furthest */
  126. setT *ridges; /* explicit ridges for nonsimplicial facets.
  127. for simplicial facets, neighbors defines ridge */
  128. setT *neighbors; /* neighbors of the facet. If simplicial, the kth
  129. neighbor is opposite the kth vertex, and the first
  130. neighbor is the horizon facet for the first vertex*/
  131. setT *outsideset; /* set of points outside this facet
  132. if non-empty, last point is furthest
  133. if NARROWhull, includes coplanars for partitioning*/
  134. setT *coplanarset; /* set of points coplanar with this facet
  135. > min_vertex and <= facet->max_outside
  136. a point is assigned to the furthest facet
  137. if non-empty, last point is furthest away */
  138. unsigned visitid; /* visit_id, for visiting all neighbors,
  139. all uses are independent */
  140. unsigned id; /* unique identifier from qh facet_id */
  141. unsigned nummerge:9; /* number of merges */
  142. #define qh_MAXnummerge 511 /* 2^9-1 */
  143. flagT newfacet:1; /* True if facet on qh newfacet_list (new or merged) */
  144. flagT visible:1; /* True if visible facet (will be deleted) */
  145. flagT toporient:1; /* True if created with top orientation
  146. after merging, use ridge orientation */
  147. flagT simplicial:1;/* True if simplicial facet, ->ridges may be implicit */
  148. flagT seen:1; /* used to perform operations only once, like visitid */
  149. flagT flipped:1; /* True if facet is flipped */
  150. flagT upperdelaunay:1; /* True if facet is upper envelope of Delaunay triangulation */
  151. flagT notfurthest:1; /* True if last point of outsideset is not furthest*/
  152. /* flags primarily for output */
  153. flagT good:1; /* True if a facet marked good for output */
  154. flagT isarea:1; /* True if facet->f.area is defined */
  155. /* flags for merging */
  156. flagT dupridge:1; /* True if duplicate ridge in facet */
  157. flagT mergeridge:1; /* True if facet or neighbor contains a qh_MERGEridge
  158. ->normal defined (also defined for mergeridge2) */
  159. flagT mergeridge2:1; /* True if neighbor contains a qh_MERGEridge (mark_dupridges */
  160. flagT coplanar:1; /* True if horizon facet is coplanar at last use */
  161. flagT mergehorizon:1; /* True if will merge into horizon (->coplanar) */
  162. flagT cycledone:1;/* True if mergecycle_all already done */
  163. flagT tested:1; /* True if facet convexity has been tested (false after merge */
  164. flagT keepcentrum:1; /* True if keep old centrum after a merge */
  165. flagT newmerge:1; /* True if facet is newly merged for reducevertices */
  166. flagT degenerate:1; /* True if facet is degenerate (degen_mergeset) */
  167. flagT redundant:1; /* True if facet is redundant (degen_mergeset) */
  168. };
  169. /*----------------------------------------------
  170. -ridgeT- specifies a ridge
  171. a ridge is hull_dim-1 simplex between two neighboring facets. If the
  172. facets are non-simplicial, there may be more than one ridge between
  173. two facets. E.G. a 4-d hypercube has two triangles between each pair
  174. of neighboring facets.
  175. */
  176. struct ridgeT {
  177. setT *vertices; /* vertices belonging to this ridge, inverse sorted by id
  178. NULL if a degen ridge (matchsame) */
  179. facetT *top; /* top facet this ridge is part of */
  180. facetT *bottom; /* bottom facet this ridge is part of */
  181. unsigned id:24; /* unique identifier, =>room for 8 flags */
  182. flagT seen:1; /* used to perform operations only once */
  183. flagT tested:1; /* True when ridge is tested for convexity */
  184. flagT nonconvex:1; /* True if getmergeset detected a non-convex neighbor
  185. only one ridge between neighbors may have nonconvex */
  186. };
  187. /* ----------------------------------------------
  188. -vertexT- specifies a vertex
  189. */
  190. struct vertexT {
  191. vertexT *next; /* next vertex in vertex_list */
  192. vertexT *previous; /* previous vertex in vertex_list */
  193. pointT *point; /* hull_dim coordinates (coordT) */
  194. setT *neighbors; /* neighboring facets of vertex, qh_vertexneighbors()
  195. inits in io.c or after first merge */
  196. unsigned visitid; /* for use with qh vertex_visit */
  197. unsigned id:24; /* unique identifier, =>room for 8 flags */
  198. flagT seen:1; /* used to perform operations only once */
  199. flagT seen2:1; /* another seen flag */
  200. flagT delridge:1; /* vertex was part of a deleted ridge */
  201. flagT deleted:1; /* true if vertex on qh del_vertices */
  202. flagT newlist:1; /* true if vertex on qh newvertex_list */
  203. };
  204. /* ======= -global variables -qh ============================
  205. all global variables for qhull are in qh, qhmem, and qhstat
  206. qhmem is defined in mem.h and qhstat is defined in stat.h
  207. access to qh_qh is via the "qh" macro. See qh_QHpointer in user.h
  208. */
  209. typedef struct qhT qhT;
  210. #if qh_QHpointer
  211. #define qh qh_qh->
  212. extern qhT *qh_qh; /* allocated in global.c */
  213. #else
  214. #define qh qh_qh.
  215. extern qhT qh_qh;
  216. #endif
  217. extern char qh_version[]; /* defined in unix.c */
  218. struct qhT {
  219. /*-user flags */
  220. boolT ALLpoints; /* true 'Qs' if search all points for initial simplex */
  221. boolT ANGLEmerge; /* true 'Qa' if sort potential merges by angle */
  222. boolT APPROXhull; /* true 'Wn' if MINoutside set */
  223. realT MINoutside; /* 'Wn' min. distance for an outside point */
  224. boolT ATinfinity; /* true 'Qz' if point num_points-1 is "at-infinity"
  225. for improving precision in Delaunay triangulations */
  226. boolT AVOIDold; /* true 'Q4' if avoid old->new merges */
  227. boolT BESToutside; /* true 'Qf' if partition points into best outsideset */
  228. boolT CDDinput; /* true 'Pc' if input uses CDD format (1.0/offset first) */
  229. boolT CDDoutput; /* true 'PC' if print normals in CDD format (offset first) */
  230. boolT CHECKfrequently; /* true 'Tc' if checking frequently */
  231. realT premerge_cos; /* 'A-n' cos_max when pre merging */
  232. realT postmerge_cos; /* 'An' cos_max when post merging */
  233. boolT DELAUNAY; /* true 'd' if computing DELAUNAY triangulation */
  234. boolT DOintersections; /* true 'Gh' if print hyperplane intersections */
  235. int DROPdim; /* drops dim 'GDn' for 4-d -> 3-d output */
  236. boolT FORCEoutput; /* true 'Po' if forcing output despite degeneracies */
  237. int GOODpoint; /* 1+n for 'QGn', good facet if visible/not(-) from point n*/
  238. pointT *GOODpointp; /* the actual point */
  239. boolT GOODthreshold; /* true if qh lower_threshold/upper_threshold defined
  240. false if qh SPLITthreshold */
  241. int GOODvertex; /* 1+n, good facet if vertex for point n */
  242. pointT *GOODvertexp; /* the actual point */
  243. boolT HALFspace; /* true 'Hn,n,n' if half-space intersection */
  244. int IStracing; /* trace execution, 0=none, 1=least, 4=most, -1=events */
  245. int KEEParea; /* 'PAn' number of largest facets to keep */
  246. boolT KEEPcoplanar; /* true 'Qc' if keeping nearest facet for coplanar points */
  247. boolT KEEPinside; /* true 'Qi' if keeping nearest facet for inside points
  248. set automatically if 'd Qc' */
  249. int KEEPmerge; /* 'PMn' number of facets to keep with most merges */
  250. realT KEEPminArea; /* 'PFn' minimum facet area to keep */
  251. boolT MERGEexact; /* true 'Qx' if exact merges (coplanar, degen, dupridge, flipped) */
  252. boolT MERGEindependent; /* true 'Q2' if merging independent sets */
  253. boolT MERGING; /* true if exact-, pre- or post-merging, with angle and centrum tests */
  254. realT premerge_centrum; /* 'C-n' centrum_radius when pre merging. Default is round-off */
  255. realT postmerge_centrum; /* 'Cn' centrum_radius when post merging. Default is round-off */
  256. boolT MERGEvertices; /* true 'Q3' if merging redundant vertices */
  257. realT MAXcoplanar; /* 'Un' max distance below a facet to be coplanar*/
  258. realT MINvisible; /* 'Vn' min. distance for a facet to be visible */
  259. boolT NOnearinside; /* true 'Q8' if ignore near-inside points when partitioning */
  260. boolT NOpremerge; /* true 'Q0' if no defaults for C-0 or Qx */
  261. boolT ONLYgood; /* true 'Qg' if process points with good visible or horizon facets */
  262. boolT ONLYmax; /* true 'Qm' if only process points that increase max_outside */
  263. boolT PICKfurthest; /* true 'Q9' if process furthest of furthest points*/
  264. boolT POINTSmalloc; /* true if qh first_point/num_points allocated */
  265. boolT POSTmerge; /* true if merging after buildhull (Cn or An) */
  266. boolT PREmerge; /* true if merging during buildhull (C-n or A-n) */
  267. /* NOTE: some of these names are similar to qh_PRINT names */
  268. boolT PRINTcentrums; /* true 'Gc' if printing centrums */
  269. boolT PRINTcoplanar; /* true 'Gp' if printing coplanar points */
  270. int PRINTdim; /* print dimension for Geomview output */
  271. boolT PRINTdots; /* true 'Ga' if printing all points as dots */
  272. boolT PRINTgood; /* true 'Pg' if printing good facets */
  273. boolT PRINTinner; /* true 'Gi' if printing inner planes */
  274. boolT PRINTneighbors; /* true 'PG' if printing neighbors of good facets */
  275. boolT PRINTnoplanes; /* true 'Gn' if printing no planes */
  276. boolT PRINToptions1st; /* true 'FO' if printing options to stderr */
  277. boolT PRINTouter; /* true 'Go' if printing outer planes */
  278. boolT PRINTprecision; /* false 'Pp' if not reporting precision problems */
  279. qh_PRINT PRINTout[qh_PRINTEND]; /* list of output formats to print */
  280. boolT PRINTridges; /* true 'Gr' if print ridges */
  281. boolT PRINTspheres; /* true 'Gv' if print vertices as spheres */
  282. boolT PRINTstatistics; /* true 'Ts' if printing statistics to stderr */
  283. boolT PRINTsummary; /* true 's' if printing summary to stderr */
  284. boolT PRINTtransparent; /* true 'Gt' if print transparent outer ridges */
  285. boolT PROJECTdelaunay; /* true if DELAUNAY, no readpoints() and
  286. need projectinput() for Delaunay in qh_init_B */
  287. int PROJECTinput; /* number of projected dimensions 'bn:0Bn:0' */
  288. boolT QUICKhelp; /* true if quick help message for degen input */
  289. boolT RANDOMdist; /* true if randomly change distplane and setfacetplane */
  290. realT RANDOMfactor; /* maximum perturbation */
  291. realT RANDOMa; /* qh_randomfactor is randr * RANDOMa + RANDOMb */
  292. realT RANDOMb;
  293. boolT RANDOMoutside; /* true if select a random outside point */
  294. int REPORTfreq; /* buildtracing reports every n facets */
  295. int REPORTfreq2; /* tracemerging reports every REPORTfreq/2 facets */
  296. int ROTATErandom; /* 'QRn' seed, 0 time, >= rotate input */
  297. boolT SCALEinput; /* true 'Qbk' if scaling input */
  298. boolT SCALElast; /* true 'Qbb' if scale last coord to max prev coord */
  299. boolT SETroundoff; /* true 'E' if qh DISTround is predefined */
  300. boolT SKIPcheckmax; /* true 'Q5' if skip qh_check_maxout */
  301. boolT SKIPconvex; /* true 'Q6' if skip convexity testing during pre-merge */
  302. boolT SPLITthresholds; /* true if upper_/lower_threshold defines a region
  303. used only for printing (not for qh ONLYgood) */
  304. int STOPcone; /* 'TCn' 1+n for stopping after cone for point n*/
  305. int STOPpoint; /* 'TVn' 'TV-n' 1+n for stopping after/before(-)
  306. adding point n */
  307. boolT TESTvneighbors; /* true 'Qv' if test vertex neighbors at end */
  308. int TRACElevel; /* 'Tn' conditional IStracing level */
  309. int TRACEpoint; /* 'TPn' start tracing when point n is a vertex */
  310. realT TRACEdist; /* 'TWn' start tracing when merge distance too big */
  311. int TRACEmerge; /* 'TMn' start tracing before this merge */
  312. boolT UPPERdelaunay; /* true 'Qu' if computing furthest-site Delaunay */
  313. boolT VERIFYoutput; /* true 'Tv' if verify output at end of qhull */
  314. boolT VIRTUALmemory; /* true 'Q7' if depth-first processing in buildhull */
  315. boolT VORONOI; /* true 'v' if computing Voronoi diagram */
  316. /* -input constants */
  317. realT AREAfactor; /* 1/(hull_dim-1)! for converting det's to area */
  318. boolT DOcheckmax; /* true if calling qh_check_maxout (qh_initqhull_globals) */
  319. char *feasible_string; /* feasible point 'Hn,n,n' for half-space intersection */
  320. coordT *feasible_point; /* as coordinates, both malloc'd */
  321. boolT GETarea; /* true if need to compute facet areas in io.c */
  322. boolT KEEPnearinside; /* true if near-inside points in coplanarset */
  323. int input_dim; /* dimension of input, set by initbuffers */
  324. int num_points; /* number of input points */
  325. pointT *first_point; /* array of input points, malloc'd */
  326. int hull_dim; /* dimension of hull, set by initbuffers */
  327. char qhull_command[256];/* command line that invoked this program */
  328. char rbox_command[256]; /* command line that produced the input points */
  329. char qhull_options[512]; /* descriptive list of options */
  330. int qhull_optionslen; /* length of line */
  331. boolT VERTEXneighbors; /* true if maintaining vertex neighbors */
  332. boolT ZEROcentrum; /* true if 'C-0' or 'C-0 Qx'. sets ZEROall_ok */
  333. realT *upper_threshold; /* don't print if facet->normal[k]>=upper_threshold[k]
  334. must set either GOODthreshold or SPLITthreshold
  335. if Delaunay, default is 0.0 for upper envelope */
  336. realT *lower_threshold; /* don't print if facet->normal[k] <=lower_threshold[k] */
  337. realT *upper_bound; /* scale point[k] to new upper bound */
  338. realT *lower_bound; /* scale point[k] to new lower bound
  339. project if both upper_ and lower_bound == 0 */
  340. /* -precision constants, computed in qh_maxmin */
  341. realT ANGLEround; /* max round off error for angles */
  342. realT centrum_radius; /* max centrum radius for convexity (roundoff added) */
  343. realT cos_max; /* max cosine for convexity (roundoff added) */
  344. realT DISTround; /* max round off error for distances, 'E' overrides */
  345. realT MAXlowcoord; /* max coordinate in all but last dimension */
  346. realT maxmaxcoord; /* max coordinate in any dimension */
  347. realT MINdenom_1; /* min. abs. value for 1/x */
  348. realT MINdenom; /* use divzero if denominator < MINdenom */
  349. realT MINdenom_1_2; /* min. abs. val for 1/x that allows normalization */
  350. realT MINdenom_2; /* use divzero if denominator < MINdenom_2 */
  351. realT MINnorm; /* min. norm for not redoing qh_sethyperplane_det */
  352. boolT NARROWhull; /* set in qh_initialhull if angle < qh_MAXnarrow */
  353. realT *NEARzero; /* hull_dim array for near zero in gausselim */
  354. realT NEARinside; /* keep points for qh_check_maxout if close to facet */
  355. realT ONEmerge; /* max distance for merging simplicial facets */
  356. realT WIDEfacet; /* size of wide facet for skipping ridge in
  357. area computation and locking centrum */
  358. /* -internal constants */
  359. char qhull[sizeof("qhull")]; /* for checking ownership */
  360. void *old_stat; /* pointer to saved qh_qhstat, qh_save_qhull */
  361. jmp_buf errexit; /* exit label for qh_errexit, defined by setjmp() */
  362. char jmpXtra[40]; /* extra bytes in case jmp_buf is defined wrong */
  363. FILE *fin; /* pointer to input file, init by qh_meminit */
  364. FILE *fout; /* pointer to output file */
  365. FILE *ferr; /* pointer to error file */
  366. pointT *interior_point; /* center point of the initial simplex*/
  367. int normal_size; /* size in bytes for facet normals and point coords*/
  368. int center_size; /* size in bytes for Voronoi centers */
  369. int TEMPsize; /* size for small, temporary sets (in quick mem) */
  370. /* -list of all facets, from facet_list to facet_tail, see qh_appendfacet */
  371. facetT *facet_list; /* first facet */
  372. facetT *facet_tail; /* end of facet_list (dummy facet) */
  373. facetT *facet_next; /* next facet for buildhull()
  374. all previous facets do not have outside sets*/
  375. facetT *newfacet_list; /* list of new facets to end of facet_list */
  376. facetT *visible_list; /* list of visible facets preceeding newfacet_list,
  377. facet->visible set */
  378. int num_visible; /* current number of visible facets */
  379. unsigned tracefacet_id; /* set at init, then can print whenever */
  380. facetT *tracefacet; /* set in newfacet/mergefacet, undone in delfacet*/
  381. unsigned tracevertex_id; /* set at buildtracing, can print whenever */
  382. vertexT *tracevertex; /* set in newvertex, undone in delvertex*/
  383. vertexT *vertex_list; /* list of all vertices, to vertex_tail */
  384. vertexT *vertex_tail;
  385. vertexT *newvertex_list; /* list of vertices in newfacet_list, to vertex_tail
  386. all vertices have 'newlist' set */
  387. int num_facets; /* number of facets in facet_list
  388. includes visble faces (num_visible) */
  389. int num_vertices; /* number of vertices in facet_list */
  390. int num_outside; /* number of points in outsidesets (for tracing) */
  391. int num_good; /* number of good facets (after findgood_all) */
  392. unsigned facet_id; /* id of next, new facet from newfacet() */
  393. unsigned ridge_id; /* id of next, new ridge from newridge() */
  394. unsigned vertex_id; /* id of next, new vertex from newvertex() */
  395. /* -variables */
  396. unsigned long hulltime; /* ignore time to set up input and randomize */
  397. /* use unsigned to avoid wrap-around errors */
  398. qh_CENTER CENTERtype; /* current type of facet->center, qh_CENTER */
  399. int furthest_id; /* pointid of furthest point, for tracing */
  400. facetT *GOODclosest; /* closest facet to GOODthreshold in qh_findgood */
  401. realT max_outside; /* maximum distance from a point to a facet,
  402. before roundoff, not simplicial vertices
  403. actual outer plane is +DISTround and
  404. computed outer plane is +2*DISTround */
  405. realT max_vertex; /* maximum distance (>0) from vertex to a facet,
  406. before roundoff, not simplicial vertices */
  407. realT min_vertex; /* minimum distance (<0) from vertex to a facet,
  408. before roundoff, not simplicial vertices
  409. actual inner plane is -DISTround and
  410. computed inner plane is -2*DISTround */
  411. boolT NEWfacets; /* true while visible facets invalid due to new or merge
  412. from makecone/attachnewfacets to deletevisible */
  413. boolT findbestnew; /* true if partitioning calls qh_findbestnew */
  414. boolT findbest_notsharp; /* true if new facets are at least 90 degrees */
  415. boolT NOerrexit; /* true if qh_errexit is not available */
  416. realT PRINTcradius; /* radius for printing centrums */
  417. realT PRINTradius; /* radius for printing vertex spheres and points */
  418. boolT POSTmerging; /* true when post merging */
  419. int printoutvar; /* temporary variable for qh_printbegin, etc. */
  420. int printoutnum; /* number of facets printed */
  421. boolT QHULLfinished; /* True after qhull() is finished */
  422. realT totarea; /* total facet area computed by qh_getarea */
  423. realT totvol; /* total volume computed by qh_getarea */
  424. unsigned int visit_id; /* unique id for searching neighborhoods, */
  425. unsigned int vertex_visit; /* unique id for searching vertices */
  426. boolT ZEROall_ok; /* True if qh_checkzero always succeeds */
  427. boolT WAScoplanar; /* True if qh_partitioncoplanar (qh_check_maxout) */
  428. /* -sets */
  429. setT *facet_mergeset; /* temporary set of merges to be done */
  430. setT *degen_mergeset; /* temporary set of degenerate and redundant merges */
  431. setT *initial_points; /* initial simplex for buildhull() */
  432. setT *hash_table; /* hash table for matching ridges in qh_matchfacets
  433. size is setsize() */
  434. int num_hashentries; /* current number of hashentries */
  435. setT *other_points; /* additional points (first is qh interior_point) */
  436. setT *del_vertices; /* vertices to partition and delete with visible
  437. facets. Have deleted set for checkfacet */
  438. /* -buffers */
  439. coordT *gm_matrix; /* (dim+1)Xdim matrix for geom.c */
  440. coordT **gm_row; /* array of gm_matrix rows */
  441. char* line; /* malloc'd input line of maxline+1 chars */
  442. int maxline;
  443. coordT *half_space; /* malloc'd input array for halfspace (qh normal_size+coordT) */
  444. coordT *temp_malloc; /* malloc'd input array for points */
  445. /* -statics */
  446. boolT ERREXITcalled; /* true during errexit (prevents duplicate calls */
  447. boolT firstcentrum; /* for qh_printcentrum */
  448. unsigned lastreport; /* for qh_buildtracing */
  449. int mergereport; /* for qh_tracemerging */
  450. boolT old_randomdist; /* save RANDOMdist when io, tracing, or statistics */
  451. int ridgeoutnum; /* number of ridges in 4OFF output */
  452. void *old_qhstat; /* for saving qh_qhstat in save_qhull() */
  453. setT *old_tempstack; /* for saving qhmem.tempstack in save_qhull */
  454. setT *searchset; /* set of facets for searching in qh_findbest() */
  455. int rand_seed; /* for qh_rand/qh_srand */
  456. };
  457. /* =========== -macros- =========================
  458. -otherfacet_(ridge, facet) return neighboring facet for a ridge in facet
  459. -getid_(p) return id or -1 if NULL
  460. */
  461. #define otherfacet_(ridge, facet) \
  462. (((ridge)->top == (facet)) ? (ridge)->bottom : (ridge)->top)
  463. #define getid_(p) ((p) ? (p)->id : -1)
  464. /* ---------------------------------------------
  465. -FORALL and FOREACH macros
  466. These all iterate using a variable of the same name, e.g. FORALLfacets
  467. and FOREACHfacet_ uses 'facet' declared by 'facetT *facet'. The macros
  468. may use auxiliary variables as indicated.
  469. -FORALLfacets iterate over all facets in facetlist
  470. -FORALLpoint_(points, num) iterate over num points (uses 'pointT *pointtemp')
  471. -FORALLvertices iterate over all vertices in vertex_list
  472. -FOREACHfacet_(facets) iterate over facet set (uses 'facetT **facetp')
  473. -FOREACHneighbor_(facet) iterate over facet->neighbors (uses 'facetT **neighborp')
  474. -FOREACHpoint_(points) iterate over point set (uses 'pointT **pointp')
  475. -FOREACHridge_(ridges) iterate over ridge set (uses 'ridgeT **ridgep')
  476. -FOREACHvertex_(vertice) iterate over vertex set (uses 'vertexT **vertexp')
  477. -FOREACHneighbor_(vertex) iterate over neighboring facets to vertex
  478. -FOREACHfacet_i_(facets) iterate over facets by facet_i and facet_n
  479. -FOREACHneighbor_i_(facet) iterate over facet->neighbors by neighbor_i, neighbor_n
  480. -FOREACHvertex_i_(vertices) iterate over vertices by vertex_i, vertex_n
  481. -FOREACHpoint_i_(points) iterate over points by point_i, point_n
  482. -FOREACHridge_i_(ridges) iterate over ridges by ridge_i, ridge_n
  483. -FOREACHneighbor_i_(vertex) iterate over vertex->neighbors by neighbor_i, neighbor_n
  484. WARNING: nested loops can't use the same variable (define another FOREACH)
  485. WARNING: strange behavior if don't fully brace when nested (including
  486. intervening blocks, e.g. FOREACH...{ if () FOREACH...} )
  487. poly.h defines other FOREACH/FORALL macros
  488. set.h defines FOREACHsetelement and contains additional notes
  489. */
  490. #define FORALLfacets for (facet=qh facet_list;facet && facet->next;facet=facet->next)
  491. #define FORALLpoints FORALLpoint_(qh first_point, qh num_points)
  492. #define FORALLvertices for (vertex=qh vertex_list;vertex && vertex->next;vertex= vertex->next)
  493. #define FORALLpoint_(points, num) for(point= (points), \
  494. pointtemp= (points)+qh hull_dim*(num); point < pointtemp; point += qh hull_dim)
  495. #define FOREACHfacet_(facets) FOREACHsetelement_(facetT, facets, facet)
  496. #define FOREACHneighbor_(facet) FOREACHsetelement_(facetT, facet->neighbors, neighbor)
  497. #define FOREACHpoint_(points) FOREACHsetelement_(pointT, points, point)
  498. #define FOREACHridge_(ridges) FOREACHsetelement_(ridgeT, ridges, ridge)
  499. #define FOREACHvertex_(vertices) FOREACHsetelement_(vertexT, vertices,vertex)
  500. #define FOREACHfacet_i_(facets) FOREACHsetelement_i_(facetT, facets, facet)
  501. #define FOREACHneighbor_i_(facet) FOREACHsetelement_i_(facetT, facet->neighbors, neighbor)
  502. #define FOREACHpoint_i_(points) FOREACHsetelement_i_(pointT, points, point)
  503. #define FOREACHridge_i_(ridges) FOREACHsetelement_i_(ridgeT, ridges, ridge)
  504. #define FOREACHvertex_i_(vertices) FOREACHsetelement_i_(vertexT, vertices,vertex)
  505. /* ======= -functions ===========
  506. see corresponding .c file for definitions
  507. Qhull functions (see qhull.c and qhull_a.h)
  508. -qhull construct the convex hull of a set of points
  509. -addpoint add point to hull (must be above facet)
  510. -printsummary print summary about the output
  511. User redefinable functions (see user.c)
  512. -errexit return exitcode to system after an error
  513. -errprint print erroneous facets, ridge, and vertex
  514. -printfacetlist print all fields for a list of facets
  515. -user_memsizes define up to 10 additional quick allocation sizes
  516. Geometric functions (see geom.c and geom.h for other useful functions)
  517. -gram_schmidt implements Gram-Schmidt orthogonalization by rows
  518. -projectinput project input along one or more dimensions + Delaunay projection
  519. -randommatrix generate a random dimXdim matrix in range (-1,1)
  520. -rotatepoints rotate numpoints points by a row matrix
  521. -scaleinput scale input to new lowbound and highbound
  522. -setdelaunay project points to paraboloid for Delaunay triangulation
  523. -sethalfspace_all generate dual for halfspace intersection with feasible point
  524. Global init/free functions (see global.c and qhull_a.h)
  525. -freeqhull free memory used by qhull
  526. -init_A called before error handling initialized
  527. -init_B called after points are defined
  528. -initflags set flags and initialized constants from command line
  529. -restore_qhull restores a saved qhull
  530. -save_qhull saves qhull for later restoring
  531. Input/output functions (see io.c and io.h)
  532. -dfacet print facet by id
  533. -dvertex print vertex by id
  534. -printsummary print summary about the output
  535. -produce_output prints out the result of qhull in desired format
  536. -readpoints read points from input
  537. Polyhedron functions (see poly.c or poly2.c)
  538. -check_output check output data structure according to user flags
  539. -check_points verify that all points are inside the hull
  540. -findbestfacet find facet or best facet below a point
  541. -setvoronoi_all compute Voronoi centers for all facets
  542. -nearvertex return nearest vertex to point
  543. -point return point for a point id, or NULL if unknown
  544. -pointid return id for a point, or -1 if not known
  545. -facetvertices returns temporary set of vertices in a set of facets
  546. -pointfacet return temporary set of facets indexed by point id
  547. -pointvertex return temporary set of vertices indexed by point id
  548. Statistics functions (see stat.c)
  549. -printallstatistics print all statistics
  550. other functions no longer needed by most users
  551. -init_qhull_command build qhull_command from argc/argv
  552. -initqhull_buffers initialize global memory buffers
  553. -initqhull_globals initialize globals
  554. -initqhull_mem initialize mem.c for qhull
  555. -initqhull_start start initialization of qhull
  556. -initthresholds set thresholds for printing and scaling from command line
  557. -findbest find visible facet for a point starting at a facet
  558. -findbestnew find best newfacet for point
  559. */
  560. /********* -qhull.c prototypes (duplicated from qhull_a.h) **********************/
  561. void qh_qhull (void);
  562. boolT qh_addpoint (pointT *furthest, facetT *facet, boolT checkdist);
  563. void qh_printsummary(FILE *fp);
  564. /********* -user.c prototypes (alphabetical) **********************/
  565. void qh_errexit(int exitcode, facetT *facet, ridgeT *ridge);
  566. void qh_errprint(char* string, facetT *atfacet, facetT *otherfacet, ridgeT *atridge, vertexT *atvertex);
  567. void qh_printfacetlist(facetT *facetlist, setT *facets, boolT printall);
  568. void qh_user_memsizes (void);
  569. /***** -geom.c/geom2.c prototypes (duplicated from geom.h) ****************/
  570. facetT *qh_findbest (pointT *point, facetT *facet,
  571. boolT bestoutside, boolT newfacets, boolT noupper,
  572. realT *dist, boolT *isoutside, int *numpart);
  573. facetT *qh_findbestnew (pointT *point, facetT *startfacet,
  574. realT *dist, boolT *isoutside, int *numpart);
  575. boolT qh_gram_schmidt(int dim, realT **rows);
  576. void qh_printsummary(FILE *fp);
  577. void qh_projectinput (void);
  578. void qh_randommatrix (realT *buffer, int dim, realT **row);
  579. void qh_rotateinput (realT **rows);
  580. void qh_scaleinput (void);
  581. void qh_setdelaunay (int dim, int count, pointT *points);
  582. coordT *qh_sethalfspace_all (int dim, int count, coordT *halfspaces, pointT *feasible);
  583. /***** -global.c prototypes (alphabetical) ***********************/
  584. void qh_freebuffers (void);
  585. void qh_freeqhull (boolT allmem);
  586. void qh_init_A (FILE *infile, FILE *outfile, FILE *errfile, int argc, char *argv[]);
  587. void qh_init_B (coordT *points, int numpoints, int dim, boolT ismalloc);
  588. void qh_init_qhull_command (int argc, char *argv[]);
  589. void qh_initbuffers (coordT *points, int numpoints, int dim, boolT ismalloc);
  590. void qh_initflags (char *command);
  591. void qh_initqhull_buffers (void);
  592. void qh_initqhull_globals (coordT *points, int numpoints, int dim, boolT ismalloc);
  593. void qh_initqhull_mem (void);
  594. void qh_initqhull_start (FILE *infile, FILE *outfile, FILE *errfile);
  595. void qh_initthresholds (char *command);
  596. #if qh_QHpointer
  597. void qh_restore_qhull (qhT **oldqh);
  598. qhT *qh_save_qhull (void);
  599. #endif
  600. /***** -io.c prototypes (duplicated from io.h) ***********************/
  601. void dfacet( unsigned id);
  602. void dvertex( unsigned id);
  603. void qh_produce_output(void);
  604. coordT *qh_readpoints(int *numpoints, int *dimension, boolT *ismalloc);
  605. /********* -mem.c prototypes (duplicated from mem.h) **********************/
  606. void qh_meminit (FILE *ferr);
  607. void qh_memfreeshort (int *curlong, int *totlong);
  608. /********* -poly.c/poly2.c prototypes (duplicated from poly.h) **********************/
  609. void qh_check_output (void);
  610. void qh_check_points (void);
  611. setT *qh_facetvertices (facetT *facetlist, setT *facets, boolT allfacets);
  612. facetT *qh_findbestfacet (pointT *point, boolT bestoutside,
  613. realT *bestdist, boolT *isoutside);
  614. vertexT *qh_nearvertex (facetT *facet, pointT *point, realT *bestdistp);
  615. pointT *qh_point (int id);
  616. setT *qh_pointfacet (void /*qh.facet_list*/);
  617. int qh_pointid (pointT *point);
  618. setT *qh_pointvertex (void /*qh.facet_list*/);
  619. void qh_setvoronoi_all (void);
  620. /********* -stat.c prototypes (duplicated from stat.h) **********************/
  621. void qh_printallstatistics (FILE *fp, char *string);
  622. #endif /* qhDEFqhull */