| 1 |
/*<html><pre> -<a href="qh-qhull.htm" |
| 2 |
>-------------------------------</a><a name="TOP">-</a> |
| 3 |
|
| 4 |
qhull.c |
| 5 |
Quickhull algorithm for convex hulls |
| 6 |
|
| 7 |
qhull() and top-level routines |
| 8 |
|
| 9 |
see qh-qhull.htm, qhull.h, unix.c |
| 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 |
/*============= functions in alphabetic order after qhull() =======*/ |
| 19 |
|
| 20 |
/*-<a href="qh-qhull.htm#TOC" |
| 21 |
>-------------------------------</a><a name="qhull">-</a> |
| 22 |
|
| 23 |
qh_qhull() |
| 24 |
compute DIM3 convex hull of qh.num_points starting at qh.first_point |
| 25 |
qh contains all global options and variables |
| 26 |
|
| 27 |
returns: |
| 28 |
returns polyhedron |
| 29 |
qh.facet_list, qh.num_facets, qh.vertex_list, qh.num_vertices, |
| 30 |
|
| 31 |
returns global variables |
| 32 |
qh.hulltime, qh.max_outside, qh.interior_point, qh.max_vertex, qh.min_vertex |
| 33 |
|
| 34 |
returns precision constants |
| 35 |
qh.ANGLEround, centrum_radius, cos_max, DISTround, MAXabs_coord, ONEmerge |
| 36 |
|
| 37 |
notes: |
| 38 |
unless needed for output |
| 39 |
qh.max_vertex and qh.min_vertex are max/min due to merges |
| 40 |
|
| 41 |
see: |
| 42 |
to add individual points to either qh.num_points |
| 43 |
please use qh_addpoint() |
| 44 |
|
| 45 |
if qh.GETarea |
| 46 |
qh_produceoutput() returns qh.totarea and qh.totvol via qh_getarea() |
| 47 |
|
| 48 |
design: |
| 49 |
record starting time |
| 50 |
initialize hull and partition points |
| 51 |
build convex hull |
| 52 |
unless early termination |
| 53 |
update facet->maxoutside for vertices, coplanar, and near-inside points |
| 54 |
error if temporary sets exist |
| 55 |
record end time |
| 56 |
*/ |
| 57 |
|
| 58 |
void qh_qhull (void) { |
| 59 |
int numoutside; |
| 60 |
|
| 61 |
qh hulltime= qh_CPUclock; |
| 62 |
if (qh RERUN || qh JOGGLEmax < REALmax/2) |
| 63 |
qh_build_withrestart(); |
| 64 |
else { |
| 65 |
qh_initbuild(); |
| 66 |
qh_buildhull(); |
| 67 |
} |
| 68 |
if (!qh STOPpoint && !qh STOPcone) { |
| 69 |
if (qh ZEROall_ok && !qh TESTvneighbors && qh MERGEexact) |
| 70 |
qh_checkzero( qh_ALL); |
| 71 |
if (qh ZEROall_ok && !qh TESTvneighbors && !qh WAScoplanar) { |
| 72 |
trace2((qh ferr, "qh_qhull: all facets are clearly convex and no coplanar points. Post-merging and check of maxout not needed.\n")); |
| 73 |
qh DOcheckmax= False; |
| 74 |
}else { |
| 75 |
if (qh MERGEexact || (qh hull_dim > qh_DIMreduceBuild && qh PREmerge)) |
| 76 |
qh_postmerge ("First post-merge", qh premerge_centrum, qh premerge_cos, |
| 77 |
(qh POSTmerge ? False : qh TESTvneighbors)); |
| 78 |
else if (!qh POSTmerge && qh TESTvneighbors) |
| 79 |
qh_postmerge ("For testing vertex neighbors", qh premerge_centrum, |
| 80 |
qh premerge_cos, True); |
| 81 |
if (qh POSTmerge) |
| 82 |
qh_postmerge ("For post-merging", qh postmerge_centrum, |
| 83 |
qh postmerge_cos, qh TESTvneighbors); |
| 84 |
if (qh visible_list == qh facet_list) { /* i.e., merging done */ |
| 85 |
qh findbestnew= True; |
| 86 |
qh_partitionvisible (/*visible_list, newfacet_list*/ !qh_ALL, &numoutside); |
| 87 |
qh findbestnew= False; |
| 88 |
qh_deletevisible (/*qh visible_list*/); |
| 89 |
qh_resetlists (False, qh_RESETvisible /*qh visible_list newvertex_list newfacet_list */); |
| 90 |
} |
| 91 |
} |
| 92 |
if (qh DOcheckmax){ |
| 93 |
if (qh REPORTfreq) { |
| 94 |
qh_buildtracing (NULL, NULL); |
| 95 |
fprintf (qh ferr, "\nTesting all coplanar points.\n"); |
| 96 |
} |
| 97 |
qh_check_maxout(); |
| 98 |
} |
| 99 |
if (qh KEEPnearinside && !qh maxoutdone) |
| 100 |
qh_nearcoplanar(); |
| 101 |
} |
| 102 |
if (qh_setsize ((setT*)qhmem.tempstack) != 0) { |
| 103 |
fprintf (qh ferr, "qhull internal error (qh_qhull): temporary sets not empty (%d)\n", |
| 104 |
qh_setsize ((setT*)qhmem.tempstack)); |
| 105 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
| 106 |
} |
| 107 |
qh hulltime= qh_CPUclock - qh hulltime; |
| 108 |
qh QHULLfinished= True; |
| 109 |
trace1((qh ferr, "qh_qhull: algorithm completed\n")); |
| 110 |
} /* qhull */ |
| 111 |
|
| 112 |
/*-<a href="qh-qhull.htm#TOC" |
| 113 |
>-------------------------------</a><a name="addpoint">-</a> |
| 114 |
|
| 115 |
qh_addpoint( furthest, facet, checkdist ) |
| 116 |
add point (usually furthest point) above facet to hull |
| 117 |
if checkdist, |
| 118 |
check that point is above facet. |
| 119 |
if point is not outside of the hull, uses qh_partitioncoplanar() |
| 120 |
assumes that facet is defined by qh_findbestfacet() |
| 121 |
else if facet specified, |
| 122 |
assumes that point is above facet (major damage if below) |
| 123 |
for Delaunay triangulations, |
| 124 |
Use qh_setdelaunay() to lift point to paraboloid and scale by 'Qbb' if needed |
| 125 |
Do not use options 'Qbk', 'QBk', or 'QbB' since they scale the coordinates. |
| 126 |
|
| 127 |
returns: |
| 128 |
returns False if user requested an early termination |
| 129 |
qh.visible_list, newfacet_list, delvertex_list, NEWfacets may be defined |
| 130 |
updates qh.facet_list, qh.num_facets, qh.vertex_list, qh.num_vertices |
| 131 |
clear qh.maxoutdone (will need to call qh_check_maxout() for facet->maxoutside) |
| 132 |
if unknown point, adds a pointer to qh.other_points |
| 133 |
do not deallocate the point's coordinates |
| 134 |
|
| 135 |
notes: |
| 136 |
assumes point is near its best facet and not at a local minimum of a lens |
| 137 |
distributions. Use qh_findbestfacet to avoid this case. |
| 138 |
uses qh.visible_list, qh.newfacet_list, qh.delvertex_list, qh.NEWfacets |
| 139 |
|
| 140 |
see also: |
| 141 |
qh_triangulate() -- triangulate non-simplicial facets |
| 142 |
|
| 143 |
design: |
| 144 |
check point in qh.first_point/.num_points |
| 145 |
if checkdist |
| 146 |
if point not above facet |
| 147 |
partition coplanar point |
| 148 |
exit |
| 149 |
exit if pre STOPpoint requested |
| 150 |
find horizon and visible facets for point |
| 151 |
make new facets for point to horizon |
| 152 |
make hyperplanes for point |
| 153 |
compute balance statistics |
| 154 |
match neighboring new facets |
| 155 |
update vertex neighbors and delete interior vertices |
| 156 |
exit if STOPcone requested |
| 157 |
merge non-convex new facets |
| 158 |
if merge found, many merges, or 'Qf' |
| 159 |
please use qh_findbestnew() instead of qh_findbest() |
| 160 |
partition outside points from visible facets |
| 161 |
delete visible facets |
| 162 |
check polyhedron if requested |
| 163 |
exit if post STOPpoint requested |
| 164 |
reset working lists of facets and vertices |
| 165 |
*/ |
| 166 |
boolT qh_addpoint (pointT *furthest, facetT *facet, boolT checkdist) { |
| 167 |
int goodvisible, goodhorizon; |
| 168 |
vertexT *vertex; |
| 169 |
facetT *newfacet; |
| 170 |
realT dist, newbalance, pbalance; |
| 171 |
boolT isoutside= False; |
| 172 |
int numpart, numpoints, numnew, firstnew; |
| 173 |
|
| 174 |
qh maxoutdone= False; |
| 175 |
if (qh_pointid (furthest) == -1) |
| 176 |
qh_setappend (&qh other_points, furthest); |
| 177 |
if (!facet) { |
| 178 |
fprintf (qh ferr, "qh_addpoint: NULL facet. Need to call qh_findbestfacet first\n"); |
| 179 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
| 180 |
} |
| 181 |
if (checkdist) { |
| 182 |
facet= qh_findbest (furthest, facet, !qh_ALL, !qh_ISnewfacets, !qh_NOupper, |
| 183 |
&dist, &isoutside, &numpart); |
| 184 |
zzadd_(Zpartition, numpart); |
| 185 |
if (!isoutside) { |
| 186 |
zinc_(Znotmax); /* last point of outsideset is no longer furthest. */ |
| 187 |
facet->notfurthest= True; |
| 188 |
qh_partitioncoplanar (furthest, facet, &dist); |
| 189 |
return True; |
| 190 |
} |
| 191 |
} |
| 192 |
qh_buildtracing (furthest, facet); |
| 193 |
if (qh STOPpoint < 0 && qh furthest_id == -qh STOPpoint-1) { |
| 194 |
facet->notfurthest= True; |
| 195 |
return False; |
| 196 |
} |
| 197 |
qh_findhorizon (furthest, facet, &goodvisible, &goodhorizon); |
| 198 |
if (qh ONLYgood && !(goodvisible+goodhorizon) && !qh GOODclosest) { |
| 199 |
zinc_(Znotgood); |
| 200 |
facet->notfurthest= True; |
| 201 |
/* last point of outsideset is no longer furthest. This is ok |
| 202 |
since all points of the outside are likely to be bad */ |
| 203 |
qh_resetlists (False, qh_RESETvisible /*qh visible_list newvertex_list newfacet_list */); |
| 204 |
return True; |
| 205 |
} |
| 206 |
zzinc_(Zprocessed); |
| 207 |
firstnew= qh facet_id; |
| 208 |
vertex= qh_makenewfacets (furthest /*visible_list, attaches if !ONLYgood */); |
| 209 |
qh_makenewplanes (/* newfacet_list */); |
| 210 |
numnew= qh facet_id - firstnew; |
| 211 |
newbalance= numnew - (realT) (qh num_facets-qh num_visible) |
| 212 |
* qh hull_dim/qh num_vertices; |
| 213 |
wadd_(Wnewbalance, newbalance); |
| 214 |
wadd_(Wnewbalance2, newbalance * newbalance); |
| 215 |
if (qh ONLYgood |
| 216 |
&& !qh_findgood (qh newfacet_list, goodhorizon) && !qh GOODclosest) { |
| 217 |
FORALLnew_facets |
| 218 |
qh_delfacet (newfacet); |
| 219 |
qh_delvertex (vertex); |
| 220 |
qh_resetlists (True, qh_RESETvisible /*qh visible_list newvertex_list newfacet_list */); |
| 221 |
zinc_(Znotgoodnew); |
| 222 |
facet->notfurthest= True; |
| 223 |
return True; |
| 224 |
} |
| 225 |
if (qh ONLYgood) |
| 226 |
qh_attachnewfacets(/*visible_list*/); |
| 227 |
qh_matchnewfacets(); |
| 228 |
qh_updatevertices(); |
| 229 |
if (qh STOPcone && qh furthest_id == qh STOPcone-1) { |
| 230 |
facet->notfurthest= True; |
| 231 |
return False; /* visible_list etc. still defined */ |
| 232 |
} |
| 233 |
qh findbestnew= False; |
| 234 |
if (qh PREmerge || qh MERGEexact) { |
| 235 |
qh_premerge (vertex, qh premerge_centrum, qh premerge_cos); |
| 236 |
if (qh_USEfindbestnew) |
| 237 |
qh findbestnew= True; |
| 238 |
else { |
| 239 |
FORALLnew_facets { |
| 240 |
if (!newfacet->simplicial) { |
| 241 |
qh findbestnew= True; /* use qh_findbestnew instead of qh_findbest*/ |
| 242 |
break; |
| 243 |
} |
| 244 |
} |
| 245 |
} |
| 246 |
}else if (qh BESToutside) |
| 247 |
qh findbestnew= True; |
| 248 |
qh_partitionvisible (/*visible_list, newfacet_list*/ !qh_ALL, &numpoints); |
| 249 |
qh findbestnew= False; |
| 250 |
qh findbest_notsharp= False; |
| 251 |
zinc_(Zpbalance); |
| 252 |
pbalance= numpoints - (realT) qh hull_dim /* assumes all points extreme */ |
| 253 |
* (qh num_points - qh num_vertices)/qh num_vertices; |
| 254 |
wadd_(Wpbalance, pbalance); |
| 255 |
wadd_(Wpbalance2, pbalance * pbalance); |
| 256 |
qh_deletevisible (/*qh visible_list*/); |
| 257 |
zmax_(Zmaxvertex, qh num_vertices); |
| 258 |
qh NEWfacets= False; |
| 259 |
if (qh IStracing >= 4) { |
| 260 |
if (qh num_facets < 2000) |
| 261 |
qh_printlists(); |
| 262 |
qh_printfacetlist (qh newfacet_list, NULL, True); |
| 263 |
qh_checkpolygon (qh facet_list); |
| 264 |
}else if (qh CHECKfrequently) { |
| 265 |
if (qh num_facets < 50) |
| 266 |
qh_checkpolygon (qh facet_list); |
| 267 |
else |
| 268 |
qh_checkpolygon (qh newfacet_list); |
| 269 |
} |
| 270 |
if (qh STOPpoint > 0 && qh furthest_id == qh STOPpoint-1) |
| 271 |
return False; |
| 272 |
qh_resetlists (True, qh_RESETvisible /*qh visible_list newvertex_list newfacet_list */); |
| 273 |
/* qh_triangulate(); to test qh.TRInormals */ |
| 274 |
trace2((qh ferr, "qh_addpoint: added p%d new facets %d new balance %2.2g point balance %2.2g\n", qh_pointid (furthest), numnew, newbalance, pbalance)); |
| 275 |
return True; |
| 276 |
} /* addpoint */ |
| 277 |
|
| 278 |
/*-<a href="qh-qhull.htm#TOC" |
| 279 |
>-------------------------------</a><a name="build_withrestart">-</a> |
| 280 |
|
| 281 |
qh_build_withrestart() |
| 282 |
allow restarts due to qh.JOGGLEmax while calling qh_buildhull() |
| 283 |
qh.FIRSTpoint/qh.NUMpoints is point array |
| 284 |
it may be moved by qh_joggleinput() |
| 285 |
*/ |
| 286 |
void qh_build_withrestart (void) { |
| 287 |
int restart; |
| 288 |
|
| 289 |
qh ALLOWrestart= True; |
| 290 |
while (True) { |
| 291 |
restart= setjmp (qh restartexit); /* simple statement for CRAY J916 */ |
| 292 |
if (restart) { /* only from qh_precision() */ |
| 293 |
zzinc_(Zretry); |
| 294 |
wmax_(Wretrymax, qh JOGGLEmax); |
| 295 |
qh ERREXITcalled= False; |
| 296 |
qh STOPcone= True; /* if break, prevents normal output */ |
| 297 |
} |
| 298 |
if (!qh RERUN && qh JOGGLEmax < REALmax/2) { |
| 299 |
if (qh build_cnt > qh_JOGGLEmaxretry) { |
| 300 |
fprintf(qh ferr, "\n\ |
| 301 |
qhull precision error: %d attempts to construct a convex hull\n\ |
| 302 |
with joggled input. Increase joggle above 'QJ%2.2g'\n\ |
| 303 |
or modify qh_JOGGLE... parameters in user.h\n", |
| 304 |
qh build_cnt, qh JOGGLEmax); |
| 305 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
| 306 |
} |
| 307 |
if (qh build_cnt && !restart) |
| 308 |
break; |
| 309 |
}else if (qh build_cnt && qh build_cnt >= qh RERUN) |
| 310 |
break; |
| 311 |
qh STOPcone= False; |
| 312 |
qh_freebuild (True); /* first call is a nop */ |
| 313 |
qh build_cnt++; |
| 314 |
if (!qh qhull_optionsiz) |
| 315 |
qh qhull_optionsiz= strlen (qh qhull_options); |
| 316 |
else { |
| 317 |
qh qhull_options [qh qhull_optionsiz]= '\0'; |
| 318 |
qh qhull_optionlen= 80; |
| 319 |
} |
| 320 |
qh_option("_run", &qh build_cnt, NULL); |
| 321 |
if (qh build_cnt == qh RERUN) { |
| 322 |
qh IStracing= qh TRACElastrun; /* duplicated from qh_initqhull_globals */ |
| 323 |
if (qh TRACEpoint != -1 || qh TRACEdist < REALmax/2 || qh TRACEmerge) { |
| 324 |
qh TRACElevel= (qh IStracing? qh IStracing : 3); |
| 325 |
qh IStracing= 0; |
| 326 |
} |
| 327 |
qhmem.IStracing= qh IStracing; |
| 328 |
} |
| 329 |
if (qh JOGGLEmax < REALmax/2) |
| 330 |
qh_joggleinput(); |
| 331 |
qh_initbuild(); |
| 332 |
qh_buildhull(); |
| 333 |
if (qh JOGGLEmax < REALmax/2 && !qh MERGING) |
| 334 |
qh_checkconvex (qh facet_list, qh_ALGORITHMfault); |
| 335 |
} |
| 336 |
qh ALLOWrestart= False; |
| 337 |
} /* qh_build_withrestart */ |
| 338 |
|
| 339 |
/*-<a href="qh-qhull.htm#TOC" |
| 340 |
>-------------------------------</a><a name="buildhull">-</a> |
| 341 |
|
| 342 |
qh_buildhull() |
| 343 |
construct a convex hull by adding outside points one at a time |
| 344 |
|
| 345 |
returns: |
| 346 |
|
| 347 |
notes: |
| 348 |
may be called multiple times |
| 349 |
checks facet and vertex lists for incorrect flags |
| 350 |
to recover from STOPcone, call qh_deletevisible and qh_resetlists |
| 351 |
|
| 352 |
design: |
| 353 |
check visible facet and newfacet flags |
| 354 |
check newlist vertex flags and qh.STOPcone/STOPpoint |
| 355 |
for each facet with a furthest outside point |
| 356 |
add point to facet |
| 357 |
exit if qh.STOPcone or qh.STOPpoint requested |
| 358 |
if qh.NARROWhull for initial simplex |
| 359 |
partition remaining outside points to coplanar sets |
| 360 |
*/ |
| 361 |
void qh_buildhull(void) { |
| 362 |
facetT *facet; |
| 363 |
pointT *furthest; |
| 364 |
vertexT *vertex; |
| 365 |
int id; |
| 366 |
|
| 367 |
trace1((qh ferr, "qh_buildhull: start build hull\n")); |
| 368 |
FORALLfacets { |
| 369 |
if (facet->visible || facet->newfacet) { |
| 370 |
fprintf (qh ferr, "qhull internal error (qh_buildhull): visible or new facet f%d in facet list\n", |
| 371 |
facet->id); |
| 372 |
qh_errexit (qh_ERRqhull, facet, NULL); |
| 373 |
} |
| 374 |
} |
| 375 |
FORALLvertices { |
| 376 |
if (vertex->newlist) { |
| 377 |
fprintf (qh ferr, "qhull internal error (qh_buildhull): new vertex f%d in vertex list\n", |
| 378 |
vertex->id); |
| 379 |
qh_errprint ("ERRONEOUS", NULL, NULL, NULL, vertex); |
| 380 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
| 381 |
} |
| 382 |
id= qh_pointid (vertex->point); |
| 383 |
if ((qh STOPpoint>0 && id == qh STOPpoint-1) || |
| 384 |
(qh STOPpoint<0 && id == -qh STOPpoint-1) || |
| 385 |
(qh STOPcone>0 && id == qh STOPcone-1)) { |
| 386 |
trace1((qh ferr,"qh_buildhull: stop point or cone P%d in initial hull\n", id)); |
| 387 |
return; |
| 388 |
} |
| 389 |
} |
| 390 |
qh facet_next= qh facet_list; /* advance facet when processed */ |
| 391 |
while ((furthest= qh_nextfurthest (&facet))) { |
| 392 |
qh num_outside--; /* if ONLYmax, furthest may not be outside */ |
| 393 |
if (!qh_addpoint (furthest, facet, qh ONLYmax)) |
| 394 |
break; |
| 395 |
} |
| 396 |
if (qh NARROWhull) /* move points from outsideset to coplanarset */ |
| 397 |
qh_outcoplanar( /* facet_list */ ); |
| 398 |
if (qh num_outside && !furthest) { |
| 399 |
fprintf (qh ferr, "qhull internal error (qh_buildhull): %d outside points were never processed.\n", qh num_outside); |
| 400 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
| 401 |
} |
| 402 |
trace1((qh ferr, "qh_buildhull: completed the hull construction\n")); |
| 403 |
} /* buildhull */ |
| 404 |
|
| 405 |
|
| 406 |
/*-<a href="qh-qhull.htm#TOC" |
| 407 |
>-------------------------------</a><a name="buildtracing">-</a> |
| 408 |
|
| 409 |
qh_buildtracing( furthest, facet ) |
| 410 |
trace an iteration of qh_buildhull() for furthest point and facet |
| 411 |
if !furthest, prints progress message |
| 412 |
|
| 413 |
returns: |
| 414 |
tracks progress with qh.lastreport |
| 415 |
updates qh.furthest_id (-3 if furthest is NULL) |
| 416 |
also resets visit_id, vertext_visit on wrap around |
| 417 |
|
| 418 |
see: |
| 419 |
qh_tracemerging() |
| 420 |
|
| 421 |
design: |
| 422 |
if !furthest |
| 423 |
print progress message |
| 424 |
exit |
| 425 |
if 'TFn' iteration |
| 426 |
print progress message |
| 427 |
else if tracing |
| 428 |
trace furthest point and facet |
| 429 |
reset qh.visit_id and qh.vertex_visit if overflow may occur |
| 430 |
set qh.furthest_id for tracing |
| 431 |
*/ |
| 432 |
void qh_buildtracing (pointT *furthest, facetT *facet) { |
| 433 |
realT dist= 0; |
| 434 |
float cpu; |
| 435 |
int total, furthestid; |
| 436 |
time_t timedata; |
| 437 |
struct tm *tp; |
| 438 |
vertexT *vertex; |
| 439 |
|
| 440 |
qh old_randomdist= qh RANDOMdist; |
| 441 |
qh RANDOMdist= False; |
| 442 |
if (!furthest) { |
| 443 |
time (&timedata); |
| 444 |
tp= localtime (&timedata); |
| 445 |
cpu= qh_CPUclock - qh hulltime; |
| 446 |
cpu /= qh_SECticks; |
| 447 |
total= zzval_(Ztotmerge) - zzval_(Zcyclehorizon) + zzval_(Zcyclefacettot); |
| 448 |
fprintf (qh ferr, "\n\ |
| 449 |
At %02d:%02d:%02d & %2.5g CPU secs, qhull has created %d facets and merged %d.\n\ |
| 450 |
The current hull contains %d facets and %d vertices. Last point was p%d\n", |
| 451 |
tp->tm_hour, tp->tm_min, tp->tm_sec, cpu, qh facet_id -1, |
| 452 |
total, qh num_facets, qh num_vertices, qh furthest_id); |
| 453 |
return; |
| 454 |
} |
| 455 |
furthestid= qh_pointid (furthest); |
| 456 |
if (qh TRACEpoint == furthestid) { |
| 457 |
qh IStracing= qh TRACElevel; |
| 458 |
qhmem.IStracing= qh TRACElevel; |
| 459 |
}else if (qh TRACEpoint != -1 && qh TRACEdist < REALmax/2) { |
| 460 |
qh IStracing= 0; |
| 461 |
qhmem.IStracing= 0; |
| 462 |
} |
| 463 |
if (qh REPORTfreq && (qh facet_id-1 > qh lastreport+qh REPORTfreq)) { |
| 464 |
qh lastreport= qh facet_id-1; |
| 465 |
time (&timedata); |
| 466 |
tp= localtime (&timedata); |
| 467 |
cpu= qh_CPUclock - qh hulltime; |
| 468 |
cpu /= qh_SECticks; |
| 469 |
total= zzval_(Ztotmerge) - zzval_(Zcyclehorizon) + zzval_(Zcyclefacettot); |
| 470 |
zinc_(Zdistio); |
| 471 |
qh_distplane (furthest, facet, &dist); |
| 472 |
fprintf (qh ferr, "\n\ |
| 473 |
At %02d:%02d:%02d & %2.5g CPU secs, qhull has created %d facets and merged %d.\n\ |
| 474 |
The current hull contains %d facets and %d vertices. There are %d\n\ |
| 475 |
outside points. Next is point p%d (v%d), %2.2g above f%d.\n", |
| 476 |
tp->tm_hour, tp->tm_min, tp->tm_sec, cpu, qh facet_id -1, |
| 477 |
total, qh num_facets, qh num_vertices, qh num_outside+1, |
| 478 |
furthestid, qh vertex_id, dist, getid_(facet)); |
| 479 |
}else if (qh IStracing >=1) { |
| 480 |
cpu= qh_CPUclock - qh hulltime; |
| 481 |
cpu /= qh_SECticks; |
| 482 |
qh_distplane (furthest, facet, &dist); |
| 483 |
fprintf (qh ferr, "qh_addpoint: add p%d (v%d) to hull of %d facets (%2.2g above f%d) and %d outside at %4.4g CPU secs. Previous was p%d.\n", |
| 484 |
furthestid, qh vertex_id, qh num_facets, dist, |
| 485 |
getid_(facet), qh num_outside+1, cpu, qh furthest_id); |
| 486 |
} |
| 487 |
if (qh visit_id > (unsigned) INT_MAX) { |
| 488 |
qh visit_id= 0; |
| 489 |
FORALLfacets |
| 490 |
facet->visitid= qh visit_id; |
| 491 |
} |
| 492 |
if (qh vertex_visit > (unsigned) INT_MAX) { |
| 493 |
qh vertex_visit= 0; |
| 494 |
FORALLvertices |
| 495 |
vertex->visitid= qh vertex_visit; |
| 496 |
} |
| 497 |
qh furthest_id= furthestid; |
| 498 |
qh RANDOMdist= qh old_randomdist; |
| 499 |
} /* buildtracing */ |
| 500 |
|
| 501 |
/*-<a href="qh-qhull.htm#TOC" |
| 502 |
>-------------------------------</a><a name="errexit2">-</a> |
| 503 |
|
| 504 |
qh_errexit2( exitcode, facet, otherfacet ) |
| 505 |
return exitcode to system after an error |
| 506 |
report two facets |
| 507 |
|
| 508 |
returns: |
| 509 |
assumes exitcode non-zero |
| 510 |
|
| 511 |
see: |
| 512 |
normally use qh_errexit() in user.c (reports a facet and a ridge) |
| 513 |
*/ |
| 514 |
void qh_errexit2(int exitcode, facetT *facet, facetT *otherfacet) { |
| 515 |
|
| 516 |
qh_errprint("ERRONEOUS", facet, otherfacet, NULL, NULL); |
| 517 |
qh_errexit (exitcode, NULL, NULL); |
| 518 |
} /* errexit2 */ |
| 519 |
|
| 520 |
|
| 521 |
/*-<a href="qh-qhull.htm#TOC" |
| 522 |
>-------------------------------</a><a name="findhorizon">-</a> |
| 523 |
|
| 524 |
qh_findhorizon( point, facet, goodvisible, goodhorizon ) |
| 525 |
given a visible facet, find the point's horizon and visible facets |
| 526 |
for all facets, !facet-visible |
| 527 |
|
| 528 |
returns: |
| 529 |
returns qh.visible_list/num_visible with all visible facets |
| 530 |
marks visible facets with ->visible |
| 531 |
updates count of good visible and good horizon facets |
| 532 |
updates qh.max_outside, qh.max_vertex, facet->maxoutside |
| 533 |
|
| 534 |
see: |
| 535 |
similar to qh_delpoint() |
| 536 |
|
| 537 |
design: |
| 538 |
move facet to qh.visible_list at end of qh.facet_list |
| 539 |
for all visible facets |
| 540 |
for each unvisited neighbor of a visible facet |
| 541 |
compute distance of point to neighbor |
| 542 |
if point above neighbor |
| 543 |
move neighbor to end of qh.visible_list |
| 544 |
else if point is coplanar with neighbor |
| 545 |
update qh.max_outside, qh.max_vertex, neighbor->maxoutside |
| 546 |
mark neighbor coplanar (will create a samecycle later) |
| 547 |
update horizon statistics |
| 548 |
*/ |
| 549 |
void qh_findhorizon(pointT *point, facetT *facet, int *goodvisible, int *goodhorizon) { |
| 550 |
facetT *neighbor, **neighborp, *visible; |
| 551 |
int numhorizon= 0, coplanar= 0; |
| 552 |
realT dist; |
| 553 |
|
| 554 |
trace1((qh ferr,"qh_findhorizon: find horizon for point p%d facet f%d\n",qh_pointid(point),facet->id)); |
| 555 |
*goodvisible= *goodhorizon= 0; |
| 556 |
zinc_(Ztotvisible); |
| 557 |
qh_removefacet(facet); /* visible_list at end of qh facet_list */ |
| 558 |
qh_appendfacet(facet); |
| 559 |
qh num_visible= 1; |
| 560 |
if (facet->good) |
| 561 |
(*goodvisible)++; |
| 562 |
qh visible_list= facet; |
| 563 |
facet->visible= True; |
| 564 |
facet->f.replace= NULL; |
| 565 |
if (qh IStracing >=4) |
| 566 |
qh_errprint ("visible", facet, NULL, NULL, NULL); |
| 567 |
qh visit_id++; |
| 568 |
FORALLvisible_facets { |
| 569 |
if (visible->tricoplanar && !qh TRInormals) { |
| 570 |
fprintf (qh ferr, "qh_findhorizon: does not work for tricoplanar facets. Use option 'Q11'\n"); |
| 571 |
qh_errexit (qh_ERRqhull, visible, NULL); |
| 572 |
} |
| 573 |
visible->visitid= qh visit_id; |
| 574 |
FOREACHneighbor_(visible) { |
| 575 |
if (neighbor->visitid == qh visit_id) |
| 576 |
continue; |
| 577 |
neighbor->visitid= qh visit_id; |
| 578 |
zzinc_(Znumvisibility); |
| 579 |
qh_distplane(point, neighbor, &dist); |
| 580 |
if (dist > qh MINvisible) { |
| 581 |
zinc_(Ztotvisible); |
| 582 |
qh_removefacet(neighbor); /* append to end of qh visible_list */ |
| 583 |
qh_appendfacet(neighbor); |
| 584 |
neighbor->visible= True; |
| 585 |
neighbor->f.replace= NULL; |
| 586 |
qh num_visible++; |
| 587 |
if (neighbor->good) |
| 588 |
(*goodvisible)++; |
| 589 |
if (qh IStracing >=4) |
| 590 |
qh_errprint ("visible", neighbor, NULL, NULL, NULL); |
| 591 |
}else { |
| 592 |
if (dist > - qh MAXcoplanar) { |
| 593 |
neighbor->coplanar= True; |
| 594 |
zzinc_(Zcoplanarhorizon); |
| 595 |
qh_precision ("coplanar horizon"); |
| 596 |
coplanar++; |
| 597 |
if (qh MERGING) { |
| 598 |
if (dist > 0) { |
| 599 |
maximize_(qh max_outside, dist); |
| 600 |
maximize_(qh max_vertex, dist); |
| 601 |
#if qh_MAXoutside |
| 602 |
maximize_(neighbor->maxoutside, dist); |
| 603 |
#endif |
| 604 |
}else |
| 605 |
minimize_(qh min_vertex, dist); /* due to merge later */ |
| 606 |
} |
| 607 |
trace2((qh ferr, "qh_findhorizon: point p%d is coplanar to horizon f%d, dist=%2.7g < qh MINvisible (%2.7g)\n", qh_pointid(point), neighbor->id, dist, qh MINvisible)); |
| 608 |
}else |
| 609 |
neighbor->coplanar= False; |
| 610 |
zinc_(Ztothorizon); |
| 611 |
numhorizon++; |
| 612 |
if (neighbor->good) |
| 613 |
(*goodhorizon)++; |
| 614 |
if (qh IStracing >=4) |
| 615 |
qh_errprint ("horizon", neighbor, NULL, NULL, NULL); |
| 616 |
} |
| 617 |
} |
| 618 |
} |
| 619 |
if (!numhorizon) { |
| 620 |
qh_precision ("empty horizon"); |
| 621 |
fprintf(qh ferr, "qhull precision error (qh_findhorizon): empty horizon\n\ |
| 622 |
Point p%d was above all facets.\n", qh_pointid(point)); |
| 623 |
qh_printfacetlist (qh facet_list, NULL, True); |
| 624 |
qh_errexit(qh_ERRprec, NULL, NULL); |
| 625 |
} |
| 626 |
trace1((qh ferr, "qh_findhorizon: %d horizon facets (good %d), %d visible (good %d), %d coplanar\n", numhorizon, *goodhorizon, qh num_visible, *goodvisible, coplanar)); |
| 627 |
if (qh IStracing >= 4 && qh num_facets < 50) |
| 628 |
qh_printlists (); |
| 629 |
} /* findhorizon */ |
| 630 |
|
| 631 |
/*-<a href="qh-qhull.htm#TOC" |
| 632 |
>-------------------------------</a><a name="nextfurthest">-</a> |
| 633 |
|
| 634 |
qh_nextfurthest( visible ) |
| 635 |
returns next furthest point and visible facet for qh_addpoint() |
| 636 |
starts search at qh.facet_next |
| 637 |
|
| 638 |
returns: |
| 639 |
removes furthest point from outside set |
| 640 |
NULL if none available |
| 641 |
advances qh.facet_next over facets with empty outside sets |
| 642 |
|
| 643 |
design: |
| 644 |
for each facet from qh.facet_next |
| 645 |
if empty outside set |
| 646 |
advance qh.facet_next |
| 647 |
else if qh.NARROWhull |
| 648 |
determine furthest outside point |
| 649 |
if furthest point is not outside |
| 650 |
advance qh.facet_next (point will be coplanar) |
| 651 |
remove furthest point from outside set |
| 652 |
*/ |
| 653 |
pointT *qh_nextfurthest (facetT **visible) { |
| 654 |
facetT *facet; |
| 655 |
int size, index; |
| 656 |
realT randr, dist; |
| 657 |
pointT *furthest; |
| 658 |
|
| 659 |
while ((facet= qh facet_next) != qh facet_tail) { |
| 660 |
if (!facet->outsideset) { |
| 661 |
qh facet_next= facet->next; |
| 662 |
continue; |
| 663 |
} |
| 664 |
SETreturnsize_(facet->outsideset, size); |
| 665 |
if (!size) { |
| 666 |
qh_setfree (&facet->outsideset); |
| 667 |
qh facet_next= facet->next; |
| 668 |
continue; |
| 669 |
} |
| 670 |
if (qh NARROWhull) { |
| 671 |
if (facet->notfurthest) |
| 672 |
qh_furthestout (facet); |
| 673 |
furthest= (pointT*)qh_setlast (facet->outsideset); |
| 674 |
#if qh_COMPUTEfurthest |
| 675 |
qh_distplane (furthest, facet, &dist); |
| 676 |
zinc_(Zcomputefurthest); |
| 677 |
#else |
| 678 |
dist= facet->furthestdist; |
| 679 |
#endif |
| 680 |
if (dist < qh MINoutside) { /* remainder of outside set is coplanar for qh_outcoplanar */ |
| 681 |
qh facet_next= facet->next; |
| 682 |
continue; |
| 683 |
} |
| 684 |
} |
| 685 |
if (!qh RANDOMoutside && !qh VIRTUALmemory) { |
| 686 |
if (qh PICKfurthest) { |
| 687 |
qh_furthestnext (/* qh facet_list */); |
| 688 |
facet= qh facet_next; |
| 689 |
} |
| 690 |
*visible= facet; |
| 691 |
return ((pointT*)qh_setdellast (facet->outsideset)); |
| 692 |
} |
| 693 |
if (qh RANDOMoutside) { |
| 694 |
int outcoplanar = 0; |
| 695 |
if (qh NARROWhull) { |
| 696 |
FORALLfacets { |
| 697 |
if (facet == qh facet_next) |
| 698 |
break; |
| 699 |
if (facet->outsideset) |
| 700 |
outcoplanar += qh_setsize( facet->outsideset); |
| 701 |
} |
| 702 |
} |
| 703 |
randr= qh_RANDOMint; |
| 704 |
randr= randr/(qh_RANDOMmax+1); |
| 705 |
index= (int)floor((qh num_outside - outcoplanar) * randr); |
| 706 |
FORALLfacet_(qh facet_next) { |
| 707 |
if (facet->outsideset) { |
| 708 |
SETreturnsize_(facet->outsideset, size); |
| 709 |
if (!size) |
| 710 |
qh_setfree (&facet->outsideset); |
| 711 |
else if (size > index) { |
| 712 |
*visible= facet; |
| 713 |
return ((pointT*)qh_setdelnth (facet->outsideset, index)); |
| 714 |
}else |
| 715 |
index -= size; |
| 716 |
} |
| 717 |
} |
| 718 |
fprintf (qh ferr, "qhull internal error (qh_nextfurthest): num_outside %d is too low\nby at least %d, or a random real %g >= 1.0\n", |
| 719 |
qh num_outside, index+1, randr); |
| 720 |
qh_errexit (qh_ERRqhull, NULL, NULL); |
| 721 |
}else { /* VIRTUALmemory */ |
| 722 |
facet= qh facet_tail->previous; |
| 723 |
if (!(furthest= (pointT*)qh_setdellast(facet->outsideset))) { |
| 724 |
if (facet->outsideset) |
| 725 |
qh_setfree (&facet->outsideset); |
| 726 |
qh_removefacet (facet); |
| 727 |
qh_prependfacet (facet, &qh facet_list); |
| 728 |
continue; |
| 729 |
} |
| 730 |
*visible= facet; |
| 731 |
return furthest; |
| 732 |
} |
| 733 |
} |
| 734 |
return NULL; |
| 735 |
} /* nextfurthest */ |
| 736 |
|
| 737 |
/*-<a href="qh-qhull.htm#TOC" |
| 738 |
>-------------------------------</a><a name="partitionall">-</a> |
| 739 |
|
| 740 |
qh_partitionall( vertices, points, numpoints ) |
| 741 |
partitions all points in points/numpoints to the outsidesets of facets |
| 742 |
vertices= vertices in qh.facet_list (not partitioned) |
| 743 |
|
| 744 |
returns: |
| 745 |
builds facet->outsideset |
| 746 |
does not partition qh.GOODpoint |
| 747 |
if qh.ONLYgood && !qh.MERGING, |
| 748 |
does not partition qh.GOODvertex |
| 749 |
|
| 750 |
notes: |
| 751 |
faster if qh.facet_list sorted by anticipated size of outside set |
| 752 |
|
| 753 |
design: |
| 754 |
initialize pointset with all points |
| 755 |
remove vertices from pointset |
| 756 |
remove qh.GOODpointp from pointset (unless it's qh.STOPcone or qh.STOPpoint) |
| 757 |
for all facets |
| 758 |
for all remaining points in pointset |
| 759 |
compute distance from point to facet |
| 760 |
if point is outside facet |
| 761 |
remove point from pointset (by not reappending) |
| 762 |
update bestpoint |
| 763 |
append point or old bestpoint to facet's outside set |
| 764 |
append bestpoint to facet's outside set (furthest) |
| 765 |
for all points remaining in pointset |
| 766 |
partition point into facets' outside sets and coplanar sets |
| 767 |
*/ |
| 768 |
void qh_partitionall(setT *vertices, pointT *points, int numpoints){ |
| 769 |
setT *pointset; |
| 770 |
vertexT *vertex, **vertexp; |
| 771 |
pointT *point, **pointp, *bestpoint; |
| 772 |
int size, point_i, point_n, point_end, remaining, i, id; |
| 773 |
facetT *facet; |
| 774 |
realT bestdist= -REALmax, dist, distoutside; |
| 775 |
|
| 776 |
trace1((qh ferr, "qh_partitionall: partition all points into outside sets\n")); |
| 777 |
pointset= qh_settemp (numpoints); |
| 778 |
qh num_outside= 0; |
| 779 |
pointp= SETaddr_(pointset, pointT); |
| 780 |
for (i=numpoints, point= points; i--; point += qh hull_dim) |
| 781 |
*(pointp++)= point; |
| 782 |
qh_settruncate (pointset, numpoints); |
| 783 |
FOREACHvertex_(vertices) { |
| 784 |
if ((id= qh_pointid(vertex->point)) >= 0) |
| 785 |
SETelem_(pointset, id)= NULL; |
| 786 |
} |
| 787 |
id= qh_pointid (qh GOODpointp); |
| 788 |
if (id >=0 && qh STOPcone-1 != id && -qh STOPpoint-1 != id) |
| 789 |
SETelem_(pointset, id)= NULL; |
| 790 |
if (qh GOODvertexp && qh ONLYgood && !qh MERGING) { /* matches qhull()*/ |
| 791 |
if ((id= qh_pointid(qh GOODvertexp)) >= 0) |
| 792 |
SETelem_(pointset, id)= NULL; |
| 793 |
} |
| 794 |
if (!qh BESToutside) { /* matches conditional for qh_partitionpoint below */ |
| 795 |
distoutside= qh_DISToutside; /* multiple of qh.MINoutside & qh.max_outside, see user.h */ |
| 796 |
zval_(Ztotpartition)= qh num_points - qh hull_dim - 1; /*misses GOOD... */ |
| 797 |
remaining= qh num_facets; |
| 798 |
point_end= numpoints; |
| 799 |
FORALLfacets { |
| 800 |
size= point_end/(remaining--) + 100; |
| 801 |
facet->outsideset= qh_setnew (size); |
| 802 |
bestpoint= NULL; |
| 803 |
point_end= 0; |
| 804 |
FOREACHpoint_i_(pointset) { |
| 805 |
if (point) { |
| 806 |
zzinc_(Zpartitionall); |
| 807 |
qh_distplane (point, facet, &dist); |
| 808 |
if (dist < distoutside) |
| 809 |
SETelem_(pointset, point_end++)= point; |
| 810 |
else { |
| 811 |
qh num_outside++; |
| 812 |
if (!bestpoint) { |
| 813 |
bestpoint= point; |
| 814 |
bestdist= dist; |
| 815 |
}else if (dist > bestdist) { |
| 816 |
qh_setappend (&facet->outsideset, bestpoint); |
| 817 |
bestpoint= point; |
| 818 |
bestdist= dist; |
| 819 |
}else |
| 820 |
qh_setappend (&facet->outsideset, point); |
| 821 |
} |
| 822 |
} |
| 823 |
} |
| 824 |
if (bestpoint) { |
| 825 |
qh_setappend (&facet->outsideset, bestpoint); |
| 826 |
#if !qh_COMPUTEfurthest |
| 827 |
facet->furthestdist= bestdist; |
| 828 |
#endif |
| 829 |
}else |
| 830 |
qh_setfree (&facet->outsideset); |
| 831 |
qh_settruncate (pointset, point_end); |
| 832 |
} |
| 833 |
} |
| 834 |
/* if !qh BESToutside, pointset contains points not assigned to outsideset */ |
| 835 |
if (qh BESToutside || qh MERGING || qh KEEPcoplanar || qh KEEPinside) { |
| 836 |
qh findbestnew= True; |
| 837 |
FOREACHpoint_i_(pointset) { |
| 838 |
if (point) |
| 839 |
qh_partitionpoint(point, qh facet_list); |
| 840 |
} |
| 841 |
qh findbestnew= False; |
| 842 |
} |
| 843 |
zzadd_(Zpartitionall, zzval_(Zpartition)); |
| 844 |
zzval_(Zpartition)= 0; |
| 845 |
qh_settempfree(&pointset); |
| 846 |
if (qh IStracing >= 4) |
| 847 |
qh_printfacetlist (qh facet_list, NULL, True); |
| 848 |
} /* partitionall */ |
| 849 |
|
| 850 |
|
| 851 |
/*-<a href="qh-qhull.htm#TOC" |
| 852 |
>-------------------------------</a><a name="partitioncoplanar">-</a> |
| 853 |
|
| 854 |
qh_partitioncoplanar( point, facet, dist ) |
| 855 |
partition coplanar point to a facet |
| 856 |
dist is distance from point to facet |
| 857 |
if dist NULL, |
| 858 |
searches for bestfacet and does nothing if inside |
| 859 |
if qh.findbestnew set, |
| 860 |
searches new facets instead of using qh_findbest() |
| 861 |
|
| 862 |
returns: |
| 863 |
qh.max_ouside updated |
| 864 |
if qh.KEEPcoplanar or qh.KEEPinside |
| 865 |
point assigned to best coplanarset |
| 866 |
|
| 867 |
notes: |
| 868 |
facet->maxoutside is updated at end by qh_check_maxout |
| 869 |
|
| 870 |
design: |
| 871 |
if dist undefined |
| 872 |
find best facet for point |
| 873 |
if point sufficiently below facet (depends on qh.NEARinside and qh.KEEPinside) |
| 874 |
exit |
| 875 |
if keeping coplanar/nearinside/inside points |
| 876 |
if point is above furthest coplanar point |
| 877 |
append point to coplanar set (it is the new furthest) |
| 878 |
update qh.max_outside |
| 879 |
else |
| 880 |
append point one before end of coplanar set |
| 881 |
else if point is clearly outside of qh.max_outside and bestfacet->coplanarset |
| 882 |
and bestfacet is more than perpendicular to facet |
| 883 |
repartition the point using qh_findbest() -- it may be put on an outsideset |
| 884 |
else |
| 885 |
update qh.max_outside |
| 886 |
*/ |
| 887 |
void qh_partitioncoplanar (pointT *point, facetT *facet, realT *dist) { |
| 888 |
facetT *bestfacet; |
| 889 |
pointT *oldfurthest; |
| 890 |
realT bestdist, dist2, angle; |
| 891 |
int numpart= 0, oldfindbest; |
| 892 |
boolT isoutside; |
| 893 |
|
| 894 |
qh WAScoplanar= True; |
| 895 |
if (!dist) { |
| 896 |
if (qh findbestnew) |
| 897 |
bestfacet= qh_findbestnew (point, facet, &bestdist, qh_ALL, &isoutside, &numpart); |
| 898 |
else |
| 899 |
bestfacet= qh_findbest (point, facet, qh_ALL, !qh_ISnewfacets, qh DELAUNAY, |
| 900 |
&bestdist, &isoutside, &numpart); |
| 901 |
zinc_(Ztotpartcoplanar); |
| 902 |
zzadd_(Zpartcoplanar, numpart); |
| 903 |
if (!qh DELAUNAY && !qh KEEPinside) { /* for 'd', bestdist skips upperDelaunay facets */ |
| 904 |
if (qh KEEPnearinside) { |
| 905 |
if (bestdist < -qh NEARinside) { |
| 906 |
zinc_(Zcoplanarinside); |
| 907 |
trace4((qh ferr, "qh_partitioncoplanar: point p%d is more than near-inside facet f%d dist %2.2g findbestnew %d\n", qh_pointid(point), bestfacet->id, bestdist, qh findbestnew)); |
| 908 |
return; |
| 909 |
} |
| 910 |
}else if (bestdist < -qh MAXcoplanar) { |
| 911 |
trace4((qh ferr, "qh_partitioncoplanar: point p%d is inside facet f%d dist %2.2g findbestnew %d\n", qh_pointid(point), bestfacet->id, bestdist, qh findbestnew)); |
| 912 |
zinc_(Zcoplanarinside); |
| 913 |
return; |
| 914 |
} |
| 915 |
} |
| 916 |
}else { |
| 917 |
bestfacet= facet; |
| 918 |
bestdist= *dist; |
| 919 |
} |
| 920 |
if (bestdist > qh max_outside) { |
| 921 |
if (!dist && facet != bestfacet) { |
| 922 |
zinc_(Zpartangle); |
| 923 |
angle= qh_getangle(facet->normal, bestfacet->normal); |
| 924 |
if (angle < 0) { |
| 925 |
/* typically due to deleted vertex and coplanar facets, e.g., |
| 926 |
RBOX 1000 s Z1 G1e-13 t1001185205 | QHULL Tv */ |
| 927 |
zinc_(Zpartflip); |
| 928 |
trace2((qh ferr, "qh_partitioncoplanar: repartition point p%d from f%d. It is above flipped facet f%d dist %2.2g\n", qh_pointid(point), facet->id, bestfacet->id, bestdist)); |
| 929 |
oldfindbest= qh findbestnew; |
| 930 |
qh findbestnew= False; |
| 931 |
qh_partitionpoint(point, bestfacet); |
| 932 |
qh findbestnew= oldfindbest; |
| 933 |
return; |
| 934 |
} |
| 935 |
} |
| 936 |
qh max_outside= bestdist; |
| 937 |
if (bestdist > qh TRACEdist) { |
| 938 |
fprintf (qh ferr, "qh_partitioncoplanar: ====== p%d from f%d increases max_outside to %2.2g of f%d last p%d\n", |
| 939 |
qh_pointid(point), facet->id, bestdist, bestfacet->id, qh furthest_id); |
| 940 |
qh_errprint ("DISTANT", facet, bestfacet, NULL, NULL); |
| 941 |
} |
| 942 |
} |
| 943 |
if (qh KEEPcoplanar + qh KEEPinside + qh KEEPnearinside) { |
| 944 |
oldfurthest= (pointT*)qh_setlast (bestfacet->coplanarset); |
| 945 |
if (oldfurthest) { |
| 946 |
zinc_(Zcomputefurthest); |
| 947 |
qh_distplane (oldfurthest, bestfacet, &dist2); |
| 948 |
} |
| 949 |
if (!oldfurthest || dist2 < bestdist) |
| 950 |
qh_setappend(&bestfacet->coplanarset, point); |
| 951 |
else |
| 952 |
qh_setappend2ndlast(&bestfacet->coplanarset, point); |
| 953 |
} |
| 954 |
trace4((qh ferr, "qh_partitioncoplanar: point p%d is coplanar with facet f%d (or inside) dist %2.2g\n", qh_pointid(point), bestfacet->id, bestdist)); |
| 955 |
} /* partitioncoplanar */ |
| 956 |
|
| 957 |
/*-<a href="qh-qhull.htm#TOC" |
| 958 |
>-------------------------------</a><a name="partitionpoint">-</a> |
| 959 |
|
| 960 |
qh_partitionpoint( point, facet ) |
| 961 |
assigns point to an outside set, coplanar set, or inside set (i.e., dropt) |
| 962 |
if qh.findbestnew |
| 963 |
uses qh_findbestnew() to search all new facets |
| 964 |
else |
| 965 |
uses qh_findbest() |
| 966 |
|
| 967 |
notes: |
| 968 |
after qh_distplane(), this and qh_findbest() are most expensive in 3-d |
| 969 |
|
| 970 |
design: |
| 971 |
find best facet for point |
| 972 |
(either exhaustive search of new facets or directed search from facet) |
| 973 |
if qh.NARROWhull |
| 974 |
retain coplanar and nearinside points as outside points |
| 975 |
if point is outside bestfacet |
| 976 |
if point above furthest point for bestfacet |
| 977 |
append point to outside set (it becomes the new furthest) |
| 978 |
if outside set was empty |
| 979 |
move bestfacet to end of qh.facet_list (i.e., after qh.facet_next) |
| 980 |
update bestfacet->furthestdist |
| 981 |
else |
| 982 |
append point one before end of outside set |
| 983 |
else if point is coplanar to bestfacet |
| 984 |
if keeping coplanar points or need to update qh.max_outside |
| 985 |
partition coplanar point into bestfacet |
| 986 |
else if near-inside point |
| 987 |
partition as coplanar point into bestfacet |
| 988 |
else is an inside point |
| 989 |
if keeping inside points |
| 990 |
partition as coplanar point into bestfacet |
| 991 |
*/ |
| 992 |
void qh_partitionpoint (pointT *point, facetT *facet) { |
| 993 |
realT bestdist; |
| 994 |
boolT isoutside; |
| 995 |
facetT *bestfacet; |
| 996 |
int numpart; |
| 997 |
#if qh_COMPUTEfurthest |
| 998 |
realT dist; |
| 999 |
#endif |
| 1000 |
|
| 1001 |
if (qh findbestnew) |
| 1002 |
bestfacet= qh_findbestnew (point, facet, &bestdist, qh BESToutside, &isoutside, &numpart); |
| 1003 |
else |
| 1004 |
bestfacet= qh_findbest (point, facet, qh BESToutside, qh_ISnewfacets, !qh_NOupper, |
| 1005 |
&bestdist, &isoutside, &numpart); |
| 1006 |
zinc_(Ztotpartition); |
| 1007 |
zzadd_(Zpartition, numpart); |
| 1008 |
if (qh NARROWhull) { |
| 1009 |
if (qh DELAUNAY && !isoutside && bestdist >= -qh MAXcoplanar) |
| 1010 |
qh_precision ("nearly incident point (narrow hull)"); |
| 1011 |
if (qh KEEPnearinside) { |
| 1012 |
if (bestdist >= -qh NEARinside) |
| 1013 |
isoutside= True; |
| 1014 |
}else if (bestdist >= -qh MAXcoplanar) |
| 1015 |
isoutside= True; |
| 1016 |
} |
| 1017 |
|
| 1018 |
if (isoutside) { |
| 1019 |
if (!bestfacet->outsideset |
| 1020 |
|| !qh_setlast (bestfacet->outsideset)) { |
| 1021 |
qh_setappend(&(bestfacet->outsideset), point); |
| 1022 |
if (!bestfacet->newfacet) { |
| 1023 |
qh_removefacet (bestfacet); /* make sure it's after qh facet_next */ |
| 1024 |
qh_appendfacet (bestfacet); |
| 1025 |
} |
| 1026 |
#if !qh_COMPUTEfurthest |
| 1027 |
bestfacet->furthestdist= bestdist; |
| 1028 |
#endif |
| 1029 |
}else { |
| 1030 |
#if qh_COMPUTEfurthest |
| 1031 |
zinc_(Zcomputefurthest); |
| 1032 |
qh_distplane (oldfurthest, bestfacet, &dist); |
| 1033 |
if (dist < bestdist) |
| 1034 |
qh_setappend(&(bestfacet->outsideset), point); |
| 1035 |
else |
| 1036 |
qh_setappend2ndlast(&(bestfacet->outsideset), point); |
| 1037 |
#else |
| 1038 |
if (bestfacet->furthestdist < bestdist) { |
| 1039 |
qh_setappend(&(bestfacet->outsideset), point); |
| 1040 |
bestfacet->furthestdist= bestdist; |
| 1041 |
}else |
| 1042 |
qh_setappend2ndlast(&(bestfacet->outsideset), point); |
| 1043 |
#endif |
| 1044 |
} |
| 1045 |
qh num_outside++; |
| 1046 |
trace4((qh ferr, "qh_partitionpoint: point p%d is outside facet f%d new? %d(or narrowhull)\n", qh_pointid(point), bestfacet->id, bestfacet->newfacet)); |
| 1047 |
}else if (qh DELAUNAY || bestdist >= -qh MAXcoplanar) { /* for 'd', bestdist skips upperDelaunay facets */ |
| 1048 |
zzinc_(Zcoplanarpart); |
| 1049 |
if (qh DELAUNAY) |
| 1050 |
qh_precision ("nearly incident point"); |
| 1051 |
if ((qh KEEPcoplanar + qh KEEPnearinside) || bestdist > qh max_outside) |
| 1052 |
qh_partitioncoplanar (point, bestfacet, &bestdist); |
| 1053 |
else { |
| 1054 |
trace4((qh ferr, "qh_partitionpoint: point p%d is coplanar to facet f%d (dropped)\n", qh_pointid(point), bestfacet->id)); |
| 1055 |
} |
| 1056 |
}else if (qh KEEPnearinside && bestdist > -qh NEARinside) { |
| 1057 |
zinc_(Zpartnear); |
| 1058 |
qh_partitioncoplanar (point, bestfacet, &bestdist); |
| 1059 |
}else { |
| 1060 |
zinc_(Zpartinside); |
| 1061 |
trace4((qh ferr, "qh_partitionpoint: point p%d is inside all facets, closest to f%d dist %2.2g\n", qh_pointid(point), bestfacet->id, bestdist)); |
| 1062 |
if (qh KEEPinside) |
| 1063 |
qh_partitioncoplanar (point, bestfacet, &bestdist); |
| 1064 |
} |
| 1065 |
} /* partitionpoint */ |
| 1066 |
|
| 1067 |
/*-<a href="qh-qhull.htm#TOC" |
| 1068 |
>-------------------------------</a><a name="partitionvisible">-</a> |
| 1069 |
|
| 1070 |
qh_partitionvisible( allpoints, numoutside ) |
| 1071 |
partitions points in visible facets to qh.newfacet_list |
| 1072 |
qh.visible_list= visible facets |
| 1073 |
for visible facets |
| 1074 |
1st neighbor (if any) points to a horizon facet or a new facet |
| 1075 |
if allpoints (not used), |
| 1076 |
repartitions coplanar points |
| 1077 |
|
| 1078 |
returns: |
| 1079 |
updates outside sets and coplanar sets of qh.newfacet_list |
| 1080 |
updates qh.num_outside (count of outside points) |
| 1081 |
|
| 1082 |
notes: |
| 1083 |
qh.findbest_notsharp should be clear (extra work if set) |
| 1084 |
|
| 1085 |
design: |
| 1086 |
for all visible facets with outside set or coplanar set |
| 1087 |
select a newfacet for visible facet |
| 1088 |
if outside set |
| 1089 |
partition outside set into new facets |
| 1090 |
if coplanar set and keeping coplanar/near-inside/inside points |
| 1091 |
if allpoints |
| 1092 |
partition coplanar set into new facets, may be assigned outside |
| 1093 |
else |
| 1094 |
partition coplanar set into coplanar sets of new facets |
| 1095 |
for each deleted vertex |
| 1096 |
if allpoints |
| 1097 |
partition vertex into new facets, may be assigned outside |
| 1098 |
else |
| 1099 |
partition vertex into coplanar sets of new facets |
| 1100 |
*/ |
| 1101 |
void qh_partitionvisible(/*visible_list*/ boolT allpoints, int *numoutside) { |
| 1102 |
facetT *visible, *newfacet; |
| 1103 |
pointT *point, **pointp; |
| 1104 |
int coplanar=0, size; |
| 1105 |
unsigned count; |
| 1106 |
vertexT *vertex, **vertexp; |
| 1107 |
|
| 1108 |
if (qh ONLYmax) |
| 1109 |
maximize_(qh MINoutside, qh max_vertex); |
| 1110 |
*numoutside= 0; |
| 1111 |
FORALLvisible_facets { |
| 1112 |
if (!visible->outsideset && !visible->coplanarset) |
| 1113 |
continue; |
| 1114 |
newfacet= visible->f.replace; |
| 1115 |
count= 0; |
| 1116 |
while (newfacet && newfacet->visible) { |
| 1117 |
newfacet= newfacet->f.replace; |
| 1118 |
if (count++ > qh facet_id) |
| 1119 |
qh_infiniteloop (visible); |
| 1120 |
} |
| 1121 |
if (!newfacet) |
| 1122 |
newfacet= qh newfacet_list; |
| 1123 |
if (newfacet == qh facet_tail) { |
| 1124 |
fprintf (qh ferr, "qhull precision error (qh_partitionvisible): all new facets deleted as\n degenerate facets. Can not continue.\n"); |
| 1125 |
qh_errexit (qh_ERRprec, NULL, NULL); |
| 1126 |
} |
| 1127 |
if (visible->outsideset) { |
| 1128 |
size= qh_setsize (visible->outsideset); |
| 1129 |
*numoutside += size; |
| 1130 |
qh num_outside -= size; |
| 1131 |
FOREACHpoint_(visible->outsideset) |
| 1132 |
qh_partitionpoint (point, newfacet); |
| 1133 |
} |
| 1134 |
if (visible->coplanarset && (qh KEEPcoplanar + qh KEEPinside + qh KEEPnearinside)) { |
| 1135 |
size= qh_setsize (visible->coplanarset); |
| 1136 |
coplanar += size; |
| 1137 |
FOREACHpoint_(visible->coplanarset) { |
| 1138 |
if (allpoints) /* not used */ |
| 1139 |
qh_partitionpoint (point, newfacet); |
| 1140 |
else |
| 1141 |
qh_partitioncoplanar (point, newfacet, NULL); |
| 1142 |
} |
| 1143 |
} |
| 1144 |
} |
| 1145 |
FOREACHvertex_(qh del_vertices) { |
| 1146 |
if (vertex->point) { |
| 1147 |
if (allpoints) /* not used */ |
| 1148 |
qh_partitionpoint (vertex->point, qh newfacet_list); |
| 1149 |
else |
| 1150 |
qh_partitioncoplanar (vertex->point, qh newfacet_list, NULL); |
| 1151 |
} |
| 1152 |
} |
| 1153 |
trace1((qh ferr,"qh_partitionvisible: partitioned %d points from outsidesets and %d points from coplanarsets\n", *numoutside, coplanar)); |
| 1154 |
} /* partitionvisible */ |
| 1155 |
|
| 1156 |
|
| 1157 |
|
| 1158 |
/*-<a href="qh-qhull.htm#TOC" |
| 1159 |
>-------------------------------</a><a name="precision">-</a> |
| 1160 |
|
| 1161 |
qh_precision( reason ) |
| 1162 |
restart on precision errors if not merging and if 'QJn' |
| 1163 |
*/ |
| 1164 |
void qh_precision (char *reason) { |
| 1165 |
|
| 1166 |
if (qh ALLOWrestart && !qh PREmerge && !qh MERGEexact) { |
| 1167 |
if (qh JOGGLEmax < REALmax/2) { |
| 1168 |
trace0((qh ferr, "qh_precision: qhull restart because of %s\n", reason)); |
| 1169 |
longjmp(qh restartexit, qh_ERRprec); |
| 1170 |
} |
| 1171 |
} |
| 1172 |
} /* qh_precision */ |
| 1173 |
|
| 1174 |
/*-<a href="qh-qhull.htm#TOC" |
| 1175 |
>-------------------------------</a><a name="printsummary">-</a> |
| 1176 |
|
| 1177 |
qh_printsummary( fp ) |
| 1178 |
prints summary to fp |
| 1179 |
|
| 1180 |
notes: |
| 1181 |
not in io.c so that user_eg.c can prevent io.c from loading |
| 1182 |
qh_printsummary and qh_countfacets must match counts |
| 1183 |
|
| 1184 |
design: |
| 1185 |
determine number of points, vertices, and coplanar points |
| 1186 |
print summary |
| 1187 |
*/ |
| 1188 |
void qh_printsummary(FILE *fp) { |
| 1189 |
realT ratio, outerplane, innerplane; |
| 1190 |
float cpu; |
| 1191 |
int size, id, nummerged, numvertices, numcoplanars= 0, nonsimplicial=0; |
| 1192 |
int goodused; |
| 1193 |
facetT *facet; |
| 1194 |
char *s; |
| 1195 |
int numdel= zzval_(Zdelvertextot); |
| 1196 |
int numtricoplanars= 0; |
| 1197 |
|
| 1198 |
size= qh num_points + qh_setsize (qh other_points); |
| 1199 |
numvertices= qh num_vertices - qh_setsize (qh del_vertices); |
| 1200 |
id= qh_pointid (qh GOODpointp); |
| 1201 |
FORALLfacets { |
| 1202 |
if (facet->coplanarset) |
| 1203 |
numcoplanars += qh_setsize( facet->coplanarset); |
| 1204 |
if (facet->good) { |
| 1205 |
if (facet->simplicial) { |
| 1206 |
if (facet->keepcentrum && facet->tricoplanar) |
| 1207 |
numtricoplanars++; |
| 1208 |
}else if (qh_setsize(facet->vertices) != qh hull_dim) |
| 1209 |
nonsimplicial++; |
| 1210 |
} |
| 1211 |
} |
| 1212 |
if (id >=0 && qh STOPcone-1 != id && -qh STOPpoint-1 != id) |
| 1213 |
size--; |
| 1214 |
if (qh STOPcone || qh STOPpoint) |
| 1215 |
fprintf (fp, "\nAt a premature exit due to 'TVn', 'TCn', 'TRn', or precision error."); |
| 1216 |
if (qh UPPERdelaunay) |
| 1217 |
goodused= qh GOODvertex + qh GOODpoint + qh SPLITthresholds; |
| 1218 |
else if (qh DELAUNAY) |
| 1219 |
goodused= qh GOODvertex + qh GOODpoint + qh GOODthreshold; |
| 1220 |
else |
| 1221 |
goodused= qh num_good; |
| 1222 |
nummerged= zzval_(Ztotmerge) - zzval_(Zcyclehorizon) + zzval_(Zcyclefacettot); |
| 1223 |
if (qh VORONOI) { |
| 1224 |
if (qh UPPERdelaunay) |
| 1225 |
fprintf (fp, "\n\ |
| 1226 |
Furthest-site Voronoi vertices by the convex hull of %d points in %d-d:\n\n", size, qh hull_dim); |
| 1227 |
else |
| 1228 |
fprintf (fp, "\n\ |
| 1229 |
Voronoi diagram by the convex hull of %d points in %d-d:\n\n", size, qh hull_dim); |
| 1230 |
fprintf(fp, " Number of Voronoi regions%s: %d\n", |
| 1231 |
qh ATinfinity ? " and at-infinity" : "", numvertices); |
| 1232 |
if (numdel) |
| 1233 |
fprintf(fp, " Total number of deleted points due to merging: %d\n", numdel); |
| 1234 |
if (numcoplanars - numdel > 0) |
| 1235 |
fprintf(fp, " Number of nearly incident points: %d\n", numcoplanars - numdel); |
| 1236 |
else if (size - numvertices - numdel > 0) |
| 1237 |
fprintf(fp, " Total number of nearly incident points: %d\n", size - numvertices - numdel); |
| 1238 |
fprintf(fp, " Number of%s Voronoi vertices: %d\n", |
| 1239 |
goodused ? " 'good'" : "", qh num_good); |
| 1240 |
if (nonsimplicial) |
| 1241 |
fprintf(fp, " Number of%s non-simplicial Voronoi vertices: %d\n", |
| 1242 |
goodused ? " 'good'" : "", nonsimplicial); |
| 1243 |
}else if (qh DELAUNAY) { |
| 1244 |
if (qh UPPERdelaunay) |
| 1245 |
fprintf (fp, "\n\ |
| 1246 |
Furthest-site Delaunay triangulation by the convex hull of %d points in %d-d:\n\n", size, qh hull_dim); |
| 1247 |
else |
| 1248 |
fprintf (fp, "\n\ |
| 1249 |
Delaunay triangulation by the convex hull of %d points in %d-d:\n\n", size, qh hull_dim); |
| 1250 |
fprintf(fp, " Number of input sites%s: %d\n", |
| 1251 |
qh ATinfinity ? " and at-infinity" : "", numvertices); |
| 1252 |
if (numdel) |
| 1253 |
fprintf(fp, " Total number of deleted points due to merging: %d\n", numdel); |
| 1254 |
if (numcoplanars - numdel > 0) |
| 1255 |
fprintf(fp, " Number of nearly incident points: %d\n", numcoplanars - numdel); |
| 1256 |
else if (size - numvertices - numdel > 0) |
| 1257 |
fprintf(fp, " Total number of nearly incident points: %d\n", size - numvertices - numdel); |
| 1258 |
fprintf(fp, " Number of%s Delaunay regions: %d\n", |
| 1259 |
goodused ? " 'good'" : "", qh num_good); |
| 1260 |
if (nonsimplicial) |
| 1261 |
fprintf(fp, " Number of%s non-simplicial Delaunay regions: %d\n", |
| 1262 |
goodused ? " 'good'" : "", nonsimplicial); |
| 1263 |
}else if (qh HALFspace) { |
| 1264 |
fprintf (fp, "\n\ |
| 1265 |
Halfspace intersection by the convex hull of %d points in %d-d:\n\n", size, qh hull_dim); |
| 1266 |
fprintf(fp, " Number of halfspaces: %d\n", size); |
| 1267 |
fprintf(fp, " Number of non-redundant halfspaces: %d\n", numvertices); |
| 1268 |
if (numcoplanars) { |
| 1269 |
if (qh KEEPinside && qh KEEPcoplanar) |
| 1270 |
s= "similar and redundant"; |
| 1271 |
else if (qh KEEPinside) |
| 1272 |
s= "redundant"; |
| 1273 |
else |
| 1274 |
s= "similar"; |
| 1275 |
fprintf(fp, " Number of %s halfspaces: %d\n", s, numcoplanars); |
| 1276 |
} |
| 1277 |
fprintf(fp, " Number of intersection points: %d\n", qh num_facets - qh num_visible); |
| 1278 |
if (goodused) |
| 1279 |
fprintf(fp, " Number of 'good' intersection points: %d\n", qh num_good); |
| 1280 |
if (nonsimplicial) |
| 1281 |
fprintf(fp, " Number of%s non-simplicial intersection points: %d\n", |
| 1282 |
goodused ? " 'good'" : "", nonsimplicial); |
| 1283 |
}else { |
| 1284 |
fprintf (fp, "\n\ |
| 1285 |
Convex hull of %d points in %d-d:\n\n", size, qh hull_dim); |
| 1286 |
fprintf(fp, " Number of vertices: %d\n", numvertices); |
| 1287 |
if (numcoplanars) { |
| 1288 |
if (qh KEEPinside && qh KEEPcoplanar) |
| 1289 |
s= "coplanar and interior"; |
| 1290 |
else if (qh KEEPinside) |
| 1291 |
s= "interior"; |
| 1292 |
else |
| 1293 |
s= "coplanar"; |
| 1294 |
fprintf(fp, " Number of %s points: %d\n", s, numcoplanars); |
| 1295 |
} |
| 1296 |
fprintf(fp, " Number of facets: %d\n", qh num_facets - qh num_visible); |
| 1297 |
if (goodused) |
| 1298 |
fprintf(fp, " Number of 'good' facets: %d\n", qh num_good); |
| 1299 |
if (nonsimplicial) |
| 1300 |
fprintf(fp, " Number of%s non-simplicial facets: %d\n", |
| 1301 |
goodused ? " 'good'" : "", nonsimplicial); |
| 1302 |
} |
| 1303 |
if (numtricoplanars) |
| 1304 |
fprintf(fp, " Number of triangulated facets: %d\n", numtricoplanars); |
| 1305 |
fprintf(fp, "\nStatistics for: %s | %s", |
| 1306 |
qh rbox_command, qh qhull_command); |
| 1307 |
if (qh ROTATErandom != INT_MIN) |
| 1308 |
fprintf(fp, " QR%d\n\n", qh ROTATErandom); |
| 1309 |
else |
| 1310 |
fprintf(fp, "\n\n"); |
| 1311 |
fprintf(fp, " Number of points processed: %d\n", zzval_(Zprocessed)); |
| 1312 |
fprintf(fp, " Number of hyperplanes created: %d\n", zzval_(Zsetplane)); |
| 1313 |
if (qh DELAUNAY) |
| 1314 |
fprintf(fp, " Number of facets in hull: %d\n", qh num_facets - qh num_visible); |
| 1315 |
fprintf(fp, " Number of distance tests for qhull: %d\n", zzval_(Zpartition)+ |
| 1316 |
zzval_(Zpartitionall)+zzval_(Znumvisibility)+zzval_(Zpartcoplanar)); |
| 1317 |
#if 0 |
| 1318 |
/* NOTE: must print before printstatistics() */ |
| 1319 |
{realT stddev, ave; |
| 1320 |
fprintf(fp, " average new facet balance: %2.2g\n", |
| 1321 |
wval_(Wnewbalance)/zval_(Zprocessed)); |
| 1322 |
stddev= qh_stddev (zval_(Zprocessed), wval_(Wnewbalance), |
| 1323 |
wval_(Wnewbalance2), &ave); |
| 1324 |
fprintf(fp, " new facet standard deviation: %2.2g\n", stddev); |
| 1325 |
fprintf(fp, " average partition balance: %2.2g\n", |
| 1326 |
wval_(Wpbalance)/zval_(Zpbalance)); |
| 1327 |
stddev= qh_stddev (zval_(Zpbalance), wval_(Wpbalance), |
| 1328 |
wval_(Wpbalance2), &ave); |
| 1329 |
fprintf(fp, " partition standard deviation: %2.2g\n", stddev); |
| 1330 |
} |
| 1331 |
#endif |
| 1332 |
if (nummerged) { |
| 1333 |
fprintf(fp," Number of distance tests for merging: %d\n",zzval_(Zbestdist)+ |
| 1334 |
zzval_(Zcentrumtests)+zzval_(Zdistconvex)+zzval_(Zdistcheck)+ |
| 1335 |
zzval_(Zdistzero)); |
| 1336 |
fprintf(fp," Number of distance tests for checking: %d\n",zzval_(Zcheckpart)); |
| 1337 |
fprintf(fp," Number of merged facets: %d\n", nummerged); |
| 1338 |
} |
| 1339 |
if (!qh RANDOMoutside && qh QHULLfinished) { |
| 1340 |
cpu= qh hulltime; |
| 1341 |
cpu /= qh_SECticks; |
| 1342 |
wval_(Wcpu)= cpu; |
| 1343 |
fprintf (fp, " CPU seconds to compute hull (after input): %2.4g\n", cpu); |
| 1344 |
} |
| 1345 |
if (qh RERUN) { |
| 1346 |
if (!qh PREmerge && !qh MERGEexact) |
| 1347 |
fprintf(fp, " Percentage of runs with precision errors: %4.1f\n", |
| 1348 |
zzval_(Zretry)*100.0/qh build_cnt); /* careful of order */ |
| 1349 |
}else if (qh JOGGLEmax < REALmax/2) { |
| 1350 |
if (zzval_(Zretry)) |
| 1351 |
fprintf(fp, " After %d retries, input joggled by: %2.2g\n", |
| 1352 |
zzval_(Zretry), qh JOGGLEmax); |
| 1353 |
else |
| 1354 |
fprintf(fp, " Input joggled by: %2.2g\n", qh JOGGLEmax); |
| 1355 |
} |
| 1356 |
if (qh totarea != 0.0) |
| 1357 |
fprintf(fp, " %s facet area: %2.8g\n", |
| 1358 |
zzval_(Ztotmerge) ? "Approximate" : "Total", qh totarea); |
| 1359 |
if (qh totvol != 0.0) |
| 1360 |
fprintf(fp, " %s volume: %2.8g\n", |
| 1361 |
zzval_(Ztotmerge) ? "Approximate" : "Total", qh totvol); |
| 1362 |
if (qh MERGING) { |
| 1363 |
qh_outerinner (NULL, &outerplane, &innerplane); |
| 1364 |
if (outerplane > 2 * qh DISTround) { |
| 1365 |
fprintf(fp, " Maximum distance of %spoint above facet: %2.2g", |
| 1366 |
(qh QHULLfinished ? "" : "merged "), outerplane); |
| 1367 |
ratio= outerplane/(qh ONEmerge + qh DISTround); |
| 1368 |
/* don't report ratio if MINoutside is large */ |
| 1369 |
if (ratio > 0.05 && 2* qh ONEmerge > qh MINoutside && qh JOGGLEmax > REALmax/2) |
| 1370 |
fprintf (fp, " (%.1fx)\n", ratio); |
| 1371 |
else |
| 1372 |
fprintf (fp, "\n"); |
| 1373 |
} |
| 1374 |
if (innerplane < -2 * qh DISTround) { |
| 1375 |
fprintf(fp, " Maximum distance of %svertex below facet: %2.2g", |
| 1376 |
(qh QHULLfinished ? "" : "merged "), innerplane); |
| 1377 |
ratio= -innerplane/(qh ONEmerge+qh DISTround); |
| 1378 |
if (ratio > 0.05 && qh JOGGLEmax > REALmax/2) |
| 1379 |
fprintf (fp, " (%.1fx)\n", ratio); |
| 1380 |
else |
| 1381 |
fprintf (fp, "\n"); |
| 1382 |
} |
| 1383 |
} |
| 1384 |
fprintf(fp, "\n"); |
| 1385 |
} /* printsummary */ |
| 1386 |
|
| 1387 |
|