|
||||
File indexing completed on 2025-01-18 10:01:54
0001 /*<html><pre> -<a href="qh-user.htm" 0002 >-------------------------------</a><a name="TOP">-</a> 0003 0004 user.h 0005 user redefinable constants 0006 0007 for each source file, user.h is included first 0008 0009 see qh-user.htm. see COPYING for copyright information. 0010 0011 See user.c for sample code. 0012 0013 before reading any code, review libqhull.h for data structure definitions 0014 and the "qh" macro. 0015 0016 Sections: 0017 ============= qhull library constants ====================== 0018 ============= data types and configuration macros ========== 0019 ============= performance related constants ================ 0020 ============= memory constants ============================= 0021 ============= joggle constants ============================= 0022 ============= conditional compilation ====================== 0023 ============= merge constants ============================== 0024 ============= Microsoft DevStudio ========================== 0025 0026 Code flags -- 0027 NOerrors -- the code does not call qh_errexit() 0028 WARN64 -- the code may be incompatible with 64-bit pointers 0029 0030 */ 0031 0032 #include <float.h> 0033 #include <limits.h> 0034 #include <time.h> 0035 0036 #ifndef qhDEFuser 0037 #define qhDEFuser 1 0038 0039 /* Derived from Qt's corelib/global/qglobal.h */ 0040 #if !defined(SAG_COM) && !defined(__CYGWIN__) && (defined(WIN64) || defined(_WIN64) || defined(__WIN64__) || defined(WIN32) || defined(_WIN32) || defined(__WIN32__) || defined(__NT__)) 0041 # define QHULL_OS_WIN 0042 #elif defined(__MWERKS__) && defined(__INTEL__) /* Metrowerks discontinued before the release of Intel Macs */ 0043 # define QHULL_OS_WIN 0044 #endif 0045 0046 /*============================================================*/ 0047 /*============= qhull library constants ======================*/ 0048 /*============================================================*/ 0049 0050 /*-<a href="qh-user.htm#TOC" 0051 >--------------------------------</a><a name="filenamelen">-</a> 0052 0053 FILENAMElen -- max length for TI and TO filenames 0054 0055 */ 0056 0057 #define qh_FILENAMElen 500 0058 0059 /*-<a href="qh-user.htm#TOC" 0060 >--------------------------------</a><a name="msgcode">-</a> 0061 0062 msgcode -- Unique message codes for qh_fprintf 0063 0064 If add new messages, assign these values and increment in user.h and user.h 0065 See QhullError.h for 10000 error codes. 0066 Cannot use '0031' since it would be octal 0067 0068 def counters = [31/32/33/38, 1067, 2113, 3079, 4097, 5006, 0069 6429, 7027/7028/7035/7068/7070/7102, 8163, 9428, 10000, 11034] 0070 0071 See: qh_ERR* [libqhull.h] 0072 */ 0073 0074 #define MSG_TRACE0 0 /* always include if logging ('Tn') */ 0075 #define MSG_TRACE1 1000 0076 #define MSG_TRACE2 2000 0077 #define MSG_TRACE3 3000 0078 #define MSG_TRACE4 4000 0079 #define MSG_TRACE5 5000 0080 #define MSG_ERROR 6000 /* errors written to qh.ferr */ 0081 #define MSG_WARNING 7000 0082 #define MSG_STDERR 8000 /* log messages Written to qh.ferr */ 0083 #define MSG_OUTPUT 9000 0084 #define MSG_QHULL_ERROR 10000 /* errors thrown by QhullError.cpp (QHULLlastError is in QhullError.h) */ 0085 #define MSG_FIX 11000 /* Document as 'QH11... FIX: ...' */ 0086 #define MSG_MAXLEN 3000 /* qh_printhelp_degenerate() in user.c */ 0087 0088 0089 /*-<a href="qh-user.htm#TOC" 0090 >--------------------------------</a><a name="qh_OPTIONline">-</a> 0091 0092 qh_OPTIONline -- max length of an option line 'FO' 0093 */ 0094 #define qh_OPTIONline 80 0095 0096 /*============================================================*/ 0097 /*============= data types and configuration macros ==========*/ 0098 /*============================================================*/ 0099 0100 /*-<a href="qh-user.htm#TOC" 0101 >--------------------------------</a><a name="realT">-</a> 0102 0103 realT 0104 set the size of floating point numbers 0105 0106 qh_REALdigits 0107 maximimum number of significant digits 0108 0109 qh_REAL_1, qh_REAL_2n, qh_REAL_3n 0110 format strings for printf 0111 0112 qh_REALmax, qh_REALmin 0113 maximum and minimum (near zero) values 0114 0115 qh_REALepsilon 0116 machine roundoff. Maximum roundoff error for addition and multiplication. 0117 0118 notes: 0119 Select whether to store floating point numbers in single precision (float) 0120 or double precision (double). 0121 0122 Use 'float' to save about 8% in time and 25% in space. This is particularly 0123 helpful if high-d where convex hulls are space limited. Using 'float' also 0124 reduces the printed size of Qhull's output since numbers have 8 digits of 0125 precision. 0126 0127 Use 'double' when greater arithmetic precision is needed. This is needed 0128 for Delaunay triangulations and Voronoi diagrams when you are not merging 0129 facets. 0130 0131 If 'double' gives insufficient precision, your data probably includes 0132 degeneracies. If so you should use facet merging (done by default) 0133 or exact arithmetic (see imprecision section of manual, qh-impre.htm). 0134 You may also use option 'Po' to force output despite precision errors. 0135 0136 You may use 'long double', but many format statements need to be changed 0137 and you may need a 'long double' square root routine. S. Grundmann 0138 (sg@eeiwzb.et.tu-dresden.de) has done this. He reports that the code runs 0139 much slower with little gain in precision. 0140 0141 WARNING: on some machines, int f(){realT a= REALmax;return (a == REALmax);} 0142 returns False. Use (a > REALmax/2) instead of (a == REALmax). 0143 0144 REALfloat = 1 all numbers are 'float' type 0145 = 0 all numbers are 'double' type 0146 */ 0147 #define REALfloat 0 0148 0149 #if (REALfloat == 1) 0150 #define realT float 0151 #define REALmax FLT_MAX 0152 #define REALmin FLT_MIN 0153 #define REALepsilon FLT_EPSILON 0154 #define qh_REALdigits 8 /* maximum number of significant digits */ 0155 #define qh_REAL_1 "%6.8g " 0156 #define qh_REAL_2n "%6.8g %6.8g\n" 0157 #define qh_REAL_3n "%6.8g %6.8g %6.8g\n" 0158 0159 #elif (REALfloat == 0) 0160 #define realT double 0161 #define REALmax DBL_MAX 0162 #define REALmin DBL_MIN 0163 #define REALepsilon DBL_EPSILON 0164 #define qh_REALdigits 16 /* maximum number of significant digits */ 0165 #define qh_REAL_1 "%6.16g " 0166 #define qh_REAL_2n "%6.16g %6.16g\n" 0167 #define qh_REAL_3n "%6.16g %6.16g %6.16g\n" 0168 0169 #else 0170 #error unknown float option 0171 #endif 0172 0173 /*-<a href="qh-user.htm#TOC" 0174 >--------------------------------</a><a name="countT">-</a> 0175 0176 countT 0177 The type for counts and identifiers (e.g., the number of points, vertex identifiers) 0178 Currently used by C++ code-only. Decided against using it for setT because most sets are small. 0179 0180 Defined as 'int' for C-code compatibility and QH11026 0181 0182 QH11026 FIX: countT may be defined as a 'unsigned int', but several code issues need to be solved first. See countT in Changes.txt 0183 */ 0184 0185 #ifndef DEFcountT 0186 #define DEFcountT 1 0187 typedef int countT; 0188 #endif 0189 #define COUNTmax INT_MAX 0190 0191 /*-<a href="qh-user.htm#TOC" 0192 >--------------------------------</a><a name="qh_POINTSmax">-</a> 0193 0194 qh_POINTSmax 0195 Maximum number of points for qh.num_points and point allocation in qh_readpoints 0196 */ 0197 #define qh_POINTSmax (INT_MAX-16) 0198 0199 /*-<a href="qh-user.htm#TOC" 0200 >--------------------------------</a><a name="CPUclock">-</a> 0201 0202 qh_CPUclock 0203 define the clock() function for reporting the total time spent by Qhull 0204 returns CPU ticks as a 'long int' 0205 qh_CPUclock is only used for reporting the total time spent by Qhull 0206 0207 qh_SECticks 0208 the number of clock ticks per second 0209 0210 notes: 0211 looks for CLOCKS_PER_SEC, CLOCKS_PER_SECOND, or assumes microseconds 0212 to define a custom clock, set qh_CLOCKtype to 0 0213 0214 if your system does not use clock() to return CPU ticks, replace 0215 qh_CPUclock with the corresponding function. It is converted 0216 to 'unsigned long' to prevent wrap-around during long runs. By default, 0217 <time.h> defines clock_t as 'long' 0218 0219 Set qh_CLOCKtype to 0220 0221 1 for CLOCKS_PER_SEC, CLOCKS_PER_SECOND, or microsecond 0222 Note: may fail if more than 1 hour elapsed time 0223 0224 2 use qh_clock() with POSIX times() (see global.c) 0225 */ 0226 #define qh_CLOCKtype 1 /* change to the desired number */ 0227 0228 #if (qh_CLOCKtype == 1) 0229 0230 #if defined(CLOCKS_PER_SECOND) 0231 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock, may be converted to approximate double */ 0232 #define qh_SECticks CLOCKS_PER_SECOND 0233 0234 #elif defined(CLOCKS_PER_SEC) 0235 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock, may be converted to approximate double */ 0236 #define qh_SECticks CLOCKS_PER_SEC 0237 0238 #elif defined(CLK_TCK) 0239 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock, may be converted to approximate double */ 0240 #define qh_SECticks CLK_TCK 0241 0242 #else 0243 #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock, may be converted to approximate double */ 0244 #define qh_SECticks 1E6 0245 #endif 0246 0247 #elif (qh_CLOCKtype == 2) 0248 #define qh_CPUclock qh_clock() /* return CPU clock, may be converted to approximate double */ 0249 #define qh_SECticks 100 0250 0251 #else /* qh_CLOCKtype == ? */ 0252 #error unknown clock option 0253 #endif 0254 0255 /*-<a href="qh-user.htm#TOC" 0256 >--------------------------------</a><a name="RANDOM">-</a> 0257 0258 qh_RANDOMtype, qh_RANDOMmax, qh_RANDOMseed 0259 define random number generator 0260 0261 qh_RANDOMint generates a random integer between 0 and qh_RANDOMmax. 0262 qh_RANDOMseed sets the random number seed for qh_RANDOMint 0263 0264 Set qh_RANDOMtype (default 5) to: 0265 1 for random() with 31 bits (UCB) 0266 2 for rand() with RAND_MAX or 15 bits (system 5) 0267 3 for rand() with 31 bits (Sun) 0268 4 for lrand48() with 31 bits (Solaris) 0269 5 for qh_rand() with 31 bits (included with Qhull) 0270 0271 notes: 0272 Random numbers are used by rbox to generate point sets. Random 0273 numbers are used by Qhull to rotate the input ('QRn' option), 0274 simulate a randomized algorithm ('Qr' option), and to simulate 0275 roundoff errors ('Rn' option). 0276 0277 Random number generators differ between systems. Most systems provide 0278 rand() but the period varies. The period of rand() is not critical 0279 since qhull does not normally use random numbers. 0280 0281 The default generator is Park & Miller's minimal standard random 0282 number generator [CACM 31:1195 '88]. It is included with Qhull. 0283 0284 If qh_RANDOMmax is wrong, qhull will report a warning and Geomview 0285 output will likely be invisible. 0286 */ 0287 #define qh_RANDOMtype 5 /* *** change to the desired number *** */ 0288 0289 #if (qh_RANDOMtype == 1) 0290 #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, random()/MAX */ 0291 #define qh_RANDOMint random() 0292 #define qh_RANDOMseed_(seed) srandom(seed); 0293 0294 #elif (qh_RANDOMtype == 2) 0295 #ifdef RAND_MAX 0296 #define qh_RANDOMmax ((realT)RAND_MAX) 0297 #else 0298 #define qh_RANDOMmax ((realT)32767) /* 15 bits (System 5) */ 0299 #endif 0300 #define qh_RANDOMint rand() 0301 #define qh_RANDOMseed_(seed) srand((unsigned int)seed); 0302 0303 #elif (qh_RANDOMtype == 3) 0304 #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, Sun */ 0305 #define qh_RANDOMint rand() 0306 #define qh_RANDOMseed_(seed) srand((unsigned int)seed); 0307 0308 #elif (qh_RANDOMtype == 4) 0309 #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, lrand38()/MAX */ 0310 #define qh_RANDOMint lrand48() 0311 #define qh_RANDOMseed_(seed) srand48(seed); 0312 0313 #elif (qh_RANDOMtype == 5) 0314 #define qh_RANDOMmax ((realT)2147483646UL) /* 31 bits, qh_rand/MAX */ 0315 #define qh_RANDOMint qh_rand( ) 0316 #define qh_RANDOMseed_(seed) qh_srand(seed); 0317 /* unlike rand(), never returns 0 */ 0318 0319 #else 0320 #error: unknown random option 0321 #endif 0322 0323 /*-<a href="qh-user.htm#TOC" 0324 >--------------------------------</a><a name="ORIENTclock">-</a> 0325 0326 qh_ORIENTclock 0327 0 for inward pointing normals by Geomview convention 0328 */ 0329 #define qh_ORIENTclock 0 0330 0331 /*-<a href="qh-user.htm#TOC" 0332 >--------------------------------</a><a name="RANDOMdist">-</a> 0333 0334 qh_RANDOMdist 0335 define for random perturbation of qh_distplane and qh_setfacetplane (qh.RANDOMdist, 'QRn') 0336 0337 For testing qh.DISTround. Qhull should not depend on computations always producing the same roundoff error 0338 0339 #define qh_RANDOMdist 1e-13 0340 */ 0341 0342 /*============================================================*/ 0343 /*============= joggle constants =============================*/ 0344 /*============================================================*/ 0345 0346 /*-<a href="qh-user.htm#TOC" 0347 >--------------------------------</a><a name="JOGGLEdefault">-</a> 0348 0349 qh_JOGGLEdefault 0350 default qh.JOGGLEmax is qh.DISTround * qh_JOGGLEdefault 0351 0352 notes: 0353 rbox s r 100 | qhull QJ1e-15 QR0 generates 90% faults at distround 7e-16 0354 rbox s r 100 | qhull QJ1e-14 QR0 generates 70% faults 0355 rbox s r 100 | qhull QJ1e-13 QR0 generates 35% faults 0356 rbox s r 100 | qhull QJ1e-12 QR0 generates 8% faults 0357 rbox s r 100 | qhull QJ1e-11 QR0 generates 1% faults 0358 rbox s r 100 | qhull QJ1e-10 QR0 generates 0% faults 0359 rbox 1000 W0 | qhull QJ1e-12 QR0 generates 86% faults 0360 rbox 1000 W0 | qhull QJ1e-11 QR0 generates 20% faults 0361 rbox 1000 W0 | qhull QJ1e-10 QR0 generates 2% faults 0362 the later have about 20 points per facet, each of which may interfere 0363 0364 pick a value large enough to avoid retries on most inputs 0365 */ 0366 #define qh_JOGGLEdefault 30000.0 0367 0368 /*-<a href="qh-user.htm#TOC" 0369 >--------------------------------</a><a name="JOGGLEincrease">-</a> 0370 0371 qh_JOGGLEincrease 0372 factor to increase qh.JOGGLEmax on qh_JOGGLEretry or qh_JOGGLEagain 0373 */ 0374 #define qh_JOGGLEincrease 10.0 0375 0376 /*-<a href="qh-user.htm#TOC" 0377 >--------------------------------</a><a name="JOGGLEretry">-</a> 0378 0379 qh_JOGGLEretry 0380 if ZZretry = qh_JOGGLEretry, increase qh.JOGGLEmax 0381 0382 notes: 0383 try twice at the original value in case of bad luck the first time 0384 */ 0385 #define qh_JOGGLEretry 2 0386 0387 /*-<a href="qh-user.htm#TOC" 0388 >--------------------------------</a><a name="JOGGLEagain">-</a> 0389 0390 qh_JOGGLEagain 0391 every following qh_JOGGLEagain, increase qh.JOGGLEmax 0392 0393 notes: 0394 1 is OK since it's already failed qh_JOGGLEretry times 0395 */ 0396 #define qh_JOGGLEagain 1 0397 0398 /*-<a href="qh-user.htm#TOC" 0399 >--------------------------------</a><a name="JOGGLEmaxincrease">-</a> 0400 0401 qh_JOGGLEmaxincrease 0402 maximum qh.JOGGLEmax due to qh_JOGGLEincrease 0403 relative to qh.MAXwidth 0404 0405 notes: 0406 qh.joggleinput will retry at this value until qh_JOGGLEmaxretry 0407 */ 0408 #define qh_JOGGLEmaxincrease 1e-2 0409 0410 /*-<a href="qh-user.htm#TOC" 0411 >--------------------------------</a><a name="JOGGLEmaxretry">-</a> 0412 0413 qh_JOGGLEmaxretry 0414 stop after qh_JOGGLEmaxretry attempts 0415 */ 0416 #define qh_JOGGLEmaxretry 50 0417 0418 /*============================================================*/ 0419 /*============= performance related constants ================*/ 0420 /*============================================================*/ 0421 0422 /*-<a href="qh-user.htm#TOC" 0423 >--------------------------------</a><a name="HASHfactor">-</a> 0424 0425 qh_HASHfactor 0426 total hash slots / used hash slots. Must be at least 1.1. 0427 0428 notes: 0429 =2 for at worst 50% occupancy for qh.hash_table and normally 25% occupancy 0430 */ 0431 #define qh_HASHfactor 2 0432 0433 /*-<a href="qh-user.htm#TOC" 0434 >--------------------------------</a><a name="VERIFYdirect">-</a> 0435 0436 qh_VERIFYdirect 0437 with 'Tv' verify all points against all facets if op count is smaller 0438 0439 notes: 0440 if greater, calls qh_check_bestdist() instead 0441 */ 0442 #define qh_VERIFYdirect 1000000 0443 0444 /*-<a href="qh-user.htm#TOC" 0445 >--------------------------------</a><a name="INITIALsearch">-</a> 0446 0447 qh_INITIALsearch 0448 if qh_INITIALmax, search points up to this dimension 0449 */ 0450 #define qh_INITIALsearch 6 0451 0452 /*-<a href="qh-user.htm#TOC" 0453 >--------------------------------</a><a name="INITIALmax">-</a> 0454 0455 qh_INITIALmax 0456 if dim >= qh_INITIALmax, use min/max coordinate points for initial simplex 0457 0458 notes: 0459 from points with non-zero determinants 0460 use option 'Qs' to override (much slower) 0461 */ 0462 #define qh_INITIALmax 8 0463 0464 /*============================================================*/ 0465 /*============= memory constants =============================*/ 0466 /*============================================================*/ 0467 0468 /*-<a href="qh-user.htm#TOC" 0469 >--------------------------------</a><a name="MEMalign">-</a> 0470 0471 qh_MEMalign 0472 memory alignment for qh_meminitbuffers() in global.c 0473 0474 notes: 0475 to avoid bus errors, memory allocation must consider alignment requirements. 0476 malloc() automatically takes care of alignment. Since mem.c manages 0477 its own memory, we need to explicitly specify alignment in 0478 qh_meminitbuffers(). 0479 0480 A safe choice is sizeof(double). sizeof(float) may be used if doubles 0481 do not occur in data structures and pointers are the same size. Be careful 0482 of machines (e.g., DEC Alpha) with large pointers. 0483 0484 If using gcc, best alignment is [fmax_() is defined in geom.h] 0485 #define qh_MEMalign fmax_(__alignof__(realT),__alignof__(void *)) 0486 */ 0487 #define qh_MEMalign ((int)(fmax_(sizeof(realT), sizeof(void *)))) 0488 0489 /*-<a href="qh-user.htm#TOC" 0490 >--------------------------------</a><a name="MEMbufsize">-</a> 0491 0492 qh_MEMbufsize 0493 size of additional memory buffers 0494 0495 notes: 0496 used for qh_meminitbuffers() in global.c 0497 */ 0498 #define qh_MEMbufsize 0x10000 /* allocate 64K memory buffers */ 0499 0500 /*-<a href="qh-user.htm#TOC" 0501 >--------------------------------</a><a name="MEMinitbuf">-</a> 0502 0503 qh_MEMinitbuf 0504 size of initial memory buffer 0505 0506 notes: 0507 use for qh_meminitbuffers() in global.c 0508 */ 0509 #define qh_MEMinitbuf 0x20000 /* initially allocate 128K buffer */ 0510 0511 /*-<a href="qh-user.htm#TOC" 0512 >--------------------------------</a><a name="INFINITE">-</a> 0513 0514 qh_INFINITE 0515 on output, indicates Voronoi center at infinity 0516 */ 0517 #define qh_INFINITE -10.101 0518 0519 /*-<a href="qh-user.htm#TOC" 0520 >--------------------------------</a><a name="DEFAULTbox">-</a> 0521 0522 qh_DEFAULTbox 0523 default box size (Geomview expects 0.5) 0524 0525 qh_DEFAULTbox 0526 default box size for integer coorindate (rbox only) 0527 */ 0528 #define qh_DEFAULTbox 0.5 0529 #define qh_DEFAULTzbox 1e6 0530 0531 /*============================================================*/ 0532 /*============= conditional compilation ======================*/ 0533 /*============================================================*/ 0534 0535 /*-<a href="qh-user.htm#TOC" 0536 >--------------------------------</a><a name="compiler">-</a> 0537 0538 __cplusplus 0539 defined by C++ compilers 0540 0541 __MSC_VER 0542 defined by Microsoft Visual C++ 0543 0544 __MWERKS__ && __INTEL__ 0545 defined by Metrowerks when compiling for Windows (not Intel-based Macintosh) 0546 0547 __MWERKS__ && __POWERPC__ 0548 defined by Metrowerks when compiling for PowerPC-based Macintosh 0549 0550 __STDC__ 0551 defined for strict ANSI C 0552 */ 0553 0554 /*-<a href="qh-user.htm#TOC" 0555 >--------------------------------</a><a name="COMPUTEfurthest">-</a> 0556 0557 qh_COMPUTEfurthest 0558 compute furthest distance to an outside point instead of storing it with the facet 0559 =1 to compute furthest 0560 0561 notes: 0562 computing furthest saves memory but costs time 0563 about 40% more distance tests for partitioning 0564 removes facet->furthestdist 0565 */ 0566 #define qh_COMPUTEfurthest 0 0567 0568 /*-<a href="qh-user.htm#TOC" 0569 >--------------------------------</a><a name="KEEPstatistics">-</a> 0570 0571 qh_KEEPstatistics 0572 =0 removes most of statistic gathering and reporting 0573 0574 notes: 0575 if 0, code size is reduced by about 4%. 0576 */ 0577 #define qh_KEEPstatistics 1 0578 0579 /*-<a href="qh-user.htm#TOC" 0580 >--------------------------------</a><a name="MAXoutside">-</a> 0581 0582 qh_MAXoutside 0583 record outer plane for each facet 0584 =1 to record facet->maxoutside 0585 0586 notes: 0587 this takes a realT per facet and slightly slows down qhull 0588 it produces better outer planes for geomview output 0589 */ 0590 #define qh_MAXoutside 1 0591 0592 /*-<a href="qh-user.htm#TOC" 0593 >--------------------------------</a><a name="NOmerge">-</a> 0594 0595 qh_NOmerge 0596 disables facet merging if defined 0597 For MSVC compiles, use qhull-exports-nomerge.def instead of qhull-exports.def 0598 0599 notes: 0600 This saves about 25% space, 30% space in combination with qh_NOtrace, 0601 and 36% with qh_NOtrace and qh_KEEPstatistics 0 0602 0603 Unless option 'Q0' is used 0604 qh_NOmerge sets 'QJ' to avoid precision errors 0605 0606 see: 0607 <a href="mem.h#NOmem">qh_NOmem</a> in mem.h 0608 0609 see user.c/user_eg.c for removing io.o 0610 0611 #define qh_NOmerge 0612 */ 0613 0614 /*-<a href="qh-user.htm#TOC" 0615 >--------------------------------</a><a name="NOtrace">-</a> 0616 0617 qh_NOtrace 0618 no tracing if defined 0619 disables 'Tn', 'TMn', 'TPn' and 'TWn' 0620 override with 'Qw' for qh_addpoint tracing and various other items 0621 0622 notes: 0623 This saves about 15% space. 0624 Removes all traceN((...)) code and substantial sections of qh.IStracing code 0625 0626 #define qh_NOtrace 0627 */ 0628 0629 /*-<a href="qh-user.htm#TOC" 0630 >--------------------------------</a><a name="QHpointer">-</a> 0631 0632 qh_QHpointer 0633 access global data with pointer or static structure 0634 0635 qh_QHpointer = 1 access globals via a pointer to allocated memory 0636 enables qh_saveqhull() and qh_restoreqhull() 0637 [2010, gcc] costs about 4% in time and 4% in space 0638 [2003, msvc] costs about 8% in time and 2% in space 0639 0640 = 0 qh_qh and qh_qhstat are static data structures 0641 only one instance of qhull() can be active at a time 0642 default value 0643 0644 For msvc, define qh_QHpointer_dllimport or qh_dllimport for qh_qh as __declspec(dllimport) [libqhull.h] 0645 0646 notes: 0647 [jan'16] qh_QHpointer is deprecated for Qhull. Use libqhull_r instead. 0648 all global variables for qhull are in qh, qhmem, and qhstat 0649 qh is defined in libqhull.h 0650 qhmem is defined in mem.h 0651 qhstat is defined in stat.h 0652 0653 #define qh_QHpointer 1 0654 */ 0655 0656 #ifdef qh_QHpointer 0657 #ifdef qh_dllimport 0658 #error QH6207 Qhull error: Define qh_QHpointer_dllimport instead of qh_dllimport with qh_QHpointer 0659 #endif 0660 #else 0661 #define qh_QHpointer 0 0662 #ifdef qh_QHpointer_dllimport 0663 #error QH6234 Qhull error: Define qh_dllimport instead of qh_QHpointer_dllimport when qh_QHpointer is not defined 0664 #endif 0665 #endif 0666 #if 0 /* sample code */ 0667 qhT *oldqhA, *oldqhB; 0668 0669 exitcode= qh_new_qhull(dim, numpoints, points, ismalloc, 0670 flags, outfile, errfile); 0671 /* use results from first call to qh_new_qhull */ 0672 oldqhA= qh_save_qhull(); 0673 exitcode= qh_new_qhull(dimB, numpointsB, pointsB, ismalloc, 0674 flags, outfile, errfile); 0675 /* use results from second call to qh_new_qhull */ 0676 oldqhB= qh_save_qhull(); 0677 qh_restore_qhull(&oldqhA); 0678 /* use results from first call to qh_new_qhull */ 0679 qh_freeqhull(qh_ALL); /* frees all memory used by first call */ 0680 qh_restore_qhull(&oldqhB); 0681 /* use results from second call to qh_new_qhull */ 0682 qh_freeqhull(!qh_ALL); /* frees long memory used by second call */ 0683 qh_memfreeshort(&curlong, &totlong); /* frees short memory and memory allocator */ 0684 #endif 0685 0686 /*-<a href="qh-user.htm#TOC" 0687 >--------------------------------</a><a name="QUICKhelp">-</a> 0688 0689 qh_QUICKhelp 0690 =1 to use abbreviated help messages, e.g., for degenerate inputs 0691 */ 0692 #define qh_QUICKhelp 0 0693 0694 /*============================================================*/ 0695 /*============= merge constants ==============================*/ 0696 /*============================================================*/ 0697 /* 0698 These constants effect facet merging. You probably will not need 0699 to modify them. They effect the performance of facet merging. 0700 */ 0701 0702 /*-<a href="qh-user.htm#TOC" 0703 >--------------------------------</a><a name="BESTcentrum">-</a> 0704 0705 qh_BESTcentrum 0706 if > 2*dim+n vertices, qh_findbestneighbor() tests centrums (faster) 0707 else, qh_findbestneighbor() tests all vertices (much better merges) 0708 0709 qh_BESTcentrum2 0710 if qh_BESTcentrum2 * DIM3 + BESTcentrum < #vertices tests centrums 0711 */ 0712 #define qh_BESTcentrum 20 0713 #define qh_BESTcentrum2 2 0714 0715 /*-<a href="qh-user.htm#TOC" 0716 >--------------------------------</a><a name="BESTnonconvex">-</a> 0717 0718 qh_BESTnonconvex 0719 if > dim+n neighbors, qh_findbestneighbor() tests nonconvex ridges. 0720 0721 notes: 0722 It is needed because qh_findbestneighbor is slow for large facets 0723 */ 0724 #define qh_BESTnonconvex 15 0725 0726 /*-<a href="qh-user.htm#TOC" 0727 >--------------------------------</a><a name="COPLANARratio">-</a> 0728 0729 qh_COPLANARratio 0730 for 3-d+ merging, qh.MINvisible is n*premerge_centrum 0731 0732 notes: 0733 for non-merging, it's DISTround 0734 */ 0735 #define qh_COPLANARratio 3 0736 0737 /*-<a href="qh-user.htm#TOC" 0738 >--------------------------------</a><a name="DIMmergeVertex">-</a> 0739 0740 qh_DIMmergeVertex 0741 max dimension for vertex merging (it is not effective in high-d) 0742 */ 0743 #define qh_DIMmergeVertex 6 0744 0745 /*-<a href="qh-user.htm#TOC" 0746 >--------------------------------</a><a name="DIMreduceBuild">-</a> 0747 0748 qh_DIMreduceBuild 0749 max dimension for vertex reduction during build (slow in high-d) 0750 */ 0751 #define qh_DIMreduceBuild 5 0752 0753 /*-<a href="qh-user.htm#TOC" 0754 >--------------------------------</a><a name="DISToutside">-</a> 0755 0756 qh_DISToutside 0757 When is a point clearly outside of a facet? 0758 Stops search in qh_findbestnew or qh_partitionall 0759 qh_findbest uses qh.MINoutside since since it is only called if no merges. 0760 0761 notes: 0762 'Qf' always searches for best facet 0763 if !qh.MERGING, same as qh.MINoutside. 0764 if qh_USEfindbestnew, increase value since neighboring facets may be ill-behaved 0765 [Note: Zdelvertextot occurs normally with interior points] 0766 RBOX 1000 s Z1 G1e-13 t1001188774 | QHULL Tv 0767 When there is a sharp edge, need to move points to a 0768 clearly good facet; otherwise may be lost in another partitioning. 0769 if too big then O(n^2) behavior for partitioning in cone 0770 if very small then important points not processed 0771 Needed in qh_partitionall for 0772 RBOX 1000 s Z1 G1e-13 t1001032651 | QHULL Tv 0773 Needed in qh_findbestnew for many instances of 0774 RBOX 1000 s Z1 G1e-13 t | QHULL Tv 0775 0776 See: 0777 qh_DISToutside -- when is a point clearly outside of a facet 0778 qh_SEARCHdist -- when is facet coplanar with the best facet? 0779 qh_USEfindbestnew -- when to use qh_findbestnew for qh_partitionpoint() 0780 */ 0781 #define qh_DISToutside ((qh_USEfindbestnew ? 2 : 1) * \ 0782 fmax_((qh MERGING ? 2 : 1)*qh MINoutside, qh max_outside)) 0783 0784 /*-<a href="qh-user.htm#TOC" 0785 >--------------------------------</a><a name="MAXcheckpoint">-</a> 0786 0787 qh_MAXcheckpoint 0788 Report up to qh_MAXcheckpoint errors per facet in qh_check_point ('Tv') 0789 */ 0790 #define qh_MAXcheckpoint 10 0791 0792 /*-<a href="qh-user.htm#TOC" 0793 >--------------------------------</a><a name="MAXcoplanarcentrum">-</a> 0794 0795 qh_MAXcoplanarcentrum 0796 if pre-merging with qh.MERGEexact ('Qx') and f.nummerge > qh_MAXcoplanarcentrum 0797 use f.maxoutside instead of qh.centrum_radius for coplanarity testing 0798 0799 notes: 0800 see qh_test_nonsimplicial_merges 0801 with qh.MERGEexact, a coplanar ridge is ignored until post-merging 0802 otherwise a large facet with many merges may take all the facets 0803 */ 0804 #define qh_MAXcoplanarcentrum 10 0805 0806 /*-<a href="qh-user.htm#TOC" 0807 >--------------------------------</a><a name="MAXnewcentrum">-</a> 0808 0809 qh_MAXnewcentrum 0810 if <= dim+n vertices (n approximates the number of merges), 0811 reset the centrum in qh_updatetested() and qh_mergecycle_facets() 0812 0813 notes: 0814 needed to reduce cost and because centrums may move too much if 0815 many vertices in high-d 0816 */ 0817 #define qh_MAXnewcentrum 5 0818 0819 /*-<a href="qh-user.htm#TOC" 0820 >--------------------------------</a><a name="MAXnewmerges">-</a> 0821 0822 qh_MAXnewmerges 0823 if >n newmerges, qh_merge_nonconvex() calls qh_reducevertices_centrums. 0824 0825 notes: 0826 It is needed because postmerge can merge many facets at once 0827 */ 0828 #define qh_MAXnewmerges 2 0829 0830 /*-<a href="qh-user.htm#TOC" 0831 >--------------------------------</a><a name="RATIOconcavehorizon">-</a> 0832 0833 qh_RATIOconcavehorizon 0834 ratio of horizon vertex distance to max_outside for concave, twisted new facets in qh_test_nonsimplicial_merge 0835 if too small, end up with vertices far below merged facets 0836 */ 0837 #define qh_RATIOconcavehorizon 20.0 0838 0839 /*-<a href="qh-user.htm#TOC" 0840 >--------------------------------</a><a name="RATIOconvexmerge">-</a> 0841 0842 qh_RATIOconvexmerge 0843 ratio of vertex distance to qh.min_vertex for clearly convex new facets in qh_test_nonsimplicial_merge 0844 0845 notes: 0846 must be convex for MRGtwisted 0847 */ 0848 #define qh_RATIOconvexmerge 10.0 0849 0850 /*-<a href="qh-user.htm#TOC" 0851 >--------------------------------</a><a name="RATIOcoplanarapex">-</a> 0852 0853 qh_RATIOcoplanarapex 0854 ratio of best distance for coplanar apex vs. vertex merge in qh_getpinchedmerges 0855 0856 notes: 0857 A coplanar apex always works, while a vertex merge may fail 0858 */ 0859 #define qh_RATIOcoplanarapex 3.0 0860 0861 /*-<a href="qh-user.htm#TOC" 0862 >--------------------------------</a><a name="RATIOcoplanaroutside">-</a> 0863 0864 qh_RATIOcoplanaroutside 0865 qh.MAXoutside ratio to repartition a coplanar point in qh_partitioncoplanar and qh_check_maxout 0866 0867 notes: 0868 combines several tests, see qh_partitioncoplanar 0869 0870 */ 0871 #define qh_RATIOcoplanaroutside 30.0 0872 0873 /*-<a href="qh-user.htm#TOC" 0874 >--------------------------------</a><a name="RATIOmaxsimplex">-</a> 0875 0876 qh_RATIOmaxsimplex 0877 ratio of max determinate to estimated determinate for searching all points in qh_maxsimplex 0878 0879 notes: 0880 As each point is added to the simplex, the max determinate is should approximate the previous determinate * qh.MAXwidth 0881 If maxdet is significantly less, the simplex may not be full-dimensional. 0882 If so, all points are searched, stopping at 10 times qh_RATIOmaxsimplex 0883 */ 0884 #define qh_RATIOmaxsimplex 1.0e-3 0885 0886 /*-<a href="qh-user.htm#TOC" 0887 >--------------------------------</a><a name="RATIOnearinside">-</a> 0888 0889 qh_RATIOnearinside 0890 ratio of qh.NEARinside to qh.ONEmerge for retaining inside points for 0891 qh_check_maxout(). 0892 0893 notes: 0894 This is overkill since do not know the correct value. 0895 It effects whether 'Qc' reports all coplanar points 0896 Not used for 'd' since non-extreme points are coplanar, nearly incident points 0897 */ 0898 #define qh_RATIOnearinside 5 0899 0900 /*-<a href="qh-user.htm#TOC" 0901 >--------------------------------</a><a name="RATIOpinchedsubridge">-</a> 0902 0903 qh_RATIOpinchedsubridge 0904 ratio to qh.ONEmerge to accept vertices in qh_findbest_pinchedvertex 0905 skips search of neighboring vertices 0906 facet width may increase by this ratio 0907 */ 0908 #define qh_RATIOpinchedsubridge 10.0 0909 0910 /*-<a href="qh-user.htm#TOC" 0911 >--------------------------------</a><a name="RATIOtrypinched">-</a> 0912 0913 qh_RATIOtrypinched 0914 ratio to qh.ONEmerge to try qh_getpinchedmerges in qh_buildcone_mergepinched 0915 otherwise a duplicate ridge will increase facet width by this amount 0916 */ 0917 #define qh_RATIOtrypinched 4.0 0918 0919 /*-<a href="qh-user.htm#TOC" 0920 >--------------------------------</a><a name="RATIOtwisted">-</a> 0921 0922 qh_RATIOtwisted 0923 maximum ratio to qh.ONEmerge to merge twisted facets in qh_merge_twisted 0924 */ 0925 #define qh_RATIOtwisted 20.0 0926 0927 /*-<a href="qh-user.htm#TOC" 0928 >--------------------------------</a><a name="SEARCHdist">-</a> 0929 0930 qh_SEARCHdist 0931 When is a facet coplanar with the best facet? 0932 qh_findbesthorizon: all coplanar facets of the best facet need to be searched. 0933 increases minsearch if ischeckmax and more than 100 neighbors (is_5x_minsearch) 0934 See: 0935 qh_DISToutside -- when is a point clearly outside of a facet 0936 qh_SEARCHdist -- when is facet coplanar with the best facet? 0937 qh_USEfindbestnew -- when to use qh_findbestnew for qh_partitionpoint() 0938 */ 0939 #define qh_SEARCHdist ((qh_USEfindbestnew ? 2 : 1) * \ 0940 (qh max_outside + 2 * qh DISTround + fmax_( qh MINvisible, qh MAXcoplanar))); 0941 0942 /*-<a href="qh-user.htm#TOC" 0943 >--------------------------------</a><a name="USEfindbestnew">-</a> 0944 0945 qh_USEfindbestnew 0946 Always use qh_findbestnew for qh_partitionpoint, otherwise use 0947 qh_findbestnew if merged new facet or sharpnewfacets. 0948 0949 See: 0950 qh_DISToutside -- when is a point clearly outside of a facet 0951 qh_SEARCHdist -- when is facet coplanar with the best facet? 0952 qh_USEfindbestnew -- when to use qh_findbestnew for qh_partitionpoint() 0953 */ 0954 #define qh_USEfindbestnew (zzval_(Ztotmerge) > 50) 0955 0956 /*-<a href="qh-user.htm#TOC" 0957 >--------------------------------</a><a name="MAXnarrow">-</a> 0958 0959 qh_MAXnarrow 0960 max. cosine in initial hull that sets qh.NARROWhull 0961 0962 notes: 0963 If qh.NARROWhull, the initial partition does not make 0964 coplanar points. If narrow, a coplanar point can be 0965 coplanar to two facets of opposite orientations and 0966 distant from the exact convex hull. 0967 0968 Conservative estimate. Don't actually see problems until it is -1.0 0969 */ 0970 #define qh_MAXnarrow -0.99999999 0971 0972 /*-<a href="qh-user.htm#TOC" 0973 >--------------------------------</a><a name="WARNnarrow">-</a> 0974 0975 qh_WARNnarrow 0976 max. cosine in initial hull to warn about qh.NARROWhull 0977 0978 notes: 0979 this is a conservative estimate. 0980 Don't actually see problems until it is -1.0. See qh-impre.htm 0981 */ 0982 #define qh_WARNnarrow -0.999999999999999 0983 0984 /*-<a href="qh-user.htm#TOC" 0985 >--------------------------------</a><a name="WIDEcoplanar">-</a> 0986 0987 qh_WIDEcoplanar 0988 n*MAXcoplanar or n*MINvisible for a WIDEfacet 0989 0990 if vertex is further than qh.WIDEfacet from the hyperplane 0991 then its ridges are not counted in computing the area, and 0992 the facet's centrum is frozen. 0993 0994 notes: 0995 qh.WIDEfacet= max(qh.MAXoutside,qh_WIDEcoplanar*qh.MAXcoplanar, 0996 qh_WIDEcoplanar * qh.MINvisible); 0997 */ 0998 #define qh_WIDEcoplanar 6 0999 1000 /*-<a href="qh-user.htm#TOC" 1001 >--------------------------------</a><a name="WIDEduplicate">-</a> 1002 1003 qh_WIDEduplicate 1004 merge ratio for errexit from qh_forcedmerges due to duplicate ridge 1005 Override with option Q12-allow-wide 1006 1007 Notes: 1008 Merging a duplicate ridge can lead to very wide facets. 1009 */ 1010 #define qh_WIDEduplicate 100 1011 1012 /*-<a href="qh-user.htm#TOC" 1013 >--------------------------------</a><a name="WIDEdupridge">-</a> 1014 1015 qh_WIDEdupridge 1016 Merge ratio for selecting a forced dupridge merge 1017 1018 Notes: 1019 Merging a dupridge can lead to very wide facets. 1020 */ 1021 #define qh_WIDEdupridge 50 1022 1023 /*-<a href="qh-user.htm#TOC" 1024 >--------------------------------</a><a name="WIDEmaxoutside">-</a> 1025 1026 qh_WIDEmaxoutside 1027 Precision ratio for maximum increase for qh.max_outside in qh_check_maxout 1028 Precision errors while constructing the hull, may lead to very wide facets when checked in qh_check_maxout 1029 Nearly incident points in 4-d and higher is the most likely culprit 1030 Skip qh_check_maxout with 'Q5' (no-check-outer) 1031 Do not error with option 'Q12' (allow-wide) 1032 Do not warn with options 'Q12 Pp' 1033 */ 1034 #define qh_WIDEmaxoutside 100 1035 1036 /*-<a href="qh-user.htm#TOC" 1037 >--------------------------------</a><a name="WIDEmaxoutside2">-</a> 1038 1039 qh_WIDEmaxoutside2 1040 Precision ratio for maximum qh.max_outside in qh_check_maxout 1041 Skip qh_check_maxout with 'Q5' no-check-outer 1042 Do not error with option 'Q12' allow-wide 1043 */ 1044 #define qh_WIDEmaxoutside2 (10*qh_WIDEmaxoutside) 1045 1046 1047 /*-<a href="qh-user.htm#TOC" 1048 >--------------------------------</a><a name="WIDEpinched">-</a> 1049 1050 qh_WIDEpinched 1051 Merge ratio for distance between pinched vertices compared to current facet width for qh_getpinchedmerges and qh_next_vertexmerge 1052 Reports warning and merges duplicate ridges instead 1053 Enable these attempts with option Q14 merge-pinched-vertices 1054 1055 notes: 1056 Merging pinched vertices should prevent duplicate ridges (see qh_WIDEduplicate) 1057 Merging the duplicate ridges may be better than merging the pinched vertices 1058 Found up to 45x ratio for qh_pointdist -- for ((i=1; i<20; i++)); do rbox 175 C1,6e-13 t | qhull d T4 2>&1 | tee x.1 | grep -E 'QH|non-simplicial|Statis|pinched'; done 1059 Actual distance to facets is a third to a tenth of the qh_pointdist (T1) 1060 */ 1061 #define qh_WIDEpinched 100 1062 1063 /*-<a href="qh-user.htm#TOC" 1064 >--------------------------------</a><a name="ZEROdelaunay">-</a> 1065 1066 qh_ZEROdelaunay 1067 a zero Delaunay facet occurs for input sites coplanar with their convex hull 1068 the last normal coefficient of a zero Delaunay facet is within 1069 qh_ZEROdelaunay * qh.ANGLEround of 0 1070 1071 notes: 1072 qh_ZEROdelaunay does not allow for joggled input ('QJ'). 1073 1074 You can avoid zero Delaunay facets by surrounding the input with a box. 1075 1076 Use option 'PDk:-n' to explicitly define zero Delaunay facets 1077 k= dimension of input sites (e.g., 3 for 3-d Delaunay triangulation) 1078 n= the cutoff for zero Delaunay facets (e.g., 'PD3:-1e-12') 1079 */ 1080 #define qh_ZEROdelaunay 2 1081 1082 /*============================================================*/ 1083 /*============= Microsoft DevStudio ==========================*/ 1084 /*============================================================*/ 1085 1086 /* 1087 Finding Memory Leaks Using the CRT Library 1088 https://msdn.microsoft.com/en-us/library/x98tx3cf(v=vs.100).aspx 1089 1090 Reports enabled in qh_lib_check for Debug window and stderr 1091 1092 From 2005=>msvcr80d, 2010=>msvcr100d, 2012=>msvcr110d 1093 1094 Watch: {,,msvcr80d.dll}_crtBreakAlloc Value from {n} in the leak report 1095 _CrtSetBreakAlloc(689); // qh_lib_check() [global.c] 1096 1097 Examples 1098 http://free-cad.sourceforge.net/SrcDocu/d2/d7f/MemDebug_8cpp_source.html 1099 https://github.com/illlust/Game/blob/master/library/MemoryLeak.cpp 1100 */ 1101 #if 0 /* off (0) by default for QHULL_CRTDBG */ 1102 #define QHULL_CRTDBG 1103 #endif 1104 1105 #if defined(_MSC_VER) && defined(_DEBUG) && defined(QHULL_CRTDBG) 1106 #define _CRTDBG_MAP_ALLOC 1107 #include <stdlib.h> 1108 #include <crtdbg.h> 1109 #endif 1110 1111 #endif /* qh_DEFuser */
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |