Doxygen Source Code Documentation
global.c File Reference
#include "qhull_a.h"Go to the source code of this file.
Functions | |
| void | qh_appendprint (qh_PRINT format) |
| void | qh_checkflags (char *command, char *hiddenflags) |
| unsigned long | qh_clock (void) |
| void | qh_freebuffers (void) |
| void | qh_freebuild (boolT allmem) |
| void | qh_freeqhull (boolT allmem) |
| void | qh_init_A (FILE *infile, FILE *outfile, FILE *errfile, int argc, char *argv[]) |
| void | qh_init_B (coordT *points, int numpoints, int dim, boolT ismalloc) |
| void | qh_init_qhull_command (int argc, char *argv[]) |
| void | qh_initflags (char *command) |
| void | qh_initqhull_buffers (void) |
| void | qh_initqhull_globals (coordT *points, int numpoints, int dim, boolT ismalloc) |
| void | qh_initqhull_mem (void) |
| void | qh_initqhull_start (FILE *infile, FILE *outfile, FILE *errfile) |
| void | qh_initthresholds (char *command) |
| void | qh_option (char *option, int *i, realT *r) |
| double | qh_strtod (const char *s, char **endp) |
| int | qh_strtol (const char *s, char **endp) |
Variables | |
| qhT | qh_qh |
Function Documentation
|
|
Definition at line 34 of file global.c. References format, i, qh, qh_PRINT, and qh_PRINTEND. Referenced by qh_initflags().
|
|
||||||||||||
|
Definition at line 61 of file global.c. References key, qh, qh_errexit(), qh_ERRinput, and qh_strtod().
00061 {
00062 char *s= command, *t, *chkerr, key, opt, prevopt;
00063 char chkkey[]= " ";
00064 char chkopt[]= " ";
00065 char chkopt2[]= " ";
00066
00067 if (*hiddenflags != ' ' || hiddenflags[strlen(hiddenflags)-1] != ' ') {
00068 fprintf(qh ferr, "qhull error (qh_checkflags): hiddenflags must start and end with a space: \"%s\"", hiddenflags);
00069 qh_errexit(qh_ERRinput, NULL, NULL);
00070 }
00071 if (strpbrk(hiddenflags, ",\n\r\t")) {
00072 fprintf(qh ferr, "qhull error (qh_checkflags): hiddenflags contains commas, newlines, or tabs: \"%s\"", hiddenflags);
00073 qh_errexit(qh_ERRinput, NULL, NULL);
00074 }
00075 while (*s && !isspace(*s)) /* skip program name */
00076 s++;
00077 while (*s) {
00078 while (*s && isspace(*s))
00079 s++;
00080 if (*s == '-')
00081 s++;
00082 if (!*s)
00083 break;
00084 key = *s++;
00085 chkerr = NULL;
00086 if (key == '\'') { /* TO 'file name' */
00087 t= strchr(s, '\'');
00088 if (!t) {
00089 fprintf(qh ferr, "qhull error (qh_checkflags): missing the 2nd single-quote for:\n%s\n", s-1);
00090 qh_errexit(qh_ERRinput, NULL, NULL);
00091 }
00092 s= t+1;
00093 continue;
00094 }
00095 chkkey[1]= key;
00096 if (strstr(hiddenflags, chkkey)) {
00097 chkerr= chkkey;
00098 }else if (isupper(key)) {
00099 opt= ' ';
00100 prevopt= ' ';
00101 chkopt[1]= key;
00102 chkopt2[1]= key;
00103 while (!chkerr && *s && !isspace(*s)) {
00104 opt= *s++;
00105 if (isalpha(opt)) {
00106 chkopt[2]= opt;
00107 if (strstr(hiddenflags, chkopt))
00108 chkerr= chkopt;
00109 if (prevopt != ' ') {
00110 chkopt2[2]= prevopt;
00111 chkopt2[3]= opt;
00112 if (strstr(hiddenflags, chkopt2))
00113 chkerr= chkopt2;
00114 }
00115 }else if (key == 'Q' && isdigit(opt) && prevopt != 'b'
00116 && (prevopt == ' ' || islower(prevopt))) {
00117 chkopt[2]= opt;
00118 if (strstr(hiddenflags, chkopt))
00119 chkerr= chkopt;
00120 }else {
00121 qh_strtod (s-1, &t);
00122 if (s < t)
00123 s= t;
00124 }
00125 prevopt= opt;
00126 }
00127 }
00128 if (chkerr) {
00129 *chkerr= '\'';
00130 chkerr[strlen(chkerr)-1]= '\'';
00131 fprintf(qh ferr, "qhull error: option %s is not used with this program.\n It may be used with qhull.\n", chkerr);
00132 qh_errexit(qh_ERRinput, NULL, NULL);
00133 }
00134 }
00135 } /* checkflags */
|
|
|
Definition at line 148 of file global.c. References qh, qh_errexit(), qh_ERRqhull, and qh_SECticks.
00148 {
00149
00150 #if (qh_CLOCKtype == 2)
00151 struct tms time;
00152 static long clktck; /* initialized first call */
00153 double ratio, cpu;
00154 unsigned long ticks;
00155
00156 if (!clktck) {
00157 if ((clktck= sysconf (_SC_CLK_TCK)) < 0) {
00158 fprintf (qh ferr, "qhull internal error (qh_clock): sysconf() failed. Use qh_CLOCKtype 1 in user.h\n");
00159 qh_errexit (qh_ERRqhull, NULL, NULL);
00160 }
00161 }
00162 if (times (&time) == -1) {
00163 fprintf (qh ferr, "qhull internal error (qh_clock): times() failed. Use qh_CLOCKtype 1 in user.h\n");
00164 qh_errexit (qh_ERRqhull, NULL, NULL);
00165 }
00166 ratio= qh_SECticks / (double)clktck;
00167 ticks= time.tms_utime * ratio;
00168 return ticks;
00169 #else
00170 fprintf (qh ferr, "qhull internal error (qh_clock): use qh_CLOCKtype 2 in user.h\n");
00171 qh_errexit (qh_ERRqhull, NULL, NULL); /* never returns */
00172 return 0;
00173 #endif
00174 } /* clock */
|
|
|
Definition at line 185 of file global.c. References coordT, free, qh, qh_memfree(), qh_setfree(), realT, and trace5. Referenced by qh_freeqhull().
00185 {
00186
00187 trace5((qh ferr, "qh_freebuffers: freeing up global memory buffers\n"));
00188 /* allocated by qh_initqhull_buffers */
00189 qh_memfree (qh NEARzero, qh hull_dim * sizeof(realT));
00190 qh_memfree (qh lower_threshold, (qh input_dim+1) * sizeof(realT));
00191 qh_memfree (qh upper_threshold, (qh input_dim+1) * sizeof(realT));
00192 qh_memfree (qh lower_bound, (qh input_dim+1) * sizeof(realT));
00193 qh_memfree (qh upper_bound, (qh input_dim+1) * sizeof(realT));
00194 qh_memfree (qh gm_matrix, (qh hull_dim+1) * qh hull_dim * sizeof(coordT));
00195 qh_memfree (qh gm_row, (qh hull_dim+1) * sizeof(coordT *));
00196 qh NEARzero= qh lower_threshold= qh upper_threshold= NULL;
00197 qh lower_bound= qh upper_bound= NULL;
00198 qh gm_matrix= NULL;
00199 qh gm_row= NULL;
00200 qh_setfree (&qh other_points);
00201 qh_setfree (&qh del_vertices);
00202 qh_setfree (&qh searchset);
00203 if (qh line) /* allocated by qh_readinput, freed if no error */
00204 free (qh line);
00205 if (qh half_space)
00206 free (qh half_space);
00207 if (qh temp_malloc)
00208 free (qh temp_malloc);
00209 if (qh feasible_point) /* allocated by qh_readfeasible */
00210 free (qh feasible_point);
00211 if (qh feasible_string) /* allocated by qh_initflags */
00212 free (qh feasible_string);
00213 qh line= qh feasible_string= NULL;
00214 qh half_space= qh feasible_point= qh temp_malloc= NULL;
00215 /* usually allocated by qh_readinput */
00216 if (qh first_point && qh POINTSmalloc) {
00217 free(qh first_point);
00218 qh first_point= NULL;
00219 }
00220 if (qh input_points && qh input_malloc) { /* set by qh_joggleinput */
00221 free (qh input_points);
00222 qh input_points= NULL;
00223 }
00224 trace5((qh ferr, "qh_freebuffers: finished\n"));
00225 } /* freebuffers */
|
|
|
Definition at line 249 of file global.c. References boolT, facetT::coplanarset, FORALLfacets, FORALLvertices, FOREACHmerge_, FOREACHridge_, vertexT::neighbors, facetT::neighbors, vertexT::next, facetT::next, otherfacet_, facetT::outsideset, qh, qh_ASnone, qh_clearcenters(), qh_delfacet(), qh_delvertex(), qh_memfree(), qh_setfree(), qh_setfreelong(), qh_settempfree_all(), qh_settruncate(), facetT::ridges, ridgeT::seen, facetT::simplicial, trace1, ridgeT::vertices, facetT::vertices, and facetT::visible. Referenced by qh_build_withrestart(), and qh_freeqhull().
00249 {
00250 facetT *facet;
00251 vertexT *vertex;
00252 ridgeT *ridge, **ridgep;
00253 mergeT *merge, **mergep;
00254
00255 trace1((qh ferr, "qh_freebuild: free memory from qh_inithull and qh_buildhull\n"));
00256 if (qh del_vertices)
00257 qh_settruncate (qh del_vertices, 0);
00258 if (allmem) {
00259 qh_clearcenters (qh_ASnone);
00260 while ((vertex= qh vertex_list)) {
00261 if (vertex->next)
00262 qh_delvertex (vertex);
00263 else {
00264 qh_memfree (vertex, sizeof(vertexT));
00265 qh newvertex_list= qh vertex_list= NULL;
00266 }
00267 }
00268 }else if (qh VERTEXneighbors) {
00269 FORALLvertices
00270 qh_setfreelong (&(vertex->neighbors));
00271 }
00272 qh VERTEXneighbors= False;
00273 qh GOODclosest= NULL;
00274 if (allmem) {
00275 FORALLfacets {
00276 FOREACHridge_(facet->ridges)
00277 ridge->seen= False;
00278 }
00279 FORALLfacets {
00280 if (facet->visible) {
00281 FOREACHridge_(facet->ridges) {
00282 if (!otherfacet_(ridge, facet)->visible)
00283 ridge->seen= True; /* an unattached ridge */
00284 }
00285 }
00286 }
00287 while ((facet= qh facet_list)) {
00288 FOREACHridge_(facet->ridges) {
00289 if (ridge->seen) {
00290 qh_setfree(&(ridge->vertices));
00291 qh_memfree(ridge, sizeof(ridgeT));
00292 }else
00293 ridge->seen= True;
00294 }
00295 qh_setfree (&(facet->outsideset));
00296 qh_setfree (&(facet->coplanarset));
00297 qh_setfree (&(facet->neighbors));
00298 qh_setfree (&(facet->ridges));
00299 qh_setfree (&(facet->vertices));
00300 if (facet->next)
00301 qh_delfacet (facet);
00302 else {
00303 qh_memfree (facet, sizeof(facetT));
00304 qh visible_list= qh newfacet_list= qh facet_list= NULL;
00305 }
00306 }
00307 }else {
00308 FORALLfacets {
00309 qh_setfreelong (&(facet->outsideset));
00310 qh_setfreelong (&(facet->coplanarset));
00311 if (!facet->simplicial) {
00312 qh_setfreelong (&(facet->neighbors));
00313 qh_setfreelong (&(facet->ridges));
00314 qh_setfreelong (&(facet->vertices));
00315 }
00316 }
00317 }
00318 qh_setfree (&(qh hash_table));
00319 qh_memfree (qh interior_point, qh normal_size);
00320 qh interior_point= NULL;
00321 FOREACHmerge_(qh facet_mergeset) /* usually empty */
00322 qh_memfree (merge, sizeof(mergeT));
00323 qh facet_mergeset= NULL; /* temp set */
00324 qh degen_mergeset= NULL; /* temp set */
00325 qh_settempfree_all();
00326 } /* freebuild */
|
|
|
Definition at line 344 of file global.c. References boolT, free, qh, qh_freebuffers(), qh_freebuild(), qh_freestatistics(), and trace1. Referenced by main().
00344 {
00345
00346 trace1((qh ferr, "qh_freeqhull: free global memory\n"));
00347 qh NOerrexit= True; /* no more setjmp since called at exit */
00348 qh_freebuild (allmem);
00349 qh_freebuffers();
00350 qh_freestatistics();
00351 #if qh_QHpointer
00352 free (qh_qh);
00353 qh_qh= NULL;
00354 #else
00355 memset((char *)&qh_qh, 0, sizeof(qhT));
00356 qh NOerrexit= True;
00357 #endif
00358 } /* freeqhull */
|
|
||||||||||||||||||||||||
|
Definition at line 378 of file global.c. References argc, qh_init_qhull_command(), qh_initqhull_start(), and qh_meminit(). Referenced by main().
00378 {
00379 qh_meminit (errfile);
00380 qh_initqhull_start (infile, outfile, errfile);
00381 qh_init_qhull_command (argc, argv);
00382 } /* init_A */
|
|
||||||||||||||||||||
|
Definition at line 425 of file global.c. References boolT, coordT, qh, qh_gram_schmidt(), qh_initqhull_buffers(), qh_initqhull_globals(), qh_initqhull_mem(), qh_initthresholds(), qh_projectinput(), qh_randommatrix(), qh_rotateinput(), and qh_scaleinput(). Referenced by main(), and qh_new_qhull().
00425 {
00426 qh_initqhull_globals (points, numpoints, dim, ismalloc);
00427 if (qhmem.LASTsize == 0)
00428 qh_initqhull_mem();
00429 /* mem.c and qset.c are initialized */
00430 qh_initqhull_buffers();
00431 qh_initthresholds (qh qhull_command);
00432 if (qh PROJECTinput || (qh DELAUNAY && qh PROJECTdelaunay))
00433 qh_projectinput();
00434 if (qh SCALEinput)
00435 qh_scaleinput();
00436 if (qh ROTATErandom >= 0) {
00437 qh_randommatrix (qh gm_matrix, qh hull_dim, qh gm_row);
00438 if (qh DELAUNAY) {
00439 int k, lastk= qh hull_dim-1;
00440 for (k= 0; k < lastk; k++) {
00441 qh gm_row[k][lastk]= 0.0;
00442 qh gm_row[lastk][k]= 0.0;
00443 }
00444 qh gm_row[lastk][lastk]= 1.0;
00445 }
00446 qh_gram_schmidt (qh hull_dim, qh gm_row);
00447 qh_rotateinput (qh gm_row);
00448 }
00449 } /* init_B */
|
|
||||||||||||
|
Definition at line 465 of file global.c. Referenced by qh_init_A().
00465 {
00466 int i;
00467 char *s;
00468
00469 if (argc) {
00470 if ((s= strrchr( argv[0], '\\'))) /* Borland gives full path */
00471 strcpy (qh qhull_command, s+1);
00472 else
00473 strcpy (qh qhull_command, argv[0]);
00474 if ((s= strstr (qh qhull_command, ".EXE"))
00475 || (s= strstr (qh qhull_command, ".exe")))
00476 *s= '\0';
00477 }
00478 for (i=1; i < argc; i++) {
00479 if (strlen (qh qhull_command) + strlen(argv[i]) + 1 < sizeof(qh qhull_command)) {
00480 strcat (qh qhull_command, " ");
00481 strcat (qh qhull_command, argv[i]);
00482 }else {
00483 fprintf (qh ferr, "qhull input error: more than %d characters in command line\n",
00484 (int)sizeof(qh qhull_command));
00485 exit (1); /* can not use qh_errexit */
00486 }
00487 }
00488 } /* init_qhull_command */
|
|
|
Definition at line 521 of file global.c. References boolT, calloc, fout, i, key, MERGING, qh, qh_appendprint(), qh_DEFAULTbox, qh_errexit(), qh_ERRinput, qh_ERRmem, qh_option(), qh_PRINTarea, qh_PRINTaverage, qh_PRINTcentrums, qh_PRINTcoplanars, qh_PRINTextremes, qh_PRINTfacets, qh_PRINTfacets_xridge, qh_PRINTgeom, qh_PRINTids, qh_PRINTincidences, qh_PRINTinner, qh_PRINTmathematica, qh_PRINTmerges, qh_PRINTneighbors, qh_PRINTnormals, qh_PRINToff, qh_PRINToptions, qh_PRINTouter, qh_PRINTpointintersect, qh_PRINTpointnearest, qh_PRINTpoints, qh_PRINTqhull, qh_PRINTsize, qh_PRINTsummary, qh_PRINTtriangles, qh_PRINTvertices, qh_PRINTvneighbors, qh_strtod(), qh_strtol(), r, realT, and trace2. Referenced by main(), and qh_new_qhull().
00521 {
00522 int k, i, lastproject;
00523 char *s= command, *t, *prev_s, *start, key;
00524 boolT isgeom= False, wasproject;
00525 realT r;
00526
00527 if (command != &qh qhull_command[0]) {
00528 *qh qhull_command= '\0';
00529 strncat( qh qhull_command, command, sizeof( qh qhull_command));
00530 }
00531 while (*s && !isspace(*s)) /* skip program name */
00532 s++;
00533 while (*s) {
00534 while (*s && isspace(*s))
00535 s++;
00536 if (*s == '-')
00537 s++;
00538 if (!*s)
00539 break;
00540 prev_s= s;
00541 switch (*s++) {
00542 case 'd':
00543 qh_option ("delaunay", NULL, NULL);
00544 qh DELAUNAY= True;
00545 break;
00546 case 'f':
00547 qh_option ("facets", NULL, NULL);
00548 qh_appendprint (qh_PRINTfacets);
00549 break;
00550 case 'i':
00551 qh_option ("incidence", NULL, NULL);
00552 qh_appendprint (qh_PRINTincidences);
00553 break;
00554 case 'm':
00555 qh_option ("mathematica", NULL, NULL);
00556 qh_appendprint (qh_PRINTmathematica);
00557 break;
00558 case 'n':
00559 qh_option ("normals", NULL, NULL);
00560 qh_appendprint (qh_PRINTnormals);
00561 break;
00562 case 'o':
00563 qh_option ("offFile", NULL, NULL);
00564 qh_appendprint (qh_PRINToff);
00565 break;
00566 case 'p':
00567 qh_option ("points", NULL, NULL);
00568 qh_appendprint (qh_PRINTpoints);
00569 break;
00570 case 's':
00571 qh_option ("summary", NULL, NULL);
00572 qh PRINTsummary= True;
00573 break;
00574 case 'v':
00575 qh_option ("voronoi", NULL, NULL);
00576 qh VORONOI= True;
00577 qh DELAUNAY= True;
00578 break;
00579 case 'A':
00580 if (!isdigit(*s) && *s != '.' && *s != '-')
00581 fprintf(qh ferr, "qhull warning: no maximum cosine angle given for option 'An'. Ignored.\n");
00582 else {
00583 if (*s == '-') {
00584 qh premerge_cos= -qh_strtod (s, &s);
00585 qh_option ("Angle-premerge-", NULL, &qh premerge_cos);
00586 qh PREmerge= True;
00587 }else {
00588 qh postmerge_cos= qh_strtod (s, &s);
00589 qh_option ("Angle-postmerge", NULL, &qh postmerge_cos);
00590 qh POSTmerge= True;
00591 }
00592 qh MERGING= True;
00593 }
00594 break;
00595 case 'C':
00596 if (!isdigit(*s) && *s != '.' && *s != '-')
00597 fprintf(qh ferr, "qhull warning: no centrum radius given for option 'Cn'. Ignored.\n");
00598 else {
00599 if (*s == '-') {
00600 qh premerge_centrum= -qh_strtod (s, &s);
00601 qh_option ("Centrum-premerge-", NULL, &qh premerge_centrum);
00602 qh PREmerge= True;
00603 }else {
00604 qh postmerge_centrum= qh_strtod (s, &s);
00605 qh_option ("Centrum-postmerge", NULL, &qh postmerge_centrum);
00606 qh POSTmerge= True;
00607 }
00608 qh MERGING= True;
00609 }
00610 break;
00611 case 'E':
00612 if (*s == '-')
00613 fprintf(qh ferr, "qhull warning: negative maximum roundoff given for option 'An'. Ignored.\n");
00614 else if (!isdigit(*s))
00615 fprintf(qh ferr, "qhull warning: no maximum roundoff given for option 'En'. Ignored.\n");
00616 else {
00617 qh DISTround= qh_strtod (s, &s);
00618 qh_option ("Distance-roundoff", NULL, &qh DISTround);
00619 qh SETroundoff= True;
00620 }
00621 break;
00622 case 'H':
00623 start= s;
00624 qh HALFspace= True;
00625 qh_strtod (s, &t);
00626 while (t > s) {
00627 if (*t && !isspace (*t)) {
00628 if (*t == ',')
00629 t++;
00630 else
00631 fprintf (qh ferr, "qhull warning: origin for Halfspace intersection should be 'Hn,n,n,...'\n");
00632 }
00633 s= t;
00634 qh_strtod (s, &t);
00635 }
00636 if (start < t) {
00637 if (!(qh feasible_string= (char*)calloc (t-start+1, 1))) {
00638 fprintf(qh ferr, "qhull error: insufficient memory for 'Hn,n,n'\n");
00639 qh_errexit(qh_ERRmem, NULL, NULL);
00640 }
00641 strncpy (qh feasible_string, start, t-start);
00642 qh_option ("Halfspace-about", NULL, NULL);
00643 qh_option (qh feasible_string, NULL, NULL);
00644 }else
00645 qh_option ("Halfspace", NULL, NULL);
00646 break;
00647 case 'R':
00648 if (!isdigit(*s))
00649 fprintf(qh ferr, "qhull warning: missing random perturbation for option 'Rn'. Ignored\n");
00650 else {
00651 qh RANDOMfactor= qh_strtod (s, &s);
00652 qh_option ("Random_perturb", NULL, &qh RANDOMfactor);
00653 qh RANDOMdist= True;
00654 }
00655 break;
00656 case 'V':
00657 if (!isdigit(*s) && *s != '-')
00658 fprintf(qh ferr, "qhull warning: missing visible distance for option 'Vn'. Ignored\n");
00659 else {
00660 qh MINvisible= qh_strtod (s, &s);
00661 qh_option ("Visible", NULL, &qh MINvisible);
00662 }
00663 break;
00664 case 'U':
00665 if (!isdigit(*s) && *s != '-')
00666 fprintf(qh ferr, "qhull warning: missing coplanar distance for option 'Un'. Ignored\n");
00667 else {
00668 qh MAXcoplanar= qh_strtod (s, &s);
00669 qh_option ("U-coplanar", NULL, &qh MAXcoplanar);
00670 }
00671 break;
00672 case 'W':
00673 if (*s == '-')
00674 fprintf(qh ferr, "qhull warning: negative outside width for option 'Wn'. Ignored.\n");
00675 else if (!isdigit(*s))
00676 fprintf(qh ferr, "qhull warning: missing outside width for option 'Wn'. Ignored\n");
00677 else {
00678 qh MINoutside= qh_strtod (s, &s);
00679 qh_option ("W-outside", NULL, &qh MINoutside);
00680 qh APPROXhull= True;
00681 }
00682 break;
00683 /************ sub menus ***************/
00684 case 'F':
00685 while (*s && !isspace(*s)) {
00686 switch(*s++) {
00687 case 'a':
00688 qh_option ("Farea", NULL, NULL);
00689 qh_appendprint (qh_PRINTarea);
00690 qh GETarea= True;
00691 break;
00692 case 'A':
00693 qh_option ("FArea-total", NULL, NULL);
00694 qh GETarea= True;
00695 break;
00696 case 'c':
00697 qh_option ("Fcoplanars", NULL, NULL);
00698 qh_appendprint (qh_PRINTcoplanars);
00699 break;
00700 case 'C':
00701 qh_option ("FCentrums", NULL, NULL);
00702 qh_appendprint (qh_PRINTcentrums);
00703 break;
00704 case 'd':
00705 qh_option ("Fd-cdd-in", NULL, NULL);
00706 qh CDDinput= True;
00707 break;
00708 case 'D':
00709 qh_option ("FD-cdd-out", NULL, NULL);
00710 qh CDDoutput= True;
00711 break;
00712 case 'F':
00713 qh_option ("FFacets-xridge", NULL, NULL);
00714 qh_appendprint (qh_PRINTfacets_xridge);
00715 break;
00716 case 'i':
00717 qh_option ("Finner", NULL, NULL);
00718 qh_appendprint (qh_PRINTinner);
00719 break;
00720 case 'I':
00721 qh_option ("FIDs", NULL, NULL);
00722 qh_appendprint (qh_PRINTids);
00723 break;
00724 case 'm':
00725 qh_option ("Fmerges", NULL, NULL);
00726 qh_appendprint (qh_PRINTmerges);
00727 break;
00728 case 'n':
00729 qh_option ("Fneighbors", NULL, NULL);
00730 qh_appendprint (qh_PRINTneighbors);
00731 break;
00732 case 'N':
00733 qh_option ("FNeighbors-vertex", NULL, NULL);
00734 qh_appendprint (qh_PRINTvneighbors);
00735 break;
00736 case 'o':
00737 qh_option ("Fouter", NULL, NULL);
00738 qh_appendprint (qh_PRINTouter);
00739 break;
00740 case 'O':
00741 if (qh PRINToptions1st) {
00742 qh_option ("FOptions", NULL, NULL);
00743 qh_appendprint (qh_PRINToptions);
00744 }else
00745 qh PRINToptions1st= True;
00746 break;
00747 case 'p':
00748 qh_option ("Fpoint-intersect", NULL, NULL);
00749 qh_appendprint (qh_PRINTpointintersect);
00750 break;
00751 case 'P':
00752 qh_option ("FPoint-nearest", NULL, NULL);
00753 qh_appendprint (qh_PRINTpointnearest);
00754 break;
00755 case 'Q':
00756 qh_option ("FQhull", NULL, NULL);
00757 qh_appendprint (qh_PRINTqhull);
00758 break;
00759 case 's':
00760 qh_option ("Fsummary", NULL, NULL);
00761 qh_appendprint (qh_PRINTsummary);
00762 break;
00763 case 'S':
00764 qh_option ("FSize", NULL, NULL);
00765 qh_appendprint (qh_PRINTsize);
00766 qh GETarea= True;
00767 break;
00768 case 't':
00769 qh_option ("Ftriangles", NULL, NULL);
00770 qh_appendprint (qh_PRINTtriangles);
00771 break;
00772 case 'v':
00773 /* option set in qh_initqhull_globals */
00774 qh_appendprint (qh_PRINTvertices);
00775 break;
00776 case 'V':
00777 qh_option ("FVertex-average", NULL, NULL);
00778 qh_appendprint (qh_PRINTaverage);
00779 break;
00780 case 'x':
00781 qh_option ("Fxtremes", NULL, NULL);
00782 qh_appendprint (qh_PRINTextremes);
00783 break;
00784 default:
00785 s--;
00786 fprintf (qh ferr, "qhull warning: unknown 'F' output option %c, rest ignored\n", (int)s[0]);
00787 while (*++s && !isspace(*s));
00788 break;
00789 }
00790 }
00791 break;
00792 case 'G':
00793 isgeom= True;
00794 qh_appendprint (qh_PRINTgeom);
00795 while (*s && !isspace(*s)) {
00796 switch(*s++) {
00797 case 'a':
00798 qh_option ("Gall-points", NULL, NULL);
00799 qh PRINTdots= True;
00800 break;
00801 case 'c':
00802 qh_option ("Gcentrums", NULL, NULL);
00803 qh PRINTcentrums= True;
00804 break;
00805 case 'h':
00806 qh_option ("Gintersections", NULL, NULL);
00807 qh DOintersections= True;
00808 break;
00809 case 'i':
00810 qh_option ("Ginner", NULL, NULL);
00811 qh PRINTinner= True;
00812 break;
00813 case 'n':
00814 qh_option ("Gno-planes", NULL, NULL);
00815 qh PRINTnoplanes= True;
00816 break;
00817 case 'o':
00818 qh_option ("Gouter", NULL, NULL);
00819 qh PRINTouter= True;
00820 break;
00821 case 'p':
00822 qh_option ("Gpoints", NULL, NULL);
00823 qh PRINTcoplanar= True;
00824 break;
00825 case 'r':
00826 qh_option ("Gridges", NULL, NULL);
00827 qh PRINTridges= True;
00828 break;
00829 case 't':
00830 qh_option ("Gtransparent", NULL, NULL);
00831 qh PRINTtransparent= True;
00832 break;
00833 case 'v':
00834 qh_option ("Gvertices", NULL, NULL);
00835 qh PRINTspheres= True;
00836 break;
00837 case 'D':
00838 if (!isdigit (*s))
00839 fprintf (qh ferr, "qhull input error: missing dimension for option 'GDn'\n");
00840 else {
00841 if (qh DROPdim >= 0)
00842 fprintf (qh ferr, "qhull warning: can only drop one dimension. Previous 'GD%d' ignored\n",
00843 qh DROPdim);
00844 qh DROPdim= qh_strtol (s, &s);
00845 qh_option ("GDrop-dim", &qh DROPdim, NULL);
00846 }
00847 break;
00848 default:
00849 s--;
00850 fprintf (qh ferr, "qhull warning: unknown 'G' print option %c, rest ignored\n", (int)s[0]);
00851 while (*++s && !isspace(*s));
00852 break;
00853 }
00854 }
00855 break;
00856 case 'P':
00857 while (*s && !isspace(*s)) {
00858 switch(*s++) {
00859 case 'd': case 'D': /* see qh_initthresholds() */
00860 key= s[-1];
00861 i= qh_strtol (s, &s);
00862 r= 0;
00863 if (*s == ':') {
00864 s++;
00865 r= qh_strtod (s, &s);
00866 }
00867 if (key == 'd')
00868 qh_option ("Pdrop-facets-dim-less", &i, &r);
00869 else
00870 qh_option ("PDrop-facets-dim-more", &i, &r);
00871 break;
00872 case 'g':
00873 qh_option ("Pgood-facets", NULL, NULL);
00874 qh PRINTgood= True;
00875 break;
00876 case 'G':
00877 qh_option ("PGood-facet-neighbors", NULL, NULL);
00878 qh PRINTneighbors= True;
00879 break;
00880 case 'o':
00881 qh_option ("Poutput-forced", NULL, NULL);
00882 qh FORCEoutput= True;
00883 break;
00884 case 'p':
00885 qh_option ("Pprecision-ignore", NULL, NULL);
00886 qh PRINTprecision= False;
00887 break;
00888 case 'A':
00889 if (!isdigit (*s))
00890 fprintf (qh ferr, "qhull input error: missing facet count for keep area option 'PAn'\n");
00891 else {
00892 qh KEEParea= qh_strtol (s, &s);
00893 qh_option ("PArea-keep", &qh KEEParea, NULL);
00894 qh GETarea= True;
00895 }
00896 break;
00897 case 'F':
00898 if (!isdigit (*s))
00899 fprintf (qh ferr, "qhull input error: missing facet area for option 'PFn'\n");
00900 else {
00901 qh KEEPminArea= qh_strtod (s, &s);
00902 qh_option ("PFacet-area-keep", NULL, &qh KEEPminArea);
00903 qh GETarea= True;
00904 }
00905 break;
00906 case 'M':
00907 if (!isdigit (*s))
00908 fprintf (qh ferr, "qhull input error: missing merge count for option 'PMn'\n");
00909 else {
00910 qh KEEPmerge= qh_strtol (s, &s);
00911 qh_option ("PMerge-keep", &qh KEEPmerge, NULL);
00912 }
00913 break;
00914 default:
00915 s--;
00916 fprintf (qh ferr, "qhull warning: unknown 'P' print option %c, rest ignored\n", (int)s[0]);
00917 while (*++s && !isspace(*s));
00918 break;
00919 }
00920 }
00921 break;
00922 case 'Q':
00923 lastproject= -1;
00924 while (*s && !isspace(*s)) {
00925 switch(*s++) {
00926 case 'b': case 'B': /* handled by qh_initthresholds */
00927 key= s[-1];
00928 if (key == 'b' && *s == 'B') {
00929 s++;
00930 r= qh_DEFAULTbox;
00931 qh SCALEinput= True;
00932 qh_option ("QbBound-unit-box", NULL, &r);
00933 break;
00934 }
00935 if (key == 'b' && *s == 'b') {
00936 s++;
00937 qh SCALElast= True;
00938 qh_option ("Qbbound-last", NULL, NULL);
00939 break;
00940 }
00941 k= qh_strtol (s, &s);
00942 r= 0.0;
00943 wasproject= False;
00944 if (*s == ':') {
00945 s++;
00946 if ((r= qh_strtod(s, &s)) == 0.0) {
00947 t= s; /* need true dimension for memory allocation */
00948 while (*t && !isspace(*t)) {
00949 if (toupper(*t++) == 'B'
00950 && k == qh_strtol (t, &t)
00951 && *t++ == ':'
00952 && qh_strtod(t, &t) == 0.0) {
00953 qh PROJECTinput++;
00954 trace2((qh ferr, "qh_initflags: project dimension %d\n", k));
00955 qh_option ("Qb-project-dim", &k, NULL);
00956 wasproject= True;
00957 lastproject= k;
00958 break;
00959 }
00960 }
00961 }
00962 }
00963 if (!wasproject) {
00964 if (lastproject == k && r == 0.0)
00965 lastproject= -1; /* doesn't catch all possible sequences */
00966 else if (key == 'b') {
00967 qh SCALEinput= True;
00968 if (r == 0.0)
00969 r= -qh_DEFAULTbox;
00970 qh_option ("Qbound-dim-low", &k, &r);
00971 }else {
00972 qh SCALEinput= True;
00973 if (r == 0.0)
00974 r= qh_DEFAULTbox;
00975 qh_option ("QBound-dim-high", &k, &r);
00976 }
00977 }
00978 break;
00979 case 'c':
00980 qh_option ("Qcoplanar-keep", NULL, NULL);
00981 qh KEEPcoplanar= True;
00982 break;
00983 case 'f':
00984 qh_option ("Qfurthest-outside", NULL, NULL);
00985 qh BESToutside= True;
00986 break;
00987 case 'g':
00988 qh_option ("Qgood-facets-only", NULL, NULL);
00989 qh ONLYgood= True;
00990 break;
00991 case 'i':
00992 qh_option ("Qinterior-keep", NULL, NULL);
00993 qh KEEPinside= True;
00994 break;
00995 case 'm':
00996 qh_option ("Qmax-outside-only", NULL, NULL);
00997 qh ONLYmax= True;
00998 break;
00999 case 'r':
01000 qh_option ("Qrandom-outside", NULL, NULL);
01001 qh RANDOMoutside= True;
01002 break;
01003 case 's':
01004 qh_option ("Qsearch-initial-simplex", NULL, NULL);
01005 qh ALLpoints= True;
01006 break;
01007 case 'u':
01008 qh_option ("QupperDelaunay", NULL, NULL);
01009 qh UPPERdelaunay= True;
01010 break;
01011 case 'v':
01012 qh_option ("Qvertex-neighbors-convex", NULL, NULL);
01013 qh TESTvneighbors= True;
01014 break;
01015 case 'x':
01016 qh_option ("Qxact-merge", NULL, NULL);
01017 qh MERGEexact= True;
01018 break;
01019 case 'z':
01020 qh_option ("Qz-infinity-point", NULL, NULL);
01021 qh ATinfinity= True;
01022 break;
01023 case '0':
01024 qh_option ("Q0-no-premerge", NULL, NULL);
01025 qh NOpremerge= True;
01026 break;
01027 case '1':
01028 qh_option ("Q1-no-angle-sort", NULL, NULL);
01029 qh ANGLEmerge= False;
01030 goto LABELcheckdigit;
01031 break; /* no warnings */
01032 case '2':
01033 qh_option ("Q2-no-merge-independent", NULL, NULL);
01034 qh MERGEindependent= False;
01035 goto LABELcheckdigit;
01036 break; /* no warnings */
01037 case '3':
01038 qh_option ("Q3-no-merge-vertices", NULL, NULL);
01039 qh MERGEvertices= False;
01040 LABELcheckdigit:
01041 if (isdigit(*s))
01042 fprintf (qh ferr, "qhull warning: can not follow '1', '2', or '3' with a digit. '%c' skipped.\n",
01043 *s++);
01044 break;
01045 case '4':
01046 qh_option ("Q4-avoid-old-into-new", NULL, NULL);
01047 qh AVOIDold= True;
01048 break;
01049 case '5':
01050 qh_option ("Q5-no-check-outer", NULL, NULL);
01051 qh SKIPcheckmax= True;
01052 break;
01053 case '6':
01054 qh_option ("Q6-no-concave-merge", NULL, NULL);
01055 qh SKIPconvex= True;
01056 break;
01057 case '7':
01058 qh_option ("Q7-no-breadth-first", NULL, NULL);
01059 qh VIRTUALmemory= True;
01060 break;
01061 case '8':
01062 qh_option ("Q8-no-near-inside", NULL, NULL);
01063 qh NOnearinside= True;
01064 break;
01065 case '9':
01066 qh_option ("Q9-pick-furthest", NULL, NULL);
01067 qh PICKfurthest= True;
01068 break;
01069 case 'G':
01070 i= qh_strtol (s, &t);
01071 if (qh GOODpoint)
01072 fprintf (qh ferr, "qhull warning: good point already defined for option 'QGn'. Ignored\n");
01073 else if (s == t)
01074 fprintf (qh ferr, "qhull warning: missing good point id for option 'QGn'. Ignored\n");
01075 else if (i < 0 || *s == '-') {
01076 qh GOODpoint= i-1;
01077 qh_option ("QGood-if-dont-see-point", &i, NULL);
01078 }else {
01079 qh GOODpoint= i+1;
01080 qh_option ("QGood-if-see-point", &i, NULL);
01081 }
01082 s= t;
01083 break;
01084 case 'J':
01085 if (!isdigit(*s) && *s != '-')
01086 qh JOGGLEmax= 0.0;
01087 else {
01088 qh JOGGLEmax= (realT) qh_strtod (s, &s);
01089 qh_option ("QJoggle", NULL, &qh JOGGLEmax);
01090 }
01091 break;
01092 case 'R':
01093 if (!isdigit(*s) && *s != '-')
01094 fprintf (qh ferr, "qhull warning: missing random seed for option 'QRn'. Ignored\n");
01095 else {
01096 qh ROTATErandom= i= qh_strtol(s, &s);
01097 if (i > 0)
01098 qh_option ("QRotate-id", &i, NULL );
01099 else if (i < -1)
01100 qh_option ("QRandom-seed", &i, NULL );
01101 }
01102 break;
01103 case 'V':
01104 i= qh_strtol (s, &t);
01105 if (qh GOODvertex)
01106 fprintf (qh ferr, "qhull warning: good vertex already defined for option 'QVn'. Ignored\n");
01107 else if (s == t)
01108 fprintf (qh ferr, "qhull warning: no good point id given for option 'QVn'. Ignored\n");
01109 else if (i < 0) {
01110 qh GOODvertex= i - 1;
01111 qh_option ("QV-good-facets-not-point", &i, NULL);
01112 }else {
01113 qh_option ("QV-good-facets-point", &i, NULL);
01114 qh GOODvertex= i + 1;
01115 }
01116 s= t;
01117 break;
01118 default:
01119 s--;
01120 fprintf (qh ferr, "qhull warning: unknown 'Q' qhull option %c, rest ignored\n", (int)s[0]);
01121 while (*++s && !isspace(*s));
01122 break;
01123 }
01124 }
01125 break;
01126 case 'T':
01127 while (*s && !isspace(*s)) {
01128 if (isdigit(*s) || *s == '-')
01129 qh IStracing= qh_strtol(s, &s);
01130 else switch(*s++) {
01131 case 'c':
01132 qh_option ("Tcheck-frequently", NULL, NULL);
01133 qh CHECKfrequently= True;
01134 break;
01135 case 's':
01136 qh_option ("Tstatistics", NULL, NULL);
01137 qh PRINTstatistics= True;
01138 break;
01139 case 'v':
01140 qh_option ("Tverify", NULL, NULL);
01141 qh VERIFYoutput= True;
01142 break;
01143 case 'z':
01144 if (!qh fout)
01145 fprintf (qh ferr, "qhull warning: output file undefined (stdout). Option 'Tz' ignored.\n");
01146 else {
01147 qh_option ("Tz-stdout", NULL, NULL);
01148 qh ferr= qh fout;
01149 qhmem.ferr= qh fout;
01150 }
01151 break;
01152 case 'C':
01153 if (!isdigit(*s))
01154 fprintf (qh ferr, "qhull warning: missing point id for cone for trace option 'TCn'. Ignored\n");
01155 else {
01156 i= qh_strtol (s, &s);
01157 qh_option ("TCone-stop", &i, NULL);
01158 qh STOPcone= i + 1;
01159 }
01160 break;
01161 case 'F':
01162 if (!isdigit(*s))
01163 fprintf (qh ferr, "qhull warning: missing frequency count for trace option 'TFn'. Ignored\n");
01164 else {
01165 qh REPORTfreq= qh_strtol (s, &s);
01166 qh_option ("TFacet-log", &qh REPORTfreq, NULL);
01167 qh REPORTfreq2= qh REPORTfreq/2; /* for tracemerging() */
01168 }
01169 break;
01170 case 'O':
01171 if (s[0] != ' ' || s[1] == '\"' || isspace (s[1])) {
01172 s++;
01173 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");
01174 }else { /* not a procedure because of qh_option (filename, NULL, NULL); */
01175 char filename[500], *t= filename;
01176 boolT isquote= False;
01177
01178 s++;
01179 if (*s == '\'') {
01180 isquote= True;
01181 s++;
01182 }
01183 while (*s) {
01184 if (t - filename >= sizeof (filename)-2) {
01185 fprintf (qh ferr, "qhull error: filename for 'TO' too long.\n");
01186 qh_errexit (qh_ERRinput, NULL, NULL);
01187 }
01188 if (isquote) {
01189 if (*s == '\'') {
01190 s++;
01191 isquote= False;
01192 break;
01193 }
01194 }else if (isspace (*s))
01195 break;
01196 *(t++)= *s++;
01197 }
01198 *t= '\0';
01199 if (isquote)
01200 fprintf (qh ferr, "qhull error: missing end quote for option 'TO'. Rest of line ignored.\n");
01201 else if (!freopen (filename, "w", stdout)) {
01202 fprintf (qh ferr, "qhull error: could not open file \"%s\".", filename);
01203 qh_errexit (qh_ERRinput, NULL, NULL);
01204 }else {
01205 qh_option ("TOutput-file", NULL, NULL);
01206 qh_option (filename, NULL, NULL);
01207 }
01208 }
01209 break;
01210 case 'P':
01211 if (!isdigit(*s))
01212 fprintf (qh ferr, "qhull warning: missing point id for trace option 'TPn'. Ignored\n");
01213 else {
01214 qh TRACEpoint= qh_strtol (s, &s);
01215 qh_option ("Trace-point", &qh TRACEpoint, NULL);
01216 }
01217 break;
01218 case 'M':
01219 if (!isdigit(*s))
01220 fprintf (qh ferr, "qhull warning: missing merge id for trace option 'TMn'. Ignored\n");
01221 else {
01222 qh TRACEmerge= qh_strtol (s, &s);
01223 qh_option ("Trace-merge", &qh TRACEmerge, NULL);
01224 }
01225 break;
01226 case 'R':
01227 if (!isdigit(*s))
01228 fprintf (qh ferr, "qhull warning: missing rerun count for trace option 'TRn'. Ignored\n");
01229 else {
01230 qh RERUN= qh_strtol (s, &s);
01231 qh_option ("TRerun", &qh RERUN, NULL);
01232 }
01233 break;
01234 case 'V':
01235 i= qh_strtol (s, &t);
01236 if (s == t)
01237 fprintf (qh ferr, "qhull warning: missing furthest point id for trace option 'TVn'. Ignored\n");
01238 else if (i < 0) {
01239 qh STOPpoint= i - 1;
01240 qh_option ("TV-stop-before-point", &i, NULL);
01241 }else {
01242 qh STOPpoint= i + 1;
01243 qh_option ("TV-stop-after-point", &i, NULL);
01244 }
01245 s= t;
01246 break;
01247 case 'W':
01248 if (!isdigit(*s))
01249 fprintf (qh ferr, "qhull warning: missing max width for trace option 'TWn'. Ignored\n");
01250 else {
01251 qh TRACEdist= (realT) qh_strtod (s, &s);
01252 qh_option ("TWide-trace", NULL, &qh TRACEdist);
01253 }
01254 break;
01255 default:
01256 s--;
01257 fprintf (qh ferr, "qhull warning: unknown 'T' trace option %c, rest ignored\n", (int)s[0]);
01258 while (*++s && !isspace(*s));
01259 break;
01260 }
01261 }
01262 break;
01263 default:
01264 fprintf (qh ferr, "qhull warning: unknown flag %c (%x)\n", (int)s[-1],
01265 (int)s[-1]);
01266 break;
01267 }
01268 if (s-1 == prev_s && *s && !isspace(*s)) {
01269 fprintf (qh ferr, "qhull warning: missing space after flag %c (%x); reserved for menu. Skipped.\n",
01270 (int)*prev_s, (int)*prev_s);
01271 while (*s && !isspace(*s))
01272 s++;
01273 }
01274 }
01275 if (isgeom && !qh FORCEoutput && qh PRINTout[1])
01276 fprintf (qh ferr, "qhull warning: additional output formats are not compatible with Geomview\n");
01277 /* set derived values in qh_initqhull_globals */
01278 } /* initflags */
|
|
|
Definition at line 1290 of file global.c. References coordT, qh, qh_memalloc(), qh_setnew(), REALmax, realT, and SETelemsize. Referenced by qh_init_B().
01290 {
01291 int k;
01292
01293 qh TEMPsize= (qhmem.LASTsize - sizeof (setT))/SETelemsize;
01294 if (qh TEMPsize <= 0 || qh TEMPsize > qhmem.LASTsize)
01295 qh TEMPsize= 8; /* e.g., if qh_NOmem */
01296 qh other_points= qh_setnew (qh TEMPsize);
01297 qh del_vertices= qh_setnew (qh TEMPsize);
01298 qh searchset= qh_setnew (qh TEMPsize);
01299 qh NEARzero= (realT *)qh_memalloc(qh hull_dim * sizeof(realT));
01300 qh lower_threshold= (realT *)qh_memalloc((qh input_dim+1) * sizeof(realT));
01301 qh upper_threshold= (realT *)qh_memalloc((qh input_dim+1) * sizeof(realT));
01302 qh lower_bound= (realT *)qh_memalloc((qh input_dim+1) * sizeof(realT));
01303 qh upper_bound= (realT *)qh_memalloc((qh input_dim+1) * sizeof(realT));
01304 for(k= qh input_dim+1; k--; ) {
01305 qh lower_threshold[k]= -REALmax;
01306 qh upper_threshold[k]= REALmax;
01307 qh lower_bound[k]= -REALmax;
01308 qh upper_bound[k]= REALmax;
01309 }
01310 qh gm_matrix= (coordT *)qh_memalloc((qh hull_dim+1) * qh hull_dim * sizeof(coordT));
01311 qh gm_row= (coordT **)qh_memalloc((qh hull_dim+1) * sizeof(coordT *));
01312 } /* initqhull_buffers */
|
|
||||||||||||||||||||
|
Definition at line 1345 of file global.c. References boolT, coordT, i, MERGING, num_points, qh, qh_AScentrum, qh_ASvoronoi, qh_DIMmergeVertex, qh_errexit(), qh_ERRinput, qh_ERRqhull, qh_HASHfactor, qh_option(), qh_PRINTcentrums, qh_PRINTcoplanars, qh_PRINTEND, qh_PRINTgeom, qh_PRINTmathematica, qh_PRINTpointintersect, qh_PRINTpointnearest, qh_PRINTtriangles, qh_PRINTvertices, qh_RANDOMint, qh_RANDOMmax, qh_RANDOMseed_, REALepsilon, REALmax, realT, seed, trace0, and trace2. Referenced by qh_init_B().
01345 {
01346 int seed, pointsneeded, extra= 0, i, randi, k;
01347 boolT printgeom= False, printmath= False, printcoplanar= False;
01348 realT randr;
01349 realT factorial;
01350
01351 time_t timedata;
01352
01353 trace0((qh ferr, "qh_initqhull_globals: for %s | %s\n", qh rbox_command,
01354 qh qhull_command));
01355 qh POINTSmalloc= ismalloc;
01356 qh first_point= points;
01357 qh num_points= numpoints;
01358 qh hull_dim= qh input_dim= dim;
01359 if (!qh NOpremerge && !qh MERGEexact && !qh PREmerge && qh JOGGLEmax > REALmax/2) {
01360 qh MERGING= True;
01361 if (qh hull_dim <= 4) {
01362 qh PREmerge= True;
01363 qh_option ("_pre-merge", NULL, NULL);
01364 }else {
01365 qh MERGEexact= True;
01366 qh_option ("Qxact_merge", NULL, NULL);
01367 }
01368 }else if (qh MERGEexact)
01369 qh MERGING= True;
01370 if (!qh NOpremerge && qh JOGGLEmax > REALmax/2) {
01371 #ifdef qh_NOmerge
01372 qh JOGGLEmax= 0.0;
01373 #endif
01374 }
01375 if (qh JOGGLEmax < REALmax/2 && qh DELAUNAY && !qh SCALEinput && !qh SCALElast) {
01376 qh SCALElast= True;
01377 qh_option ("Qbbound-last-qj", NULL, NULL);
01378 }
01379 if (qh MERGING && !qh POSTmerge && qh premerge_cos > REALmax/2
01380 && qh premerge_centrum == 0) {
01381 qh ZEROcentrum= True;
01382 qh ZEROall_ok= True;
01383 qh_option ("_zero-centrum", NULL, NULL);
01384 }
01385 if (qh JOGGLEmax < REALmax/2 && REALepsilon > 2e-8 && qh PRINTprecision)
01386 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",
01387 REALepsilon);
01388 #ifdef qh_NOmerge
01389 if (qh MERGING) {
01390 fprintf (qh ferr, "qhull input error: merging not installed (qh_NOmerge + 'Qx', 'Cn' or 'An')\n");
01391 qh_errexit (qh_ERRinput, NULL, NULL);
01392 }
01393 #endif
01394 if (!(qh PRINTgood || qh PRINTneighbors)) {
01395 if (qh KEEParea || qh KEEPminArea < REALmax/2 || qh KEEPmerge || qh DELAUNAY
01396 || (!qh ONLYgood && (qh GOODvertex || qh GOODpoint))) {
01397 qh PRINTgood= True;
01398 qh_option ("Pgood", NULL, NULL);
01399 }
01400 }
01401 if (qh DELAUNAY && qh KEEPcoplanar && !qh KEEPinside) {
01402 qh KEEPinside= True;
01403 qh_option ("Qinterior-keep", NULL, NULL);
01404 }
01405 if (qh DELAUNAY && qh HALFspace) {
01406 fprintf (qh ferr, "qhull input error: can not use Delaunay ('d') or Voronoi ('v') with halfspace intersection ('H')\n");
01407 qh_errexit (qh_ERRinput, NULL, NULL);
01408 }
01409 if (!qh DELAUNAY && (qh UPPERdelaunay || qh ATinfinity)) {
01410 fprintf (qh ferr, "qhull input error: use upper-Delaunay ('Qu') or infinity-point ('Qz') with Delaunay ('d') or Voronoi ('v')\n");
01411 qh_errexit (qh_ERRinput, NULL, NULL);
01412 }
01413 if (qh UPPERdelaunay && qh ATinfinity) {
01414 fprintf (qh ferr, "qhull input error: can not use infinity-point ('Qz') with upper-Delaunay ('Qu')\n");
01415 qh_errexit (qh_ERRinput, NULL, NULL);
01416 }
01417 if (qh SCALElast && !qh DELAUNAY && qh PRINTprecision)
01418 fprintf (qh ferr, "qhull input warning: option 'Qbb' (scale-last-coordinate) is normally used with 'd' or 'v'\n");
01419 qh DOcheckmax= (!qh FORCEoutput && !qh SKIPcheckmax && qh MERGING );
01420 qh KEEPnearinside= (qh DOcheckmax && !(qh KEEPinside && qh KEEPcoplanar)
01421 && !qh NOnearinside);
01422 if (qh MERGING)
01423 qh CENTERtype= qh_AScentrum;
01424 else if (qh VORONOI)
01425 qh CENTERtype= qh_ASvoronoi;
01426 if (qh TESTvneighbors && !qh MERGING) {
01427 fprintf(qh ferr, "qhull input error: test vertex neighbors ('Qv') needs a merge option\n");
01428 qh_errexit (qh_ERRinput, NULL ,NULL);
01429 }
01430 if (qh PROJECTinput || (qh DELAUNAY && qh PROJECTdelaunay)) {
01431 qh hull_dim -= qh PROJECTinput;
01432 if (qh DELAUNAY) {
01433 qh hull_dim++;
01434 extra= 1;
01435 }
01436 }
01437 if (qh hull_dim <= 1) {
01438 fprintf(qh ferr, "qhull error: dimension %d must be > 1\n", qh hull_dim);
01439 qh_errexit (qh_ERRinput, NULL, NULL);
01440 }
01441 for (k= 2, factorial=1.0; k < qh hull_dim; k++)
01442 factorial *= k;
01443 qh AREAfactor= 1.0 / factorial;
01444 trace2((qh ferr, "qh_initqhull_globals: initialize globals. dim %d numpoints %d malloc? %d projected %d to hull_dim %d\n",
01445 dim, numpoints, ismalloc, qh PROJECTinput, qh hull_dim));
01446 qh normal_size= qh hull_dim * sizeof(coordT);
01447 qh center_size= qh normal_size - sizeof(coordT);
01448 pointsneeded= qh hull_dim+1;
01449 if (qh hull_dim > qh_DIMmergeVertex) {
01450 qh MERGEvertices= False;
01451 qh_option ("Q3-no-merge-vertices-dim-high", NULL, NULL);
01452 }
01453 if (qh GOODpoint)
01454 pointsneeded++;
01455 #ifdef qh_NOtrace
01456 if (qh IStracing) {
01457 fprintf (qh ferr, "qhull input error: tracing is not installed (qh_NOtrace in user.h)");
01458 qh_errexit (qh_ERRqhull, NULL, NULL);
01459 }
01460 #endif
01461 if (qh RERUN > 1) {
01462 qh TRACElastrun= qh IStracing; /* qh_build_withrestart duplicates next conditional */
01463 if (qh IStracing != -1)
01464 qh IStracing= 0;
01465 }else if (qh TRACEpoint != -1 || qh TRACEdist < REALmax/2 || qh TRACEmerge) {
01466 qh TRACElevel= (qh IStracing? qh IStracing : 3);
01467 qh IStracing= 0;
01468 }
01469 if (qh ROTATErandom == 0 || qh ROTATErandom == -1) {
01470 seed= time (&timedata);
01471 if (qh ROTATErandom == -1) {
01472 seed= -seed;
01473 qh_option ("QRandom-seed", &seed, NULL );
01474 }else
01475 qh_option ("QRotate-random", &seed, NULL);
01476 qh ROTATErandom= seed;
01477 }
01478 seed= qh ROTATErandom;
01479 if (seed == INT_MIN) /* default value */
01480 seed= 1;
01481 else if (seed < 0)
01482 seed= -seed;
01483 qh_RANDOMseed_(seed);
01484 randr= 0.0;
01485 for (i= 1000; i--; ) {
01486 randi= qh_RANDOMint;
01487 randr += randi;
01488 if (randi > qh_RANDOMmax) {
01489 fprintf (qh ferr, "\
01490 qhull configuration error (qh_RANDOMmax in user.h):\n\
01491 random integer %d > qh_RANDOMmax (%.8g)\n",
01492 randi, qh_RANDOMmax);
01493 qh_errexit (qh_ERRinput, NULL, NULL);
01494 }
01495 }
01496 qh_RANDOMseed_(seed);
01497 randr = randr/1000;
01498 if (randr < qh_RANDOMmax/10
01499 || randr > qh_RANDOMmax * 5)
01500 fprintf (qh ferr, "\
01501 qhull configuration warning (qh_RANDOMmax in user.h):\n\
01502 average of 1000 random integers (%.2g) is much different than expected (%.2g).\n\
01503 Is qh_RANDOMmax (%g) wrong?\n",
01504 randr, qh_RANDOMmax/2.0, qh_RANDOMmax);
01505 qh RANDOMa= 2.0 * qh RANDOMfactor/qh_RANDOMmax;
01506 qh RANDOMb= 1.0 - qh RANDOMfactor;
01507 if (qh_HASHfactor < 1.1) {
01508 fprintf(qh ferr, "qhull internal error (qh_initqhull_globals): qh_HASHfactor %d must be at least 1.1. Qhull uses linear hash probing\n",
01509 qh_HASHfactor);
01510 qh_errexit (qh_ERRqhull, NULL, NULL);
01511 }
01512 if (numpoints+extra < pointsneeded) {
01513 fprintf(qh ferr,"qhull input error: not enough points (%d) to construct initial simplex (need %d)\n",
01514 numpoints, pointsneeded);
01515 qh_errexit(qh_ERRinput, NULL, NULL);
01516 }
01517 if (qh PRINTtransparent) {
01518 if (qh hull_dim != 4 || !qh DELAUNAY || qh VORONOI || qh DROPdim >= 0) {
01519 fprintf(qh ferr,"qhull input error: transparent Delaunay ('Gt') needs 3-d Delaunay ('d') w/o 'GDn'\n");
01520 qh_errexit(qh_ERRinput, NULL, NULL);
01521 }
01522 qh DROPdim = 3;
01523 qh PRINTridges = True;
01524 }
01525 for (i= qh_PRINTEND; i--; ) {
01526 if (qh PRINTout[i] == qh_PRINTgeom)
01527 printgeom= True;
01528 else if (qh PRINTout[i] == qh_PRINTmathematica)
01529 printmath= True;
01530 else if (qh PRINTout[i] == qh_PRINTcoplanars)
01531 printcoplanar= True;
01532 else if (qh PRINTout[i] == qh_PRINTpointnearest)
01533 printcoplanar= True;
01534 else if (qh PRINTout[i] == qh_PRINTpointintersect && !qh HALFspace) {
01535 fprintf (qh ferr, "qhull input error: option 'Fp' is only used for \nhalfspace intersection ('Hn,n,n').\n");
01536 qh_errexit (qh_ERRinput, NULL, NULL);
01537 }else if (qh PRINTout[i] == qh_PRINTtriangles && (qh HALFspace || qh VORONOI)) {
01538 fprintf (qh ferr, "qhull input error: option 'Ft' is not available for Voronoi vertices or halfspace intersection\n");
01539 qh_errexit (qh_ERRinput, NULL, NULL);
01540 }else if (qh PRINTout[i] == qh_PRINTcentrums && qh VORONOI) {
01541 fprintf (qh ferr, "qhull input error: option 'FC' is not available for Voronoi vertices ('v')\n");
01542 qh_errexit (qh_ERRinput, NULL, NULL);
01543 }else if (qh PRINTout[i] == qh_PRINTvertices) {
01544 if (qh VORONOI)
01545 qh_option ("Fvoronoi", NULL, NULL);
01546 else
01547 qh_option ("Fvertices", NULL, NULL);
01548 }
01549 }
01550 if (printcoplanar && qh DELAUNAY && qh JOGGLEmax < REALmax/2) {
01551 if (qh PRINTprecision)
01552 fprintf (qh ferr, "qhull input warning: 'QJ' (joggle) will usually prevent coincident input sites for options 'Fc' and 'FP'\n");
01553 }
01554 if (!qh KEEPcoplanar && !qh KEEPinside && !qh ONLYgood) {
01555 if ((qh PRINTcoplanar && qh PRINTspheres) || printcoplanar) {
01556 qh KEEPcoplanar = True;
01557 qh_option ("Qcoplanar", NULL, NULL);
01558 }
01559 }
01560 if (printmath && (qh hull_dim > 3 || qh VORONOI)) {
01561 fprintf (qh ferr, "qhull input error: Mathematica output is only available for 2-d and 3-d convex hulls and 2-d Delaunay triangulations\n");
01562 qh_errexit (qh_ERRinput, NULL, NULL);
01563 }
01564 if (printgeom) {
01565 if (qh hull_dim > 4) {
01566 fprintf (qh ferr, "qhull input error: Geomview output is only available for 2-d, 3-d and 4-d\n");
01567 qh_errexit (qh_ERRinput, NULL, NULL);
01568 }
01569 if (qh PRINTnoplanes && !(qh PRINTcoplanar + qh PRINTcentrums
01570 + qh PRINTdots + qh PRINTspheres + qh DOintersections + qh PRINTridges)) {
01571 fprintf (qh ferr, "qhull input error: no output specified for Geomview\n");
01572 qh_errexit (qh_ERRinput, NULL, NULL);
01573 }
01574 if (qh VORONOI && (qh hull_dim > 3 || qh DROPdim >= 0)) {
01575 fprintf (qh ferr, "qhull input error: Geomview output for Voronoi diagrams only for 2-d\n");
01576 qh_errexit (qh_ERRinput, NULL, NULL);
01577 }
01578 /* can not warn about furthest-site Geomview output: no lower_threshold */
01579 if (qh hull_dim == 4 && qh DROPdim == -1 &&
01580 (qh PRINTcoplanar || qh PRINTspheres || qh PRINTcentrums)) {
01581 fprintf (qh ferr, "qhull input warning: coplanars, vertices, and centrums output not\n\
01582 available for 4-d output (ignored). Could use 'GDn' instead.\n");
01583 qh PRINTcoplanar= qh PRINTspheres= qh PRINTcentrums= False;
01584 }
01585 }
01586 qh PRINTdim= qh hull_dim;
01587 if (qh DROPdim >=0) { /* after Geomview checks */
01588 if (qh DROPdim < qh hull_dim) {
01589 qh PRINTdim--;
01590 if (!printgeom || qh hull_dim < 3)
01591 fprintf (qh ferr, "qhull input warning: drop dimension 'GD%d' is only available for 3-d/4-d Geomview\n", qh DROPdim);
01592 }else
01593 qh DROPdim= -1;
01594 }else if (qh VORONOI) {
01595 qh DROPdim= qh hull_dim-1;
01596 qh PRINTdim= qh hull_dim-1;
01597 }
01598 } /* initqhull_globals */
|
|
|
Definition at line 1618 of file global.c. References i, MERGING, qh, qh_MEMalign, qh_MEMbufsize, qh_MEMinitbuf, qh_meminitbuffers(), qh_memsetup(), qh_memsize(), qh_user_memsizes(), and SETelemsize. Referenced by qh_init_B().
01618 {
01619 int numsizes;
01620 int i;
01621
01622 numsizes= 8+10;
01623 qh_meminitbuffers (qh IStracing, qh_MEMalign, numsizes,
01624 qh_MEMbufsize,qh_MEMinitbuf);
01625 qh_memsize(sizeof(vertexT));
01626 if (qh MERGING) {
01627 qh_memsize(sizeof(ridgeT));
01628 qh_memsize(sizeof(mergeT));
01629 }
01630 qh_memsize(sizeof(facetT));
01631 i= sizeof(setT) + (qh hull_dim - 1) * SETelemsize; /* ridge.vertices */
01632 qh_memsize(i);
01633 qh_memsize(qh normal_size); /* normal */
01634 i += SETelemsize; /* facet.vertices, .ridges, .neighbors */
01635 qh_memsize(i);
01636 qh_user_memsizes();
01637 qh_memsetup();
01638 } /* initqhull_mem */
|
|
||||||||||||||||
|
Definition at line 1650 of file global.c. References fmax_, fout, malloc, qh, qh_CPUclock, qh_ERRmem, qh_initstatistics(), qh_RANDOMseed_, REALmax, and REALmin. Referenced by qh_init_A(), and qh_new_qhull().
01650 {
01651
01652 qh_CPUclock; /* start the clock */
01653 #if qh_QHpointer
01654 if (!(qh_qh= (qhT *)malloc (sizeof(qhT)))) {
01655 fprintf (errfile, "qhull error (qh_initqhull_globals): insufficient memory\n");
01656 exit (qh_ERRmem); /* no error handler */
01657 }
01658 memset((char *)qh_qh, 0, sizeof(qhT)); /* every field is 0, FALSE, NULL */
01659 #else
01660 memset((char *)&qh_qh, 0, sizeof(qhT));
01661 #endif
01662 strcat (qh qhull, "qhull");
01663 qh_initstatistics();
01664 qh ANGLEmerge= True;
01665 qh DROPdim= -1;
01666 qh ferr= errfile;
01667 qh fin= infile;
01668 qh fout= outfile;
01669 qh furthest_id= -1;
01670 qh JOGGLEmax= REALmax;
01671 qh KEEPminArea = REALmax;
01672 qh last_low= REALmax;
01673 qh last_high= REALmax;
01674 qh last_newhigh= REALmax;
01675 qh max_outside= 0.0;
01676 qh max_vertex= 0.0;
01677 qh MAXabs_coord= 0.0;
01678 qh MAXsumcoord= 0.0;
01679 qh MAXwidth= -REALmax;
01680 qh MERGEindependent= True;
01681 qh MINdenom_1= fmax_(1.0/REALmax, REALmin); /* used by qh_scalepoints */
01682 qh MINoutside= 0.0;
01683 qh MINvisible= REALmax;
01684 qh MAXcoplanar= REALmax;
01685 qh outside_err= REALmax;
01686 qh premerge_centrum= 0.0;
01687 qh premerge_cos= REALmax;
01688 qh PRINTprecision= True;
01689 qh PRINTradius= 0.0;
01690 qh postmerge_cos= REALmax;
01691 qh postmerge_centrum= 0.0;
01692 qh ROTATErandom= INT_MIN;
01693 qh MERGEvertices= True;
01694 qh totarea= 0.0;
01695 qh totvol= 0.0;
01696 qh TRACEdist= REALmax;
01697 qh TRACEpoint= -1; /* recompile to trace a point */
01698 qh tracefacet_id= UINT_MAX; /* recompile to trace a facet */
01699 qh tracevertex_id= UINT_MAX; /* recompile to trace a vertex */
01700 qh_RANDOMseed_(1);
01701 } /* initqhull_start */
|
|
|
Definition at line 1723 of file global.c. References key, qh, qh_DEFAULTbox, qh_strtod(), qh_strtol(), REALmax, and realT. Referenced by qh_init_B().
01723 {
01724 realT value;
01725 int index, maxdim, k;
01726 char *s= command;
01727 char key;
01728
01729 maxdim= qh input_dim;
01730 if (qh DELAUNAY && (qh PROJECTdelaunay || qh PROJECTinput))
01731 maxdim++;
01732 while (*s) {
01733 if (*s == '-')
01734 s++;
01735 if (*s == 'P') {
01736 s++;
01737 while (*s && !isspace(key= *s++)) {
01738 if (key == 'd' || key == 'D') {
01739 if (!isdigit(*s)) {
01740 fprintf(qh ferr, "qhull warning: no dimension given for Print option '%c' at: %s. Ignored\n",
01741 key, s-1);
01742 continue;
01743 }
01744 index= qh_strtol (s, &s);
01745 if (index >= qh hull_dim) {
01746 fprintf(qh ferr, "qhull warning: dimension %d for Print option '%c' is >= %d. Ignored\n",
01747 index, key, qh hull_dim);
01748 continue;
01749 }
01750 if (*s == ':') {
01751 s++;
01752 value= qh_strtod(s, &s);
01753 if (fabs((double)value) > 1.0) {
01754 fprintf(qh ferr, "qhull warning: value %2.4g for Print option %c is > +1 or < -1. Ignored\n",
01755 value, key);
01756 continue;
01757 }
01758 }else
01759 value= 0.0;
01760 if (key == 'd')
01761 qh lower_threshold[index]= value;
01762 else
01763 qh upper_threshold[index]= value;
01764 }
01765 }
01766 }else if (*s == 'Q') {
01767 s++;
01768 while (*s && !isspace(key= *s++)) {
01769 if (key == 'b' && *s == 'B') {
01770 s++;
01771 for (k=maxdim; k--; ) {
01772 qh lower_bound[k]= -qh_DEFAULTbox;
01773 qh upper_bound[k]= qh_DEFAULTbox;
01774 }
01775 }else if (key == 'b' && *s == 'b')
01776 s++;
01777 else if (key == 'b' || key == 'B') {
01778 if (!isdigit(*s)) {
01779 fprintf(qh ferr, "qhull warning: no dimension given for Qhull option %c. Ignored\n",
01780 key);
01781 continue;
01782 }
01783 index= qh_strtol (s, &s);
01784 if (index >= maxdim) {
01785 fprintf(qh ferr, "qhull warning: dimension %d for Qhull option %c is >= %d. Ignored\n",
01786 index, key, maxdim);
01787 continue;
01788 }
01789 if (*s == ':') {
01790 s++;
01791 value= qh_strtod(s, &s);
01792 }else if (key == 'b')
01793 value= -qh_DEFAULTbox;
01794 else
01795 value= qh_DEFAULTbox;
01796 if (key == 'b')
01797 qh lower_bound[index]= value;
01798 else
01799 qh upper_bound[index]= value;
01800 }
01801 }
01802 }else {
01803 while (*s && !isspace (*s))
01804 s++;
01805 }
01806 while (isspace (*s))
01807 s++;
01808 }
01809 for (k= qh hull_dim; k--; ) {
01810 if (qh lower_threshold[k] > -REALmax/2) {
01811 qh GOODthreshold= True;
01812 if (qh upper_threshold[k] < REALmax/2) {
01813 qh SPLITthresholds= True;
01814 qh GOODthreshold= False;
01815 break;
01816 }
01817 }else if (qh upper_threshold[k] < REALmax/2)
01818 qh GOODthreshold= True;
01819 }
01820 } /* initthresholds */
|
|
||||||||||||||||
|
Definition at line 1832 of file global.c. References i, maximize_, qh, r, and realT. Referenced by qh_build_withrestart(), qh_detroundoff(), qh_initflags(), qh_initialhull(), qh_initqhull_globals(), and qh_joggleinput().
01832 {
01833 char buf[200];
01834 int len, maxlen;
01835
01836 sprintf (buf, " %s", option);
01837 if (i)
01838 sprintf (buf+strlen(buf), " %d", *i);
01839 if (r)
01840 sprintf (buf+strlen(buf), " %2.2g", *r);
01841 len= strlen(buf);
01842 qh qhull_optionlen += len;
01843 maxlen= sizeof (qh qhull_options) - len -1;
01844 maximize_(maxlen, 0);
01845 if (qh qhull_optionlen >= 80 && maxlen > 0) {
01846 qh qhull_optionlen= len;
01847 strncat (qh qhull_options, "\n", maxlen--);
01848 }
01849 strncat (qh qhull_options, buf, maxlen);
01850 } /* option */
|
|
||||||||||||
|
Definition at line 1942 of file global.c. References strtod(). Referenced by qh_checkflags(), qh_initflags(), qh_initthresholds(), qh_readfeasible(), qh_readpoints(), and qh_setfeasible().
01942 {
01943 double result;
01944
01945 result= strtod (s, endp);
01946 if (s < (*endp) && (*endp)[-1] == ' ')
01947 (*endp)--;
01948 return result;
01949 } /* strtod */
|
|
||||||||||||
|
Definition at line 1951 of file global.c. Referenced by qh_initflags(), qh_initthresholds(), and qh_readpoints().
01951 {
01952 int result;
01953
01954 result= (int) strtol (s, endp, 10);
01955 if (s< (*endp) && (*endp)[-1] == ' ')
01956 (*endp)--;
01957 return result;
01958 } /* strtol */
|
Variable Documentation
|
|
|