00001 /*<html><pre> -<a href="qh-c.htm#user" 00002 >-------------------------------</a><a name="TOP">-</a> 00003 00004 user.h 00005 user redefinable constants 00006 00007 see qh-c.htm. see COPYING for copyright information. 00008 00009 before reading any code, review qhull.h for data structure definitions and 00010 the "qh" macro. 00011 */ 00012 00013 #ifndef qhDEFuser 00014 #define qhDEFuser 1 00015 00016 /*============= data types and configuration macros ==========*/ 00017 00018 /*-<a href="qh-c.htm#user" 00019 >--------------------------------</a><a name="realT">-</a> 00020 00021 realT 00022 set the size of floating point numbers 00023 00024 qh_REALdigits 00025 maximimum number of significant digits 00026 00027 qh_REAL_1, qh_REAL_2n, qh_REAL_3n 00028 format strings for printf 00029 00030 qh_REALmax, qh_REALmin 00031 maximum and minimum (near zero) values 00032 00033 qh_REALepsilon 00034 machine roundoff 00035 00036 notes: 00037 Select whether to store floating point numbers in single precision (float) 00038 or double precision (double). 00039 00040 Use 'float' to save about 8% in time and 25% in space. This is particularly 00041 help if high-d where convex hulls are space limited. Using 'float' also 00042 reduces the printed size of Qhull's output since numbers have 8 digits of 00043 precision. 00044 00045 Use 'double' when greater arithmetic precision is needed. This is needed 00046 for Delaunay triangulations and Voronoi diagrams when you are not merging 00047 facets. 00048 00049 If 'double' gives insufficient precision, your data probably includes 00050 degeneracies. If so you should use facet merging (done by default) 00051 or exact arithmetic (see imprecision section of manual, qh-impre.htm). 00052 You may also use option 'Po' to force output despite precision errors. 00053 00054 You may use 'long double', but many format statements need to be changed 00055 and you may need a 'long double' square root routine. S. Grundmann 00056 (sg@eeiwzb.et.tu-dresden.de) has done this. He reports that the code runs 00057 much slower with little gain in precision. 00058 00059 WARNING: on some machines, int f(){realT a= REALmax;return (a == REALmax);} 00060 returns False. Use (a > REALmax/2) instead of (a == REALmax). 00061 00062 REALfloat = 1 all numbers are 'float' type 00063 = 0 all numbers are 'double' type 00064 */ 00065 #define REALfloat 0 00066 00067 #if (REALfloat == 1) 00068 #define realT float 00069 #define REALmax FLT_MAX 00070 #define REALmin FLT_MIN 00071 #define REALepsilon FLT_EPSILON 00072 #define qh_REALdigits 8 /* maximum number of significant digits */ 00073 #define qh_REAL_1 "%6.8g " 00074 #define qh_REAL_2n "%6.8g %6.8g\n" 00075 #define qh_REAL_3n "%6.8g %6.8g %6.8g\n" 00076 00077 #elif (REALfloat == 0) 00078 #define realT double 00079 #define REALmax DBL_MAX 00080 #define REALmin DBL_MIN 00081 #define REALepsilon DBL_EPSILON 00082 #define qh_REALdigits 16 /* maximum number of significant digits */ 00083 #define qh_REAL_1 "%6.16g " 00084 #define qh_REAL_2n "%6.16g %6.16g\n" 00085 #define qh_REAL_3n "%6.16g %6.16g %6.16g\n" 00086 00087 #else 00088 #error unknown float option 00089 #endif 00090 00091 /*-<a href="qh-c.htm#user" 00092 >--------------------------------</a><a name="CPUclock">-</a> 00093 00094 qh_CPUclock 00095 define the clock() function for reporting the total time spent by Qhull 00096 returns CPU ticks as a 'long int' 00097 qh_CPUclock is only used for reporting the total time spent by Qhull 00098 00099 qh_SECticks 00100 the number of clock ticks per second 00101 00102 notes: 00103 looks for CLOCKS_PER_SEC, CLOCKS_PER_SECOND, or assumes microseconds 00104 to define a custom clock, set qh_CLOCKtype to 0 00105 00106 if your system does not use clock() to return CPU ticks, replace 00107 qh_CPUclock with the corresponding function. It is converted 00108 to unsigned long to prevent wrap-around during long runs. 00109 00110 00111 Set qh_CLOCKtype to 00112 00113 1 for CLOCKS_PER_SEC, CLOCKS_PER_SECOND, or microsecond 00114 */ 00115 #define qh_CLOCKtype 1 /* change to the desired number */ 00116 00117 #if (qh_CLOCKtype == 1) 00118 00119 #if defined (CLOCKS_PER_SECOND) 00120 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock */ 00121 #define qh_SECticks CLOCKS_PER_SECOND 00122 00123 #elif defined (CLOCKS_PER_SEC) 00124 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock */ 00125 #define qh_SECticks CLOCKS_PER_SEC 00126 00127 #elif defined (CLK_TCK) 00128 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock */ 00129 #define qh_SECticks CLK_TCK 00130 00131 #else 00132 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock */ 00133 #define qh_SECticks 1E6 00134 #endif 00135 00136 #else /* qh_CLOCKtype == ? */ 00137 #error unknown clock option 00138 #endif 00139 00140 /*-<a href="qh-c.htm#user" 00141 >--------------------------------</a><a name="RANDOM">-</a> 00142 00143 qh_RANDOMtype, qh_RANDOMmax, qh_RANDOMseed 00144 define random number generator 00145 00146 qh_RANDOMint generates a random integer between 0 and qh_RANDOMmax. 00147 qh_RANDOMseed sets the random number seed for qh_RANDOMint 00148 00149 Set qh_RANDOMtype (default 5) to: 00150 1 for random() with 31 bits (UCB) 00151 2 for rand() with RAND_MAX or 15 bits (system 5) 00152 3 for rand() with 31 bits (Sun) 00153 4 for lrand48() with 31 bits (Solaris) 00154 5 for qh_rand() with 31 bits (included with Qhull) 00155 00156 notes: 00157 Random numbers are used by rbox to generate point sets. Random 00158 numbers are used by Qhull to rotate the input ('QRn' option), 00159 simulate a randomized algorithm ('Qr' option), and to simulate 00160 roundoff errors ('Rn' option). 00161 00162 Random number generators differ between systems. Most systems provide 00163 rand() but the period varies. The period of rand() is not critical 00164 since qhull does not normally use random numbers. 00165 00166 The default generator is Park & Miller's minimal standard random 00167 number generator [CACM 31:1195 '88]. It is included with Qhull. 00168 00169 If qh_RANDOMmax is wrong, qhull will report a warning and Geomview 00170 output will likely be invisible. 00171 */ 00172 #define qh_RANDOMtype 5 /* *** change to the desired number *** */ 00173 00174 #if (qh_RANDOMtype == 1) 00175 #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, random()/MAX */ 00176 #define qh_RANDOMint random() 00177 #define qh_RANDOMseed_(seed) srandom(seed); 00178 00179 #elif (qh_RANDOMtype == 2) 00180 #ifdef RAND_MAX 00181 #define qh_RANDOMmax ((realT)RAND_MAX) 00182 #else 00183 #define qh_RANDOMmax ((realT)32767) /* 15 bits (System 5) */ 00184 #endif 00185 #define qh_RANDOMint rand() 00186 #define qh_RANDOMseed_(seed) srand((unsigned)seed); 00187 00188 #elif (qh_RANDOMtype == 3) 00189 #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, Sun */ 00190 #define qh_RANDOMint rand() 00191 #define qh_RANDOMseed_(seed) srand((unsigned)seed); 00192 00193 #elif (qh_RANDOMtype == 4) 00194 #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, lrand38()/MAX */ 00195 #define qh_RANDOMint lrand48() 00196 #define qh_RANDOMseed_(seed) srand48(seed); 00197 00198 #elif (qh_RANDOMtype == 5) 00199 #define qh_RANDOMmax ((realT)2147483646UL) /* 31 bits, qh_rand/MAX */ 00200 #define qh_RANDOMint qh_rand() 00201 #define qh_RANDOMseed_(seed) qh_srand(seed); 00202 /* unlike rand(), never returns 0 */ 00203 00204 #else 00205 #error: unknown random option 00206 #endif 00207 00208 /*-<a href="qh-c.htm#user" 00209 >--------------------------------</a><a name="ORIENTclock">-</a> 00210 00211 qh_ORIENTclock 00212 0 for inward pointing normals by Geomview convention 00213 */ 00214 #define qh_ORIENTclock 0 00215 00216 00217 /*========= performance related constants =========*/ 00218 00219 /*-<a href="qh-c.htm#user" 00220 >--------------------------------</a><a name="HASHfactor">-</a> 00221 00222 qh_HASHfactor 00223 total hash slots / used hash slots. Must be at least 1.1. 00224 00225 notes: 00226 =2 for at worst 50% occupancy for qh hash_table and normally 25% occupancy 00227 */ 00228 #define qh_HASHfactor 2 00229 00230 /*-<a href="qh-c.htm#user" 00231 >--------------------------------</a><a name="VERIFYdirect">-</a> 00232 00233 qh_VERIFYdirect 00234 with 'Tv' verify all points against all facets if op count is smaller 00235 00236 notes: 00237 if greater, calls qh_check_bestdist() instead 00238 */ 00239 #define qh_VERIFYdirect 1000000 00240 00241 /*-<a href="qh-c.htm#user" 00242 >--------------------------------</a><a name="INITIALsearch">-</a> 00243 00244 qh_INITIALsearch 00245 if qh_INITIALmax, search points up to this dimension 00246 */ 00247 #define qh_INITIALsearch 6 00248 00249 /*-<a href="qh-c.htm#user" 00250 >--------------------------------</a><a name="INITIALmax">-</a> 00251 00252 qh_INITIALmax 00253 if dim >= qh_INITIALmax, use min/max coordinate points for initial simplex 00254 00255 notes: 00256 from points with non-zero determinants 00257 use option 'Qs' to override (much slower) 00258 */ 00259 #define qh_INITIALmax 8 00260 00261 /*-<a href="qh-c.htm#user" 00262 >--------------------------------</a><a name="JOGGLEdefault">-</a> 00263 00264 qh_JOGGLEdefault 00265 default qh.JOGGLEmax is qh.DISTround * qh_JOGGLEdefault 00266 00267 notes: 00268 rbox s r 100 | qhull QJ1e-15 QR0 generates 90% faults at distround 7e-16 00269 rbox s r 100 | qhull QJ1e-14 QR0 generates 70% faults 00270 rbox s r 100 | qhull QJ1e-13 QR0 generates 35% faults 00271 rbox s r 100 | qhull QJ1e-12 QR0 generates 8% faults 00272 rbox s r 100 | qhull QJ1e-11 QR0 generates 1% faults 00273 rbox s r 100 | qhull QJ1e-10 QR0 generates 0% faults 00274 rbox 1000 W0 | qhull QJ1e-12 QR0 generates 86% faults 00275 rbox 1000 W0 | qhull QJ1e-11 QR0 generates 20% faults 00276 rbox 1000 W0 | qhull QJ1e-10 QR0 generates 2% faults 00277 the later have about 20 points per facet, each of which may interfere 00278 00279 pick a value large enough to avoid retries on most inputs 00280 */ 00281 #define qh_JOGGLEdefault 30000.0 00282 00283 /*-<a href="qh-c.htm#user" 00284 >--------------------------------</a><a name="JOGGLEincrease">-</a> 00285 00286 qh_JOGGLEincrease 00287 factor to increase qh.JOGGLEmax on qh_JOGGLEretry or qh_JOGGLEagain 00288 */ 00289 #define qh_JOGGLEincrease 10.0 00290 00291 /*-<a href="qh-c.htm#user" 00292 >--------------------------------</a><a name="JOGGLEretry">-</a> 00293 00294 qh_JOGGLEretry 00295 if ZZretry = qh_JOGGLEretry, increase qh.JOGGLEmax 00296 00297 notes: 00298 try twice at the original value in case of bad luck the first time 00299 */ 00300 #define qh_JOGGLEretry 2 00301 00302 /*-<a href="qh-c.htm#user" 00303 >--------------------------------</a><a name="JOGGLEagain">-</a> 00304 00305 qh_JOGGLEagain 00306 every following qh_JOGGLEagain, increase qh.JOGGLEmax 00307 00308 notes: 00309 1 is OK since it's already failed qh_JOGGLEretry times 00310 */ 00311 #define qh_JOGGLEagain 1 00312 00313 /*-<a href="qh-c.htm#user" 00314 >--------------------------------</a><a name="JOGGLEmaxincrease">-</a> 00315 00316 qh_JOGGLEmaxincrease 00317 maximum qh.JOGGLEmax due to qh_JOGGLEincrease 00318 relative to qh.MAXwidth 00319 00320 notes: 00321 qh.joggleinput will retry at this value until qh_JOGGLEmaxretry 00322 */ 00323 #define qh_JOGGLEmaxincrease 1e-2 00324 00325 /*-<a href="qh-c.htm#user" 00326 >--------------------------------</a><a name="JOGGLEmaxretry">-</a> 00327 00328 qh_JOGGLEmaxretry 00329 stop after qh_JOGGLEmaxretry attempts 00330 */ 00331 #define qh_JOGGLEmaxretry 100 00332 00333 /*========= memory constants =========*/ 00334 00335 /*-<a href="qh-c.htm#user" 00336 >--------------------------------</a><a name="MEMalign">-</a> 00337 00338 qh_MEMalign 00339 memory alignment for qh_meminitbuffers() in global.c 00340 00341 notes: 00342 to avoid bus errors, memory allocation must consider alignment requirements. 00343 malloc() automatically takes care of alignment. Since mem.c manages 00344 its own memory, we need to explicitly specify alignment in 00345 qh_meminitbuffers(). 00346 00347 A safe choice is sizeof(double). sizeof(float) may be used if doubles 00348 do not occur in data structures and pointers are the same size. Be careful 00349 of machines (e.g., DEC Alpha) with large pointers. 00350 00351 If using gcc, best alignment is 00352 #define qh_MEMalign fmax_(__alignof__(realT),__alignof__(void *)) 00353 */ 00354 #define qh_MEMalign fmax_(sizeof(realT), sizeof(void *)) 00355 00356 /*-<a href="qh-c.htm#user" 00357 >--------------------------------</a><a name="MEMbufsize">-</a> 00358 00359 qh_MEMbufsize 00360 size of additional memory buffers 00361 00362 notes: 00363 used for qh_meminitbuffers() in global.c 00364 */ 00365 #define qh_MEMbufsize 0x10000 /* allocate 64K memory buffers */ 00366 00367 /*-<a href="qh-c.htm#user" 00368 >--------------------------------</a><a name="MEMinitbuf">-</a> 00369 00370 qh_MEMinitbuf 00371 size of initial memory buffer 00372 00373 notes: 00374 use for qh_meminitbuffers() in global.c 00375 */ 00376 #define qh_MEMinitbuf 0x20000 /* initially allocate 128K buffer */ 00377 00378 /*-<a href="qh-c.htm#user" 00379 >--------------------------------</a><a name="INFINITE">-</a> 00380 00381 qh_INFINITE 00382 on output, indicates Voronoi center at infinity 00383 */ 00384 #define qh_INFINITE -10.101 00385 00386 /*-<a href="qh-c.htm#user" 00387 >--------------------------------</a><a name="DEFAULTbox">-</a> 00388 00389 qh_DEFAULTbox 00390 default box size (Geomview expects 0.5) 00391 */ 00392 #define qh_DEFAULTbox 0.5 00393 00394 /*======= conditional compilation ============================*/ 00395 00396 /*-<a href="qh-c.htm#user" 00397 >--------------------------------</a><a name="compiler">-</a> 00398 00399 __cplusplus 00400 defined by C++ compilers 00401 00402 __MSC_VER 00403 defined by Microsoft Visual C++ 00404 00405 __MWERKS__ && __POWERPC__ 00406 defined by Metrowerks when compiling for the Power Macintosh 00407 00408 __STDC__ 00409 defined for strict ANSI C 00410 */ 00411 00412 /*-<a href="qh-c.htm#user" 00413 >--------------------------------</a><a name="COMPUTEfurthest">-</a> 00414 00415 qh_COMPUTEfurthest 00416 compute furthest distance to an outside point instead of storing it with the facet 00417 =1 to compute furthest 00418 00419 notes: 00420 computing furthest saves memory but costs time 00421 about 40% more distance tests for partitioning 00422 removes facet->furthestdist 00423 */ 00424 #define qh_COMPUTEfurthest 0 00425 00426 /*-<a href="qh-c.htm#user" 00427 >--------------------------------</a><a name="KEEPstatistics">-</a> 00428 00429 qh_KEEPstatistics 00430 =0 removes most of statistic gathering and reporting 00431 00432 notes: 00433 if 0, code size is reduced by about 4%. 00434 */ 00435 #define qh_KEEPstatistics 1 00436 00437 /*-<a href="qh-c.htm#user" 00438 >--------------------------------</a><a name="MAXoutside">-</a> 00439 00440 qh_MAXoutside 00441 record outer plane for each facet 00442 =1 to record facet->maxoutside 00443 00444 notes: 00445 this takes a realT per facet and slightly slows down qhull 00446 it produces better outer planes for geomview output 00447 */ 00448 #define qh_MAXoutside 1 00449 00450 /*-<a href="qh-c.htm#user" 00451 >--------------------------------</a><a name="NOmerge">-</a> 00452 00453 qh_NOmerge 00454 disables facet merging if defined 00455 00456 notes: 00457 This saves about 10% space. 00458 00459 Unless 'Q0' 00460 qh_NOmerge sets 'QJ' to avoid precision errors 00461 00462 #define qh_NOmerge 00463 00464 see: 00465 <a href="mem.h#NOmem">qh_NOmem</a> in mem.c 00466 00467 see user.c/user_eg.c for removing io.o 00468 */ 00469 00470 /*-<a href="qh-c.htm#user" 00471 >--------------------------------</a><a name="NOtrace">-</a> 00472 00473 qh_NOtrace 00474 no tracing if defined 00475 00476 notes: 00477 This saves about 5% space. 00478 00479 #define qh_NOtrace 00480 */ 00481 00482 /*-<a href="qh-c.htm#user" 00483 >--------------------------------</a><a name="QHpointer">-</a> 00484 00485 qh_QHpointer 00486 access global data with pointer or static structure 00487 00488 qh_QHpointer = 1 access globals via a pointer to allocated memory 00489 enables qh_saveqhull() and qh_restoreqhull() 00490 costs about 8% in time and 2% in space 00491 00492 = 0 qh_qh and qh_qhstat are static data structures 00493 only one instance of qhull() can be active at a time 00494 default value 00495 00496 notes: 00497 all global variables for qhull are in qh, qhmem, and qhstat 00498 qh is defined in qhull.h 00499 qhmem is defined in mem.h 00500 qhstat is defined in stat.h 00501 */ 00502 #define qh_QHpointer 0 00503 00504 /*-<a href="qh-c.htm#user" 00505 >--------------------------------</a><a name="QUICKhelp">-</a> 00506 00507 qh_QUICKhelp 00508 =1 to use abbreviated help messages, e.g., for degenerate inputs 00509 */ 00510 #define qh_QUICKhelp 0 00511 00512 /* ============ -merge constants- ==================== 00513 00514 These constants effect facet merging. You probably will not need 00515 to modify these. They effect the performance of facet merging. 00516 */ 00517 00518 /*-<a href="qh-c.htm#user" 00519 >--------------------------------</a><a name="DIMmergeVertex">-</a> 00520 00521 qh_DIMmergeVertex 00522 max dimension for vertex merging (it is not effective in high-d) 00523 */ 00524 #define qh_DIMmergeVertex 6 00525 00526 /*-<a href="qh-c.htm#user" 00527 >--------------------------------</a><a name="DIMreduceBuild">-</a> 00528 00529 qh_DIMreduceBuild 00530 max dimension for vertex reduction during build (slow in high-d) 00531 */ 00532 #define qh_DIMreduceBuild 5 00533 00534 /*-<a href="qh-c.htm#user" 00535 >--------------------------------</a><a name="BESTcentrum">-</a> 00536 00537 qh_BESTcentrum 00538 if > 2*dim+n vertices, qh_findbestneighbor() tests centrums (faster) 00539 else, qh_findbestneighbor() tests all vertices (much better merges) 00540 00541 qh_BESTcentrum2 00542 if qh_BESTcentrum2 * DIM3 + BESTcentrum < #vertices tests centrums 00543 */ 00544 #define qh_BESTcentrum 20 00545 #define qh_BESTcentrum2 2 00546 00547 /*-<a href="qh-c.htm#user" 00548 >--------------------------------</a><a name="BESTnonconvex">-</a> 00549 00550 qh_BESTnonconvex 00551 if > dim+n neighbors, qh_findbestneighbor() tests nonconvex ridges. 00552 00553 notes: 00554 It is needed because qh_findbestneighbor is slow for large facets 00555 */ 00556 #define qh_BESTnonconvex 15 00557 00558 /*-<a href="qh-c.htm#user" 00559 >--------------------------------</a><a name="MAXnewmerges">-</a> 00560 00561 qh_MAXnewmerges 00562 if >n newmerges, qh_merge_nonconvex() calls qh_reducevertices_centrums. 00563 00564 notes: 00565 It is needed because postmerge can merge many facets at once 00566 */ 00567 #define qh_MAXnewmerges 2 00568 00569 /*-<a href="qh-c.htm#user" 00570 >--------------------------------</a><a name="MAXnewcentrum">-</a> 00571 00572 qh_MAXnewcentrum 00573 if <= dim+n vertices (n approximates the number of merges), 00574 reset the centrum in qh_updatetested() and qh_mergecycle_facets() 00575 00576 notes: 00577 needed to reduce cost and because centrums may move too much if 00578 many vertices in high-d 00579 */ 00580 #define qh_MAXnewcentrum 5 00581 00582 /*-<a href="qh-c.htm#user" 00583 >--------------------------------</a><a name="COPLANARratio">-</a> 00584 00585 qh_COPLANARratio 00586 for 3-d+ merging, qh.MINvisible is n*premerge_centrum 00587 00588 notes: 00589 for non-merging, it's DISTround 00590 */ 00591 #define qh_COPLANARratio 3 00592 00593 /*-<a href="qh-c.htm#user" 00594 >--------------------------------</a><a name="DISToutside">-</a> 00595 00596 qh_DISToutside 00597 minimum distance when merging to stop search in qh_findbestnew 00598 or qh_partitionall 00599 00600 notes: 00601 if too big then O(n^2) behavior for partitioning in cone 00602 if very small then important points not processed 00603 */ 00604 #define qh_DISToutside fmax_(4*qh MINoutside, 2*qh max_outside) 00605 00606 /*-<a href="qh-c.htm#user" 00607 >--------------------------------</a><a name="RATIOnearinside">-</a> 00608 00609 qh_RATIOnearinside 00610 ratio of qh.NEARinside to qh.ONEmerge for retaining inside points for 00611 qh_check_maxout(). 00612 00613 notes: 00614 This is overkill since do not know the correct value. 00615 It effects whether 'Qc' reports all coplanar points 00616 */ 00617 #define qh_RATIOnearinside 5 00618 00619 /*-<a href="qh-c.htm#user" 00620 >--------------------------------</a><a name="USEfindbestnew">-</a> 00621 00622 qh_USEfindbestnew 00623 cut-off between qh_findbest/qh_findbestnew in #merged facets. 00624 */ 00625 #define qh_USEfindbestnew 50 00626 00627 /*-<a href="qh-c.htm#user" 00628 >--------------------------------</a><a name="WIDEcoplanar">-</a> 00629 00630 qh_WIDEcoplanar 00631 n*MAXcoplanar or n*MINvisible for a WIDEfacet 00632 00633 if vertex is further than qh.WIDEfacet from the hyperplane 00634 then its ridges are not counted in computing the area, and 00635 the facet's centrum is frozen. 00636 00637 notes: 00638 qh.WIDEfacet= max(qh.MAXoutside,qh_WIDEcoplanar*qh.MAXcoplanar, 00639 qh_WIDEcoplanar * qh.MINvisible); 00640 */ 00641 #define qh_WIDEcoplanar 6 00642 00643 /*-<a href="qh-c.htm#user" 00644 >--------------------------------</a><a name="MAXnarrow">-</a> 00645 00646 qh_MAXnarrow 00647 max. cosine in initial hull that sets qh.NARROWhull 00648 00649 notes: 00650 If qh.NARROWhull, the initial partition does not make 00651 coplanar points. If narrow, a coplanar point can be 00652 coplanar to two facets of opposite orientations and 00653 distant from the exact convex hull. 00654 00655 set to -cos( 0.2 ) for about 1:10 ratio between coplanar point and distance to exact hull 00656 */ 00657 #define qh_MAXnarrow -0.98 00658 00659 /*-<a href="qh-c.htm#user" 00660 >--------------------------------</a><a name="WARNnarrow">-</a> 00661 00662 qh_WARNnarrow 00663 max. cosine in initial hull to warn about qh.NARROWhull 00664 00665 notes: 00666 this is a conservative estimate. 00667 Don't actually see problems until much narrower. 00668 */ 00669 #define qh_WARNnarrow -0.9999 00670 00671 #endif /* qh_DEFuser */ 00672 00673 00674