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