pyhull package¶
Subpackages¶
Submodules¶
pyhull.convex_hull module¶
This module implements a ConvexHull class.
-
class
ConvexHull
(points, joggle=False)[source]¶ Bases:
object
Convex hull for a set of points.
Initializes a ConvexHull from points.
- Args:
- points ([[float]]): All the points as a sequence of sequences.
- e.g., [[-0.5, -0.5], [-0.5, 0.5], [0.5, -0.5], [0.5, 0.5]]
- joggle (bool): Use qhull option to joggle inputs until simplical
- result is obtained instead of merging facets.
-
simplices
¶ Returns the simplices of the convex hull.
pyhull.delaunay module¶
This module implements a DelaunayTri class representing a Delaunay triangulation of a set of points.
-
class
DelaunayTri
(points, joggle=False)[source]¶ Bases:
object
Delaunay triangulation for a set of points.
Initializes a DelaunayTri from points.
- Args:
- points ([[float]]): All the points as a sequence of sequences.
- e.g., [[-0.5, -0.5], [-0.5, 0.5], [0.5, -0.5], [0.5, 0.5]]
- joggle (bool): Use qhull option to joggle inputs until simplical
- result is obtained instead of merging facets.
-
simplices
¶ Returns the simplices of the triangulation.
pyhull.halfspace module¶
This module defines classes for computing halfspace intersections
-
class
Halfspace
(normal, offset)[source]¶ Bases:
object
A halfspace defined by dot(normal, coords) + offset <= 0
Initializes a Halfspace.
- Args:
- normal: vector normal to hyperplane offset: offset of hyperplane from origin
-
static
from_hyperplane
(basis, origin, point, internal=True)[source]¶ Returns a Halfspace defined by a list of vectors parallel to the bounding hyperplane.
- Args:
- basis: basis for the hyperplane (array with vector rows) origin: point on the hyperplane point: point not on the hyperplane internal: whether point is inside the halfspace
-
class
HalfspaceIntersection
(halfspaces, interior_point)[source]¶ Bases:
object
Uses qhalf to calculate the vertex representation of the intersection of a set of halfspaces
-
facets_by_halfspace
¶ Returns a list of vertex indices for each halfspace e.g: facets_by_halfspace[0] is the list of indices ov vertices incident to halfspace 0
-
facets_by_vertex
¶ Returns a list of non-redundant halfspace indices for each vertex e.g: facets_by_vertex[0] is the list of indices of halfspaces incident to vertex 0
-
vertices
¶ Returns the vertices of the halfspace intersection
-
pyhull.simplex module¶
This module defines a class representing an arbitrary Simplex in arbitrary dimensional space.
-
class
Simplex
(coords)[source]¶ Bases:
object
A generalized simplex object. See http://en.wikipedia.org/wiki/Simplex.
Initializes a Simplex from coordinates.
- Args:
- coords ([[float]]): Coords of the vertices of the simplex. E.g.,
- [[1, 2, 3], [2, 4, 5], [6, 7, 8]].
-
coords
¶ Returns a copy of the vertex coordinates in the simplex.
-
in_simplex
(point, tolerance=1e-08)[source]¶ Checks if a point is in the simplex using the standard barycentric coordinate system algorithm.
Taking an arbitrary vertex as an origin, we compute the basis for the simplex from this origin by subtracting all other vertices from the origin. We then project the point into this coordinate system and determine the linear decomposition coefficients in this coordinate system. If the coeffs satisfy that all coeffs >= 0, the composition is in the facet.
- Args:
- point ([float]): Point to test tolerance (float): Tolerance to test if point is in simplex.
pyhull.voronoi module¶
This module implements a VoronoiTess class representing a Voronoi tessellation of a set of points.
-
class
VoronoiTess
(points, add_bounding_box=False)[source]¶ Bases:
object
Voronoi tessellation for a set of points.
Initializes a VoronoiTess from points.
- Args:
- points ([[float]]): All the points as a sequence of sequences.
- e.g., [[-0.5, -0.5], [-0.5, 0.5], [0.5, -0.5], [0.5, 0.5]]
- add_bounding_box (bool): If True, a hypercube corresponding to
- the extremes of each coordinate will be added to the list of points.
Module contents¶
pyhull is a Python wrapper to Qhull (http://www.qhull.org/) for the computation of the convex hull, Delaunay triangulation and Voronoi diagram.
Useful low-level functions are implemented for direct import in the base package and can be called as pyhull.qconvex, pyhull.qdelauany, etc.
-
qconvex
(options, points)[source]¶ Similar to qconvex command in command-line qhull.
- Args:
- options:
- An option string. Up to two options separated by spaces are supported. See Qhull’s qconvex help for info. Typically used options are: s - summary of results (default) i - vertices incident to each facet n - normals with offsets p - vertex coordinates (includes coplanar points if ‘Qc’)
- points:
- Sequence of points as input.
- Returns:
- Output as a list of strings. E.g., [‘4’, ‘0 2’, ‘1 0’, ‘2 3’, ‘3 1’]
-
qdelaunay
(options, points)[source]¶ Similar to qdelaunay command in command-line qhull.
- Args:
- options:
- An options string. Up to two options separated by spaces are supported. See Qhull’s qdelaunay help for info. Typically used options are: s - summary of results (default) i - vertices incident to each Delaunay region o - OFF format (shows the points lifted to a paraboloid)
- points:
- Sequence of points as input.
- Returns:
- Output as a list of strings. E.g., [‘4’, ‘2 4 0’, ‘4 1 0’, ‘3 4 2’, ‘4 3 1’]
-
qhalf
(options, halfspaces, interior_point)[source]¶ Similar to qvoronoi command in command-line qhull.
- Args:
- option:
- An options string. Up to two options separated by spaces are supported. See Qhull’s qhalf help for info. Typically used options are: Fp
- halfspaces:
- List of Halfspaces as input.
- interior_point:
- An interior point (see qhalf documentation)
- Returns:
- Output as a list of strings. E.g., [‘3’, ‘4’, ‘ 1 1 0 ‘, ‘ 1 -1 2 ‘, ‘ -1 1 2 ‘, ‘ 1 1 2 ‘]
-
qhull_cmd
(cmd, options, points)[source]¶ Generalized helper method to perform a qhull based command.
- Args:
- cmd:
- Command to perform. Supported commands are qconvex, qdelaunay and qvoronoi.
- options:
- Options to be provided for qhull command. See specific methods for info on supported options. Up to two options separated by spaces are supported.
- points:
- Sequence of points as input to qhull command.
- Returns:
- Output as a list of strings. E.g., [‘4’, ‘0 2’, ‘1 0’, ‘2 3 ‘, ‘3 1’]
-
qvoronoi
(options, points)[source]¶ Similar to qvoronoi command in command-line qhull.
- Args:
- option:
- An options string. Up to two options separated by spaces are supported. See Qhull’s qvoronoi help for info. Typically used options are: s - summary of results p - Voronoi vertices o - OFF file format (dim, Voronoi vertices, and Voronoi regions)
- points:
- Sequence of points as input.
- Returns:
- Output as a list of strings. E.g., [‘2’, ‘5 5 1’, ‘-10.101 -10.101’, ‘0 -0.5’, ‘-0.5 0’, ‘0.5 0’, ‘0 0.5’, ‘3 2 0 1’, ‘3 4 0 2’, ‘3 3 0 1’, ‘3 4 0 3’, ‘4 4 2 1 3’]