A program for convex hulls

Summary

Hull is an ANSI C program that computes the convex hull of a point set in general (but small!) dimension. The input is a list of points, and the output is a list of facets of the convex hull of the points, each facet presented as a list of its vertices. (The facets are assumed to be simplices, such as triangles in 3d; this is enforced by tiebreaking, giving a triangulation of a facet by "placing".) The program can also compute Delaunay triangulations and alpha shapes, and volumes of Voronoi regions. The program uses exact arithmetic when possible, with a moderate speed penalty. (Typically a factor of 2 or 3 for Delaunay triangulation, less for convex hulls). Output in postscript and OFF format for geomview is supported.

This document is an html version of the manpage, with links to papers describing the algorithm.

At long last, I've made the code a bit more portable: it compiles under gcc with the provided makefile, and I've included Visual Studio 7 solution and project files to build it on Windows (as a console application). I'm told it also compiles and runs fine under Solaris. clarkson@research.bell-labs.com .


Source Code


Examples of use


Synopsis

hull -d -f<format> -A -aa<alpha> -af<format> -oN -ov -s<seed> -r -m<multiplier> -X<debug_file> -i<input_file> -oF<output_file>


Description

The input_file (default stdin) is a sequence of points (also called sites), separated by \n; a d-dimensional point is specified as a group of d floats separated by whitespace (other than \n).

The output_file (default stdout) gives d-tuples of the input points, where a point is given as an integer i if it was the i'th point in the input_file. If the convex hull lies in a flat (affine subspace) of dimension d', the output will comprise a list of d'-tuples, the vertices of the convex hull relative to that flat.

The output tuples represent the facets of the convex hull of the input set. A facet which is not a simplex is output implicitly as the collection of simplices of its triangulation.

The output for Delaunay triangulation includes a "point at infinity", numbered -1; facets including it correspond to facets of the convex hull of the sites.

Some chatter will appear on debug_file (default stderr).

The coordinates of the input points are rounded to integers. With -m<multiplier>, the coordinates are multiplied by multiplier before this rounding.

The program attempts to use a method that gives exact answers to numerical tests; in some circumstances, this may fail, and with some warnings, the program continues, using approximate arithmetic.

More detail on the options:

-d
compute the Delaunay triangulation of the input set.
-f<format>
give the main output (convex hull or Delaunay triangulation) in output <format>, which is by default the list of vertex numbers described above, or ps, for postscript output of planar points, or off, for OFF output of 3d points.
-aa<alpha>
compute the alpha shape using parameter < alpha >: the output is the set of all d-tuples of sites such that there is a ball of radius < alpha > containing those sites on its bounding sphere, and containing no other sites. A Delaunay triangulation computation is implied by this option and by -A.
-af<format>
output the alpha shape in format <format>, as in option -f.
-A
compute the alpha shape of the input, finding the smallest alpha so that the sites are all contained in the alpha-shape.
-r
randomly shuffle the input points; this may speed up the program.
-s<seed>
randomly shuffle using < seed > for the random number generator.
-oN
do not produce main output (convex hull or Delaunay triangulation). If you want an alpha shape only, you need this to turn off the Delaunay output.
-ov
Give volumes of Voronoi regions of input sites, and in general d'-volumes of d'-dimensional Voronoi cells. Implies -d.


Bugs/Portability

If the convex hull is a single point, the algorithm will fail to report it. All other degeneracies should be handled. Tie-breaking is done so that all reported facets are simplices. If the input points are degenerate, some hull facets may be; for example, some Delaunay simplices may have zero volume. Determining non-simplicial facets or deleting zero-volume Delaunauy simplices could be done in post-processing (not implemented).

The file rand.c includes calls to pseudo-random number generators;

No simplices are deleted; the only way to free storage is to free it all using free_hull_storage.


Author

Ken Clarkson, (clarkson@research.bell-labs.com),
using an earlier version written by Susan Dorward, who is not to blame.


Algorithms

Combinatorial

The program implements an incremental algorithm, which can be randomized if desired; the algorithm is roughly as described in Section 3 (page 7) of

K. L.Clarkson, K.Mehlhorn, and R.Seidel. Four results on randomized incremental constructions. Comp. Geom.: Theory and Applications, pages 185--121, 1993. Preliminary version in Proc. Symp. Theor. Aspects of Comp. Sci., 1992.

The code is an implementation of the "all-visibilities" search procedure. The "anti-origin" point of that discussion is represented as a NULL point in the code, so that simplex s represents a facet iff s->peak.vert is NULL.

Numerical

The numerical code, for computing normals to facets and so on, is an implementation of

K. L.Clarkson. Safe and effective determinant evaluation. In Proc. 31st IEEE Symposium on Foundations of Computer Science, pages 387--395, Pittsburgh, PA, October 1992.

The conditioning procedure of that paper is used, and applied to finding a normal vector to a facet of the convex hull. The visibility test is then a dot product. Vectors a_1, a_2... are obtained from points p_0, p_1,...by translating p_0 to the origin: a_i = p_i-p_0, for i=1..cdim.


Function call

The main function to call is

simplex *build_convex_hull(gsitef *get_s, site_n *site_numm, short dim, short vdd)
returns pointer to root simplex of triangulation of convex hull (or Delaunay triangulation if vdd==1).

see hull.h for gsitef typedef etc.

get_s
returns next site each call; hull construction stops when NULL returned;
site_numm
returns number of site when given site;
dim
dimension of point set;
vdd
if (vdd) then return Delaunay triangulation, otherwise convex hull;