Free cookie consent management tool by TermsFeed Policy Generator

source: branches/HeuristicLab.Analysis.AlgorithmBehavior/qhull-2012.1/src/libqhull/poly2.c @ 10635

Last change on this file since 10635 was 10207, checked in by ascheibe, 11 years ago

#1886 added a unit test for volume calculation and the qhull library

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