Up: Home page for Qhull (local)
Up: Qhull Wiki and
FAQ (local)
Up: Qhull manual: contents
To: Imprecision in Qhull
To: Programs
Options
Output
Formats
Geomview
Print
Qhull
Precision
Trace
Functions (local)
To: FAQ: contents
If your question does not appear here, see:
Qhull is a general dimension code for computing convex hulls, Delaunay triangulations, halfspace intersections about a point, Voronoi diagrams, furthest-site Delaunay triangulations, and furthest-site Voronoi diagrams. These structures have applications in science, engineering, statistics, and mathematics. For a detailed introduction, see O'Rourke ['94], Computational Geometry in C.
There are separate programs for each application of Qhull. These programs disable experimental and inappropriate options. If you prefer, you may use Qhull directly. All programs run the same code.
Version 2019.1 adds an experimental option for vertex merging of nearly adjacent vertices ('Q14'). It may resolve topological issues such as "dupridges" with more than two facet neighbors.
Version 2015.1 introduced the reentrant library. It should be used for all code that calls Qhull. The 'qhull' program is built with the reentrant library.
Version 3.1 added triangulated output ('Qt'). It may be used for Delaunay triangulations instead of using joggled input ('QJ').
Brad Barber, Arlington MA, 2019/02/11
Copyright © 1998-2020 C.B. Barber
Within each category, the most recently asked questions are first.
Qhull is a console program. You will first need a command window (i.e., a "command prompt"). You can double click on 'eg\Qhull-go.bat'.
- Type 'qconvex', 'qdelaunay', 'qhalf', 'qvoronoi, 'qhull', and 'rbox' for a synopsis of each program.
- Type 'rbox c D2 | qconvex s i' to compute the convex hull of a square.
- Type 'rbox c D2 | qconvex s i TO results.txt' to write the results to the file 'results.txt'. A summary is still printed on the the console.
- Type 'rbox c D2' to see the input format for qconvex.
- Type 'qconvex < data.txt s i TO results.txt' to read input data from 'data.txt'.
- If you want to enter data by hand, type 'qconvex s i TO results.txt' to read input data from the console. Type in the numbers and end with a ctrl-D.
If you regularly use Qhull on a Windows host, install a bash shell such as
- Git for Windows (wiki, based on MSYS2) -- Git for Windows v2.21 requires arguments for 'qhull', otherwise it waits for stdin. Use 'qhull --help' for a usage note instead of 'qhull'.
- MSYS2 (wiki)
- Cygwin
If you use Windows XP or Windows 8, you may use
Qhull takes its data from standard input (stdin). For example, create a file named 'data.txt' with the following contents:
2 #sample 2-d input 5 #number of points 1 2 #coordinates of points -1.1 3 3 2.2 4 5 -10 -10Then call qconvex with 'qconvex < data.txt'. It will print a summary of the convex hull. Use 'qconvex < data.txt o' to print the vertices and edges. See also input format.
You can generate sample data with rbox. For example, 'rbox 10' generates 10 random points in 3-d. Use a pipe ('|') to run rbox and qhull together, e.g.,
rbox c | qconvex o
computes the convex hull of a cube.
First read:
- Introduction to Qhull
- When to use Qhull
- qconvex -- convex hull
- qdelaunay -- Delaunay triangulation
- qhalf -- half-space intersection about a point
- qvoronoi -- Voronoi diagram
- Rbox, for sample inputs
- Examples of Qhull
Look at Qhull's on-line documentation:
- 'rbox' lists all of the options for generating point sets
- 'qconvex --help' gives a synopsis of qconvex and its options
- 'qconvex -' lists all of the options for qconvex
- 'qconvex .' gives a concise list of options
- 'qdelaunay', 'qhalf', 'qvoronoi', and 'qhull' also have a synopsis and options
Then try out the Qhull programs on small examples.
- 'rbox c' -- lists the vertices of a cube
- 'rbox c D2 | qconvex' -- is the convex hull of a square
- 'rbox c D2 | qconvex o' -- lists the vertices and facets of a square
- 'rbox c | qconvex' -- is the convex hull of a cube
- 'rbox c | qconvex o' -- lists the vertices and facets of a cube
- 'rbox c | qconvex Qt o' -- triangulates the cube
- 'rbox c | qconvex QJ o' -- joggles the input and triangulates the cube
- 'rbox c D4 | qconvex' -- is the convex hull of a hypercube
- 'rbox 6 s D2 t | qconvex p Fx' -- is the convex hull of 6 random, cocircular points. Option 'p' lists the points while option 'Fx' lists the vertices in order.
- 'rbox d D2 c G2 | qdelaunay' -- is the Delaunay triangulation of a diamond and a square. The diamond's vertices are cocircular.
- 'rbox d D2 c G2 | qdelaunay o' -- lists the input sites projected to a paraboloid and the Delaunay regions. The region with 4 vertices is the diamond.
- 'rbox d D2 c G2 | qdelaunay o Qt' -- the cocircular diamond is triangulated as two Delaunay regions.
- 'rbox d D2 c G2 | qdelaunay o QJ' -- the input is joggled and the diamond is triangulated.
- 'rbox d D2 c G2 | qvoronoi o' -- is the Voronoi regions for a diamond and a square. The Voronoi vertex for the diamond is the origin (0,0). Unbounded regions are represented by the first vertex (-10.101 -10.101)
- 'rbox d D2 c G2 | qvoronoi Fv' -- shows the Voronoi diagram for the previous example. Each line is one edge of the diagram. The first number is 4, the next two numbers list a pair of input sites, and the last two numbers list the corresponding pair of Voronoi vertices.
- 'rbox d D2 c G2 | qvoronoi o Qt' -- the cocircular Delaunay region is triangulated. Instead of one Voronoi vertex for the diamond, there are two Voronoi vertices (0,0) and (0,0).
Install Geomview if you are running SGI Irix, Solaris, SunOS, Linux, HP, IBM RS/6000, DEC Alpha, or Next. You can then visualize the output of Qhull. Qhull comes with Geomview examples.
Then try Qhull with a small example of your application. Work out the results by hand. Then experiment with Qhull's options to find the ones that you need.
You will need to decide how Qhull should handle precision problems. It can triangulate the output ('Qt'), joggle the input ('QJ'), or merge facets (the default).
- With triangulated output, Qhull merges facets and triangulates the result.
- With joggle, Qhull produces simplicial (i.e., triangular) output by joggling the input. After joggle, no points are cocircular or cospherical.
- With facet merging, Qhull produces a better approximation than joggle, nor does it modify the input.
- See Merged facets or joggled input.
Use option 'FS' or 'FA'. The area is the area of the surface of the convex hull, while the volume is the total volume of the convex hull.For example,
rbox 10 | qconvex FS 0 2 2.192915621644613 0.2027867899638665 rbox 10 | qconvex FA Convex hull of 10 points in 3-d: Number of vertices: 10 Number of facets: 16 Statistics for: RBOX 10 | QCONVEX FA Number of points processed: 10 Number of hyperplanes created: 28 Number of distance tests for qhull: 44 CPU seconds to compute hull (after input): 0 Total facet area: 2.1929156 Total volume: 0.20278679In 2-d, the convex hull is a polygon. Its surface is the edges of a polygon. So in 2-d, the 'area' is the length of the polygon's edges, while the 'volume' is the area of the polygon.
For example the convex hull of a square,
rbox c D2 | qconvex FS 0 2 4 1 rbox c D2 | qconvex FA Convex hull of 4 points in 2-d: Number of vertices: 4 Number of facets: 4 Statistics for: rbox c D2 | qconvex FA Number of points processed: 4 Number of hyperplanes created: 6 Number of distance tests for qhull: 5 CPU seconds to compute hull (after input): 0 Total facet area: 4 Total volume: 1
Options 'i' (in 4-D and higher) and 'Ft' (in 3-D and higher) use "extra" points for non-simplicial facets (e.g., a face of a cube or hypercube). These points are not part of the convex hull. Options 'i' and 'Ft' triangulate non-simplicial facets using the facet's centrum.
For example, Qhull reports the following for one facet of the convex hull of a hypercube. The facets of a 4-D hypercube are 3-d cubes. Option 'Pd0:0.5' returns the facet along the positive-x axis. Point 17 represents the centrum of this facet. The facet's vertices are eight points: point 8 to point 15
rbox c D4 | qconvex i Pd0:0.5 12 17 13 14 15 17 13 12 14 17 11 13 15 17 14 11 15 17 10 11 14 17 14 12 8 17 12 13 8 17 10 14 8 17 11 10 8 17 13 9 8 17 9 11 8 17 11 9 13 rbox c D4 | qconvex Fx Pd0:0.5 8 8 9 10 11 12 13 14 15The 4-d hypercube has 16 vertices; so point "17" was added by qconvex. Qhull adds the point in order to report a simplicial decomposition of the facet. The point corresponds to the "centrum" which Qhull computes to test for convexity.
Triangulate the output ('Qt') to avoid the extra points. Since the hypercube is 4-d, each simplicial facet is a tetrahedron.
C:\qhull3.1>rbox c D4 | qconvex i Pd0:0.5 Qt 9 9 13 14 15 12 9 13 14 9 11 13 15 11 9 14 15 9 10 11 14 12 9 14 8 9 12 13 8 9 10 14 8 10 9 11 8Use the 'Fv' option to print the vertices of simplicial and non-simplicial facets. For example, here is the same hypercube facet with option 'Fv' instead of 'i':
C:\qhull>rbox c D4 | qconvex Pd0:0.5 Fv 1 8 9 10 12 11 13 14 15 8The coordinates of the extra point are printed with the 'Ft' option. For centrums, option 'Ft' uses indices one less than option 'i'. In this case, point 16 represents the centrum of the facet.
rbox c D4 | qconvex Pd0:0.5 Ft 4 17 12 3 -0.5 -0.5 -0.5 -0.5 -0.5 -0.5 -0.5 0.5 -0.5 -0.5 0.5 -0.5 -0.5 -0.5 0.5 0.5 -0.5 0.5 -0.5 -0.5 -0.5 0.5 -0.5 0.5 -0.5 0.5 0.5 -0.5 -0.5 0.5 0.5 0.5 0.5 -0.5 -0.5 -0.5 0.5 -0.5 -0.5 0.5 0.5 -0.5 0.5 -0.5 0.5 -0.5 0.5 0.5 0.5 0.5 -0.5 -0.5 0.5 0.5 -0.5 0.5 0.5 0.5 0.5 -0.5 0.5 0.5 0.5 0.5 0.5 0 0 0 4 16 13 14 15 4 16 13 12 14 4 16 11 13 15 4 16 14 11 15 4 16 10 11 14 4 16 14 12 8 4 16 12 13 8 4 16 10 14 8 4 16 11 10 8 4 16 13 9 8 4 16 9 11 8 4 16 11 9 13
There's no direct way. You can use option 'FP' to report the distance to the nearest vertex for coplanar input points. Select the minimum distance for a duplicated vertex, and locate all input sites less than this distance.
For Delaunay triangulations, all coplanar points are nearly incident to a vertex. If you want a report of coincident input sites, do not use option 'QJ'. By adding a small random quantity to each input coordinate, it prevents coincident input sites.
Nearly flat triangles occur when boundary points are nearly collinear or coplanar. They also occur for nearly coincident points. Both events can easily occur when using joggle. For example (rbox 10 W0 D2 | qdelaunay QJ Fa) lists the areas of the Delaunay triangles of 10 points on the boundary of a square. Some of these triangles are nearly flat. This occurs when one point is joggled inside of two other points. In this case, nearly flat triangles do not occur with triangulated output (rbox 10 W0 D2 | qdelaunay Qt Fa).
Another example, (rbox c P0 P0 D2 | qdelaunay QJ Fa), computes the areas of the Delaunay triangles for the unit square and two instances of the origin. Four of the triangles have an area of 0.25 while two have an area of 2.0e-11. The later are due to the duplicated origin. With triangulated output (rbox c P0 P0 D2 | qdelaunay Qt Fa) there are four triangles of equal area.
Nearly flat triangles also occur without using joggle. For example, (rbox c P0 P0,0.4999999999 | qdelaunay Fa), computes the areas of the Delaunay triangles for the unit square, a nearly collinear point, and the origin. One triangle has an area of 3.3e-11.
Unfortunately, none of Qhull's merging options remove nearly flat Delaunay triangles due to nearly collinear or coplanar boundary points. The merging options concern the empty circumsphere property of Delaunay triangles. This is independent of the area of the Delaunay triangles. Qhull does handle nearly coincident points.
If you are calling Qhull from a program, you can merge slivers into an adjacent facet. In d dimensions with simplicial facets (e.g., from 'Qt'), each facet has d+1 neighbors. Each neighbor shares d vertices of the facet's d+1 vertices. Let the other vertex be the opposite vertex. For each neighboring facet, if its circumsphere includes the opposite.vertex, the two facets can be merged. [M. Treacy]
You can handle collinear or coplanar boundary points by enclosing the points in a box. For example, (rbox c P0 P0,0.4999999999 c G1 | qdelaunay Fa), surrounds the previous points with [(1,1), (1,-1), (-1,-1), (-1, 1)]. Its Delaunay triangulation does not include a nearly flat triangle. The box also simplifies the graphical output from Qhull.
Without joggle, Qhull lists coincident points as "coplanar" points. For example, (rbox c P0 P0 D2 | qdelaunay Fa), ignores the duplicated origin and lists four triangles of size 0.25. Use 'Fc' to list the coincident points (e.g., rbox c P0 P0 D2 | qdelaunay Fc).
There is no easy way to determine coincident points with joggle. Joggle removes all coincident, cocircular, and cospherical points before running Qhull. Instead use facet merging (the default) or triangulated output ('Qt').
A similar question is "How do I mesh a volume from a set of triangulated surface points?"
This is an instance of the constrained Delaunay Triangulation problem. Qhull does not handle constraints. The boundary of the Delaunay triangulation is always convex. But if the input set contains enough points, the triangulation will include the boundary. The number of points needed depends on the input.
Shewchuk has developed a theory of constrained Delaunay triangulations. See his paper at the 1998 Computational Geometry Conference. Using these ideas, constraints could be added to Qhull. They would have many applications.
There is a large literature on mesh generation and many commercial offerings. For pointers see Owen's International Meshing Roundtable and Schneiders' Finite Element Mesh Generation page.
Yes for convex objects, no for non-convex objects. For non-convex objects, it triangulates the concavities. Unless the object has many points on its surface, triangles may cross the surface.
For points in general position, a 3-d Delaunay triangulation generates tetrahedron. Each face of a tetrahedron is a triangle. For example, the 3-d Delaunay triangulation of random points on the surface of a cube, is a cellular structure of tetrahedron.
Use triangulated output ('qdelaunay Qt i') or joggled input ('qdelaunay QJ i') to generate the Delaunay triangulation. Option 'i' reports each tetrahedron. The triangles are every combination of 3 vertices. Each triangle is a "ridge" of the Delaunay triangulation.
For example,
rbox 10 | qdelaunay Qt i 14 9 5 8 7 0 9 8 7 5 3 8 7 3 0 8 7 5 4 8 1 4 6 8 1 2 9 5 8 4 2 5 8 4 2 9 5 6 2 4 8 9 2 0 8 2 6 0 8 2 4 9 1 2 6 4 1is the Delaunay triangulation of 10 random points. Ridge 9-5-8 occurs twice. Once for tetrahedron 9 5 8 7 and the other for tetrahedron 2 9 5 8.
You can also use the Qhull library to generate the triangles. See 'How do I visit the ridges of a Delaunay triangulation?'
For 3-d Delaunay triangulations with cospherical input sites, use triangulated output ('Qt') or joggled input ('QJ'). Otherwise option 'i' will triangulate non-simplicial facets with the facet's centrum.
If you want non-simplicial output for cospherical sites, use option 'Fv' or 'o'. For option 'o', ignore the last coordinate. It is the lifted coordinate for the corresponding convex hull in 4-d.
The following example is a cube inside a tetrahedron. The 8-vertex facet is the cube. Ignore the last coordinates.
C:\qhull>rbox r y c G0.1 | qdelaunay Fv 4 12 20 44 0.5 0 0 0.3055555555555555 0 0.5 0 0.3055555555555555 0 0 0.5 0.3055555555555555 -0.5 -0.5 -0.5 0.9999999999999999 -0.1 -0.1 -0.1 -6.938893903907228e-018 -0.1 -0.1 0.1 -6.938893903907228e-018 -0.1 0.1 -0.1 -6.938893903907228e-018 -0.1 0.1 0.1 -6.938893903907228e-018 0.1 -0.1 -0.1 -6.938893903907228e-018 0.1 -0.1 0.1 -6.938893903907228e-018 0.1 0.1 -0.1 -6.938893903907228e-018 0.1 0.1 0.1 -6.938893903907228e-018 4 2 11 1 0 4 10 1 0 3 4 11 10 1 0 4 2 9 0 3 4 9 11 2 0 4 7 2 1 3 4 11 7 2 1 4 8 10 0 3 4 9 8 0 3 5 8 9 10 11 0 4 10 6 1 3 4 6 7 1 3 5 6 8 10 4 3 5 6 7 10 11 1 4 5 9 2 3 4 7 5 2 3 5 5 8 9 4 3 5 5 6 7 4 3 8 5 6 8 7 9 10 11 4 5 5 7 9 11 2If you want simplicial output use options 'Qt i' or 'QJ i', e.g.,
rbox r y c G0.1 | qdelaunay Qt i 31 2 11 1 0 11 10 1 0 9 11 2 0 11 7 2 1 8 10 0 3 9 8 0 3 10 6 1 3 6 7 1 3 5 9 2 3 7 5 2 3 9 8 10 11 8 10 11 0 9 8 11 0 6 8 10 4 8 6 10 3 6 8 4 3 6 7 10 11 10 6 11 1 6 7 11 1 8 5 4 3 5 8 9 3 5 6 4 3 6 5 7 3 5 9 10 11 8 5 9 10 7 5 10 11 5 6 7 10 8 5 10 4 5 6 10 4 5 9 11 2 7 5 11 2
To compute the Delaunay triangles indexed by the indices of the input sites, use
rbox 10 D2 | qdelaunay Qt i
To compute the Voronoi vertices and the Voronoi region for each input site, use
rbox 10 D2 | qvoronoi o
To compute each edge ("ridge") of the Voronoi diagram for each pair of adjacent input sites, use
rbox 10 D2 | qvoronoi Fv
To compute the area and volume of the Voronoi region for input site 5 (site 0 is the first one), use
rbox 10 D2 | qvoronoi QV5 p | qconvex s FS
To compute the lines ("hyperplanes") that define the Voronoi region for input site 5, use
orrbox 10 D2 | qvoronoi QV5 p | qconvex n
rbox 10 D2 | qvoronoi QV5 Fi Fo
To list the extreme points of the input sites use
rbox 10 D2 | qdelaunay Fx
You will get the same point ids with
rbox 10 D2 | qconvex Fx
No. This is an immense structure. A triangulation of 19, 16-d points has 43 simplices. If you add one point at a time, the triangulation increased as follows: 43, 189, 523, 1289, 2830, 6071, 11410, 20487. The last triangulation for 26 points used 13 megabytes of memory. When Qhull uses virtual memory, it becomes too slow to use.
For each Voronoi region, compute the convex hull of the region's Voronoi vertices. The volume of each convex hull is the volume of the corresponding Vornoi region.
For example, to compute the volume of the bounded Voronoi region about [0,0,0]: output the origin's Voronoi vertices and compute the volume of their convex hull. The last number from option 'FS' is the volume.
rbox P0 10 | qvoronoi QV0 p | qhull FS 0 2 1.448134756744281 0.1067973560800857For another example, see How do I get the triangles for a 2-d Delaunay triangulation and the vertices of its Voronoi diagram?
This approach is slow if you are using the command line. A faster approcach is to call Qhull from a program. The fastest method is Clarkson's hull program. It computes the volume for all Voronoi regions.
An unbounded Voronoi region does not have a volume.
Use option 'Fi' to list each bisector (i.e. Delaunay ridge). Then compute the minimum distance for each Voronoi vertex.There's other ways to get the same information. Let me know if you find a better method.
Consider a square,
C:\qhull>rbox c D2 2 RBOX c D2 4 -0.5 -0.5 -0.5 0.5 0.5 -0.5 0.5 0.5There's two ways to compute the Voronoi diagram: with facet merging or with joggle. With facet merging, the result is:
C:\qhull>rbox c D2 | qvoronoi Qz Voronoi diagram by the convex hull of 5 points in 3-d: Number of Voronoi regions and at-infinity: 5 Number of Voronoi vertices: 1 Number of facets in hull: 5 Statistics for: RBOX c D2 | QVORONOI Qz Number of points processed: 5 Number of hyperplanes created: 7 Number of distance tests for qhull: 8 Number of merged facets: 1 Number of distance tests for merging: 29 CPU seconds to compute hull (after input): 0 C:\qhull>rbox c D2 | qvoronoi Qz o 2 2 5 1 -10.101 -10.101 0 0 2 0 1 2 0 1 2 0 1 2 0 1 0 C:\qhull>rbox c D2 | qvoronoi Qz Fv 4 4 0 1 0 1 4 0 2 0 1 4 1 3 0 1 4 2 3 0 1There is one Voronoi vertex at the origin and rays from the origin along each of the coordinate axes. The last line '4 2 3 0 1' means that there is a ray that bisects input points #2 and #3 from infinity (vertex 0) to the origin (vertex 1). Option 'Qz' adds an artificial point since the input is cocircular. Coordinates -10.101 indicate the vertex at infinity.
With triangulated output, the Voronoi vertex is duplicated:
C:\qhull3.1>rbox c D2 | qvoronoi Qt Qz Voronoi diagram by the convex hull of 5 points in 3-d: Number of Voronoi regions and at-infinity: 5 Number of Voronoi vertices: 2 Number of triangulated facets: 1 Statistics for: RBOX c D2 | QVORONOI Qt Qz Number of points processed: 5 Number of hyperplanes created: 7 Number of facets in hull: 6 Number of distance tests for qhull: 8 Number of distance tests for merging: 33 Number of distance tests for checking: 30 Number of merged facets: 1 CPU seconds to compute hull (after input): 0.05 C:\qhull3.1>rbox c D2 | qvoronoi Qt Qz o 2 3 5 1 -10.101 -10.101 0 0 0 0 3 2 0 1 2 1 0 2 2 0 3 2 0 1 0 C:\qhull3.1>rbox c D2 | qvoronoi Qt Qz Fv 4 4 0 2 0 2 4 0 1 0 1 4 1 3 0 1 4 2 3 0 2With joggle, the input is no longer cocircular and the Voronoi vertex is split into two:
C:\qhull>rbox c D2 | qvoronoi Qt Qz C:\qhull>rbox c D2 | qvoronoi QJ o 2 3 4 1 -10.101 -10.101 -4.71511718558304e-012 -1.775812830118184e-011 9.020340030474472e-012 -4.02267108512433e-012 2 0 1 3 2 1 0 3 2 0 1 2 2 0 C:\qhull>rbox c D2 | qvoronoi QJ Fv 5 4 0 2 0 1 4 0 1 0 1 4 1 2 1 2 4 1 3 0 2 4 2 3 0 2Note that the Voronoi diagram includes the same rays as before plus a short edge between the two vertices.
Three-dimensional terrain data can be approximated with cospherical points. The Delaunay triangulation of cospherical points is the same as their convex hull. If the points lie on the unit sphere, the facet normals are the Voronoi vertices [via S. Fortune].
For example, consider the points {[1,0,0], [-1,0,0], [0,1,0], ...}. Their convex hull is:
rbox d G1 | qconvex o 3 6 8 12 0 0 -1 0 0 1 0 -1 0 0 1 0 -1 0 0 1 0 0 3 3 1 4 3 1 3 5 3 0 3 4 3 3 0 5 3 2 1 5 3 1 2 4 3 2 0 4 3 0 2 5The facet normals are:
rbox d G1 | qconvex n 4 8 -0.5773502691896258 0.5773502691896258 0.5773502691896258 -0.5773502691896258 0.5773502691896258 0.5773502691896258 0.5773502691896258 -0.5773502691896258 -0.5773502691896258 0.5773502691896258 -0.5773502691896258 -0.5773502691896258 0.5773502691896258 0.5773502691896258 -0.5773502691896258 -0.5773502691896258 0.5773502691896258 -0.5773502691896258 0.5773502691896258 -0.5773502691896258 -0.5773502691896258 -0.5773502691896258 0.5773502691896258 -0.5773502691896258 -0.5773502691896258 -0.5773502691896258 -0.5773502691896258 -0.5773502691896258 0.5773502691896258 -0.5773502691896258 -0.5773502691896258 -0.5773502691896258If you drop the offset from each line (the last number), each line is the Voronoi vertex for the corresponding facet. The neighboring facets for each point define the Voronoi region for each point. For example:
rbox d G1 | qconvex FN 6 4 7 3 2 6 4 5 0 1 4 4 7 4 5 6 4 3 1 0 2 4 6 2 0 5 4 7 3 1 4The Voronoi vertices {7, 3, 2, 6} define the Voronoi region for point 0. Point 0 is [0,0,-1]. Its Voronoi vertices are
-0.5773502691896258 0.5773502691896258 -0.5773502691896258 0.5773502691896258 0.5773502691896258 -0.5773502691896258 -0.5773502691896258 -0.5773502691896258 -0.5773502691896258 0.5773502691896258 -0.5773502691896258 -0.5773502691896258In this case, the Voronoi vertices are oriented, but in general they are unordered.
By taking the dual of the Delaunay triangulation, you can construct the Voronoi diagram. For cospherical points, the convex hull vertices for each facet, define the input sites for each Voronoi vertex. In 3-d, the input sites are oriented. For example:
rbox d G1 | qconvex i 8 3 1 4 1 3 5 0 3 4 3 0 5 2 1 5 1 2 4 2 0 4 0 2 5The convex hull vertices for facet 0 are {3, 1, 4}. So Voronoi vertex 0 (i.e., [-0.577, 0.577, 0.577]) is the Voronoi vertex for input sites {3, 1, 4} (i.e., {[0,1,0], [0,0,1], [-1,0,0]}).
Use 'Fo' to compute the separating hyperplanes for unbounded Voronoi regions. The corresponding ray goes to infinity from the Voronoi vertices. The midpoint between input sites replaces the Voronoi vertex at infinity. Alternatively, if you enclose the input sites in a large enough box, the outermost bounded regions will represent the unbounded regions of the original points.
If you do not box the input sites, you can identify the unbounded regions. They list '0' as a vertex. Vertex 0 represents "infinity". Each unbounded ray includes vertex 0 in option 'Fv. See Voronoi graphics and Voronoi notes.
Qhull may be used to help select a simplex that approximates a data set. It will take experimentation. Geomview will help to visualize the results. This task may be difficult to do in 5-d and higher. Use rbox options 'x' and 'y' to produce random distributions within a simplex. Your methods work if you can recover the simplex.
Use Qhull's precision options to get a first approximation to the hull, say with 10 to 50 facets. For example, try 'C0.05' to remove small facets after constructing the hull. Use 'W0.05' to ignore points within 0.05 of a facet. Use 'PA5' to print the five largest facets by area.
Then use other methods to fit a simplex to this data. Remove outlying vertices with few nearby points. Look for large facets in different quadrants. You can use option 'Pd0d1d2' to print all the facets in a quadrant.
In 4-d and higher, use the outer planes (option 'Fo' or 'facet->maxoutside') since the hyperplane of an approximate facet may be below many of the input points.
For example, consider fitting a cube to 1000 uniformly random points in the unit cube. In this case, the first try was good:
rbox 1000 | qconvex W0.05 C0.05 PA6 Fo 4 6 0.35715408374381 0.08706467018177928 -0.9299788727015564 -0.5985514741284483 0.995841591359023 -0.02512604712761577 0.08756829720435189 -0.5258834069202866 0.02448099521570909 -0.02685210459017302 0.9993396046151313 -0.5158104982631999 -0.9990223929415094 -0.01261133513150079 0.04236994958247349 -0.509218270408407 -0.0128069014364698 -0.9998380680115362 0.01264203427283151 -0.5002512653670584 0.01120895057872914 0.01803671994177704 -0.9997744926535512 -0.5056824072956361
Qhull computes the halfspace intersection about a point. The point must be inside all of the halfspaces. Given a point, a duality turns a halfspace intersection problem into a convex hull problem.
Use linear programming if you do not know a point in the interior of the halfspaces. See the notes for qhalf. You will need a linear programming code. This may require a fair amount of work to implement.
MATLAB
Z. You of MathWorks added qhull to MATLAB 6. See functions convhulln, delaunayn, griddatan, tsearchn, and voronoin. V. Brumberg update MATLAB R14 for Qhull 2003.1 and triangulated output.
Engwirda wrote mesh2d for unstructured mesh generation in MATLAB. It is based on the iterative method of Persson and generally results in better quality meshes than delaunay refinement.
Mathematica and Maple
See qh-math for a Delaunay interface to Mathematica. It includes projects for CodeWarrior on the Macintosh and Visual C++ on Win32 PCs.
The following sample code may produce fewer ridges than expected:facetT *facetp; ridgeT *ridge, **ridgep; FORALLfacets { printf("facet f%d\n", facet->id); FOREACHridge_(facet->ridges) { printf(" ridge r%d between f%d and f%d\n", ridge->id, ridge->top->id, ridge->bottom->id); } }Qhull does not create ridges for simplicial facets. Instead it computes ridges from facet->neighbors. To make ridges for a simplicial facet, use qh_makeridges() in merge.c. Use facet->visit_id to visit each ridge once (instead of twice). For example,
facetT *facet, *neighbor; ridgeT *ridge, **ridgep; qh visit_id++; FORALLfacets { printf("facet f%d\n", facet->id); qh_makeridges(facet); facet->visitId= qh visit_id; FOREACHridge_(facet->ridges) { neighbor= otherfacet_(ridge, visible); if (neighbor->visitid != qh visit_id) printf(" ridge r%d between f%d and f%d\n", ridge->id, ridge->top->id, ridge->bottom->id); } }
You may call Qhull from a program. Please use the reentrant Qhull library (libqhullstatic_r.a, libqhull_r.so, or qhull_r.dll). See user_eg.c and "Qhull-template" in user_r.c for examples.. See Qhull code for an introduction to Qhull's reentrant library and its C++ interface.
Hint: Start with a small example for which you know the answer.
Qhull uses a general-dimension data structure. The size depends on the dimension. Use option 'Ts' to print out the memory statistics [e.g., 'rbox D2 10 | qconvex Ts'].
Qhull's data structures use many pointers. For 64-bit code, pointers are twice the size of integers. For 64-bit code, Qhull uses 50% more memory. It there is not enough memory in the computer's level 1 and level 2 caches, Qhull will run slower as it retrieves data from main memory. A future version of Qhull will include memory and performance improvements for 64-bit code.
The Qhull library may be used to construct convex hulls and Delaunay triangulations one point at a time. It may not be used for deleting points or moving points.
Qhull is designed for batch processing. Neither Clarkson's randomized incremental algorithm nor Qhull are designed for on-line operation. For many applications, it is better to reconstruct the convex hull or Delaunay triangulation from scratch for each new point.
With random point sets and on-line processing, Clarkson's algorithm should run faster than Qhull. Clarkson uses the intermediate facets to reject new, interior points, while Qhull, when used on-line, visits every facet to reject such points. If used on-line for n points, Clarkson may take O(n) times as much memory as the average off-line case, while Qhull's space requirement does not change.
If you triangulate the output before adding all the points (option 'Qt' and procedure qh_triangulate), you must set option 'Q11'. It duplicates the normals of triangulated facets and recomputes the centrums. This should be avoided for regular use since triangulated facets are not clearly convex with their neighbors. It appears to work most of the time, but fails for cases that Qhull normally handles well [see the test call to qh_triangulate in qh_addpoint].
To visit the ridges of a Delaunay triangulation, visit each facet. Each ridge will appear twice since it belongs to two facets. In pseudo-code:
for each facet of the triangulation if the facet is Delaunay (i.e., part of the lower convex hull) for each ridge of the facet if the ridge's neighboring facet has not been visited ... process a ridge of the Delaunay triangulation ...In undebugged, C code:
qh visit_id++; FORALLfacets_(facetlist) if (!facet->upperdelaunay) { facet->visitid= qh visit_id; qh_makeridges(facet); FOREACHridge_(facet->ridges) { neighbor= otherfacet_(ridge, facet); if (neighbor->visitid != qh visit_id) { /* Print ridge here with facet-id and neighbor-id */ /*fprintf(fp, "f%d\tf%d\t",facet->id,neighbor->ID);*/ FOREACHvertex_(ridge->vertices) fprintf(fp,"%d ",qh_pointid (vertex->point) ); qh_printfacetNvertex_simplicial (fp, facet, format); fprintf(fp," "); if(neighbor->upperdelaunay) fprintf(fp," -1 -1 -1 -1 "); else qh_printfacetNvertex_simplicial (fp, neighbor, format); fprintf(fp,"\n"); } } } }
Qhull constructs a Delaunay triangulation by lifting the input sites to a paraboloid. The Delaunay triangulation corresponds to the lower convex hull of the lifted points. To visit each facet of the lower convex hull, use:
facetT *facet; ... FORALLfacets { if (!facet->upperdelaunay) { ... only facets for Delaunay regions ... } }
A point is outside of a facet if it is clearly outside the facet's outer plane. The outer plane is defined by an offset (facet->maxoutside) from the facet's hyperplane.
facetT *facet; pointT *point; realT dist; ... qh_distplane(point, facet, &dist); if (dist > facet->maxoutside + 2 * qh DISTround) { /* point is clearly outside of facet */ }A point is inside of a facet if it is clearly inside the facet's inner plane. The inner plane is computed as the maximum distance of a vertex to the facet. It may be computed for an individual facet, or you may use the maximum over all facets. For example:
facetT *facet; pointT *point; realT dist; ... qh_distplane(point, facet, &dist); if (dist < qh min_vertex - 2 * qh DISTround) { /* point is clearly inside of facet */ }Both tests include two qh.DISTrounds because the computation of the furthest point from a facet may be off by qh.DISTround and the computation of the current distance to the facet may be off by qh.DISTround.
See Locate facet with qh_findbestfacet. For Delaunay triangulations, qh_findbestfacet returns the Delaunay triangle or adjacent triangle that contains the point.
Use qh_findbestfacet(). For example,
coordT point[ DIM ]; boolT isoutside; realT bestdist; facetT *facet; ... set coordinates for point ... facet= qh_findbestfacet (point, qh_ALL, &bestdist, &isoutside); /* 'facet' or an adjacent facet is the closest facet to 'point' */qh_findbestfacet() performs a directed search for the facet furthest below the point. If the point lies inside this facet, qh_findbestfacet() performs an exhaustive search of all facets. An exhaustive search may be needed because a facet on the far side of a lens-shaped distribution may be closer to a point than all of the facet's neighbors. The exhaustive search may be skipped for spherical distributions.
Also see, 'How do I find the Delaunay triangle that is closest to a point?'
A Delaunay triangulation subdivides the plane, or in general dimension, subdivides space. Given a point, how do you determine the subdivision containing the point? Or, given a set of points, how do you determine the subdivision containing each point of the set? Efficiency is important -- an exhaustive search of the subdivision is too slow.
First compute the Delaunay triangle with qh_new_qhull() in user_r.c or Qhull::runQhull(). Lift the point to the paraboloid by summing the squares of the coordinates. Use qh_findbestfacet [poly2_r.c] to find the closest Delaunay facet or adjacent facet. Determine the closest vertex to find the corresponding Voronoi region. Do not use options 'Qbb', 'QbB', 'Qbk:n', or 'QBk:n' since these scale the last coordinate. Optimizations of qh_findbestfacet() should be possible for Delaunay triangulations.
You first need to lift the point to the paraboloid (i.e., the last coordinate is the sum of the squares of the point's coordinates). The routine, qh_setdelaunay() [geom2.c], lifts an array of points to the paraboloid. The following excerpt is from findclosest() in user_eg.c.
coordT point[ DIM + 1]; /* one extra coordinate for lifting the point */ boolT isoutside; realT bestdist; facetT *facet; ... set coordinates for point[] ... qh_setdelaunay (DIM+1, 1, point); facet= qh_findbestfacet (point, qh_ALL, &bestdist, &isoutside); /* 'facet' or an adjacent facet is the closest Delaunay triangle to 'point' */The returned facet either contains the point, or an adjacent facet contains the point, or it is the closest Delaunay triangle along the convex hull of the input set.
Point location is an active research area in Computational Geometry. For a practical approach, see Mucke, et al, "Fast randomized point location without preprocessing in two- and three-dimensional Delaunay triangulations," Computational Geometry '96, p. 274-283, May 1996. For an introduction to planar point location see [O'Rourke '93]. Also see, 'How do I find the facet that is closest to a point?'
To locate the closest Voronoi region, determine the closest vertex of the closest Delaunay triangle.
realT dist, bestdist= REALmax; vertexT *bestvertex= NULL, *vertex, **vertexp; /* 'facet' is the closest Delaunay triangle to 'point' */ FOREACHvertex_( facet->vertices ) { dist= qh_pointdist( point, vertex->point, DIM ); if (dist < bestdist) { bestdist= dist; bestvertex= vertex; } } /* 'bestvertex' represents the Voronoi region closest to 'point'. The corresponding input site is 'bestvertex->point' */
To list the vertices (i.e., extreme points) of the convex hull use
vertexT *vertex; FORALLvertices { ... // vertex->point is the coordinates of the vertex // qh_pointid(vertex->point) is the point ID of the vertex ... }
Compare the output from your program with the output from the Qhull program. Use option 'T1' or 'T4' to trace what Qhull is doing. Prepare a small example for which you know the output. Run the example through the Qhull program and your code. Compare the trace outputs. If you do everything right, the two trace outputs should be almost the same. The trace output will also guide you to the functions that you need to review.
Qhull orients simplicial facets, and prints oriented output for 'i', 'Ft', and other options. The orientation depends on both the vertex order and the flag facet->toporient.
Qhull does not orient non-simplicial facets. Instead it orients the facet's ridges. These are printed with the 'Qt' and 'Ft' option. The facet's hyperplane is oriented.
Up:Home page for Qhull (local)
Up: Qhull Wiki and
FAQ (local)
Up: Qhull manual: contents
To: Imprecision in Qhull
To: Programs
Options
Output
Formats
Geomview
Print
Qhull
Precision
Trace
Functions (local)
To: FAQ: contents
Comments to: qhull@qhull.org
Created:
Sept. 25, 1995 --- Last modified: see top