Up: Home page for Qhull
Up: Qhull manual: Table of Contents
To: ProgramsOptionsOutputFormatsGeomviewPrintQhullPrecisionTrace
To: Qhull internals: Table of Contents (please wait while loading)
Dn: Qhull functions, macros, and data structures


[4-d cube] Qhull internals

This section discusses the internals of Qhull.

Copyright © 1995-2003 The Geometry Center, Minneapolis MN


»Qhull internals: Table of Contents


»Performance of Qhull

Empirically, Qhull's performance is balanced in the sense that the average case happens on average. This may always be true if the precision of the input is limited to at most O(log n) bits. Empirically, the maximum number of vertices occurs at the end of constructing the hull.

Let n be the number of input points, v be the number of output vertices, and f_v be the maximum number of facets for a convex hull of v vertices. If both conditions hold, Qhull runs in O(n log v) in 2-d and 3-d and O(n f_v/v) otherwise. The function f_v increases rapidly with dimension. It is O(v^floor(d/2) / floor(d/2)!).

The time complexity for merging is unknown. Options 'C-0' and 'Qx' (defaults) handle precision problems due to floating-point arithmetic. They are optimized for simplicial outputs.

When running large data sets, you should monitor Qhull's performance with the 'TFn' option. The time per facet is approximately constant. In high-d with many merged facets, the size of the ridge sets grows rapidly. For example the product of 8-d simplices contains 18 facets and 500,000 ridges. This will increase the time needed per facet.

As dimension increases, the number of facets and ridges in a convex hull grows rapidly for the same number of vertices. For example, the convex hull of 300 cospherical points in 6-d has 30,000 facets.

If Qhull appears to stop processing facets, check the memory usage of Qhull. If more than 5-10% of Qhull is in virtual memory, its performance will degrade rapidly.

When building hulls in 20-d and higher, you can follow the progress of Qhull with option 'T1'. It reports each major event in processing a point.

To reduce memory requirements, recompile Qhull for single-precision reals (REALfloat in user.h). Single-precision does not work with joggle ('QJ'). Check qh_MEMalign in user.h and the match between free list sizes and data structure sizes (see the end of the statistics report from 'Ts'). If free list sizes do not match, you may be able to use a smaller qh_MEMalign. Setting qh_COMPUTEfurthest saves a small amount of memory, as does clearing qh_MAXoutside (both in user.h).

Shewchuk is working on a 3-d version of his triangle program. It is optimized for 3-d simplicial Delaunay triangulation and uses less memory than Qhull.

To reduce the size of the Qhull executable, consider qh_NOtrace and qh_KEEPstatistics 0 in user.h. By changing user.c you can also remove the input/output code in io.c. If you don't need facet merging, then version 1.01 of Qhull is much smaller. It contains some bugs that prevent Qhull from initializing in simple test cases. It is slower in high dimensions.

The precision options, 'Vn', 'Wn', 'Un'. 'A-n', 'C-n', 'An', 'Cn', and 'Qx', may have large effects on Qhull performance. You will need to experiment to find the best combination for your application.

The verify option ('Tv') checks every point after the hull is complete. If facet merging is used, it checks that every point is inside every facet. This can take a very long time if there are many points and many facets. You can interrupt the verify without losing your output. If facet merging is not used and there are many points and facets, Qhull uses a directed search instead of an exhaustive search. This should be fast enough for most point sets. Directed search is not used for facet merging because directed search was already used for updating the facets' outer planes.

The check-frequently option ('Tc') becomes expensive as the dimension increases. The verify option ('Tv') performs many of the same checks before outputting the results.

Options 'Q0' (no pre-merging), 'Q3' (no checks for redundant vertices), 'Q5' (no updates for outer planes), and 'Q8' (no near-interior points) increase Qhull's speed. The corresponding operations may not be needed in your application.

In 2-d and 3-d, a partial hull may be faster to produce. Option 'QgGn' only builds facets visible to point n. Option 'QgVn' only builds facets that contain point n. In higher-dimensions, this does not reduce the number of facets.

User.h includes a number of performance-related constants. Changes may improve Qhull performance on your data sets. To understand their effect on performance, you will need to read the corresponding code.

GNU gprof reports that the dominate cost for 3-d convex hull of cosperical points is qh_distplane(), mainly called from qh_findbestnew(). The dominate cost for 3-d Delaunay triangulation is creating new facets in qh_addpoint(), while qh_distplane() remains the most expensive function.

»Calling Qhull from your program

Warning: Qhull was not designed for calling from other programs. There is neither API nor Qhull classes. It can be done, but it takes work and head scratching. You will need to understand the data structures and read the code. Most users will find it easier to call Qhull as an external command.

For examples of calling Qhull, see GNU Octave's computational geometry code, and Qhull's user_eg.c, user_eg2.c, user.c, and qhull_interface.cpp. Qhull's programs use the same library: unix.c, qconvex.c, qdelaun.c, qhalf.c, and qvoronoi.c.

The BGL Boost Graph Library [aka GGCL] provides C++ classes for graph data structures and algorithms [Dr. Dobb's 9/00 p. 29-38; OOPSLA '99 p. 399-414]. It is modelled after the Standard Template Library. It would provide a good interface to Qhull. If you are interested in adapting BGL to Qhull, please contact bradb@qhull.org.

See Qhull functions, macros, and data structures for internal documentation of Qhull. The documentation provides an overview and index. To use the library you will need to read and understand the code. For most users, it is better to write data to a file, call the qhull program, and read the results from the output file.

When you read the code, be aware of the macros "qh" and "qhstat", e.g., "qh hull_dim". They are defined in qhull.h. They allow the global data structures to be pre-allocated (faster access) or dynamically allocated (allows multiple copies).

Qhull's Makefile produces a library, libqhull.a, for inclusion in your programs. First review qhull.h. This defines the data structures used by Qhull and provides prototypes for the top-level functions. Most users will only need qhull.h in their programs. For example, the Qhull program is defined with qhull.h and unix.c. To access all functions, use qhull_a.h. Include the file with "#include <qhull/qhull_a.h>". This avoids potential name conflicts.

If you use the Qhull library, you are on your own as far as bugs go. Start with small examples for which you know the output. If you get a bug, try to duplicate it with the Qhull program. The 'Tc' option will catch many problems as they occur. When an error occurs, use 'T4 TPn' to trace from the last point added to the hull. Compare your trace with the trace output from the Qhull program.

Errors in the Qhull library are more likely than errors in the Qhull program. These are usually due to feature interactions that do not occur in the Qhull program. Please report all errors that you find in the Qhull library. Please include suggestions for improvement.

»sets and quick memory allocation

You can use mem.c and qset.c individually. Mem.c implements quick-fit memory allocation. It is faster than malloc/free in applications that allocate and deallocate lots of memory.

Qset.c implements sets and related collections. It's the inner loop of Qhull, so speed is more important than abstraction. Set iteration is particularly fast. qset.c just includes the functions needed for Qhull.

»Delaunay triangulations and point indices

Here some unchecked code to print the point indices of each Delaunay triangle. Use option 'QJ' if you want to avoid non-simplicial facets. Note that upper Delaunay regions are skipped. These facets correspond to the furthest-site Delaunay triangulation.

  facetT *facet;
  vertexT *vertex, **vertexp;
  
  FORALLfacets {
    if (!facet->upperdelaunay) {
      printf ("%d", qh_setsize (facet->vertices);
      FOREACHvertex_(facet->vertices)
        printf (" %d", qh_pointid (vertex->point));
      printf ("\n");
    }
  }

»locate a facet with qh_findbestfacet()

The routine qh_findbestfacet in poly2.c is particularly useful. It uses a directed search to locate the facet that is furthest below a point. For Delaunay triangulations, this facet is the Delaunay triangle that contains the lifted point. For convex hulls, the distance of a point to the convex hull is either the distance to this facet or the distance to a subface of the facet.

Warning: If triangulated output ('Qt') and the best facet is triangulated, qh_findbestfacet() returns one of the corresponding 'tricoplanar' facets. The actual best facet may be a different tricoplanar facet.

See qh_nearvertex() in poly2.c for sample code to visit each tricoplanar facet. To identify the correct tricoplanar facet, see Devillers, et. al., ['01] and Mucke, et al ['96]. If you implement this test in general dimension, please notify qhull@qhull.org.

qh_findbestfacet performs an exhaustive search if its directed search returns a facet that is above the point. This occurs when the point is inside the hull or if the curvature of the convex hull is less than the curvature of a sphere centered at the point (e.g., a point near a lens-shaped convex hull). When the later occurs, the distance function is bimodal and a directed search may return a facet on the far side of the convex hull.

Algorithms that retain the previously constructed hulls usually avoid an exhaustive search for the best facet. You may use a hierarchical decomposition of the convex hull [Dobkin and Kirkpatrick '90].

To use qh_findbestfacet with Delaunay triangulations, lift the point to a paraboloid by summing the squares of its coordinates (see qh_setdelaunay in geom2.c). Do not scale the input with options 'Qbk', 'QBk', 'QbB' or 'Qbb'. See Mucke, et al ['96] for a good point location algorithm.

The intersection of a ray with the convex hull may be found by locating the facet closest to a distant point on the ray. Intersecting the ray with the facet's hyperplane gives a new point to test.

»on-line construction with qh_addpoint()

The Qhull library may be used for the on-line construction of convex hulls, Delaunay triangulations, and halfspace intersections about a point. It may be slower than implementations that retain intermediate convex hulls (e.g., Clarkson's hull program). These implementations always use a directed search. For the on-line construction of convex hulls and halfspace intersections, Qhull may use an exhaustive search (qh_findbestfacet).

You may use qh_findbestfacet and qh_addpoint (qhull.c) to add a point to a convex hull. Do not modify the point's coordinates since qh_addpoint does not make a copy of the coordinates. For Delaunay triangulations, you need to lift the point to a paraboloid by summing the squares of the coordinates (see qh_setdelaunay in geom2.c). Do not scale the input with options 'Qbk', 'QBk', 'QbB' or 'Qbb'. Do not deallocate the point's coordinates. You need to provide a facet that is below the point (qh_findbestfacet).

You can not delete points. Another limitation is that Qhull uses the initial set of points to determine the maximum roundoff error (via the upper and lower bounds for each coordinate).

For many applications, it is better to rebuild the hull from scratch for each new point. This is especially true if the point set is small or if many points are added at a time.

Calling qh_addpoint from your program may be slower than recomputing the convex hull with qh_qhull. This is especially true if the added points are not appended to the qh first_point array. In this case, Qhull must search a set to determine a point's ID. [R. Weber]

See user_eg.c for examples of the on-line construction of convex hulls, Delaunay triangulations, and halfspace intersections. The outline is:

initialize qhull with an initial set of points
qh_qhull(); 

for each additional point p
   append p to the end of the point array or allocate p separately
   lift p to the paraboloid by calling qh_setdelaunay
   facet= qh_findbestfacet (p, !qh_ALL, &bestdist, &isoutside);
   if (isoutside) 
      if (!qh_addpoint (point, facet, False))
         break;  /* user requested an early exit with 'TVn' or 'TCn' */
   
call qh_check_maxout() to compute outer planes
terminate qhull

»Constrained Delaunay triangulation

With a fair amount of work, Qhull is suitable for constrained Delaunay triangulation. See Shewchuk, ACM Symposium on Computational Geometry, Minneapolis 1998.

Here's a quick way to add a constraint to a Delaunay triangulation: subdivide the constraint into pieces shorter than the minimum feature separation. You will need an independent check of the constraint in the output since the minimum feature separation may be incorrect. [H. Geron]

»Tricoplanar facets and option 'Qt'

Option 'Qt' triangulates non-simplicial facets (e.g., a square facet in 3-d or a cubical facet in 4-d). All facets share the same apex (i.e., the first vertex in facet->vertices). For each triangulated facet, Qhull sets facet->tricoplanar true and copies facet->center, facet->normal, facet->offset, and facet->maxoutside. One of the facets owns facet->normal; its facet->keepcentrum is true. If facet->isarea is false, facet->triowner points to the owning facet.

Qhull sets facet->degenerate if the facet's vertices belong to the same ridge of the non-simplicial facet.

To visit each tricoplanar facet of a non-simplicial facet, either visit all neighbors of the apex or recursively visit all neighbors of a tricoplanar facet. The tricoplanar facets will have the same facet->center.

See qh_detvridge for an example of ignoring tricoplanar facets.

»Voronoi vertices of a region

The following code iterates over all Voronoi vertices for each Voronoi region. Qhull computes Voronoi vertices from the convex hull that corresponds to a Delaunay triangulation. An input site corresponds to a vertex of the convex hull and a Voronoi vertex corresponds to an adjacent facet. A facet is "upperdelaunay" if it corresponds to a Voronoi vertex "at-infinity". Qhull uses qh_printvoronoi in io.c for 'qvoronoi o'

/* please review this code for correctness */
qh_setvoronoi_all(); 
FORALLvertices {
   site_id = qh_pointid (vertex->point);
   if (qh hull_dim == 3)
      qh_order_vertexneighbors(vertex);
   infinity_seen = 0;
   FOREACHneighbor_(vertex) {
      if (neighbor->upperdelaunay) {
        if (!infinity_seen) {
          infinity_seen = 1;
          ... process a Voronoi vertex "at infinity" ...
        }
      }else {
        voronoi_vertex = neighbor->center;
        ... your code goes here ...
      }
   }
}

»Voronoi vertices of a ridge

Qhull uses qh_printvdiagram() in io.c to print the ridges of a Voronoi diagram for option 'Fv'. The helper function qh_eachvoronoi() does the real work. It calls the callback 'printvridge' for each ridge of the Voronoi diagram.

You may call qh_printvdiagram2(), qh_eachvoronoi(), or qh_eachvoronoi_all() with your own function. If you do not need the total number of ridges, you can skip the first call to qh_printvdiagram2(). See qh_printvridge() and qh_printvnorm() in io.c for examples.

»vertex neighbors of a vertex

To visit all of the vertices that share an edge with a vertex:

For non-simplicial facets, the ridges form a simplicial decomposition of the (d-2)-faces between each pair of facets -- if you need 1-faces, you probably need to generate the full face graph of the convex hull.

»Enhancements to Qhull

There are many ways in which Qhull can be improved.

Here is a partial list:
 - fix finddelaunay() in user_eg.c for tricoplanar facets
 - write a BGL, C++ interface to Qhull 
     http://www.boost.org/libs/graph/doc/table_of_contents.html
 - change qh_save_qhull to swap the qhT structure instead of using pointers
 - change error handling and tracing to be independent of 'qh ferr'
 - determine the maximum width for a given set of parameters
 - prove that directed search locates all coplanar facets
 - in high-d merging, can a loop of facets become disconnected?
 - find a way to improve inner hulls in 5-d and higher
 - determine the best policy for facet visibility ('Vn')
 - determine the limitations of 'Qg'

Precision improvements:
 - For 'Qt', resolve cross-linked, butterfly ridges.  
     May allow retriangulation in qh_addpoint().
 - for Delaunay triangulations ('d' or 'v') under joggled input ('QJ'),
     remove vertical facets whose lowest vertex may be coplanar with convex hull 
 - review use of 'Qbb' with 'd QJ'.  Is MAXabs_coord better than MAXwidth?
 - check Sugihara and Iri's better in-sphere test [Canadian 
     Conf. on Comp. Geo., 1989; Univ. of Tokyo RMI 89-05]
 - replace centrum with center of mass and facet area
 - handle numeric overflow in qh_normalize and elsewhere
 - merge flipped facets into non-flipped neighbors.
     currently they merge into best neighbor (appears ok)
 - determine min norm for Cramer's rule (qh_sethyperplane_det).  It looks high.
 - improve facet width for very narrow distributions

New features:
 - implement Matlab's tsearch() using Qhull
 - compute volume of Voronoi regions.  You need to determine the dual face
   graph in all dimensions [see Clarkson's hull program]
 - compute alpha shapes [see Clarkson's hull program]
 - implement deletion of Delaunay vertices 
      see Devillers, ACM Symposium on Computational Geometry, Minneapolis 1999.
 - compute largest empty circle [see O'Rourke, chapter 5.5.3] [Hase]
 - list redundant (i.e., coincident) vertices [Spitz]
 - implement Mucke, et al, ['96] for point location in Delaunay triangulations
 - implement convex hull of moving points
 - implement constrained Delaunay diagrams
      see Shewchuk, ACM Symposium on Computational Geometry, Minneapolis 1998.
 - estimate outer volume of hull 
 - automatically determine lower dimensional hulls
 - allow "color" data for input points
      need to insert a coordinate for Delaunay triangulations

Input/output improvements:
 - Support the VTK Visualization Toolkit, http://www.kitware.com/vtk.html
 - generate output data array for Qhull library [Gautier]
 - need improved DOS window with screen fonts, scrollbar, cut/paste
 - generate Geomview output for Voronoi ridges and unbounded rays
 - generate Geomview output for halfspace intersection
 - generate Geomview display of furthest-site Voronoi diagram
 - use 'GDn' to view 5-d facets in 4-d
 - convert Geomview output for other 3-d viewers
 - add interactive output option to avoid recomputing a hull
 - orient vertex neighbors for 'Fv' in 3-d and 2-d
 - track total number of ridges for summary and logging

Performance improvements:
 - optimize Qhull for 2-d Delaunay triangulations
 -   use O'Rourke's '94 vertex->duplicate_edge
 -   add bucketing 
 -   better to specialize all of the code (ca. 2-3x faster w/o merging)
 - use updated LU decomposition to speed up hyperplane construction
 -        [Gill et al. 1974, Math. Comp. 28:505-35]
 - construct hyperplanes from the corresponding horizon/visible facets
 - for merging in high d, do not use vertex->neighbors

Please let us know about your applications and improvements.


Up: Home page for Qhull
Up: Qhull manual: Table of Contents
To: ProgramsOptionsOutputFormatsGeomviewPrintQhullPrecisionTrace
To: Qhull internals: Table of Contents
Dn: Qhull functions, macros, and data structures


The Geometry Center Home Page

Comments to: qhull@qhull.org
Created: Sept. 25, 1995 --- Last modified: see changes.txt