|
||||
File indexing completed on 2025-01-18 10:01:53
0001 /*<html><pre> -<a href="qh-qhull.htm" 0002 >-------------------------------</a><a name="TOP">-</a> 0003 0004 libqhull.h 0005 user-level header file for using qhull.a library 0006 0007 see qh-qhull.htm, qhull_a.h 0008 0009 Copyright (c) 1993-2020 The Geometry Center. 0010 $Id: //main/2019/qhull/src/libqhull/libqhull.h#11 $$Change: 2958 $ 0011 $DateTime: 2020/05/26 16:17:49 $$Author: bbarber $ 0012 0013 NOTE: access to qh_qh is via the 'qh' macro. This allows 0014 qh_qh to be either a pointer or a structure. An example 0015 of using qh is "qh.DROPdim" which accesses the DROPdim 0016 field of qh_qh. Similarly, access to qh_qhstat is via 0017 the 'qhstat' macro. 0018 0019 includes function prototypes for libqhull.c, geom.c, global.c, io.c, user.c 0020 0021 use mem.h for mem.c 0022 use qset.h for qset.c 0023 0024 see unix.c for an example of using libqhull.h 0025 0026 recompile qhull if you change this file 0027 */ 0028 0029 #ifndef qhDEFlibqhull 0030 #define qhDEFlibqhull 1 0031 0032 /*=========================== -included files ==============*/ 0033 0034 /* user.h first for QHULL_CRTDBG */ 0035 #include "user.h" /* user definable constants (e.g., realT). */ 0036 0037 #include "mem.h" /* Needed for qhT in libqhull_r.h. Here for compatibility */ 0038 #include "qset.h" /* Needed for QHULL_LIB_CHECK */ 0039 /* include stat.h after defining boolT. Needed for qhT in libqhull_r.h. Here for compatibility */ 0040 0041 #include <setjmp.h> 0042 #include <float.h> 0043 #include <limits.h> 0044 #include <time.h> 0045 #include <stdio.h> 0046 0047 #if defined(__MWERKS__) && defined(__POWERPC__) 0048 #include <SIOUX.h> 0049 #include <Files.h> 0050 #include <Desk.h> 0051 #endif 0052 0053 #ifndef __STDC__ 0054 #ifndef __cplusplus 0055 #if !defined(_MSC_VER) 0056 #error Neither __STDC__ nor __cplusplus is defined. Please use strict ANSI C or C++ to compile 0057 #error Qhull. You may need to turn off compiler extensions in your project configuration. If 0058 #error your compiler is a standard C compiler, you can delete this warning from libqhull.h 0059 #endif 0060 #endif 0061 #endif 0062 0063 /*============ constants and basic types ====================*/ 0064 0065 extern const char qh_version[]; /* defined in global.c */ 0066 extern const char qh_version2[]; /* defined in global.c */ 0067 0068 /*-<a href="qh-geom.htm#TOC" 0069 >--------------------------------</a><a name="coordT">-</a> 0070 0071 coordT 0072 coordinates and coefficients are stored as realT (i.e., double) 0073 0074 notes: 0075 Qhull works well if realT is 'float'. If so joggle (QJ) is not effective. 0076 0077 Could use 'float' for data and 'double' for calculations (realT vs. coordT) 0078 This requires many type casts, and adjusted error bounds. 0079 Also C compilers may do expressions in double anyway. 0080 */ 0081 #define coordT realT 0082 0083 /*-<a href="qh-geom.htm#TOC" 0084 >--------------------------------</a><a name="pointT">-</a> 0085 0086 pointT 0087 a point is an array of coordinates, usually qh.hull_dim 0088 qh_pointid returns 0089 qh_IDnone if point==0 or qh is undefined 0090 qh_IDinterior for qh.interior_point 0091 qh_IDunknown if point is neither in qh.first_point... nor qh.other_points 0092 0093 notes: 0094 qh.STOPcone and qh.STOPpoint assume that qh_IDunknown==-1 (other negative numbers indicate points) 0095 qh_IDunknown is also returned by getid_() for unknown facet, ridge, or vertex 0096 */ 0097 #define pointT coordT 0098 typedef enum 0099 { 0100 qh_IDnone= -3, qh_IDinterior= -2, qh_IDunknown= -1 0101 } 0102 qh_pointT; 0103 0104 /*-<a href="qh-qhull.htm#TOC" 0105 >--------------------------------</a><a name="flagT">-</a> 0106 0107 flagT 0108 Boolean flag as a bit 0109 */ 0110 #define flagT unsigned int 0111 0112 /*-<a href="qh-qhull.htm#TOC" 0113 >--------------------------------</a><a name="boolT">-</a> 0114 0115 boolT 0116 boolean value, either True or False 0117 0118 notes: 0119 needed for portability 0120 Use qh_False/qh_True as synonyms 0121 */ 0122 #define boolT unsigned int 0123 #ifdef False 0124 #undef False 0125 #endif 0126 #ifdef True 0127 #undef True 0128 #endif 0129 #define False 0 0130 #define True 1 0131 #define qh_False 0 0132 #define qh_True 1 0133 0134 #include "stat.h" /* needs boolT */ 0135 0136 /*-<a href="qh-qhull.htm#TOC" 0137 >--------------------------------</a><a name="CENTERtype">-</a> 0138 0139 qh_CENTER 0140 to distinguish facet->center 0141 */ 0142 typedef enum 0143 { 0144 qh_ASnone= 0, /* If not MERGING and not VORONOI */ 0145 qh_ASvoronoi, /* Set by qh_clearcenters on qh_prepare_output, or if not MERGING and VORONOI */ 0146 qh_AScentrum /* If MERGING (assumed during merging) */ 0147 } 0148 qh_CENTER; 0149 0150 /*-<a href="qh-qhull.htm#TOC" 0151 >--------------------------------</a><a name="qh_PRINT">-</a> 0152 0153 qh_PRINT 0154 output formats for printing (qh.PRINTout). 0155 'Fa' 'FV' 'Fc' 'FC' 0156 0157 0158 notes: 0159 some of these names are similar to qhT names. The similar names are only 0160 used in switch statements in qh_printbegin() etc. 0161 */ 0162 typedef enum {qh_PRINTnone= 0, 0163 qh_PRINTarea, qh_PRINTaverage, /* 'Fa' 'FV' 'Fc' 'FC' */ 0164 qh_PRINTcoplanars, qh_PRINTcentrums, 0165 qh_PRINTfacets, qh_PRINTfacets_xridge, /* 'f' 'FF' 'G' 'FI' 'Fi' 'Fn' */ 0166 qh_PRINTgeom, qh_PRINTids, qh_PRINTinner, qh_PRINTneighbors, 0167 qh_PRINTnormals, qh_PRINTouter, qh_PRINTmaple, /* 'n' 'Fo' 'i' 'm' 'Fm' 'FM', 'o' */ 0168 qh_PRINTincidences, qh_PRINTmathematica, qh_PRINTmerges, qh_PRINToff, 0169 qh_PRINToptions, qh_PRINTpointintersect, /* 'FO' 'Fp' 'FP' 'p' 'FQ' 'FS' */ 0170 qh_PRINTpointnearest, qh_PRINTpoints, qh_PRINTqhull, qh_PRINTsize, 0171 qh_PRINTsummary, qh_PRINTtriangles, /* 'Fs' 'Ft' 'Fv' 'FN' 'Fx' */ 0172 qh_PRINTvertices, qh_PRINTvneighbors, qh_PRINTextremes, 0173 qh_PRINTEND} qh_PRINT; 0174 0175 /*-<a href="qh-qhull.htm#TOC" 0176 >--------------------------------</a><a name="qh_ALL">-</a> 0177 0178 qh_ALL 0179 argument flag for selecting everything 0180 */ 0181 #define qh_ALL True 0182 #define qh_NOupper True /* argument for qh_findbest */ 0183 #define qh_IScheckmax True /* argument for qh_findbesthorizon */ 0184 #define qh_ISnewfacets True /* argument for qh_findbest */ 0185 #define qh_RESETvisible True /* argument for qh_resetlists */ 0186 0187 /*-<a href="qh-qhull.htm#TOC" 0188 >--------------------------------</a><a name="qh_ERR">-</a> 0189 0190 qh_ERR... 0191 Qhull exit status codes, for indicating errors 0192 See: MSG_ERROR (6000) and MSG_WARNING (7000) [user.h] 0193 */ 0194 #define qh_ERRnone 0 /* no error occurred during qhull */ 0195 #define qh_ERRinput 1 /* input inconsistency */ 0196 #define qh_ERRsingular 2 /* singular input data, calls qh_printhelp_singular */ 0197 #define qh_ERRprec 3 /* precision error, calls qh_printhelp_degenerate */ 0198 #define qh_ERRmem 4 /* insufficient memory, matches mem.h */ 0199 #define qh_ERRqhull 5 /* internal error detected, matches mem.h, calls qh_printhelp_internal */ 0200 #define qh_ERRother 6 /* other error detected */ 0201 #define qh_ERRtopology 7 /* topology error, maybe due to nearly adjacent vertices, calls qh_printhelp_topology */ 0202 #define qh_ERRwide 8 /* wide facet error, maybe due to nearly adjacent vertices, calls qh_printhelp_wide */ 0203 #define qh_ERRdebug 9 /* qh_errexit from debugging code */ 0204 0205 /*-<a href="qh-qhull.htm#TOC" 0206 >--------------------------------</a><a name="qh_FILEstderr">-</a> 0207 0208 qh_FILEstderr 0209 Fake stderr to distinguish error output from normal output 0210 For C++ interface. Must redefine qh_fprintf_qhull 0211 */ 0212 #define qh_FILEstderr ((FILE *)1) 0213 0214 /* ============ -structures- ==================== 0215 each of the following structures is defined by a typedef 0216 all realT and coordT fields occur at the beginning of a structure 0217 (otherwise space may be wasted due to alignment) 0218 define all flags together and pack into 32-bit number 0219 0220 DEFsetT is likewise defined in mem.h and qset.h 0221 */ 0222 0223 typedef struct vertexT vertexT; 0224 typedef struct ridgeT ridgeT; 0225 typedef struct facetT facetT; 0226 typedef struct qhT qhT; 0227 #ifndef DEFsetT 0228 #define DEFsetT 1 0229 typedef struct setT setT; /* defined in qset.h */ 0230 #endif 0231 0232 /*-<a href="qh-poly.htm#TOC" 0233 >--------------------------------</a><a name="facetT">-</a> 0234 0235 facetT 0236 defines a facet 0237 0238 notes: 0239 qhull() generates the hull as a list of facets. 0240 0241 topological information: 0242 f.previous,next doubly-linked list of facets, next is always defined 0243 f.vertices set of vertices 0244 f.ridges set of ridges 0245 f.neighbors set of neighbors 0246 f.toporient True if facet has top-orientation (else bottom) 0247 0248 geometric information: 0249 f.offset,normal hyperplane equation 0250 f.maxoutside offset to outer plane -- all points inside 0251 f.center centrum for testing convexity or Voronoi center for output 0252 f.simplicial True if facet is simplicial 0253 f.flipped True if facet does not include qh.interior_point 0254 0255 for constructing hull: 0256 f.visible True if facet on list of visible facets (will be deleted) 0257 f.newfacet True if facet on list of newly created facets 0258 f.coplanarset set of points coplanar with this facet 0259 (includes near-inside points for later testing) 0260 f.outsideset set of points outside of this facet 0261 f.furthestdist distance to furthest point of outside set 0262 f.visitid marks visited facets during a loop 0263 f.replace replacement facet for to-be-deleted, visible facets 0264 f.samecycle,newcycle cycle of facets for merging into horizon facet 0265 0266 see below for other flags and fields 0267 */ 0268 struct facetT { 0269 #if !qh_COMPUTEfurthest 0270 coordT furthestdist;/* distance to furthest point of outsideset */ 0271 #endif 0272 #if qh_MAXoutside 0273 coordT maxoutside; /* max computed distance of point to facet 0274 Before QHULLfinished this is an approximation 0275 since maxdist not always set for qh_mergefacet 0276 Actual outer plane is +DISTround and 0277 computed outer plane is +2*DISTround. 0278 Initial maxoutside is qh.DISTround, otherwise distance tests need to account for DISTround */ 0279 #endif 0280 coordT offset; /* exact offset of hyperplane from origin */ 0281 coordT *normal; /* normal of hyperplane, hull_dim coefficients */ 0282 /* if f.tricoplanar, shared with a neighbor */ 0283 union { /* in order of testing */ 0284 realT area; /* area of facet, only in io.c if f.isarea */ 0285 facetT *replace; /* replacement facet for qh.NEWfacets with f.visible 0286 NULL if qh_mergedegen_redundant, interior, or !NEWfacets */ 0287 facetT *samecycle; /* cycle of facets from the same visible/horizon intersection, 0288 if ->newfacet */ 0289 facetT *newcycle; /* in horizon facet, current samecycle of new facets */ 0290 facetT *trivisible; /* visible facet for ->tricoplanar facets during qh_triangulate() */ 0291 facetT *triowner; /* owner facet for ->tricoplanar, !isarea facets w/ ->keepcentrum */ 0292 }f; 0293 coordT *center; /* set according to qh.CENTERtype */ 0294 /* qh_ASnone: no center (not MERGING) */ 0295 /* qh_AScentrum: centrum for testing convexity (qh_getcentrum) */ 0296 /* assumed qh_AScentrum while merging */ 0297 /* qh_ASvoronoi: Voronoi center (qh_facetcenter) */ 0298 /* after constructing the hull, it may be changed (qh_clearcenter) */ 0299 /* if tricoplanar and !keepcentrum, shared with a neighbor */ 0300 facetT *previous; /* previous facet in the facet_list or NULL, for C++ interface */ 0301 facetT *next; /* next facet in the facet_list or facet_tail */ 0302 setT *vertices; /* vertices for this facet, inverse sorted by ID 0303 if simplicial, 1st vertex was apex/furthest 0304 qh_reduce_vertices removes extraneous vertices via qh_remove_extravertices 0305 if f.visible, vertices may be on qh.del_vertices */ 0306 setT *ridges; /* explicit ridges for nonsimplicial facets or nonsimplicial neighbors. 0307 For simplicial facets, neighbors define the ridges 0308 qh_makeridges() converts simplicial facets by creating ridges prior to merging 0309 If qh.NEWtentative, new facets have horizon ridge, but not vice versa 0310 if f.visible && qh.NEWfacets, ridges is empty */ 0311 setT *neighbors; /* neighbors of the facet. Neighbors may be f.visible 0312 If simplicial, the kth neighbor is opposite the kth vertex and the 0313 first neighbor is the horizon facet for the first vertex. 0314 dupridges marked by qh_DUPLICATEridge (0x01) and qh_MERGEridge (0x02) 0315 if f.visible && qh.NEWfacets, neighbors is empty */ 0316 setT *outsideset; /* set of points outside this facet 0317 if non-empty, last point is furthest 0318 if NARROWhull, includes coplanars (less than qh.MINoutside) for partitioning*/ 0319 setT *coplanarset; /* set of points coplanar with this facet 0320 >= qh.min_vertex and <= facet->max_outside 0321 a point is assigned to the furthest facet 0322 if non-empty, last point is furthest away */ 0323 unsigned int visitid; /* visit_id, for visiting all neighbors, 0324 all uses are independent */ 0325 unsigned int id; /* unique identifier from qh.facet_id, 1..qh.facet_id, 0 is sentinel, printed as 'f%d' */ 0326 unsigned int nummerge:9; /* number of merges */ 0327 #define qh_MAXnummerge 511 /* 2^9-1 */ 0328 /* 23 flags (at most 23 due to nummerge), printed by "flags:" in io.c */ 0329 flagT tricoplanar:1; /* True if TRIangulate and simplicial and coplanar with a neighbor */ 0330 /* all tricoplanars share the same apex */ 0331 /* all tricoplanars share the same ->center, ->normal, ->offset, ->maxoutside */ 0332 /* ->keepcentrum is true for the owner. It has the ->coplanareset */ 0333 /* if ->degenerate, does not span facet (one logical ridge) */ 0334 /* during qh_triangulate, f.trivisible points to original facet */ 0335 flagT newfacet:1; /* True if facet on qh.newfacet_list (new/qh.first_newfacet or merged) */ 0336 flagT visible:1; /* True if visible facet (will be deleted) */ 0337 flagT toporient:1; /* True if created with top orientation 0338 after merging, use ridge orientation */ 0339 flagT simplicial:1;/* True if simplicial facet, ->ridges may be implicit */ 0340 flagT seen:1; /* used to perform operations only once, like visitid */ 0341 flagT seen2:1; /* used to perform operations only once, like visitid */ 0342 flagT flipped:1; /* True if facet is flipped */ 0343 flagT upperdelaunay:1; /* True if facet is upper envelope of Delaunay triangulation */ 0344 flagT notfurthest:1; /* True if last point of outsideset is not furthest */ 0345 0346 /*-------- flags primarily for output ---------*/ 0347 flagT good:1; /* True if a facet marked good for output */ 0348 flagT isarea:1; /* True if facet->f.area is defined */ 0349 0350 /*-------- flags for merging ------------------*/ 0351 flagT dupridge:1; /* True if facet has one or more dupridge in a new facet (qh_matchneighbor), 0352 a dupridge has a subridge shared by more than one new facet */ 0353 flagT mergeridge:1; /* True if facet or neighbor has a qh_MERGEridge (qh_mark_dupridges) 0354 ->normal defined for mergeridge and mergeridge2 */ 0355 flagT mergeridge2:1; /* True if neighbor has a qh_MERGEridge (qh_mark_dupridges) */ 0356 flagT coplanarhorizon:1; /* True if horizon facet is coplanar at last use */ 0357 flagT mergehorizon:1; /* True if will merge into horizon (its first neighbor w/ f.coplanarhorizon). */ 0358 flagT cycledone:1;/* True if mergecycle_all already done */ 0359 flagT tested:1; /* True if facet convexity has been tested (false after merge */ 0360 flagT keepcentrum:1; /* True if keep old centrum after a merge, or marks owner for ->tricoplanar 0361 Set by qh_updatetested if more than qh_MAXnewcentrum extra vertices 0362 Set by qh_mergefacet if |maxdist| > qh.WIDEfacet */ 0363 flagT newmerge:1; /* True if facet is newly merged for reducevertices */ 0364 flagT degenerate:1; /* True if facet is degenerate (degen_mergeset or ->tricoplanar) */ 0365 flagT redundant:1; /* True if facet is redundant (degen_mergeset) 0366 Maybe merge degenerate and redundant to gain another flag */ 0367 }; 0368 0369 0370 /*-<a href="qh-poly.htm#TOC" 0371 >--------------------------------</a><a name="ridgeT">-</a> 0372 0373 ridgeT 0374 defines a ridge 0375 0376 notes: 0377 a ridge is hull_dim-1 simplex between two neighboring facets. If the 0378 facets are non-simplicial, there may be more than one ridge between 0379 two facets. E.G. a 4-d hypercube has two triangles between each pair 0380 of neighboring facets. 0381 0382 topological information: 0383 vertices a set of vertices 0384 top,bottom neighboring facets with orientation 0385 0386 geometric information: 0387 tested True if ridge is clearly convex 0388 nonconvex True if ridge is non-convex 0389 */ 0390 struct ridgeT { 0391 setT *vertices; /* vertices belonging to this ridge, inverse sorted by ID 0392 NULL if a degen ridge (matchsame) */ 0393 facetT *top; /* top facet for this ridge */ 0394 facetT *bottom; /* bottom facet for this ridge 0395 ridge oriented by odd/even vertex order and top/bottom */ 0396 unsigned int id; /* unique identifier. Same size as vertex_id, printed as 'r%d' */ 0397 flagT seen:1; /* used to perform operations only once */ 0398 flagT tested:1; /* True when ridge is tested for convexity by centrum or opposite vertices */ 0399 flagT nonconvex:1; /* True if getmergeset detected a non-convex neighbor 0400 only one ridge between neighbors may have nonconvex */ 0401 flagT mergevertex:1; /* True if pending qh_appendvertexmerge due to 0402 qh_maybe_duplicateridge or qh_maybe_duplicateridges 0403 disables check for duplicate vertices in qh_checkfacet */ 0404 flagT mergevertex2:1; /* True if qh_drop_mergevertex of MRGvertices, printed but not used */ 0405 flagT simplicialtop:1; /* True if top was simplicial (original vertices) */ 0406 flagT simplicialbot:1; /* True if bottom was simplicial (original vertices) 0407 use qh_test_centrum_merge if top and bot, need to retest since centrum may change */ 0408 }; 0409 0410 /*-<a href="qh-poly.htm#TOC" 0411 >--------------------------------</a><a name="vertexT">-</a> 0412 0413 vertexT 0414 defines a vertex 0415 0416 topological information: 0417 next,previous doubly-linked list of all vertices 0418 neighbors set of adjacent facets (only if qh.VERTEXneighbors) 0419 0420 geometric information: 0421 point array of DIM3 coordinates 0422 */ 0423 struct vertexT { 0424 vertexT *next; /* next vertex in vertex_list or vertex_tail */ 0425 vertexT *previous; /* previous vertex in vertex_list or NULL, for C++ interface */ 0426 pointT *point; /* hull_dim coordinates (coordT) */ 0427 setT *neighbors; /* neighboring facets of vertex, qh_vertexneighbors() 0428 initialized in io.c or after first merge 0429 qh_update_vertices for qh_addpoint or qh_triangulate 0430 updated by merges 0431 qh_order_vertexneighbors by 2-d (orientation) 3-d (adjacency), n-d (f.visitid,id) */ 0432 unsigned int id; /* unique identifier, 1..qh.vertex_id, 0 for sentinel, printed as 'r%d' */ 0433 unsigned int visitid; /* for use with qh.vertex_visit, size must match */ 0434 flagT seen:1; /* used to perform operations only once */ 0435 flagT seen2:1; /* another seen flag */ 0436 flagT deleted:1; /* vertex will be deleted via qh.del_vertices */ 0437 flagT delridge:1; /* vertex belonged to a deleted ridge, cleared by qh_reducevertices */ 0438 flagT newfacet:1; /* true if vertex is in a new facet 0439 vertex is on qh.newvertex_list and it has a facet on qh.newfacet_list 0440 or vertex is on qh.newvertex_list due to qh_newvertices while merging 0441 cleared by qh_resetlists */ 0442 flagT partitioned:1; /* true if deleted vertex has been partitioned */ 0443 }; 0444 0445 /*======= -global variables -qh ============================*/ 0446 0447 /*-<a href="qh-globa.htm#TOC" 0448 >--------------------------------</a><a name="qh">-</a> 0449 0450 qhT 0451 All global variables for qhull are in qhT, qhmemT, and qhstatT 0452 0453 notes: 0454 qhmem is defined in mem.h, qhstat is defined in stat.h, qhrbox is defined in rboxpoints.h 0455 Access to qh_qh is via the "qh" macro. See qh_QHpointer in user.h 0456 0457 All global variables for qhull are in qh, qhmem, and qhstat 0458 qh must be unique for each instance of qhull 0459 qhstat may be shared between qhull instances. 0460 qhmem may be shared across multiple instances of Qhull. 0461 Rbox uses global variables rbox_inuse and rbox, but does not persist data across calls. 0462 0463 Qhull is not multi-threaded. Global state could be stored in thread-local storage. 0464 0465 QHULL_LIB_CHECK checks that a program and the corresponding 0466 qhull library were built with the same type of header files. 0467 0468 QHULL_LIB_TYPE is QHULL_NON_REENTRANT, QHULL_QH_POINTER, or QHULL_REENTRANT 0469 */ 0470 0471 #define QHULL_NON_REENTRANT 0 0472 #define QHULL_QH_POINTER 1 0473 #define QHULL_REENTRANT 2 0474 0475 #ifdef qh_QHpointer_dllimport 0476 #define qh qh_qh-> 0477 __declspec(dllimport) extern qhT *qh_qh; /* allocated in global.c */ 0478 #define QHULL_LIB_TYPE QHULL_QH_POINTER 0479 0480 #elif qh_QHpointer 0481 #define qh qh_qh-> 0482 extern qhT *qh_qh; /* allocated in global.c */ 0483 #define QHULL_LIB_TYPE QHULL_QH_POINTER 0484 0485 #elif defined(qh_dllimport) 0486 #define qh qh_qh. 0487 __declspec(dllimport) extern qhT qh_qh; /* allocated in global.c */ 0488 #define QHULL_LIB_TYPE QHULL_NON_REENTRANT 0489 0490 #else 0491 #define qh qh_qh. 0492 extern qhT qh_qh; 0493 #define QHULL_LIB_TYPE QHULL_NON_REENTRANT 0494 #endif 0495 0496 #define QHULL_LIB_CHECK qh_lib_check(QHULL_LIB_TYPE, sizeof(qhT), sizeof(vertexT), sizeof(ridgeT), sizeof(facetT), sizeof(setT), sizeof(qhmemT)); 0497 #define QHULL_LIB_CHECK_RBOX qh_lib_check(QHULL_LIB_TYPE, sizeof(qhT), sizeof(vertexT), sizeof(ridgeT), sizeof(facetT), 0, 0); 0498 0499 struct qhT { 0500 0501 /*-<a href="qh-globa.htm#TOC" 0502 >--------------------------------</a><a name="qh-const">-</a> 0503 0504 qh constants 0505 configuration flags and constants for Qhull 0506 0507 notes: 0508 The user configures Qhull by defining flags. They are 0509 copied into qh by qh_setflags(). qh-quick.htm#options defines the flags. 0510 */ 0511 boolT ALLpoints; /* true 'Qs' if search all points for initial simplex */ 0512 boolT ALLOWshort; /* true 'Qa' allow input with fewer or more points than coordinates */ 0513 boolT ALLOWwarning; /* true 'Qw' if allow option warnings */ 0514 boolT ALLOWwide; /* true 'Q12' if allow wide facets and wide dupridges, c.f. qh_WIDEmaxoutside */ 0515 boolT ANGLEmerge; /* true 'Q1' if sort potential merges by type/angle instead of type/distance */ 0516 boolT APPROXhull; /* true 'Wn' if MINoutside set */ 0517 realT MINoutside; /* Minimum distance for an outside point ('Wn' or 2*qh.MINvisible) */ 0518 boolT ANNOTATEoutput; /* true 'Ta' if annotate output with message codes */ 0519 boolT ATinfinity; /* true 'Qz' if point num_points-1 is "at-infinity" 0520 for improving precision in Delaunay triangulations */ 0521 boolT AVOIDold; /* true 'Q4' if avoid old->new merges */ 0522 boolT BESToutside; /* true 'Qf' if partition points into best outsideset */ 0523 boolT CDDinput; /* true 'Pc' if input uses CDD format (1.0/offset first) */ 0524 boolT CDDoutput; /* true 'PC' if print normals in CDD format (offset first) */ 0525 boolT CHECKduplicates; /* true 'Q15' if qh_maybe_duplicateridges after each qh_mergefacet */ 0526 boolT CHECKfrequently; /* true 'Tc' if checking frequently */ 0527 realT premerge_cos; /* 'A-n' cos_max when pre merging */ 0528 realT postmerge_cos; /* 'An' cos_max when post merging */ 0529 boolT DELAUNAY; /* true 'd' or 'v' if computing DELAUNAY triangulation */ 0530 boolT DOintersections; /* true 'Gh' if print hyperplane intersections */ 0531 int DROPdim; /* drops dim 'GDn' for 4-d -> 3-d output */ 0532 boolT FLUSHprint; /* true 'Tf' if flush after qh_fprintf for segfaults */ 0533 boolT FORCEoutput; /* true 'Po' if forcing output despite degeneracies */ 0534 int GOODpoint; /* 'QGn' or 'QG-n' (n+1, n-1), good facet if visible from point n (or not) */ 0535 pointT *GOODpointp; /* the actual point */ 0536 boolT GOODthreshold; /* true 'Pd/PD' if qh.lower_threshold/upper_threshold defined 0537 set if qh.UPPERdelaunay (qh_initbuild) 0538 false if qh.SPLITthreshold */ 0539 int GOODvertex; /* 'QVn' or 'QV-n' (n+1, n-1), good facet if vertex for point n (or not) */ 0540 pointT *GOODvertexp; /* the actual point */ 0541 boolT HALFspace; /* true 'Hn,n,n' if halfspace intersection */ 0542 boolT ISqhullQh; /* Set by Qhull.cpp on initialization */ 0543 int IStracing; /* 'Tn' trace execution, 0=none, 1=least, 4=most, -1=events */ 0544 int KEEParea; /* 'PAn' number of largest facets to keep */ 0545 boolT KEEPcoplanar; /* true 'Qc' if keeping nearest facet for coplanar points */ 0546 boolT KEEPinside; /* true 'Qi' if keeping nearest facet for inside points 0547 set automatically if 'd Qc' */ 0548 int KEEPmerge; /* 'PMn' number of facets to keep with most merges */ 0549 realT KEEPminArea; /* 'PFn' minimum facet area to keep */ 0550 realT MAXcoplanar; /* 'Un' max distance below a facet to be coplanar*/ 0551 int MAXwide; /* 'QWn' max ratio for wide facet, otherwise error unless Q12-allow-wide */ 0552 boolT MERGEexact; /* true 'Qx' if exact merges (concave, degen, dupridge, flipped) 0553 tested by qh_checkzero and qh_test_*_merge */ 0554 boolT MERGEindependent; /* true if merging independent sets of coplanar facets. 'Q2' disables */ 0555 boolT MERGING; /* true if exact-, pre- or post-merging, with angle and centrum tests */ 0556 realT premerge_centrum; /* 'C-n' centrum_radius when pre merging. Default is round-off */ 0557 realT postmerge_centrum; /* 'Cn' centrum_radius when post merging. Default is round-off */ 0558 boolT MERGEpinched; /* true 'Q14' if merging pinched vertices due to dupridge */ 0559 boolT MERGEvertices; /* true if merging redundant vertices, 'Q3' disables or qh.hull_dim > qh_DIMmergeVertex */ 0560 realT MINvisible; /* 'Vn' min. distance for a facet to be visible */ 0561 boolT NOnarrow; /* true 'Q10' if no special processing for narrow distributions */ 0562 boolT NOnearinside; /* true 'Q8' if ignore near-inside points when partitioning, qh_check_points may fail */ 0563 boolT NOpremerge; /* true 'Q0' if no defaults for C-0 or Qx */ 0564 boolT ONLYgood; /* true 'Qg' if process points with good visible or horizon facets */ 0565 boolT ONLYmax; /* true 'Qm' if only process points that increase max_outside */ 0566 boolT PICKfurthest; /* true 'Q9' if process furthest of furthest points*/ 0567 boolT POSTmerge; /* true if merging after buildhull ('Cn' or 'An') */ 0568 boolT PREmerge; /* true if merging during buildhull ('C-n' or 'A-n') */ 0569 /* NOTE: some of these names are similar to qh_PRINT names */ 0570 boolT PRINTcentrums; /* true 'Gc' if printing centrums */ 0571 boolT PRINTcoplanar; /* true 'Gp' if printing coplanar points */ 0572 int PRINTdim; /* print dimension for Geomview output */ 0573 boolT PRINTdots; /* true 'Ga' if printing all points as dots */ 0574 boolT PRINTgood; /* true 'Pg' if printing good facets 0575 PGood set if 'd', 'PAn', 'PFn', 'PMn', 'QGn', 'QG-n', 'QVn', or 'QV-n' */ 0576 boolT PRINTinner; /* true 'Gi' if printing inner planes */ 0577 boolT PRINTneighbors; /* true 'PG' if printing neighbors of good facets */ 0578 boolT PRINTnoplanes; /* true 'Gn' if printing no planes */ 0579 boolT PRINToptions1st; /* true 'FO' if printing options to stderr */ 0580 boolT PRINTouter; /* true 'Go' if printing outer planes */ 0581 boolT PRINTprecision; /* false 'Pp' if not reporting precision problems */ 0582 qh_PRINT PRINTout[qh_PRINTEND]; /* list of output formats to print */ 0583 boolT PRINTridges; /* true 'Gr' if print ridges */ 0584 boolT PRINTspheres; /* true 'Gv' if print vertices as spheres */ 0585 boolT PRINTstatistics; /* true 'Ts' if printing statistics to stderr */ 0586 boolT PRINTsummary; /* true 's' if printing summary to stderr */ 0587 boolT PRINTtransparent; /* true 'Gt' if print transparent outer ridges */ 0588 boolT PROJECTdelaunay; /* true if DELAUNAY, no readpoints() and 0589 need projectinput() for Delaunay in qh_init_B */ 0590 int PROJECTinput; /* number of projected dimensions 'bn:0Bn:0' */ 0591 boolT RANDOMdist; /* true 'Rn' if randomly change distplane and setfacetplane */ 0592 realT RANDOMfactor; /* maximum random perturbation */ 0593 realT RANDOMa; /* qh_randomfactor is randr * RANDOMa + RANDOMb */ 0594 realT RANDOMb; 0595 boolT RANDOMoutside; /* true 'Qr' if select a random outside point */ 0596 int REPORTfreq; /* 'TFn' buildtracing reports every n facets */ 0597 int REPORTfreq2; /* tracemerging reports every REPORTfreq/2 facets */ 0598 int RERUN; /* 'TRn' rerun qhull n times (qh.build_cnt) */ 0599 int ROTATErandom; /* 'QRn' n<-1 random seed, n==-1 time is seed, n==0 random rotation by time, n>0 rotate input */ 0600 boolT SCALEinput; /* true 'Qbk' if scaling input */ 0601 boolT SCALElast; /* true 'Qbb' if scale last coord to max prev coord */ 0602 boolT SETroundoff; /* true 'En' if qh.DISTround is predefined */ 0603 boolT SKIPcheckmax; /* true 'Q5' if skip qh_check_maxout, qh_check_points may fail */ 0604 boolT SKIPconvex; /* true 'Q6' if skip convexity testing during pre-merge */ 0605 boolT SPLITthresholds; /* true 'Pd/PD' if upper_/lower_threshold defines a region 0606 else qh.GOODthresholds 0607 set if qh.DELAUNAY (qh_initbuild) 0608 used only for printing (!for qh.ONLYgood) */ 0609 int STOPadd; /* 'TAn' 1+n for stop after adding n vertices */ 0610 int STOPcone; /* 'TCn' 1+n for stopping after cone for point n */ 0611 /* also used by qh_build_withresart for err exit*/ 0612 int STOPpoint; /* 'TVn' 'TV-n' 1+n for stopping after/before(-) 0613 adding point n */ 0614 int TESTpoints; /* 'QTn' num of test points after qh.num_points. Test points always coplanar. */ 0615 boolT TESTvneighbors; /* true 'Qv' if test vertex neighbors at end */ 0616 int TRACElevel; /* 'Tn' conditional IStracing level */ 0617 int TRACElastrun; /* qh.TRACElevel applies to last qh.RERUN */ 0618 int TRACEpoint; /* 'TPn' start tracing when point n is a vertex, use qh_IDunknown (-1) after qh_buildhull and qh_postmerge */ 0619 realT TRACEdist; /* 'TWn' start tracing when merge distance too big */ 0620 int TRACEmerge; /* 'TMn' start tracing before this merge */ 0621 boolT TRIangulate; /* true 'Qt' if triangulate non-simplicial facets */ 0622 boolT TRInormals; /* true 'Q11' if triangulate duplicates ->normal and ->center (sets Qt) */ 0623 boolT UPPERdelaunay; /* true 'Qu' if computing furthest-site Delaunay */ 0624 boolT USEstdout; /* true 'Tz' if using stdout instead of stderr */ 0625 boolT VERIFYoutput; /* true 'Tv' if verify output at end of qhull */ 0626 boolT VIRTUALmemory; /* true 'Q7' if depth-first processing in buildhull */ 0627 boolT VORONOI; /* true 'v' if computing Voronoi diagram, also sets qh.DELAUNAY */ 0628 0629 /*--------input constants ---------*/ 0630 realT AREAfactor; /* 1/(hull_dim-1)! for converting det's to area */ 0631 boolT DOcheckmax; /* true if calling qh_check_maxout (!qh.SKIPcheckmax && qh.MERGING) */ 0632 char *feasible_string; /* feasible point 'Hn,n,n' for halfspace intersection */ 0633 coordT *feasible_point; /* as coordinates, both malloc'd */ 0634 boolT GETarea; /* true 'Fa', 'FA', 'FS', 'PAn', 'PFn' if compute facet area/Voronoi volume in io.c */ 0635 boolT KEEPnearinside; /* true if near-inside points in coplanarset */ 0636 int hull_dim; /* dimension of hull, set by initbuffers */ 0637 int input_dim; /* dimension of input, set by initbuffers */ 0638 int num_points; /* number of input points */ 0639 pointT *first_point; /* array of input points, see POINTSmalloc */ 0640 boolT POINTSmalloc; /* true if qh.first_point/num_points allocated */ 0641 pointT *input_points; /* copy of original qh.first_point for input points for qh_joggleinput */ 0642 boolT input_malloc; /* true if qh.input_points malloc'd */ 0643 char qhull_command[256];/* command line that invoked this program */ 0644 int qhull_commandsiz2; /* size of qhull_command at qh_clear_outputflags */ 0645 char rbox_command[256]; /* command line that produced the input points */ 0646 char qhull_options[512];/* descriptive list of options */ 0647 int qhull_optionlen; /* length of last line */ 0648 int qhull_optionsiz; /* size of qhull_options at qh_build_withrestart */ 0649 int qhull_optionsiz2; /* size of qhull_options at qh_clear_outputflags */ 0650 int run_id; /* non-zero, random identifier for this instance of qhull */ 0651 boolT VERTEXneighbors; /* true if maintaining vertex neighbors */ 0652 boolT ZEROcentrum; /* true if 'C-0' or 'C-0 Qx' and not post-merging or 'A-n'. Sets ZEROall_ok */ 0653 realT *upper_threshold; /* don't print if facet->normal[k]>=upper_threshold[k] 0654 must set either GOODthreshold or SPLITthreshold 0655 if qh.DELAUNAY, default is 0.0 for upper envelope (qh_initbuild) */ 0656 realT *lower_threshold; /* don't print if facet->normal[k] <=lower_threshold[k] */ 0657 realT *upper_bound; /* scale point[k] to new upper bound */ 0658 realT *lower_bound; /* scale point[k] to new lower bound 0659 project if both upper_ and lower_bound == 0 */ 0660 0661 /*-<a href="qh-globa.htm#TOC" 0662 >--------------------------------</a><a name="qh-prec">-</a> 0663 0664 qh precision constants 0665 precision constants for Qhull 0666 0667 notes: 0668 qh_detroundoff [geom2.c] computes the maximum roundoff error for distance 0669 and other computations. It also sets default values for the 0670 qh constants above. 0671 */ 0672 realT ANGLEround; /* max round off error for angles */ 0673 realT centrum_radius; /* max centrum radius for convexity ('Cn' + 2*qh.DISTround) */ 0674 realT cos_max; /* max cosine for convexity (roundoff added) */ 0675 realT DISTround; /* max round off error for distances, qh.SETroundoff ('En') overrides qh_distround */ 0676 realT MAXabs_coord; /* max absolute coordinate */ 0677 realT MAXlastcoord; /* max last coordinate for qh_scalelast */ 0678 realT MAXoutside; /* max target for qh.max_outside/f.maxoutside, base for qh_RATIO... 0679 recomputed at qh_addpoint, unrelated to qh_MAXoutside */ 0680 realT MAXsumcoord; /* max sum of coordinates */ 0681 realT MAXwidth; /* max rectilinear width of point coordinates */ 0682 realT MINdenom_1; /* min. abs. value for 1/x */ 0683 realT MINdenom; /* use divzero if denominator < MINdenom */ 0684 realT MINdenom_1_2; /* min. abs. val for 1/x that allows normalization */ 0685 realT MINdenom_2; /* use divzero if denominator < MINdenom_2 */ 0686 realT MINlastcoord; /* min. last coordinate for qh_scalelast */ 0687 realT *NEARzero; /* hull_dim array for near zero in gausselim */ 0688 realT NEARinside; /* keep points for qh_check_maxout if close to facet */ 0689 realT ONEmerge; /* max distance for merging simplicial facets */ 0690 realT outside_err; /* application's epsilon for coplanar points 0691 qh_check_bestdist() qh_check_points() reports error if point outside */ 0692 realT WIDEfacet; /* size of wide facet for skipping ridge in 0693 area computation and locking centrum */ 0694 boolT NARROWhull; /* set in qh_initialhull if angle < qh_MAXnarrow */ 0695 0696 /*-<a href="qh-globa.htm#TOC" 0697 >--------------------------------</a><a name="qh-codetern">-</a> 0698 0699 qh internal constants 0700 internal constants for Qhull 0701 */ 0702 char qhull[sizeof("qhull")]; /* "qhull" for checking ownership while debugging */ 0703 jmp_buf errexit; /* exit label for qh_errexit, defined by setjmp() and NOerrexit */ 0704 char jmpXtra[40]; /* extra bytes in case jmp_buf is defined wrong by compiler */ 0705 jmp_buf restartexit; /* restart label for qh_errexit, defined by setjmp() and ALLOWrestart */ 0706 char jmpXtra2[40]; /* extra bytes in case jmp_buf is defined wrong by compiler*/ 0707 FILE * fin; /* pointer to input file, init by qh_initqhull_start2 */ 0708 FILE * fout; /* pointer to output file */ 0709 FILE * ferr; /* pointer to error file */ 0710 pointT *interior_point; /* center point of the initial simplex*/ 0711 int normal_size; /* size in bytes for facet normals and point coords */ 0712 int center_size; /* size in bytes for Voronoi centers */ 0713 int TEMPsize; /* size for small, temporary sets (in quick mem) */ 0714 0715 /*-<a href="qh-globa.htm#TOC" 0716 >--------------------------------</a><a name="qh-lists">-</a> 0717 0718 qh facet and vertex lists 0719 defines lists of facets, new facets, visible facets, vertices, and 0720 new vertices. Includes counts, next ids, and trace ids. 0721 see: 0722 qh_resetlists() 0723 */ 0724 facetT *facet_list; /* first facet */ 0725 facetT *facet_tail; /* end of facet_list (dummy facet with id 0 and next==NULL) */ 0726 facetT *facet_next; /* next facet for buildhull() 0727 previous facets do not have outside sets 0728 NARROWhull: previous facets may have coplanar outside sets for qh_outcoplanar */ 0729 facetT *newfacet_list; /* list of new facets to end of facet_list 0730 qh_postmerge sets newfacet_list to facet_list */ 0731 facetT *visible_list; /* list of visible facets preceding newfacet_list, 0732 end of visible list if !facet->visible, same as newfacet_list 0733 qh_findhorizon sets visible_list at end of facet_list 0734 qh_willdelete prepends to visible_list 0735 qh_triangulate appends mirror facets to visible_list at end of facet_list 0736 qh_postmerge sets visible_list to facet_list 0737 qh_deletevisible deletes the visible facets */ 0738 int num_visible; /* current number of visible facets */ 0739 unsigned int tracefacet_id; /* set at init, then can print whenever */ 0740 facetT *tracefacet; /* set in newfacet/mergefacet, undone in delfacet and qh_errexit */ 0741 unsigned int traceridge_id; /* set at init, then can print whenever */ 0742 ridgeT *traceridge; /* set in newridge, undone in delridge, errexit, errexit2, makenew_nonsimplicial, mergecycle_ridges */ 0743 unsigned int tracevertex_id; /* set at buildtracing, can print whenever */ 0744 vertexT *tracevertex; /* set in newvertex, undone in delvertex and qh_errexit */ 0745 vertexT *vertex_list; /* list of all vertices, to vertex_tail */ 0746 vertexT *vertex_tail; /* end of vertex_list (dummy vertex with ID 0, next NULL) */ 0747 vertexT *newvertex_list; /* list of vertices in newfacet_list, to vertex_tail 0748 all vertices have 'newfacet' set */ 0749 int num_facets; /* number of facets in facet_list 0750 includes visible faces (num_visible) */ 0751 int num_vertices; /* number of vertices in facet_list */ 0752 int num_outside; /* number of points in outsidesets (for tracing and RANDOMoutside) 0753 includes coplanar outsideset points for NARROWhull/qh_outcoplanar() */ 0754 int num_good; /* number of good facets (after qh_findgood_all or qh_markkeep) */ 0755 unsigned int facet_id; /* ID of next, new facet from newfacet() */ 0756 unsigned int ridge_id; /* ID of next, new ridge from newridge() */ 0757 unsigned int vertex_id; /* ID of next, new vertex from newvertex() */ 0758 unsigned int first_newfacet; /* ID of first_newfacet for qh_buildcone, or 0 if none */ 0759 0760 /*-<a href="qh-globa.htm#TOC" 0761 >--------------------------------</a><a name="qh-var">-</a> 0762 0763 qh global variables 0764 defines minimum and maximum distances, next visit ids, several flags, 0765 and other global variables. 0766 initialize in qh_initbuild or qh_maxmin if used in qh_buildhull 0767 */ 0768 unsigned long hulltime; /* ignore time to set up input and randomize */ 0769 /* use 'unsigned long' to avoid wrap-around errors */ 0770 boolT ALLOWrestart; /* true if qh_joggle_restart can use qh.restartexit */ 0771 int build_cnt; /* number of calls to qh_initbuild */ 0772 qh_CENTER CENTERtype; /* current type of facet->center, qh_CENTER */ 0773 int furthest_id; /* pointid of furthest point, for tracing */ 0774 int last_errcode; /* last errcode from qh_fprintf, reset in qh_build_withrestart */ 0775 facetT *GOODclosest; /* closest facet to GOODthreshold in qh_findgood */ 0776 pointT *coplanar_apex; /* last apex declared a coplanar point by qh_getpinchedmerges, prevents infinite loop */ 0777 boolT hasAreaVolume; /* true if totarea, totvol was defined by qh_getarea */ 0778 boolT hasTriangulation; /* true if triangulation created by qh_triangulate */ 0779 boolT isRenameVertex; /* true during qh_merge_pinchedvertices, disables duplicate ridge vertices in qh_checkfacet */ 0780 realT JOGGLEmax; /* set 'QJn' if randomly joggle input. 'QJ'/'QJ0.0' sets default (qh_detjoggle) */ 0781 boolT maxoutdone; /* set qh_check_maxout(), cleared by qh_addpoint() */ 0782 realT max_outside; /* maximum distance from a point to a facet, 0783 before roundoff, not simplicial vertices 0784 actual outer plane is +DISTround and 0785 computed outer plane is +2*DISTround */ 0786 realT max_vertex; /* maximum distance (>0) from vertex to a facet, 0787 before roundoff, due to a merge */ 0788 realT min_vertex; /* minimum distance (<0) from vertex to a facet, 0789 before roundoff, due to a merge 0790 if qh.JOGGLEmax, qh_makenewplanes sets it 0791 recomputed if qh.DOcheckmax, default -qh.DISTround */ 0792 boolT NEWfacets; /* true while visible facets invalid due to new or merge 0793 from qh_makecone/qh_attachnewfacets to qh_resetlists */ 0794 boolT NEWtentative; /* true while new facets are tentative due to !qh.IGNOREpinched or qh.ONLYgood 0795 from qh_makecone to qh_attachnewfacets */ 0796 boolT findbestnew; /* true if partitioning calls qh_findbestnew */ 0797 boolT findbest_notsharp; /* true if new facets are at least 90 degrees */ 0798 boolT NOerrexit; /* true if qh.errexit is not available, cleared after setjmp. See qh.ERREXITcalled */ 0799 realT PRINTcradius; /* radius for printing centrums */ 0800 realT PRINTradius; /* radius for printing vertex spheres and points */ 0801 boolT POSTmerging; /* true when post merging */ 0802 int printoutvar; /* temporary variable for qh_printbegin, etc. */ 0803 int printoutnum; /* number of facets printed */ 0804 unsigned int repart_facetid; /* previous facetid to prevent recursive qh_partitioncoplanar+qh_partitionpoint */ 0805 int retry_addpoint; /* number of retries of qh_addpoint due to merging pinched vertices */ 0806 boolT QHULLfinished; /* True after qhull() is finished */ 0807 realT totarea; /* 'FA': total facet area computed by qh_getarea, hasAreaVolume */ 0808 realT totvol; /* 'FA': total volume computed by qh_getarea, hasAreaVolume */ 0809 unsigned int visit_id; /* unique ID for searching neighborhoods, */ 0810 unsigned int vertex_visit; /* unique ID for searching vertices, reset with qh_buildtracing */ 0811 boolT WAScoplanar; /* True if qh_partitioncoplanar (qh_check_maxout) */ 0812 boolT ZEROall_ok; /* True if qh_checkzero always succeeds */ 0813 0814 /*-<a href="qh-globa.htm#TOC" 0815 >--------------------------------</a><a name="qh-set">-</a> 0816 0817 qh global sets 0818 defines sets for merging, initial simplex, hashing, extra input points, 0819 and deleted vertices 0820 */ 0821 setT *facet_mergeset; /* temporary set of merges to be done */ 0822 setT *degen_mergeset; /* temporary set of degenerate and redundant merges */ 0823 setT *vertex_mergeset; /* temporary set of vertex merges */ 0824 setT *hash_table; /* hash table for matching ridges in qh_matchfacets 0825 size is setsize() */ 0826 setT *other_points; /* additional points */ 0827 setT *del_vertices; /* vertices to partition and delete with visible 0828 facets. v.deleted is set for checkfacet */ 0829 0830 /*-<a href="qh-globa.htm#TOC" 0831 >--------------------------------</a><a name="qh-buf">-</a> 0832 0833 qh global buffers 0834 defines buffers for maxtrix operations, input, and error messages 0835 */ 0836 coordT *gm_matrix; /* (dim+1)Xdim matrix for geom.c */ 0837 coordT **gm_row; /* array of gm_matrix rows */ 0838 char* line; /* malloc'd input line of maxline+1 chars */ 0839 int maxline; 0840 coordT *half_space; /* malloc'd input array for halfspace (qh.normal_size+coordT) */ 0841 coordT *temp_malloc; /* malloc'd input array for points */ 0842 0843 /*-<a href="qh-globa.htm#TOC" 0844 >--------------------------------</a><a name="qh-static">-</a> 0845 0846 qh static variables 0847 defines static variables for individual functions 0848 0849 notes: 0850 do not use 'static' within a function. Multiple instances of qhull 0851 may exist. 0852 0853 do not assume zero initialization, 'QPn' may cause a restart 0854 */ 0855 boolT ERREXITcalled; /* true during qh_errexit (prevents duplicate calls). see qh.NOerrexit */ 0856 boolT firstcentrum; /* for qh_printcentrum */ 0857 boolT old_randomdist; /* save RANDOMdist flag during io, tracing, or statistics */ 0858 setT *coplanarfacetset; /* set of coplanar facets for searching qh_findbesthorizon() */ 0859 realT last_low; /* qh_scalelast parameters for qh_setdelaunay */ 0860 realT last_high; 0861 realT last_newhigh; 0862 realT lastcpu; /* for qh_buildtracing */ 0863 int lastfacets; /* last qh.num_facets */ 0864 int lastmerges; /* last zzval_(Ztotmerge) */ 0865 int lastplanes; /* last zzval_(Zsetplane) */ 0866 int lastdist; /* last zzval_(Zdistplane) */ 0867 unsigned int lastreport; /* last qh.facet_id */ 0868 int mergereport; /* for qh_tracemerging */ 0869 qhstatT *old_qhstat; /* for saving qh_qhstat in save_qhull() and UsingLibQhull. Free with qh_free() */ 0870 setT *old_tempstack; /* for saving qhmem.tempstack in save_qhull */ 0871 int ridgeoutnum; /* number of ridges for 4OFF output (qh_printbegin,etc) */ 0872 }; 0873 0874 /*=========== -macros- =========================*/ 0875 0876 /*-<a href="qh-poly.htm#TOC" 0877 >--------------------------------</a><a name="otherfacet_">-</a> 0878 0879 otherfacet_(ridge, facet) 0880 return neighboring facet for a ridge in facet 0881 */ 0882 #define otherfacet_(ridge, facet) \ 0883 (((ridge)->top == (facet)) ? (ridge)->bottom : (ridge)->top) 0884 0885 /*-<a href="qh-poly.htm#TOC" 0886 >--------------------------------</a><a name="getid_">-</a> 0887 0888 getid_(p) 0889 return int ID for facet, ridge, or vertex 0890 return qh_IDunknown(-1) if NULL 0891 return 0 if facet_tail or vertex_tail 0892 */ 0893 #define getid_(p) ((p) ? (int)((p)->id) : qh_IDunknown) 0894 0895 /*============== FORALL macros ===================*/ 0896 0897 /*-<a href="qh-poly.htm#TOC" 0898 >--------------------------------</a><a name="FORALLfacets">-</a> 0899 0900 FORALLfacets { ... } 0901 assign 'facet' to each facet in qh.facet_list 0902 0903 notes: 0904 uses 'facetT *facet;' 0905 assumes last facet is a sentinel 0906 0907 see: 0908 FORALLfacet_( facetlist ) 0909 */ 0910 #define FORALLfacets for (facet=qh facet_list;facet && facet->next;facet=facet->next) 0911 0912 /*-<a href="qh-poly.htm#TOC" 0913 >--------------------------------</a><a name="FORALLpoints">-</a> 0914 0915 FORALLpoints { ... } 0916 assign 'point' to each point in qh.first_point, qh.num_points 0917 0918 notes: 0919 0920 declare: 0921 coordT *point, *pointtemp; 0922 */ 0923 #define FORALLpoints FORALLpoint_(qh first_point, qh num_points) 0924 0925 /*-<a href="qh-poly.htm#TOC" 0926 >--------------------------------</a><a name="FORALLpoint_">-</a> 0927 0928 FORALLpoint_(points, num) { ... } 0929 assign 'point' to each point in points array of num points 0930 0931 declare: 0932 coordT *point, *pointtemp; 0933 */ 0934 #define FORALLpoint_(points, num) for (point=(points), \ 0935 pointtemp= (points)+qh hull_dim*(num); point < pointtemp; point += qh hull_dim) 0936 0937 /*-<a href="qh-poly.htm#TOC" 0938 >--------------------------------</a><a name="FORALLvertices">-</a> 0939 0940 FORALLvertices { ... } 0941 assign 'vertex' to each vertex in qh.vertex_list 0942 0943 declare: 0944 vertexT *vertex; 0945 0946 notes: 0947 assumes qh.vertex_list terminated by NULL or a sentinel (v.next==NULL) 0948 */ 0949 #define FORALLvertices for (vertex=qh vertex_list;vertex && vertex->next;vertex= vertex->next) 0950 0951 /*-<a href="qh-poly.htm#TOC" 0952 >--------------------------------</a><a name="FOREACHfacet_">-</a> 0953 0954 FOREACHfacet_( facets ) { ... } 0955 assign 'facet' to each facet in facets 0956 0957 declare: 0958 facetT *facet, **facetp; 0959 0960 notes: 0961 assumes set is not modified 0962 0963 see: 0964 <a href="qset.h#FOREACHsetelement_">FOREACHsetelement_</a> 0965 */ 0966 #define FOREACHfacet_(facets) FOREACHsetelement_(facetT, facets, facet) 0967 0968 /*-<a href="qh-poly.htm#TOC" 0969 >--------------------------------</a><a name="FOREACHneighbor_">-</a> 0970 0971 FOREACHneighbor_( facet ) { ... } 0972 assign 'neighbor' to each neighbor in facet->neighbors 0973 0974 FOREACHneighbor_( vertex ) { ... } 0975 assign 'neighbor' to each neighbor in vertex->neighbors 0976 0977 declare: 0978 facetT *neighbor, **neighborp; 0979 0980 notes: 0981 assumes set is not modified 0982 0983 see: 0984 <a href="qset.h#FOREACHsetelement_">FOREACHsetelement_</a> 0985 */ 0986 #define FOREACHneighbor_(facet) FOREACHsetelement_(facetT, facet->neighbors, neighbor) 0987 0988 /*-<a href="qh-poly.htm#TOC" 0989 >--------------------------------</a><a name="FOREACHpoint_">-</a> 0990 0991 FOREACHpoint_( points ) { ... } 0992 assign 'point' to each point in points set 0993 0994 declare: 0995 pointT *point, **pointp; 0996 0997 notes: 0998 assumes set is not modified 0999 1000 see: 1001 <a href="qset.h#FOREACHsetelement_">FOREACHsetelement_</a> 1002 */ 1003 #define FOREACHpoint_(points) FOREACHsetelement_(pointT, points, point) 1004 1005 /*-<a href="qh-poly.htm#TOC" 1006 >--------------------------------</a><a name="FOREACHridge_">-</a> 1007 1008 FOREACHridge_( ridges ) { ... } 1009 assign 'ridge' to each ridge in ridges set 1010 1011 declare: 1012 ridgeT *ridge, **ridgep; 1013 1014 notes: 1015 assumes set is not modified 1016 1017 see: 1018 <a href="qset.h#FOREACHsetelement_">FOREACHsetelement_</a> 1019 */ 1020 #define FOREACHridge_(ridges) FOREACHsetelement_(ridgeT, ridges, ridge) 1021 1022 /*-<a href="qh-poly.htm#TOC" 1023 >--------------------------------</a><a name="FOREACHvertex_">-</a> 1024 1025 FOREACHvertex_( vertices ) { ... } 1026 assign 'vertex' to each vertex in vertices set 1027 1028 declare: 1029 vertexT *vertex, **vertexp; 1030 1031 notes: 1032 assumes set is not modified 1033 1034 see: 1035 <a href="qset.h#FOREACHsetelement_">FOREACHsetelement_</a> 1036 */ 1037 #define FOREACHvertex_(vertices) FOREACHsetelement_(vertexT, vertices,vertex) 1038 1039 /*-<a href="qh-poly.htm#TOC" 1040 >--------------------------------</a><a name="FOREACHfacet_i_">-</a> 1041 1042 FOREACHfacet_i_( facets ) { ... } 1043 assign 'facet' and 'facet_i' for each facet in facets set 1044 1045 declare: 1046 facetT *facet; 1047 int facet_n, facet_i; 1048 1049 see: 1050 <a href="qset.h#FOREACHsetelement_i_">FOREACHsetelement_i_</a> 1051 */ 1052 #define FOREACHfacet_i_(facets) FOREACHsetelement_i_(facetT, facets, facet) 1053 1054 /*-<a href="qh-poly.htm#TOC" 1055 >--------------------------------</a><a name="FOREACHneighbor_i_">-</a> 1056 1057 FOREACHneighbor_i_( facet ) { ... } 1058 assign 'neighbor' and 'neighbor_i' for each neighbor in facet->neighbors 1059 1060 declare: 1061 facetT *neighbor; 1062 int neighbor_n, neighbor_i; 1063 1064 notes: 1065 see <a href="qset.h#FOREACHsetelement_i_">FOREACHsetelement_i_</a> 1066 for facet neighbors of vertex, need to define a new macro 1067 */ 1068 #define FOREACHneighbor_i_(facet) FOREACHsetelement_i_(facetT, facet->neighbors, neighbor) 1069 1070 /*-<a href="qh-poly.htm#TOC" 1071 >--------------------------------</a><a name="FOREACHpoint_i_">-</a> 1072 1073 FOREACHpoint_i_( points ) { ... } 1074 assign 'point' and 'point_i' for each point in points set 1075 1076 declare: 1077 pointT *point; 1078 int point_n, point_i; 1079 1080 see: 1081 <a href="qset.h#FOREACHsetelement_i_">FOREACHsetelement_i_</a> 1082 */ 1083 #define FOREACHpoint_i_(points) FOREACHsetelement_i_(pointT, points, point) 1084 1085 /*-<a href="qh-poly.htm#TOC" 1086 >--------------------------------</a><a name="FOREACHridge_i_">-</a> 1087 1088 FOREACHridge_i_( ridges ) { ... } 1089 assign 'ridge' and 'ridge_i' for each ridge in ridges set 1090 1091 declare: 1092 ridgeT *ridge; 1093 int ridge_n, ridge_i; 1094 1095 see: 1096 <a href="qset.h#FOREACHsetelement_i_">FOREACHsetelement_i_</a> 1097 */ 1098 #define FOREACHridge_i_(ridges) FOREACHsetelement_i_(ridgeT, ridges, ridge) 1099 1100 /*-<a href="qh-poly.htm#TOC" 1101 >--------------------------------</a><a name="FOREACHvertex_i_">-</a> 1102 1103 FOREACHvertex_i_( vertices ) { ... } 1104 assign 'vertex' and 'vertex_i' for each vertex in vertices set 1105 1106 declare: 1107 vertexT *vertex; 1108 int vertex_n, vertex_i; 1109 1110 see: 1111 <a href="qset.h#FOREACHsetelement_i_">FOREACHsetelement_i_</a> 1112 */ 1113 #define FOREACHvertex_i_(vertices) FOREACHsetelement_i_(vertexT, vertices, vertex) 1114 1115 /********* -libqhull.c prototypes (duplicated from qhull_a.h) **********************/ 1116 1117 void qh_qhull(void); 1118 boolT qh_addpoint(pointT *furthest, facetT *facet, boolT checkdist); 1119 void qh_errexit2(int exitcode, facetT *facet, facetT *otherfacet); 1120 void qh_printsummary(FILE *fp); 1121 1122 /********* -user.c prototypes (alphabetical) **********************/ 1123 1124 void qh_errexit(int exitcode, facetT *facet, ridgeT *ridge); 1125 void qh_errprint(const char* string, facetT *atfacet, facetT *otherfacet, ridgeT *atridge, vertexT *atvertex); 1126 int qh_new_qhull(int dim, int numpoints, coordT *points, boolT ismalloc, 1127 char *qhull_cmd, FILE *outfile, FILE *errfile); 1128 void qh_printfacetlist(facetT *facetlist, setT *facets, boolT printall); 1129 void qh_printhelp_degenerate(FILE *fp); 1130 void qh_printhelp_internal(FILE *fp); 1131 void qh_printhelp_narrowhull(FILE *fp, realT minangle); 1132 void qh_printhelp_singular(FILE *fp); 1133 void qh_printhelp_topology(FILE *fp); 1134 void qh_printhelp_wide(FILE *fp); 1135 void qh_user_memsizes(void); 1136 1137 /********* -usermem.c prototypes (alphabetical) **********************/ 1138 void qh_exit(int exitcode); 1139 void qh_fprintf_stderr(int msgcode, const char *fmt, ... ); 1140 void qh_free(void *mem); 1141 void *qh_malloc(size_t size); 1142 1143 /********* -userprintf.c and userprintf_rbox.c prototypes **********************/ 1144 void qh_fprintf(FILE *fp, int msgcode, const char *fmt, ... ); 1145 void qh_fprintf_rbox(FILE *fp, int msgcode, const char *fmt, ... ); 1146 1147 /***** -geom.c/geom2.c/random.c prototypes (duplicated from geom.h, random.h) ****************/ 1148 1149 facetT *qh_findbest(pointT *point, facetT *startfacet, 1150 boolT bestoutside, boolT newfacets, boolT noupper, 1151 realT *dist, boolT *isoutside, int *numpart); 1152 facetT *qh_findbestnew(pointT *point, facetT *startfacet, 1153 realT *dist, boolT bestoutside, boolT *isoutside, int *numpart); 1154 boolT qh_gram_schmidt(int dim, realT **rows); 1155 void qh_outerinner(facetT *facet, realT *outerplane, realT *innerplane); 1156 void qh_printsummary(FILE *fp); 1157 void qh_projectinput(void); 1158 void qh_randommatrix(realT *buffer, int dim, realT **row); 1159 void qh_rotateinput(realT **rows); 1160 void qh_scaleinput(void); 1161 void qh_setdelaunay(int dim, int count, pointT *points); 1162 coordT *qh_sethalfspace_all(int dim, int count, coordT *halfspaces, pointT *feasible); 1163 1164 /***** -global.c prototypes (alphabetical) ***********************/ 1165 1166 unsigned long qh_clock(void); 1167 void qh_checkflags(char *command, char *hiddenflags); 1168 void qh_clear_outputflags(void); 1169 void qh_freebuffers(void); 1170 void qh_freeqhull(boolT allmem); 1171 void qh_freeqhull2(boolT allmem); 1172 void qh_init_A(FILE *infile, FILE *outfile, FILE *errfile, int argc, char *argv[]); 1173 void qh_init_B(coordT *points, int numpoints, int dim, boolT ismalloc); 1174 void qh_init_qhull_command(int argc, char *argv[]); 1175 void qh_initbuffers(coordT *points, int numpoints, int dim, boolT ismalloc); 1176 void qh_initflags(char *command); 1177 void qh_initqhull_buffers(void); 1178 void qh_initqhull_globals(coordT *points, int numpoints, int dim, boolT ismalloc); 1179 void qh_initqhull_mem(void); 1180 void qh_initqhull_outputflags(void); 1181 void qh_initqhull_start(FILE *infile, FILE *outfile, FILE *errfile); 1182 void qh_initqhull_start2(FILE *infile, FILE *outfile, FILE *errfile); 1183 void qh_initthresholds(char *command); 1184 void qh_lib_check(int qhullLibraryType, int qhTsize, int vertexTsize, int ridgeTsize, int facetTsize, int setTsize, int qhmemTsize); 1185 void qh_option(const char *option, int *i, realT *r); 1186 #if qh_QHpointer 1187 void qh_restore_qhull(qhT **oldqh); 1188 qhT *qh_save_qhull(void); 1189 #endif 1190 1191 /***** -io.c prototypes (duplicated from io.h) ***********************/ 1192 1193 void qh_dfacet(unsigned int id); 1194 void qh_dvertex(unsigned int id); 1195 void qh_printneighborhood(FILE *fp, qh_PRINT format, facetT *facetA, facetT *facetB, boolT printall); 1196 void qh_produce_output(void); 1197 coordT *qh_readpoints(int *numpoints, int *dimension, boolT *ismalloc); 1198 1199 1200 /********* -mem.c prototypes (duplicated from mem.h) **********************/ 1201 1202 void qh_meminit(FILE *ferr); 1203 void qh_memfreeshort(int *curlong, int *totlong); 1204 1205 /********* -poly.c/poly2.c prototypes (duplicated from poly.h) **********************/ 1206 1207 void qh_check_output(void); 1208 void qh_check_points(void); 1209 setT *qh_facetvertices(facetT *facetlist, setT *facets, boolT allfacets); 1210 facetT *qh_findbestfacet(pointT *point, boolT bestoutside, 1211 realT *bestdist, boolT *isoutside); 1212 vertexT *qh_nearvertex(facetT *facet, pointT *point, realT *bestdistp); 1213 pointT *qh_point(int id); 1214 setT *qh_pointfacet(void /* qh.facet_list */); 1215 int qh_pointid(pointT *point); 1216 setT *qh_pointvertex(void /* qh.facet_list */); 1217 void qh_setvoronoi_all(void); 1218 void qh_triangulate(void /* qh.facet_list */); 1219 1220 /********* -rboxlib.c prototypes **********************/ 1221 int qh_rboxpoints(FILE* fout, FILE* ferr, char* rbox_command); 1222 void qh_errexit_rbox(int exitcode); 1223 1224 /********* -stat.c prototypes (duplicated from stat.h) **********************/ 1225 1226 void qh_collectstatistics(void); 1227 void qh_printallstatistics(FILE *fp, const char *string); 1228 1229 #endif /* qhDEFlibqhull */
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |