Module math

Class QuickHull3D

java.lang.Object
de.grogra.math.convexhull.QuickHull3D

public class QuickHull3D extends Object
Computes the convex hull of a set of three dimensional points.

The algorithm is a three dimensional implementation of Quickhull, as described in Barber, Dobkin, and Huhdanpaa, ``The Quickhull Algorithm for Convex Hulls'' (ACM Transactions on Mathematical Software, Vol. 22, No. 4, December 1996), and has a complexity of O(n log(n)) with respect to the number of points. A well-known C implementation of Quickhull that works for arbitrary dimensions is provided by qhull.

A hull is constructed by providing a set of points to either a constructor or a build method. After the hull is built, its vertices and faces can be retrieved using getVertices and getFaces. A typical usage might look like this:

 // x y z coordinates of 6 points
 Point3d[] points = new Point3d[] { new Point3d(0.0, 0.0, 0.0),
                new Point3d(1.0, 0.5, 0.0), new Point3d(2.0, 0.0, 0.0),
                new Point3d(0.5, 0.5, 0.5), new Point3d(0.0, 0.0, 2.0),
                new Point3d(0.1, 0.2, 0.3), new Point3d(0.0, 2.0, 0.0), };
 
 QuickHull3D hull = new QuickHull3D();
 hull.build(points);
 
 System.out.println("Vertices:");
 Point3d[] vertices = hull.getVertices();
 for (int i = 0; i < vertices.length; i++) {
        Point3d pnt = vertices[i];
        System.out.println(pnt.x + " " + pnt.y + " " + pnt.z);
 }
 
 System.out.println("Faces:");
 int[][] faceIndices = hull.getFaces();
 for (int i = 0; i < faceIndices.length; i++) {
        for (int k = 0; k < faceIndices[i].length; k++) {
                System.out.print(faceIndices[i][k] + " ");
        }
        System.out.println("");
 }
 
As a convenience, there are also build and getVertex methods which pass point information using an array of doubles.

Robustness

Because this algorithm uses floating point arithmetic, it is potentially vulnerable to errors arising from numerical imprecision. We address this problem in the same way as qhull, by merging faces whose edges are not clearly convex. A face is convex if its edges are convex, and an edge is convex if the centroid of each adjacent plane is clearly below the plane of the other face. The centroid is considered below a plane if its distance to the plane is less than the negative of a distance tolerance. This tolerance represents the smallest distance that can be reliably computed within the available numeric precision. It is normally computed automatically from the point data, although an application may set this tolerance explicitly.

Numerical problems are more likely to arise in situations where data points lie on or within the faces or edges of the convex hull. We have tested QuickHull3D for such situations by computing the convex hull of a random point set, then adding additional randomly chosen points which lie very close to the hull vertices and edges, and computing the convex hull again. The hull is deemed correct if check returns true. These tests have been successful for a large number of trials and so we are confident that QuickHull3D is reasonably robust.

Merged Faces

The merging of faces means that the faces returned by QuickHull3D may be convex polygons instead of triangles. If triangles are desired, the application may triangulate the faces, but it should be noted that this may result in triangles which are very small or thin and hence difficult to perform reliable convexity tests on. In other words, triangulating a merged face is likely to restore the numerical problems which the merging process removed. Hence is it possible that, after triangulation, check will fail (the same behavior is observed with triangulated output from qhull).

Degenerate Input

It is assumed that the input points are non-degenerate in that they are not coincident, colinear, or colplanar, and thus the convex hull has a non-zero volume. If the input points are detected to be degenerate within the distance tolerance, an IllegalArgumentException will be thrown.
Author:
John E. Lloyd, Fall 2004
  • Field Summary

    Fields
    Modifier and Type
    Field
    Description
    static final double
    Specifies that the distance tolerance should be computed automatically from the input point data.
    protected double
     
    static final int
    Specifies that (on output) vertex indices for a face should be listed in clockwise order.
    protected boolean
     
    protected double
     
    protected Vector<de.grogra.math.convexhull.Face>
     
    protected int
     
    protected Vector<de.grogra.math.convexhull.HalfEdge>
     
    static final int
    Specifies that (on output) the vertex indices for a face should be numbered starting from 1.
    static final int
    Specifies that (on output) the vertex indices for a face should be numbered starting from 0.
    protected int
     
    protected int
     
    protected int
     
    static final int
    Specifies that (on output) the vertex indices for a face should be numbered with respect to the original input points.
    protected de.grogra.math.convexhull.Vertex[]
     
    protected double
     
    protected int[]
     
  • Constructor Summary

    Constructors
    Constructor
    Description
    Creates an empty convex hull object.
    QuickHull3D(double[] coords)
    Creates a convex hull object and initializes it to the convex hull of a set of points whose coordinates are given by an array of doubles.
    QuickHull3D(Point3d[] points)
    Creates a convex hull object and initializes it to the convex hull of a set of points.
  • Method Summary

    Modifier and Type
    Method
    Description
    protected void
    addNewFaces(de.grogra.math.convexhull.FaceList newFaces, de.grogra.math.convexhull.Vertex eyeVtx, Vector<de.grogra.math.convexhull.HalfEdge> horizon)
     
    protected void
    addPointToHull(de.grogra.math.convexhull.Vertex eyeVtx)
     
    void
    build(double[] coords)
    Constructs the convex hull of a set of points whose coordinates are given by an array of doubles.
    void
    build(double[] coords, int nump)
    Constructs the convex hull of a set of points whose coordinates are given by an array of doubles.
    void
    build(Point3d[] points)
    Constructs the convex hull of a set of points.
    void
    build(Point3d[] points, int nump)
    Constructs the convex hull of a set of points.
    protected void
     
    protected void
    calculateHorizon(Point3d eyePnt, de.grogra.math.convexhull.HalfEdge edge0, de.grogra.math.convexhull.Face face, Vector<de.grogra.math.convexhull.HalfEdge> horizon)
     
    boolean
    Checks the correctness of the hull using the distance tolerance returned by getDistanceTolerance; see check(PrintStream,double) for details.
    boolean
    check(PrintStream ps, double tol)
    Checks the correctness of the hull.
    protected boolean
    checkFaceConvexity(de.grogra.math.convexhull.Face face, double tol, PrintStream ps)
     
    protected boolean
    checkFaces(double tol, PrintStream ps)
     
    protected void
     
    protected void
    Creates the initial simplex from which the hull will be built.
    protected void
    deleteFacePoints(de.grogra.math.convexhull.Face face, de.grogra.math.convexhull.Face absorbingFace)
     
    Calculates the centroid of this convex hull.
    boolean
    Returns true if debugging is enabled.
    double
    Returns the distance tolerance that was used for the most recently computed hull.
    static double
    getElement(Tuple3d v, int i)
    Gets a single element of this vector.
    double
    Returns the explicit distance tolerance.
    void
    getFaceIndices(int[] indices, de.grogra.math.convexhull.Face face, int flags)
     
    int[][]
    Returns the faces associated with this hull.
    int[][]
    getFaces(int indexFlags)
    Returns the faces associated with this hull.
    Vector<de.grogra.math.convexhull.Face>
     
    int
    Returns the number of faces in this hull.
    int
    Returns the number of vertices in this hull.
    Returns the set of base points of this hull.
    double
    Calculates the area of the surface generated by this convex hull.
    int[]
    Returns an array specifing the index of each hull vertex with respect to the original input points.
    Returns the vertex points in this hull.
    int
    getVertices(double[] coords)
    Returns the coordinates of the vertex points of this hull.
    double
    Calculates the volume of this convex hull.
    protected void
    initBuffers(int nump)
     
    protected de.grogra.math.convexhull.Vertex
     
    protected double
    oppFaceDistance(de.grogra.math.convexhull.HalfEdge he)
     
    void
    Prints the vertices and faces of this hull to the stream ps.
    void
    print(PrintStream ps, int indexFlags)
    Prints the vertices and faces of this hull to the stream ps.
    protected void
     
    protected void
    resolveUnclaimedPoints(de.grogra.math.convexhull.FaceList newFaces)
     
    void
    setDebug(boolean enable)
    Enables the printing of debugging diagnostics.
    void
    Sets an explicit distance tolerance for convexity tests.
    protected void
    setFromQhull(double[] coords, int nump, boolean triangulate)
     
    protected void
    setHull(double[] coords, int nump, int[][] faceIndices, int numf)
     
    protected void
    setPoints(double[] coords, int nump)
     
    protected void
    setPoints(Point3d[] pnts, int nump)
     
    void
    Triangulates any non-triangular hull faces.

    Methods inherited from class java.lang.Object

    clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
  • Field Details

    • CLOCKWISE

      public static final int CLOCKWISE
      Specifies that (on output) vertex indices for a face should be listed in clockwise order.
      See Also:
    • INDEXED_FROM_ONE

      public static final int INDEXED_FROM_ONE
      Specifies that (on output) the vertex indices for a face should be numbered starting from 1.
      See Also:
    • INDEXED_FROM_ZERO

      public static final int INDEXED_FROM_ZERO
      Specifies that (on output) the vertex indices for a face should be numbered starting from 0.
      See Also:
    • POINT_RELATIVE

      public static final int POINT_RELATIVE
      Specifies that (on output) the vertex indices for a face should be numbered with respect to the original input points.
      See Also:
    • AUTOMATIC_TOLERANCE

      public static final double AUTOMATIC_TOLERANCE
      Specifies that the distance tolerance should be computed automatically from the input point data.
      See Also:
    • findIndex

      protected int findIndex
    • charLength

      protected double charLength
    • debug

      protected boolean debug
    • pointBuffer

      protected de.grogra.math.convexhull.Vertex[] pointBuffer
    • vertexPointIndices

      protected int[] vertexPointIndices
    • faces

      protected Vector<de.grogra.math.convexhull.Face> faces
    • horizon

      protected Vector<de.grogra.math.convexhull.HalfEdge> horizon
    • numVertices

      protected int numVertices
    • numFaces

      protected int numFaces
    • numPoints

      protected int numPoints
    • explicitTolerance

      protected double explicitTolerance
    • tolerance

      protected double tolerance
  • Constructor Details

    • QuickHull3D

      public QuickHull3D()
      Creates an empty convex hull object.
    • QuickHull3D

      public QuickHull3D(double[] coords) throws IllegalArgumentException
      Creates a convex hull object and initializes it to the convex hull of a set of points whose coordinates are given by an array of doubles.
      Parameters:
      coords - x, y, and z coordinates of each input point. The length of this array will be three times the the number of input points.
      Throws:
      IllegalArgumentException - the number of input points is less than four, or the points appear to be coincident, colinear, or coplanar.
    • QuickHull3D

      public QuickHull3D(Point3d[] points) throws IllegalArgumentException
      Creates a convex hull object and initializes it to the convex hull of a set of points.
      Parameters:
      points - input points.
      Throws:
      IllegalArgumentException - the number of input points is less than four, or the points appear to be coincident, colinear, or coplanar.
  • Method Details

    • addNewFaces

      protected void addNewFaces(de.grogra.math.convexhull.FaceList newFaces, de.grogra.math.convexhull.Vertex eyeVtx, Vector<de.grogra.math.convexhull.HalfEdge> horizon)
    • addPointToHull

      protected void addPointToHull(de.grogra.math.convexhull.Vertex eyeVtx)
    • build

      public void build(double[] coords) throws IllegalArgumentException
      Constructs the convex hull of a set of points whose coordinates are given by an array of doubles.
      Parameters:
      coords - x, y, and z coordinates of each input point. The length of this array will be three times the number of input points.
      Throws:
      IllegalArgumentException - the number of input points is less than four, or the points appear to be coincident, colinear, or coplanar.
    • build

      public void build(double[] coords, int nump) throws IllegalArgumentException
      Constructs the convex hull of a set of points whose coordinates are given by an array of doubles.
      Parameters:
      coords - x, y, and z coordinates of each input point. The length of this array must be at least three times nump.
      nump - number of input points
      Throws:
      IllegalArgumentException - the number of input points is less than four or greater than 1/3 the length of coords, or the points appear to be coincident, colinear, or coplanar.
    • build

      public void build(Point3d[] points) throws IllegalArgumentException
      Constructs the convex hull of a set of points.
      Parameters:
      points - input points
      Throws:
      IllegalArgumentException - the number of input points is less than four, or the points appear to be coincident, colinear, or coplanar.
    • build

      public void build(Point3d[] points, int nump) throws IllegalArgumentException
      Constructs the convex hull of a set of points.
      Parameters:
      points - input points
      nump - number of input points
      Throws:
      IllegalArgumentException - the number of input points is less than four or greater then the length of points, or the points appear to be coincident, colinear, or coplanar.
    • buildHull

      protected void buildHull()
    • calculateHorizon

      protected void calculateHorizon(Point3d eyePnt, de.grogra.math.convexhull.HalfEdge edge0, de.grogra.math.convexhull.Face face, Vector<de.grogra.math.convexhull.HalfEdge> horizon)
    • check

      public boolean check(PrintStream ps)
      Checks the correctness of the hull using the distance tolerance returned by getDistanceTolerance; see check(PrintStream,double) for details.
      Parameters:
      ps - print stream for diagnostic messages; may be set to null if no messages are desired.
      Returns:
      true if the hull is valid
      See Also:
    • check

      public boolean check(PrintStream ps, double tol)
      Checks the correctness of the hull. This is done by making sure that no faces are non-convex and that no points are outside any face. These tests are performed using the distance tolerance tol. Faces are considered non-convex if any edge is non-convex, and an edge is non-convex if the centroid of either adjoining face is more than tol above the plane of the other face. Similarly, a point is considered outside a face if its distance to that face's plane is more than 10 times tol.

      If the hull has been triangulated, then this routine may fail if some of the resulting triangles are very small or thin.

      Parameters:
      ps - print stream for diagnostic messages; may be set to null if no messages are desired.
      tol - distance tolerance
      Returns:
      true if the hull is valid
      See Also:
    • checkFaceConvexity

      protected boolean checkFaceConvexity(de.grogra.math.convexhull.Face face, double tol, PrintStream ps)
    • checkFaces

      protected boolean checkFaces(double tol, PrintStream ps)
    • computeMaxAndMin

      protected void computeMaxAndMin()
    • getElement

      public static double getElement(Tuple3d v, int i)
      Gets a single element of this vector. Elements 0, 1, and 2 correspond to x, y, and z.
      Parameters:
      i - element index
      Returns:
      element value throws ArrayIndexOutOfBoundsException if i is not in the range 0 to 2.
    • createInitialSimplex

      protected void createInitialSimplex() throws IllegalArgumentException
      Creates the initial simplex from which the hull will be built.
      Throws:
      IllegalArgumentException
    • deleteFacePoints

      protected void deleteFacePoints(de.grogra.math.convexhull.Face face, de.grogra.math.convexhull.Face absorbingFace)
    • getDebug

      public boolean getDebug()
      Returns true if debugging is enabled.
      Returns:
      true is debugging is enabled
      See Also:
    • getDistanceTolerance

      public double getDistanceTolerance()
      Returns the distance tolerance that was used for the most recently computed hull. The distance tolerance is used to determine when faces are unambiguously convex with respect to each other, and when points are unambiguously above or below a face plane, in the presence of numerical imprecision. Normally, this tolerance is computed automatically for each set of input points, but it can be set explicitly by the application.
      Returns:
      distance tolerance
      See Also:
    • getSurfaceArea

      public double getSurfaceArea()
      Calculates the area of the surface generated by this convex hull.
      Returns:
      surface area
    • getExplicitDistanceTolerance

      public double getExplicitDistanceTolerance()
      Returns the explicit distance tolerance.
      Returns:
      explicit tolerance
      See Also:
    • getFaceIndices

      public void getFaceIndices(int[] indices, de.grogra.math.convexhull.Face face, int flags)
    • getFaces

      public int[][] getFaces()
      Returns the faces associated with this hull.

      Each face is represented by an integer array which gives the indices of the vertices. These indices are numbered relative to the hull vertices, are zero-based, and are arranged counter-clockwise. More control over the index format can be obtained using getFaces(indexFlags).

      Returns:
      array of integer arrays, giving the vertex indices for each face.
      See Also:
    • getFaces

      public int[][] getFaces(int indexFlags)
      Returns the faces associated with this hull.

      Each face is represented by an integer array which gives the indices of the vertices. By default, these indices are numbered with respect to the hull vertices (as opposed to the input points), are zero-based, and are arranged counter-clockwise. However, this can be changed by setting POINT_RELATIVE, INDEXED_FROM_ONE, or CLOCKWISE in the indexFlags parameter.

      Parameters:
      indexFlags - specifies index characteristics (0 results in the default)
      Returns:
      array of integer arrays, giving the vertex indices for each face.
      See Also:
    • getFacesVector

      public Vector<de.grogra.math.convexhull.Face> getFacesVector()
    • getNumFaces

      public int getNumFaces()
      Returns the number of faces in this hull.
      Returns:
      number of faces
    • getNumVertices

      public int getNumVertices()
      Returns the number of vertices in this hull.
      Returns:
      number of vertices
    • getVertexPointIndices

      public int[] getVertexPointIndices()
      Returns an array specifing the index of each hull vertex with respect to the original input points.
      Returns:
      vertex indices with respect to the original points
    • getVertices

      public Point3d[] getVertices()
      Returns the vertex points in this hull.
      Returns:
      array of vertex points
      See Also:
    • getPointBuffer

      public Point3d[] getPointBuffer()
      Returns the set of base points of this hull.
      Returns:
      array of points
    • getVertices

      public int getVertices(double[] coords)
      Returns the coordinates of the vertex points of this hull.
      Parameters:
      coords - returns the x, y, z coordinates of each vertex. This length of this array must be at least three times the number of vertices.
      Returns:
      the number of vertices
      See Also:
    • initBuffers

      protected void initBuffers(int nump)
    • nextPointToAdd

      protected de.grogra.math.convexhull.Vertex nextPointToAdd()
    • getCentroid

      public Point3d getCentroid()
      Calculates the centroid of this convex hull. based on: http://stackoverflow.com/questions/9325303/centroid-of-convex-polyhedron/9325812#9325812 The centroid = SUM(pos*volume)/SUM(volume) when you split the part into finite volumes each with a location pos and volume value volume. This is precisely the calculation done for finding the center of gravity of a composite part.
      Returns:
      point3d of the centroid
    • getVolume

      public double getVolume()
      Calculates the volume of this convex hull. based on: http://stackoverflow.com/questions/9325303/centroid-of-convex-polyhedron/9325812#9325812
      Returns:
      volume
    • oppFaceDistance

      protected double oppFaceDistance(de.grogra.math.convexhull.HalfEdge he)
    • print

      public void print(PrintStream ps)
      Prints the vertices and faces of this hull to the stream ps.

      This is done using the Alias Wavefront .obj file format, with the vertices printed first (each preceding by the letter v), followed by the vertex indices for each face (each preceded by the letter f).

      The face indices are numbered with respect to the hull vertices (as opposed to the input points), with a lowest index of 1, and are arranged counter-clockwise. More control over the index format can be obtained using print(ps,indexFlags).

      Parameters:
      ps - stream used for printing
      See Also:
    • print

      public void print(PrintStream ps, int indexFlags)
      Prints the vertices and faces of this hull to the stream ps.

      This is done using the Alias Wavefront .obj file format, with the vertices printed first (each preceding by the letter v), followed by the vertex indices for each face (each preceded by the letter f).

      By default, the face indices are numbered with respect to the hull vertices (as opposed to the input points), with a lowest index of 1, and are arranged counter-clockwise. However, this can be changed by setting POINT_RELATIVE, INDEXED_FROM_ZERO, or CLOCKWISE in the indexFlags parameter.

      Parameters:
      ps - stream used for printing
      indexFlags - specifies index characteristics (0 results in the default).
      See Also:
    • reindexFacesAndVertices

      protected void reindexFacesAndVertices()
    • resolveUnclaimedPoints

      protected void resolveUnclaimedPoints(de.grogra.math.convexhull.FaceList newFaces)
    • setDebug

      public void setDebug(boolean enable)
      Enables the printing of debugging diagnostics.
      Parameters:
      enable - if true, enables debugging
    • setExplicitDistanceTolerance

      public void setExplicitDistanceTolerance(double tol)
      Sets an explicit distance tolerance for convexity tests. If AUTOMATIC_TOLERANCE is specified (the default), then the tolerance will be computed automatically from the point data.
      Parameters:
      tol - explicit tolerance
      See Also:
    • setFromQhull

      protected void setFromQhull(double[] coords, int nump, boolean triangulate)
    • setHull

      protected void setHull(double[] coords, int nump, int[][] faceIndices, int numf)
    • setPoints

      protected void setPoints(double[] coords, int nump)
    • setPoints

      protected void setPoints(Point3d[] pnts, int nump)
    • triangulate

      public void triangulate()
      Triangulates any non-triangular hull faces. In some cases, due to precision issues, the resulting triangles may be very thin or small, and hence appear to be non-convex (this same limitation is present in qhull).