_geometric_slerp.py 8.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256
  1. __all__ = ['geometric_slerp']
  2. import warnings
  3. from typing import TYPE_CHECKING
  4. import numpy as np
  5. from scipy.spatial.distance import euclidean
  6. if TYPE_CHECKING:
  7. import numpy.typing as npt
  8. def _geometric_slerp(start, end, t):
  9. # create an orthogonal basis using QR decomposition
  10. basis = np.vstack([start, end])
  11. Q, R = np.linalg.qr(basis.T)
  12. signs = 2 * (np.diag(R) >= 0) - 1
  13. Q = Q.T * signs.T[:, np.newaxis]
  14. R = R.T * signs.T[:, np.newaxis]
  15. # calculate the angle between `start` and `end`
  16. c = np.dot(start, end)
  17. s = np.linalg.det(R)
  18. omega = np.arctan2(s, c)
  19. # interpolate
  20. start, end = Q
  21. s = np.sin(t * omega)
  22. c = np.cos(t * omega)
  23. return start * c[:, np.newaxis] + end * s[:, np.newaxis]
  24. def geometric_slerp(
  25. start: "npt.ArrayLike",
  26. end: "npt.ArrayLike",
  27. t: "npt.ArrayLike",
  28. tol: float = 1e-7,
  29. ) -> np.ndarray:
  30. """
  31. Geometric spherical linear interpolation.
  32. The interpolation occurs along a unit-radius
  33. great circle arc in arbitrary dimensional space.
  34. Parameters
  35. ----------
  36. start : (n_dimensions, ) array-like
  37. Single n-dimensional input coordinate in a 1-D array-like
  38. object. `n` must be greater than 1.
  39. end : (n_dimensions, ) array-like
  40. Single n-dimensional input coordinate in a 1-D array-like
  41. object. `n` must be greater than 1.
  42. t : float or (n_points,) 1D array-like
  43. A float or 1D array-like of doubles representing interpolation
  44. parameters. A common approach is to generate the array
  45. with ``np.linspace(0, 1, n_pts)`` for linearly spaced points.
  46. Ascending, descending, and scrambled orders are permitted.
  47. .. versionchanged:: 1.17.0
  48. Extrapolation is permitted, allowing values below `0`
  49. and above `1`.
  50. tol : float
  51. The absolute tolerance for determining if the start and end
  52. coordinates are antipodes.
  53. Returns
  54. -------
  55. result : (t.size, D)
  56. An array of doubles containing the interpolated
  57. spherical path and including start and
  58. end when 0 and 1 t are used. The
  59. interpolated values should correspond to the
  60. same sort order provided in the t array. The result
  61. may be 1-dimensional if ``t`` is a float.
  62. Raises
  63. ------
  64. ValueError
  65. If ``start`` or ``end`` are not on the
  66. unit n-sphere, or for a variety of degenerate conditions.
  67. See Also
  68. --------
  69. scipy.spatial.transform.Slerp : 3-D Slerp that works with quaternions
  70. Notes
  71. -----
  72. The implementation is based on the mathematical formula provided in [1]_,
  73. and the first known presentation of this algorithm, derived from study of
  74. 4-D geometry, is credited to Glenn Davis in a footnote of the original
  75. quaternion Slerp publication by Ken Shoemake [2]_.
  76. .. versionadded:: 1.5.0
  77. References
  78. ----------
  79. .. [1] https://en.wikipedia.org/wiki/Slerp#Geometric_Slerp
  80. .. [2] Ken Shoemake (1985) Animating rotation with quaternion curves.
  81. ACM SIGGRAPH Computer Graphics, 19(3): 245-254.
  82. Examples
  83. --------
  84. Interpolate four linearly-spaced values on the circumference of
  85. a circle spanning 90 degrees:
  86. >>> import numpy as np
  87. >>> from scipy.spatial import geometric_slerp
  88. >>> import matplotlib.pyplot as plt
  89. >>> fig = plt.figure()
  90. >>> ax = fig.add_subplot(111)
  91. >>> start = np.array([1, 0])
  92. >>> end = np.array([0, 1])
  93. >>> t_vals = np.linspace(0, 1, 4)
  94. >>> result = geometric_slerp(start,
  95. ... end,
  96. ... t_vals)
  97. The interpolated results should be at 30 degree intervals
  98. recognizable on the unit circle:
  99. >>> ax.scatter(result[...,0], result[...,1], c='k')
  100. >>> circle = plt.Circle((0, 0), 1, color='grey')
  101. >>> ax.add_artist(circle)
  102. >>> ax.set_aspect('equal')
  103. >>> plt.show()
  104. Attempting to interpolate between antipodes on a circle is
  105. ambiguous because there are two possible paths, and on a
  106. sphere there are infinite possible paths on the geodesic surface.
  107. Nonetheless, one of the ambiguous paths is returned along
  108. with a warning:
  109. >>> import warnings
  110. >>> opposite_pole = np.array([-1, 0])
  111. >>> with warnings.catch_warnings():
  112. ... warnings.simplefilter("ignore", UserWarning)
  113. ... geometric_slerp(start,
  114. ... opposite_pole,
  115. ... t_vals)
  116. array([[ 1.00000000e+00, 0.00000000e+00],
  117. [ 5.00000000e-01, 8.66025404e-01],
  118. [-5.00000000e-01, 8.66025404e-01],
  119. [-1.00000000e+00, 1.22464680e-16]])
  120. Extend the original example to a sphere and plot interpolation
  121. points in 3D:
  122. >>> from mpl_toolkits.mplot3d import proj3d
  123. >>> fig = plt.figure()
  124. >>> ax = fig.add_subplot(111, projection='3d')
  125. Plot the unit sphere for reference (optional):
  126. >>> u = np.linspace(0, 2 * np.pi, 100)
  127. >>> v = np.linspace(0, np.pi, 100)
  128. >>> x = np.outer(np.cos(u), np.sin(v))
  129. >>> y = np.outer(np.sin(u), np.sin(v))
  130. >>> z = np.outer(np.ones(np.size(u)), np.cos(v))
  131. >>> ax.plot_surface(x, y, z, color='y', alpha=0.1)
  132. Interpolating over a larger number of points
  133. may provide the appearance of a smooth curve on
  134. the surface of the sphere, which is also useful
  135. for discretized integration calculations on a
  136. sphere surface:
  137. >>> start = np.array([1, 0, 0])
  138. >>> end = np.array([0, 0, 1])
  139. >>> t_vals = np.linspace(0, 1, 200)
  140. >>> result = geometric_slerp(start,
  141. ... end,
  142. ... t_vals)
  143. >>> ax.plot(result[...,0],
  144. ... result[...,1],
  145. ... result[...,2],
  146. ... c='k')
  147. >>> plt.show()
  148. It is also possible to perform extrapolations outside
  149. the interpolation interval by using interpolation parameter
  150. values below 0 or above 1. For example, the above example
  151. may be adjusted to extrapolate to an antipodal position:
  152. >>> fig = plt.figure()
  153. >>> ax = fig.add_subplot(111, projection='3d')
  154. >>> ax.plot_surface(x, y, z, color='y', alpha=0.1)
  155. >>> start = np.array([1, 0, 0])
  156. >>> end = np.array([0, 0, 1])
  157. >>> t_vals = np.linspace(0, 2, 400)
  158. >>> result = geometric_slerp(start,
  159. ... end,
  160. ... t_vals)
  161. >>> ax.plot(result[...,0], result[...,1], result[...,2], c='k')
  162. >>> plt.show()
  163. """
  164. start = np.asarray(start, dtype=np.float64)
  165. end = np.asarray(end, dtype=np.float64)
  166. t = np.asarray(t)
  167. if t.ndim > 1:
  168. raise ValueError("The interpolation parameter "
  169. "value must be one dimensional.")
  170. if start.ndim != 1 or end.ndim != 1:
  171. raise ValueError("Start and end coordinates "
  172. "must be one-dimensional")
  173. if start.size != end.size:
  174. raise ValueError("The dimensions of start and "
  175. "end must match (have same size)")
  176. if start.size < 2 or end.size < 2:
  177. raise ValueError("The start and end coordinates must "
  178. "both be in at least two-dimensional "
  179. "space")
  180. if np.array_equal(start, end):
  181. return np.linspace(start, start, t.size)
  182. # for points that violate equation for n-sphere
  183. for coord in [start, end]:
  184. if not np.allclose(np.linalg.norm(coord), 1.0,
  185. rtol=1e-9,
  186. atol=0):
  187. raise ValueError("start and end are not"
  188. " on a unit n-sphere")
  189. if not isinstance(tol, float):
  190. raise ValueError("tol must be a float")
  191. else:
  192. tol = np.fabs(tol)
  193. coord_dist = euclidean(start, end)
  194. # diameter of 2 within tolerance means antipodes, which is a problem
  195. # for all unit n-spheres (even the 0-sphere would have an ambiguous path)
  196. if np.allclose(coord_dist, 2.0, rtol=0, atol=tol):
  197. warnings.warn("start and end are antipodes "
  198. "using the specified tolerance; "
  199. "this may cause ambiguous slerp paths",
  200. stacklevel=2)
  201. t = np.asarray(t, dtype=np.float64)
  202. if t.size == 0:
  203. return np.empty((0, start.size))
  204. if t.ndim == 0:
  205. return _geometric_slerp(start,
  206. end,
  207. np.atleast_1d(t)).ravel()
  208. else:
  209. return _geometric_slerp(start,
  210. end,
  211. t)