Up: Home page for Qhull (local)
Up: Qhull manual: contents
To: ProgramsOptionsOutputFormatsGeomviewPrintQhullPrecisionTraceFunctions (local)
To: Imprecision in Qhull
To: Qhull code: contents


[4-d cube] Qhull code

This section discusses the code for Qhull.

Copyright © 1995-2020 C.B. Barber


»Qhull code: contents


»Reentrant Qhull

Qhull-2015 introduces reentrant Qhull (libqhull_r). Reentrant Qhull uses a qhT* argument instead of global data structures. The qhT* pointer is the first argument to most Qhull routines. It allows multiple instances of Qhull to run at the same time. It simplifies the C++ interface to Qhull.

New code should be written with libqhull_r. Existing users of libqhull should consider converting to libqhull_r. Although libqhull will be supported indefinitely, improvements may not be implemented. Reentrant qhull is 1-2% slower than non-reentrant qhull.

Note: Reentrant Qhull is not thread safe. Do not invoke Qhull routines with the same qhT* pointer from multiple threads.

»How to convert code to reentrant Qhull

C++ users need to convert to libqhull_r. The new C++ interface does a better, but not perfect, job of hiding Qhull's C data structures. The previous C++ interface was unusual due to Qhull's global data structures.

All other users should consider converting to libqhull_r. The conversion is straight forward. Most of the changes may be made with global search and replace. The resulting files may be checked via eg/make-qhull_qh.sh. It performs the inverse mapping for comparison with non-reentrant code.

For example, even though the original conversion of libqhull to libqhull_r required thousands of changes, the first run of reentrant Qhull (unix_r.c) produced the same output, and nearly the same log files, as the original, non-reentrant Qhull (unix.c). The original conversion was made without the help of eg/make-qhull_qh.sh. Conversion errors almost always produce compiler errors.

Suggestions to help with conversion.

Suggestions for updating src/libqhull/* from src/libqhull_r/*:

»Qhull on 64-bit computers

Qhull compiles for 64-bit hosts. Since the size of a pointer on a 64-bit host is double the size on a 32-bit host, memory consumption increases about 50% for simplicial facets and up-to 100% for non-simplicial facets. If the convex hull does not fit in the computer's level 1 and level 2 cache memory, Qhull will run slower as it retrieves data from main memory.

If your data fits in 32-bits, run Qhull as 32-bit code. It will use less memory and run faster.

You can check memory consumption with option Ts. It includes the size of each data structure:

For Qhull 2015, the maximum identifier for ridges, vertices, and facets was increased from 24-bits to 32-bits. This allows for larger convex hulls, but may increase the size of the corresponding data structures. The sizes for Qhull 2012.1 were

»Calling Qhull from C++ programs

Qhull's C++ interface allows you to explore the results of running Qhull. It provides access to Qhull's data structures. Most of the classes derive from the corresponding qhull data structure. For example, QhullFacet is an instance of Qhull's facetT.

You can retain most of the data in Qhull and use the C++ interface to explore its results. Each object contains a reference to Qhull's data structure (via QhullQh), making the C++ representation less memory efficient.

Besides using the C++ interface, you can also use libqhull_r directly. For example, the FOREACHfacet_(...) macro will visit each facet in turn.

The C++ interface to Qhull is incomplete. You may need to extend the interface. If so, you will need to understand Qhull's data structures and read the code.

The C++ interface is not documented. You will need to read the code and review user_eg3 and Qhull's test program qhulltest. Please consider documenting the C++ interface with Doxygen or another javadoc-style processor.

user_eg3 demonstrates the C++ interface. For example, user_eg3 eg-100 prints the facets generated by Qhull.

    RboxPoints rbox;
    rbox.appendPoints("100");
    Qhull qhull;
    qhull.runQhull(rbox, "");
    cout << qhull.facetList();

The C++ interface for RboxPoints redefines the fprintf() calls in rboxlib.c. Instead of writing its output to stdout, RboxPoints appends the output to a std::vector.

The same technique may be used for calling Qhull from C++. The class QhullUser provides a starting point. See user_eg3 eg-fifo for a demonstration of Voronoi diagrams.

Since the C++ interface uses reentrant Qhull, multiple threads may run Qhull at the same time. Each thread is one run of Qhull.

Do not have two threads accessing the same Qhull instance. Qhull is not thread-safe.

Qhull 2015 used reentrant Qhull for its C++ interface. If you used the C++ interface from qhull 2012.1, you may need to adjust how you initialize and use the Qhull classes. See How to convert code to reentrant Qhull.

»CoordinateIterator

A CoordinateIterator or ConstCoordinateIterator [RboxPoints.cpp] is a std::vector<realT>::iterator for Rbox and Qhull coordinates. It is the result type of RboxPoints.coordinates().

Qhull does not use CoordinateIterator for its data structures. A point in Qhull is an array of reals instead of a std::vector. See QhullPoint.

»Qhull

Qhull is the top-level class for running Qhull. It initializes Qhull, runs the computation, and records errors. It provides access to the global data structure QhullQh, Qhull's facets, and vertices.

»QhullError

QhullError is derived from std::exception. It reports errors from Qhull and captures the output to stderr.

If error handling is not set up, Qhull exits with a code from 1 to 5. The codes are defined by qh_ERR* in libqhull_r.h. The exit is via qh_exit() in usermem_r.c. The C++ interface does not report the captured output in QhullError. Call Qhull::setErrorStream to send output to cerr instead.

»QhullFacet

A QhullFacet is a facet of the convex hull, a region of the Delaunay triangulation, a vertex of a Voronoi diagram, or an intersection of the halfspace intersection about a point. A QhullFacet has a set of QhullVertex, a set of QhullRidge, and a set of neighboring QhullFacets.

»QhullFacetList

A QhullFacetList is a linked list of QhullFacet. The result of Qhull.runQhull is a QhullFacetList stored in QhullQh.

»QhullFacetSet

A QhullFacetSet is a QhullSet of QhullFacet. QhullFacetSet may be ordered or unordered. The neighboring facets of a QhullFacet is a QhullFacetSet. The neighbors of a QhullFacet is a QhullFacetSet. The neighbors are ordered for simplicial facets, matching the opposite vertex of the facet.

»QhullIterator

QhullIterator contains macros for defining Java-style iterator templates from a STL-style iterator template.

»QhullLinkedList

A QhullLinkedLIst is a template for linked lists with next and previous pointers. QhullFacetList and QhullVertexList are QhullLinkedLists.

»QhullPoint

A QhullPoint is an array of point coordinates, typically doubles. The length of the array is QhullQh.hull_dim. The identifier of a QhullPoint is its 0-based index from QhullQh.first_point followed by QhullQh.other_points.

»QhullPointSet

A QhullPointSet is a QhullSet of QhullPoint. The QhullPointSet of a QhullFacet is its coplanar points.

»QhullQh

QhullQh is the root of Qhull's data structure. It contains initialized constants, sets, buffers, and variables. It contains an array and a set of QhullPoint, a list of QhullFacet, and a list of QhullVertex. The points are the input to Qhull. The facets and vertices are the result of running Qhull.

Qhull's functions access QhullQh through the global variable, qh_qh. The global data structures, qh_stat and qh_mem, record statistics and manage memory respectively.

»QhullRidge

A QhullRidge represents the edge between two QhullFacet's. It is always simplicial with qh.hull_dim-1 QhullVertex)'s.

»QhullRidgeSet

A QhullRidgeSet is a QhullSet of QhullRidge. Each QhullFacet contains a QhullRidgeSet.

»QhullSet

A QhullSet is a set of pointers to objects. QhullSets may be ordered or unordered. They are the core data structure for Qhull.

»QhullVertex

A QhullVertex is a vertex of the convex hull. A simplicial QhullFacet has qh.hull_dim-1 vertices. A QhullVertex contains a QhullPoint. It may list its neighboring QhullFacet's.

»QhullVertexList

A QhullVertexList is a QhullLinkedList of QhullVertex. The global data structure, QhullQh contains a QhullVertexList of all the vertices.

»QhullVertexSet

A QhullVertexSet is a QhullSet of QhullVertex. The QhullVertexSet of a QhullFacet is the vertices of the facet. It is ordered for simplicial facets and unordered for non-simplicial facets.

»QhullUser

QhullUser is a framework for capturing data from any Qhull output. The C++ demonstration program, user_eg3, uses QhullUser for its Voronoi 'eg-fifo' output.

»RboxPoints

RboxPoints is a std::vector of point coordinates (QhullPoint). Its iterator is CoordinateIterator.

RboxPoints.appendPoints() appends points from a variety of distributions such as uniformly distributed within a cube and random points on a sphere. It can also append a cube's vertices or specific points.

»Cpp questions for Qhull

Developing C++ code requires many conventions, idioms, and technical details. The following questions have either mystified the author or do not have a clear answer. See also C++ and Perl Guidelines and 'QH110nn FIX' notes in the code. Please add notes to Qhull Wiki.

»Calling Qhull from C programs

Warning: Qhull was not designed for calling from C programs. You may find the C++ interface easier to use. 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_r.c, user_eg2_r.c, and user_r.c. To see how Qhull calls its library, read unix_r.c, qconvex.c, qdelaun.c, qhalf.c, and qvoronoi.c. The '*_r.c' files are reentrant, otherwise they are non-reentrant. Either version may be used. New code should use reentrant Qhull.

See Functions (local) 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.

If you use non-reentrant Qhull, be aware of the macros "qh" and "qhstat", e.g., "qh hull_dim". They are defined in libqhull.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_r.a, for inclusion in your programs. First review libqhull_r.h. This defines the data structures used by Qhull and provides prototypes for the top-level functions. Most users will only need libqhull_r.h in their programs. For example, the Qhull program is defined with libqhull_r.h and unix_r.c. To access all functions, use qhull_ra.h. Include the file with "#include <libqhull_r/qhull_ra.h>". This avoids potential name conflicts.

Qhull provides build/qhull.pc.in for pkg-config support and CMakeLists.txt for CMake. Using back-ticks, you can compile your C program with Qhull. For example:

  	gcc `pkg-config --cflags --libs qhull_r` -o my_app my_app.c

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.

»How to avoid exit(), fprintf(), stderr, and stdout

Qhull sends output to qh.fout and errors, log messages, and summaries to qh.ferr. qh.fout is normally stdout and qh.ferr is stderr. qh.fout may be redefined by option 'TO' or the caller. qh.ferr may be redirected to qh.fout by option 'Tz'.

Qhull does not use stderr, stdout, fprintf(), or exit() directly.

Qhull reports errors via qh_errexit() by writting a message to qh.ferr and invoking longjmp(). This returns the caller to the corresponding setjmp() (c.f., QH_TRY_ in QhullQh.h). If qh_errexit() is not available, Qhull functions call qh_exit(). qh_exit() normally calls exit(), but may be redefined by the user. An example is libqhullcpp/usermem_r-cpp.cpp. It redefines qh_exit() as a 'throw'.

If qh_meminit() or qh_new_qhull() is called with ferr==NULL, then they set ferr to stderr. Otherwise the Qhull libraries use qh->ferr and qh->qhmem.ferr for error output.

If an error occurs before qh->ferr is initialized, Qhull invokes qh_fprintf_stderr(). The user may redefine this function along with qh_exit(), qh_malloc(), and qh_free().

The Qhull libraries write output via qh_fprintf() [userprintf_r.c]. Otherwise, the Qhull libraries do not use stdout, fprintf(), or printf(). Like qh_exit(), the user may redefine qh_fprintf().

»sets and quick memory allocation

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

qset_r.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_r.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_r.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 either the Delaunay triangle or a neighbor of the Delaunay triangle that contains the lifted point. Qhull determines the Delaunay triangulation by projecting the input sites to a paraboloid. The convex hull matches the Delaunay triangulation at the input sites, but does not match along the edges. See this image by F. Drielsma. A point is green or yellow depending upon the facet returned by qh_findbestfacet. For points near an edge, the circumcircles overlap and the adjacent facet may be returned.

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 was triangulated, qh_findbestfacet() returns one of the corresponding 'tricoplanar' facets. The actual best facet may be a different tricoplanar facet from the same set of facets.

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_r.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 (libqhull.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.

»How to debug Qhull

Qhull continually checks its execution, so most errors will stop Qhull with an error message. Additional checking occurs for verified output ('Tv'), check frequently ('Tc'), check for duplicate ridges ('Q15'), and tracing at level 4 ('T4').

If Qhull detects an error, it writes a descriptive error message to stderr, and exits with an exit status code (see following). The C++ interface captures the message in Qhull::qhullMessage. If Qhull::setErrorStream was called, it writes the error message to Qhull::errorStream.

If a Qhull segfault occurs, turn on tracing with option 'T4' and flush output (qh_fprintf) with option 'Tf'. See core dumps and segfaults.

If Qhull never finishes, is Qhull running slow or was there an infinite loop?

If a Qhull error occurs, try to simplify the problem.

If the segfault, infinite loop, or internal error was due to Qhull, please report the error to 'bradb@shore.net. Please include the input data (i.e., point coordinates) that triggered the error.

»Qhull errors

Qhull errors start with 'QH6...' and Qhull warnings start with 'QH7...'. The error message and error code are arguments to a qh_fprintf call. After printing the error message, Qhull exits with an exit status code. The exit status code indicates the type of error:

»Qhull infinite loops

Except for list traversals, most loops in Qhull are limited by a count or the size of set. Linked lists of facets and vertices terminate with a sentinel whose next element is NULL. If a programming error inserts a link to a previous facet or vertex, an infinite loop occurs on the next traversal. Qhull periodically checks and corrects its linked lists for infinite loops (qh_checklists).

»Qhull trace options

Qhull's trace options are the key to debugging Qhull. They describe an execution of Qhull at various levels of detail, with various options to control what is traced.

These options select when tracing starts or stops. It limits the amount of tracing, especially in high dimensions.

Additional logging by facet id (fnnn), ridge id (rnnn) or vertex id (vnnn), may be enabled by setting qh.tracefacet_id, qh.traceridge_id, or qh.tracevertex_id in global_r.c/qh_initqhull_start2.

»Qhull core dumps and segfaults

If a segfault occurs, use option 'Tf' to flush output after every qh_fprintf. Logging will be significantly slower than normal.

The following debugging plan usually identifies the error

  1. Trace execution at level 1 with flush after each qh_fprintf and output to stdout ('T1 Tf Tz').
  2. Repeat at level 4 after the last qh_addpoint (QH1049, 'TPn'). Add line numbers to the log by piping the output through 'grep -n .'.
  3. Identify the location of the failure using a build of Qhull with debug symbols.
  4. Use the debugger to find relevant facet ids, ridge ids, and vertex ids. These identifiers will appear in the level 4 log.

»eg/qtest.sh for intermittent errors and logging

For intermittent errors, use 'rbox' to generate random test cases, and eg/qtest.sh to invoke multiple runs of qhull. When a failing case is found, rerun eg/qtest.sh with the test case identifier. It produces qhull.log and the corresponding reduced log, qhull-step.log. These logs include line numbers generated by 'grep -n .'

qtest.sh provides the following options

»Memory errors

Qhull checks memory usage before exiting. To locate memory that is not freed ("QH7079 qhull internal warning (main): did not free ..."):
  1. Run qhull with memory tracing 'T5'.
    See 'Level 5' in Qhull trace options (above)
  2. Sort lines that start with 'qh_mem'. It matches qh_memalloc with the corresponding qh_memfree.
  3. For long allocations, sort lines that contain -- qh_mem.*long:
  4. Replace -- qh_mem.*alloc.*\nqh_mem.*free.* -- with 'Match' (Textpad supports internal newlines in match expressions).
  5. Sort by column 25 (n...). It shows unallocated actions. Long allocations are in execution order. Short and quick allocations are in execution order.
  6. For example: qh_mem 0000000000537440 n 10053 alloc long: 128 bytes (tot 484800 cnt 1209)
  7. To check quick vs. long allocations -- grep "qh_mem .*alloc " qhull.log | sed -e 's/^.*long/long/' -e 's/^.*short/short/' -e 's/^.*quick/quick/' -e 's/bytes.*/bytes/' | sort | uniq -c >x.1

Option 'Ts' reports numerous memory statistics.

»Qhull debugging tips

»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. The default options 'C-0' (2-d, 3-d, 4-d) and 'Qx' (5-d and higher) 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.

Additional detail is provided by QH1049 in the level-1 trace ('T1'). For each qh_addpoint, it provides vertex id, facet count, outside point count, CPU time for the previous point, deltas for facets/hyperplanes/distplanes, and the number of retries due to merged pinched vertices. For example:

[QH1049]qh_addpoint: add p260(v176) to hull of 286 facets (1.4e-12 above f830) and 2 outside at 1.192 CPU secs. Previous p710(v175) delta 0.007 CPU, 2 facets, 3 hyperplanes, 443 distplanes, 0 retries

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.

»eg/q_benchmark for optimizing Qhull

eg/q_benchmark and eg/qtest.sh make multiple runs of Qhull for testing, benchmarking, and debugging. They help test and analyze intermittent errors, performance issues, and precision issues. Each release updates eg/q_benchmark-ok.txt.

Qhull 2019.1 is 15% larger than Qhull 2015.2 due to enhanced error reporting, tracing, and facet merging. The increased code size may increase startup times.

Qhull is single threaded. Gcc's gprof works well for profiling Qhull performance.

»Enhancements to Qhull

There are many ways in which Qhull can be improved.

Top Suggestions
 - Document the C++ interface using Doxygen
 - Construct the full Voronoi Diagram using the C++ interface.  See "To do" in Changes.txt
 - Optimize for 64-bit code
   Custom facetT for simplicial facets
   32-bit indices for facets and vertices
 - Bulk qh_addpoint with a custom point partition
 - Full-dimensional flats
   Add points at right angles like 'Qz'
   Ignore added facets in output (cf. f.upperDelaunay and f.good)
 - Per-vertex joggle
   Joggle by random flip of low-order and adjustable-order bits in mantissa
   Allows consistent triangulations across distributed partitions
   Detect integer input data and automatically translate to the origin
 - Develop a theory for merging Qhull's non-simplicial facets
   A merge creates constraints on subsequent merges, what are these constraints?
   Identify topological errors in qh_findbest_test (merge_r.c)f
   Prevent duplicate ridges (Q15-check-duplicates) or facets with the same vertices
   Preserve facet-ridge orientation for nonsimplicial facets (ridge top/bot)
   Produce conformant triangulations for nonsimplicial facets (option 'Qt', QH2088)
   Should vertex merge account for facet orientation?
   Rationalize the merge options qh_RATIOtwisted, qh_WIDEdupridge, etc.
   Should wide merges be proportional to qh.ONEmerge or f.maxoutside?
   Can dupridges be avoided with geometric and topological constraints?
   Review coplanar tests across sharp ridges (coplanar horizon, qh_test_appendmerge, qh_check_maxout)
 - Improve Qhull's computations, particularly qh_setfacetplane for hyperplanes
   Toronto, N., McCarthy, J., "Practically accurate floating-point math,", Computing in
   Science & Engineering, IEEE, July/August 2014, p. 80-95.
 - Octave creates endpoints for unbounded ridges, for drawing Delaunay/Voronoi diagrams [M. Voss]
 - Option to select bounded Voronoi regions [A. Uzunovic]
 - Review Qhull performance.  qh_next_vertexmerge and qh_check_maxout are slower than expected
   Compare to Peterka et al and Li and Snoeyink, particularly 64-bit vs. 32-bit
 - Use Gaussian distribution for random cospherical points in rbox
 - Implement dimension reduction via Johnson-Lindenstrauss flattening
 - Implement bulk qh_addpoint via a subset of the facets, perhaps a box about qh.interior_point
   Allow qh_triangulate to run after each increment [coldfix, scipy #4974]
 - Write incremental addPoint with bucketed inputs and facet search (CGAL)
 - Compute hyperplanes in parallel (cf. qh_setfactplane)
 - Create Voronoi volumes and topology in parallel
 - Implement Delaunay to Voronoi tesselation [Peterka et al, 2014, www.mcs.anl.gov/papers/P5154-0614.pdf]
 - Integrate 4dview with Geomview 1.9.5
 - Use coverage testing to expand Qhull's test programs
 - Add RPM and Debian builds to Qhull (CMake's CPackRMP and CPackDeb).
 - Create a mirror/archive web site for old and new Qhull builds
 - Constrain delaunay triangulations via Shewchuk's algorithm (ACM Symp. Comp. Geo., 1998)

-----------
To do for a furture version of the C++ interface
 - Document C++ using Doxygen conventions (//! and //!<)
 - Add defineAs() to each object
 - Add Qtest::toString() functions for QhullPoint and others.  QByteArray and qstrdup()
 - Add toQVector() for Qt container support.  QVector is preferred over QList
 - Add mutable Java-style iterators for containers.  Limited due to memory-allocation issues.
 - Should Qhull manage the output formats for doubles?  QH11010 FIX: user_r.h defines qh_REAL_1 as %6.8g
 - Allocate memory for QhullSet using Qhull.qhmem.  Create default constructors for QhullVertexSet etc.  Also mid() etc.
 - Add interior point for automatic translation?
 - Write a program with concurrent Qhull
 - Write QhullStat and QhullStat_test
 - Add QList and vector instance of facetT*, etc.
 - Generalize QhullPointSetIterator
 - qh-code.htm: Document changes to C++ interface.
      Organize C++ documentation into collection classes, etc.
 - Review all C++ classes and C++ tests
 - QhullVertexSet uses QhullSetBase::referenceSetT() to free its memory.   Probably needed elsewhere
 - The Boost Graph Library provides C++ classes for graph data structures. It may help
   enhance Qhull's C++ interface [Dr. Dobb's 9/00 p. 29-38; OOPSLA 99 p. 399-414].

[Sep 2020] Suggestions
- In qh_mergefacet, check for nearly adjacent vertices via hash of high-order bits
  If close enough, merge a pair of vertices instead of merging the facets
  Is there an inexpensive way to identify slivers?  It may help select which hyperplane to keep
- qhulltest: Add check for missing Qt dlls.  Qhulltest exits without an error message
- Add C++ documentation for 'calling Qhull'
- Add examples to rbox.htm

[May 2020] Suggestions
- Check that the midpoint for Voronoi option 'Fo' is not a Voronoi vertex (rbox c D2 P0 | qvoronoi Fo)
- How to detect degenerate hyperplanes for Voronoi option 'Fo' and 'Fi'?
  qh_sethyperplane_gauss reports nearzero for axis parallel hyperplanes.
- Add a 'Tv' test for Voronoi option 'Fo' that does not use midpoints

[May 2019] Suggestions
------------
Precision issues
- Improve handling of data with positive, integer coordinates, particularly for Delaunay triangulation
  eg Sterratt's github issue #25
  Add a warning that available precision is reduced
  Add an option to automatically translate the data to the origin
- Review qh.MAXcoplanar ('Un'), it varies by dimension compared to qh.ONEmerge

Topology issues
- Need theory for facet merging, vertex merging, and topological errors
- Does qh_triangulate produce a consistent orientation if qh_renamevertex is not called?

Facet and vertex merging
- Reduce the overhead of qh.NEWtentative ('Q14') and improve the speed of facet and vertex merging
- Review MRGconcavecoplanar and early out for isconcave in qh_test_nonsimplicial_merge
- Review user_r.h ratios and precision constants for merging
  Pre-compute derived precision values (e.g., qh_detmaxoutside)
- Use a fixed target instead of a relative wide-max ratio.
  Why should qh.MAXoutside increase when qh.max_outside increases dramatically
  Why should a slow but steady increase in qh.max_outside be OK?
  Define an option to specify wide-max ratio -- 100x is borderline, bad cases can produce 400x,
- Add statistics for dupridge matching in qh_matchneighbor and qh_matchdupridge.  Report as a "precision problem"
- Collect statistics for MRGdegen and MRGredundant
- In qh_all_merges, why is isreduce set when qh.POSTmerging && qh.hull_dim >= 4?
- In qh_forcedmerges and qh_initmergesets, remove assumption that qh.facet_mergeset is the last temporary set
- Review comparisons for qh_compare_anglemerge and qh_compare_facetmerge (after developing a theory)

Other
- Add a version flag to 'rbox' (cf. global_r.c/qh_version).  Currently, the release date is part of its help prompt.
- Review effect of qh.GOODclosest on qh_buildcone_onlygood ('Qg', QH11030 FIX).  qh_findgood preserves old value if didn't find a good facet.  See qh_findgood_all for disabling
- Review the rules for -REALmax -- they look inconsistent.
  Why REALmax/2 and -REALmax/2?  The comments say 'avoid underflow'.  When was it introduced?
- Review comment in qh_matchnewfacets -- "do not allocate memory after qh.hash_table (need to free it cleanly)"
- Chase circular dependencies when compiling qhulltest with Microsoft Devstudio
  Warning MSB8017 A circular dependency has been detected while executing custom build commands for item "moc\Coordinates_test.moc". This may cause incremental build to work incorrectly.        qhulltest-64    C:\Program Files (x86)\Microsoft Visual Studio\2017\Professional\Common7\IDE\VC\VCTargets\Microsoft.CppCommon.targets   209
- Add 'design:' documentation for poly2_r.c/qh_triangulate
  Consider splitting up
- Geomview for 4-d merge is difficult to understand.  Not able to see the convexity of the edges
- Review memory sizes (mem_r.c/qh_memsize) and quick allocations for 64-bit code
- Review Qhull's collection API conventions, http://www.qhull.org/road/road-faq/xml/qhull-cpp.xml
  See http://gnuwin32.sourceforge.net/packages.html and https://google-styleguide.googlecode.com/svn/trunk/cppguide.html

[Jan 2019] Suggestions
- Optimize poly_r.c/qh_update_vertexneighbors for qh_triangulate. qh_setunique and qh_setcompact are slow
- The size of maxpoints in qh_initialvertices/qh_maxsimplex should be d+3 unique points to help avoid QH6154
- Review coordT vs. realT.  Should parameters and variables be coordT when they are distances or coordinates?
  'coordT' is defined as 'realT'
  Having computations as 'double' with coordinates stored as 'float' requires many type conversions
  Expressions are often computed as 'double' anyway
  Source code sometimes uses 'coordT' and sometimes 'realT'
- Need a separate, hash check for duplicate ridge vertices in a facet list -- faster than current qh_checkfacet
- Add joggle for 'almost incident' vertices (order of 100), may clean up Qt as well, projected to hyperplane
- Consider using r.mergevertex2 to optimize qh_postmerge
- Report two facets with same ridge vertices, opposite orientation (topology error)
  add warning (see QH7084) for duplicate ridge with opposite orientation (only two facets in the list)
- Check 'qh_NOmerge' compiler option

[Jan 2016] Suggestions
------------
 - Add a post-merge pass for Delaunay slivers.  Merge into a neighbor with a circumsphere that includes the opposite point. [M. Treacy]
 - Option to add a bounding box for Delaunay triangulations, e,g., nearly coincident points
 - Rescale output to match 'QbB' on input [J. Metz, 1/30/2014 12:21p]
 - Run through valgrind
 - Notes to compgeom on conformant triangulation and Voronoi volume
 - Implement weighted Delaunay triangulation and weighted Voronoi diagram [A. Liebscher]
   e.g., Sugihara, "Three-dimensional convex hull as a fruitful source of diagrams," Theoretical Computer Science, 2000, 235:325-337
 - testqset: test qh_setdelnth and move-to-front
 - Makefile: Re-review gcc/g++ warnings.  OK in 2011.
 - Break up -Wextra into its components or figure out how to override -Wunused-but-set-variable
   unused-but-set-variable is reporting incorrectly.  All instances are annotated.

 - Can countT be defined as 'int', 'unsigned int', or 64-bit int?
   countT is currently defined as 'int' in qset_r.h
   Vertex ID and ridge ID perhaps should be countT, They are currently 'unsigned'
   Check use of 'int' vs. countT in all cpp code
   Check use of 'int' vs. countT in all c code
   qset_r.h defines countT -- duplicates code in user_r.h -- need to add to qset.h/user.h
   countT -1 used as a flag in Coordinates.mid(), QhullFacet->id()
   Also QhullPoints indexOf and lastIndexOf
   Also QhullPointSet indexOf and lastIndexOf
   Coordinates.indexOf assumes countT is signed (from end)
   Coordinates.lastIndexOf assumes countT is signed (from end)
   All error messages with countT are wrong, convert to int?
   RboxPoints.qh_fprintf_rbox, etc. message 9393 assumes countT but may be int, va_arg(args, countT);  Need to split

[Jan 2010] Suggestions
 - Generate vcproj from qtpro files
   cd qtpro && qmake -spec win32-msvc2005 -tp vc -recursive
   sed -i 's/C\:\/bash\/local\/qhull\/qtpro\///' qhull-all.sln
   Change qhullcpp to libqhull.dll
   Allow both builds on same host (keep /tmp separate)
 - C++ class for access to statistics, accumulate vs. add
 - Add dialog box to RoadError-- a virtual function?
 - Option 'Gt' does not make visible all facets of the mesh example, rbox 32 M1,0,1 | qhull d Gt
 - Merge small volume boundary cells into unbounded regions [Dominik Szczerba]
 - Postmerge with merge options
 - Add modify operators and MutablePointCoordinateIterator to PointCoordinates
 - Fix option Qt for conformant triangulations of merged facets
 - Investigate flipped facet -- rbox 100 s D3 t1263080158 | qhull R1e-3 Tcv Qc
 - Add doc comments to c++ code
 - Measure performance of Qhull, seconds per point by dimension
 - Report potential wraparound of 64-bit ints -- e.g., a large set or points

Documentation
- Qhull::addPoint().  Problems with qh_findbestfacet and otherpoints see
   qh-code.htm#inc on-line construction with qh_addpoint()
- How to handle 64-bit possible loss of data.  WARN64, ptr_intT, size_t/int
- Show custom of qh_fprintf
- cat x.x | grep 'qh_mem ' | sort | awk '{ print $2; }' | uniq -c | grep -vE ' (2|4|6|8|10|12|14|16|20|64|162)[^0-9]'
- qtpro/qhulltest contains .pro and Makefile.  Remove Makefiles by setting shadow directory to ../../tmp/projectname
- Rules for use of qh_qh and multi processes
    UsingQhull
    errorIfAnotherUser
    ~QhullPoints() needs ownership of qh_qh
    Does !qh_pointer work?
    When is qh_qh required?  Minimize the time.
   qhmem, qhstat.ferr
   qhull_inuse==1 when qhull globals active [not useful?]
   rbox_inuse==1 when rbox globals active
   - Multithreaded -- call largest dimension for infinityPoint() and origin()
 - Better documentation for qhmem totshort, freesize, etc.
 - how to change .h, .c, and .cpp to text/html.  OK in Opera
 - QhullVertex.dimension() is not quite correct, epensive
 - Check globalAngleEpsilon
 - Deprecate save_qhull()

[Dec 2003] 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:
 - GPU hardware acceleration, particularly for qh_setplane [D. Reese]
 - 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 meSrging)
 - 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 (local)
Up: Qhull manual: contents
To: ProgramsOptionsOutputFormatsGeomviewPrintQhullPrecisionTraceFunctions (local)
To: Imprecision in Qhull
To: Qhull code: contents


The Geometry Center Home Page

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