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


hpp-fcl
Author(s):
autogenerated on Fri Jan 26 2024 03:46:15