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