ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/QuickHull/poly2.c
Revision: 1138
Committed: Tue May 29 22:51:00 2007 UTC (17 years, 11 months ago) by chuckv
Content type: text/plain
File size: 105936 byte(s)
Log Message:
Addded QuickHull to cvs.

File Contents

# Content
1 /*<html><pre> -<a href="qh-poly.htm"
2 >-------------------------------</a><a name="TOP">-</a>
3
4 poly2.c
5 implements polygons and simplices
6
7 see qh-poly.htm, poly.h and qhull.h
8
9 frequently used code is in poly.c
10
11 copyright (c) 1993-2003, The Geometry Center
12 */
13
14 #include "QuickHull/qhull_a.h"
15
16 /*======== functions in alphabetical order ==========*/
17
18 /*-<a href="qh-poly.htm#TOC"
19 >-------------------------------</a><a name="addhash">-</a>
20
21 qh_addhash( newelem, hashtable, hashsize, hash )
22 add newelem to linear hash table at hash if not already there
23 */
24 void qh_addhash (void* newelem, setT *hashtable, int hashsize, unsigned hash) {
25 int scan;
26 void *elem;
27
28 for (scan= (int)hash; (elem= SETelem_(hashtable, scan));
29 scan= (++scan >= hashsize ? 0 : scan)) {
30 if (elem == newelem)
31 break;
32 }
33 /* loop terminates because qh_HASHfactor >= 1.1 by qh_initbuffers */
34 if (!elem)
35 SETelem_(hashtable, scan)= newelem;
36 } /* addhash */
37
38 /*-<a href="qh-poly.htm#TOC"
39 >-------------------------------</a><a name="check_bestdist">-</a>
40
41 qh_check_bestdist()
42 check that all points are within max_outside of the nearest facet
43 if qh.ONLYgood,
44 ignores !good facets
45
46 see:
47 qh_check_maxout(), qh_outerinner()
48
49 notes:
50 only called from qh_check_points()
51 seldom used since qh.MERGING is almost always set
52 if notverified>0 at end of routine
53 some points were well inside the hull. If the hull contains
54 a lens-shaped component, these points were not verified. Use
55 options 'Qi Tv' to verify all points. (Exhaustive check also verifies)
56
57 design:
58 determine facet for each point (if any)
59 for each point
60 start with the assigned facet or with the first facet
61 find the best facet for the point and check all coplanar facets
62 error if point is outside of facet
63 */
64 void qh_check_bestdist (void) {
65 boolT waserror= False, unassigned;
66 facetT *facet, *bestfacet, *errfacet1= NULL, *errfacet2= NULL;
67 facetT *facetlist;
68 realT dist, maxoutside, maxdist= -REALmax;
69 pointT *point;
70 int numpart= 0, facet_i, facet_n, notgood= 0, notverified= 0;
71 setT *facets;
72
73 trace1((qh ferr, "qh_check_bestdist: check points below nearest facet. Facet_list f%d\n", qh facet_list->id));
74 maxoutside= qh_maxouter();
75 maxoutside += qh DISTround;
76 /* one more qh.DISTround for check computation */
77 trace1((qh ferr, "qh_check_bestdist: check that all points are within %2.2g of best facet\n", maxoutside));
78 facets= qh_pointfacet (/*qh facet_list*/);
79 if (!qh_QUICKhelp && qh PRINTprecision)
80 fprintf (qh ferr, "\n\
81 qhull output completed. Verifying that %d points are\n\
82 below %2.2g of the nearest %sfacet.\n",
83 qh_setsize(facets), maxoutside, (qh ONLYgood ? "good " : ""));
84 FOREACHfacet_i_(facets) { /* for each point with facet assignment */
85 if (facet)
86 unassigned= False;
87 else {
88 unassigned= True;
89 facet= qh facet_list;
90 }
91 point= qh_point(facet_i);
92 if (point == qh GOODpointp)
93 continue;
94 qh_distplane(point, facet, &dist);
95 numpart++;
96 bestfacet= qh_findbesthorizon (!qh_IScheckmax, point, facet, qh_NOupper, &dist, &numpart);
97 /* occurs after statistics reported */
98 maximize_(maxdist, dist);
99 if (dist > maxoutside) {
100 if (qh ONLYgood && !bestfacet->good
101 && !((bestfacet= qh_findgooddist (point, bestfacet, &dist, &facetlist))
102 && dist > maxoutside))
103 notgood++;
104 else {
105 waserror= True;
106 fprintf(qh ferr, "qhull precision error: point p%d is outside facet f%d, distance= %6.8g maxoutside= %6.8g\n",
107 facet_i, bestfacet->id, dist, maxoutside);
108 if (errfacet1 != bestfacet) {
109 errfacet2= errfacet1;
110 errfacet1= bestfacet;
111 }
112 }
113 }else if (unassigned && dist < -qh MAXcoplanar)
114 notverified++;
115 }
116 qh_settempfree (&facets);
117 if (notverified && !qh DELAUNAY && !qh_QUICKhelp && qh PRINTprecision)
118 fprintf(qh ferr, "\n%d points were well inside the hull. If the hull contains\n\
119 a lens-shaped component, these points were not verified. Use\n\
120 options 'Qci Tv' to verify all points.\n", notverified);
121 if (maxdist > qh outside_err) {
122 fprintf( qh ferr, "qhull precision error (qh_check_bestdist): a coplanar point is %6.2g from convex hull. The maximum value (qh.outside_err) is %6.2g\n",
123 maxdist, qh outside_err);
124 qh_errexit2 (qh_ERRprec, errfacet1, errfacet2);
125 }else if (waserror && qh outside_err > REALmax/2)
126 qh_errexit2 (qh_ERRprec, errfacet1, errfacet2);
127 else if (waserror)
128 ; /* the error was logged to qh.ferr but does not effect the output */
129 trace0((qh ferr, "qh_check_bestdist: max distance outside %2.2g\n", maxdist));
130 } /* check_bestdist */
131
132 /*-<a href="qh-poly.htm#TOC"
133 >-------------------------------</a><a name="check_maxout">-</a>
134
135 qh_check_maxout()
136 updates qh.max_outside by checking all points against bestfacet
137 if qh.ONLYgood, ignores !good facets
138
139 returns:
140 updates facet->maxoutside via qh_findbesthorizon()
141 sets qh.maxoutdone
142 if printing qh.min_vertex (qh_outerinner),
143 it is updated to the current vertices
144 removes inside/coplanar points from coplanarset as needed
145
146 notes:
147 defines coplanar as min_vertex instead of MAXcoplanar
148 may not need to check near-inside points because of qh.MAXcoplanar
149 and qh.KEEPnearinside (before it was -DISTround)
150
151 see also:
152 qh_check_bestdist()
153
154 design:
155 if qh.min_vertex is needed
156 for all neighbors of all vertices
157 test distance from vertex to neighbor
158 determine facet for each point (if any)
159 for each point with an assigned facet
160 find the best facet for the point and check all coplanar facets
161 (updates outer planes)
162 remove near-inside points from coplanar sets
163 */
164 #ifndef qh_NOmerge
165 void qh_check_maxout (void) {
166 facetT *facet, *bestfacet, *neighbor, **neighborp, *facetlist;
167 realT dist, maxoutside, minvertex, old_maxoutside;
168 pointT *point;
169 int numpart= 0, facet_i, facet_n, notgood= 0;
170 setT *facets, *vertices;
171 vertexT *vertex;
172
173 trace1((qh ferr, "qh_check_maxout: check and update maxoutside for each facet.\n"));
174 maxoutside= minvertex= 0;
175 if (qh VERTEXneighbors
176 && (qh PRINTsummary || qh KEEPinside || qh KEEPcoplanar
177 || qh TRACElevel || qh PRINTstatistics
178 || qh PRINTout[0] == qh_PRINTsummary || qh PRINTout[0] == qh_PRINTnone)) {
179 trace1((qh ferr, "qh_check_maxout: determine actual maxoutside and minvertex\n"));
180 vertices= qh_pointvertex (/*qh facet_list*/);
181 FORALLvertices {
182 FOREACHneighbor_(vertex) {
183 zinc_(Zdistvertex); /* distance also computed by main loop below */
184 qh_distplane (vertex->point, neighbor, &dist);
185 minimize_(minvertex, dist);
186 if (-dist > qh TRACEdist || dist > qh TRACEdist
187 || neighbor == qh tracefacet || vertex == qh tracevertex)
188 fprintf (qh ferr, "qh_check_maxout: p%d (v%d) is %.2g from f%d\n",
189 qh_pointid (vertex->point), vertex->id, dist, neighbor->id);
190 }
191 }
192 if (qh MERGING) {
193 wmin_(Wminvertex, qh min_vertex);
194 }
195 qh min_vertex= minvertex;
196 qh_settempfree (&vertices);
197 }
198 facets= qh_pointfacet (/*qh facet_list*/);
199 do {
200 old_maxoutside= fmax_(qh max_outside, maxoutside);
201 FOREACHfacet_i_(facets) { /* for each point with facet assignment */
202 if (facet) {
203 point= qh_point(facet_i);
204 if (point == qh GOODpointp)
205 continue;
206 zinc_(Ztotcheck);
207 qh_distplane(point, facet, &dist);
208 numpart++;
209 bestfacet= qh_findbesthorizon (qh_IScheckmax, point, facet, !qh_NOupper, &dist, &numpart);
210 if (bestfacet && dist > maxoutside) {
211 if (qh ONLYgood && !bestfacet->good
212 && !((bestfacet= qh_findgooddist (point, bestfacet, &dist, &facetlist))
213 && dist > maxoutside))
214 notgood++;
215 else
216 maxoutside= dist;
217 }
218 if (dist > qh TRACEdist || (bestfacet && bestfacet == qh tracefacet))
219 fprintf (qh ferr, "qh_check_maxout: p%d is %.2g above f%d\n",
220 qh_pointid (point), dist, bestfacet->id);
221 }
222 }
223 }while
224 (maxoutside > 2*old_maxoutside);
225 /* if qh.maxoutside increases substantially, qh_SEARCHdist is not valid
226 e.g., RBOX 5000 s Z1 G1e-13 t1001200614 | qhull */
227 zzadd_(Zcheckpart, numpart);
228 qh_settempfree (&facets);
229 wval_(Wmaxout)= maxoutside - qh max_outside;
230 wmax_(Wmaxoutside, qh max_outside);
231 qh max_outside= maxoutside;
232 qh_nearcoplanar (/*qh.facet_list*/);
233 qh maxoutdone= True;
234 trace1((qh ferr, "qh_check_maxout: maxoutside %2.2g, min_vertex %2.2g, outside of not good %d\n", maxoutside, qh min_vertex, notgood));
235 } /* check_maxout */
236 #else /* qh_NOmerge */
237 void qh_check_maxout (void) {
238 }
239 #endif
240
241 /*-<a href="qh-poly.htm#TOC"
242 >-------------------------------</a><a name="check_output">-</a>
243
244 qh_check_output()
245 performs the checks at the end of qhull algorithm
246 Maybe called after voronoi output. Will recompute otherwise centrums are Voronoi centers instead
247 */
248 void qh_check_output (void) {
249 int i;
250
251 if (qh STOPcone)
252 return;
253 if (qh VERIFYoutput | qh IStracing | qh CHECKfrequently) {
254 qh_checkpolygon (qh facet_list);
255 qh_checkflipped_all (qh facet_list);
256 qh_checkconvex (qh facet_list, qh_ALGORITHMfault);
257 }else if (!qh MERGING && qh_newstats (qhstat precision, &i)) {
258 qh_checkflipped_all (qh facet_list);
259 qh_checkconvex (qh facet_list, qh_ALGORITHMfault);
260 }
261 } /* check_output */
262
263
264
265 /*-<a href="qh-poly.htm#TOC"
266 >-------------------------------</a><a name="check_point">-</a>
267
268 qh_check_point( point, facet, maxoutside, maxdist, errfacet1, errfacet2 )
269 check that point is less than maxoutside from facet
270 */
271 void qh_check_point (pointT *point, facetT *facet, realT *maxoutside, realT *maxdist, facetT **errfacet1, facetT **errfacet2) {
272 realT dist;
273
274 /* occurs after statistics reported */
275 qh_distplane(point, facet, &dist);
276 if (dist > *maxoutside) {
277 if (*errfacet1 != facet) {
278 *errfacet2= *errfacet1;
279 *errfacet1= facet;
280 }
281 fprintf(qh ferr, "qhull precision error: point p%d is outside facet f%d, distance= %6.8g maxoutside= %6.8g\n",
282 qh_pointid(point), facet->id, dist, *maxoutside);
283 }
284 maximize_(*maxdist, dist);
285 } /* qh_check_point */
286
287
288 /*-<a href="qh-poly.htm#TOC"
289 >-------------------------------</a><a name="check_points">-</a>
290
291 qh_check_points()
292 checks that all points are inside all facets
293
294 notes:
295 if many points and qh_check_maxout not called (i.e., !qh.MERGING),
296 calls qh_findbesthorizon (seldom done).
297 ignores flipped facets
298 maxoutside includes 2 qh.DISTrounds
299 one qh.DISTround for the computed distances in qh_check_points
300 qh_printafacet and qh_printsummary needs only one qh.DISTround
301 the computation for qh.VERIFYdirect does not account for qh.other_points
302
303 design:
304 if many points
305 please use qh_check_bestdist()
306 else
307 for all facets
308 for all points
309 check that point is inside facet
310 */
311 void qh_check_points (void) {
312 facetT *facet, *errfacet1= NULL, *errfacet2= NULL;
313 realT total, maxoutside, maxdist= -REALmax;
314 pointT *point, **pointp, *pointtemp;
315 boolT testouter;
316
317 maxoutside= qh_maxouter();
318 maxoutside += qh DISTround;
319 /* one more qh.DISTround for check computation */
320 trace1((qh ferr, "qh_check_points: check all points below %2.2g of all facet planes\n", maxoutside));
321 if (qh num_good) /* miss counts other_points and !good facets */
322 total= (float) qh num_good * qh num_points;
323 else
324 total= (float) qh num_facets * qh num_points;
325 if (total >= qh_VERIFYdirect && !qh maxoutdone) {
326 if (!qh_QUICKhelp && qh SKIPcheckmax && qh MERGING)
327 fprintf (qh ferr, "\n\
328 qhull input warning: merging without checking outer planes ('Q5' or 'Po').\n\
329 Verify may report that a point is outside of a facet.\n");
330 qh_check_bestdist();
331 }else {
332 if (qh_MAXoutside && qh maxoutdone)
333 testouter= True;
334 else
335 testouter= False;
336 if (!qh_QUICKhelp) {
337 if (qh MERGEexact)
338 fprintf (qh ferr, "\n\
339 qhull input warning: exact merge ('Qx'). Verify may report that a point\n\
340 is outside of a facet. See qh-optq.htm#Qx\n");
341 else if (qh SKIPcheckmax || qh NOnearinside)
342 fprintf (qh ferr, "\n\
343 qhull input warning: no outer plane check ('Q5') or no processing of\n\
344 near-inside points ('Q8'). Verify may report that a point is outside\n\
345 of a facet.\n");
346 }
347 if (qh PRINTprecision) {
348 if (testouter)
349 fprintf (qh ferr, "\n\
350 Output completed. Verifying that all points are below outer planes of\n\
351 all %sfacets. Will make %2.0f distance computations.\n",
352 (qh ONLYgood ? "good " : ""), total);
353 else
354 fprintf (qh ferr, "\n\
355 Output completed. Verifying that all points are below %2.2g of\n\
356 all %sfacets. Will make %2.0f distance computations.\n",
357 maxoutside, (qh ONLYgood ? "good " : ""), total);
358 }
359 FORALLfacets {
360 if (!facet->good && qh ONLYgood)
361 continue;
362 if (facet->flipped)
363 continue;
364 if (!facet->normal) {
365 fprintf( qh ferr, "qhull warning (qh_check_points): missing normal for facet f%d\n", facet->id);
366 continue;
367 }
368 if (testouter) {
369 #if qh_MAXoutside
370 maxoutside= facet->maxoutside + 2* qh DISTround;
371 /* one DISTround to actual point and another to computed point */
372 #endif
373 }
374 FORALLpoints {
375 if (point != qh GOODpointp)
376 qh_check_point (point, facet, &maxoutside, &maxdist, &errfacet1, &errfacet2);
377 }
378 FOREACHpoint_(qh other_points) {
379 if (point != qh GOODpointp)
380 qh_check_point (point, facet, &maxoutside, &maxdist, &errfacet1, &errfacet2);
381 }
382 }
383 if (maxdist > qh outside_err) {
384 fprintf( qh ferr, "qhull precision error (qh_check_points): a coplanar point is %6.2g from convex hull. The maximum value (qh.outside_err) is %6.2g\n",
385 maxdist, qh outside_err );
386 qh_errexit2( qh_ERRprec, errfacet1, errfacet2 );
387 }else if (errfacet1 && qh outside_err > REALmax/2)
388 qh_errexit2( qh_ERRprec, errfacet1, errfacet2 );
389 else if (errfacet1)
390 ; /* the error was logged to qh.ferr but does not effect the output */
391 trace0((qh ferr, "qh_check_points: max distance outside %2.2g\n", maxdist));
392 }
393 } /* check_points */
394
395
396 /*-<a href="qh-poly.htm#TOC"
397 >-------------------------------</a><a name="checkconvex">-</a>
398
399 qh_checkconvex( facetlist, fault )
400 check that each ridge in facetlist is convex
401 fault = qh_DATAfault if reporting errors
402 = qh_ALGORITHMfault otherwise
403
404 returns:
405 counts Zconcaveridges and Zcoplanarridges
406 errors if concaveridge or if merging an coplanar ridge
407
408 note:
409 if not merging,
410 tests vertices for neighboring simplicial facets
411 else if ZEROcentrum,
412 tests vertices for neighboring simplicial facets
413 else
414 tests centrums of neighboring facets
415
416 design:
417 for all facets
418 report flipped facets
419 if ZEROcentrum and simplicial neighbors
420 test vertices for neighboring simplicial facets
421 else
422 test centrum against all neighbors
423 */
424 void qh_checkconvex(facetT *facetlist, int fault) {
425 facetT *facet, *neighbor, **neighborp, *errfacet1=NULL, *errfacet2=NULL;
426 vertexT *vertex;
427 realT dist;
428 pointT *centrum;
429 boolT waserror= False, centrum_warning= False, tempcentrum= False, allsimplicial;
430 int neighbor_i;
431
432 trace1((qh ferr, "qh_checkconvex: check all ridges are convex\n"));
433 if (!qh RERUN) {
434 zzval_(Zconcaveridges)= 0;
435 zzval_(Zcoplanarridges)= 0;
436 }
437 FORALLfacet_(facetlist) {
438 if (facet->flipped) {
439 qh_precision ("flipped facet");
440 fprintf (qh ferr, "qhull precision error: f%d is flipped (interior point is outside)\n",
441 facet->id);
442 errfacet1= facet;
443 waserror= True;
444 continue;
445 }
446 if (qh MERGING && (!qh ZEROcentrum || !facet->simplicial || facet->tricoplanar))
447 allsimplicial= False;
448 else {
449 allsimplicial= True;
450 neighbor_i= 0;
451 FOREACHneighbor_(facet) {
452 vertex= SETelemt_(facet->vertices, neighbor_i++, vertexT);
453 if (!neighbor->simplicial || neighbor->tricoplanar) {
454 allsimplicial= False;
455 continue;
456 }
457 qh_distplane (vertex->point, neighbor, &dist);
458 if (dist > -qh DISTround) {
459 if (fault == qh_DATAfault) {
460 qh_precision ("coplanar or concave ridge");
461 fprintf (qh ferr, "qhull precision error: initial simplex is not convex. Distance=%.2g\n", dist);
462 qh_errexit(qh_ERRsingular, NULL, NULL);
463 }
464 if (dist > qh DISTround) {
465 zzinc_(Zconcaveridges);
466 qh_precision ("concave ridge");
467 fprintf (qh ferr, "qhull precision error: f%d is concave to f%d, since p%d (v%d) is %6.4g above\n",
468 facet->id, neighbor->id, qh_pointid(vertex->point), vertex->id, dist);
469 errfacet1= facet;
470 errfacet2= neighbor;
471 waserror= True;
472 }else if (qh ZEROcentrum) {
473 if (dist > 0) { /* qh_checkzero checks that dist < - qh DISTround */
474 zzinc_(Zcoplanarridges);
475 qh_precision ("coplanar ridge");
476 fprintf (qh ferr, "qhull precision error: f%d is clearly not convex to f%d, since p%d (v%d) is %6.4g above\n",
477 facet->id, neighbor->id, qh_pointid(vertex->point), vertex->id, dist);
478 errfacet1= facet;
479 errfacet2= neighbor;
480 waserror= True;
481 }
482 }else {
483 zzinc_(Zcoplanarridges);
484 qh_precision ("coplanar ridge");
485 trace0((qh ferr, "qhull precision error: f%d may be coplanar to f%d, since p%d (v%d) is within %6.4g during p%d\n", facet->id, neighbor->id, qh_pointid(vertex->point), vertex->id, dist, qh furthest_id));
486 }
487 }
488 }
489 }
490 if (!allsimplicial) {
491 if (qh CENTERtype == qh_AScentrum) {
492 if (!facet->center)
493 facet->center= qh_getcentrum (facet);
494 centrum= facet->center;
495 }else {
496 if (!centrum_warning && (!facet->simplicial || facet->tricoplanar)) {
497 centrum_warning= True;
498 fprintf (qh ferr, "qhull note: recomputing centrums for convexity test. This may lead to false, precision errors.\n");
499 }
500 centrum= qh_getcentrum(facet);
501 tempcentrum= True;
502 }
503 FOREACHneighbor_(facet) {
504 if (qh ZEROcentrum && facet->simplicial && neighbor->simplicial)
505 continue;
506 if (facet->tricoplanar || neighbor->tricoplanar)
507 continue;
508 zzinc_(Zdistconvex);
509 qh_distplane (centrum, neighbor, &dist);
510 if (dist > qh DISTround) {
511 zzinc_(Zconcaveridges);
512 qh_precision ("concave ridge");
513 fprintf (qh ferr, "qhull precision error: f%d is concave to f%d. Centrum of f%d is %6.4g above f%d\n",
514 facet->id, neighbor->id, facet->id, dist, neighbor->id);
515 errfacet1= facet;
516 errfacet2= neighbor;
517 waserror= True;
518 }else if (dist >= 0.0) { /* if arithmetic always rounds the same,
519 can test against centrum radius instead */
520 zzinc_(Zcoplanarridges);
521 qh_precision ("coplanar ridge");
522 fprintf (qh ferr, "qhull precision error: f%d is coplanar or concave to f%d. Centrum of f%d is %6.4g above f%d\n",
523 facet->id, neighbor->id, facet->id, dist, neighbor->id);
524 errfacet1= facet;
525 errfacet2= neighbor;
526 waserror= True;
527 }
528 }
529 if (tempcentrum)
530 qh_memfree(centrum, qh normal_size);
531 }
532 }
533 if (waserror && !qh FORCEoutput)
534 qh_errexit2 (qh_ERRprec, errfacet1, errfacet2);
535 } /* checkconvex */
536
537
538 /*-<a href="qh-poly.htm#TOC"
539 >-------------------------------</a><a name="checkfacet">-</a>
540
541 qh_checkfacet( facet, newmerge, waserror )
542 checks for consistency errors in facet
543 newmerge set if from merge.c
544
545 returns:
546 sets waserror if any error occurs
547
548 checks:
549 vertex ids are inverse sorted
550 unless newmerge, at least hull_dim neighbors and vertices (exactly if simplicial)
551 if non-simplicial, at least as many ridges as neighbors
552 neighbors are not duplicated
553 ridges are not duplicated
554 in 3-d, ridges=verticies
555 (qh.hull_dim-1) ridge vertices
556 neighbors are reciprocated
557 ridge neighbors are facet neighbors and a ridge for every neighbor
558 simplicial neighbors match facetintersect
559 vertex intersection matches vertices of common ridges
560 vertex neighbors and facet vertices agree
561 all ridges have distinct vertex sets
562
563 notes:
564 uses neighbor->seen
565
566 design:
567 check sets
568 check vertices
569 check sizes of neighbors and vertices
570 check for qh_MERGEridge and qh_DUPLICATEridge flags
571 check neighbor set
572 check ridge set
573 check ridges, neighbors, and vertices
574 */
575 void qh_checkfacet(facetT *facet, boolT newmerge, boolT *waserrorp) {
576 facetT *neighbor, **neighborp, *errother=NULL;
577 ridgeT *ridge, **ridgep, *errridge= NULL, *ridge2;
578 vertexT *vertex, **vertexp;
579 unsigned previousid= INT_MAX;
580 int numneighbors, numvertices, numridges=0, numRvertices=0;
581 boolT waserror= False;
582 int skipA, skipB, ridge_i, ridge_n, i;
583 setT *intersection;
584
585 if (facet->visible) {
586 fprintf (qh ferr, "qhull internal error (qh_checkfacet): facet f%d is on the visible_list\n",
587 facet->id);
588 qh_errexit (qh_ERRqhull, facet, NULL);
589 }
590 if (!facet->normal) {
591 fprintf (qh ferr, "qhull internal error (qh_checkfacet): facet f%d does not have a normal\n",
592 facet->id);
593 waserror= True;
594 }
595 qh_setcheck (facet->vertices, "vertices for f", facet->id);
596 qh_setcheck (facet->ridges, "ridges for f", facet->id);
597 qh_setcheck (facet->outsideset, "outsideset for f", facet->id);
598 qh_setcheck (facet->coplanarset, "coplanarset for f", facet->id);
599 qh_setcheck (facet->neighbors, "neighbors for f", facet->id);
600 FOREACHvertex_(facet->vertices) {
601 if (vertex->deleted) {
602 fprintf(qh ferr, "qhull internal error (qh_checkfacet): deleted vertex v%d in f%d\n", vertex->id, facet->id);
603 qh_errprint ("ERRONEOUS", NULL, NULL, NULL, vertex);
604 waserror= True;
605 }
606 if (vertex->id >= previousid) {
607 fprintf(qh ferr, "qhull internal error (qh_checkfacet): vertices of f%d are not in descending id order at v%d\n", facet->id, vertex->id);
608 waserror= True;
609 break;
610 }
611 previousid= vertex->id;
612 }
613 numneighbors= qh_setsize(facet->neighbors);
614 numvertices= qh_setsize(facet->vertices);
615 numridges= qh_setsize(facet->ridges);
616 if (facet->simplicial) {
617 if (numvertices+numneighbors != 2*qh hull_dim
618 && !facet->degenerate && !facet->redundant) {
619 fprintf(qh ferr, "qhull internal error (qh_checkfacet): for simplicial facet f%d, #vertices %d + #neighbors %d != 2*qh hull_dim\n",
620 facet->id, numvertices, numneighbors);
621 qh_setprint (qh ferr, "", facet->neighbors);
622 waserror= True;
623 }
624 }else { /* non-simplicial */
625 if (!newmerge
626 &&(numvertices < qh hull_dim || numneighbors < qh hull_dim)
627 && !facet->degenerate && !facet->redundant) {
628 fprintf(qh ferr, "qhull internal error (qh_checkfacet): for facet f%d, #vertices %d or #neighbors %d < qh hull_dim\n",
629 facet->id, numvertices, numneighbors);
630 waserror= True;
631 }
632 /* in 3-d, can get a vertex twice in an edge list, e.g., RBOX 1000 s W1e-13 t995849315 D2 | QHULL d Tc Tv TP624 TW1e-13 T4 */
633 if (numridges < numneighbors
634 ||(qh hull_dim == 3 && numvertices > numridges && !qh NEWfacets)
635 ||(qh hull_dim == 2 && numridges + numvertices + numneighbors != 6)) {
636 if (!facet->degenerate && !facet->redundant) {
637 fprintf(qh ferr, "qhull internal error (qh_checkfacet): for facet f%d, #ridges %d < #neighbors %d or (3-d) > #vertices %d or (2-d) not all 2\n",
638 facet->id, numridges, numneighbors, numvertices);
639 waserror= True;
640 }
641 }
642 }
643 FOREACHneighbor_(facet) {
644 if (neighbor == qh_MERGEridge || neighbor == qh_DUPLICATEridge) {
645 fprintf(qh ferr, "qhull internal error (qh_checkfacet): facet f%d still has a MERGE or DUP neighbor\n", facet->id);
646 qh_errexit (qh_ERRqhull, facet, NULL);
647 }
648 neighbor->seen= True;
649 }
650 FOREACHneighbor_(facet) {
651 if (!qh_setin(neighbor->neighbors, facet)) {
652 fprintf(qh ferr, "qhull internal error (qh_checkfacet): facet f%d has neighbor f%d, but f%d does not have neighbor f%d\n",
653 facet->id, neighbor->id, neighbor->id, facet->id);
654 errother= neighbor;
655 waserror= True;
656 }
657 if (!neighbor->seen) {
658 fprintf(qh ferr, "qhull internal error (qh_checkfacet): facet f%d has a duplicate neighbor f%d\n",
659 facet->id, neighbor->id);
660 errother= neighbor;
661 waserror= True;
662 }
663 neighbor->seen= False;
664 }
665 FOREACHridge_(facet->ridges) {
666 qh_setcheck (ridge->vertices, "vertices for r", ridge->id);
667 ridge->seen= False;
668 }
669 FOREACHridge_(facet->ridges) {
670 if (ridge->seen) {
671 fprintf(qh ferr, "qhull internal error (qh_checkfacet): facet f%d has a duplicate ridge r%d\n",
672 facet->id, ridge->id);
673 errridge= ridge;
674 waserror= True;
675 }
676 ridge->seen= True;
677 numRvertices= qh_setsize(ridge->vertices);
678 if (numRvertices != qh hull_dim - 1) {
679 fprintf(qh ferr, "qhull internal error (qh_checkfacet): ridge between f%d and f%d has %d vertices\n",
680 ridge->top->id, ridge->bottom->id, numRvertices);
681 errridge= ridge;
682 waserror= True;
683 }
684 neighbor= otherfacet_(ridge, facet);
685 neighbor->seen= True;
686 if (!qh_setin(facet->neighbors, neighbor)) {
687 fprintf(qh ferr, "qhull internal error (qh_checkfacet): for facet f%d, neighbor f%d of ridge r%d not in facet\n",
688 facet->id, neighbor->id, ridge->id);
689 errridge= ridge;
690 waserror= True;
691 }
692 }
693 if (!facet->simplicial) {
694 FOREACHneighbor_(facet) {
695 if (!neighbor->seen) {
696 fprintf(qh ferr, "qhull internal error (qh_checkfacet): facet f%d does not have a ridge for neighbor f%d\n",
697 facet->id, neighbor->id);
698 errother= neighbor;
699 waserror= True;
700 }
701 intersection= qh_vertexintersect_new(facet->vertices, neighbor->vertices);
702 qh_settemppush (intersection);
703 FOREACHvertex_(facet->vertices) {
704 vertex->seen= False;
705 vertex->seen2= False;
706 }
707 FOREACHvertex_(intersection)
708 vertex->seen= True;
709 FOREACHridge_(facet->ridges) {
710 if (neighbor != otherfacet_(ridge, facet))
711 continue;
712 FOREACHvertex_(ridge->vertices) {
713 if (!vertex->seen) {
714 fprintf (qh ferr, "qhull internal error (qh_checkfacet): vertex v%d in r%d not in f%d intersect f%d\n",
715 vertex->id, ridge->id, facet->id, neighbor->id);
716 qh_errexit (qh_ERRqhull, facet, ridge);
717 }
718 vertex->seen2= True;
719 }
720 }
721 if (!newmerge) {
722 FOREACHvertex_(intersection) {
723 if (!vertex->seen2) {
724 if (qh IStracing >=3 || !qh MERGING) {
725 fprintf (qh ferr, "qhull precision error (qh_checkfacet): vertex v%d in f%d intersect f%d but\n\
726 not in a ridge. This is ok under merging. Last point was p%d\n",
727 vertex->id, facet->id, neighbor->id, qh furthest_id);
728 if (!qh FORCEoutput && !qh MERGING) {
729 qh_errprint ("ERRONEOUS", facet, neighbor, NULL, vertex);
730 if (!qh MERGING)
731 qh_errexit (qh_ERRqhull, NULL, NULL);
732 }
733 }
734 }
735 }
736 }
737 qh_settempfree (&intersection);
738 }
739 }else { /* simplicial */
740 FOREACHneighbor_(facet) {
741 if (neighbor->simplicial) {
742 skipA= SETindex_(facet->neighbors, neighbor);
743 skipB= qh_setindex (neighbor->neighbors, facet);
744 if (!qh_setequal_skip (facet->vertices, skipA, neighbor->vertices, skipB)) {
745 fprintf (qh ferr, "qhull internal error (qh_checkfacet): facet f%d skip %d and neighbor f%d skip %d do not match \n",
746 facet->id, skipA, neighbor->id, skipB);
747 errother= neighbor;
748 waserror= True;
749 }
750 }
751 }
752 }
753 if (qh hull_dim < 5 && (qh IStracing > 2 || qh CHECKfrequently)) {
754 FOREACHridge_i_(facet->ridges) { /* expensive */
755 for (i= ridge_i+1; i < ridge_n; i++) {
756 ridge2= SETelemt_(facet->ridges, i, ridgeT);
757 if (qh_setequal (ridge->vertices, ridge2->vertices)) {
758 fprintf (qh ferr, "qh_checkfacet: ridges r%d and r%d have the same vertices\n",
759 ridge->id, ridge2->id);
760 errridge= ridge;
761 waserror= True;
762 }
763 }
764 }
765 }
766 if (waserror) {
767 qh_errprint("ERRONEOUS", facet, errother, errridge, NULL);
768 *waserrorp= True;
769 }
770 } /* checkfacet */
771
772
773 /*-<a href="qh-poly.htm#TOC"
774 >-------------------------------</a><a name="checkflipped_all">-</a>
775
776 qh_checkflipped_all( facetlist )
777 checks orientation of facets in list against interior point
778 */
779 void qh_checkflipped_all (facetT *facetlist) {
780 facetT *facet;
781 boolT waserror= False;
782 realT dist;
783
784 if (facetlist == qh facet_list)
785 zzval_(Zflippedfacets)= 0;
786 FORALLfacet_(facetlist) {
787 if (facet->normal && !qh_checkflipped (facet, &dist, !qh_ALL)) {
788 fprintf(qh ferr, "qhull precision error: facet f%d is flipped, distance= %6.12g\n",
789 facet->id, dist);
790 if (!qh FORCEoutput) {
791 qh_errprint("ERRONEOUS", facet, NULL, NULL, NULL);
792 waserror= True;
793 }
794 }
795 }
796 if (waserror) {
797 fprintf (qh ferr, "\n\
798 A flipped facet occurs when its distance to the interior point is\n\
799 greater than %2.2g, the maximum roundoff error.\n", -qh DISTround);
800 qh_errexit(qh_ERRprec, NULL, NULL);
801 }
802 } /* checkflipped_all */
803
804 /*-<a href="qh-poly.htm#TOC"
805 >-------------------------------</a><a name="checkpolygon">-</a>
806
807 qh_checkpolygon( facetlist )
808 checks the correctness of the structure
809
810 notes:
811 call with either qh.facet_list or qh.newfacet_list
812 checks num_facets and num_vertices if qh.facet_list
813
814 design:
815 for each facet
816 checks facet and outside set
817 initializes vertexlist
818 for each facet
819 checks vertex set
820 if checking all facets (qh.facetlist)
821 check facet count
822 if qh.VERTEXneighbors
823 check vertex neighbors and count
824 check vertex count
825 */
826 void qh_checkpolygon(facetT *facetlist) {
827 facetT *facet;
828 vertexT *vertex, **vertexp, *vertexlist;
829 int numfacets= 0, numvertices= 0, numridges= 0;
830 int totvneighbors= 0, totvertices= 0;
831 boolT waserror= False, nextseen= False, visibleseen= False;
832
833 trace1((qh ferr, "qh_checkpolygon: check all facets from f%d\n", facetlist->id));
834 if (facetlist != qh facet_list || qh ONLYgood)
835 nextseen= True;
836 FORALLfacet_(facetlist) {
837 if (facet == qh visible_list)
838 visibleseen= True;
839 if (!facet->visible) {
840 if (!nextseen) {
841 if (facet == qh facet_next)
842 nextseen= True;
843 else if (qh_setsize (facet->outsideset)) {
844 if (!qh NARROWhull
845 #if !qh_COMPUTEfurthest
846 || facet->furthestdist >= qh MINoutside
847 #endif
848 ) {
849 fprintf (qh ferr, "qhull internal error (qh_checkpolygon): f%d has outside points before qh facet_next\n",
850 facet->id);
851 qh_errexit (qh_ERRqhull, facet, NULL);
852 }
853 }
854 }
855 numfacets++;
856 qh_checkfacet(facet, False, &waserror);
857 }
858 }
859 if (qh visible_list && !visibleseen && facetlist == qh facet_list) {
860 fprintf (qh ferr, "qhull internal error (qh_checkpolygon): visible list f%d no longer on facet list\n", qh visible_list->id);
861 qh_printlists();
862 qh_errexit (qh_ERRqhull, qh visible_list, NULL);
863 }
864 if (facetlist == qh facet_list)
865 vertexlist= qh vertex_list;
866 else if (facetlist == qh newfacet_list)
867 vertexlist= qh newvertex_list;
868 else
869 vertexlist= NULL;
870 FORALLvertex_(vertexlist) {
871 vertex->seen= False;
872 vertex->visitid= 0;
873 }
874 FORALLfacet_(facetlist) {
875 if (facet->visible)
876 continue;
877 if (facet->simplicial)
878 numridges += qh hull_dim;
879 else
880 numridges += qh_setsize (facet->ridges);
881 FOREACHvertex_(facet->vertices) {
882 vertex->visitid++;
883 if (!vertex->seen) {
884 vertex->seen= True;
885 numvertices++;
886 if (qh_pointid (vertex->point) == -1) {
887 fprintf (qh ferr, "qhull internal error (qh_checkpolygon): unknown point %p for vertex v%d first_point %p\n",
888 vertex->point, vertex->id, qh first_point);
889 waserror= True;
890 }
891 }
892 }
893 }
894 qh vertex_visit += numfacets;
895 if (facetlist == qh facet_list) {
896 if (numfacets != qh num_facets - qh num_visible) {
897 fprintf(qh ferr, "qhull internal error (qh_checkpolygon): actual number of facets is %d, cumulative facet count is %d - %d visible facets\n",
898 numfacets, qh num_facets, qh num_visible);
899 waserror= True;
900 }
901 qh vertex_visit++;
902 if (qh VERTEXneighbors) {
903 FORALLvertices {
904 qh_setcheck (vertex->neighbors, "neighbors for v", vertex->id);
905 if (vertex->deleted)
906 continue;
907 totvneighbors += qh_setsize (vertex->neighbors);
908 }
909 FORALLfacet_(facetlist)
910 totvertices += qh_setsize (facet->vertices);
911 if (totvneighbors != totvertices) {
912 fprintf(qh ferr, "qhull internal error (qh_checkpolygon): vertex neighbors inconsistent. Totvneighbors %d, totvertices %d\n",
913 totvneighbors, totvertices);
914 waserror= True;
915 }
916 }
917 if (numvertices != qh num_vertices - qh_setsize(qh del_vertices)) {
918 fprintf(qh ferr, "qhull internal error (qh_checkpolygon): actual number of vertices is %d, cumulative vertex count is %d\n",
919 numvertices, qh num_vertices - qh_setsize(qh del_vertices));
920 waserror= True;
921 }
922 if (qh hull_dim == 2 && numvertices != numfacets) {
923 fprintf (qh ferr, "qhull internal error (qh_checkpolygon): #vertices %d != #facets %d\n",
924 numvertices, numfacets);
925 waserror= True;
926 }
927 if (qh hull_dim == 3 && numvertices + numfacets - numridges/2 != 2) {
928 fprintf (qh ferr, "qhull warning: #vertices %d + #facets %d - #edges %d != 2\n\
929 A vertex appears twice in a edge list. May occur during merging.",
930 numvertices, numfacets, numridges/2);
931 /* occurs if lots of merging and a vertex ends up twice in an edge list. e.g., RBOX 1000 s W1e-13 t995849315 D2 | QHULL d Tc Tv */
932 }
933 }
934 if (waserror)
935 qh_errexit(qh_ERRqhull, NULL, NULL);
936 } /* checkpolygon */
937
938
939 /*-<a href="qh-poly.htm#TOC"
940 >-------------------------------</a><a name="checkvertex">-</a>
941
942 qh_checkvertex( vertex )
943 check vertex for consistency
944 checks vertex->neighbors
945
946 notes:
947 neighbors checked efficiently in checkpolygon
948 */
949 void qh_checkvertex (vertexT *vertex) {
950 boolT waserror= False;
951 facetT *neighbor, **neighborp, *errfacet=NULL;
952
953 if (qh_pointid (vertex->point) == -1) {
954 fprintf (qh ferr, "qhull internal error (qh_checkvertex): unknown point id %p\n", vertex->point);
955 waserror= True;
956 }
957 if (vertex->id >= qh vertex_id) {
958 fprintf (qh ferr, "qhull internal error (qh_checkvertex): unknown vertex id %d\n", vertex->id);
959 waserror= True;
960 }
961 if (!waserror && !vertex->deleted) {
962 if (qh_setsize (vertex->neighbors)) {
963 FOREACHneighbor_(vertex) {
964 if (!qh_setin (neighbor->vertices, vertex)) {
965 fprintf (qh ferr, "qhull internal error (qh_checkvertex): neighbor f%d does not contain v%d\n", neighbor->id, vertex->id);
966 errfacet= neighbor;
967 waserror= True;
968 }
969 }
970 }
971 }
972 if (waserror) {
973 qh_errprint ("ERRONEOUS", NULL, NULL, NULL, vertex);
974 qh_errexit (qh_ERRqhull, errfacet, NULL);
975 }
976 } /* checkvertex */
977
978 /*-<a href="qh-poly.htm#TOC"
979 >-------------------------------</a><a name="clearcenters">-</a>
980
981 qh_clearcenters( type )
982 clear old data from facet->center
983
984 notes:
985 sets new centertype
986 nop if CENTERtype is the same
987 */
988 void qh_clearcenters (qh_CENTER type) {
989 facetT *facet;
990
991 if (qh CENTERtype != type) {
992 FORALLfacets {
993 if (qh CENTERtype == qh_ASvoronoi){
994 if (facet->center) {
995 qh_memfree (facet->center, qh center_size);
996 facet->center= NULL;
997 }
998 }else /* qh CENTERtype == qh_AScentrum */ {
999 if (facet->center) {
1000 qh_memfree (facet->center, qh normal_size);
1001 facet->center= NULL;
1002 }
1003 }
1004 }
1005 qh CENTERtype= type;
1006 }
1007 trace2((qh ferr, "qh_clearcenters: switched to center type %d\n", type));
1008 } /* clearcenters */
1009
1010 /*-<a href="qh-poly.htm#TOC"
1011 >-------------------------------</a><a name="createsimplex">-</a>
1012
1013 qh_createsimplex( vertices )
1014 creates a simplex from a set of vertices
1015
1016 returns:
1017 initializes qh.facet_list to the simplex
1018 initializes qh.newfacet_list, .facet_tail
1019 initializes qh.vertex_list, .newvertex_list, .vertex_tail
1020
1021 design:
1022 initializes lists
1023 for each vertex
1024 create a new facet
1025 for each new facet
1026 create its neighbor set
1027 */
1028 void qh_createsimplex(setT *vertices) {
1029 facetT *facet= NULL, *newfacet;
1030 boolT toporient= True;
1031 int vertex_i, vertex_n, nth;
1032 setT *newfacets= qh_settemp (qh hull_dim+1);
1033 vertexT *vertex;
1034
1035 qh facet_list= qh newfacet_list= qh facet_tail= qh_newfacet();
1036 qh num_facets= qh num_vertices= qh num_visible= 0;
1037 qh vertex_list= qh newvertex_list= qh vertex_tail= qh_newvertex(NULL);
1038 FOREACHvertex_i_(vertices) {
1039 newfacet= qh_newfacet();
1040 newfacet->vertices= qh_setnew_delnthsorted (vertices, vertex_n,
1041 vertex_i, 0);
1042 newfacet->toporient= toporient;
1043 qh_appendfacet(newfacet);
1044 newfacet->newfacet= True;
1045 qh_appendvertex (vertex);
1046 qh_setappend (&newfacets, newfacet);
1047 toporient ^= True;
1048 }
1049 FORALLnew_facets {
1050 nth= 0;
1051 FORALLfacet_(qh newfacet_list) {
1052 if (facet != newfacet)
1053 SETelem_(newfacet->neighbors, nth++)= facet;
1054 }
1055 qh_settruncate (newfacet->neighbors, qh hull_dim);
1056 }
1057 qh_settempfree (&newfacets);
1058 trace1((qh ferr, "qh_createsimplex: created simplex\n"));
1059 } /* createsimplex */
1060
1061 /*-<a href="qh-poly.htm#TOC"
1062 >-------------------------------</a><a name="delridge">-</a>
1063
1064 qh_delridge( ridge )
1065 deletes ridge from data structures it belongs to
1066 frees up its memory
1067
1068 notes:
1069 in merge.c, caller sets vertex->delridge for each vertex
1070 ridges also freed in qh_freeqhull
1071 */
1072 void qh_delridge(ridgeT *ridge) {
1073 void **freelistp; /* used !qh_NOmem */
1074
1075 qh_setdel(ridge->top->ridges, ridge);
1076 qh_setdel(ridge->bottom->ridges, ridge);
1077 qh_setfree(&(ridge->vertices));
1078 qh_memfree_(ridge, sizeof(ridgeT), freelistp);
1079 } /* delridge */
1080
1081
1082 /*-<a href="qh-poly.htm#TOC"
1083 >-------------------------------</a><a name="delvertex">-</a>
1084
1085 qh_delvertex( vertex )
1086 deletes a vertex and frees its memory
1087
1088 notes:
1089 assumes vertex->adjacencies have been updated if needed
1090 unlinks from vertex_list
1091 */
1092 void qh_delvertex (vertexT *vertex) {
1093
1094 if (vertex == qh tracevertex)
1095 qh tracevertex= NULL;
1096 qh_removevertex (vertex);
1097 qh_setfree (&vertex->neighbors);
1098 qh_memfree(vertex, sizeof(vertexT));
1099 } /* delvertex */
1100
1101
1102 /*-<a href="qh-poly.htm#TOC"
1103 >-------------------------------</a><a name="facet3vertex">-</a>
1104
1105 qh_facet3vertex( )
1106 return temporary set of 3-d vertices in qh_ORIENTclock order
1107
1108 design:
1109 if simplicial facet
1110 build set from facet->vertices with facet->toporient
1111 else
1112 for each ridge in order
1113 build set from ridge's vertices
1114 */
1115 setT *qh_facet3vertex (facetT *facet) {
1116 ridgeT *ridge, *firstridge;
1117 vertexT *vertex;
1118 int cntvertices, cntprojected=0;
1119 setT *vertices;
1120
1121 cntvertices= qh_setsize(facet->vertices);
1122 vertices= qh_settemp (cntvertices);
1123 if (facet->simplicial) {
1124 if (cntvertices != 3) {
1125 fprintf (qh ferr, "qhull internal error (qh_facet3vertex): only %d vertices for simplicial facet f%d\n",
1126 cntvertices, facet->id);
1127 qh_errexit(qh_ERRqhull, facet, NULL);
1128 }
1129 qh_setappend (&vertices, SETfirst_(facet->vertices));
1130 if (facet->toporient ^ qh_ORIENTclock)
1131 qh_setappend (&vertices, SETsecond_(facet->vertices));
1132 else
1133 qh_setaddnth (&vertices, 0, SETsecond_(facet->vertices));
1134 qh_setappend (&vertices, SETelem_(facet->vertices, 2));
1135 }else {
1136 ridge= firstridge= SETfirstt_(facet->ridges, ridgeT); /* no infinite */
1137 while ((ridge= qh_nextridge3d (ridge, facet, &vertex))) {
1138 qh_setappend (&vertices, vertex);
1139 if (++cntprojected > cntvertices || ridge == firstridge)
1140 break;
1141 }
1142 if (!ridge || cntprojected != cntvertices) {
1143 fprintf (qh ferr, "qhull internal error (qh_facet3vertex): ridges for facet %d don't match up. got at least %d\n",
1144 facet->id, cntprojected);
1145 qh_errexit(qh_ERRqhull, facet, ridge);
1146 }
1147 }
1148 return vertices;
1149 } /* facet3vertex */
1150
1151 /*-<a href="qh-poly.htm#TOC"
1152 >-------------------------------</a><a name="findbestfacet">-</a>
1153
1154 qh_findbestfacet( point, bestoutside, bestdist, isoutside )
1155 find facet that is furthest below a point
1156
1157 for Delaunay triangulations,
1158 Please use qh_setdelaunay() to lift point to paraboloid and scale by 'Qbb' if needed
1159 Do not use options 'Qbk', 'QBk', or 'QbB' since they scale the coordinates.
1160
1161 returns:
1162 if bestoutside is set (e.g., qh_ALL)
1163 returns best facet that is not upperdelaunay
1164 if Delaunay and inside, point is outside circumsphere of bestfacet
1165 else
1166 returns first facet below point
1167 if point is inside, returns nearest, !upperdelaunay facet
1168 distance to facet
1169 isoutside set if outside of facet
1170
1171 notes:
1172 For tricoplanar facets, this finds one of the tricoplanar facets closest
1173 to the point. For Delaunay triangulations, the point may be inside a
1174 different tricoplanar facet. See <a href="../html/qh-in.htm#findfacet">locate a facet with qh_findbestfacet()</a>
1175
1176 If inside, qh_findbestfacet performs an exhaustive search
1177 this may be too conservative. Sometimes it is clearly required.
1178
1179 qh_findbestfacet is not used by qhull.
1180 uses qh.visit_id and qh.coplanarset
1181
1182 see:
1183 <a href="geom.c#findbest">qh_findbest</a>
1184 */
1185 facetT *qh_findbestfacet (pointT *point, boolT bestoutside,
1186 realT *bestdist, boolT *isoutside) {
1187 facetT *bestfacet= NULL;
1188 int numpart, totpart= 0;
1189
1190 bestfacet= qh_findbest (point, qh facet_list,
1191 bestoutside, !qh_ISnewfacets, bestoutside /* qh_NOupper */,
1192 bestdist, isoutside, &totpart);
1193 if (*bestdist < -qh DISTround) {
1194 bestfacet= qh_findfacet_all (point, bestdist, isoutside, &numpart);
1195 totpart += numpart;
1196 if ((isoutside && bestoutside)
1197 || (!isoutside && bestfacet->upperdelaunay)) {
1198 bestfacet= qh_findbest (point, bestfacet,
1199 bestoutside, False, bestoutside,
1200 bestdist, isoutside, &totpart);
1201 totpart += numpart;
1202 }
1203 }
1204 trace3((qh ferr, "qh_findbestfacet: f%d dist %2.2g isoutside %d totpart %d\n", bestfacet->id, *bestdist, *isoutside, totpart));
1205 return bestfacet;
1206 } /* findbestfacet */
1207
1208 /*-<a href="qh-poly.htm#TOC"
1209 >-------------------------------</a><a name="findbestlower">-</a>
1210
1211 qh_findbestlower( facet, point, bestdist, numpart )
1212 returns best non-upper, non-flipped neighbor of facet for point
1213 if needed, searches vertex neighbors
1214
1215 returns:
1216 returns bestdist and updates numpart
1217
1218 notes:
1219 if Delaunay and inside, point is outside of circumsphere of bestfacet
1220 called by qh_findbest() for points above an upperdelaunay facet
1221
1222 */
1223 facetT *qh_findbestlower (facetT *upperfacet, pointT *point, realT *bestdistp, int *numpart) {
1224 facetT *neighbor, **neighborp, *bestfacet= NULL;
1225 realT bestdist= -REALmax/2 /* avoid underflow */;
1226 realT dist;
1227 vertexT *vertex;
1228
1229 zinc_(Zbestlower);
1230 FOREACHneighbor_(upperfacet) {
1231 if (neighbor->upperdelaunay || neighbor->flipped)
1232 continue;
1233 (*numpart)++;
1234 qh_distplane (point, neighbor, &dist);
1235 if (dist > bestdist) {
1236 bestfacet= neighbor;
1237 bestdist= dist;
1238 }
1239 }
1240 if (!bestfacet) {
1241 zinc_(Zbestlowerv);
1242 /* rarely called, numpart does not count nearvertex computations */
1243 vertex= qh_nearvertex (upperfacet, point, &dist);
1244 qh_vertexneighbors();
1245 FOREACHneighbor_(vertex) {
1246 if (neighbor->upperdelaunay || neighbor->flipped)
1247 continue;
1248 (*numpart)++;
1249 qh_distplane (point, neighbor, &dist);
1250 if (dist > bestdist) {
1251 bestfacet= neighbor;
1252 bestdist= dist;
1253 }
1254 }
1255 }
1256 if (!bestfacet) {
1257 fprintf(qh ferr, "\n\
1258 qh_findbestlower: all neighbors of facet %d are flipped or upper Delaunay.\n\
1259 Please report this error to qhull_bug@qhull.org with the input and all of the output.\n",
1260 upperfacet->id);
1261 qh_errexit (qh_ERRqhull, upperfacet, NULL);
1262 }
1263 *bestdistp= bestdist;
1264 trace3((qh ferr, "qh_findbestlower: f%d dist %2.2g for f%d p%d\n", bestfacet->id, bestdist, upperfacet->id, qh_pointid(point)));
1265 return bestfacet;
1266 } /* findbestlower */
1267
1268 /*-<a href="qh-poly.htm#TOC"
1269 >-------------------------------</a><a name="findfacet_all">-</a>
1270
1271 qh_findfacet_all( point, bestdist, isoutside, numpart )
1272 exhaustive search for facet below a point
1273
1274 for Delaunay triangulations,
1275 Please use qh_setdelaunay() to lift point to paraboloid and scale by 'Qbb' if needed
1276 Do not use options 'Qbk', 'QBk', or 'QbB' since they scale the coordinates.
1277
1278 returns:
1279 returns first facet below point
1280 if point is inside,
1281 returns nearest facet
1282 distance to facet
1283 isoutside if point is outside of the hull
1284 number of distance tests
1285 */
1286 facetT *qh_findfacet_all (pointT *point, realT *bestdist, boolT *isoutside,
1287 int *numpart) {
1288 facetT *bestfacet= NULL, *facet;
1289 realT dist;
1290 int totpart= 0;
1291
1292 *bestdist= REALmin;
1293 *isoutside= False;
1294 FORALLfacets {
1295 if (facet->flipped || !facet->normal)
1296 continue;
1297 totpart++;
1298 qh_distplane (point, facet, &dist);
1299 if (dist > *bestdist) {
1300 *bestdist= dist;
1301 bestfacet= facet;
1302 if (dist > qh MINoutside) {
1303 *isoutside= True;
1304 break;
1305 }
1306 }
1307 }
1308 *numpart= totpart;
1309 trace3((qh ferr, "qh_findfacet_all: f%d dist %2.2g isoutside %d totpart %d\n", getid_(bestfacet), *bestdist, *isoutside, totpart));
1310 return bestfacet;
1311 } /* findfacet_all */
1312
1313 /*-<a href="qh-poly.htm#TOC"
1314 >-------------------------------</a><a name="findgood">-</a>
1315
1316 qh_findgood( facetlist, goodhorizon )
1317 identify good facets for qh.PRINTgood
1318 if qh.GOODvertex>0
1319 facet includes point as vertex
1320 if !match, returns goodhorizon
1321 inactive if qh.MERGING
1322 if qh.GOODpoint
1323 facet is visible or coplanar (>0) or not visible (<0)
1324 if qh.GOODthreshold
1325 facet->normal matches threshold
1326 if !goodhorizon and !match,
1327 selects facet with closest angle
1328 sets GOODclosest
1329
1330 returns:
1331 number of new, good facets found
1332 determines facet->good
1333 may update qh.GOODclosest
1334
1335 notes:
1336 qh_findgood_all further reduces the good region
1337
1338 design:
1339 count good facets
1340 mark good facets for qh.GOODpoint
1341 mark good facets for qh.GOODthreshold
1342 if necessary
1343 update qh.GOODclosest
1344 */
1345 int qh_findgood (facetT *facetlist, int goodhorizon) {
1346 facetT *facet, *bestfacet= NULL;
1347 realT angle, bestangle= REALmax, dist;
1348 int numgood=0;
1349
1350 FORALLfacet_(facetlist) {
1351 if (facet->good)
1352 numgood++;
1353 }
1354 if (qh GOODvertex>0 && !qh MERGING) {
1355 FORALLfacet_(facetlist) {
1356 if (!qh_isvertex (qh GOODvertexp, facet->vertices)) {
1357 facet->good= False;
1358 numgood--;
1359 }
1360 }
1361 }
1362 if (qh GOODpoint && numgood) {
1363 FORALLfacet_(facetlist) {
1364 if (facet->good && facet->normal) {
1365 zinc_(Zdistgood);
1366 qh_distplane (qh GOODpointp, facet, &dist);
1367 if ((qh GOODpoint > 0) ^ (dist > 0.0)) {
1368 facet->good= False;
1369 numgood--;
1370 }
1371 }
1372 }
1373 }
1374 if (qh GOODthreshold && (numgood || goodhorizon || qh GOODclosest)) {
1375 FORALLfacet_(facetlist) {
1376 if (facet->good && facet->normal) {
1377 if (!qh_inthresholds (facet->normal, &angle)) {
1378 facet->good= False;
1379 numgood--;
1380 if (angle < bestangle) {
1381 bestangle= angle;
1382 bestfacet= facet;
1383 }
1384 }
1385 }
1386 }
1387 if (!numgood && (!goodhorizon || qh GOODclosest)) {
1388 if (qh GOODclosest) {
1389 if (qh GOODclosest->visible)
1390 qh GOODclosest= NULL;
1391 else {
1392 qh_inthresholds (qh GOODclosest->normal, &angle);
1393 if (angle < bestangle)
1394 bestfacet= qh GOODclosest;
1395 }
1396 }
1397 if (bestfacet && bestfacet != qh GOODclosest) {
1398 if (qh GOODclosest)
1399 qh GOODclosest->good= False;
1400 qh GOODclosest= bestfacet;
1401 bestfacet->good= True;
1402 numgood++;
1403 trace2((qh ferr, "qh_findgood: f%d is closest (%2.2g) to thresholds\n", bestfacet->id, bestangle));
1404 return numgood;
1405 }
1406 }else if (qh GOODclosest) { /* numgood > 0 */
1407 qh GOODclosest->good= False;
1408 qh GOODclosest= NULL;
1409 }
1410 }
1411 zadd_(Zgoodfacet, numgood);
1412 trace2((qh ferr, "qh_findgood: found %d good facets with %d good horizon\n", numgood, goodhorizon));
1413 if (!numgood && qh GOODvertex>0 && !qh MERGING)
1414 return goodhorizon;
1415 return numgood;
1416 } /* findgood */
1417
1418 /*-<a href="qh-poly.htm#TOC"
1419 >-------------------------------</a><a name="findgood_all">-</a>
1420
1421 qh_findgood_all( facetlist )
1422 apply other constraints for good facets (used by qh.PRINTgood)
1423 if qh.GOODvertex
1424 facet includes (>0) or doesn't include (<0) point as vertex
1425 if last good facet and ONLYgood, prints warning and continues
1426 if qh.SPLITthresholds
1427 facet->normal matches threshold, or if none, the closest one
1428 calls qh_findgood
1429 nop if good not used
1430
1431 returns:
1432 clears facet->good if not good
1433 sets qh.num_good
1434
1435 notes:
1436 this is like qh_findgood but more restrictive
1437
1438 design:
1439 uses qh_findgood to mark good facets
1440 marks facets for qh.GOODvertex
1441 marks facets for qh.SPLITthreholds
1442 */
1443 void qh_findgood_all (facetT *facetlist) {
1444 facetT *facet, *bestfacet=NULL;
1445 realT angle, bestangle= REALmax;
1446 int numgood=0, startgood;
1447
1448 if (!qh GOODvertex && !qh GOODthreshold && !qh GOODpoint
1449 && !qh SPLITthresholds)
1450 return;
1451 if (!qh ONLYgood)
1452 qh_findgood (qh facet_list, 0);
1453 FORALLfacet_(facetlist) {
1454 if (facet->good)
1455 numgood++;
1456 }
1457 if (qh GOODvertex <0 || (qh GOODvertex > 0 && qh MERGING)) {
1458 FORALLfacet_(facetlist) {
1459 if (facet->good && ((qh GOODvertex > 0) ^ !!qh_isvertex (qh GOODvertexp, facet->vertices))) {
1460 if (!--numgood) {
1461 if (qh ONLYgood) {
1462 fprintf (qh ferr, "qhull warning: good vertex p%d does not match last good facet f%d. Ignored.\n",
1463 qh_pointid(qh GOODvertexp), facet->id);
1464 return;
1465 }else if (qh GOODvertex > 0)
1466 fprintf (qh ferr, "qhull warning: point p%d is not a vertex ('QV%d').\n",
1467 qh GOODvertex-1, qh GOODvertex-1);
1468 else
1469 fprintf (qh ferr, "qhull warning: point p%d is a vertex for every facet ('QV-%d').\n",
1470 -qh GOODvertex - 1, -qh GOODvertex - 1);
1471 }
1472 facet->good= False;
1473 }
1474 }
1475 }
1476 startgood= numgood;
1477 if (qh SPLITthresholds) {
1478 FORALLfacet_(facetlist) {
1479 if (facet->good) {
1480 if (!qh_inthresholds (facet->normal, &angle)) {
1481 facet->good= False;
1482 numgood--;
1483 if (angle < bestangle) {
1484 bestangle= angle;
1485 bestfacet= facet;
1486 }
1487 }
1488 }
1489 }
1490 if (!numgood && bestfacet) {
1491 bestfacet->good= True;
1492 numgood++;
1493 trace0((qh ferr, "qh_findgood_all: f%d is closest (%2.2g) to thresholds\n", bestfacet->id, bestangle));
1494 return;
1495 }
1496 }
1497 qh num_good= numgood;
1498 trace0((qh ferr, "qh_findgood_all: %d good facets remain out of %d facets\n", numgood, startgood));
1499 } /* findgood_all */
1500
1501 /*-<a href="qh-poly.htm#TOC"
1502 >-------------------------------</a><a name="furthestnext">-</a>
1503
1504 qh_furthestnext()
1505 set qh.facet_next to facet with furthest of all furthest points
1506 searches all facets on qh.facet_list
1507
1508 notes:
1509 this may help avoid precision problems
1510 */
1511 void qh_furthestnext (void /* qh facet_list */) {
1512 facetT *facet, *bestfacet= NULL;
1513 realT dist, bestdist= -REALmax;
1514
1515 FORALLfacets {
1516 if (facet->outsideset) {
1517 #if qh_COMPUTEfurthest
1518 pointT *furthest;
1519 furthest= (pointT*)qh_setlast (facet->outsideset);
1520 zinc_(Zcomputefurthest);
1521 qh_distplane (furthest, facet, &dist);
1522 #else
1523 dist= facet->furthestdist;
1524 #endif
1525 if (dist > bestdist) {
1526 bestfacet= facet;
1527 bestdist= dist;
1528 }
1529 }
1530 }
1531 if (bestfacet) {
1532 qh_removefacet (bestfacet);
1533 qh_prependfacet (bestfacet, &qh facet_next);
1534 trace1((qh ferr, "qh_furthestnext: made f%d next facet (dist %.2g)\n", bestfacet->id, bestdist));
1535 }
1536 } /* furthestnext */
1537
1538 /*-<a href="qh-poly.htm#TOC"
1539 >-------------------------------</a><a name="furthestout">-</a>
1540
1541 qh_furthestout( facet )
1542 make furthest outside point the last point of outsideset
1543
1544 returns:
1545 updates facet->outsideset
1546 clears facet->notfurthest
1547 sets facet->furthestdist
1548
1549 design:
1550 determine best point of outsideset
1551 make it the last point of outsideset
1552 */
1553 void qh_furthestout (facetT *facet) {
1554 pointT *point, **pointp, *bestpoint= NULL;
1555 realT dist, bestdist= -REALmax;
1556
1557 FOREACHpoint_(facet->outsideset) {
1558 qh_distplane (point, facet, &dist);
1559 zinc_(Zcomputefurthest);
1560 if (dist > bestdist) {
1561 bestpoint= point;
1562 bestdist= dist;
1563 }
1564 }
1565 if (bestpoint) {
1566 qh_setdel (facet->outsideset, point);
1567 qh_setappend (&facet->outsideset, point);
1568 #if !qh_COMPUTEfurthest
1569 facet->furthestdist= bestdist;
1570 #endif
1571 }
1572 facet->notfurthest= False;
1573 trace3((qh ferr, "qh_furthestout: p%d is furthest outside point of f%d\n", qh_pointid (point), facet->id));
1574 } /* furthestout */
1575
1576
1577 /*-<a href="qh-qhull.htm#TOC"
1578 >-------------------------------</a><a name="infiniteloop">-</a>
1579
1580 qh_infiniteloop( facet )
1581 report infinite loop error due to facet
1582 */
1583 void qh_infiniteloop (facetT *facet) {
1584
1585 fprintf (qh ferr, "qhull internal error (qh_infiniteloop): potential infinite loop detected\n");
1586 qh_errexit (qh_ERRqhull, facet, NULL);
1587 } /* qh_infiniteloop */
1588
1589 /*-<a href="qh-poly.htm#TOC"
1590 >-------------------------------</a><a name="initbuild">-</a>
1591
1592 qh_initbuild()
1593 initialize hull and outside sets with point array
1594 qh.FIRSTpoint/qh.NUMpoints is point array
1595 if qh.GOODpoint
1596 adds qh.GOODpoint to initial hull
1597
1598 returns:
1599 qh_facetlist with initial hull
1600 points partioned into outside sets, coplanar sets, or inside
1601 initializes qh.GOODpointp, qh.GOODvertexp,
1602
1603 design:
1604 initialize global variables used during qh_buildhull
1605 determine precision constants and points with max/min coordinate values
1606 if qh.SCALElast, scale last coordinate (for 'd')
1607 build initial simplex
1608 partition input points into facets of initial simplex
1609 set up lists
1610 if qh.ONLYgood
1611 check consistency
1612 add qh.GOODvertex if defined
1613 */
1614 void qh_initbuild( void) {
1615 setT *maxpoints, *vertices;
1616 facetT *facet;
1617 int i, numpart;
1618 realT dist;
1619 boolT isoutside;
1620
1621 qh furthest_id= -1;
1622 qh lastreport= 0;
1623 qh facet_id= qh vertex_id= qh ridge_id= 0;
1624 qh visit_id= qh vertex_visit= 0;
1625 qh maxoutdone= False;
1626
1627 if (qh GOODpoint > 0)
1628 qh GOODpointp= qh_point (qh GOODpoint-1);
1629 else if (qh GOODpoint < 0)
1630 qh GOODpointp= qh_point (-qh GOODpoint-1);
1631 if (qh GOODvertex > 0)
1632 qh GOODvertexp= qh_point (qh GOODvertex-1);
1633 else if (qh GOODvertex < 0)
1634 qh GOODvertexp= qh_point (-qh GOODvertex-1);
1635 if ((qh GOODpoint
1636 && (qh GOODpointp < qh first_point /* also catches !GOODpointp */
1637 || qh GOODpointp > qh_point (qh num_points-1)))
1638 || (qh GOODvertex
1639 && (qh GOODvertexp < qh first_point /* also catches !GOODvertexp */
1640 || qh GOODvertexp > qh_point (qh num_points-1)))) {
1641 fprintf (qh ferr, "qhull input error: either QGn or QVn point is > p%d\n",
1642 qh num_points-1);
1643 qh_errexit (qh_ERRinput, NULL, NULL);
1644 }
1645 maxpoints= qh_maxmin(qh first_point, qh num_points, qh hull_dim);
1646 if (qh SCALElast)
1647 qh_scalelast (qh first_point, qh num_points, qh hull_dim,
1648 qh MINlastcoord, qh MAXlastcoord, qh MAXwidth);
1649 qh_detroundoff();
1650 if (qh DELAUNAY && qh upper_threshold[qh hull_dim-1] > REALmax/2
1651 && qh lower_threshold[qh hull_dim-1] < -REALmax/2) {
1652 for (i= qh_PRINTEND; i--; ) {
1653 if (qh PRINTout[i] == qh_PRINTgeom && qh DROPdim < 0
1654 && !qh GOODthreshold && !qh SPLITthresholds)
1655 break; /* in this case, don't set upper_threshold */
1656 }
1657 if (i < 0) {
1658 if (qh UPPERdelaunay) { /* matches qh.upperdelaunay in qh_setfacetplane */
1659 qh lower_threshold[qh hull_dim-1]= qh ANGLEround * qh_ZEROdelaunay;
1660 qh GOODthreshold= True;
1661 }else {
1662 qh upper_threshold[qh hull_dim-1]= -qh ANGLEround * qh_ZEROdelaunay;
1663 if (!qh GOODthreshold)
1664 qh SPLITthresholds= True; /* build upper-convex hull even if Qg */
1665 /* qh_initqhull_globals errors if Qg without Pdk/etc. */
1666 }
1667 }
1668 }
1669 vertices= qh_initialvertices(qh hull_dim, maxpoints, qh first_point, qh num_points);
1670 qh_initialhull (vertices); /* initial qh facet_list */
1671 qh_partitionall (vertices, qh first_point, qh num_points);
1672 if (qh PRINToptions1st || qh TRACElevel || qh IStracing) {
1673 if (qh TRACElevel || qh IStracing)
1674 fprintf (qh ferr, "\nTrace level %d for %s | %s\n",
1675 qh IStracing ? qh IStracing : qh TRACElevel, qh rbox_command, qh qhull_command);
1676 fprintf (qh ferr, "Options selected for Qhull %s:\n%s\n", qh_version, qh qhull_options);
1677 }
1678 qh_resetlists (False, qh_RESETvisible /*qh visible_list newvertex_list newfacet_list */);
1679 qh facet_next= qh facet_list;
1680 qh_furthestnext (/* qh facet_list */);
1681 if (qh PREmerge) {
1682 qh cos_max= qh premerge_cos;
1683 qh centrum_radius= qh premerge_centrum;
1684 }
1685 if (qh ONLYgood) {
1686 if (qh GOODvertex > 0 && qh MERGING) {
1687 fprintf (qh ferr, "qhull input error: 'Qg QVn' (only good vertex) does not work with merging.\nUse 'QJ' to joggle the input or 'Q0' to turn off merging.\n");
1688 qh_errexit (qh_ERRinput, NULL, NULL);
1689 }
1690 if (!(qh GOODthreshold || qh GOODpoint
1691 || (!qh MERGEexact && !qh PREmerge && qh GOODvertexp))) {
1692 fprintf (qh ferr, "qhull input error: 'Qg' (ONLYgood) needs a good threshold ('Pd0D0'), a\n\
1693 good point (QGn or QG-n), or a good vertex with 'QJ' or 'Q0' (QVn).\n");
1694 qh_errexit (qh_ERRinput, NULL, NULL);
1695 }
1696 if (qh GOODvertex > 0 && !qh MERGING /* matches qh_partitionall */
1697 && !qh_isvertex (qh GOODvertexp, vertices)) {
1698 facet= qh_findbestnew (qh GOODvertexp, qh facet_list,
1699 &dist, !qh_ALL, &isoutside, &numpart);
1700 zadd_(Zdistgood, numpart);
1701 if (!isoutside) {
1702 fprintf (qh ferr, "qhull input error: point for QV%d is inside initial simplex. It can not be made a vertex.\n",
1703 qh_pointid(qh GOODvertexp));
1704 qh_errexit (qh_ERRinput, NULL, NULL);
1705 }
1706 if (!qh_addpoint (qh GOODvertexp, facet, False)) {
1707 qh_settempfree(&vertices);
1708 qh_settempfree(&maxpoints);
1709 return;
1710 }
1711 }
1712 qh_findgood (qh facet_list, 0);
1713 }
1714 qh_settempfree(&vertices);
1715 qh_settempfree(&maxpoints);
1716 trace1((qh ferr, "qh_initbuild: initial hull created and points partitioned\n"));
1717 } /* initbuild */
1718
1719 /*-<a href="qh-poly.htm#TOC"
1720 >-------------------------------</a><a name="initialhull">-</a>
1721
1722 qh_initialhull( vertices )
1723 constructs the initial hull as a DIM3 simplex of vertices
1724
1725 design:
1726 creates a simplex (initializes lists)
1727 determines orientation of simplex
1728 sets hyperplanes for facets
1729 doubles checks orientation (in case of axis-parallel facets with Gaussian elimination)
1730 checks for flipped facets and qh.NARROWhull
1731 checks the result
1732 */
1733 void qh_initialhull(setT *vertices) {
1734 facetT *facet, *firstfacet, *neighbor, **neighborp;
1735 realT dist, angle, minangle= REALmax;
1736 #ifndef qh_NOtrace
1737 int k;
1738 #endif
1739
1740 qh_createsimplex(vertices); /* qh facet_list */
1741 qh_resetlists (False, qh_RESETvisible);
1742 qh facet_next= qh facet_list; /* advance facet when processed */
1743 qh interior_point= qh_getcenter(vertices);
1744 firstfacet= qh facet_list;
1745 qh_setfacetplane(firstfacet);
1746 zinc_(Znumvisibility); /* needs to be in printsummary */
1747 qh_distplane(qh interior_point, firstfacet, &dist);
1748 if (dist > 0) {
1749 FORALLfacets
1750 facet->toporient ^= True;
1751 }
1752 FORALLfacets
1753 qh_setfacetplane(facet);
1754 FORALLfacets {
1755 if (!qh_checkflipped (facet, NULL, qh_ALL)) {/* due to axis-parallel facet */
1756 trace1((qh ferr, "qh_initialhull: initial orientation incorrect. Correct all facets\n"));
1757 facet->flipped= False;
1758 FORALLfacets {
1759 facet->toporient ^= True;
1760 qh_orientoutside (facet);
1761 }
1762 break;
1763 }
1764 }
1765 FORALLfacets {
1766 if (!qh_checkflipped (facet, NULL, !qh_ALL)) { /* can happen with 'R0.1' */
1767 qh_precision ("initial facet is coplanar with interior point");
1768 fprintf (qh ferr, "qhull precision error: initial facet %d is coplanar with the interior point\n",
1769 facet->id);
1770 qh_errexit (qh_ERRsingular, facet, NULL);
1771 }
1772 FOREACHneighbor_(facet) {
1773 angle= qh_getangle (facet->normal, neighbor->normal);
1774 minimize_( minangle, angle);
1775 }
1776 }
1777 if (minangle < qh_MAXnarrow && !qh NOnarrow) {
1778 realT diff= 1.0 + minangle;
1779
1780 qh NARROWhull= True;
1781 qh_option ("_narrow-hull", NULL, &diff);
1782 if (minangle < qh_WARNnarrow && !qh RERUN && qh PRINTprecision)
1783 fprintf (qh ferr, "qhull precision warning: \n\
1784 The initial hull is narrow (cosine of min. angle is %.16f).\n\
1785 A coplanar point may lead to a wide facet. Options 'QbB' (scale to unit box)\n\
1786 or 'Qbb' (scale last coordinate) may remove this warning. Use 'Pp' to skip\n\
1787 this warning. See 'Limitations' in qh-impre.htm.\n",
1788 -minangle); /* convert from angle between normals to angle between facets */
1789 }
1790 zzval_(Zprocessed)= qh hull_dim+1;
1791 qh_checkpolygon (qh facet_list);
1792 qh_checkconvex(qh facet_list, qh_DATAfault);
1793 #ifndef qh_NOtrace
1794 if (qh IStracing >= 1) {
1795 fprintf(qh ferr, "qh_initialhull: simplex constructed, interior point:");
1796 for (k=0; k < qh hull_dim; k++)
1797 fprintf (qh ferr, " %6.4g", qh interior_point[k]);
1798 fprintf (qh ferr, "\n");
1799 }
1800 #endif
1801 } /* initialhull */
1802
1803 /*-<a href="qh-poly.htm#TOC"
1804 >-------------------------------</a><a name="initialvertices">-</a>
1805
1806 qh_initialvertices( dim, maxpoints, points, numpoints )
1807 determines a non-singular set of initial vertices
1808 maxpoints may include duplicate points
1809
1810 returns:
1811 temporary set of dim+1 vertices in descending order by vertex id
1812 if qh.RANDOMoutside && !qh.ALLpoints
1813 picks random points
1814 if dim >= qh_INITIALmax,
1815 uses min/max x and max points with non-zero determinants
1816
1817 notes:
1818 unless qh.ALLpoints,
1819 uses maxpoints as long as determinate is non-zero
1820 */
1821 setT *qh_initialvertices(int dim, setT *maxpoints, pointT *points, int numpoints) {
1822 pointT *point, **pointp;
1823 setT *vertices, *simplex, *tested;
1824 realT randr;
1825 int index, point_i, point_n, k;
1826 boolT nearzero= False;
1827
1828 vertices= qh_settemp (dim + 1);
1829 simplex= qh_settemp (dim+1);
1830 if (qh ALLpoints)
1831 qh_maxsimplex (dim, NULL, points, numpoints, &simplex);
1832 else if (qh RANDOMoutside) {
1833 while (qh_setsize (simplex) != dim+1) {
1834 randr= qh_RANDOMint;
1835 randr= randr/(qh_RANDOMmax+1);
1836 index= (int)floor(qh num_points * randr);
1837 while (qh_setin (simplex, qh_point (index))) {
1838 index++; /* in case qh_RANDOMint always returns the same value */
1839 index= index < qh num_points ? index : 0;
1840 }
1841 qh_setappend (&simplex, qh_point (index));
1842 }
1843 }else if (qh hull_dim >= qh_INITIALmax) {
1844 tested= qh_settemp (dim+1);
1845 qh_setappend (&simplex, SETfirst_(maxpoints)); /* max and min X coord */
1846 qh_setappend (&simplex, SETsecond_(maxpoints));
1847 qh_maxsimplex (fmin_(qh_INITIALsearch, dim), maxpoints, points, numpoints, &simplex);
1848 k= qh_setsize (simplex);
1849 FOREACHpoint_i_(maxpoints) {
1850 if (point_i & 0x1) { /* first pick up max. coord. points */
1851 if (!qh_setin (simplex, point) && !qh_setin (tested, point)){
1852 qh_detsimplex(point, simplex, k, &nearzero);
1853 if (nearzero)
1854 qh_setappend (&tested, point);
1855 else {
1856 qh_setappend (&simplex, point);
1857 if (++k == dim) /* use search for last point */
1858 break;
1859 }
1860 }
1861 }
1862 }
1863 while (k != dim && (point= (pointT*)qh_setdellast (maxpoints))) {
1864 if (!qh_setin (simplex, point) && !qh_setin (tested, point)){
1865 qh_detsimplex (point, simplex, k, &nearzero);
1866 if (nearzero)
1867 qh_setappend (&tested, point);
1868 else {
1869 qh_setappend (&simplex, point);
1870 k++;
1871 }
1872 }
1873 }
1874 index= 0;
1875 while (k != dim && (point= qh_point (index++))) {
1876 if (!qh_setin (simplex, point) && !qh_setin (tested, point)){
1877 qh_detsimplex (point, simplex, k, &nearzero);
1878 if (!nearzero){
1879 qh_setappend (&simplex, point);
1880 k++;
1881 }
1882 }
1883 }
1884 qh_settempfree (&tested);
1885 qh_maxsimplex (dim, maxpoints, points, numpoints, &simplex);
1886 }else
1887 qh_maxsimplex (dim, maxpoints, points, numpoints, &simplex);
1888 FOREACHpoint_(simplex)
1889 qh_setaddnth (&vertices, 0, qh_newvertex(point)); /* descending order */
1890 qh_settempfree (&simplex);
1891 return vertices;
1892 } /* initialvertices */
1893
1894
1895 /*-<a href="qh-poly.htm#TOC"
1896 >-------------------------------</a><a name="isvertex">-</a>
1897
1898 qh_isvertex( )
1899 returns vertex if point is in vertex set, else returns NULL
1900
1901 notes:
1902 for qh.GOODvertex
1903 */
1904 vertexT *qh_isvertex (pointT *point, setT *vertices) {
1905 vertexT *vertex, **vertexp;
1906
1907 FOREACHvertex_(vertices) {
1908 if (vertex->point == point)
1909 return vertex;
1910 }
1911 return NULL;
1912 } /* isvertex */
1913
1914 /*-<a href="qh-poly.htm#TOC"
1915 >-------------------------------</a><a name="makenewfacets">-</a>
1916
1917 qh_makenewfacets( point )
1918 make new facets from point and qh.visible_list
1919
1920 returns:
1921 qh.newfacet_list= list of new facets with hyperplanes and ->newfacet
1922 qh.newvertex_list= list of vertices in new facets with ->newlist set
1923
1924 if (qh.ONLYgood)
1925 newfacets reference horizon facets, but not vice versa
1926 ridges reference non-simplicial horizon ridges, but not vice versa
1927 does not change existing facets
1928 else
1929 sets qh.NEWfacets
1930 new facets attached to horizon facets and ridges
1931 for visible facets,
1932 visible->r.replace is corresponding new facet
1933
1934 see also:
1935 qh_makenewplanes() -- make hyperplanes for facets
1936 qh_attachnewfacets() -- attachnewfacets if not done here (qh ONLYgood)
1937 qh_matchnewfacets() -- match up neighbors
1938 qh_updatevertices() -- update vertex neighbors and delvertices
1939 qh_deletevisible() -- delete visible facets
1940 qh_checkpolygon() --check the result
1941 qh_triangulate() -- triangulate a non-simplicial facet
1942
1943 design:
1944 for each visible facet
1945 make new facets to its horizon facets
1946 update its f.replace
1947 clear its neighbor set
1948 */
1949 vertexT *qh_makenewfacets (pointT *point /*visible_list*/) {
1950 facetT *visible, *newfacet= NULL, *newfacet2= NULL, *neighbor, **neighborp;
1951 vertexT *apex;
1952 int numnew=0;
1953
1954 qh newfacet_list= qh facet_tail;
1955 qh newvertex_list= qh vertex_tail;
1956 apex= qh_newvertex(point);
1957 qh_appendvertex (apex);
1958 qh visit_id++;
1959 if (!qh ONLYgood)
1960 qh NEWfacets= True;
1961 FORALLvisible_facets {
1962 FOREACHneighbor_(visible)
1963 neighbor->seen= False;
1964 if (visible->ridges) {
1965 visible->visitid= qh visit_id;
1966 newfacet2= qh_makenew_nonsimplicial (visible, apex, &numnew);
1967 }
1968 if (visible->simplicial)
1969 newfacet= qh_makenew_simplicial (visible, apex, &numnew);
1970 if (!qh ONLYgood) {
1971 if (newfacet2) /* newfacet is null if all ridges defined */
1972 newfacet= newfacet2;
1973 if (newfacet)
1974 visible->f.replace= newfacet;
1975 else
1976 zinc_(Zinsidevisible);
1977 SETfirst_(visible->neighbors)= NULL;
1978 }
1979 }
1980 trace1((qh ferr, "qh_makenewfacets: created %d new facets from point p%d to horizon\n", numnew, qh_pointid(point)));
1981 if (qh IStracing >= 4)
1982 qh_printfacetlist (qh newfacet_list, NULL, qh_ALL);
1983 return apex;
1984 } /* makenewfacets */
1985
1986 /*-<a href="qh-poly.htm#TOC"
1987 >-------------------------------</a><a name="matchduplicates">-</a>
1988
1989 qh_matchduplicates( atfacet, atskip, hashsize, hashcount )
1990 match duplicate ridges in qh.hash_table for atfacet/atskip
1991 duplicates marked with ->dupridge and qh_DUPLICATEridge
1992
1993 returns:
1994 picks match with worst merge (min distance apart)
1995 updates hashcount
1996
1997 see also:
1998 qh_matchneighbor
1999
2000 notes:
2001
2002 design:
2003 compute hash value for atfacet and atskip
2004 repeat twice -- once to make best matches, once to match the rest
2005 for each possible facet in qh.hash_table
2006 if it is a matching facet and pass 2
2007 make match
2008 unless tricoplanar, mark match for merging (qh_MERGEridge)
2009 [e.g., tricoplanar RBOX s 1000 t993602376 | QHULL C-1e-3 d Qbb FA Qt]
2010 if it is a matching facet and pass 1
2011 test if this is a better match
2012 if pass 1,
2013 make best match (it will not be merged)
2014 */
2015 #ifndef qh_NOmerge
2016 void qh_matchduplicates (facetT *atfacet, int atskip, int hashsize, int *hashcount) {
2017 boolT same, ismatch;
2018 int hash, scan;
2019 facetT *facet, *newfacet, *maxmatch= NULL, *maxmatch2= NULL, *nextfacet;
2020 int skip, newskip, nextskip= 0, maxskip= 0, maxskip2= 0, makematch;
2021 realT maxdist= -REALmax, mindist, dist2, low, high;
2022
2023 hash= (int)qh_gethash (hashsize, atfacet->vertices, qh hull_dim, 1,
2024 SETelem_(atfacet->vertices, atskip));
2025 trace2((qh ferr, "qh_matchduplicates: find duplicate matches for f%d skip %d hash %d hashcount %d\n", atfacet->id, atskip, hash, *hashcount));
2026 for (makematch= 0; makematch < 2; makematch++) {
2027 qh visit_id++;
2028 for (newfacet= atfacet, newskip= atskip; newfacet; newfacet= nextfacet, newskip= nextskip) {
2029 zinc_(Zhashlookup);
2030 nextfacet= NULL;
2031 newfacet->visitid= qh visit_id;
2032 for (scan= hash; (facet= SETelemt_(qh hash_table, scan, facetT));
2033 scan= (++scan >= hashsize ? 0 : scan)) {
2034 if (!facet->dupridge || facet->visitid == qh visit_id)
2035 continue;
2036 zinc_(Zhashtests);
2037 if (qh_matchvertices (1, newfacet->vertices, newskip, facet->vertices, &skip, &same)) {
2038 ismatch= (same == (newfacet->toporient ^ facet->toporient));
2039 if (SETelemt_(facet->neighbors, skip, facetT) != qh_DUPLICATEridge) {
2040 if (!makematch) {
2041 fprintf (qh ferr, "qhull internal error (qh_matchduplicates): missing dupridge at f%d skip %d for new f%d skip %d hash %d\n",
2042 facet->id, skip, newfacet->id, newskip, hash);
2043 qh_errexit2 (qh_ERRqhull, facet, newfacet);
2044 }
2045 }else if (ismatch && makematch) {
2046 if (SETelemt_(newfacet->neighbors, newskip, facetT) == qh_DUPLICATEridge) {
2047 SETelem_(facet->neighbors, skip)= newfacet;
2048 if (newfacet->tricoplanar)
2049 SETelem_(newfacet->neighbors, newskip)= facet;
2050 else
2051 SETelem_(newfacet->neighbors, newskip)= qh_MERGEridge;
2052 *hashcount -= 2; /* removed two unmatched facets */
2053 trace4((qh ferr, "qh_matchduplicates: duplicate f%d skip %d matched with new f%d skip %d merge\n", facet->id, skip, newfacet->id, newskip));
2054 }
2055 }else if (ismatch) {
2056 mindist= qh_getdistance (facet, newfacet, &low, &high);
2057 dist2= qh_getdistance (newfacet, facet, &low, &high);
2058 minimize_(mindist, dist2);
2059 if (mindist > maxdist) {
2060 maxdist= mindist;
2061 maxmatch= facet;
2062 maxskip= skip;
2063 maxmatch2= newfacet;
2064 maxskip2= newskip;
2065 }
2066 trace3((qh ferr, "qh_matchduplicates: duplicate f%d skip %d new f%d skip %d at dist %2.2g, max is now f%d f%d\n", facet->id, skip, newfacet->id, newskip, mindist, maxmatch->id, maxmatch2->id));
2067 }else { /* !ismatch */
2068 nextfacet= facet;
2069 nextskip= skip;
2070 }
2071 }
2072 if (makematch && !facet
2073 && SETelemt_(facet->neighbors, skip, facetT) == qh_DUPLICATEridge) {
2074 fprintf (qh ferr, "qhull internal error (qh_matchduplicates): no MERGEridge match for duplicate f%d skip %d at hash %d\n",
2075 newfacet->id, newskip, hash);
2076 qh_errexit (qh_ERRqhull, newfacet, NULL);
2077 }
2078 }
2079 } /* end of for each new facet at hash */
2080 if (!makematch) {
2081 if (!maxmatch) {
2082 fprintf (qh ferr, "qhull internal error (qh_matchduplicates): no maximum match at duplicate f%d skip %d at hash %d\n",
2083 atfacet->id, atskip, hash);
2084 qh_errexit (qh_ERRqhull, atfacet, NULL);
2085 }
2086 SETelem_(maxmatch->neighbors, maxskip)= maxmatch2;
2087 SETelem_(maxmatch2->neighbors, maxskip2)= maxmatch;
2088 *hashcount -= 2; /* removed two unmatched facets */
2089 zzinc_(Zmultiridge);
2090 trace0((qh ferr, "qh_matchduplicates: duplicate f%d skip %d matched with new f%d skip %d keep\n", maxmatch->id, maxskip, maxmatch2->id, maxskip2));
2091 qh_precision ("ridge with multiple neighbors");
2092 if (qh IStracing >= 4)
2093 qh_errprint ("DUPLICATED/MATCH", maxmatch, maxmatch2, NULL, NULL);
2094 }
2095 }
2096 } /* matchduplicates */
2097
2098 /*-<a href="qh-poly.htm#TOC"
2099 >-------------------------------</a><a name="nearcoplanar">-</a>
2100
2101 qh_nearcoplanar()
2102 for all facets, remove near-inside points from facet->coplanarset</li>
2103 coplanar points defined by innerplane from qh_outerinner()
2104
2105 returns:
2106 if qh KEEPcoplanar && !qh KEEPinside
2107 facet->coplanarset only contains coplanar points
2108 if qh.JOGGLEmax
2109 drops inner plane by another qh.JOGGLEmax diagonal since a
2110 vertex could shift out while a coplanar point shifts in
2111
2112 notes:
2113 used for qh.PREmerge and qh.JOGGLEmax
2114 must agree with computation of qh.NEARcoplanar in qh_detroundoff()
2115 design:
2116 if not keeping coplanar or inside points
2117 free all coplanar sets
2118 else if not keeping both coplanar and inside points
2119 remove !coplanar or !inside points from coplanar sets
2120 */
2121 void qh_nearcoplanar ( void /* qh.facet_list */) {
2122 facetT *facet;
2123 pointT *point, **pointp;
2124 int numpart;
2125 realT dist, innerplane;
2126
2127 if (!qh KEEPcoplanar && !qh KEEPinside) {
2128 FORALLfacets {
2129 if (facet->coplanarset)
2130 qh_setfree( &facet->coplanarset);
2131 }
2132 }else if (!qh KEEPcoplanar || !qh KEEPinside) {
2133 qh_outerinner (NULL, NULL, &innerplane);
2134 if (qh JOGGLEmax < REALmax/2)
2135 innerplane -= qh JOGGLEmax * sqrt (qh hull_dim);
2136 numpart= 0;
2137 FORALLfacets {
2138 if (facet->coplanarset) {
2139 FOREACHpoint_(facet->coplanarset) {
2140 numpart++;
2141 qh_distplane (point, facet, &dist);
2142 if (dist < innerplane) {
2143 if (!qh KEEPinside)
2144 SETref_(point)= NULL;
2145 }else if (!qh KEEPcoplanar)
2146 SETref_(point)= NULL;
2147 }
2148 qh_setcompact (facet->coplanarset);
2149 }
2150 }
2151 zzadd_(Zcheckpart, numpart);
2152 }
2153 } /* nearcoplanar */
2154
2155 /*-<a href="qh-poly.htm#TOC"
2156 >-------------------------------</a><a name="nearvertex">-</a>
2157
2158 qh_nearvertex( facet, point, bestdist )
2159 return nearest vertex in facet to point
2160
2161 returns:
2162 vertex and its distance
2163
2164 notes:
2165 if qh.DELAUNAY
2166 distance is measured in the input set
2167 searches neighboring tricoplanar facets (requires vertexneighbors)
2168 Slow implementation. Recomputes vertex set for each point.
2169 The vertex set could be stored in the qh.keepcentrum facet.
2170 */
2171 vertexT *qh_nearvertex (facetT *facet, pointT *point, realT *bestdistp) {
2172 realT bestdist= REALmax, dist;
2173 vertexT *bestvertex= NULL, *vertex, **vertexp, *apex;
2174 coordT *center;
2175 facetT *neighbor, **neighborp;
2176 setT *vertices;
2177 int dim= qh hull_dim;
2178
2179 if (qh DELAUNAY)
2180 dim--;
2181 if (facet->tricoplanar) {
2182 if (!qh VERTEXneighbors || !facet->center) {
2183 fprintf(qh ferr, "qhull internal error (qh_nearvertex): qh.VERTEXneighbors and facet->center required for tricoplanar facets\n");
2184 qh_errexit(qh_ERRqhull, facet, NULL);
2185 }
2186 vertices= qh_settemp (qh TEMPsize);
2187 apex= SETfirst_(facet->vertices);
2188 center= facet->center;
2189 FOREACHneighbor_(apex) {
2190 if (neighbor->center == center) {
2191 FOREACHvertex_(neighbor->vertices)
2192 qh_setappend(&vertices, vertex);
2193 }
2194 }
2195 }else
2196 vertices= facet->vertices;
2197 FOREACHvertex_(vertices) {
2198 dist= qh_pointdist (vertex->point, point, -dim);
2199 if (dist < bestdist) {
2200 bestdist= dist;
2201 bestvertex= vertex;
2202 }
2203 }
2204 if (facet->tricoplanar)
2205 qh_settempfree (&vertices);
2206 *bestdistp= sqrt (bestdist);
2207 trace3((qh ferr, "qh_nearvertex: v%d dist %2.2g for f%d p%d\n", bestvertex->id, *bestdistp, facet->id, qh_pointid(point)));
2208 return bestvertex;
2209 } /* nearvertex */
2210
2211 /*-<a href="qh-poly.htm#TOC"
2212 >-------------------------------</a><a name="newhashtable">-</a>
2213
2214 qh_newhashtable( newsize )
2215 returns size of qh.hash_table of at least newsize slots
2216
2217 notes:
2218 assumes qh.hash_table is NULL
2219 qh_HASHfactor determines the number of extra slots
2220 size is not divisible by 2, 3, or 5
2221 */
2222 int qh_newhashtable(int newsize) {
2223 int size;
2224
2225 size= ((newsize+1)*qh_HASHfactor) | 0x1; /* odd number */
2226 while (True) {
2227 if ((size%3) && (size%5))
2228 break;
2229 size += 2;
2230 /* loop terminates because there is an infinite number of primes */
2231 }
2232 qh hash_table= qh_setnew (size);
2233 qh_setzero (qh hash_table, 0, size);
2234 return size;
2235 } /* newhashtable */
2236
2237 /*-<a href="qh-poly.htm#TOC"
2238 >-------------------------------</a><a name="newvertex">-</a>
2239
2240 qh_newvertex( point )
2241 returns a new vertex for point
2242 */
2243 vertexT *qh_newvertex(pointT *point) {
2244 vertexT *vertex;
2245
2246 zinc_(Ztotvertices);
2247 vertex= (vertexT *)qh_memalloc(sizeof(vertexT));
2248 memset ((char *) vertex, 0, sizeof (vertexT));
2249 if (qh vertex_id == 0xFFFFFF) {
2250 fprintf(qh ferr, "qhull input error: more than %d vertices. ID field overflows and two vertices\n\
2251 may have the same identifier. Vertices not sorted correctly.\n", 0xFFFFFF);
2252 qh_errexit(qh_ERRinput, NULL, NULL);
2253 }
2254 if (qh vertex_id == qh tracevertex_id)
2255 qh tracevertex= vertex;
2256 vertex->id= qh vertex_id++;
2257 vertex->point= point;
2258 trace4((qh ferr, "qh_newvertex: vertex p%d (v%d) created\n", qh_pointid(vertex->point), vertex->id));
2259 return (vertex);
2260 } /* newvertex */
2261
2262 /*-<a href="qh-poly.htm#TOC"
2263 >-------------------------------</a><a name="nextridge3d">-</a>
2264
2265 qh_nextridge3d( atridge, facet, vertex )
2266 return next ridge and vertex for a 3d facet
2267
2268 notes:
2269 in qh_ORIENTclock order
2270 this is a O(n^2) implementation to trace all ridges
2271 be sure to stop on any 2nd visit
2272
2273 design:
2274 for each ridge
2275 exit if it is the ridge after atridge
2276 */
2277 ridgeT *qh_nextridge3d (ridgeT *atridge, facetT *facet, vertexT **vertexp) {
2278 vertexT *atvertex, *vertex, *othervertex;
2279 ridgeT *ridge, **ridgep;
2280
2281 if ((atridge->top == facet) ^ qh_ORIENTclock)
2282 atvertex= SETsecondt_(atridge->vertices, vertexT);
2283 else
2284 atvertex= SETfirstt_(atridge->vertices, vertexT);
2285 FOREACHridge_(facet->ridges) {
2286 if (ridge == atridge)
2287 continue;
2288 if ((ridge->top == facet) ^ qh_ORIENTclock) {
2289 othervertex= SETsecondt_(ridge->vertices, vertexT);
2290 vertex= SETfirstt_(ridge->vertices, vertexT);
2291 }else {
2292 vertex= SETsecondt_(ridge->vertices, vertexT);
2293 othervertex= SETfirstt_(ridge->vertices, vertexT);
2294 }
2295 if (vertex == atvertex) {
2296 if (vertexp)
2297 *vertexp= othervertex;
2298 return ridge;
2299 }
2300 }
2301 return NULL;
2302 } /* nextridge3d */
2303 #else /* qh_NOmerge */
2304 void qh_matchduplicates (facetT *atfacet, int atskip, int hashsize, int *hashcount) {
2305 }
2306 ridgeT *qh_nextridge3d (ridgeT *atridge, facetT *facet, vertexT **vertexp) {
2307
2308 return NULL;
2309 }
2310 #endif /* qh_NOmerge */
2311
2312 /*-<a href="qh-poly.htm#TOC"
2313 >-------------------------------</a><a name="outcoplanar">-</a>
2314
2315 qh_outcoplanar()
2316 move points from all facets' outsidesets to their coplanarsets
2317
2318 notes:
2319 for post-processing under qh.NARROWhull
2320
2321 design:
2322 for each facet
2323 for each outside point for facet
2324 partition point into coplanar set
2325 */
2326 void qh_outcoplanar (void /* facet_list */) {
2327 pointT *point, **pointp;
2328 facetT *facet;
2329 realT dist;
2330
2331 trace1((qh ferr, "qh_outcoplanar: move outsideset to coplanarset for qh NARROWhull\n"));
2332 FORALLfacets {
2333 FOREACHpoint_(facet->outsideset) {
2334 qh num_outside--;
2335 if (qh KEEPcoplanar || qh KEEPnearinside) {
2336 qh_distplane (point, facet, &dist);
2337 zinc_(Zpartition);
2338 qh_partitioncoplanar (point, facet, &dist);
2339 }
2340 }
2341 qh_setfree (&facet->outsideset);
2342 }
2343 } /* outcoplanar */
2344
2345 /*-<a href="qh-poly.htm#TOC"
2346 >-------------------------------</a><a name="point">-</a>
2347
2348 qh_point( id )
2349 return point for a point id, or NULL if unknown
2350
2351 alternative code:
2352 return ((pointT *)((unsigned long)qh.first_point
2353 + (unsigned long)((id)*qh.normal_size)));
2354 */
2355 pointT *qh_point (int id) {
2356
2357 if (id < 0)
2358 return NULL;
2359 if (id < qh num_points)
2360 return qh first_point + id * qh hull_dim;
2361 id -= qh num_points;
2362 if (id < qh_setsize (qh other_points))
2363 return SETelemt_(qh other_points, id, pointT);
2364 return NULL;
2365 } /* point */
2366
2367 /*-<a href="qh-poly.htm#TOC"
2368 >-------------------------------</a><a name="point_add">-</a>
2369
2370 qh_point_add( set, point, elem )
2371 stores elem at set[point.id]
2372
2373 returns:
2374 access function for qh_pointfacet and qh_pointvertex
2375
2376 notes:
2377 checks point.id
2378 */
2379 void qh_point_add (setT *set, pointT *point, void *elem) {
2380 int id, size;
2381
2382 SETreturnsize_(set, size);
2383 if ((id= qh_pointid(point)) < 0)
2384 fprintf (qh ferr, "qhull internal warning (point_add): unknown point %p id %d\n",
2385 point, id);
2386 else if (id >= size) {
2387 fprintf (qh ferr, "qhull internal errror (point_add): point p%d is out of bounds (%d)\n",
2388 id, size);
2389 qh_errexit (qh_ERRqhull, NULL, NULL);
2390 }else
2391 SETelem_(set, id)= elem;
2392 } /* point_add */
2393
2394
2395 /*-<a href="qh-poly.htm#TOC"
2396 >-------------------------------</a><a name="pointfacet">-</a>
2397
2398 qh_pointfacet()
2399 return temporary set of facet for each point
2400 the set is indexed by point id
2401
2402 notes:
2403 vertices assigned to one of the facets
2404 coplanarset assigned to the facet
2405 outside set assigned to the facet
2406 NULL if no facet for point (inside)
2407 includes qh.GOODpointp
2408
2409 access:
2410 FOREACHfacet_i_(facets) { ... }
2411 SETelem_(facets, i)
2412
2413 design:
2414 for each facet
2415 add each vertex
2416 add each coplanar point
2417 add each outside point
2418 */
2419 setT *qh_pointfacet (void /*qh facet_list*/) {
2420 int numpoints= qh num_points + qh_setsize (qh other_points);
2421 setT *facets;
2422 facetT *facet;
2423 vertexT *vertex, **vertexp;
2424 pointT *point, **pointp;
2425
2426 facets= qh_settemp (numpoints);
2427 qh_setzero (facets, 0, numpoints);
2428 qh vertex_visit++;
2429 FORALLfacets {
2430 FOREACHvertex_(facet->vertices) {
2431 if (vertex->visitid != qh vertex_visit) {
2432 vertex->visitid= qh vertex_visit;
2433 qh_point_add (facets, vertex->point, facet);
2434 }
2435 }
2436 FOREACHpoint_(facet->coplanarset)
2437 qh_point_add (facets, point, facet);
2438 FOREACHpoint_(facet->outsideset)
2439 qh_point_add (facets, point, facet);
2440 }
2441 return facets;
2442 } /* pointfacet */
2443
2444 /*-<a href="qh-poly.htm#TOC"
2445 >-------------------------------</a><a name="pointvertex">-</a>
2446
2447 qh_pointvertex( )
2448 return temporary set of vertices indexed by point id
2449 entry is NULL if no vertex for a point
2450 this will include qh.GOODpointp
2451
2452 access:
2453 FOREACHvertex_i_(vertices) { ... }
2454 SETelem_(vertices, i)
2455 */
2456 setT *qh_pointvertex (void /*qh facet_list*/) {
2457 int numpoints= qh num_points + qh_setsize (qh other_points);
2458 setT *vertices;
2459 vertexT *vertex;
2460
2461 vertices= qh_settemp (numpoints);
2462 qh_setzero (vertices, 0, numpoints);
2463 FORALLvertices
2464 qh_point_add (vertices, vertex->point, vertex);
2465 return vertices;
2466 } /* pointvertex */
2467
2468
2469 /*-<a href="qh-poly.htm#TOC"
2470 >-------------------------------</a><a name="prependfacet">-</a>
2471
2472 qh_prependfacet( facet, facetlist )
2473 prepend facet to the start of a facetlist
2474
2475 returns:
2476 increments qh.numfacets
2477 updates facetlist, qh.facet_list, facet_next
2478
2479 notes:
2480 be careful of prepending since it can lose a pointer.
2481 e.g., can lose _next by deleting and then prepending before _next
2482 */
2483 void qh_prependfacet(facetT *facet, facetT **facetlist) {
2484 facetT *prevfacet, *list;
2485
2486
2487 trace4((qh ferr, "qh_prependfacet: prepend f%d before f%d\n", facet->id, getid_(*facetlist)));
2488 if (!*facetlist)
2489 (*facetlist)= qh facet_tail;
2490 list= *facetlist;
2491 prevfacet= list->previous;
2492 facet->previous= prevfacet;
2493 if (prevfacet)
2494 prevfacet->next= facet;
2495 list->previous= facet;
2496 facet->next= *facetlist;
2497 if (qh facet_list == list) /* this may change *facetlist */
2498 qh facet_list= facet;
2499 if (qh facet_next == list)
2500 qh facet_next= facet;
2501 *facetlist= facet;
2502 qh num_facets++;
2503 } /* prependfacet */
2504
2505
2506 /*-<a href="qh-poly.htm#TOC"
2507 >-------------------------------</a><a name="printhashtable">-</a>
2508
2509 qh_printhashtable( fp )
2510 print hash table to fp
2511
2512 notes:
2513 not in I/O to avoid bringing io.c in
2514
2515 design:
2516 for each hash entry
2517 if defined
2518 if unmatched or will merge (NULL, qh_MERGEridge, qh_DUPLICATEridge)
2519 print entry and neighbors
2520 */
2521 void qh_printhashtable(FILE *fp) {
2522 facetT *facet, *neighbor;
2523 int id, facet_i, facet_n, neighbor_i= 0, neighbor_n= 0;
2524 vertexT *vertex, **vertexp;
2525
2526 FOREACHfacet_i_(qh hash_table) {
2527 if (facet) {
2528 FOREACHneighbor_i_(facet) {
2529 if (!neighbor || neighbor == qh_MERGEridge || neighbor == qh_DUPLICATEridge)
2530 break;
2531 }
2532 if (neighbor_i == neighbor_n)
2533 continue;
2534 fprintf (fp, "hash %d f%d ", facet_i, facet->id);
2535 FOREACHvertex_(facet->vertices)
2536 fprintf (fp, "v%d ", vertex->id);
2537 fprintf (fp, "\n neighbors:");
2538 FOREACHneighbor_i_(facet) {
2539 if (neighbor == qh_MERGEridge)
2540 id= -3;
2541 else if (neighbor == qh_DUPLICATEridge)
2542 id= -2;
2543 else
2544 id= getid_(neighbor);
2545 fprintf (fp, " %d", id);
2546 }
2547 fprintf (fp, "\n");
2548 }
2549 }
2550 } /* printhashtable */
2551
2552
2553 /*-<a href="qh-poly.htm#TOC"
2554 >-------------------------------</a><a name="printlists">-</a>
2555
2556 qh_printlists( fp )
2557 print out facet and vertex list for debugging (without 'f/v' tags)
2558 */
2559 void qh_printlists (void) {
2560 facetT *facet;
2561 vertexT *vertex;
2562 int count= 0;
2563
2564 fprintf (qh ferr, "qh_printlists: facets:");
2565 FORALLfacets {
2566 if (++count % 100 == 0)
2567 fprintf (qh ferr, "\n ");
2568 fprintf (qh ferr, " %d", facet->id);
2569 }
2570 fprintf (qh ferr, "\n new facets %d visible facets %d next facet for qh_addpoint %d\n vertices (new %d):",
2571 getid_(qh newfacet_list), getid_(qh visible_list), getid_(qh facet_next),
2572 getid_(qh newvertex_list));
2573 count = 0;
2574 FORALLvertices {
2575 if (++count % 100 == 0)
2576 fprintf (qh ferr, "\n ");
2577 fprintf (qh ferr, " %d", vertex->id);
2578 }
2579 fprintf (qh ferr, "\n");
2580 } /* printlists */
2581
2582 /*-<a href="qh-poly.htm#TOC"
2583 >-------------------------------</a><a name="resetlists">-</a>
2584
2585 qh_resetlists( stats, qh_RESETvisible )
2586 reset newvertex_list, newfacet_list, visible_list
2587 if stats,
2588 maintains statistics
2589
2590 returns:
2591 visible_list is empty if qh_deletevisible was called
2592 */
2593 void qh_resetlists (boolT stats, boolT resetVisible /*qh newvertex_list newfacet_list visible_list*/) {
2594 vertexT *vertex;
2595 facetT *newfacet, *visible;
2596 int totnew=0, totver=0;
2597
2598 if (stats) {
2599 FORALLvertex_(qh newvertex_list)
2600 totver++;
2601 FORALLnew_facets
2602 totnew++;
2603 zadd_(Zvisvertextot, totver);
2604 zmax_(Zvisvertexmax, totver);
2605 zadd_(Znewfacettot, totnew);
2606 zmax_(Znewfacetmax, totnew);
2607 }
2608 FORALLvertex_(qh newvertex_list)
2609 vertex->newlist= False;
2610 qh newvertex_list= NULL;
2611 FORALLnew_facets
2612 newfacet->newfacet= False;
2613 qh newfacet_list= NULL;
2614 if (resetVisible) {
2615 FORALLvisible_facets {
2616 visible->f.replace= NULL;
2617 visible->visible= False;
2618 }
2619 qh num_visible= 0;
2620 }
2621 qh visible_list= NULL; /* may still have visible facets via qh_triangulate */
2622 qh NEWfacets= False;
2623 } /* resetlists */
2624
2625 /*-<a href="qh-poly.htm#TOC"
2626 >-------------------------------</a><a name="setvoronoi_all">-</a>
2627
2628 qh_setvoronoi_all()
2629 compute Voronoi centers for all facets
2630 includes upperDelaunay facets if qh.UPPERdelaunay ('Qu')
2631
2632 returns:
2633 facet->center is the Voronoi center
2634
2635 notes:
2636 this is unused/untested code
2637 please email bradb@shore.net if this works ok for you
2638
2639 please use:
2640 FORALLvertices {...} to locate the vertex for a point.
2641 FOREACHneighbor_(vertex) {...} to visit the Voronoi centers for a Voronoi cell.
2642 */
2643 void qh_setvoronoi_all (void) {
2644 facetT *facet;
2645
2646 qh_clearcenters (qh_ASvoronoi);
2647 qh_vertexneighbors();
2648
2649 FORALLfacets {
2650 if (!facet->normal || !facet->upperdelaunay || qh UPPERdelaunay) {
2651 if (!facet->center)
2652 facet->center= qh_facetcenter (facet->vertices);
2653 }
2654 }
2655 } /* setvoronoi_all */
2656
2657 #ifndef qh_NOmerge
2658
2659 /*-<a href="qh-poly.htm#TOC"
2660 >-------------------------------</a><a name="triangulate">-</a>
2661
2662 qh_triangulate()
2663 triangulate non-simplicial facets on qh.facet_list,
2664 if qh.CENTERtype=qh_ASvoronoi, sets Voronoi centers of non-simplicial facets
2665
2666 returns:
2667 all facets simplicial
2668 each tricoplanar facet has ->f.triowner == owner of ->center,normal,etc.
2669
2670 notes:
2671 call after qh_check_output since may switch to Voronoi centers
2672 Output may overwrite ->f.triowner with ->f.area
2673 */
2674 void qh_triangulate (void /*qh facet_list*/) {
2675 facetT *facet, *nextfacet, *owner;
2676 int onlygood= qh ONLYgood;
2677 facetT *neighbor, *visible= NULL, *facet1, *facet2, *new_facet_list= NULL;
2678 facetT *orig_neighbor= NULL, *otherfacet;
2679 vertexT *new_vertex_list= NULL;
2680 mergeT *merge;
2681 mergeType mergetype;
2682 int neighbor_i, neighbor_n;
2683
2684 trace1((qh ferr, "qh_triangulate: triangulate non-simplicial facets\n"));
2685 if (qh hull_dim == 2)
2686 return;
2687 if (qh VORONOI) { /* otherwise lose Voronoi centers [could rebuild vertex set from tricoplanar] */
2688 qh_clearcenters (qh_ASvoronoi);
2689 qh_vertexneighbors();
2690 }
2691 qh ONLYgood= False; /* for makenew_nonsimplicial */
2692 qh visit_id++;
2693 qh NEWfacets= True;
2694 qh degen_mergeset= qh_settemp (qh TEMPsize);
2695 qh newvertex_list= qh vertex_tail;
2696 for (facet= qh facet_list; facet && facet->next; facet= nextfacet) { /* non-simplicial facets moved to end */
2697 nextfacet= facet->next;
2698 if (facet->visible || facet->simplicial)
2699 continue;
2700 /* triangulate all non-simplicial facets, otherwise merging does not work, e.g., RBOX c P-0.1 P+0.1 P+0.1 D3 | QHULL d Qt Tv */
2701 if (!new_facet_list)
2702 new_facet_list= facet; /* will be moved to end */
2703 qh_triangulate_facet (facet, &new_vertex_list);
2704 }
2705 trace2((qh ferr, "qh_triangulate: delete null facets from f%d -- apex same as second vertex\n", getid_(new_facet_list)));
2706 for (facet= new_facet_list; facet && facet->next; facet= nextfacet) { /* null facets moved to end */
2707 nextfacet= facet->next;
2708 if (facet->visible)
2709 continue;
2710 if (facet->ridges) {
2711 if (qh_setsize(facet->ridges) > 0) {
2712 fprintf( qh ferr, "qhull error (qh_triangulate): ridges still defined for f%d\n", facet->id);
2713 qh_errexit (qh_ERRqhull, facet, NULL);
2714 }
2715 qh_setfree (&facet->ridges);
2716 }
2717 if (SETfirst_(facet->vertices) == SETsecond_(facet->vertices)) {
2718 zinc_(Ztrinull);
2719 qh_triangulate_null (facet);
2720 }
2721 }
2722 trace2((qh ferr, "qh_triangulate: delete %d or more mirror facets -- same vertices and neighbors\n", qh_setsize(qh degen_mergeset)));
2723 qh visible_list= qh facet_tail;
2724 while ((merge= (mergeT*)qh_setdellast (qh degen_mergeset))) {
2725 facet1= merge->facet1;
2726 facet2= merge->facet2;
2727 mergetype= merge->type;
2728 qh_memfree (merge, sizeof(mergeT));
2729 if (mergetype == MRGmirror) {
2730 zinc_(Ztrimirror);
2731 qh_triangulate_mirror (facet1, facet2);
2732 }
2733 }
2734 qh_settempfree(&qh degen_mergeset);
2735 trace2((qh ferr, "qh_triangulate: update neighbor lists for vertices from v%d\n", getid_(new_vertex_list)));
2736 qh newvertex_list= new_vertex_list; /* all vertices of new facets */
2737 qh visible_list= NULL;
2738 qh_updatevertices(/*qh newvertex_list, empty newfacet_list and visible_list*/);
2739 qh_resetlists (False, !qh_RESETvisible /*qh newvertex_list, empty newfacet_list and visible_list*/);
2740
2741 trace2((qh ferr, "qh_triangulate: identify degenerate tricoplanar facets from f%d\n", getid_(new_facet_list)));
2742 trace2((qh ferr, "qh_triangulate: and replace facet->f.triowner with tricoplanar facets that own center, normal, etc.\n"));
2743 FORALLfacet_(new_facet_list) {
2744 if (facet->tricoplanar && !facet->visible) {
2745 FOREACHneighbor_i_(facet) {
2746 if (neighbor_i == 0) { /* first iteration */
2747 if (neighbor->tricoplanar)
2748 orig_neighbor= neighbor->f.triowner;
2749 else
2750 orig_neighbor= neighbor;
2751 }else {
2752 if (neighbor->tricoplanar)
2753 otherfacet= neighbor->f.triowner;
2754 else
2755 otherfacet= neighbor;
2756 if (orig_neighbor == otherfacet) {
2757 zinc_(Ztridegen);
2758 facet->degenerate= True;
2759 break;
2760 }
2761 }
2762 }
2763 }
2764 }
2765
2766 trace2((qh ferr, "qh_triangulate: delete visible facets -- non-simplicial, null, and mirrored facets\n"));
2767 owner= NULL;
2768 visible= NULL;
2769 for (facet= new_facet_list; facet && facet->next; facet= nextfacet) { /* may delete facet */
2770 nextfacet= facet->next;
2771 if (facet->visible) {
2772 if (facet->tricoplanar) { /* a null or mirrored facet */
2773 qh_delfacet(facet);
2774 qh num_visible--;
2775 }else { /* a non-simplicial facet followed by its tricoplanars */
2776 if (visible && !owner) {
2777 /* RBOX 200 s D5 t1001471447 | QHULL Qt C-0.01 Qx Qc Tv Qt -- f4483 had 6 vertices/neighbors and 8 ridges */
2778 trace2((qh ferr, "qh_triangulate: all tricoplanar facets degenerate for non-simplicial facet f%d\n",
2779 visible->id));
2780 qh_delfacet(visible);
2781 qh num_visible--;
2782 }
2783 visible= facet;
2784 owner= NULL;
2785 }
2786 }else if (facet->tricoplanar) {
2787 if (facet->f.triowner != visible) {
2788 fprintf( qh ferr, "qhull error (qh_triangulate): tricoplanar facet f%d not owned by its visible, non-simplicial facet f%d\n", facet->id, getid_(visible));
2789 qh_errexit2 (qh_ERRqhull, facet, visible);
2790 }
2791 if (owner)
2792 facet->f.triowner= owner;
2793 else if (!facet->degenerate) {
2794 owner= facet;
2795 nextfacet= visible->next; /* rescan tricoplanar facets with owner */
2796 facet->keepcentrum= True; /* one facet owns ->normal, etc. */
2797 facet->coplanarset= visible->coplanarset;
2798 facet->outsideset= visible->outsideset;
2799 visible->coplanarset= NULL;
2800 visible->outsideset= NULL;
2801 if (!qh TRInormals) { /* center and normal copied to tricoplanar facets */
2802 visible->center= NULL;
2803 visible->normal= NULL;
2804 }
2805 qh_delfacet(visible);
2806 qh num_visible--;
2807 }
2808 }
2809 }
2810 if (visible && !owner) {
2811 trace2((qh ferr, "qh_triangulate: all tricoplanar facets degenerate for last non-simplicial facet f%d\n", visible->id));
2812 qh_delfacet(visible);
2813 qh num_visible--;
2814 }
2815 qh NEWfacets= False;
2816 qh ONLYgood= onlygood; /* restore value */
2817 if (qh CHECKfrequently)
2818 qh_checkpolygon (qh facet_list);
2819 } /* triangulate */
2820
2821
2822 /*-<a href="qh-poly.htm#TOC"
2823 >-------------------------------</a><a name="triangulate_facet">-</a>
2824
2825 qh_triangulate_facet (facetA)
2826 triangulate a non-simplicial facet
2827 if qh.CENTERtype=qh_ASvoronoi, sets its Voronoi center
2828 returns:
2829 qh.newfacet_list == simplicial facets
2830 facet->tricoplanar set and ->keepcentrum false
2831 facet->degenerate set if duplicated apex
2832 facet->f.trivisible set to facetA
2833 facet->center copied from facetA (created if qh_ASvoronoi)
2834 qh_eachvoronoi, qh_detvridge, qh_detvridge3 assume centers copied
2835 facet->normal,offset,maxoutside copied from facetA
2836
2837 notes:
2838 qh_makenew_nonsimplicial uses neighbor->seen for the same
2839
2840 see also:
2841 qh_addpoint() -- add a point
2842 qh_makenewfacets() -- construct a cone of facets for a new vertex
2843
2844 design:
2845 if qh_ASvoronoi,
2846 compute Voronoi center (facet->center)
2847 select first vertex (highest ID to preserve ID ordering of ->vertices)
2848 triangulate from vertex to ridges
2849 copy facet->center, normal, offset
2850 update vertex neighbors
2851 */
2852 void qh_triangulate_facet (facetT *facetA, vertexT **first_vertex) {
2853 facetT *newfacet;
2854 facetT *neighbor, **neighborp;
2855 vertexT *apex;
2856 int numnew=0;
2857
2858 trace3((qh ferr, "qh_triangulate_facet: triangulate facet f%d\n", facetA->id));
2859
2860 if (qh IStracing >= 4)
2861 qh_printfacet (qh ferr, facetA);
2862 FOREACHneighbor_(facetA) {
2863 neighbor->seen= False;
2864 neighbor->coplanar= False;
2865 }
2866 if (qh CENTERtype == qh_ASvoronoi && !facetA->center /* matches upperdelaunay in qh_setfacetplane() */
2867 && fabs_(facetA->normal[qh hull_dim -1]) >= qh ANGLEround * qh_ZEROdelaunay) {
2868 facetA->center= qh_facetcenter (facetA->vertices);
2869 }
2870 qh_willdelete (facetA, NULL);
2871 qh newfacet_list= qh facet_tail;
2872 facetA->visitid= qh visit_id;
2873 apex= SETfirst_(facetA->vertices);
2874 qh_makenew_nonsimplicial (facetA, apex, &numnew);
2875 SETfirst_(facetA->neighbors)= NULL;
2876 FORALLnew_facets {
2877 newfacet->tricoplanar= True;
2878 newfacet->f.trivisible= facetA;
2879 newfacet->degenerate= False;
2880 newfacet->upperdelaunay= facetA->upperdelaunay;
2881 newfacet->good= facetA->good;
2882 if (qh TRInormals) {
2883 newfacet->keepcentrum= True;
2884 newfacet->normal= qh_copypoints (facetA->normal, 1, qh hull_dim);
2885 if (qh CENTERtype == qh_AScentrum)
2886 newfacet->center= qh_getcentrum (newfacet);
2887 else
2888 newfacet->center= qh_copypoints (facetA->center, 1, qh hull_dim);
2889 }else {
2890 newfacet->keepcentrum= False;
2891 newfacet->normal= facetA->normal;
2892 newfacet->center= facetA->center;
2893 }
2894 newfacet->offset= facetA->offset;
2895 #if qh_MAXoutside
2896 newfacet->maxoutside= facetA->maxoutside;
2897 #endif
2898 }
2899 qh_matchnewfacets(/*qh newfacet_list*/);
2900 zinc_(Ztricoplanar);
2901 zadd_(Ztricoplanartot, numnew);
2902 zmax_(Ztricoplanarmax, numnew);
2903 qh visible_list= NULL;
2904 if (!(*first_vertex))
2905 (*first_vertex)= qh newvertex_list;
2906 qh newvertex_list= NULL;
2907 qh_updatevertices(/*qh newfacet_list, empty visible_list and newvertex_list*/);
2908 qh_resetlists (False, !qh_RESETvisible /*qh newfacet_list, empty visible_list and newvertex_list*/);
2909 } /* triangulate_facet */
2910
2911 /*-<a href="qh-poly.htm#TOC"
2912 >-------------------------------</a><a name="triangulate_link">-</a>
2913
2914 qh_triangulate_link (oldfacetA, facetA, oldfacetB, facetB)
2915 relink facetA to facetB via oldfacets
2916 returns:
2917 adds mirror facets to qh degen_mergeset (4-d and up only)
2918 design:
2919 if they are already neighbors, the opposing neighbors become MRGmirror facets
2920 */
2921 void qh_triangulate_link (facetT *oldfacetA, facetT *facetA, facetT *oldfacetB, facetT *facetB) {
2922 int errmirror= False;
2923
2924 trace3((qh ferr, "qh_triangulate_link: relink old facets f%d and f%d between neighbors f%d and f%d\n", oldfacetA->id, oldfacetB->id, facetA->id, facetB->id));
2925 if (qh_setin (facetA->neighbors, facetB)) {
2926 if (!qh_setin (facetB->neighbors, facetA))
2927 errmirror= True;
2928 else
2929 qh_appendmergeset (facetA, facetB, MRGmirror, NULL);
2930 }else if (qh_setin (facetB->neighbors, facetA))
2931 errmirror= True;
2932 if (errmirror) {
2933 fprintf( qh ferr, "qhull error (qh_triangulate_link): mirror facets f%d and f%d do not match for old facets f%d and f%d\n",
2934 facetA->id, facetB->id, oldfacetA->id, oldfacetB->id);
2935 qh_errexit2 (qh_ERRqhull, facetA, facetB);
2936 }
2937 qh_setreplace (facetB->neighbors, oldfacetB, facetA);
2938 qh_setreplace (facetA->neighbors, oldfacetA, facetB);
2939 } /* triangulate_link */
2940
2941 /*-<a href="qh-poly.htm#TOC"
2942 >-------------------------------</a><a name="triangulate_mirror">-</a>
2943
2944 qh_triangulate_mirror (facetA, facetB)
2945 delete mirrored facets from qh_triangulate_null() and qh_triangulate_mirror
2946 a mirrored facet shares the same vertices of a logical ridge
2947 design:
2948 since a null facet duplicates the first two vertices, the opposing neighbors absorb the null facet
2949 if they are already neighbors, the opposing neighbors become MRGmirror facets
2950 */
2951 void qh_triangulate_mirror (facetT *facetA, facetT *facetB) {
2952 facetT *neighbor, *neighborB;
2953 int neighbor_i, neighbor_n;
2954
2955 trace3((qh ferr, "qh_triangulate_mirror: delete mirrored facets f%d and f%d\n", facetA->id, facetB->id));
2956 FOREACHneighbor_i_(facetA) {
2957 neighborB= SETelemt_(facetB->neighbors, neighbor_i, facetT);
2958 if (neighbor == neighborB)
2959 continue; /* occurs twice */
2960 qh_triangulate_link (facetA, neighbor, facetB, neighborB);
2961 }
2962 qh_willdelete (facetA, NULL);
2963 qh_willdelete (facetB, NULL);
2964 } /* triangulate_mirror */
2965
2966 /*-<a href="qh-poly.htm#TOC"
2967 >-------------------------------</a><a name="triangulate_null">-</a>
2968
2969 qh_triangulate_null (facetA)
2970 remove null facetA from qh_triangulate_facet()
2971 a null facet has vertex #1 (apex) == vertex #2
2972 returns:
2973 adds facetA to ->visible for deletion after qh_updatevertices
2974 qh degen_mergeset contains mirror facets (4-d and up only)
2975 design:
2976 since a null facet duplicates the first two vertices, the opposing neighbors absorb the null facet
2977 if they are already neighbors, the opposing neighbors become MRGmirror facets
2978 */
2979 void qh_triangulate_null (facetT *facetA) {
2980 facetT *neighbor, *otherfacet;
2981
2982 trace3((qh ferr, "qh_triangulate_null: delete null facet f%d\n", facetA->id));
2983 neighbor= SETfirst_(facetA->neighbors);
2984 otherfacet= SETsecond_(facetA->neighbors);
2985 qh_triangulate_link (facetA, neighbor, facetA, otherfacet);
2986 qh_willdelete (facetA, NULL);
2987 } /* triangulate_null */
2988
2989 #else /* qh_NOmerge */
2990 void qh_triangulate (void) {
2991 }
2992 #endif /* qh_NOmerge */
2993
2994 /*-<a href="qh-poly.htm#TOC"
2995 >-------------------------------</a><a name="vertexintersect">-</a>
2996
2997 qh_vertexintersect( vertexsetA, vertexsetB )
2998 intersects two vertex sets (inverse id ordered)
2999 vertexsetA is a temporary set at the top of qhmem.tempstack
3000
3001 returns:
3002 replaces vertexsetA with the intersection
3003
3004 notes:
3005 could overwrite vertexsetA if currently too slow
3006 */
3007 void qh_vertexintersect(setT **vertexsetA,setT *vertexsetB) {
3008 setT *intersection;
3009
3010 intersection= qh_vertexintersect_new (*vertexsetA, vertexsetB);
3011 qh_settempfree (vertexsetA);
3012 *vertexsetA= intersection;
3013 qh_settemppush (intersection);
3014 } /* vertexintersect */
3015
3016 /*-<a href="qh-poly.htm#TOC"
3017 >-------------------------------</a><a name="vertexintersect_new">-</a>
3018
3019 qh_vertexintersect_new( )
3020 intersects two vertex sets (inverse id ordered)
3021
3022 returns:
3023 a new set
3024 */
3025 setT *qh_vertexintersect_new (setT *vertexsetA,setT *vertexsetB) {
3026 setT *intersection= qh_setnew (qh hull_dim - 1);
3027 vertexT **vertexA= SETaddr_(vertexsetA, vertexT);
3028 vertexT **vertexB= SETaddr_(vertexsetB, vertexT);
3029
3030 while (*vertexA && *vertexB) {
3031 if (*vertexA == *vertexB) {
3032 qh_setappend(&intersection, *vertexA);
3033 vertexA++; vertexB++;
3034 }else {
3035 if ((*vertexA)->id > (*vertexB)->id)
3036 vertexA++;
3037 else
3038 vertexB++;
3039 }
3040 }
3041 return intersection;
3042 } /* vertexintersect_new */
3043
3044 /*-<a href="qh-poly.htm#TOC"
3045 >-------------------------------</a><a name="vertexneighbors">-</a>
3046
3047 qh_vertexneighbors()
3048 for each vertex in qh.facet_list,
3049 determine its neighboring facets
3050
3051 returns:
3052 sets qh.VERTEXneighbors
3053 nop if qh.VERTEXneighbors already set
3054 qh_addpoint() will maintain them
3055
3056 notes:
3057 assumes all vertex->neighbors are NULL
3058
3059 design:
3060 for each facet
3061 for each vertex
3062 append facet to vertex->neighbors
3063 */
3064 void qh_vertexneighbors (void /*qh facet_list*/) {
3065 facetT *facet;
3066 vertexT *vertex, **vertexp;
3067
3068 if (qh VERTEXneighbors)
3069 return;
3070 trace1((qh ferr, "qh_vertexneighbors: determing neighboring facets for each vertex\n"));
3071 qh vertex_visit++;
3072 FORALLfacets {
3073 if (facet->visible)
3074 continue;
3075 FOREACHvertex_(facet->vertices) {
3076 if (vertex->visitid != qh vertex_visit) {
3077 vertex->visitid= qh vertex_visit;
3078 vertex->neighbors= qh_setnew (qh hull_dim);
3079 }
3080 qh_setappend (&vertex->neighbors, facet);
3081 }
3082 }
3083 qh VERTEXneighbors= True;
3084 } /* vertexneighbors */
3085
3086 /*-<a href="qh-poly.htm#TOC"
3087 >-------------------------------</a><a name="vertexsubset">-</a>
3088
3089 qh_vertexsubset( vertexsetA, vertexsetB )
3090 returns True if vertexsetA is a subset of vertexsetB
3091 assumes vertexsets are sorted
3092
3093 note:
3094 empty set is a subset of any other set
3095 */
3096 boolT qh_vertexsubset(setT *vertexsetA, setT *vertexsetB) {
3097 vertexT **vertexA= (vertexT **) SETaddr_(vertexsetA, vertexT);
3098 vertexT **vertexB= (vertexT **) SETaddr_(vertexsetB, vertexT);
3099
3100 while (True) {
3101 if (!*vertexA)
3102 return True;
3103 if (!*vertexB)
3104 return False;
3105 if ((*vertexA)->id > (*vertexB)->id)
3106 return False;
3107 if (*vertexA == *vertexB)
3108 vertexA++;
3109 vertexB++;
3110 }
3111 return False; /* avoid warnings */
3112 } /* vertexsubset */