_spherical_voronoi.py 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342
  1. """
  2. Spherical Voronoi Code
  3. .. versionadded:: 0.18.0
  4. """
  5. #
  6. # Copyright (C) Tyler Reddy, Ross Hemsley, Edd Edmondson,
  7. # Nikolai Nowaczyk, Joe Pitt-Francis, 2015.
  8. #
  9. # Distributed under the same BSD license as SciPy.
  10. #
  11. import numpy as np
  12. import scipy
  13. from . import _voronoi
  14. from scipy.spatial import cKDTree
  15. __all__ = ['SphericalVoronoi']
  16. def calculate_solid_angles(R):
  17. """Calculates the solid angles of plane triangles. Implements the method of
  18. Van Oosterom and Strackee [VanOosterom]_ with some modifications. Assumes
  19. that input points have unit norm."""
  20. # Original method uses a triple product `R1 . (R2 x R3)` for the numerator.
  21. # This is equal to the determinant of the matrix [R1 R2 R3], which can be
  22. # computed with better stability.
  23. numerator = np.linalg.det(R)
  24. denominator = 1 + (np.einsum('ij,ij->i', R[:, 0], R[:, 1]) +
  25. np.einsum('ij,ij->i', R[:, 1], R[:, 2]) +
  26. np.einsum('ij,ij->i', R[:, 2], R[:, 0]))
  27. return np.abs(2 * np.arctan2(numerator, denominator))
  28. class SphericalVoronoi:
  29. """ Voronoi diagrams on the surface of a sphere.
  30. .. versionadded:: 0.18.0
  31. Parameters
  32. ----------
  33. points : ndarray of floats, shape (npoints, ndim)
  34. Coordinates of points from which to construct a spherical
  35. Voronoi diagram.
  36. radius : float, optional
  37. Radius of the sphere (Default: 1)
  38. center : ndarray of floats, shape (ndim,)
  39. Center of sphere (Default: origin)
  40. threshold : float
  41. Threshold for detecting duplicate points and
  42. mismatches between points and sphere parameters.
  43. (Default: 1e-06)
  44. Attributes
  45. ----------
  46. points : double array of shape (npoints, ndim)
  47. the points in `ndim` dimensions to generate the Voronoi diagram from
  48. radius : double
  49. radius of the sphere
  50. center : double array of shape (ndim,)
  51. center of the sphere
  52. vertices : double array of shape (nvertices, ndim)
  53. Voronoi vertices corresponding to points
  54. regions : list of list of integers of shape (npoints, _ )
  55. the n-th entry is a list consisting of the indices
  56. of the vertices belonging to the n-th point in points
  57. Methods
  58. -------
  59. calculate_areas
  60. Calculates the areas of the Voronoi regions. For 2D point sets, the
  61. regions are circular arcs. The sum of the areas is `2 * pi * radius`.
  62. For 3D point sets, the regions are spherical polygons. The sum of the
  63. areas is `4 * pi * radius**2`.
  64. Raises
  65. ------
  66. ValueError
  67. If there are duplicates in `points`.
  68. If the provided `radius` is not consistent with `points`.
  69. Notes
  70. -----
  71. The spherical Voronoi diagram algorithm proceeds as follows. The Convex
  72. Hull of the input points (generators) is calculated, and is equivalent to
  73. their Delaunay triangulation on the surface of the sphere [Caroli]_.
  74. The Convex Hull neighbour information is then used to
  75. order the Voronoi region vertices around each generator. The latter
  76. approach is substantially less sensitive to floating point issues than
  77. angle-based methods of Voronoi region vertex sorting.
  78. Empirical assessment of spherical Voronoi algorithm performance suggests
  79. quadratic time complexity (loglinear is optimal, but algorithms are more
  80. challenging to implement).
  81. References
  82. ----------
  83. .. [Caroli] Caroli et al. Robust and Efficient Delaunay triangulations of
  84. points on or close to a sphere. Research Report RR-7004, 2009.
  85. .. [VanOosterom] Van Oosterom and Strackee. The solid angle of a plane
  86. triangle. IEEE Transactions on Biomedical Engineering,
  87. 2, 1983, pp 125--126.
  88. See Also
  89. --------
  90. Voronoi : Conventional Voronoi diagrams in N dimensions.
  91. Examples
  92. --------
  93. Do some imports and take some points on a cube:
  94. >>> import numpy as np
  95. >>> import matplotlib.pyplot as plt
  96. >>> from scipy.spatial import SphericalVoronoi, geometric_slerp
  97. >>> from mpl_toolkits.mplot3d import proj3d
  98. >>> # set input data
  99. >>> points = np.array([[0, 0, 1], [0, 0, -1], [1, 0, 0],
  100. ... [0, 1, 0], [0, -1, 0], [-1, 0, 0], ])
  101. Calculate the spherical Voronoi diagram:
  102. >>> radius = 1
  103. >>> center = np.array([0, 0, 0])
  104. >>> sv = SphericalVoronoi(points, radius, center)
  105. Generate plot:
  106. >>> # sort vertices (optional, helpful for plotting)
  107. >>> sv.sort_vertices_of_regions()
  108. >>> t_vals = np.linspace(0, 1, 2000)
  109. >>> fig = plt.figure()
  110. >>> ax = fig.add_subplot(111, projection='3d')
  111. >>> # plot the unit sphere for reference (optional)
  112. >>> u = np.linspace(0, 2 * np.pi, 100)
  113. >>> v = np.linspace(0, np.pi, 100)
  114. >>> x = np.outer(np.cos(u), np.sin(v))
  115. >>> y = np.outer(np.sin(u), np.sin(v))
  116. >>> z = np.outer(np.ones(np.size(u)), np.cos(v))
  117. >>> ax.plot_surface(x, y, z, color='y', alpha=0.1)
  118. >>> # plot generator points
  119. >>> ax.scatter(points[:, 0], points[:, 1], points[:, 2], c='b')
  120. >>> # plot Voronoi vertices
  121. >>> ax.scatter(sv.vertices[:, 0], sv.vertices[:, 1], sv.vertices[:, 2],
  122. ... c='g')
  123. >>> # indicate Voronoi regions (as Euclidean polygons)
  124. >>> for region in sv.regions:
  125. ... n = len(region)
  126. ... for i in range(n):
  127. ... start = sv.vertices[region][i]
  128. ... end = sv.vertices[region][(i + 1) % n]
  129. ... result = geometric_slerp(start, end, t_vals)
  130. ... ax.plot(result[..., 0],
  131. ... result[..., 1],
  132. ... result[..., 2],
  133. ... c='k')
  134. >>> ax.azim = 10
  135. >>> ax.elev = 40
  136. >>> _ = ax.set_xticks([])
  137. >>> _ = ax.set_yticks([])
  138. >>> _ = ax.set_zticks([])
  139. >>> fig.set_size_inches(4, 4)
  140. >>> plt.show()
  141. """
  142. def __init__(self, points, radius=1, center=None, threshold=1e-06):
  143. if radius is None:
  144. raise ValueError('`radius` is `None`. '
  145. 'Please provide a floating point number '
  146. '(i.e. `radius=1`).')
  147. self.radius = float(radius)
  148. self.points = np.array(points).astype(np.double)
  149. self._dim = self.points.shape[1]
  150. if center is None:
  151. self.center = np.zeros(self._dim)
  152. else:
  153. self.center = np.array(center, dtype=float)
  154. # test degenerate input
  155. self._rank = np.linalg.matrix_rank(self.points - self.points[0],
  156. tol=threshold * self.radius)
  157. if self._rank < self._dim:
  158. raise ValueError("Rank of input points must be at least {0}".format(self._dim))
  159. if cKDTree(self.points).query_pairs(threshold * self.radius):
  160. raise ValueError("Duplicate generators present.")
  161. radii = np.linalg.norm(self.points - self.center, axis=1)
  162. max_discrepancy = np.abs(radii - self.radius).max()
  163. if max_discrepancy >= threshold * self.radius:
  164. raise ValueError("Radius inconsistent with generators.")
  165. self._calc_vertices_regions()
  166. def _calc_vertices_regions(self):
  167. """
  168. Calculates the Voronoi vertices and regions of the generators stored
  169. in self.points. The vertices will be stored in self.vertices and the
  170. regions in self.regions.
  171. This algorithm was discussed at PyData London 2015 by
  172. Tyler Reddy, Ross Hemsley and Nikolai Nowaczyk
  173. """
  174. # get Convex Hull
  175. conv = scipy.spatial.ConvexHull(self.points)
  176. # get circumcenters of Convex Hull triangles from facet equations
  177. # for 3D input circumcenters will have shape: (2N-4, 3)
  178. self.vertices = self.radius * conv.equations[:, :-1] + self.center
  179. self._simplices = conv.simplices
  180. # calculate regions from triangulation
  181. # for 3D input simplex_indices will have shape: (2N-4,)
  182. simplex_indices = np.arange(len(self._simplices))
  183. # for 3D input tri_indices will have shape: (6N-12,)
  184. tri_indices = np.column_stack([simplex_indices] * self._dim).ravel()
  185. # for 3D input point_indices will have shape: (6N-12,)
  186. point_indices = self._simplices.ravel()
  187. # for 3D input indices will have shape: (6N-12,)
  188. indices = np.argsort(point_indices, kind='mergesort')
  189. # for 3D input flattened_groups will have shape: (6N-12,)
  190. flattened_groups = tri_indices[indices].astype(np.intp)
  191. # intervals will have shape: (N+1,)
  192. intervals = np.cumsum(np.bincount(point_indices + 1))
  193. # split flattened groups to get nested list of unsorted regions
  194. groups = [list(flattened_groups[intervals[i]:intervals[i + 1]])
  195. for i in range(len(intervals) - 1)]
  196. self.regions = groups
  197. def sort_vertices_of_regions(self):
  198. """Sort indices of the vertices to be (counter-)clockwise ordered.
  199. Raises
  200. ------
  201. TypeError
  202. If the points are not three-dimensional.
  203. Notes
  204. -----
  205. For each region in regions, it sorts the indices of the Voronoi
  206. vertices such that the resulting points are in a clockwise or
  207. counterclockwise order around the generator point.
  208. This is done as follows: Recall that the n-th region in regions
  209. surrounds the n-th generator in points and that the k-th
  210. Voronoi vertex in vertices is the circumcenter of the k-th triangle
  211. in self._simplices. For each region n, we choose the first triangle
  212. (=Voronoi vertex) in self._simplices and a vertex of that triangle
  213. not equal to the center n. These determine a unique neighbor of that
  214. triangle, which is then chosen as the second triangle. The second
  215. triangle will have a unique vertex not equal to the current vertex or
  216. the center. This determines a unique neighbor of the second triangle,
  217. which is then chosen as the third triangle and so forth. We proceed
  218. through all the triangles (=Voronoi vertices) belonging to the
  219. generator in points and obtain a sorted version of the vertices
  220. of its surrounding region.
  221. """
  222. if self._dim != 3:
  223. raise TypeError("Only supported for three-dimensional point sets")
  224. _voronoi.sort_vertices_of_regions(self._simplices, self.regions)
  225. def _calculate_areas_3d(self):
  226. self.sort_vertices_of_regions()
  227. sizes = [len(region) for region in self.regions]
  228. csizes = np.cumsum(sizes)
  229. num_regions = csizes[-1]
  230. # We create a set of triangles consisting of one point and two Voronoi
  231. # vertices. The vertices of each triangle are adjacent in the sorted
  232. # regions list.
  233. point_indices = [i for i, size in enumerate(sizes)
  234. for j in range(size)]
  235. nbrs1 = np.array([r for region in self.regions for r in region])
  236. # The calculation of nbrs2 is a vectorized version of:
  237. # np.array([r for region in self.regions for r in np.roll(region, 1)])
  238. nbrs2 = np.roll(nbrs1, 1)
  239. indices = np.roll(csizes, 1)
  240. indices[0] = 0
  241. nbrs2[indices] = nbrs1[csizes - 1]
  242. # Normalize points and vertices.
  243. pnormalized = (self.points - self.center) / self.radius
  244. vnormalized = (self.vertices - self.center) / self.radius
  245. # Create the complete set of triangles and calculate their solid angles
  246. triangles = np.hstack([pnormalized[point_indices],
  247. vnormalized[nbrs1],
  248. vnormalized[nbrs2]
  249. ]).reshape((num_regions, 3, 3))
  250. triangle_solid_angles = calculate_solid_angles(triangles)
  251. # Sum the solid angles of the triangles in each region
  252. solid_angles = np.cumsum(triangle_solid_angles)[csizes - 1]
  253. solid_angles[1:] -= solid_angles[:-1]
  254. # Get polygon areas using A = omega * r**2
  255. return solid_angles * self.radius**2
  256. def _calculate_areas_2d(self):
  257. # Find start and end points of arcs
  258. arcs = self.points[self._simplices] - self.center
  259. # Calculate the angle subtended by arcs
  260. cosine = np.einsum('ij,ij->i', arcs[:, 0], arcs[:, 1])
  261. sine = np.abs(np.linalg.det(arcs))
  262. theta = np.arctan2(sine, cosine)
  263. # Get areas using A = r * theta
  264. areas = self.radius * theta
  265. # Correct arcs which go the wrong way (single-hemisphere inputs)
  266. signs = np.sign(np.einsum('ij,ij->i', arcs[:, 0],
  267. self.vertices - self.center))
  268. indices = np.where(signs < 0)
  269. areas[indices] = 2 * np.pi * self.radius - areas[indices]
  270. return areas
  271. def calculate_areas(self):
  272. """Calculates the areas of the Voronoi regions.
  273. For 2D point sets, the regions are circular arcs. The sum of the areas
  274. is `2 * pi * radius`.
  275. For 3D point sets, the regions are spherical polygons. The sum of the
  276. areas is `4 * pi * radius**2`.
  277. .. versionadded:: 1.5.0
  278. Returns
  279. -------
  280. areas : double array of shape (npoints,)
  281. The areas of the Voronoi regions.
  282. """
  283. if self._dim == 2:
  284. return self._calculate_areas_2d()
  285. elif self._dim == 3:
  286. return self._calculate_areas_3d()
  287. else:
  288. raise TypeError("Only supported for 2D and 3D point sets")