00001 /*<html><pre> -<a href="qh-user.htm" 00002 >-------------------------------</a><a name="TOP">-</a> 00003 00004 user.h 00005 user redefinable constants 00006 00007 see qh-user.htm. see COPYING for copyright information. 00008 00009 before reading any code, review libqhull.h for data structure definitions and 00010 the "qh" macro. 00011 00012 Sections: 00013 ============= qhull library constants ====================== 00014 ============= data types and configuration macros ========== 00015 ============= performance related constants ================ 00016 ============= memory constants ============================= 00017 ============= joggle constants ============================= 00018 ============= conditional compilation ====================== 00019 ============= -merge constants- ============================ 00020 00021 Code flags -- 00022 NOerrors -- the code does not call qh_errexit() 00023 WARN64 -- the code may be incompatible with 64-bit pointers 00024 00025 */ 00026 00027 #include <time.h> 00028 00029 #ifndef qhDEFuser 00030 #define qhDEFuser 1 00031 00032 /*============================================================*/ 00033 /*============= qhull library constants ======================*/ 00034 /*============================================================*/ 00035 00036 /*-<a href="qh-user.htm#TOC" 00037 >--------------------------------</a><a name="filenamelen">-</a> 00038 00039 FILENAMElen -- max length for TI and TO filenames 00040 00041 */ 00042 00043 #define qh_FILENAMElen 500 00044 00045 /*-<a href="qh-user.htm#TOC" 00046 >--------------------------------</a><a name="msgcode">-</a> 00047 00048 msgcode -- Unique message codes for qh_fprintf 00049 00050 If add new messages, assign these values and increment. 00051 00052 def counters = [27, 1047, 2059, 3025, 4068, 5003, 00053 6241, 7079, 8143, 9410, 10000, 11026] 00054 00055 See: qh_ERR* [libqhull.h] 00056 */ 00057 00058 #define MSG_TRACE0 0 00059 #define MSG_TRACE1 1000 00060 #define MSG_TRACE2 2000 00061 #define MSG_TRACE3 3000 00062 #define MSG_TRACE4 4000 00063 #define MSG_TRACE5 5000 00064 #define MSG_ERROR 6000 /* errors written to qh.ferr */ 00065 #define MSG_WARNING 7000 00066 #define MSG_STDERR 8000 /* log messages Written to qh.ferr */ 00067 #define MSG_OUTPUT 9000 00068 #define MSG_QHULL_ERROR 10000 /* errors thrown by QhullError [QhullError.h] */ 00069 #define MSG_FIXUP 11000 /* FIXUP QH11... */ 00070 #define MSG_MAXLEN 3000 /* qh_printhelp_degenerate() in user.c */ 00071 00072 00073 /*-<a href="qh-user.htm#TOC" 00074 >--------------------------------</a><a name="qh_OPTIONline">-</a> 00075 00076 qh_OPTIONline -- max length of an option line 'FO' 00077 */ 00078 #define qh_OPTIONline 80 00079 00080 /*============================================================*/ 00081 /*============= data types and configuration macros ==========*/ 00082 /*============================================================*/ 00083 00084 /*-<a href="qh-user.htm#TOC" 00085 >--------------------------------</a><a name="realT">-</a> 00086 00087 realT 00088 set the size of floating point numbers 00089 00090 qh_REALdigits 00091 maximimum number of significant digits 00092 00093 qh_REAL_1, qh_REAL_2n, qh_REAL_3n 00094 format strings for printf 00095 00096 qh_REALmax, qh_REALmin 00097 maximum and minimum (near zero) values 00098 00099 qh_REALepsilon 00100 machine roundoff. Maximum roundoff error for addition and multiplication. 00101 00102 notes: 00103 Select whether to store floating point numbers in single precision (float) 00104 or double precision (double). 00105 00106 Use 'float' to save about 8% in time and 25% in space. This is particularly 00107 helpful if high-d where convex hulls are space limited. Using 'float' also 00108 reduces the printed size of Qhull's output since numbers have 8 digits of 00109 precision. 00110 00111 Use 'double' when greater arithmetic precision is needed. This is needed 00112 for Delaunay triangulations and Voronoi diagrams when you are not merging 00113 facets. 00114 00115 If 'double' gives insufficient precision, your data probably includes 00116 degeneracies. If so you should use facet merging (done by default) 00117 or exact arithmetic (see imprecision section of manual, qh-impre.htm). 00118 You may also use option 'Po' to force output despite precision errors. 00119 00120 You may use 'long double', but many format statements need to be changed 00121 and you may need a 'long double' square root routine. S. Grundmann 00122 (sg@eeiwzb.et.tu-dresden.de) has done this. He reports that the code runs 00123 much slower with little gain in precision. 00124 00125 WARNING: on some machines, int f(){realT a= REALmax;return (a == REALmax);} 00126 returns False. Use (a > REALmax/2) instead of (a == REALmax). 00127 00128 REALfloat = 1 all numbers are 'float' type 00129 = 0 all numbers are 'double' type 00130 */ 00131 #define REALfloat 0 00132 00133 #if (REALfloat == 1) 00134 #define realT float 00135 #define REALmax FLT_MAX 00136 #define REALmin FLT_MIN 00137 #define REALepsilon FLT_EPSILON 00138 #define qh_REALdigits 8 /* maximum number of significant digits */ 00139 #define qh_REAL_1 "%6.8g " 00140 #define qh_REAL_2n "%6.8g %6.8g\n" 00141 #define qh_REAL_3n "%6.8g %6.8g %6.8g\n" 00142 00143 #elif (REALfloat == 0) 00144 #define realT double 00145 #define REALmax DBL_MAX 00146 #define REALmin DBL_MIN 00147 #define REALepsilon DBL_EPSILON 00148 #define qh_REALdigits 16 /* maximum number of significant digits */ 00149 #define qh_REAL_1 "%6.16g " 00150 #define qh_REAL_2n "%6.16g %6.16g\n" 00151 #define qh_REAL_3n "%6.16g %6.16g %6.16g\n" 00152 00153 #else 00154 #error unknown float option 00155 #endif 00156 00157 /*-<a href="qh-user.htm#TOC" 00158 >--------------------------------</a><a name="CPUclock">-</a> 00159 00160 qh_CPUclock 00161 define the clock() function for reporting the total time spent by Qhull 00162 returns CPU ticks as a 'long int' 00163 qh_CPUclock is only used for reporting the total time spent by Qhull 00164 00165 qh_SECticks 00166 the number of clock ticks per second 00167 00168 notes: 00169 looks for CLOCKS_PER_SEC, CLOCKS_PER_SECOND, or assumes microseconds 00170 to define a custom clock, set qh_CLOCKtype to 0 00171 00172 if your system does not use clock() to return CPU ticks, replace 00173 qh_CPUclock with the corresponding function. It is converted 00174 to 'unsigned long' to prevent wrap-around during long runs. By default, 00175 <time.h> defines clock_t as 'long' 00176 00177 Set qh_CLOCKtype to 00178 00179 1 for CLOCKS_PER_SEC, CLOCKS_PER_SECOND, or microsecond 00180 Note: may fail if more than 1 hour elapsed time 00181 00182 2 use qh_clock() with POSIX times() (see global.c) 00183 */ 00184 #define qh_CLOCKtype 1 /* change to the desired number */ 00185 00186 #if (qh_CLOCKtype == 1) 00187 00188 #if defined(CLOCKS_PER_SECOND) 00189 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock */ 00190 #define qh_SECticks CLOCKS_PER_SECOND 00191 00192 #elif defined(CLOCKS_PER_SEC) 00193 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock */ 00194 #define qh_SECticks CLOCKS_PER_SEC 00195 00196 #elif defined(CLK_TCK) 00197 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock */ 00198 #define qh_SECticks CLK_TCK 00199 00200 #else 00201 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock */ 00202 #define qh_SECticks 1E6 00203 #endif 00204 00205 #elif (qh_CLOCKtype == 2) 00206 #define qh_CPUclock qh_clock() /* return CPU clock */ 00207 #define qh_SECticks 100 00208 00209 #else /* qh_CLOCKtype == ? */ 00210 #error unknown clock option 00211 #endif 00212 00213 /*-<a href="qh-user.htm#TOC" 00214 >--------------------------------</a><a name="RANDOM">-</a> 00215 00216 qh_RANDOMtype, qh_RANDOMmax, qh_RANDOMseed 00217 define random number generator 00218 00219 qh_RANDOMint generates a random integer between 0 and qh_RANDOMmax. 00220 qh_RANDOMseed sets the random number seed for qh_RANDOMint 00221 00222 Set qh_RANDOMtype (default 5) to: 00223 1 for random() with 31 bits (UCB) 00224 2 for rand() with RAND_MAX or 15 bits (system 5) 00225 3 for rand() with 31 bits (Sun) 00226 4 for lrand48() with 31 bits (Solaris) 00227 5 for qh_rand() with 31 bits (included with Qhull) 00228 00229 notes: 00230 Random numbers are used by rbox to generate point sets. Random 00231 numbers are used by Qhull to rotate the input ('QRn' option), 00232 simulate a randomized algorithm ('Qr' option), and to simulate 00233 roundoff errors ('Rn' option). 00234 00235 Random number generators differ between systems. Most systems provide 00236 rand() but the period varies. The period of rand() is not critical 00237 since qhull does not normally use random numbers. 00238 00239 The default generator is Park & Miller's minimal standard random 00240 number generator [CACM 31:1195 '88]. It is included with Qhull. 00241 00242 If qh_RANDOMmax is wrong, qhull will report a warning and Geomview 00243 output will likely be invisible. 00244 */ 00245 #define qh_RANDOMtype 5 /* *** change to the desired number *** */ 00246 00247 #if (qh_RANDOMtype == 1) 00248 #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, random()/MAX */ 00249 #define qh_RANDOMint random() 00250 #define qh_RANDOMseed_(seed) srandom(seed); 00251 00252 #elif (qh_RANDOMtype == 2) 00253 #ifdef RAND_MAX 00254 #define qh_RANDOMmax ((realT)RAND_MAX) 00255 #else 00256 #define qh_RANDOMmax ((realT)32767) /* 15 bits (System 5) */ 00257 #endif 00258 #define qh_RANDOMint rand() 00259 #define qh_RANDOMseed_(seed) srand((unsigned)seed); 00260 00261 #elif (qh_RANDOMtype == 3) 00262 #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, Sun */ 00263 #define qh_RANDOMint rand() 00264 #define qh_RANDOMseed_(seed) srand((unsigned)seed); 00265 00266 #elif (qh_RANDOMtype == 4) 00267 #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, lrand38()/MAX */ 00268 #define qh_RANDOMint lrand48() 00269 #define qh_RANDOMseed_(seed) srand48(seed); 00270 00271 #elif (qh_RANDOMtype == 5) 00272 #define qh_RANDOMmax ((realT)2147483646UL) /* 31 bits, qh_rand/MAX */ 00273 #define qh_RANDOMint qh_rand() 00274 #define qh_RANDOMseed_(seed) qh_srand(seed); 00275 /* unlike rand(), never returns 0 */ 00276 00277 #else 00278 #error: unknown random option 00279 #endif 00280 00281 /*-<a href="qh-user.htm#TOC" 00282 >--------------------------------</a><a name="ORIENTclock">-</a> 00283 00284 qh_ORIENTclock 00285 0 for inward pointing normals by Geomview convention 00286 */ 00287 #define qh_ORIENTclock 0 00288 00289 00290 /*============================================================*/ 00291 /*============= joggle constants =============================*/ 00292 /*============================================================*/ 00293 00294 /*-<a href="qh-user.htm#TOC" 00295 >--------------------------------</a><a name="JOGGLEdefault">-</a> 00296 00297 qh_JOGGLEdefault 00298 default qh.JOGGLEmax is qh.DISTround * qh_JOGGLEdefault 00299 00300 notes: 00301 rbox s r 100 | qhull QJ1e-15 QR0 generates 90% faults at distround 7e-16 00302 rbox s r 100 | qhull QJ1e-14 QR0 generates 70% faults 00303 rbox s r 100 | qhull QJ1e-13 QR0 generates 35% faults 00304 rbox s r 100 | qhull QJ1e-12 QR0 generates 8% faults 00305 rbox s r 100 | qhull QJ1e-11 QR0 generates 1% faults 00306 rbox s r 100 | qhull QJ1e-10 QR0 generates 0% faults 00307 rbox 1000 W0 | qhull QJ1e-12 QR0 generates 86% faults 00308 rbox 1000 W0 | qhull QJ1e-11 QR0 generates 20% faults 00309 rbox 1000 W0 | qhull QJ1e-10 QR0 generates 2% faults 00310 the later have about 20 points per facet, each of which may interfere 00311 00312 pick a value large enough to avoid retries on most inputs 00313 */ 00314 #define qh_JOGGLEdefault 30000.0 00315 00316 /*-<a href="qh-user.htm#TOC" 00317 >--------------------------------</a><a name="JOGGLEincrease">-</a> 00318 00319 qh_JOGGLEincrease 00320 factor to increase qh.JOGGLEmax on qh_JOGGLEretry or qh_JOGGLEagain 00321 */ 00322 #define qh_JOGGLEincrease 10.0 00323 00324 /*-<a href="qh-user.htm#TOC" 00325 >--------------------------------</a><a name="JOGGLEretry">-</a> 00326 00327 qh_JOGGLEretry 00328 if ZZretry = qh_JOGGLEretry, increase qh.JOGGLEmax 00329 00330 notes: 00331 try twice at the original value in case of bad luck the first time 00332 */ 00333 #define qh_JOGGLEretry 2 00334 00335 /*-<a href="qh-user.htm#TOC" 00336 >--------------------------------</a><a name="JOGGLEagain">-</a> 00337 00338 qh_JOGGLEagain 00339 every following qh_JOGGLEagain, increase qh.JOGGLEmax 00340 00341 notes: 00342 1 is OK since it's already failed qh_JOGGLEretry times 00343 */ 00344 #define qh_JOGGLEagain 1 00345 00346 /*-<a href="qh-user.htm#TOC" 00347 >--------------------------------</a><a name="JOGGLEmaxincrease">-</a> 00348 00349 qh_JOGGLEmaxincrease 00350 maximum qh.JOGGLEmax due to qh_JOGGLEincrease 00351 relative to qh.MAXwidth 00352 00353 notes: 00354 qh.joggleinput will retry at this value until qh_JOGGLEmaxretry 00355 */ 00356 #define qh_JOGGLEmaxincrease 1e-2 00357 00358 /*-<a href="qh-user.htm#TOC" 00359 >--------------------------------</a><a name="JOGGLEmaxretry">-</a> 00360 00361 qh_JOGGLEmaxretry 00362 stop after qh_JOGGLEmaxretry attempts 00363 */ 00364 #define qh_JOGGLEmaxretry 100 00365 00366 /*============================================================*/ 00367 /*============= performance related constants ================*/ 00368 /*============================================================*/ 00369 00370 /*-<a href="qh-user.htm#TOC" 00371 >--------------------------------</a><a name="HASHfactor">-</a> 00372 00373 qh_HASHfactor 00374 total hash slots / used hash slots. Must be at least 1.1. 00375 00376 notes: 00377 =2 for at worst 50% occupancy for qh hash_table and normally 25% occupancy 00378 */ 00379 #define qh_HASHfactor 2 00380 00381 /*-<a href="qh-user.htm#TOC" 00382 >--------------------------------</a><a name="VERIFYdirect">-</a> 00383 00384 qh_VERIFYdirect 00385 with 'Tv' verify all points against all facets if op count is smaller 00386 00387 notes: 00388 if greater, calls qh_check_bestdist() instead 00389 */ 00390 #define qh_VERIFYdirect 1000000 00391 00392 /*-<a href="qh-user.htm#TOC" 00393 >--------------------------------</a><a name="INITIALsearch">-</a> 00394 00395 qh_INITIALsearch 00396 if qh_INITIALmax, search points up to this dimension 00397 */ 00398 #define qh_INITIALsearch 6 00399 00400 /*-<a href="qh-user.htm#TOC" 00401 >--------------------------------</a><a name="INITIALmax">-</a> 00402 00403 qh_INITIALmax 00404 if dim >= qh_INITIALmax, use min/max coordinate points for initial simplex 00405 00406 notes: 00407 from points with non-zero determinants 00408 use option 'Qs' to override (much slower) 00409 */ 00410 #define qh_INITIALmax 8 00411 00412 /*============================================================*/ 00413 /*============= memory constants =============================*/ 00414 /*============================================================*/ 00415 00416 /*-<a href="qh-user.htm#TOC" 00417 >--------------------------------</a><a name="MEMalign">-</a> 00418 00419 qh_MEMalign 00420 memory alignment for qh_meminitbuffers() in global.c 00421 00422 notes: 00423 to avoid bus errors, memory allocation must consider alignment requirements. 00424 malloc() automatically takes care of alignment. Since mem.c manages 00425 its own memory, we need to explicitly specify alignment in 00426 qh_meminitbuffers(). 00427 00428 A safe choice is sizeof(double). sizeof(float) may be used if doubles 00429 do not occur in data structures and pointers are the same size. Be careful 00430 of machines (e.g., DEC Alpha) with large pointers. 00431 00432 If using gcc, best alignment is 00433 #define qh_MEMalign fmax_(__alignof__(realT),__alignof__(void *)) 00434 */ 00435 #define qh_MEMalign ((int)(fmax_(sizeof(realT), sizeof(void *)))) 00436 00437 /*-<a href="qh-user.htm#TOC" 00438 >--------------------------------</a><a name="MEMbufsize">-</a> 00439 00440 qh_MEMbufsize 00441 size of additional memory buffers 00442 00443 notes: 00444 used for qh_meminitbuffers() in global.c 00445 */ 00446 #define qh_MEMbufsize 0x10000 /* allocate 64K memory buffers */ 00447 00448 /*-<a href="qh-user.htm#TOC" 00449 >--------------------------------</a><a name="MEMinitbuf">-</a> 00450 00451 qh_MEMinitbuf 00452 size of initial memory buffer 00453 00454 notes: 00455 use for qh_meminitbuffers() in global.c 00456 */ 00457 #define qh_MEMinitbuf 0x20000 /* initially allocate 128K buffer */ 00458 00459 /*-<a href="qh-user.htm#TOC" 00460 >--------------------------------</a><a name="INFINITE">-</a> 00461 00462 qh_INFINITE 00463 on output, indicates Voronoi center at infinity 00464 */ 00465 #define qh_INFINITE -10.101 00466 00467 /*-<a href="qh-user.htm#TOC" 00468 >--------------------------------</a><a name="DEFAULTbox">-</a> 00469 00470 qh_DEFAULTbox 00471 default box size (Geomview expects 0.5) 00472 00473 qh_DEFAULTbox 00474 default box size for integer coorindate (rbox only) 00475 */ 00476 #define qh_DEFAULTbox 0.5 00477 #define qh_DEFAULTzbox 1e6 00478 00479 /*============================================================*/ 00480 /*============= conditional compilation ======================*/ 00481 /*============================================================*/ 00482 00483 /*-<a href="qh-user.htm#TOC" 00484 >--------------------------------</a><a name="compiler">-</a> 00485 00486 __cplusplus 00487 defined by C++ compilers 00488 00489 __MSC_VER 00490 defined by Microsoft Visual C++ 00491 00492 __MWERKS__ && __POWERPC__ 00493 defined by Metrowerks when compiling for the Power Macintosh 00494 00495 __STDC__ 00496 defined for strict ANSI C 00497 */ 00498 00499 /*-<a href="qh-user.htm#TOC" 00500 >--------------------------------</a><a name="COMPUTEfurthest">-</a> 00501 00502 qh_COMPUTEfurthest 00503 compute furthest distance to an outside point instead of storing it with the facet 00504 =1 to compute furthest 00505 00506 notes: 00507 computing furthest saves memory but costs time 00508 about 40% more distance tests for partitioning 00509 removes facet->furthestdist 00510 */ 00511 #define qh_COMPUTEfurthest 0 00512 00513 /*-<a href="qh-user.htm#TOC" 00514 >--------------------------------</a><a name="KEEPstatistics">-</a> 00515 00516 qh_KEEPstatistics 00517 =0 removes most of statistic gathering and reporting 00518 00519 notes: 00520 if 0, code size is reduced by about 4%. 00521 */ 00522 #define qh_KEEPstatistics 1 00523 00524 /*-<a href="qh-user.htm#TOC" 00525 >--------------------------------</a><a name="MAXoutside">-</a> 00526 00527 qh_MAXoutside 00528 record outer plane for each facet 00529 =1 to record facet->maxoutside 00530 00531 notes: 00532 this takes a realT per facet and slightly slows down qhull 00533 it produces better outer planes for geomview output 00534 */ 00535 #define qh_MAXoutside 1 00536 00537 /*-<a href="qh-user.htm#TOC" 00538 >--------------------------------</a><a name="NOmerge">-</a> 00539 00540 qh_NOmerge 00541 disables facet merging if defined 00542 00543 notes: 00544 This saves about 10% space. 00545 00546 Unless 'Q0' 00547 qh_NOmerge sets 'QJ' to avoid precision errors 00548 00549 #define qh_NOmerge 00550 00551 see: 00552 <a href="mem.h#NOmem">qh_NOmem</a> in mem.c 00553 00554 see user.c/user_eg.c for removing io.o 00555 */ 00556 00557 /*-<a href="qh-user.htm#TOC" 00558 >--------------------------------</a><a name="NOtrace">-</a> 00559 00560 qh_NOtrace 00561 no tracing if defined 00562 00563 notes: 00564 This saves about 5% space. 00565 00566 #define qh_NOtrace 00567 */ 00568 00569 /*-<a href="qh-user.htm#TOC" 00570 >--------------------------------</a><a name="QHpointer">-</a> 00571 00572 qh_QHpointer 00573 access global data with pointer or static structure 00574 00575 qh_QHpointer = 1 access globals via a pointer to allocated memory 00576 enables qh_saveqhull() and qh_restoreqhull() 00577 [2010, gcc] costs about 4% in time and 4% in space 00578 [2003, msvc] costs about 8% in time and 2% in space 00579 00580 = 0 qh_qh and qh_qhstat are static data structures 00581 only one instance of qhull() can be active at a time 00582 default value 00583 00584 qh_QHpointer_dllimport and qh_dllimport define qh_qh as __declspec(dllimport) [libqhull.h] 00585 It is required for msvc-2005. It is not needed for gcc. 00586 00587 notes: 00588 all global variables for qhull are in qh, qhmem, and qhstat 00589 qh is defined in libqhull.h 00590 qhmem is defined in mem.h 00591 qhstat is defined in stat.h 00592 C++ build defines qh_QHpointer [libqhullp.pro, libqhullcpp.pro] 00593 00594 see: 00595 user_eg.c for an example 00596 */ 00597 #ifdef qh_QHpointer 00598 #if qh_dllimport 00599 #error QH6207 Qhull error: Use qh_QHpointer_dllimport instead of qh_dllimport with qh_QHpointer 00600 #endif 00601 #else 00602 #define qh_QHpointer 0 00603 #if qh_QHpointer_dllimport 00604 #error QH6234 Qhull error: Use qh_dllimport instead of qh_QHpointer_dllimport when qh_QHpointer is not defined 00605 #endif 00606 #endif 00607 #if 0 /* sample code */ 00608 qhT *oldqhA, *oldqhB; 00609 00610 exitcode= qh_new_qhull(dim, numpoints, points, ismalloc, 00611 flags, outfile, errfile); 00612 /* use results from first call to qh_new_qhull */ 00613 oldqhA= qh_save_qhull(); 00614 exitcode= qh_new_qhull(dimB, numpointsB, pointsB, ismalloc, 00615 flags, outfile, errfile); 00616 /* use results from second call to qh_new_qhull */ 00617 oldqhB= qh_save_qhull(); 00618 qh_restore_qhull(&oldqhA); 00619 /* use results from first call to qh_new_qhull */ 00620 qh_freeqhull(qh_ALL); /* frees all memory used by first call */ 00621 qh_restore_qhull(&oldqhB); 00622 /* use results from second call to qh_new_qhull */ 00623 qh_freeqhull(!qh_ALL); /* frees long memory used by second call */ 00624 qh_memfreeshort(&curlong, &totlong); /* frees short memory and memory allocator */ 00625 #endif 00626 00627 /*-<a href="qh-user.htm#TOC" 00628 >--------------------------------</a><a name="QUICKhelp">-</a> 00629 00630 qh_QUICKhelp 00631 =1 to use abbreviated help messages, e.g., for degenerate inputs 00632 */ 00633 #define qh_QUICKhelp 0 00634 00635 /*============================================================*/ 00636 /*============= -merge constants- ============================*/ 00637 /*============================================================*/ 00638 /* 00639 These constants effect facet merging. You probably will not need 00640 to modify them. They effect the performance of facet merging. 00641 */ 00642 00643 /*-<a href="qh-user.htm#TOC" 00644 >--------------------------------</a><a name="DIMmergeVertex">-</a> 00645 00646 qh_DIMmergeVertex 00647 max dimension for vertex merging (it is not effective in high-d) 00648 */ 00649 #define qh_DIMmergeVertex 6 00650 00651 /*-<a href="qh-user.htm#TOC" 00652 >--------------------------------</a><a name="DIMreduceBuild">-</a> 00653 00654 qh_DIMreduceBuild 00655 max dimension for vertex reduction during build (slow in high-d) 00656 */ 00657 #define qh_DIMreduceBuild 5 00658 00659 /*-<a href="qh-user.htm#TOC" 00660 >--------------------------------</a><a name="BESTcentrum">-</a> 00661 00662 qh_BESTcentrum 00663 if > 2*dim+n vertices, qh_findbestneighbor() tests centrums (faster) 00664 else, qh_findbestneighbor() tests all vertices (much better merges) 00665 00666 qh_BESTcentrum2 00667 if qh_BESTcentrum2 * DIM3 + BESTcentrum < #vertices tests centrums 00668 */ 00669 #define qh_BESTcentrum 20 00670 #define qh_BESTcentrum2 2 00671 00672 /*-<a href="qh-user.htm#TOC" 00673 >--------------------------------</a><a name="BESTnonconvex">-</a> 00674 00675 qh_BESTnonconvex 00676 if > dim+n neighbors, qh_findbestneighbor() tests nonconvex ridges. 00677 00678 notes: 00679 It is needed because qh_findbestneighbor is slow for large facets 00680 */ 00681 #define qh_BESTnonconvex 15 00682 00683 /*-<a href="qh-user.htm#TOC" 00684 >--------------------------------</a><a name="MAXnewmerges">-</a> 00685 00686 qh_MAXnewmerges 00687 if >n newmerges, qh_merge_nonconvex() calls qh_reducevertices_centrums. 00688 00689 notes: 00690 It is needed because postmerge can merge many facets at once 00691 */ 00692 #define qh_MAXnewmerges 2 00693 00694 /*-<a href="qh-user.htm#TOC" 00695 >--------------------------------</a><a name="MAXnewcentrum">-</a> 00696 00697 qh_MAXnewcentrum 00698 if <= dim+n vertices (n approximates the number of merges), 00699 reset the centrum in qh_updatetested() and qh_mergecycle_facets() 00700 00701 notes: 00702 needed to reduce cost and because centrums may move too much if 00703 many vertices in high-d 00704 */ 00705 #define qh_MAXnewcentrum 5 00706 00707 /*-<a href="qh-user.htm#TOC" 00708 >--------------------------------</a><a name="COPLANARratio">-</a> 00709 00710 qh_COPLANARratio 00711 for 3-d+ merging, qh.MINvisible is n*premerge_centrum 00712 00713 notes: 00714 for non-merging, it's DISTround 00715 */ 00716 #define qh_COPLANARratio 3 00717 00718 /*-<a href="qh-user.htm#TOC" 00719 >--------------------------------</a><a name="DISToutside">-</a> 00720 00721 qh_DISToutside 00722 When is a point clearly outside of a facet? 00723 Stops search in qh_findbestnew or qh_partitionall 00724 qh_findbest uses qh.MINoutside since since it is only called if no merges. 00725 00726 notes: 00727 'Qf' always searches for best facet 00728 if !qh.MERGING, same as qh.MINoutside. 00729 if qh_USEfindbestnew, increase value since neighboring facets may be ill-behaved 00730 [Note: Zdelvertextot occurs normally with interior points] 00731 RBOX 1000 s Z1 G1e-13 t1001188774 | QHULL Tv 00732 When there is a sharp edge, need to move points to a 00733 clearly good facet; otherwise may be lost in another partitioning. 00734 if too big then O(n^2) behavior for partitioning in cone 00735 if very small then important points not processed 00736 Needed in qh_partitionall for 00737 RBOX 1000 s Z1 G1e-13 t1001032651 | QHULL Tv 00738 Needed in qh_findbestnew for many instances of 00739 RBOX 1000 s Z1 G1e-13 t | QHULL Tv 00740 00741 See: 00742 qh_DISToutside -- when is a point clearly outside of a facet 00743 qh_SEARCHdist -- when is facet coplanar with the best facet? 00744 qh_USEfindbestnew -- when to use qh_findbestnew for qh_partitionpoint() 00745 */ 00746 #define qh_DISToutside ((qh_USEfindbestnew ? 2 : 1) * \ 00747 fmax_((qh MERGING ? 2 : 1)*qh MINoutside, qh max_outside)) 00748 00749 /*-<a href="qh-user.htm#TOC" 00750 >--------------------------------</a><a name="RATIOnearinside">-</a> 00751 00752 qh_RATIOnearinside 00753 ratio of qh.NEARinside to qh.ONEmerge for retaining inside points for 00754 qh_check_maxout(). 00755 00756 notes: 00757 This is overkill since do not know the correct value. 00758 It effects whether 'Qc' reports all coplanar points 00759 Not used for 'd' since non-extreme points are coplanar 00760 */ 00761 #define qh_RATIOnearinside 5 00762 00763 /*-<a href="qh-user.htm#TOC" 00764 >--------------------------------</a><a name="SEARCHdist">-</a> 00765 00766 qh_SEARCHdist 00767 When is a facet coplanar with the best facet? 00768 qh_findbesthorizon: all coplanar facets of the best facet need to be searched. 00769 00770 See: 00771 qh_DISToutside -- when is a point clearly outside of a facet 00772 qh_SEARCHdist -- when is facet coplanar with the best facet? 00773 qh_USEfindbestnew -- when to use qh_findbestnew for qh_partitionpoint() 00774 */ 00775 #define qh_SEARCHdist ((qh_USEfindbestnew ? 2 : 1) * \ 00776 (qh max_outside + 2 * qh DISTround + fmax_( qh MINvisible, qh MAXcoplanar))); 00777 00778 /*-<a href="qh-user.htm#TOC" 00779 >--------------------------------</a><a name="USEfindbestnew">-</a> 00780 00781 qh_USEfindbestnew 00782 Always use qh_findbestnew for qh_partitionpoint, otherwise use 00783 qh_findbestnew if merged new facet or sharpnewfacets. 00784 00785 See: 00786 qh_DISToutside -- when is a point clearly outside of a facet 00787 qh_SEARCHdist -- when is facet coplanar with the best facet? 00788 qh_USEfindbestnew -- when to use qh_findbestnew for qh_partitionpoint() 00789 */ 00790 #define qh_USEfindbestnew (zzval_(Ztotmerge) > 50) 00791 00792 /*-<a href="qh-user.htm#TOC" 00793 >--------------------------------</a><a name="WIDEcoplanar">-</a> 00794 00795 qh_WIDEcoplanar 00796 n*MAXcoplanar or n*MINvisible for a WIDEfacet 00797 00798 if vertex is further than qh.WIDEfacet from the hyperplane 00799 then its ridges are not counted in computing the area, and 00800 the facet's centrum is frozen. 00801 00802 notes: 00803 qh.WIDEfacet= max(qh.MAXoutside,qh_WIDEcoplanar*qh.MAXcoplanar, 00804 qh_WIDEcoplanar * qh.MINvisible); 00805 */ 00806 #define qh_WIDEcoplanar 6 00807 00808 /*-<a href="qh-user.htm#TOC" 00809 >--------------------------------</a><a name="MAXnarrow">-</a> 00810 00811 qh_MAXnarrow 00812 max. cosine in initial hull that sets qh.NARROWhull 00813 00814 notes: 00815 If qh.NARROWhull, the initial partition does not make 00816 coplanar points. If narrow, a coplanar point can be 00817 coplanar to two facets of opposite orientations and 00818 distant from the exact convex hull. 00819 00820 Conservative estimate. Don't actually see problems until it is -1.0 00821 */ 00822 #define qh_MAXnarrow -0.99999999 00823 00824 /*-<a href="qh-user.htm#TOC" 00825 >--------------------------------</a><a name="WARNnarrow">-</a> 00826 00827 qh_WARNnarrow 00828 max. cosine in initial hull to warn about qh.NARROWhull 00829 00830 notes: 00831 this is a conservative estimate. 00832 Don't actually see problems until it is -1.0. See qh-impre.htm 00833 */ 00834 #define qh_WARNnarrow -0.999999999999999 00835 00836 /*-<a href="qh-user.htm#TOC" 00837 >--------------------------------</a><a name="ZEROdelaunay">-</a> 00838 00839 qh_ZEROdelaunay 00840 a zero Delaunay facet occurs for input sites coplanar with their convex hull 00841 the last normal coefficient of a zero Delaunay facet is within 00842 qh_ZEROdelaunay * qh.ANGLEround of 0 00843 00844 notes: 00845 qh_ZEROdelaunay does not allow for joggled input ('QJ'). 00846 00847 You can avoid zero Delaunay facets by surrounding the input with a box. 00848 00849 Use option 'PDk:-n' to explicitly define zero Delaunay facets 00850 k= dimension of input sites (e.g., 3 for 3-d Delaunay triangulation) 00851 n= the cutoff for zero Delaunay facets (e.g., 'PD3:-1e-12') 00852 */ 00853 #define qh_ZEROdelaunay 2 00854 00855 #endif /* qh_DEFuser */ 00856 00857 00858