1 |
/*<html><pre> -<a href="qh-globa.htm" |
2 |
>-------------------------------</a><a name="TOP">-</a> |
3 |
|
4 |
global.c |
5 |
initializes all the globals of the qhull application |
6 |
|
7 |
see README |
8 |
|
9 |
see qhull.h for qh.globals and function prototypes |
10 |
|
11 |
see qhull_a.h for internal functions |
12 |
|
13 |
copyright (c) 1993-2003, The Geometry Center |
14 |
*/ |
15 |
|
16 |
#include "QuickHull/qhull_a.h" |
17 |
|
18 |
/*========= qh definition (see qhull.h) =======================*/ |
19 |
|
20 |
#if qh_QHpointer |
21 |
qhT *qh_qh= NULL; /* pointer to all global variables */ |
22 |
#else |
23 |
qhT qh_qh; /* all global variables. |
24 |
Add "= {0}" if this causes a compiler error. |
25 |
Also qh_qhstat in stat.c and qhmem in mem.c. */ |
26 |
#endif |
27 |
|
28 |
/*-<a href ="qh-globa.htm#TOC" |
29 |
>--------------------------------</a><a name="qh_version">-</a> |
30 |
|
31 |
qh_version |
32 |
version string by year and date |
33 |
|
34 |
the revision increases on code changes only |
35 |
|
36 |
notes: |
37 |
change date: Changes.txt, Announce.txt, README.txt, |
38 |
qhull.man, qhull.txt, qhull-news.html, Eudora signatures, |
39 |
change version: README.txt, qh-get.htm, File_id.diz, Makefile.txt |
40 |
change year: Copying.txt |
41 |
check download size |
42 |
recompile user_eg.c, rbox.c, qhull.c, qconvex.c, qdelaun.c qvoronoi.c, qhalf.c |
43 |
*/ |
44 |
|
45 |
char *qh_version = "2003.1 2003/12/30"; |
46 |
|
47 |
/*-<a href="qh-globa.htm#TOC" |
48 |
>-------------------------------</a><a name="appendprint">-</a> |
49 |
|
50 |
qh_appendprint( printFormat ) |
51 |
append printFormat to qh.PRINTout unless already defined |
52 |
*/ |
53 |
void qh_appendprint (qh_PRINT format) { |
54 |
int i; |
55 |
|
56 |
for (i=0; i < qh_PRINTEND; i++) { |
57 |
if (qh PRINTout[i] == format && format != qh_PRINTqhull) |
58 |
break; |
59 |
if (!qh PRINTout[i]) { |
60 |
qh PRINTout[i]= format; |
61 |
break; |
62 |
} |
63 |
} |
64 |
} /* appendprint */ |
65 |
|
66 |
/*-<a href="qh-globa.htm#TOC" |
67 |
>-------------------------------</a><a name="checkflags">-</a> |
68 |
|
69 |
qh_checkflags( commandStr, hiddenFlags ) |
70 |
errors if commandStr contains hiddenFlags |
71 |
hiddenFlags starts and ends with a space and is space deliminated (checked) |
72 |
|
73 |
notes: |
74 |
ignores first word (e.g., "qconvex i") |
75 |
please use qh_strtol/strtod since strtol/strtod may or may not skip trailing spaces |
76 |
|
77 |
see: |
78 |
qh_initflags() initializes Qhull according to commandStr |
79 |
*/ |
80 |
void qh_checkflags(char *command, char *hiddenflags) { |
81 |
char *s= command, *t, *chkerr, key, opt, prevopt; |
82 |
char chkkey[]= " "; |
83 |
char chkopt[]= " "; |
84 |
char chkopt2[]= " "; |
85 |
|
86 |
if (*hiddenflags != ' ' || hiddenflags[strlen(hiddenflags)-1] != ' ') { |
87 |
fprintf(qh ferr, "qhull error (qh_checkflags): hiddenflags must start and end with a space: \"%s\"", hiddenflags); |
88 |
qh_errexit(qh_ERRinput, NULL, NULL); |
89 |
} |
90 |
if (strpbrk(hiddenflags, ",\n\r\t")) { |
91 |
fprintf(qh ferr, "qhull error (qh_checkflags): hiddenflags contains commas, newlines, or tabs: \"%s\"", hiddenflags); |
92 |
qh_errexit(qh_ERRinput, NULL, NULL); |
93 |
} |
94 |
while (*s && !isspace(*s)) /* skip program name */ |
95 |
s++; |
96 |
while (*s) { |
97 |
while (*s && isspace(*s)) |
98 |
s++; |
99 |
if (*s == '-') |
100 |
s++; |
101 |
if (!*s) |
102 |
break; |
103 |
key = *s++; |
104 |
chkerr = NULL; |
105 |
if (key == '\'') { /* TO 'file name' */ |
106 |
t= strchr(s, '\''); |
107 |
if (!t) { |
108 |
fprintf(qh ferr, "qhull error (qh_checkflags): missing the 2nd single-quote for:\n%s\n", s-1); |
109 |
qh_errexit(qh_ERRinput, NULL, NULL); |
110 |
} |
111 |
s= t+1; |
112 |
continue; |
113 |
} |
114 |
chkkey[1]= key; |
115 |
if (strstr(hiddenflags, chkkey)) { |
116 |
chkerr= chkkey; |
117 |
}else if (isupper(key)) { |
118 |
opt= ' '; |
119 |
prevopt= ' '; |
120 |
chkopt[1]= key; |
121 |
chkopt2[1]= key; |
122 |
while (!chkerr && *s && !isspace(*s)) { |
123 |
opt= *s++; |
124 |
if (isalpha(opt)) { |
125 |
chkopt[2]= opt; |
126 |
if (strstr(hiddenflags, chkopt)) |
127 |
chkerr= chkopt; |
128 |
if (prevopt != ' ') { |
129 |
chkopt2[2]= prevopt; |
130 |
chkopt2[3]= opt; |
131 |
if (strstr(hiddenflags, chkopt2)) |
132 |
chkerr= chkopt2; |
133 |
} |
134 |
}else if (key == 'Q' && isdigit(opt) && prevopt != 'b' |
135 |
&& (prevopt == ' ' || islower(prevopt))) { |
136 |
chkopt[2]= opt; |
137 |
if (strstr(hiddenflags, chkopt)) |
138 |
chkerr= chkopt; |
139 |
}else { |
140 |
qh_strtod (s-1, &t); |
141 |
if (s < t) |
142 |
s= t; |
143 |
} |
144 |
prevopt= opt; |
145 |
} |
146 |
} |
147 |
if (chkerr) { |
148 |
*chkerr= '\''; |
149 |
chkerr[strlen(chkerr)-1]= '\''; |
150 |
fprintf(qh ferr, "qhull error: option %s is not used with this program.\n It may be used with qhull.\n", chkerr); |
151 |
qh_errexit(qh_ERRinput, NULL, NULL); |
152 |
} |
153 |
} |
154 |
} /* checkflags */ |
155 |
|
156 |
/*-<a href="qh-globa.htm#TOC" |
157 |
>-------------------------------</a><a name="clock">-</a> |
158 |
|
159 |
qh_clock() |
160 |
return user CPU time in 100ths (qh_SECtick) |
161 |
only defined for qh_CLOCKtype == 2 |
162 |
|
163 |
notes: |
164 |
please use first value to determine time 0 |
165 |
from Stevens '92 8.15 |
166 |
*/ |
167 |
unsigned long qh_clock (void) { |
168 |
|
169 |
#if (qh_CLOCKtype == 2) |
170 |
struct tms time; |
171 |
static long clktck; /* initialized first call */ |
172 |
double ratio, cpu; |
173 |
unsigned long ticks; |
174 |
|
175 |
if (!clktck) { |
176 |
if ((clktck= sysconf (_SC_CLK_TCK)) < 0) { |
177 |
fprintf (qh ferr, "qhull internal error (qh_clock): sysconf() failed. Use qh_CLOCKtype 1 in user.h\n"); |
178 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
179 |
} |
180 |
} |
181 |
if (times (&time) == -1) { |
182 |
fprintf (qh ferr, "qhull internal error (qh_clock): times() failed. Use qh_CLOCKtype 1 in user.h\n"); |
183 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
184 |
} |
185 |
ratio= qh_SECticks / (double)clktck; |
186 |
ticks= time.tms_utime * ratio; |
187 |
return ticks; |
188 |
#else |
189 |
fprintf (qh ferr, "qhull internal error (qh_clock): use qh_CLOCKtype 2 in user.h\n"); |
190 |
qh_errexit (qh_ERRqhull, NULL, NULL); /* never returns */ |
191 |
return 0; |
192 |
#endif |
193 |
} /* clock */ |
194 |
|
195 |
/*-<a href="qh-globa.htm#TOC" |
196 |
>-------------------------------</a><a name="freebuffers">-</a> |
197 |
|
198 |
qh_freebuffers() |
199 |
free up global memory buffers |
200 |
|
201 |
notes: |
202 |
must match qh_initbuffers() |
203 |
*/ |
204 |
void qh_freebuffers (void) { |
205 |
|
206 |
trace5((qh ferr, "qh_freebuffers: freeing up global memory buffers\n")); |
207 |
/* allocated by qh_initqhull_buffers */ |
208 |
qh_memfree (qh NEARzero, qh hull_dim * sizeof(realT)); |
209 |
qh_memfree (qh lower_threshold, (qh input_dim+1) * sizeof(realT)); |
210 |
qh_memfree (qh upper_threshold, (qh input_dim+1) * sizeof(realT)); |
211 |
qh_memfree (qh lower_bound, (qh input_dim+1) * sizeof(realT)); |
212 |
qh_memfree (qh upper_bound, (qh input_dim+1) * sizeof(realT)); |
213 |
qh_memfree (qh gm_matrix, (qh hull_dim+1) * qh hull_dim * sizeof(coordT)); |
214 |
qh_memfree (qh gm_row, (qh hull_dim+1) * sizeof(coordT *)); |
215 |
qh NEARzero= qh lower_threshold= qh upper_threshold= NULL; |
216 |
qh lower_bound= qh upper_bound= NULL; |
217 |
qh gm_matrix= NULL; |
218 |
qh gm_row= NULL; |
219 |
qh_setfree (&qh other_points); |
220 |
qh_setfree (&qh del_vertices); |
221 |
qh_setfree (&qh coplanarset); |
222 |
if (qh line) /* allocated by qh_readinput, freed if no error */ |
223 |
free (qh line); |
224 |
if (qh half_space) |
225 |
free (qh half_space); |
226 |
if (qh temp_malloc) |
227 |
free (qh temp_malloc); |
228 |
if (qh feasible_point) /* allocated by qh_readfeasible */ |
229 |
free (qh feasible_point); |
230 |
if (qh feasible_string) /* allocated by qh_initflags */ |
231 |
free (qh feasible_string); |
232 |
qh line= qh feasible_string= NULL; |
233 |
qh half_space= qh feasible_point= qh temp_malloc= NULL; |
234 |
/* usually allocated by qh_readinput */ |
235 |
if (qh first_point && qh POINTSmalloc) { |
236 |
free(qh first_point); |
237 |
qh first_point= NULL; |
238 |
} |
239 |
if (qh input_points && qh input_malloc) { /* set by qh_joggleinput */ |
240 |
free (qh input_points); |
241 |
qh input_points= NULL; |
242 |
} |
243 |
trace5((qh ferr, "qh_freebuffers: finished\n")); |
244 |
} /* freebuffers */ |
245 |
|
246 |
|
247 |
/*-<a href="qh-globa.htm#TOC" |
248 |
>-------------------------------</a><a name="freebuild">-</a> |
249 |
|
250 |
qh_freebuild( allmem ) |
251 |
free global memory used by qh_initbuild and qh_buildhull |
252 |
if !allmem, |
253 |
does not free short memory (freed by qh_memfreeshort) |
254 |
|
255 |
design: |
256 |
free centrums |
257 |
free each vertex |
258 |
mark unattached ridges |
259 |
for each facet |
260 |
free ridges |
261 |
free outside set, coplanar set, neighbor set, ridge set, vertex set |
262 |
free facet |
263 |
free hash table |
264 |
free interior point |
265 |
free merge set |
266 |
free temporary sets |
267 |
*/ |
268 |
void qh_freebuild (boolT allmem) { |
269 |
facetT *facet; |
270 |
vertexT *vertex; |
271 |
ridgeT *ridge, **ridgep; |
272 |
mergeT *merge, **mergep; |
273 |
|
274 |
trace1((qh ferr, "qh_freebuild: free memory from qh_inithull and qh_buildhull\n")); |
275 |
if (qh del_vertices) |
276 |
qh_settruncate (qh del_vertices, 0); |
277 |
if (allmem) { |
278 |
qh_clearcenters (qh_ASnone); |
279 |
while ((vertex= qh vertex_list)) { |
280 |
if (vertex->next) |
281 |
qh_delvertex (vertex); |
282 |
else { |
283 |
qh_memfree (vertex, sizeof(vertexT)); |
284 |
qh newvertex_list= qh vertex_list= NULL; |
285 |
} |
286 |
} |
287 |
}else if (qh VERTEXneighbors) { |
288 |
FORALLvertices |
289 |
qh_setfreelong (&(vertex->neighbors)); |
290 |
} |
291 |
qh VERTEXneighbors= False; |
292 |
qh GOODclosest= NULL; |
293 |
if (allmem) { |
294 |
FORALLfacets { |
295 |
FOREACHridge_(facet->ridges) |
296 |
ridge->seen= False; |
297 |
} |
298 |
FORALLfacets { |
299 |
if (facet->visible) { |
300 |
FOREACHridge_(facet->ridges) { |
301 |
if (!otherfacet_(ridge, facet)->visible) |
302 |
ridge->seen= True; /* an unattached ridge */ |
303 |
} |
304 |
} |
305 |
} |
306 |
while ((facet= qh facet_list)) { |
307 |
FOREACHridge_(facet->ridges) { |
308 |
if (ridge->seen) { |
309 |
qh_setfree(&(ridge->vertices)); |
310 |
qh_memfree(ridge, sizeof(ridgeT)); |
311 |
}else |
312 |
ridge->seen= True; |
313 |
} |
314 |
qh_setfree (&(facet->outsideset)); |
315 |
qh_setfree (&(facet->coplanarset)); |
316 |
qh_setfree (&(facet->neighbors)); |
317 |
qh_setfree (&(facet->ridges)); |
318 |
qh_setfree (&(facet->vertices)); |
319 |
if (facet->next) |
320 |
qh_delfacet (facet); |
321 |
else { |
322 |
qh_memfree (facet, sizeof(facetT)); |
323 |
qh visible_list= qh newfacet_list= qh facet_list= NULL; |
324 |
} |
325 |
} |
326 |
}else { |
327 |
FORALLfacets { |
328 |
qh_setfreelong (&(facet->outsideset)); |
329 |
qh_setfreelong (&(facet->coplanarset)); |
330 |
if (!facet->simplicial) { |
331 |
qh_setfreelong (&(facet->neighbors)); |
332 |
qh_setfreelong (&(facet->ridges)); |
333 |
qh_setfreelong (&(facet->vertices)); |
334 |
} |
335 |
} |
336 |
} |
337 |
qh_setfree (&(qh hash_table)); |
338 |
qh_memfree (qh interior_point, qh normal_size); |
339 |
qh interior_point= NULL; |
340 |
FOREACHmerge_(qh facet_mergeset) /* usually empty */ |
341 |
qh_memfree (merge, sizeof(mergeT)); |
342 |
qh facet_mergeset= NULL; /* temp set */ |
343 |
qh degen_mergeset= NULL; /* temp set */ |
344 |
qh_settempfree_all(); |
345 |
} /* freebuild */ |
346 |
|
347 |
/*-<a href="qh-globa.htm#TOC" |
348 |
>-------------------------------</a><a name="freeqhull">-</a> |
349 |
|
350 |
qh_freeqhull( allmem ) |
351 |
free global memory |
352 |
if !allmem, |
353 |
does not free short memory (freed by qh_memfreeshort) |
354 |
|
355 |
notes: |
356 |
sets qh.NOerrexit in case caller forgets to |
357 |
|
358 |
design: |
359 |
free global and temporary memory from qh_initbuild and qh_buildhull |
360 |
free buffers |
361 |
free statistics |
362 |
*/ |
363 |
void qh_freeqhull (boolT allmem) { |
364 |
|
365 |
trace1((qh ferr, "qh_freeqhull: free global memory\n")); |
366 |
qh NOerrexit= True; /* no more setjmp since called at exit */ |
367 |
qh_freebuild (allmem); |
368 |
qh_freebuffers(); |
369 |
qh_freestatistics(); |
370 |
#if qh_QHpointer |
371 |
free (qh_qh); |
372 |
qh_qh= NULL; |
373 |
#else |
374 |
memset((char *)&qh_qh, 0, sizeof(qhT)); |
375 |
qh NOerrexit= True; |
376 |
#endif |
377 |
} /* freeqhull */ |
378 |
|
379 |
/*-<a href="qh-globa.htm#TOC" |
380 |
>-------------------------------</a><a name="init_A">-</a> |
381 |
|
382 |
qh_init_A( infile, outfile, errfile, argc, argv ) |
383 |
initialize memory and stdio files |
384 |
convert input options to option string (qh.qhull_command) |
385 |
|
386 |
notes: |
387 |
infile may be NULL if qh_readpoints() is not called |
388 |
|
389 |
errfile should always be defined. It is used for reporting |
390 |
errors. outfile is used for output and format options. |
391 |
|
392 |
argc/argv may be 0/NULL |
393 |
|
394 |
called before error handling initialized |
395 |
qh_errexit() may not be used |
396 |
*/ |
397 |
void qh_init_A (FILE *infile, FILE *outfile, FILE *errfile, int argc, char *argv[]) { |
398 |
qh_meminit (errfile); |
399 |
qh_initqhull_start (infile, outfile, errfile); |
400 |
qh_init_qhull_command (argc, argv); |
401 |
} /* init_A */ |
402 |
|
403 |
/*-<a href="qh-globa.htm#TOC" |
404 |
>-------------------------------</a><a name="init_B">-</a> |
405 |
|
406 |
qh_init_B( points, numpoints, dim, ismalloc ) |
407 |
initialize globals for points array |
408 |
|
409 |
points has numpoints dim-dimensional points |
410 |
points[0] is the first coordinate of the first point |
411 |
points[1] is the second coordinate of the first point |
412 |
points[dim] is the first coordinate of the second point |
413 |
|
414 |
ismalloc=True |
415 |
Qhull will call free(points) on exit or input transformation |
416 |
ismalloc=False |
417 |
Qhull will allocate a new point array if needed for input transformation |
418 |
|
419 |
qh.qhull_command |
420 |
is the option string. |
421 |
It is defined by qh_init_B(), qh_qhull_command(), or qh_initflags |
422 |
|
423 |
returns: |
424 |
if qh.PROJECTinput or (qh.DELAUNAY and qh.PROJECTdelaunay) |
425 |
projects the input to a new point array |
426 |
|
427 |
if qh.DELAUNAY, |
428 |
qh.hull_dim is increased by one |
429 |
if qh.ATinfinity, |
430 |
qh_projectinput adds point-at-infinity for Delaunay tri. |
431 |
|
432 |
if qh.SCALEinput |
433 |
changes the upper and lower bounds of the input, see qh_scaleinput() |
434 |
|
435 |
if qh.ROTATEinput |
436 |
rotates the input by a random rotation, see qh_rotateinput() |
437 |
if qh.DELAUNAY |
438 |
rotates about the last coordinate |
439 |
|
440 |
notes: |
441 |
called after points are defined |
442 |
qh_errexit() may be used |
443 |
*/ |
444 |
void qh_init_B (coordT *points, int numpoints, int dim, boolT ismalloc) { |
445 |
qh_initqhull_globals (points, numpoints, dim, ismalloc); |
446 |
if (qhmem.LASTsize == 0) |
447 |
qh_initqhull_mem(); |
448 |
/* mem.c and qset.c are initialized */ |
449 |
qh_initqhull_buffers(); |
450 |
qh_initthresholds (qh qhull_command); |
451 |
if (qh PROJECTinput || (qh DELAUNAY && qh PROJECTdelaunay)) |
452 |
qh_projectinput(); |
453 |
if (qh SCALEinput) |
454 |
qh_scaleinput(); |
455 |
if (qh ROTATErandom >= 0) { |
456 |
qh_randommatrix (qh gm_matrix, qh hull_dim, qh gm_row); |
457 |
if (qh DELAUNAY) { |
458 |
int k, lastk= qh hull_dim-1; |
459 |
for (k= 0; k < lastk; k++) { |
460 |
qh gm_row[k][lastk]= 0.0; |
461 |
qh gm_row[lastk][k]= 0.0; |
462 |
} |
463 |
qh gm_row[lastk][lastk]= 1.0; |
464 |
} |
465 |
qh_gram_schmidt (qh hull_dim, qh gm_row); |
466 |
qh_rotateinput (qh gm_row); |
467 |
} |
468 |
} /* init_B */ |
469 |
|
470 |
/*-<a href="qh-globa.htm#TOC" |
471 |
>-------------------------------</a><a name="init_qhull_command">-</a> |
472 |
|
473 |
qh_init_qhull_command( argc, argv ) |
474 |
build qh.qhull_command from argc/argv |
475 |
|
476 |
returns: |
477 |
a space-deliminated string of options (just as typed) |
478 |
|
479 |
notes: |
480 |
makes option string easy to input and output |
481 |
|
482 |
argc/argv may be 0/NULL |
483 |
*/ |
484 |
void qh_init_qhull_command(int argc, char *argv[]) { |
485 |
int i; |
486 |
char *s; |
487 |
|
488 |
if (argc) { |
489 |
if ((s= strrchr( argv[0], '\\'))) /* Borland gives full path */ |
490 |
strcpy (qh qhull_command, s+1); |
491 |
else |
492 |
strcpy (qh qhull_command, argv[0]); |
493 |
if ((s= strstr (qh qhull_command, ".EXE")) |
494 |
|| (s= strstr (qh qhull_command, ".exe"))) |
495 |
*s= '\0'; |
496 |
} |
497 |
for (i=1; i < argc; i++) { |
498 |
if (strlen (qh qhull_command) + strlen(argv[i]) + 1 < sizeof(qh qhull_command)) { |
499 |
strcat (qh qhull_command, " "); |
500 |
strcat (qh qhull_command, argv[i]); |
501 |
}else { |
502 |
fprintf (qh ferr, "qhull input error: more than %d characters in command line\n", |
503 |
(int)sizeof(qh qhull_command)); |
504 |
exit (1); /* can not use qh_errexit */ |
505 |
} |
506 |
} |
507 |
} /* init_qhull_command */ |
508 |
|
509 |
/*-<a href="qh-globa.htm#TOC" |
510 |
>-------------------------------</a><a name="initflags">-</a> |
511 |
|
512 |
qh_initflags( commandStr ) |
513 |
set flags and initialized constants from commandStr |
514 |
|
515 |
returns: |
516 |
sets qh.qhull_command to command if needed |
517 |
|
518 |
notes: |
519 |
ignores first word (e.g., "qhull d") |
520 |
please use qh_strtol/strtod since strtol/strtod may or may not skip trailing spaces |
521 |
|
522 |
see: |
523 |
qh_initthresholds() continues processing of 'Pdn' and 'PDn' |
524 |
'prompt' in unix.c for documentation |
525 |
|
526 |
design: |
527 |
for each space-deliminated option group |
528 |
if top-level option |
529 |
check syntax |
530 |
append approriate option to option string |
531 |
set appropriate global variable or append printFormat to print options |
532 |
else |
533 |
for each sub-option |
534 |
check syntax |
535 |
append approriate option to option string |
536 |
set appropriate global variable or append printFormat to print options |
537 |
|
538 |
|
539 |
*/ |
540 |
void qh_initflags(char *command) { |
541 |
int k, i, lastproject; |
542 |
char *s= command, *t, *prev_s, *start, key; |
543 |
boolT isgeom= False, wasproject; |
544 |
realT r; |
545 |
|
546 |
if (command != &qh qhull_command[0]) { |
547 |
*qh qhull_command= '\0'; |
548 |
strncat( qh qhull_command, command, sizeof( qh qhull_command)); |
549 |
} |
550 |
while (*s && !isspace(*s)) /* skip program name */ |
551 |
s++; |
552 |
while (*s) { |
553 |
while (*s && isspace(*s)) |
554 |
s++; |
555 |
if (*s == '-') |
556 |
s++; |
557 |
if (!*s) |
558 |
break; |
559 |
prev_s= s; |
560 |
switch (*s++) { |
561 |
case 'd': |
562 |
qh_option ("delaunay", NULL, NULL); |
563 |
qh DELAUNAY= True; |
564 |
break; |
565 |
case 'f': |
566 |
qh_option ("facets", NULL, NULL); |
567 |
qh_appendprint (qh_PRINTfacets); |
568 |
break; |
569 |
case 'i': |
570 |
qh_option ("incidence", NULL, NULL); |
571 |
qh_appendprint (qh_PRINTincidences); |
572 |
break; |
573 |
case 'm': |
574 |
qh_option ("mathematica", NULL, NULL); |
575 |
qh_appendprint (qh_PRINTmathematica); |
576 |
break; |
577 |
case 'n': |
578 |
qh_option ("normals", NULL, NULL); |
579 |
qh_appendprint (qh_PRINTnormals); |
580 |
break; |
581 |
case 'o': |
582 |
qh_option ("offFile", NULL, NULL); |
583 |
qh_appendprint (qh_PRINToff); |
584 |
break; |
585 |
case 'p': |
586 |
qh_option ("points", NULL, NULL); |
587 |
qh_appendprint (qh_PRINTpoints); |
588 |
break; |
589 |
case 's': |
590 |
qh_option ("summary", NULL, NULL); |
591 |
qh PRINTsummary= True; |
592 |
break; |
593 |
case 'v': |
594 |
qh_option ("voronoi", NULL, NULL); |
595 |
qh VORONOI= True; |
596 |
qh DELAUNAY= True; |
597 |
break; |
598 |
case 'A': |
599 |
if (!isdigit(*s) && *s != '.' && *s != '-') |
600 |
fprintf(qh ferr, "qhull warning: no maximum cosine angle given for option 'An'. Ignored.\n"); |
601 |
else { |
602 |
if (*s == '-') { |
603 |
qh premerge_cos= -qh_strtod (s, &s); |
604 |
qh_option ("Angle-premerge-", NULL, &qh premerge_cos); |
605 |
qh PREmerge= True; |
606 |
}else { |
607 |
qh postmerge_cos= qh_strtod (s, &s); |
608 |
qh_option ("Angle-postmerge", NULL, &qh postmerge_cos); |
609 |
qh POSTmerge= True; |
610 |
} |
611 |
qh MERGING= True; |
612 |
} |
613 |
break; |
614 |
case 'C': |
615 |
if (!isdigit(*s) && *s != '.' && *s != '-') |
616 |
fprintf(qh ferr, "qhull warning: no centrum radius given for option 'Cn'. Ignored.\n"); |
617 |
else { |
618 |
if (*s == '-') { |
619 |
qh premerge_centrum= -qh_strtod (s, &s); |
620 |
qh_option ("Centrum-premerge-", NULL, &qh premerge_centrum); |
621 |
qh PREmerge= True; |
622 |
}else { |
623 |
qh postmerge_centrum= qh_strtod (s, &s); |
624 |
qh_option ("Centrum-postmerge", NULL, &qh postmerge_centrum); |
625 |
qh POSTmerge= True; |
626 |
} |
627 |
qh MERGING= True; |
628 |
} |
629 |
break; |
630 |
case 'E': |
631 |
if (*s == '-') |
632 |
fprintf(qh ferr, "qhull warning: negative maximum roundoff given for option 'An'. Ignored.\n"); |
633 |
else if (!isdigit(*s)) |
634 |
fprintf(qh ferr, "qhull warning: no maximum roundoff given for option 'En'. Ignored.\n"); |
635 |
else { |
636 |
qh DISTround= qh_strtod (s, &s); |
637 |
qh_option ("Distance-roundoff", NULL, &qh DISTround); |
638 |
qh SETroundoff= True; |
639 |
} |
640 |
break; |
641 |
case 'H': |
642 |
start= s; |
643 |
qh HALFspace= True; |
644 |
qh_strtod (s, &t); |
645 |
while (t > s) { |
646 |
if (*t && !isspace (*t)) { |
647 |
if (*t == ',') |
648 |
t++; |
649 |
else |
650 |
fprintf (qh ferr, "qhull warning: origin for Halfspace intersection should be 'Hn,n,n,...'\n"); |
651 |
} |
652 |
s= t; |
653 |
qh_strtod (s, &t); |
654 |
} |
655 |
if (start < t) { |
656 |
if (!(qh feasible_string= (char*)calloc (t-start+1, 1))) { |
657 |
fprintf(qh ferr, "qhull error: insufficient memory for 'Hn,n,n'\n"); |
658 |
qh_errexit(qh_ERRmem, NULL, NULL); |
659 |
} |
660 |
strncpy (qh feasible_string, start, t-start); |
661 |
qh_option ("Halfspace-about", NULL, NULL); |
662 |
qh_option (qh feasible_string, NULL, NULL); |
663 |
}else |
664 |
qh_option ("Halfspace", NULL, NULL); |
665 |
break; |
666 |
case 'R': |
667 |
if (!isdigit(*s)) |
668 |
fprintf(qh ferr, "qhull warning: missing random perturbation for option 'Rn'. Ignored\n"); |
669 |
else { |
670 |
qh RANDOMfactor= qh_strtod (s, &s); |
671 |
qh_option ("Random_perturb", NULL, &qh RANDOMfactor); |
672 |
qh RANDOMdist= True; |
673 |
} |
674 |
break; |
675 |
case 'V': |
676 |
if (!isdigit(*s) && *s != '-') |
677 |
fprintf(qh ferr, "qhull warning: missing visible distance for option 'Vn'. Ignored\n"); |
678 |
else { |
679 |
qh MINvisible= qh_strtod (s, &s); |
680 |
qh_option ("Visible", NULL, &qh MINvisible); |
681 |
} |
682 |
break; |
683 |
case 'U': |
684 |
if (!isdigit(*s) && *s != '-') |
685 |
fprintf(qh ferr, "qhull warning: missing coplanar distance for option 'Un'. Ignored\n"); |
686 |
else { |
687 |
qh MAXcoplanar= qh_strtod (s, &s); |
688 |
qh_option ("U-coplanar", NULL, &qh MAXcoplanar); |
689 |
} |
690 |
break; |
691 |
case 'W': |
692 |
if (*s == '-') |
693 |
fprintf(qh ferr, "qhull warning: negative outside width for option 'Wn'. Ignored.\n"); |
694 |
else if (!isdigit(*s)) |
695 |
fprintf(qh ferr, "qhull warning: missing outside width for option 'Wn'. Ignored\n"); |
696 |
else { |
697 |
qh MINoutside= qh_strtod (s, &s); |
698 |
qh_option ("W-outside", NULL, &qh MINoutside); |
699 |
qh APPROXhull= True; |
700 |
} |
701 |
break; |
702 |
/************ sub menus ***************/ |
703 |
case 'F': |
704 |
while (*s && !isspace(*s)) { |
705 |
switch(*s++) { |
706 |
case 'a': |
707 |
qh_option ("Farea", NULL, NULL); |
708 |
qh_appendprint (qh_PRINTarea); |
709 |
qh GETarea= True; |
710 |
break; |
711 |
case 'A': |
712 |
qh_option ("FArea-total", NULL, NULL); |
713 |
qh GETarea= True; |
714 |
break; |
715 |
case 'c': |
716 |
qh_option ("Fcoplanars", NULL, NULL); |
717 |
qh_appendprint (qh_PRINTcoplanars); |
718 |
break; |
719 |
case 'C': |
720 |
qh_option ("FCentrums", NULL, NULL); |
721 |
qh_appendprint (qh_PRINTcentrums); |
722 |
break; |
723 |
case 'd': |
724 |
qh_option ("Fd-cdd-in", NULL, NULL); |
725 |
qh CDDinput= True; |
726 |
break; |
727 |
case 'D': |
728 |
qh_option ("FD-cdd-out", NULL, NULL); |
729 |
qh CDDoutput= True; |
730 |
break; |
731 |
case 'F': |
732 |
qh_option ("FFacets-xridge", NULL, NULL); |
733 |
qh_appendprint (qh_PRINTfacets_xridge); |
734 |
break; |
735 |
case 'i': |
736 |
qh_option ("Finner", NULL, NULL); |
737 |
qh_appendprint (qh_PRINTinner); |
738 |
break; |
739 |
case 'I': |
740 |
qh_option ("FIDs", NULL, NULL); |
741 |
qh_appendprint (qh_PRINTids); |
742 |
break; |
743 |
case 'm': |
744 |
qh_option ("Fmerges", NULL, NULL); |
745 |
qh_appendprint (qh_PRINTmerges); |
746 |
break; |
747 |
case 'M': |
748 |
qh_option ("FMaple", NULL, NULL); |
749 |
qh_appendprint (qh_PRINTmaple); |
750 |
break; |
751 |
case 'n': |
752 |
qh_option ("Fneighbors", NULL, NULL); |
753 |
qh_appendprint (qh_PRINTneighbors); |
754 |
break; |
755 |
case 'N': |
756 |
qh_option ("FNeighbors-vertex", NULL, NULL); |
757 |
qh_appendprint (qh_PRINTvneighbors); |
758 |
break; |
759 |
case 'o': |
760 |
qh_option ("Fouter", NULL, NULL); |
761 |
qh_appendprint (qh_PRINTouter); |
762 |
break; |
763 |
case 'O': |
764 |
if (qh PRINToptions1st) { |
765 |
qh_option ("FOptions", NULL, NULL); |
766 |
qh_appendprint (qh_PRINToptions); |
767 |
}else |
768 |
qh PRINToptions1st= True; |
769 |
break; |
770 |
case 'p': |
771 |
qh_option ("Fpoint-intersect", NULL, NULL); |
772 |
qh_appendprint (qh_PRINTpointintersect); |
773 |
break; |
774 |
case 'P': |
775 |
qh_option ("FPoint-nearest", NULL, NULL); |
776 |
qh_appendprint (qh_PRINTpointnearest); |
777 |
break; |
778 |
case 'Q': |
779 |
qh_option ("FQhull", NULL, NULL); |
780 |
qh_appendprint (qh_PRINTqhull); |
781 |
break; |
782 |
case 's': |
783 |
qh_option ("Fsummary", NULL, NULL); |
784 |
qh_appendprint (qh_PRINTsummary); |
785 |
break; |
786 |
case 'S': |
787 |
qh_option ("FSize", NULL, NULL); |
788 |
qh_appendprint (qh_PRINTsize); |
789 |
qh GETarea= True; |
790 |
break; |
791 |
case 't': |
792 |
qh_option ("Ftriangles", NULL, NULL); |
793 |
qh_appendprint (qh_PRINTtriangles); |
794 |
break; |
795 |
case 'v': |
796 |
/* option set in qh_initqhull_globals */ |
797 |
qh_appendprint (qh_PRINTvertices); |
798 |
break; |
799 |
case 'V': |
800 |
qh_option ("FVertex-average", NULL, NULL); |
801 |
qh_appendprint (qh_PRINTaverage); |
802 |
break; |
803 |
case 'x': |
804 |
qh_option ("Fxtremes", NULL, NULL); |
805 |
qh_appendprint (qh_PRINTextremes); |
806 |
break; |
807 |
default: |
808 |
s--; |
809 |
fprintf (qh ferr, "qhull warning: unknown 'F' output option %c, rest ignored\n", (int)s[0]); |
810 |
while (*++s && !isspace(*s)); |
811 |
break; |
812 |
} |
813 |
} |
814 |
break; |
815 |
case 'G': |
816 |
isgeom= True; |
817 |
qh_appendprint (qh_PRINTgeom); |
818 |
while (*s && !isspace(*s)) { |
819 |
switch(*s++) { |
820 |
case 'a': |
821 |
qh_option ("Gall-points", NULL, NULL); |
822 |
qh PRINTdots= True; |
823 |
break; |
824 |
case 'c': |
825 |
qh_option ("Gcentrums", NULL, NULL); |
826 |
qh PRINTcentrums= True; |
827 |
break; |
828 |
case 'h': |
829 |
qh_option ("Gintersections", NULL, NULL); |
830 |
qh DOintersections= True; |
831 |
break; |
832 |
case 'i': |
833 |
qh_option ("Ginner", NULL, NULL); |
834 |
qh PRINTinner= True; |
835 |
break; |
836 |
case 'n': |
837 |
qh_option ("Gno-planes", NULL, NULL); |
838 |
qh PRINTnoplanes= True; |
839 |
break; |
840 |
case 'o': |
841 |
qh_option ("Gouter", NULL, NULL); |
842 |
qh PRINTouter= True; |
843 |
break; |
844 |
case 'p': |
845 |
qh_option ("Gpoints", NULL, NULL); |
846 |
qh PRINTcoplanar= True; |
847 |
break; |
848 |
case 'r': |
849 |
qh_option ("Gridges", NULL, NULL); |
850 |
qh PRINTridges= True; |
851 |
break; |
852 |
case 't': |
853 |
qh_option ("Gtransparent", NULL, NULL); |
854 |
qh PRINTtransparent= True; |
855 |
break; |
856 |
case 'v': |
857 |
qh_option ("Gvertices", NULL, NULL); |
858 |
qh PRINTspheres= True; |
859 |
break; |
860 |
case 'D': |
861 |
if (!isdigit (*s)) |
862 |
fprintf (qh ferr, "qhull input error: missing dimension for option 'GDn'\n"); |
863 |
else { |
864 |
if (qh DROPdim >= 0) |
865 |
fprintf (qh ferr, "qhull warning: can only drop one dimension. Previous 'GD%d' ignored\n", |
866 |
qh DROPdim); |
867 |
qh DROPdim= qh_strtol (s, &s); |
868 |
qh_option ("GDrop-dim", &qh DROPdim, NULL); |
869 |
} |
870 |
break; |
871 |
default: |
872 |
s--; |
873 |
fprintf (qh ferr, "qhull warning: unknown 'G' print option %c, rest ignored\n", (int)s[0]); |
874 |
while (*++s && !isspace(*s)); |
875 |
break; |
876 |
} |
877 |
} |
878 |
break; |
879 |
case 'P': |
880 |
while (*s && !isspace(*s)) { |
881 |
switch(*s++) { |
882 |
case 'd': case 'D': /* see qh_initthresholds() */ |
883 |
key= s[-1]; |
884 |
i= qh_strtol (s, &s); |
885 |
r= 0; |
886 |
if (*s == ':') { |
887 |
s++; |
888 |
r= qh_strtod (s, &s); |
889 |
} |
890 |
if (key == 'd') |
891 |
qh_option ("Pdrop-facets-dim-less", &i, &r); |
892 |
else |
893 |
qh_option ("PDrop-facets-dim-more", &i, &r); |
894 |
break; |
895 |
case 'g': |
896 |
qh_option ("Pgood-facets", NULL, NULL); |
897 |
qh PRINTgood= True; |
898 |
break; |
899 |
case 'G': |
900 |
qh_option ("PGood-facet-neighbors", NULL, NULL); |
901 |
qh PRINTneighbors= True; |
902 |
break; |
903 |
case 'o': |
904 |
qh_option ("Poutput-forced", NULL, NULL); |
905 |
qh FORCEoutput= True; |
906 |
break; |
907 |
case 'p': |
908 |
qh_option ("Pprecision-ignore", NULL, NULL); |
909 |
qh PRINTprecision= False; |
910 |
break; |
911 |
case 'A': |
912 |
if (!isdigit (*s)) |
913 |
fprintf (qh ferr, "qhull input error: missing facet count for keep area option 'PAn'\n"); |
914 |
else { |
915 |
qh KEEParea= qh_strtol (s, &s); |
916 |
qh_option ("PArea-keep", &qh KEEParea, NULL); |
917 |
qh GETarea= True; |
918 |
} |
919 |
break; |
920 |
case 'F': |
921 |
if (!isdigit (*s)) |
922 |
fprintf (qh ferr, "qhull input error: missing facet area for option 'PFn'\n"); |
923 |
else { |
924 |
qh KEEPminArea= qh_strtod (s, &s); |
925 |
qh_option ("PFacet-area-keep", NULL, &qh KEEPminArea); |
926 |
qh GETarea= True; |
927 |
} |
928 |
break; |
929 |
case 'M': |
930 |
if (!isdigit (*s)) |
931 |
fprintf (qh ferr, "qhull input error: missing merge count for option 'PMn'\n"); |
932 |
else { |
933 |
qh KEEPmerge= qh_strtol (s, &s); |
934 |
qh_option ("PMerge-keep", &qh KEEPmerge, NULL); |
935 |
} |
936 |
break; |
937 |
default: |
938 |
s--; |
939 |
fprintf (qh ferr, "qhull warning: unknown 'P' print option %c, rest ignored\n", (int)s[0]); |
940 |
while (*++s && !isspace(*s)); |
941 |
break; |
942 |
} |
943 |
} |
944 |
break; |
945 |
case 'Q': |
946 |
lastproject= -1; |
947 |
while (*s && !isspace(*s)) { |
948 |
switch(*s++) { |
949 |
case 'b': case 'B': /* handled by qh_initthresholds */ |
950 |
key= s[-1]; |
951 |
if (key == 'b' && *s == 'B') { |
952 |
s++; |
953 |
r= qh_DEFAULTbox; |
954 |
qh SCALEinput= True; |
955 |
qh_option ("QbBound-unit-box", NULL, &r); |
956 |
break; |
957 |
} |
958 |
if (key == 'b' && *s == 'b') { |
959 |
s++; |
960 |
qh SCALElast= True; |
961 |
qh_option ("Qbbound-last", NULL, NULL); |
962 |
break; |
963 |
} |
964 |
k= qh_strtol (s, &s); |
965 |
r= 0.0; |
966 |
wasproject= False; |
967 |
if (*s == ':') { |
968 |
s++; |
969 |
if ((r= qh_strtod(s, &s)) == 0.0) { |
970 |
t= s; /* need true dimension for memory allocation */ |
971 |
while (*t && !isspace(*t)) { |
972 |
if (toupper(*t++) == 'B' |
973 |
&& k == qh_strtol (t, &t) |
974 |
&& *t++ == ':' |
975 |
&& qh_strtod(t, &t) == 0.0) { |
976 |
qh PROJECTinput++; |
977 |
trace2((qh ferr, "qh_initflags: project dimension %d\n", k)); |
978 |
qh_option ("Qb-project-dim", &k, NULL); |
979 |
wasproject= True; |
980 |
lastproject= k; |
981 |
break; |
982 |
} |
983 |
} |
984 |
} |
985 |
} |
986 |
if (!wasproject) { |
987 |
if (lastproject == k && r == 0.0) |
988 |
lastproject= -1; /* doesn't catch all possible sequences */ |
989 |
else if (key == 'b') { |
990 |
qh SCALEinput= True; |
991 |
if (r == 0.0) |
992 |
r= -qh_DEFAULTbox; |
993 |
qh_option ("Qbound-dim-low", &k, &r); |
994 |
}else { |
995 |
qh SCALEinput= True; |
996 |
if (r == 0.0) |
997 |
r= qh_DEFAULTbox; |
998 |
qh_option ("QBound-dim-high", &k, &r); |
999 |
} |
1000 |
} |
1001 |
break; |
1002 |
case 'c': |
1003 |
qh_option ("Qcoplanar-keep", NULL, NULL); |
1004 |
qh KEEPcoplanar= True; |
1005 |
break; |
1006 |
case 'f': |
1007 |
qh_option ("Qfurthest-outside", NULL, NULL); |
1008 |
qh BESToutside= True; |
1009 |
break; |
1010 |
case 'g': |
1011 |
qh_option ("Qgood-facets-only", NULL, NULL); |
1012 |
qh ONLYgood= True; |
1013 |
break; |
1014 |
case 'i': |
1015 |
qh_option ("Qinterior-keep", NULL, NULL); |
1016 |
qh KEEPinside= True; |
1017 |
break; |
1018 |
case 'm': |
1019 |
qh_option ("Qmax-outside-only", NULL, NULL); |
1020 |
qh ONLYmax= True; |
1021 |
break; |
1022 |
case 'r': |
1023 |
qh_option ("Qrandom-outside", NULL, NULL); |
1024 |
qh RANDOMoutside= True; |
1025 |
break; |
1026 |
case 's': |
1027 |
qh_option ("Qsearch-initial-simplex", NULL, NULL); |
1028 |
qh ALLpoints= True; |
1029 |
break; |
1030 |
case 't': |
1031 |
qh_option ("Qtriangulate", NULL, NULL); |
1032 |
qh TRIangulate= True; |
1033 |
break; |
1034 |
case 'T': |
1035 |
qh_option ("QTestPoints", NULL, NULL); |
1036 |
if (!isdigit (*s)) |
1037 |
fprintf (qh ferr, "qhull input error: missing number of test points for option 'QTn'\n"); |
1038 |
else { |
1039 |
qh TESTpoints= qh_strtol (s, &s); |
1040 |
qh_option ("QTestPoints", &qh TESTpoints, NULL); |
1041 |
} |
1042 |
break; |
1043 |
case 'u': |
1044 |
qh_option ("QupperDelaunay", NULL, NULL); |
1045 |
qh UPPERdelaunay= True; |
1046 |
break; |
1047 |
case 'v': |
1048 |
qh_option ("Qvertex-neighbors-convex", NULL, NULL); |
1049 |
qh TESTvneighbors= True; |
1050 |
break; |
1051 |
case 'x': |
1052 |
qh_option ("Qxact-merge", NULL, NULL); |
1053 |
qh MERGEexact= True; |
1054 |
break; |
1055 |
case 'z': |
1056 |
qh_option ("Qz-infinity-point", NULL, NULL); |
1057 |
qh ATinfinity= True; |
1058 |
break; |
1059 |
case '0': |
1060 |
qh_option ("Q0-no-premerge", NULL, NULL); |
1061 |
qh NOpremerge= True; |
1062 |
break; |
1063 |
case '1': |
1064 |
if (!isdigit(*s)) { |
1065 |
qh_option ("Q1-no-angle-sort", NULL, NULL); |
1066 |
qh ANGLEmerge= False; |
1067 |
break; |
1068 |
} |
1069 |
switch(*s++) { |
1070 |
case '0': |
1071 |
qh_option ("Q10-no-narrow", NULL, NULL); |
1072 |
qh NOnarrow= True; |
1073 |
break; |
1074 |
case '1': |
1075 |
qh_option ("Q11-trinormals Qtriangulate", NULL, NULL); |
1076 |
qh TRInormals= True; |
1077 |
qh TRIangulate= True; |
1078 |
break; |
1079 |
default: |
1080 |
s--; |
1081 |
fprintf (qh ferr, "qhull warning: unknown 'Q' qhull option 1%c, rest ignored\n", (int)s[0]); |
1082 |
while (*++s && !isspace(*s)); |
1083 |
break; |
1084 |
} |
1085 |
break; |
1086 |
case '2': |
1087 |
qh_option ("Q2-no-merge-independent", NULL, NULL); |
1088 |
qh MERGEindependent= False; |
1089 |
goto LABELcheckdigit; |
1090 |
break; /* no warnings */ |
1091 |
case '3': |
1092 |
qh_option ("Q3-no-merge-vertices", NULL, NULL); |
1093 |
qh MERGEvertices= False; |
1094 |
LABELcheckdigit: |
1095 |
if (isdigit(*s)) |
1096 |
fprintf (qh ferr, "qhull warning: can not follow '1', '2', or '3' with a digit. '%c' skipped.\n", |
1097 |
*s++); |
1098 |
break; |
1099 |
case '4': |
1100 |
qh_option ("Q4-avoid-old-into-new", NULL, NULL); |
1101 |
qh AVOIDold= True; |
1102 |
break; |
1103 |
case '5': |
1104 |
qh_option ("Q5-no-check-outer", NULL, NULL); |
1105 |
qh SKIPcheckmax= True; |
1106 |
break; |
1107 |
case '6': |
1108 |
qh_option ("Q6-no-concave-merge", NULL, NULL); |
1109 |
qh SKIPconvex= True; |
1110 |
break; |
1111 |
case '7': |
1112 |
qh_option ("Q7-no-breadth-first", NULL, NULL); |
1113 |
qh VIRTUALmemory= True; |
1114 |
break; |
1115 |
case '8': |
1116 |
qh_option ("Q8-no-near-inside", NULL, NULL); |
1117 |
qh NOnearinside= True; |
1118 |
break; |
1119 |
case '9': |
1120 |
qh_option ("Q9-pick-furthest", NULL, NULL); |
1121 |
qh PICKfurthest= True; |
1122 |
break; |
1123 |
case 'G': |
1124 |
i= qh_strtol (s, &t); |
1125 |
if (qh GOODpoint) |
1126 |
fprintf (qh ferr, "qhull warning: good point already defined for option 'QGn'. Ignored\n"); |
1127 |
else if (s == t) |
1128 |
fprintf (qh ferr, "qhull warning: missing good point id for option 'QGn'. Ignored\n"); |
1129 |
else if (i < 0 || *s == '-') { |
1130 |
qh GOODpoint= i-1; |
1131 |
qh_option ("QGood-if-dont-see-point", &i, NULL); |
1132 |
}else { |
1133 |
qh GOODpoint= i+1; |
1134 |
qh_option ("QGood-if-see-point", &i, NULL); |
1135 |
} |
1136 |
s= t; |
1137 |
break; |
1138 |
case 'J': |
1139 |
if (!isdigit(*s) && *s != '-') |
1140 |
qh JOGGLEmax= 0.0; |
1141 |
else { |
1142 |
qh JOGGLEmax= (realT) qh_strtod (s, &s); |
1143 |
qh_option ("QJoggle", NULL, &qh JOGGLEmax); |
1144 |
} |
1145 |
break; |
1146 |
case 'R': |
1147 |
if (!isdigit(*s) && *s != '-') |
1148 |
fprintf (qh ferr, "qhull warning: missing random seed for option 'QRn'. Ignored\n"); |
1149 |
else { |
1150 |
qh ROTATErandom= i= qh_strtol(s, &s); |
1151 |
if (i > 0) |
1152 |
qh_option ("QRotate-id", &i, NULL ); |
1153 |
else if (i < -1) |
1154 |
qh_option ("QRandom-seed", &i, NULL ); |
1155 |
} |
1156 |
break; |
1157 |
case 'V': |
1158 |
i= qh_strtol (s, &t); |
1159 |
if (qh GOODvertex) |
1160 |
fprintf (qh ferr, "qhull warning: good vertex already defined for option 'QVn'. Ignored\n"); |
1161 |
else if (s == t) |
1162 |
fprintf (qh ferr, "qhull warning: no good point id given for option 'QVn'. Ignored\n"); |
1163 |
else if (i < 0) { |
1164 |
qh GOODvertex= i - 1; |
1165 |
qh_option ("QV-good-facets-not-point", &i, NULL); |
1166 |
}else { |
1167 |
qh_option ("QV-good-facets-point", &i, NULL); |
1168 |
qh GOODvertex= i + 1; |
1169 |
} |
1170 |
s= t; |
1171 |
break; |
1172 |
default: |
1173 |
s--; |
1174 |
fprintf (qh ferr, "qhull warning: unknown 'Q' qhull option %c, rest ignored\n", (int)s[0]); |
1175 |
while (*++s && !isspace(*s)); |
1176 |
break; |
1177 |
} |
1178 |
} |
1179 |
break; |
1180 |
case 'T': |
1181 |
while (*s && !isspace(*s)) { |
1182 |
if (isdigit(*s) || *s == '-') |
1183 |
qh IStracing= qh_strtol(s, &s); |
1184 |
else switch(*s++) { |
1185 |
case 'c': |
1186 |
qh_option ("Tcheck-frequently", NULL, NULL); |
1187 |
qh CHECKfrequently= True; |
1188 |
break; |
1189 |
case 's': |
1190 |
qh_option ("Tstatistics", NULL, NULL); |
1191 |
qh PRINTstatistics= True; |
1192 |
break; |
1193 |
case 'v': |
1194 |
qh_option ("Tverify", NULL, NULL); |
1195 |
qh VERIFYoutput= True; |
1196 |
break; |
1197 |
case 'z': |
1198 |
if (!qh fout) |
1199 |
fprintf (qh ferr, "qhull warning: output file undefined (stdout). Option 'Tz' ignored.\n"); |
1200 |
else { |
1201 |
qh_option ("Tz-stdout", NULL, NULL); |
1202 |
qh ferr= qh fout; |
1203 |
qhmem.ferr= qh fout; |
1204 |
} |
1205 |
break; |
1206 |
case 'C': |
1207 |
if (!isdigit(*s)) |
1208 |
fprintf (qh ferr, "qhull warning: missing point id for cone for trace option 'TCn'. Ignored\n"); |
1209 |
else { |
1210 |
i= qh_strtol (s, &s); |
1211 |
qh_option ("TCone-stop", &i, NULL); |
1212 |
qh STOPcone= i + 1; |
1213 |
} |
1214 |
break; |
1215 |
case 'F': |
1216 |
if (!isdigit(*s)) |
1217 |
fprintf (qh ferr, "qhull warning: missing frequency count for trace option 'TFn'. Ignored\n"); |
1218 |
else { |
1219 |
qh REPORTfreq= qh_strtol (s, &s); |
1220 |
qh_option ("TFacet-log", &qh REPORTfreq, NULL); |
1221 |
qh REPORTfreq2= qh REPORTfreq/2; /* for tracemerging() */ |
1222 |
} |
1223 |
break; |
1224 |
case 'I': |
1225 |
if (s[0] != ' ' || s[1] == '\"' || s[1] == '\'' ||isspace (s[1])) { |
1226 |
s++; |
1227 |
fprintf (qh ferr, "qhull warning: option 'TI' mistyped.\nUse 'TI', one space, file name, and space or end-of-line.\nDo not use quotes. Option 'FI' ignored.\n"); |
1228 |
}else { /* not a procedure because of qh_option (filename, NULL, NULL); */ |
1229 |
char filename[500], *t= filename; |
1230 |
|
1231 |
s++; |
1232 |
while (*s) { |
1233 |
if (t - filename >= sizeof (filename)-2) { |
1234 |
fprintf (qh ferr, "qhull error: filename for 'TI' too long.\n"); |
1235 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1236 |
} |
1237 |
if (isspace (*s)) |
1238 |
break; |
1239 |
*(t++)= *s++; |
1240 |
} |
1241 |
*t= '\0'; |
1242 |
if (!freopen (filename, "r", stdin)) { |
1243 |
fprintf (qh ferr, "qhull error: could not open file \"%s\".", filename); |
1244 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1245 |
}else { |
1246 |
qh_option ("TInput-file", NULL, NULL); |
1247 |
qh_option (filename, NULL, NULL); |
1248 |
} |
1249 |
} |
1250 |
break; |
1251 |
case 'O': |
1252 |
if (s[0] != ' ' || s[1] == '\"' || isspace (s[1])) { |
1253 |
s++; |
1254 |
fprintf (qh ferr, "qhull warning: option 'TO' mistyped.\nUse 'TO', one space, file name, and space or end-of-line.\nThe file name may be enclosed in single quotes.\nDo not use double quotes. Option 'FO' ignored.\n"); |
1255 |
}else { /* not a procedure because of qh_option (filename, NULL, NULL); */ |
1256 |
char filename[500], *t= filename; |
1257 |
boolT isquote= False; |
1258 |
|
1259 |
s++; |
1260 |
if (*s == '\'') { |
1261 |
isquote= True; |
1262 |
s++; |
1263 |
} |
1264 |
while (*s) { |
1265 |
if (t - filename >= sizeof (filename)-2) { |
1266 |
fprintf (qh ferr, "qhull error: filename for 'TO' too long.\n"); |
1267 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1268 |
} |
1269 |
if (isquote) { |
1270 |
if (*s == '\'') { |
1271 |
s++; |
1272 |
isquote= False; |
1273 |
break; |
1274 |
} |
1275 |
}else if (isspace (*s)) |
1276 |
break; |
1277 |
*(t++)= *s++; |
1278 |
} |
1279 |
*t= '\0'; |
1280 |
if (isquote) |
1281 |
fprintf (qh ferr, "qhull error: missing end quote for option 'TO'. Rest of line ignored.\n"); |
1282 |
else if (!freopen (filename, "w", stdout)) { |
1283 |
fprintf (qh ferr, "qhull error: could not open file \"%s\".", filename); |
1284 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1285 |
}else { |
1286 |
qh_option ("TOutput-file", NULL, NULL); |
1287 |
qh_option (filename, NULL, NULL); |
1288 |
} |
1289 |
} |
1290 |
break; |
1291 |
case 'P': |
1292 |
if (!isdigit(*s)) |
1293 |
fprintf (qh ferr, "qhull warning: missing point id for trace option 'TPn'. Ignored\n"); |
1294 |
else { |
1295 |
qh TRACEpoint= qh_strtol (s, &s); |
1296 |
qh_option ("Trace-point", &qh TRACEpoint, NULL); |
1297 |
} |
1298 |
break; |
1299 |
case 'M': |
1300 |
if (!isdigit(*s)) |
1301 |
fprintf (qh ferr, "qhull warning: missing merge id for trace option 'TMn'. Ignored\n"); |
1302 |
else { |
1303 |
qh TRACEmerge= qh_strtol (s, &s); |
1304 |
qh_option ("Trace-merge", &qh TRACEmerge, NULL); |
1305 |
} |
1306 |
break; |
1307 |
case 'R': |
1308 |
if (!isdigit(*s)) |
1309 |
fprintf (qh ferr, "qhull warning: missing rerun count for trace option 'TRn'. Ignored\n"); |
1310 |
else { |
1311 |
qh RERUN= qh_strtol (s, &s); |
1312 |
qh_option ("TRerun", &qh RERUN, NULL); |
1313 |
} |
1314 |
break; |
1315 |
case 'V': |
1316 |
i= qh_strtol (s, &t); |
1317 |
if (s == t) |
1318 |
fprintf (qh ferr, "qhull warning: missing furthest point id for trace option 'TVn'. Ignored\n"); |
1319 |
else if (i < 0) { |
1320 |
qh STOPpoint= i - 1; |
1321 |
qh_option ("TV-stop-before-point", &i, NULL); |
1322 |
}else { |
1323 |
qh STOPpoint= i + 1; |
1324 |
qh_option ("TV-stop-after-point", &i, NULL); |
1325 |
} |
1326 |
s= t; |
1327 |
break; |
1328 |
case 'W': |
1329 |
if (!isdigit(*s)) |
1330 |
fprintf (qh ferr, "qhull warning: missing max width for trace option 'TWn'. Ignored\n"); |
1331 |
else { |
1332 |
qh TRACEdist= (realT) qh_strtod (s, &s); |
1333 |
qh_option ("TWide-trace", NULL, &qh TRACEdist); |
1334 |
} |
1335 |
break; |
1336 |
default: |
1337 |
s--; |
1338 |
fprintf (qh ferr, "qhull warning: unknown 'T' trace option %c, rest ignored\n", (int)s[0]); |
1339 |
while (*++s && !isspace(*s)); |
1340 |
break; |
1341 |
} |
1342 |
} |
1343 |
break; |
1344 |
default: |
1345 |
fprintf (qh ferr, "qhull warning: unknown flag %c (%x)\n", (int)s[-1], |
1346 |
(int)s[-1]); |
1347 |
break; |
1348 |
} |
1349 |
if (s-1 == prev_s && *s && !isspace(*s)) { |
1350 |
fprintf (qh ferr, "qhull warning: missing space after flag %c (%x); reserved for menu. Skipped.\n", |
1351 |
(int)*prev_s, (int)*prev_s); |
1352 |
while (*s && !isspace(*s)) |
1353 |
s++; |
1354 |
} |
1355 |
} |
1356 |
if (isgeom && !qh FORCEoutput && qh PRINTout[1]) |
1357 |
fprintf (qh ferr, "qhull warning: additional output formats are not compatible with Geomview\n"); |
1358 |
/* set derived values in qh_initqhull_globals */ |
1359 |
} /* initflags */ |
1360 |
|
1361 |
|
1362 |
/*-<a href="qh-globa.htm#TOC" |
1363 |
>-------------------------------</a><a name="initqhull_buffers">-</a> |
1364 |
|
1365 |
qh_initqhull_buffers() |
1366 |
initialize global memory buffers |
1367 |
|
1368 |
notes: |
1369 |
must match qh_freebuffers() |
1370 |
*/ |
1371 |
void qh_initqhull_buffers (void) { |
1372 |
int k; |
1373 |
|
1374 |
qh TEMPsize= (qhmem.LASTsize - sizeof (setT))/SETelemsize; |
1375 |
if (qh TEMPsize <= 0 || qh TEMPsize > qhmem.LASTsize) |
1376 |
qh TEMPsize= 8; /* e.g., if qh_NOmem */ |
1377 |
qh other_points= qh_setnew (qh TEMPsize); |
1378 |
qh del_vertices= qh_setnew (qh TEMPsize); |
1379 |
qh coplanarset= qh_setnew (qh TEMPsize); |
1380 |
qh NEARzero= (realT *)qh_memalloc(qh hull_dim * sizeof(realT)); |
1381 |
qh lower_threshold= (realT *)qh_memalloc((qh input_dim+1) * sizeof(realT)); |
1382 |
qh upper_threshold= (realT *)qh_memalloc((qh input_dim+1) * sizeof(realT)); |
1383 |
qh lower_bound= (realT *)qh_memalloc((qh input_dim+1) * sizeof(realT)); |
1384 |
qh upper_bound= (realT *)qh_memalloc((qh input_dim+1) * sizeof(realT)); |
1385 |
for(k= qh input_dim+1; k--; ) { |
1386 |
qh lower_threshold[k]= -REALmax; |
1387 |
qh upper_threshold[k]= REALmax; |
1388 |
qh lower_bound[k]= -REALmax; |
1389 |
qh upper_bound[k]= REALmax; |
1390 |
} |
1391 |
qh gm_matrix= (coordT *)qh_memalloc((qh hull_dim+1) * qh hull_dim * sizeof(coordT)); |
1392 |
qh gm_row= (coordT **)qh_memalloc((qh hull_dim+1) * sizeof(coordT *)); |
1393 |
} /* initqhull_buffers */ |
1394 |
|
1395 |
/*-<a href="qh-globa.htm#TOC" |
1396 |
>-------------------------------</a><a name="initqhull_globals">-</a> |
1397 |
|
1398 |
qh_initqhull_globals( points, numpoints, dim, ismalloc ) |
1399 |
initialize globals |
1400 |
if ismalloc |
1401 |
points were malloc'd and qhull should free at end |
1402 |
|
1403 |
returns: |
1404 |
sets qh.first_point, num_points, input_dim, hull_dim and others |
1405 |
seeds random number generator (seed=1 if tracing) |
1406 |
modifies qh.hull_dim if ((qh.DELAUNAY and qh.PROJECTdelaunay) or qh.PROJECTinput) |
1407 |
adjust user flags as needed |
1408 |
also checks DIM3 dependencies and constants |
1409 |
|
1410 |
notes: |
1411 |
do not use qh_point() since an input transformation may move them elsewhere |
1412 |
|
1413 |
see: |
1414 |
qh_initqhull_start() sets default values for non-zero globals |
1415 |
|
1416 |
design: |
1417 |
initialize points array from input arguments |
1418 |
test for qh.ZEROcentrum |
1419 |
(i.e., use opposite vertex instead of cetrum for convexity testing) |
1420 |
test for qh.PRINTgood (i.e., only print 'good' facets) |
1421 |
initialize qh.CENTERtype, qh.normal_size, |
1422 |
qh.center_size, qh.TRACEpoint/level, |
1423 |
initialize and test random numbers |
1424 |
check for conflicting print output options |
1425 |
*/ |
1426 |
void qh_initqhull_globals (coordT *points, int numpoints, int dim, boolT ismalloc) { |
1427 |
int seed, pointsneeded, extra= 0, i, randi, k; |
1428 |
boolT printgeom= False, printmath= False, printcoplanar= False; |
1429 |
realT randr; |
1430 |
realT factorial; |
1431 |
|
1432 |
time_t timedata; |
1433 |
|
1434 |
trace0((qh ferr, "qh_initqhull_globals: for %s | %s\n", qh rbox_command, qh qhull_command)); |
1435 |
qh POINTSmalloc= ismalloc; |
1436 |
qh first_point= points; |
1437 |
qh num_points= numpoints; |
1438 |
qh hull_dim= qh input_dim= dim; |
1439 |
if (!qh NOpremerge && !qh MERGEexact && !qh PREmerge && qh JOGGLEmax > REALmax/2) { |
1440 |
qh MERGING= True; |
1441 |
if (qh hull_dim <= 4) { |
1442 |
qh PREmerge= True; |
1443 |
qh_option ("_pre-merge", NULL, NULL); |
1444 |
}else { |
1445 |
qh MERGEexact= True; |
1446 |
qh_option ("Qxact_merge", NULL, NULL); |
1447 |
} |
1448 |
}else if (qh MERGEexact) |
1449 |
qh MERGING= True; |
1450 |
if (!qh NOpremerge && qh JOGGLEmax > REALmax/2) { |
1451 |
#ifdef qh_NOmerge |
1452 |
qh JOGGLEmax= 0.0; |
1453 |
#endif |
1454 |
} |
1455 |
if (qh TRIangulate && qh JOGGLEmax < REALmax/2 && qh PRINTprecision) |
1456 |
fprintf(qh ferr, "qhull warning: joggle ('QJ') always produces simplicial output. Triangulated output ('Qt') does nothing.\n"); |
1457 |
if (qh JOGGLEmax < REALmax/2 && qh DELAUNAY && !qh SCALEinput && !qh SCALElast) { |
1458 |
qh SCALElast= True; |
1459 |
qh_option ("Qbbound-last-qj", NULL, NULL); |
1460 |
} |
1461 |
if (qh MERGING && !qh POSTmerge && qh premerge_cos > REALmax/2 |
1462 |
&& qh premerge_centrum == 0) { |
1463 |
qh ZEROcentrum= True; |
1464 |
qh ZEROall_ok= True; |
1465 |
qh_option ("_zero-centrum", NULL, NULL); |
1466 |
} |
1467 |
if (qh JOGGLEmax < REALmax/2 && REALepsilon > 2e-8 && qh PRINTprecision) |
1468 |
fprintf(qh ferr, "qhull warning: real epsilon, %2.2g, is probably too large for joggle ('QJn')\nRecompile with double precision reals (see user.h).\n", |
1469 |
REALepsilon); |
1470 |
#ifdef qh_NOmerge |
1471 |
if (qh MERGING) { |
1472 |
fprintf (qh ferr, "qhull input error: merging not installed (qh_NOmerge + 'Qx', 'Cn' or 'An')\n"); |
1473 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1474 |
} |
1475 |
#endif |
1476 |
if (!(qh PRINTgood || qh PRINTneighbors)) { |
1477 |
if (qh KEEParea || qh KEEPminArea < REALmax/2 || qh KEEPmerge || qh DELAUNAY |
1478 |
|| (!qh ONLYgood && (qh GOODvertex || qh GOODpoint))) { |
1479 |
qh PRINTgood= True; |
1480 |
qh_option ("Pgood", NULL, NULL); |
1481 |
} |
1482 |
} |
1483 |
if (qh DELAUNAY && qh KEEPcoplanar && !qh KEEPinside) { |
1484 |
qh KEEPinside= True; |
1485 |
qh_option ("Qinterior-keep", NULL, NULL); |
1486 |
} |
1487 |
if (qh DELAUNAY && qh HALFspace) { |
1488 |
fprintf (qh ferr, "qhull input error: can not use Delaunay ('d') or Voronoi ('v') with halfspace intersection ('H')\n"); |
1489 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1490 |
} |
1491 |
if (!qh DELAUNAY && (qh UPPERdelaunay || qh ATinfinity)) { |
1492 |
fprintf (qh ferr, "qhull input error: use upper-Delaunay ('Qu') or infinity-point ('Qz') with Delaunay ('d') or Voronoi ('v')\n"); |
1493 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1494 |
} |
1495 |
if (qh UPPERdelaunay && qh ATinfinity) { |
1496 |
fprintf (qh ferr, "qhull input error: can not use infinity-point ('Qz') with upper-Delaunay ('Qu')\n"); |
1497 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1498 |
} |
1499 |
if (qh SCALElast && !qh DELAUNAY && qh PRINTprecision) |
1500 |
fprintf (qh ferr, "qhull input warning: option 'Qbb' (scale-last-coordinate) is normally used with 'd' or 'v'\n"); |
1501 |
qh DOcheckmax= (!qh SKIPcheckmax && qh MERGING ); |
1502 |
qh KEEPnearinside= (qh DOcheckmax && !(qh KEEPinside && qh KEEPcoplanar) |
1503 |
&& !qh NOnearinside); |
1504 |
if (qh MERGING) |
1505 |
qh CENTERtype= qh_AScentrum; |
1506 |
else if (qh VORONOI) |
1507 |
qh CENTERtype= qh_ASvoronoi; |
1508 |
if (qh TESTvneighbors && !qh MERGING) { |
1509 |
fprintf(qh ferr, "qhull input error: test vertex neighbors ('Qv') needs a merge option\n"); |
1510 |
qh_errexit (qh_ERRinput, NULL ,NULL); |
1511 |
} |
1512 |
if (qh PROJECTinput || (qh DELAUNAY && qh PROJECTdelaunay)) { |
1513 |
qh hull_dim -= qh PROJECTinput; |
1514 |
if (qh DELAUNAY) { |
1515 |
qh hull_dim++; |
1516 |
extra= 1; |
1517 |
} |
1518 |
} |
1519 |
if (qh hull_dim <= 1) { |
1520 |
fprintf(qh ferr, "qhull error: dimension %d must be > 1\n", qh hull_dim); |
1521 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1522 |
} |
1523 |
for (k= 2, factorial=1.0; k < qh hull_dim; k++) |
1524 |
factorial *= k; |
1525 |
qh AREAfactor= 1.0 / factorial; |
1526 |
trace2((qh ferr, "qh_initqhull_globals: initialize globals. dim %d numpoints %d malloc? %d projected %d to hull_dim %d\n", dim, numpoints, ismalloc, qh PROJECTinput, qh hull_dim)); |
1527 |
qh normal_size= qh hull_dim * sizeof(coordT); |
1528 |
qh center_size= qh normal_size - sizeof(coordT); |
1529 |
pointsneeded= qh hull_dim+1; |
1530 |
if (qh hull_dim > qh_DIMmergeVertex) { |
1531 |
qh MERGEvertices= False; |
1532 |
qh_option ("Q3-no-merge-vertices-dim-high", NULL, NULL); |
1533 |
} |
1534 |
if (qh GOODpoint) |
1535 |
pointsneeded++; |
1536 |
#ifdef qh_NOtrace |
1537 |
if (qh IStracing) { |
1538 |
fprintf (qh ferr, "qhull input error: tracing is not installed (qh_NOtrace in user.h)"); |
1539 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
1540 |
} |
1541 |
#endif |
1542 |
if (qh RERUN > 1) { |
1543 |
qh TRACElastrun= qh IStracing; /* qh_build_withrestart duplicates next conditional */ |
1544 |
if (qh IStracing != -1) |
1545 |
qh IStracing= 0; |
1546 |
}else if (qh TRACEpoint != -1 || qh TRACEdist < REALmax/2 || qh TRACEmerge) { |
1547 |
qh TRACElevel= (qh IStracing? qh IStracing : 3); |
1548 |
qh IStracing= 0; |
1549 |
} |
1550 |
if (qh ROTATErandom == 0 || qh ROTATErandom == -1) { |
1551 |
seed= time (&timedata); |
1552 |
if (qh ROTATErandom == -1) { |
1553 |
seed= -seed; |
1554 |
qh_option ("QRandom-seed", &seed, NULL ); |
1555 |
}else |
1556 |
qh_option ("QRotate-random", &seed, NULL); |
1557 |
qh ROTATErandom= seed; |
1558 |
} |
1559 |
seed= qh ROTATErandom; |
1560 |
if (seed == INT_MIN) /* default value */ |
1561 |
seed= 1; |
1562 |
else if (seed < 0) |
1563 |
seed= -seed; |
1564 |
qh_RANDOMseed_(seed); |
1565 |
randr= 0.0; |
1566 |
for (i= 1000; i--; ) { |
1567 |
randi= qh_RANDOMint; |
1568 |
randr += randi; |
1569 |
if (randi > qh_RANDOMmax) { |
1570 |
fprintf (qh ferr, "\ |
1571 |
qhull configuration error (qh_RANDOMmax in user.h):\n\ |
1572 |
random integer %d > qh_RANDOMmax (%.8g)\n", |
1573 |
randi, qh_RANDOMmax); |
1574 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1575 |
} |
1576 |
} |
1577 |
qh_RANDOMseed_(seed); |
1578 |
randr = randr/1000; |
1579 |
if (randr < qh_RANDOMmax/10 |
1580 |
|| randr > qh_RANDOMmax * 5) |
1581 |
fprintf (qh ferr, "\ |
1582 |
qhull configuration warning (qh_RANDOMmax in user.h):\n\ |
1583 |
average of 1000 random integers (%.2g) is much different than expected (%.2g).\n\ |
1584 |
Is qh_RANDOMmax (%.2g) wrong?\n", |
1585 |
randr, qh_RANDOMmax/2.0, qh_RANDOMmax); |
1586 |
qh RANDOMa= 2.0 * qh RANDOMfactor/qh_RANDOMmax; |
1587 |
qh RANDOMb= 1.0 - qh RANDOMfactor; |
1588 |
if (qh_HASHfactor < 1.1) { |
1589 |
fprintf(qh ferr, "qhull internal error (qh_initqhull_globals): qh_HASHfactor %d must be at least 1.1. Qhull uses linear hash probing\n", |
1590 |
qh_HASHfactor); |
1591 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
1592 |
} |
1593 |
if (numpoints+extra < pointsneeded) { |
1594 |
fprintf(qh ferr,"qhull input error: not enough points (%d) to construct initial simplex (need %d)\n", |
1595 |
numpoints, pointsneeded); |
1596 |
qh_errexit(qh_ERRinput, NULL, NULL); |
1597 |
} |
1598 |
if (qh PRINTtransparent) { |
1599 |
if (qh hull_dim != 4 || !qh DELAUNAY || qh VORONOI || qh DROPdim >= 0) { |
1600 |
fprintf(qh ferr,"qhull input error: transparent Delaunay ('Gt') needs 3-d Delaunay ('d') w/o 'GDn'\n"); |
1601 |
qh_errexit(qh_ERRinput, NULL, NULL); |
1602 |
} |
1603 |
qh DROPdim = 3; |
1604 |
qh PRINTridges = True; |
1605 |
} |
1606 |
for (i= qh_PRINTEND; i--; ) { |
1607 |
if (qh PRINTout[i] == qh_PRINTgeom) |
1608 |
printgeom= True; |
1609 |
else if (qh PRINTout[i] == qh_PRINTmathematica || qh PRINTout[i] == qh_PRINTmaple) |
1610 |
printmath= True; |
1611 |
else if (qh PRINTout[i] == qh_PRINTcoplanars) |
1612 |
printcoplanar= True; |
1613 |
else if (qh PRINTout[i] == qh_PRINTpointnearest) |
1614 |
printcoplanar= True; |
1615 |
else if (qh PRINTout[i] == qh_PRINTpointintersect && !qh HALFspace) { |
1616 |
fprintf (qh ferr, "qhull input error: option 'Fp' is only used for \nhalfspace intersection ('Hn,n,n').\n"); |
1617 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1618 |
}else if (qh PRINTout[i] == qh_PRINTtriangles && (qh HALFspace || qh VORONOI)) { |
1619 |
fprintf (qh ferr, "qhull input error: option 'Ft' is not available for Voronoi vertices or halfspace intersection\n"); |
1620 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1621 |
}else if (qh PRINTout[i] == qh_PRINTcentrums && qh VORONOI) { |
1622 |
fprintf (qh ferr, "qhull input error: option 'FC' is not available for Voronoi vertices ('v')\n"); |
1623 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1624 |
}else if (qh PRINTout[i] == qh_PRINTvertices) { |
1625 |
if (qh VORONOI) |
1626 |
qh_option ("Fvoronoi", NULL, NULL); |
1627 |
else |
1628 |
qh_option ("Fvertices", NULL, NULL); |
1629 |
} |
1630 |
} |
1631 |
if (printcoplanar && qh DELAUNAY && qh JOGGLEmax < REALmax/2) { |
1632 |
if (qh PRINTprecision) |
1633 |
fprintf (qh ferr, "qhull input warning: 'QJ' (joggle) will usually prevent coincident input sites for options 'Fc' and 'FP'\n"); |
1634 |
} |
1635 |
if (!qh KEEPcoplanar && !qh KEEPinside && !qh ONLYgood) { |
1636 |
if ((qh PRINTcoplanar && qh PRINTspheres) || printcoplanar) { |
1637 |
qh KEEPcoplanar = True; |
1638 |
qh_option ("Qcoplanar", NULL, NULL); |
1639 |
} |
1640 |
} |
1641 |
if (printmath && (qh hull_dim > 3 || qh VORONOI)) { |
1642 |
fprintf (qh ferr, "qhull input error: Mathematica and Maple output is only available for 2-d and 3-d convex hulls and 2-d Delaunay triangulations\n"); |
1643 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1644 |
} |
1645 |
if (printgeom) { |
1646 |
if (qh hull_dim > 4) { |
1647 |
fprintf (qh ferr, "qhull input error: Geomview output is only available for 2-d, 3-d and 4-d\n"); |
1648 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1649 |
} |
1650 |
if (qh PRINTnoplanes && !(qh PRINTcoplanar + qh PRINTcentrums |
1651 |
+ qh PRINTdots + qh PRINTspheres + qh DOintersections + qh PRINTridges)) { |
1652 |
fprintf (qh ferr, "qhull input error: no output specified for Geomview\n"); |
1653 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1654 |
} |
1655 |
if (qh VORONOI && (qh hull_dim > 3 || qh DROPdim >= 0)) { |
1656 |
fprintf (qh ferr, "qhull input error: Geomview output for Voronoi diagrams only for 2-d\n"); |
1657 |
qh_errexit (qh_ERRinput, NULL, NULL); |
1658 |
} |
1659 |
/* can not warn about furthest-site Geomview output: no lower_threshold */ |
1660 |
if (qh hull_dim == 4 && qh DROPdim == -1 && |
1661 |
(qh PRINTcoplanar || qh PRINTspheres || qh PRINTcentrums)) { |
1662 |
fprintf (qh ferr, "qhull input warning: coplanars, vertices, and centrums output not\n\ |
1663 |
available for 4-d output (ignored). Could use 'GDn' instead.\n"); |
1664 |
qh PRINTcoplanar= qh PRINTspheres= qh PRINTcentrums= False; |
1665 |
} |
1666 |
} |
1667 |
qh PRINTdim= qh hull_dim; |
1668 |
if (qh DROPdim >=0) { /* after Geomview checks */ |
1669 |
if (qh DROPdim < qh hull_dim) { |
1670 |
qh PRINTdim--; |
1671 |
if (!printgeom || qh hull_dim < 3) |
1672 |
fprintf (qh ferr, "qhull input warning: drop dimension 'GD%d' is only available for 3-d/4-d Geomview\n", qh DROPdim); |
1673 |
}else |
1674 |
qh DROPdim= -1; |
1675 |
}else if (qh VORONOI) { |
1676 |
qh DROPdim= qh hull_dim-1; |
1677 |
qh PRINTdim= qh hull_dim-1; |
1678 |
} |
1679 |
} /* initqhull_globals */ |
1680 |
|
1681 |
/*-<a href="qh-globa.htm#TOC" |
1682 |
>-------------------------------</a><a name="initqhull_mem">-</a> |
1683 |
|
1684 |
qh_initqhull_mem( ) |
1685 |
initialize mem.c for qhull |
1686 |
qh.hull_dim and qh.normal_size determine some of the allocation sizes |
1687 |
if qh.MERGING, |
1688 |
includes ridgeT |
1689 |
calls qh_user_memsizes() to add up to 10 additional sizes for quick allocation |
1690 |
(see numsizes below) |
1691 |
|
1692 |
returns: |
1693 |
mem.c already for qh_memalloc/qh_memfree (errors if called beforehand) |
1694 |
|
1695 |
notes: |
1696 |
qh_produceoutput() prints memsizes |
1697 |
|
1698 |
*/ |
1699 |
void qh_initqhull_mem (void) { |
1700 |
int numsizes; |
1701 |
int i; |
1702 |
|
1703 |
numsizes= 8+10; |
1704 |
qh_meminitbuffers (qh IStracing, qh_MEMalign, numsizes, |
1705 |
qh_MEMbufsize,qh_MEMinitbuf); |
1706 |
qh_memsize(sizeof(vertexT)); |
1707 |
if (qh MERGING) { |
1708 |
qh_memsize(sizeof(ridgeT)); |
1709 |
qh_memsize(sizeof(mergeT)); |
1710 |
} |
1711 |
qh_memsize(sizeof(facetT)); |
1712 |
i= sizeof(setT) + (qh hull_dim - 1) * SETelemsize; /* ridge.vertices */ |
1713 |
qh_memsize(i); |
1714 |
qh_memsize(qh normal_size); /* normal */ |
1715 |
i += SETelemsize; /* facet.vertices, .ridges, .neighbors */ |
1716 |
qh_memsize(i); |
1717 |
qh_user_memsizes(); |
1718 |
qh_memsetup(); |
1719 |
} /* initqhull_mem */ |
1720 |
|
1721 |
/*-<a href="qh-globa.htm#TOC" |
1722 |
>-------------------------------</a><a name="initqhull_start">-</a> |
1723 |
|
1724 |
qh_initqhull_start( infile, outfile, errfile ) |
1725 |
start initialization of qhull |
1726 |
initialize statistics, stdio, default values for global variables |
1727 |
|
1728 |
see: |
1729 |
qh_maxmin() determines the precision constants |
1730 |
*/ |
1731 |
void qh_initqhull_start (FILE *infile, FILE *outfile, FILE *errfile) { |
1732 |
|
1733 |
qh_CPUclock; /* start the clock */ |
1734 |
#if qh_QHpointer |
1735 |
if (!(qh_qh= (qhT *)malloc (sizeof(qhT)))) { |
1736 |
fprintf (errfile, "qhull error (qh_initqhull_globals): insufficient memory\n"); |
1737 |
exit (qh_ERRmem); /* no error handler */ |
1738 |
} |
1739 |
memset((char *)qh_qh, 0, sizeof(qhT)); /* every field is 0, FALSE, NULL */ |
1740 |
#else |
1741 |
memset((char *)&qh_qh, 0, sizeof(qhT)); |
1742 |
#endif |
1743 |
strcat (qh qhull, "qhull"); |
1744 |
qh_initstatistics(); |
1745 |
qh ANGLEmerge= True; |
1746 |
qh DROPdim= -1; |
1747 |
qh ferr= errfile; |
1748 |
qh fin= infile; |
1749 |
qh fout= outfile; |
1750 |
qh furthest_id= -1; |
1751 |
qh JOGGLEmax= REALmax; |
1752 |
qh KEEPminArea = REALmax; |
1753 |
qh last_low= REALmax; |
1754 |
qh last_high= REALmax; |
1755 |
qh last_newhigh= REALmax; |
1756 |
qh max_outside= 0.0; |
1757 |
qh max_vertex= 0.0; |
1758 |
qh MAXabs_coord= 0.0; |
1759 |
qh MAXsumcoord= 0.0; |
1760 |
qh MAXwidth= -REALmax; |
1761 |
qh MERGEindependent= True; |
1762 |
qh MINdenom_1= fmax_(1.0/REALmax, REALmin); /* used by qh_scalepoints */ |
1763 |
qh MINoutside= 0.0; |
1764 |
qh MINvisible= REALmax; |
1765 |
qh MAXcoplanar= REALmax; |
1766 |
qh outside_err= REALmax; |
1767 |
qh premerge_centrum= 0.0; |
1768 |
qh premerge_cos= REALmax; |
1769 |
qh PRINTprecision= True; |
1770 |
qh PRINTradius= 0.0; |
1771 |
qh postmerge_cos= REALmax; |
1772 |
qh postmerge_centrum= 0.0; |
1773 |
qh ROTATErandom= INT_MIN; |
1774 |
qh MERGEvertices= True; |
1775 |
qh totarea= 0.0; |
1776 |
qh totvol= 0.0; |
1777 |
qh TRACEdist= REALmax; |
1778 |
qh TRACEpoint= -1; /* recompile or use 'TPn' */ |
1779 |
qh tracefacet_id= UINT_MAX; /* recompile to trace a facet */ |
1780 |
qh tracevertex_id= UINT_MAX; /* recompile to trace a vertex */ |
1781 |
qh_RANDOMseed_(1); |
1782 |
} /* initqhull_start */ |
1783 |
|
1784 |
/*-<a href="qh-globa.htm#TOC" |
1785 |
>-------------------------------</a><a name="initthresholds">-</a> |
1786 |
|
1787 |
qh_initthresholds( commandString ) |
1788 |
set thresholds for printing and scaling from commandString |
1789 |
|
1790 |
returns: |
1791 |
sets qh.GOODthreshold or qh.SPLITthreshold if 'Pd0D1' used |
1792 |
|
1793 |
see: |
1794 |
qh_initflags(), 'Qbk' 'QBk' 'Pdk' and 'PDk' |
1795 |
qh_inthresholds() |
1796 |
|
1797 |
design: |
1798 |
for each 'Pdn' or 'PDn' option |
1799 |
check syntax |
1800 |
set qh.lower_threshold or qh.upper_threshold |
1801 |
set qh.GOODthreshold if an unbounded threshold is used |
1802 |
set qh.SPLITthreshold if a bounded threshold is used |
1803 |
*/ |
1804 |
void qh_initthresholds(char *command) { |
1805 |
realT value; |
1806 |
int index, maxdim, k; |
1807 |
char *s= command; |
1808 |
char key; |
1809 |
|
1810 |
maxdim= qh input_dim; |
1811 |
if (qh DELAUNAY && (qh PROJECTdelaunay || qh PROJECTinput)) |
1812 |
maxdim++; |
1813 |
while (*s) { |
1814 |
if (*s == '-') |
1815 |
s++; |
1816 |
if (*s == 'P') { |
1817 |
s++; |
1818 |
while (*s && !isspace(key= *s++)) { |
1819 |
if (key == 'd' || key == 'D') { |
1820 |
if (!isdigit(*s)) { |
1821 |
fprintf(qh ferr, "qhull warning: no dimension given for Print option '%c' at: %s. Ignored\n", |
1822 |
key, s-1); |
1823 |
continue; |
1824 |
} |
1825 |
index= qh_strtol (s, &s); |
1826 |
if (index >= qh hull_dim) { |
1827 |
fprintf(qh ferr, "qhull warning: dimension %d for Print option '%c' is >= %d. Ignored\n", |
1828 |
index, key, qh hull_dim); |
1829 |
continue; |
1830 |
} |
1831 |
if (*s == ':') { |
1832 |
s++; |
1833 |
value= qh_strtod(s, &s); |
1834 |
if (fabs((double)value) > 1.0) { |
1835 |
fprintf(qh ferr, "qhull warning: value %2.4g for Print option %c is > +1 or < -1. Ignored\n", |
1836 |
value, key); |
1837 |
continue; |
1838 |
} |
1839 |
}else |
1840 |
value= 0.0; |
1841 |
if (key == 'd') |
1842 |
qh lower_threshold[index]= value; |
1843 |
else |
1844 |
qh upper_threshold[index]= value; |
1845 |
} |
1846 |
} |
1847 |
}else if (*s == 'Q') { |
1848 |
s++; |
1849 |
while (*s && !isspace(key= *s++)) { |
1850 |
if (key == 'b' && *s == 'B') { |
1851 |
s++; |
1852 |
for (k=maxdim; k--; ) { |
1853 |
qh lower_bound[k]= -qh_DEFAULTbox; |
1854 |
qh upper_bound[k]= qh_DEFAULTbox; |
1855 |
} |
1856 |
}else if (key == 'b' && *s == 'b') |
1857 |
s++; |
1858 |
else if (key == 'b' || key == 'B') { |
1859 |
if (!isdigit(*s)) { |
1860 |
fprintf(qh ferr, "qhull warning: no dimension given for Qhull option %c. Ignored\n", |
1861 |
key); |
1862 |
continue; |
1863 |
} |
1864 |
index= qh_strtol (s, &s); |
1865 |
if (index >= maxdim) { |
1866 |
fprintf(qh ferr, "qhull warning: dimension %d for Qhull option %c is >= %d. Ignored\n", |
1867 |
index, key, maxdim); |
1868 |
continue; |
1869 |
} |
1870 |
if (*s == ':') { |
1871 |
s++; |
1872 |
value= qh_strtod(s, &s); |
1873 |
}else if (key == 'b') |
1874 |
value= -qh_DEFAULTbox; |
1875 |
else |
1876 |
value= qh_DEFAULTbox; |
1877 |
if (key == 'b') |
1878 |
qh lower_bound[index]= value; |
1879 |
else |
1880 |
qh upper_bound[index]= value; |
1881 |
} |
1882 |
} |
1883 |
}else { |
1884 |
while (*s && !isspace (*s)) |
1885 |
s++; |
1886 |
} |
1887 |
while (isspace (*s)) |
1888 |
s++; |
1889 |
} |
1890 |
for (k= qh hull_dim; k--; ) { |
1891 |
if (qh lower_threshold[k] > -REALmax/2) { |
1892 |
qh GOODthreshold= True; |
1893 |
if (qh upper_threshold[k] < REALmax/2) { |
1894 |
qh SPLITthresholds= True; |
1895 |
qh GOODthreshold= False; |
1896 |
break; |
1897 |
} |
1898 |
}else if (qh upper_threshold[k] < REALmax/2) |
1899 |
qh GOODthreshold= True; |
1900 |
} |
1901 |
} /* initthresholds */ |
1902 |
|
1903 |
/*-<a href="qh-globa.htm#TOC" |
1904 |
>-------------------------------</a><a name="option">-</a> |
1905 |
|
1906 |
qh_option( option, intVal, realVal ) |
1907 |
add an option description to qh.qhull_options |
1908 |
|
1909 |
notes: |
1910 |
will be printed with statistics ('Ts') and errors |
1911 |
strlen(option) < 40 |
1912 |
*/ |
1913 |
void qh_option (char *option, int *i, realT *r) { |
1914 |
char buf[200]; |
1915 |
int len, maxlen; |
1916 |
|
1917 |
sprintf (buf, " %s", option); |
1918 |
if (i) |
1919 |
sprintf (buf+strlen(buf), " %d", *i); |
1920 |
if (r) |
1921 |
sprintf (buf+strlen(buf), " %2.2g", *r); |
1922 |
len= strlen(buf); |
1923 |
qh qhull_optionlen += len; |
1924 |
maxlen= sizeof (qh qhull_options) - len -1; |
1925 |
maximize_(maxlen, 0); |
1926 |
if (qh qhull_optionlen >= 80 && maxlen > 0) { |
1927 |
qh qhull_optionlen= len; |
1928 |
strncat (qh qhull_options, "\n", maxlen--); |
1929 |
} |
1930 |
strncat (qh qhull_options, buf, maxlen); |
1931 |
} /* option */ |
1932 |
|
1933 |
#if qh_QHpointer |
1934 |
/*-<a href="qh-globa.htm#TOC" |
1935 |
>-------------------------------</a><a name="restore_qhull">-</a> |
1936 |
|
1937 |
qh_restore_qhull( oldqh ) |
1938 |
restores a previously saved qhull |
1939 |
also restores qh_qhstat and qhmem.tempstack |
1940 |
|
1941 |
notes: |
1942 |
errors if current qhull hasn't been saved or freed |
1943 |
uses qhmem for error reporting |
1944 |
|
1945 |
NOTE 1998/5/11: |
1946 |
Freeing memory after qh_save_qhull and qh_restore_qhull |
1947 |
is complicated. The procedures will be redesigned. |
1948 |
|
1949 |
see: |
1950 |
qh_save_qhull() |
1951 |
*/ |
1952 |
void qh_restore_qhull (qhT **oldqh) { |
1953 |
|
1954 |
if (*oldqh && strcmp ((*oldqh)->qhull, "qhull")) { |
1955 |
fprintf (qhmem.ferr, "qhull internal error (qh_restore_qhull): %p is not a qhull data structure\n", |
1956 |
*oldqh); |
1957 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
1958 |
} |
1959 |
if (qh_qh) { |
1960 |
fprintf (qhmem.ferr, "qhull internal error (qh_restore_qhull): did not save or free existing qhull\n"); |
1961 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
1962 |
} |
1963 |
if (!*oldqh || !(*oldqh)->old_qhstat) { |
1964 |
fprintf (qhmem.ferr, "qhull internal error (qh_restore_qhull): did not previously save qhull %p\n", |
1965 |
*oldqh); |
1966 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
1967 |
} |
1968 |
qh_qh= *oldqh; |
1969 |
*oldqh= NULL; |
1970 |
qh_qhstat= qh old_qhstat; |
1971 |
qhmem.tempstack= qh old_tempstack; |
1972 |
trace1((qh ferr, "qh_restore_qhull: restored qhull from %p\n", *oldqh)); |
1973 |
} /* restore_qhull */ |
1974 |
|
1975 |
/*-<a href="qh-globa.htm#TOC" |
1976 |
>-------------------------------</a><a name="save_qhull">-</a> |
1977 |
|
1978 |
qh_save_qhull( ) |
1979 |
saves qhull for a later qh_restore_qhull |
1980 |
also saves qh_qhstat and qhmem.tempstack |
1981 |
|
1982 |
returns: |
1983 |
qh_qh=NULL |
1984 |
|
1985 |
notes: |
1986 |
need to initialize qhull or call qh_restore_qhull before continuing |
1987 |
|
1988 |
NOTE 1998/5/11: |
1989 |
Freeing memory after qh_save_qhull and qh_restore_qhull |
1990 |
is complicated. The procedures will be redesigned. |
1991 |
|
1992 |
see: |
1993 |
qh_restore_qhull() |
1994 |
*/ |
1995 |
qhT *qh_save_qhull (void) { |
1996 |
qhT *oldqh; |
1997 |
|
1998 |
trace1((qhmem.ferr, "qh_save_qhull: save qhull %p\n", qh_qh)); |
1999 |
if (!qh_qh) { |
2000 |
fprintf (qhmem.ferr, "qhull internal error (qh_save_qhull): qhull not initialized\n"); |
2001 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
2002 |
} |
2003 |
qh old_qhstat= qh_qhstat; |
2004 |
qh_qhstat= NULL; |
2005 |
qh old_tempstack= qhmem.tempstack; |
2006 |
qhmem.tempstack= NULL; |
2007 |
oldqh= qh_qh; |
2008 |
qh_qh= NULL; |
2009 |
return oldqh; |
2010 |
} /* save_qhull */ |
2011 |
|
2012 |
#endif |
2013 |
|
2014 |
/*-<a href="qh-globa.htm#TOC" |
2015 |
>-------------------------------</a><a name="strtol">-</a> |
2016 |
|
2017 |
qh_strtol( s, endp) qh_strtod( s, endp) |
2018 |
internal versions of strtol() and strtod() |
2019 |
does not skip trailing spaces |
2020 |
notes: |
2021 |
some implementations of strtol()/strtod() skip trailing spaces |
2022 |
*/ |
2023 |
double qh_strtod (const char *s, char **endp) { |
2024 |
double result; |
2025 |
|
2026 |
result= strtod (s, endp); |
2027 |
if (s < (*endp) && (*endp)[-1] == ' ') |
2028 |
(*endp)--; |
2029 |
return result; |
2030 |
} /* strtod */ |
2031 |
|
2032 |
int qh_strtol (const char *s, char **endp) { |
2033 |
int result; |
2034 |
|
2035 |
result= (int) strtol (s, endp, 10); |
2036 |
if (s< (*endp) && (*endp)[-1] == ' ') |
2037 |
(*endp)--; |
2038 |
return result; |
2039 |
} /* strtol */ |