user_r.h 36 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060
  1. /*<html><pre> -<a href="qh-user_r.htm"
  2. >-------------------------------</a><a name="TOP">-</a>
  3. user_r.h
  4. user redefinable constants
  5. for each source file, user_r.h is included first
  6. see qh-user_r.htm. see COPYING for copyright information.
  7. See user_r.c for sample code.
  8. before reading any code, review libqhull_r.h for data structure definitions
  9. Sections:
  10. ============= qhull library constants ======================
  11. ============= data types and configuration macros ==========
  12. ============= performance related constants ================
  13. ============= memory constants =============================
  14. ============= joggle constants =============================
  15. ============= conditional compilation ======================
  16. ============= merge constants ==============================
  17. ============= Microsoft DevStudio ==========================
  18. Code flags --
  19. NOerrors -- the code does not call qh_errexit()
  20. WARN64 -- the code may be incompatible with 64-bit pointers
  21. */
  22. #include <float.h>
  23. #include <limits.h>
  24. #include <time.h>
  25. #ifndef qhDEFuser
  26. #define qhDEFuser 1
  27. /* Derived from Qt's corelib/global/qglobal.h */
  28. #if !defined(SAG_COM) && !defined(__CYGWIN__) && (defined(WIN64) || defined(_WIN64) || defined(__WIN64__) || defined(WIN32) || defined(_WIN32) || defined(__WIN32__) || defined(__NT__))
  29. # define QHULL_OS_WIN
  30. #elif defined(__MWERKS__) && defined(__INTEL__) /* Metrowerks discontinued before the release of Intel Macs */
  31. # define QHULL_OS_WIN
  32. #endif
  33. /*============================================================*/
  34. /*============= qhull library constants ======================*/
  35. /*============================================================*/
  36. /*-<a href="qh-user_r.htm#TOC"
  37. >--------------------------------</a><a name="filenamelen">-</a>
  38. FILENAMElen -- max length for TI and TO filenames
  39. */
  40. #define qh_FILENAMElen 500
  41. /*-<a href="qh-user_r.htm#TOC"
  42. >--------------------------------</a><a name="msgcode">-</a>
  43. msgcode -- Unique message codes for qh_fprintf
  44. If add new messages, assign these values and increment in user.h and user_r.h
  45. See QhullError.h for 10000 error codes.
  46. Cannot use '0031' since it would be octal
  47. def counters = [31/32/33/38, 1067, 2113, 3079, 4097, 5006,
  48. 6429, 7027/7028/7035/7068/7070/7102, 8163, 9428, 10000, 11034]
  49. See: qh_ERR* [libqhull_r.h]
  50. */
  51. #define MSG_TRACE0 0 /* always include if logging ('Tn') */
  52. #define MSG_TRACE1 1000
  53. #define MSG_TRACE2 2000
  54. #define MSG_TRACE3 3000
  55. #define MSG_TRACE4 4000
  56. #define MSG_TRACE5 5000
  57. #define MSG_ERROR 6000 /* errors written to qh.ferr */
  58. #define MSG_WARNING 7000
  59. #define MSG_STDERR 8000 /* log messages Written to qh.ferr */
  60. #define MSG_OUTPUT 9000
  61. #define MSG_QHULL_ERROR 10000 /* errors thrown by QhullError.cpp (QHULLlastError is in QhullError.h) */
  62. #define MSG_FIX 11000 /* Document as 'QH11... FIX: ...' */
  63. #define MSG_MAXLEN 3000 /* qh_printhelp_degenerate() in user_r.c */
  64. /*-<a href="qh-user_r.htm#TOC"
  65. >--------------------------------</a><a name="qh_OPTIONline">-</a>
  66. qh_OPTIONline -- max length of an option line 'FO'
  67. */
  68. #define qh_OPTIONline 80
  69. /*============================================================*/
  70. /*============= data types and configuration macros ==========*/
  71. /*============================================================*/
  72. /*-<a href="qh-user_r.htm#TOC"
  73. >--------------------------------</a><a name="realT">-</a>
  74. realT
  75. set the size of floating point numbers
  76. qh_REALdigits
  77. maximimum number of significant digits
  78. qh_REAL_1, qh_REAL_2n, qh_REAL_3n
  79. format strings for printf
  80. qh_REALmax, qh_REALmin
  81. maximum and minimum (near zero) values
  82. qh_REALepsilon
  83. machine roundoff. Maximum roundoff error for addition and multiplication.
  84. notes:
  85. Select whether to store floating point numbers in single precision (float)
  86. or double precision (double).
  87. Use 'float' to save about 8% in time and 25% in space. This is particularly
  88. helpful if high-d where convex hulls are space limited. Using 'float' also
  89. reduces the printed size of Qhull's output since numbers have 8 digits of
  90. precision.
  91. Use 'double' when greater arithmetic precision is needed. This is needed
  92. for Delaunay triangulations and Voronoi diagrams when you are not merging
  93. facets.
  94. If 'double' gives insufficient precision, your data probably includes
  95. degeneracies. If so you should use facet merging (done by default)
  96. or exact arithmetic (see imprecision section of manual, qh-impre.htm).
  97. You may also use option 'Po' to force output despite precision errors.
  98. You may use 'long double', but many format statements need to be changed
  99. and you may need a 'long double' square root routine. S. Grundmann
  100. (sg@eeiwzb.et.tu-dresden.de) has done this. He reports that the code runs
  101. much slower with little gain in precision.
  102. WARNING: on some machines, int f(){realT a= REALmax;return (a == REALmax);}
  103. returns False. Use (a > REALmax/2) instead of (a == REALmax).
  104. REALfloat = 1 all numbers are 'float' type
  105. = 0 all numbers are 'double' type
  106. */
  107. #define REALfloat 0
  108. #if (REALfloat == 1)
  109. #define realT float
  110. #define REALmax FLT_MAX
  111. #define REALmin FLT_MIN
  112. #define REALepsilon FLT_EPSILON
  113. #define qh_REALdigits 8 /* maximum number of significant digits */
  114. #define qh_REAL_1 "%6.8g "
  115. #define qh_REAL_2n "%6.8g %6.8g\n"
  116. #define qh_REAL_3n "%6.8g %6.8g %6.8g\n"
  117. #elif (REALfloat == 0)
  118. #define realT double
  119. #define REALmax DBL_MAX
  120. #define REALmin DBL_MIN
  121. #define REALepsilon DBL_EPSILON
  122. #define qh_REALdigits 16 /* maximum number of significant digits */
  123. #define qh_REAL_1 "%6.16g "
  124. #define qh_REAL_2n "%6.16g %6.16g\n"
  125. #define qh_REAL_3n "%6.16g %6.16g %6.16g\n"
  126. #else
  127. #error unknown float option
  128. #endif
  129. /*-<a href="qh-user_r.htm#TOC"
  130. >--------------------------------</a><a name="countT">-</a>
  131. countT
  132. The type for counts and identifiers (e.g., the number of points, vertex identifiers)
  133. Currently used by C++ code-only. Decided against using it for setT because most sets are small.
  134. Defined as 'int' for C-code compatibility and QH11026
  135. QH11026 FIX: countT may be defined as a 'unsigned int', but several code issues need to be solved first. See countT in Changes.txt
  136. */
  137. #ifndef DEFcountT
  138. #define DEFcountT 1
  139. typedef int countT;
  140. #endif
  141. #define COUNTmax INT_MAX
  142. /*-<a href="qh-user_r.htm#TOC"
  143. >--------------------------------</a><a name="qh_POINTSmax">-</a>
  144. qh_POINTSmax
  145. Maximum number of points for qh.num_points and point allocation in qh_readpoints
  146. */
  147. #define qh_POINTSmax (INT_MAX-16)
  148. /*-<a href="qh-user_r.htm#TOC"
  149. >--------------------------------</a><a name="CPUclock">-</a>
  150. qh_CPUclock
  151. define the clock() function for reporting the total time spent by Qhull
  152. returns CPU ticks as a 'long int'
  153. qh_CPUclock is only used for reporting the total time spent by Qhull
  154. qh_SECticks
  155. the number of clock ticks per second
  156. notes:
  157. looks for CLOCKS_PER_SEC, CLOCKS_PER_SECOND, or assumes microseconds
  158. to define a custom clock, set qh_CLOCKtype to 0
  159. if your system does not use clock() to return CPU ticks, replace
  160. qh_CPUclock with the corresponding function. It is converted
  161. to 'unsigned long' to prevent wrap-around during long runs. By default,
  162. <time.h> defines clock_t as 'long'
  163. Set qh_CLOCKtype to
  164. 1 for CLOCKS_PER_SEC, CLOCKS_PER_SECOND, or microsecond
  165. Note: may fail if more than 1 hour elapsed time
  166. 2 use qh_clock() with POSIX times() (see global_r.c)
  167. */
  168. #define qh_CLOCKtype 1 /* change to the desired number */
  169. #if (qh_CLOCKtype == 1)
  170. #if defined(CLOCKS_PER_SECOND)
  171. #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock, may be converted to approximate double */
  172. #define qh_SECticks CLOCKS_PER_SECOND
  173. #elif defined(CLOCKS_PER_SEC)
  174. #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock, may be converted to approximate double */
  175. #define qh_SECticks CLOCKS_PER_SEC
  176. #elif defined(CLK_TCK)
  177. #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock, may be converted to approximate double */
  178. #define qh_SECticks CLK_TCK
  179. #else
  180. #define qh_CPUclock ((unsigned long)clock()) /* return CPU clock, may be converted to approximate double */
  181. #define qh_SECticks 1E6
  182. #endif
  183. #elif (qh_CLOCKtype == 2)
  184. #define qh_CPUclock qh_clock() /* return CPU clock, may be converted to approximate double */
  185. #define qh_SECticks 100
  186. #else /* qh_CLOCKtype == ? */
  187. #error unknown clock option
  188. #endif
  189. /*-<a href="qh-user_r.htm#TOC"
  190. >--------------------------------</a><a name="RANDOM">-</a>
  191. qh_RANDOMtype, qh_RANDOMmax, qh_RANDOMseed
  192. define random number generator
  193. qh_RANDOMint generates a random integer between 0 and qh_RANDOMmax.
  194. qh_RANDOMseed sets the random number seed for qh_RANDOMint
  195. Set qh_RANDOMtype (default 5) to:
  196. 1 for random() with 31 bits (UCB)
  197. 2 for rand() with RAND_MAX or 15 bits (system 5)
  198. 3 for rand() with 31 bits (Sun)
  199. 4 for lrand48() with 31 bits (Solaris)
  200. 5 for qh_rand(qh) with 31 bits (included with Qhull, requires 'qh')
  201. notes:
  202. Random numbers are used by rbox to generate point sets. Random
  203. numbers are used by Qhull to rotate the input ('QRn' option),
  204. simulate a randomized algorithm ('Qr' option), and to simulate
  205. roundoff errors ('Rn' option).
  206. Random number generators differ between systems. Most systems provide
  207. rand() but the period varies. The period of rand() is not critical
  208. since qhull does not normally use random numbers.
  209. The default generator is Park & Miller's minimal standard random
  210. number generator [CACM 31:1195 '88]. It is included with Qhull.
  211. If qh_RANDOMmax is wrong, qhull will report a warning and Geomview
  212. output will likely be invisible.
  213. */
  214. #define qh_RANDOMtype 5 /* *** change to the desired number *** */
  215. #if (qh_RANDOMtype == 1)
  216. #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, random()/MAX */
  217. #define qh_RANDOMint random()
  218. #define qh_RANDOMseed_(qh, seed) srandom(seed);
  219. #elif (qh_RANDOMtype == 2)
  220. #ifdef RAND_MAX
  221. #define qh_RANDOMmax ((realT)RAND_MAX)
  222. #else
  223. #define qh_RANDOMmax ((realT)32767) /* 15 bits (System 5) */
  224. #endif
  225. #define qh_RANDOMint rand()
  226. #define qh_RANDOMseed_(qh, seed) srand((unsigned int)seed);
  227. #elif (qh_RANDOMtype == 3)
  228. #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, Sun */
  229. #define qh_RANDOMint rand()
  230. #define qh_RANDOMseed_(qh, seed) srand((unsigned int)seed);
  231. #elif (qh_RANDOMtype == 4)
  232. #define qh_RANDOMmax ((realT)0x7fffffffUL) /* 31 bits, lrand38()/MAX */
  233. #define qh_RANDOMint lrand48()
  234. #define qh_RANDOMseed_(qh, seed) srand48(seed);
  235. #elif (qh_RANDOMtype == 5) /* 'qh' is an implicit parameter */
  236. #define qh_RANDOMmax ((realT)2147483646UL) /* 31 bits, qh_rand/MAX */
  237. #define qh_RANDOMint qh_rand(qh)
  238. #define qh_RANDOMseed_(qh, seed) qh_srand(qh, seed);
  239. /* unlike rand(), never returns 0 */
  240. #else
  241. #error: unknown random option
  242. #endif
  243. /*-<a href="qh-user_r.htm#TOC"
  244. >--------------------------------</a><a name="ORIENTclock">-</a>
  245. qh_ORIENTclock
  246. 0 for inward pointing normals by Geomview convention
  247. */
  248. #define qh_ORIENTclock 0
  249. /*-<a href="qh-user_r.htm#TOC"
  250. >--------------------------------</a><a name="RANDOMdist">-</a>
  251. qh_RANDOMdist
  252. define for random perturbation of qh_distplane and qh_setfacetplane (qh.RANDOMdist, 'QRn')
  253. For testing qh.DISTround. Qhull should not depend on computations always producing the same roundoff error
  254. #define qh_RANDOMdist 1e-13
  255. */
  256. /*============================================================*/
  257. /*============= joggle constants =============================*/
  258. /*============================================================*/
  259. /*-<a href="qh-user_r.htm#TOC"
  260. >--------------------------------</a><a name="JOGGLEdefault">-</a>
  261. qh_JOGGLEdefault
  262. default qh.JOGGLEmax is qh.DISTround * qh_JOGGLEdefault
  263. notes:
  264. rbox s r 100 | qhull QJ1e-15 QR0 generates 90% faults at distround 7e-16
  265. rbox s r 100 | qhull QJ1e-14 QR0 generates 70% faults
  266. rbox s r 100 | qhull QJ1e-13 QR0 generates 35% faults
  267. rbox s r 100 | qhull QJ1e-12 QR0 generates 8% faults
  268. rbox s r 100 | qhull QJ1e-11 QR0 generates 1% faults
  269. rbox s r 100 | qhull QJ1e-10 QR0 generates 0% faults
  270. rbox 1000 W0 | qhull QJ1e-12 QR0 generates 86% faults
  271. rbox 1000 W0 | qhull QJ1e-11 QR0 generates 20% faults
  272. rbox 1000 W0 | qhull QJ1e-10 QR0 generates 2% faults
  273. the later have about 20 points per facet, each of which may interfere
  274. pick a value large enough to avoid retries on most inputs
  275. */
  276. #define qh_JOGGLEdefault 30000.0
  277. /*-<a href="qh-user_r.htm#TOC"
  278. >--------------------------------</a><a name="JOGGLEincrease">-</a>
  279. qh_JOGGLEincrease
  280. factor to increase qh.JOGGLEmax on qh_JOGGLEretry or qh_JOGGLEagain
  281. */
  282. #define qh_JOGGLEincrease 10.0
  283. /*-<a href="qh-user_r.htm#TOC"
  284. >--------------------------------</a><a name="JOGGLEretry">-</a>
  285. qh_JOGGLEretry
  286. if ZZretry = qh_JOGGLEretry, increase qh.JOGGLEmax
  287. notes:
  288. try twice at the original value in case of bad luck the first time
  289. */
  290. #define qh_JOGGLEretry 2
  291. /*-<a href="qh-user_r.htm#TOC"
  292. >--------------------------------</a><a name="JOGGLEagain">-</a>
  293. qh_JOGGLEagain
  294. every following qh_JOGGLEagain, increase qh.JOGGLEmax
  295. notes:
  296. 1 is OK since it's already failed qh_JOGGLEretry times
  297. */
  298. #define qh_JOGGLEagain 1
  299. /*-<a href="qh-user_r.htm#TOC"
  300. >--------------------------------</a><a name="JOGGLEmaxincrease">-</a>
  301. qh_JOGGLEmaxincrease
  302. maximum qh.JOGGLEmax due to qh_JOGGLEincrease
  303. relative to qh.MAXwidth
  304. notes:
  305. qh.joggleinput will retry at this value until qh_JOGGLEmaxretry
  306. */
  307. #define qh_JOGGLEmaxincrease 1e-2
  308. /*-<a href="qh-user_r.htm#TOC"
  309. >--------------------------------</a><a name="JOGGLEmaxretry">-</a>
  310. qh_JOGGLEmaxretry
  311. stop after qh_JOGGLEmaxretry attempts
  312. */
  313. #define qh_JOGGLEmaxretry 50
  314. /*============================================================*/
  315. /*============= performance related constants ================*/
  316. /*============================================================*/
  317. /*-<a href="qh-user_r.htm#TOC"
  318. >--------------------------------</a><a name="HASHfactor">-</a>
  319. qh_HASHfactor
  320. total hash slots / used hash slots. Must be at least 1.1.
  321. notes:
  322. =2 for at worst 50% occupancy for qh.hash_table and normally 25% occupancy
  323. */
  324. #define qh_HASHfactor 2
  325. /*-<a href="qh-user_r.htm#TOC"
  326. >--------------------------------</a><a name="VERIFYdirect">-</a>
  327. qh_VERIFYdirect
  328. with 'Tv' verify all points against all facets if op count is smaller
  329. notes:
  330. if greater, calls qh_check_bestdist() instead
  331. */
  332. #define qh_VERIFYdirect 1000000
  333. /*-<a href="qh-user_r.htm#TOC"
  334. >--------------------------------</a><a name="INITIALsearch">-</a>
  335. qh_INITIALsearch
  336. if qh_INITIALmax, search points up to this dimension
  337. */
  338. #define qh_INITIALsearch 6
  339. /*-<a href="qh-user_r.htm#TOC"
  340. >--------------------------------</a><a name="INITIALmax">-</a>
  341. qh_INITIALmax
  342. if dim >= qh_INITIALmax, use min/max coordinate points for initial simplex
  343. notes:
  344. from points with non-zero determinants
  345. use option 'Qs' to override (much slower)
  346. */
  347. #define qh_INITIALmax 8
  348. /*============================================================*/
  349. /*============= memory constants =============================*/
  350. /*============================================================*/
  351. /*-<a href="qh-user_r.htm#TOC"
  352. >--------------------------------</a><a name="MEMalign">-</a>
  353. qh_MEMalign
  354. memory alignment for qh_meminitbuffers() in global_r.c
  355. notes:
  356. to avoid bus errors, memory allocation must consider alignment requirements.
  357. malloc() automatically takes care of alignment. Since mem_r.c manages
  358. its own memory, we need to explicitly specify alignment in
  359. qh_meminitbuffers().
  360. A safe choice is sizeof(double). sizeof(float) may be used if doubles
  361. do not occur in data structures and pointers are the same size. Be careful
  362. of machines (e.g., DEC Alpha) with large pointers.
  363. If using gcc, best alignment is [fmax_() is defined in geom_r.h]
  364. #define qh_MEMalign fmax_(__alignof__(realT),__alignof__(void *))
  365. */
  366. #define qh_MEMalign ((int)(fmax_(sizeof(realT), sizeof(void *))))
  367. /*-<a href="qh-user_r.htm#TOC"
  368. >--------------------------------</a><a name="MEMbufsize">-</a>
  369. qh_MEMbufsize
  370. size of additional memory buffers
  371. notes:
  372. used for qh_meminitbuffers() in global_r.c
  373. */
  374. #define qh_MEMbufsize 0x10000 /* allocate 64K memory buffers */
  375. /*-<a href="qh-user_r.htm#TOC"
  376. >--------------------------------</a><a name="MEMinitbuf">-</a>
  377. qh_MEMinitbuf
  378. size of initial memory buffer
  379. notes:
  380. use for qh_meminitbuffers() in global_r.c
  381. */
  382. #define qh_MEMinitbuf 0x20000 /* initially allocate 128K buffer */
  383. /*-<a href="qh-user_r.htm#TOC"
  384. >--------------------------------</a><a name="INFINITE">-</a>
  385. qh_INFINITE
  386. on output, indicates Voronoi center at infinity
  387. */
  388. #define qh_INFINITE -10.101
  389. /*-<a href="qh-user_r.htm#TOC"
  390. >--------------------------------</a><a name="DEFAULTbox">-</a>
  391. qh_DEFAULTbox
  392. default box size (Geomview expects 0.5)
  393. qh_DEFAULTbox
  394. default box size for integer coorindate (rbox only)
  395. */
  396. #define qh_DEFAULTbox 0.5
  397. #define qh_DEFAULTzbox 1e6
  398. /*============================================================*/
  399. /*============= conditional compilation ======================*/
  400. /*============================================================*/
  401. /*-<a href="qh-user_r.htm#TOC"
  402. >--------------------------------</a><a name="compiler">-</a>
  403. __cplusplus
  404. defined by C++ compilers
  405. __MSC_VER
  406. defined by Microsoft Visual C++
  407. __MWERKS__ && __INTEL__
  408. defined by Metrowerks when compiling for Windows (not Intel-based Macintosh)
  409. __MWERKS__ && __POWERPC__
  410. defined by Metrowerks when compiling for PowerPC-based Macintosh
  411. __STDC__
  412. defined for strict ANSI C
  413. */
  414. /*-<a href="qh-user_r.htm#TOC"
  415. >--------------------------------</a><a name="COMPUTEfurthest">-</a>
  416. qh_COMPUTEfurthest
  417. compute furthest distance to an outside point instead of storing it with the facet
  418. =1 to compute furthest
  419. notes:
  420. computing furthest saves memory but costs time
  421. about 40% more distance tests for partitioning
  422. removes facet->furthestdist
  423. */
  424. #define qh_COMPUTEfurthest 0
  425. /*-<a href="qh-user_r.htm#TOC"
  426. >--------------------------------</a><a name="KEEPstatistics">-</a>
  427. qh_KEEPstatistics
  428. =0 removes most of statistic gathering and reporting
  429. notes:
  430. if 0, code size is reduced by about 4%.
  431. */
  432. #define qh_KEEPstatistics 1
  433. /*-<a href="qh-user_r.htm#TOC"
  434. >--------------------------------</a><a name="MAXoutside">-</a>
  435. qh_MAXoutside
  436. record outer plane for each facet
  437. =1 to record facet->maxoutside
  438. notes:
  439. this takes a realT per facet and slightly slows down qhull
  440. it produces better outer planes for geomview output
  441. */
  442. #define qh_MAXoutside 1
  443. /*-<a href="qh-user_r.htm#TOC"
  444. >--------------------------------</a><a name="NOmerge">-</a>
  445. qh_NOmerge
  446. disables facet merging if defined
  447. For MSVC compiles, use qhull_r-exports-nomerge.def instead of qhull_r-exports.def
  448. notes:
  449. This saves about 25% space, 30% space in combination with qh_NOtrace,
  450. and 36% with qh_NOtrace and qh_KEEPstatistics 0
  451. Unless option 'Q0' is used
  452. qh_NOmerge sets 'QJ' to avoid precision errors
  453. see:
  454. <a href="mem_r.h#NOmem">qh_NOmem</a> in mem_r.h
  455. see user_r.c/user_eg.c for removing io_r.o
  456. #define qh_NOmerge
  457. */
  458. /*-<a href="qh-user_r.htm#TOC"
  459. >--------------------------------</a><a name="NOtrace">-</a>
  460. qh_NOtrace
  461. no tracing if defined
  462. disables 'Tn', 'TMn', 'TPn' and 'TWn'
  463. override with 'Qw' for qh_addpoint tracing and various other items
  464. notes:
  465. This saves about 15% space.
  466. Removes all traceN((...)) code and substantial sections of qh.IStracing code
  467. #define qh_NOtrace
  468. */
  469. #if 0 /* sample code */
  470. exitcode= qh_new_qhull(qhT *qh, dim, numpoints, points, ismalloc,
  471. flags, outfile, errfile);
  472. qh_freeqhull(qhT *qh, !qh_ALL); /* frees long memory used by second call */
  473. qh_memfreeshort(qhT *qh, &curlong, &totlong); /* frees short memory and memory allocator */
  474. #endif
  475. /*-<a href="qh-user_r.htm#TOC"
  476. >--------------------------------</a><a name="QUICKhelp">-</a>
  477. qh_QUICKhelp
  478. =1 to use abbreviated help messages, e.g., for degenerate inputs
  479. */
  480. #define qh_QUICKhelp 0
  481. /*============================================================*/
  482. /*============= merge constants ==============================*/
  483. /*============================================================*/
  484. /*
  485. These constants effect facet merging. You probably will not need
  486. to modify them. They effect the performance of facet merging.
  487. */
  488. /*-<a href="qh-user_r.htm#TOC"
  489. >--------------------------------</a><a name="BESTcentrum">-</a>
  490. qh_BESTcentrum
  491. if > 2*dim+n vertices, qh_findbestneighbor() tests centrums (faster)
  492. else, qh_findbestneighbor() tests all vertices (much better merges)
  493. qh_BESTcentrum2
  494. if qh_BESTcentrum2 * DIM3 + BESTcentrum < #vertices tests centrums
  495. */
  496. #define qh_BESTcentrum 20
  497. #define qh_BESTcentrum2 2
  498. /*-<a href="qh-user_r.htm#TOC"
  499. >--------------------------------</a><a name="BESTnonconvex">-</a>
  500. qh_BESTnonconvex
  501. if > dim+n neighbors, qh_findbestneighbor() tests nonconvex ridges.
  502. notes:
  503. It is needed because qh_findbestneighbor is slow for large facets
  504. */
  505. #define qh_BESTnonconvex 15
  506. /*-<a href="qh-user_r.htm#TOC"
  507. >--------------------------------</a><a name="COPLANARratio">-</a>
  508. qh_COPLANARratio
  509. for 3-d+ merging, qh.MINvisible is n*premerge_centrum
  510. notes:
  511. for non-merging, it's DISTround
  512. */
  513. #define qh_COPLANARratio 3
  514. /*-<a href="qh-user_r.htm#TOC"
  515. >--------------------------------</a><a name="DIMmergeVertex">-</a>
  516. qh_DIMmergeVertex
  517. max dimension for vertex merging (it is not effective in high-d)
  518. */
  519. #define qh_DIMmergeVertex 6
  520. /*-<a href="qh-user_r.htm#TOC"
  521. >--------------------------------</a><a name="DIMreduceBuild">-</a>
  522. qh_DIMreduceBuild
  523. max dimension for vertex reduction during build (slow in high-d)
  524. */
  525. #define qh_DIMreduceBuild 5
  526. /*-<a href="qh-user_r.htm#TOC"
  527. >--------------------------------</a><a name="DISToutside">-</a>
  528. qh_DISToutside
  529. When is a point clearly outside of a facet?
  530. Stops search in qh_findbestnew or qh_partitionall
  531. qh_findbest uses qh.MINoutside since since it is only called if no merges.
  532. notes:
  533. 'Qf' always searches for best facet
  534. if !qh.MERGING, same as qh.MINoutside.
  535. if qh_USEfindbestnew, increase value since neighboring facets may be ill-behaved
  536. [Note: Zdelvertextot occurs normally with interior points]
  537. RBOX 1000 s Z1 G1e-13 t1001188774 | QHULL Tv
  538. When there is a sharp edge, need to move points to a
  539. clearly good facet; otherwise may be lost in another partitioning.
  540. if too big then O(n^2) behavior for partitioning in cone
  541. if very small then important points not processed
  542. Needed in qh_partitionall for
  543. RBOX 1000 s Z1 G1e-13 t1001032651 | QHULL Tv
  544. Needed in qh_findbestnew for many instances of
  545. RBOX 1000 s Z1 G1e-13 t | QHULL Tv
  546. See:
  547. qh_DISToutside -- when is a point clearly outside of a facet
  548. qh_SEARCHdist -- when is facet coplanar with the best facet?
  549. qh_USEfindbestnew -- when to use qh_findbestnew for qh_partitionpoint()
  550. */
  551. #define qh_DISToutside ((qh_USEfindbestnew ? 2 : 1) * \
  552. fmax_((qh->MERGING ? 2 : 1)*qh->MINoutside, qh->max_outside))
  553. /*-<a href="qh-user_r.htm#TOC"
  554. >--------------------------------</a><a name="MAXcheckpoint">-</a>
  555. qh_MAXcheckpoint
  556. Report up to qh_MAXcheckpoint errors per facet in qh_check_point ('Tv')
  557. */
  558. #define qh_MAXcheckpoint 10
  559. /*-<a href="qh-user_r.htm#TOC"
  560. >--------------------------------</a><a name="MAXcoplanarcentrum">-</a>
  561. qh_MAXcoplanarcentrum
  562. if pre-merging with qh.MERGEexact ('Qx') and f.nummerge > qh_MAXcoplanarcentrum
  563. use f.maxoutside instead of qh.centrum_radius for coplanarity testing
  564. notes:
  565. see qh_test_nonsimplicial_merges
  566. with qh.MERGEexact, a coplanar ridge is ignored until post-merging
  567. otherwise a large facet with many merges may take all the facets
  568. */
  569. #define qh_MAXcoplanarcentrum 10
  570. /*-<a href="qh-user_r.htm#TOC"
  571. >--------------------------------</a><a name="MAXnewcentrum">-</a>
  572. qh_MAXnewcentrum
  573. if <= dim+n vertices (n approximates the number of merges),
  574. reset the centrum in qh_updatetested() and qh_mergecycle_facets()
  575. notes:
  576. needed to reduce cost and because centrums may move too much if
  577. many vertices in high-d
  578. */
  579. #define qh_MAXnewcentrum 5
  580. /*-<a href="qh-user_r.htm#TOC"
  581. >--------------------------------</a><a name="MAXnewmerges">-</a>
  582. qh_MAXnewmerges
  583. if >n newmerges, qh_merge_nonconvex() calls qh_reducevertices_centrums.
  584. notes:
  585. It is needed because postmerge can merge many facets at once
  586. */
  587. #define qh_MAXnewmerges 2
  588. /*-<a href="qh-user_r.htm#TOC"
  589. >--------------------------------</a><a name="RATIOconcavehorizon">-</a>
  590. qh_RATIOconcavehorizon
  591. ratio of horizon vertex distance to max_outside for concave, twisted new facets in qh_test_nonsimplicial_merge
  592. if too small, end up with vertices far below merged facets
  593. */
  594. #define qh_RATIOconcavehorizon 20.0
  595. /*-<a href="qh-user_r.htm#TOC"
  596. >--------------------------------</a><a name="RATIOconvexmerge">-</a>
  597. qh_RATIOconvexmerge
  598. ratio of vertex distance to qh.min_vertex for clearly convex new facets in qh_test_nonsimplicial_merge
  599. notes:
  600. must be convex for MRGtwisted
  601. */
  602. #define qh_RATIOconvexmerge 10.0
  603. /*-<a href="qh-user_r.htm#TOC"
  604. >--------------------------------</a><a name="RATIOcoplanarapex">-</a>
  605. qh_RATIOcoplanarapex
  606. ratio of best distance for coplanar apex vs. vertex merge in qh_getpinchedmerges
  607. notes:
  608. A coplanar apex always works, while a vertex merge may fail
  609. */
  610. #define qh_RATIOcoplanarapex 3.0
  611. /*-<a href="qh-user_r.htm#TOC"
  612. >--------------------------------</a><a name="RATIOcoplanaroutside">-</a>
  613. qh_RATIOcoplanaroutside
  614. qh.MAXoutside ratio to repartition a coplanar point in qh_partitioncoplanar and qh_check_maxout
  615. notes:
  616. combines several tests, see qh_partitioncoplanar
  617. */
  618. #define qh_RATIOcoplanaroutside 30.0
  619. /*-<a href="qh-user_r.htm#TOC"
  620. >--------------------------------</a><a name="RATIOmaxsimplex">-</a>
  621. qh_RATIOmaxsimplex
  622. ratio of max determinate to estimated determinate for searching all points in qh_maxsimplex
  623. notes:
  624. As each point is added to the simplex, the max determinate is should approximate the previous determinate * qh.MAXwidth
  625. If maxdet is significantly less, the simplex may not be full-dimensional.
  626. If so, all points are searched, stopping at 10 times qh_RATIOmaxsimplex
  627. */
  628. #define qh_RATIOmaxsimplex 1.0e-3
  629. /*-<a href="qh-user_r.htm#TOC"
  630. >--------------------------------</a><a name="RATIOnearinside">-</a>
  631. qh_RATIOnearinside
  632. ratio of qh.NEARinside to qh.ONEmerge for retaining inside points for
  633. qh_check_maxout().
  634. notes:
  635. This is overkill since do not know the correct value.
  636. It effects whether 'Qc' reports all coplanar points
  637. Not used for 'd' since non-extreme points are coplanar, nearly incident points
  638. */
  639. #define qh_RATIOnearinside 5
  640. /*-<a href="qh-user_r.htm#TOC"
  641. >--------------------------------</a><a name="RATIOpinchedsubridge">-</a>
  642. qh_RATIOpinchedsubridge
  643. ratio to qh.ONEmerge to accept vertices in qh_findbest_pinchedvertex
  644. skips search of neighboring vertices
  645. facet width may increase by this ratio
  646. */
  647. #define qh_RATIOpinchedsubridge 10.0
  648. /*-<a href="qh-user_r.htm#TOC"
  649. >--------------------------------</a><a name="RATIOtrypinched">-</a>
  650. qh_RATIOtrypinched
  651. ratio to qh.ONEmerge to try qh_getpinchedmerges in qh_buildcone_mergepinched
  652. otherwise a duplicate ridge will increase facet width by this amount
  653. */
  654. #define qh_RATIOtrypinched 4.0
  655. /*-<a href="qh-user_r.htm#TOC"
  656. >--------------------------------</a><a name="RATIOtwisted">-</a>
  657. qh_RATIOtwisted
  658. maximum ratio to qh.ONEmerge to merge twisted facets in qh_merge_twisted
  659. */
  660. #define qh_RATIOtwisted 20.0
  661. /*-<a href="qh-user_r.htm#TOC"
  662. >--------------------------------</a><a name="SEARCHdist">-</a>
  663. qh_SEARCHdist
  664. When is a facet coplanar with the best facet?
  665. qh_findbesthorizon: all coplanar facets of the best facet need to be searched.
  666. increases minsearch if ischeckmax and more than 100 neighbors (is_5x_minsearch)
  667. See:
  668. qh_DISToutside -- when is a point clearly outside of a facet
  669. qh_SEARCHdist -- when is facet coplanar with the best facet?
  670. qh_USEfindbestnew -- when to use qh_findbestnew for qh_partitionpoint()
  671. */
  672. #define qh_SEARCHdist ((qh_USEfindbestnew ? 2 : 1) * \
  673. (qh->max_outside + 2 * qh->DISTround + fmax_( qh->MINvisible, qh->MAXcoplanar)));
  674. /*-<a href="qh-user_r.htm#TOC"
  675. >--------------------------------</a><a name="USEfindbestnew">-</a>
  676. qh_USEfindbestnew
  677. Always use qh_findbestnew for qh_partitionpoint, otherwise use
  678. qh_findbestnew if merged new facet or sharpnewfacets.
  679. See:
  680. qh_DISToutside -- when is a point clearly outside of a facet
  681. qh_SEARCHdist -- when is facet coplanar with the best facet?
  682. qh_USEfindbestnew -- when to use qh_findbestnew for qh_partitionpoint()
  683. */
  684. #define qh_USEfindbestnew (zzval_(Ztotmerge) > 50)
  685. /*-<a href="qh-user_r.htm#TOC"
  686. >--------------------------------</a><a name="MAXnarrow">-</a>
  687. qh_MAXnarrow
  688. max. cosine in initial hull that sets qh.NARROWhull
  689. notes:
  690. If qh.NARROWhull, the initial partition does not make
  691. coplanar points. If narrow, a coplanar point can be
  692. coplanar to two facets of opposite orientations and
  693. distant from the exact convex hull.
  694. Conservative estimate. Don't actually see problems until it is -1.0
  695. */
  696. #define qh_MAXnarrow -0.99999999
  697. /*-<a href="qh-user_r.htm#TOC"
  698. >--------------------------------</a><a name="WARNnarrow">-</a>
  699. qh_WARNnarrow
  700. max. cosine in initial hull to warn about qh.NARROWhull
  701. notes:
  702. this is a conservative estimate.
  703. Don't actually see problems until it is -1.0. See qh-impre.htm
  704. */
  705. #define qh_WARNnarrow -0.999999999999999
  706. /*-<a href="qh-user_r.htm#TOC"
  707. >--------------------------------</a><a name="WIDEcoplanar">-</a>
  708. qh_WIDEcoplanar
  709. n*MAXcoplanar or n*MINvisible for a WIDEfacet
  710. if vertex is further than qh.WIDEfacet from the hyperplane
  711. then its ridges are not counted in computing the area, and
  712. the facet's centrum is frozen.
  713. notes:
  714. qh.WIDEfacet= max(qh.MAXoutside,qh_WIDEcoplanar*qh.MAXcoplanar,
  715. qh_WIDEcoplanar * qh.MINvisible);
  716. */
  717. #define qh_WIDEcoplanar 6
  718. /*-<a href="qh-user_r.htm#TOC"
  719. >--------------------------------</a><a name="WIDEduplicate">-</a>
  720. qh_WIDEduplicate
  721. merge ratio for errexit from qh_forcedmerges due to duplicate ridge
  722. Override with option Q12-allow-wide
  723. Notes:
  724. Merging a duplicate ridge can lead to very wide facets.
  725. */
  726. #define qh_WIDEduplicate 100
  727. /*-<a href="qh-user_r.htm#TOC"
  728. >--------------------------------</a><a name="WIDEdupridge">-</a>
  729. qh_WIDEdupridge
  730. Merge ratio for selecting a forced dupridge merge
  731. Notes:
  732. Merging a dupridge can lead to very wide facets.
  733. */
  734. #define qh_WIDEdupridge 50
  735. /*-<a href="qh-user_r.htm#TOC"
  736. >--------------------------------</a><a name="WIDEmaxoutside">-</a>
  737. qh_WIDEmaxoutside
  738. Precision ratio for maximum increase for qh.max_outside in qh_check_maxout
  739. Precision errors while constructing the hull, may lead to very wide facets when checked in qh_check_maxout
  740. Nearly incident points in 4-d and higher is the most likely culprit
  741. Skip qh_check_maxout with 'Q5' (no-check-outer)
  742. Do not error with option 'Q12' (allow-wide)
  743. Do not warn with options 'Q12 Pp'
  744. */
  745. #define qh_WIDEmaxoutside 100
  746. /*-<a href="qh-user_r.htm#TOC"
  747. >--------------------------------</a><a name="WIDEmaxoutside2">-</a>
  748. qh_WIDEmaxoutside2
  749. Precision ratio for maximum qh.max_outside in qh_check_maxout
  750. Skip qh_check_maxout with 'Q5' no-check-outer
  751. Do not error with option 'Q12' allow-wide
  752. */
  753. #define qh_WIDEmaxoutside2 (10*qh_WIDEmaxoutside)
  754. /*-<a href="qh-user_r.htm#TOC"
  755. >--------------------------------</a><a name="WIDEpinched">-</a>
  756. qh_WIDEpinched
  757. Merge ratio for distance between pinched vertices compared to current facet width for qh_getpinchedmerges and qh_next_vertexmerge
  758. Reports warning and merges duplicate ridges instead
  759. Enable these attempts with option Q14 merge-pinched-vertices
  760. notes:
  761. Merging pinched vertices should prevent duplicate ridges (see qh_WIDEduplicate)
  762. Merging the duplicate ridges may be better than merging the pinched vertices
  763. 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
  764. Actual distance to facets is a third to a tenth of the qh_pointdist (T1)
  765. */
  766. #define qh_WIDEpinched 100
  767. /*-<a href="qh-user_r.htm#TOC"
  768. >--------------------------------</a><a name="ZEROdelaunay">-</a>
  769. qh_ZEROdelaunay
  770. a zero Delaunay facet occurs for input sites coplanar with their convex hull
  771. the last normal coefficient of a zero Delaunay facet is within
  772. qh_ZEROdelaunay * qh.ANGLEround of 0
  773. notes:
  774. qh_ZEROdelaunay does not allow for joggled input ('QJ').
  775. You can avoid zero Delaunay facets by surrounding the input with a box.
  776. Use option 'PDk:-n' to explicitly define zero Delaunay facets
  777. k= dimension of input sites (e.g., 3 for 3-d Delaunay triangulation)
  778. n= the cutoff for zero Delaunay facets (e.g., 'PD3:-1e-12')
  779. */
  780. #define qh_ZEROdelaunay 2
  781. /*============================================================*/
  782. /*============= Microsoft DevStudio ==========================*/
  783. /*============================================================*/
  784. /*
  785. Finding Memory Leaks Using the CRT Library
  786. https://msdn.microsoft.com/en-us/library/x98tx3cf(v=vs.100).aspx
  787. Reports enabled in qh_lib_check for Debug window and stderr
  788. From 2005=>msvcr80d, 2010=>msvcr100d, 2012=>msvcr110d
  789. Watch: {,,msvcr80d.dll}_crtBreakAlloc Value from {n} in the leak report
  790. _CrtSetBreakAlloc(689); // qh_lib_check() [global_r.c]
  791. Examples
  792. http://free-cad.sourceforge.net/SrcDocu/d2/d7f/MemDebug_8cpp_source.html
  793. https://github.com/illlust/Game/blob/master/library/MemoryLeak.cpp
  794. */
  795. #if 0 /* off (0) by default for QHULL_CRTDBG */
  796. #define QHULL_CRTDBG
  797. #endif
  798. #if defined(_MSC_VER) && defined(_DEBUG) && defined(QHULL_CRTDBG)
  799. #define _CRTDBG_MAP_ALLOC
  800. #include <stdlib.h>
  801. #include <crtdbg.h>
  802. #endif
  803. #endif /* qh_DEFuser */