_rotation_spline.py 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460
  1. import numpy as np
  2. from scipy.linalg import solve_banded
  3. from ._rotation import Rotation
  4. def _create_skew_matrix(x):
  5. """Create skew-symmetric matrices corresponding to vectors.
  6. Parameters
  7. ----------
  8. x : ndarray, shape (n, 3)
  9. Set of vectors.
  10. Returns
  11. -------
  12. ndarray, shape (n, 3, 3)
  13. """
  14. result = np.zeros((len(x), 3, 3))
  15. result[:, 0, 1] = -x[:, 2]
  16. result[:, 0, 2] = x[:, 1]
  17. result[:, 1, 0] = x[:, 2]
  18. result[:, 1, 2] = -x[:, 0]
  19. result[:, 2, 0] = -x[:, 1]
  20. result[:, 2, 1] = x[:, 0]
  21. return result
  22. def _matrix_vector_product_of_stacks(A, b):
  23. """Compute the product of stack of matrices and vectors."""
  24. return np.einsum("ijk,ik->ij", A, b)
  25. def _angular_rate_to_rotvec_dot_matrix(rotvecs):
  26. """Compute matrices to transform angular rates to rot. vector derivatives.
  27. The matrices depend on the current attitude represented as a rotation
  28. vector.
  29. Parameters
  30. ----------
  31. rotvecs : ndarray, shape (n, 3)
  32. Set of rotation vectors.
  33. Returns
  34. -------
  35. ndarray, shape (n, 3, 3)
  36. """
  37. norm = np.linalg.norm(rotvecs, axis=1)
  38. k = np.empty_like(norm)
  39. mask = norm > 1e-4
  40. nm = norm[mask]
  41. k[mask] = (1 - 0.5 * nm / np.tan(0.5 * nm)) / nm**2
  42. mask = ~mask
  43. nm = norm[mask]
  44. k[mask] = 1/12 + 1/720 * nm**2
  45. skew = _create_skew_matrix(rotvecs)
  46. result = np.empty((len(rotvecs), 3, 3))
  47. result[:] = np.identity(3)
  48. result[:] += 0.5 * skew
  49. result[:] += k[:, None, None] * np.matmul(skew, skew)
  50. return result
  51. def _rotvec_dot_to_angular_rate_matrix(rotvecs):
  52. """Compute matrices to transform rot. vector derivatives to angular rates.
  53. The matrices depend on the current attitude represented as a rotation
  54. vector.
  55. Parameters
  56. ----------
  57. rotvecs : ndarray, shape (n, 3)
  58. Set of rotation vectors.
  59. Returns
  60. -------
  61. ndarray, shape (n, 3, 3)
  62. """
  63. norm = np.linalg.norm(rotvecs, axis=1)
  64. k1 = np.empty_like(norm)
  65. k2 = np.empty_like(norm)
  66. mask = norm > 1e-4
  67. nm = norm[mask]
  68. k1[mask] = (1 - np.cos(nm)) / nm ** 2
  69. k2[mask] = (nm - np.sin(nm)) / nm ** 3
  70. mask = ~mask
  71. nm = norm[mask]
  72. k1[mask] = 0.5 - nm ** 2 / 24
  73. k2[mask] = 1 / 6 - nm ** 2 / 120
  74. skew = _create_skew_matrix(rotvecs)
  75. result = np.empty((len(rotvecs), 3, 3))
  76. result[:] = np.identity(3)
  77. result[:] -= k1[:, None, None] * skew
  78. result[:] += k2[:, None, None] * np.matmul(skew, skew)
  79. return result
  80. def _angular_acceleration_nonlinear_term(rotvecs, rotvecs_dot):
  81. """Compute the non-linear term in angular acceleration.
  82. The angular acceleration contains a quadratic term with respect to
  83. the derivative of the rotation vector. This function computes that.
  84. Parameters
  85. ----------
  86. rotvecs : ndarray, shape (n, 3)
  87. Set of rotation vectors.
  88. rotvecs_dot : ndarray, shape (n, 3)
  89. Set of rotation vector derivatives.
  90. Returns
  91. -------
  92. ndarray, shape (n, 3)
  93. """
  94. norm = np.linalg.norm(rotvecs, axis=1)
  95. dp = np.sum(rotvecs * rotvecs_dot, axis=1)
  96. cp = np.cross(rotvecs, rotvecs_dot)
  97. ccp = np.cross(rotvecs, cp)
  98. dccp = np.cross(rotvecs_dot, cp)
  99. k1 = np.empty_like(norm)
  100. k2 = np.empty_like(norm)
  101. k3 = np.empty_like(norm)
  102. mask = norm > 1e-4
  103. nm = norm[mask]
  104. k1[mask] = (-nm * np.sin(nm) - 2 * (np.cos(nm) - 1)) / nm ** 4
  105. k2[mask] = (-2 * nm + 3 * np.sin(nm) - nm * np.cos(nm)) / nm ** 5
  106. k3[mask] = (nm - np.sin(nm)) / nm ** 3
  107. mask = ~mask
  108. nm = norm[mask]
  109. k1[mask] = 1/12 - nm ** 2 / 180
  110. k2[mask] = -1/60 + nm ** 2 / 12604
  111. k3[mask] = 1/6 - nm ** 2 / 120
  112. dp = dp[:, None]
  113. k1 = k1[:, None]
  114. k2 = k2[:, None]
  115. k3 = k3[:, None]
  116. return dp * (k1 * cp + k2 * ccp) + k3 * dccp
  117. def _compute_angular_rate(rotvecs, rotvecs_dot):
  118. """Compute angular rates given rotation vectors and its derivatives.
  119. Parameters
  120. ----------
  121. rotvecs : ndarray, shape (n, 3)
  122. Set of rotation vectors.
  123. rotvecs_dot : ndarray, shape (n, 3)
  124. Set of rotation vector derivatives.
  125. Returns
  126. -------
  127. ndarray, shape (n, 3)
  128. """
  129. return _matrix_vector_product_of_stacks(
  130. _rotvec_dot_to_angular_rate_matrix(rotvecs), rotvecs_dot)
  131. def _compute_angular_acceleration(rotvecs, rotvecs_dot, rotvecs_dot_dot):
  132. """Compute angular acceleration given rotation vector and its derivatives.
  133. Parameters
  134. ----------
  135. rotvecs : ndarray, shape (n, 3)
  136. Set of rotation vectors.
  137. rotvecs_dot : ndarray, shape (n, 3)
  138. Set of rotation vector derivatives.
  139. rotvecs_dot_dot : ndarray, shape (n, 3)
  140. Set of rotation vector second derivatives.
  141. Returns
  142. -------
  143. ndarray, shape (n, 3)
  144. """
  145. return (_compute_angular_rate(rotvecs, rotvecs_dot_dot) +
  146. _angular_acceleration_nonlinear_term(rotvecs, rotvecs_dot))
  147. def _create_block_3_diagonal_matrix(A, B, d):
  148. """Create a 3-diagonal block matrix as banded.
  149. The matrix has the following structure:
  150. DB...
  151. ADB..
  152. .ADB.
  153. ..ADB
  154. ...AD
  155. The blocks A, B and D are 3-by-3 matrices. The D matrices has the form
  156. d * I.
  157. Parameters
  158. ----------
  159. A : ndarray, shape (n, 3, 3)
  160. Stack of A blocks.
  161. B : ndarray, shape (n, 3, 3)
  162. Stack of B blocks.
  163. d : ndarray, shape (n + 1,)
  164. Values for diagonal blocks.
  165. Returns
  166. -------
  167. ndarray, shape (11, 3 * (n + 1))
  168. Matrix in the banded form as used by `scipy.linalg.solve_banded`.
  169. """
  170. ind = np.arange(3)
  171. ind_blocks = np.arange(len(A))
  172. A_i = np.empty_like(A, dtype=int)
  173. A_i[:] = ind[:, None]
  174. A_i += 3 * (1 + ind_blocks[:, None, None])
  175. A_j = np.empty_like(A, dtype=int)
  176. A_j[:] = ind
  177. A_j += 3 * ind_blocks[:, None, None]
  178. B_i = np.empty_like(B, dtype=int)
  179. B_i[:] = ind[:, None]
  180. B_i += 3 * ind_blocks[:, None, None]
  181. B_j = np.empty_like(B, dtype=int)
  182. B_j[:] = ind
  183. B_j += 3 * (1 + ind_blocks[:, None, None])
  184. diag_i = diag_j = np.arange(3 * len(d))
  185. i = np.hstack((A_i.ravel(), B_i.ravel(), diag_i))
  186. j = np.hstack((A_j.ravel(), B_j.ravel(), diag_j))
  187. values = np.hstack((A.ravel(), B.ravel(), np.repeat(d, 3)))
  188. u = 5
  189. l = 5
  190. result = np.zeros((u + l + 1, 3 * len(d)))
  191. result[u + i - j, j] = values
  192. return result
  193. class RotationSpline:
  194. """Interpolate rotations with continuous angular rate and acceleration.
  195. The rotation vectors between each consecutive orientation are cubic
  196. functions of time and it is guaranteed that angular rate and acceleration
  197. are continuous. Such interpolation are analogous to cubic spline
  198. interpolation.
  199. Refer to [1]_ for math and implementation details.
  200. Parameters
  201. ----------
  202. times : array_like, shape (N,)
  203. Times of the known rotations. At least 2 times must be specified.
  204. rotations : `Rotation` instance
  205. Rotations to perform the interpolation between. Must contain N
  206. rotations.
  207. Methods
  208. -------
  209. __call__
  210. References
  211. ----------
  212. .. [1] `Smooth Attitude Interpolation
  213. <https://github.com/scipy/scipy/files/2932755/attitude_interpolation.pdf>`_
  214. Examples
  215. --------
  216. >>> from scipy.spatial.transform import Rotation, RotationSpline
  217. >>> import numpy as np
  218. Define the sequence of times and rotations from the Euler angles:
  219. >>> times = [0, 10, 20, 40]
  220. >>> angles = [[-10, 20, 30], [0, 15, 40], [-30, 45, 30], [20, 45, 90]]
  221. >>> rotations = Rotation.from_euler('XYZ', angles, degrees=True)
  222. Create the interpolator object:
  223. >>> spline = RotationSpline(times, rotations)
  224. Interpolate the Euler angles, angular rate and acceleration:
  225. >>> angular_rate = np.rad2deg(spline(times, 1))
  226. >>> angular_acceleration = np.rad2deg(spline(times, 2))
  227. >>> times_plot = np.linspace(times[0], times[-1], 100)
  228. >>> angles_plot = spline(times_plot).as_euler('XYZ', degrees=True)
  229. >>> angular_rate_plot = np.rad2deg(spline(times_plot, 1))
  230. >>> angular_acceleration_plot = np.rad2deg(spline(times_plot, 2))
  231. On this plot you see that Euler angles are continuous and smooth:
  232. >>> import matplotlib.pyplot as plt
  233. >>> plt.plot(times_plot, angles_plot)
  234. >>> plt.plot(times, angles, 'x')
  235. >>> plt.title("Euler angles")
  236. >>> plt.show()
  237. The angular rate is also smooth:
  238. >>> plt.plot(times_plot, angular_rate_plot)
  239. >>> plt.plot(times, angular_rate, 'x')
  240. >>> plt.title("Angular rate")
  241. >>> plt.show()
  242. The angular acceleration is continuous, but not smooth. Also note that
  243. the angular acceleration is not a piecewise-linear function, because
  244. it is different from the second derivative of the rotation vector (which
  245. is a piecewise-linear function as in the cubic spline).
  246. >>> plt.plot(times_plot, angular_acceleration_plot)
  247. >>> plt.plot(times, angular_acceleration, 'x')
  248. >>> plt.title("Angular acceleration")
  249. >>> plt.show()
  250. """
  251. # Parameters for the solver for angular rate.
  252. MAX_ITER = 10
  253. TOL = 1e-9
  254. def _solve_for_angular_rates(self, dt, angular_rates, rotvecs):
  255. angular_rate_first = angular_rates[0].copy()
  256. A = _angular_rate_to_rotvec_dot_matrix(rotvecs)
  257. A_inv = _rotvec_dot_to_angular_rate_matrix(rotvecs)
  258. M = _create_block_3_diagonal_matrix(
  259. 2 * A_inv[1:-1] / dt[1:-1, None, None],
  260. 2 * A[1:-1] / dt[1:-1, None, None],
  261. 4 * (1 / dt[:-1] + 1 / dt[1:]))
  262. b0 = 6 * (rotvecs[:-1] * dt[:-1, None] ** -2 +
  263. rotvecs[1:] * dt[1:, None] ** -2)
  264. b0[0] -= 2 / dt[0] * A_inv[0].dot(angular_rate_first)
  265. b0[-1] -= 2 / dt[-1] * A[-1].dot(angular_rates[-1])
  266. for iteration in range(self.MAX_ITER):
  267. rotvecs_dot = _matrix_vector_product_of_stacks(A, angular_rates)
  268. delta_beta = _angular_acceleration_nonlinear_term(
  269. rotvecs[:-1], rotvecs_dot[:-1])
  270. b = b0 - delta_beta
  271. angular_rates_new = solve_banded((5, 5), M, b.ravel())
  272. angular_rates_new = angular_rates_new.reshape((-1, 3))
  273. delta = np.abs(angular_rates_new - angular_rates[:-1])
  274. angular_rates[:-1] = angular_rates_new
  275. if np.all(delta < self.TOL * (1 + np.abs(angular_rates_new))):
  276. break
  277. rotvecs_dot = _matrix_vector_product_of_stacks(A, angular_rates)
  278. angular_rates = np.vstack((angular_rate_first, angular_rates[:-1]))
  279. return angular_rates, rotvecs_dot
  280. def __init__(self, times, rotations):
  281. from scipy.interpolate import PPoly
  282. if rotations.single:
  283. raise ValueError("`rotations` must be a sequence of rotations.")
  284. if len(rotations) == 1:
  285. raise ValueError("`rotations` must contain at least 2 rotations.")
  286. times = np.asarray(times, dtype=float)
  287. if times.ndim != 1:
  288. raise ValueError("`times` must be 1-dimensional.")
  289. if len(times) != len(rotations):
  290. raise ValueError("Expected number of rotations to be equal to "
  291. "number of timestamps given, got {} rotations "
  292. "and {} timestamps."
  293. .format(len(rotations), len(times)))
  294. dt = np.diff(times)
  295. if np.any(dt <= 0):
  296. raise ValueError("Values in `times` must be in a strictly "
  297. "increasing order.")
  298. rotvecs = (rotations[:-1].inv() * rotations[1:]).as_rotvec()
  299. angular_rates = rotvecs / dt[:, None]
  300. if len(rotations) == 2:
  301. rotvecs_dot = angular_rates
  302. else:
  303. angular_rates, rotvecs_dot = self._solve_for_angular_rates(
  304. dt, angular_rates, rotvecs)
  305. dt = dt[:, None]
  306. coeff = np.empty((4, len(times) - 1, 3))
  307. coeff[0] = (-2 * rotvecs + dt * angular_rates
  308. + dt * rotvecs_dot) / dt ** 3
  309. coeff[1] = (3 * rotvecs - 2 * dt * angular_rates
  310. - dt * rotvecs_dot) / dt ** 2
  311. coeff[2] = angular_rates
  312. coeff[3] = 0
  313. self.times = times
  314. self.rotations = rotations
  315. self.interpolator = PPoly(coeff, times)
  316. def __call__(self, times, order=0):
  317. """Compute interpolated values.
  318. Parameters
  319. ----------
  320. times : float or array_like
  321. Times of interest.
  322. order : {0, 1, 2}, optional
  323. Order of differentiation:
  324. * 0 (default) : return Rotation
  325. * 1 : return the angular rate in rad/sec
  326. * 2 : return the angular acceleration in rad/sec/sec
  327. Returns
  328. -------
  329. Interpolated Rotation, angular rate or acceleration.
  330. """
  331. if order not in [0, 1, 2]:
  332. raise ValueError("`order` must be 0, 1 or 2.")
  333. times = np.asarray(times, dtype=float)
  334. if times.ndim > 1:
  335. raise ValueError("`times` must be at most 1-dimensional.")
  336. singe_time = times.ndim == 0
  337. times = np.atleast_1d(times)
  338. rotvecs = self.interpolator(times)
  339. if order == 0:
  340. index = np.searchsorted(self.times, times, side='right')
  341. index -= 1
  342. index[index < 0] = 0
  343. n_segments = len(self.times) - 1
  344. index[index > n_segments - 1] = n_segments - 1
  345. result = self.rotations[index] * Rotation.from_rotvec(rotvecs)
  346. elif order == 1:
  347. rotvecs_dot = self.interpolator(times, 1)
  348. result = _compute_angular_rate(rotvecs, rotvecs_dot)
  349. elif order == 2:
  350. rotvecs_dot = self.interpolator(times, 1)
  351. rotvecs_dot_dot = self.interpolator(times, 2)
  352. result = _compute_angular_acceleration(rotvecs, rotvecs_dot,
  353. rotvecs_dot_dot)
  354. else:
  355. assert False
  356. if singe_time:
  357. result = result[0]
  358. return result