transforms.py 64 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871187218731874187518761877187818791880188118821883188418851886188718881889189018911892189318941895189618971898189919001901190219031904190519061907190819091910191119121913191419151916
  1. # -*- coding: utf-8 -*-
  2. # transformations.py
  3. # Copyright (c) 2006-2015, Christoph Gohlke
  4. # Copyright (c) 2006-2015, The Regents of the University of California
  5. # Produced at the Laboratory for Fluorescence Dynamics
  6. # All rights reserved.
  7. #
  8. # Redistribution and use in source and binary forms, with or without
  9. # modification, are permitted provided that the following conditions are met:
  10. #
  11. # * Redistributions of source code must retain the above copyright
  12. # notice, this list of conditions and the following disclaimer.
  13. # * Redistributions in binary form must reproduce the above copyright
  14. # notice, this list of conditions and the following disclaimer in the
  15. # documentation and/or other materials provided with the distribution.
  16. # * Neither the name of the copyright holders nor the names of any
  17. # contributors may be used to endorse or promote products derived
  18. # from this software without specific prior written permission.
  19. #
  20. # THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
  21. # AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
  22. # IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
  23. # ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
  24. # LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
  25. # CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
  26. # SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
  27. # INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
  28. # CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
  29. # ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
  30. # POSSIBILITY OF SUCH DAMAGE.
  31. """Homogeneous Transformation Matrices and Quaternions.
  32. A library for calculating 4x4 matrices for translating, rotating, reflecting,
  33. scaling, shearing, projecting, orthogonalizing, and superimposing arrays of
  34. 3D homogeneous coordinates as well as for converting between rotation matrices,
  35. Euler angles, and quaternions. Also includes an Arcball control object and
  36. functions to decompose transformation matrices.
  37. :Author:
  38. `Christoph Gohlke <http://www.lfd.uci.edu/~gohlke/>`_
  39. :Organization:
  40. Laboratory for Fluorescence Dynamics, University of California, Irvine
  41. :Version: 2015.01.29
  42. Requirements
  43. ------------
  44. * `CPython 2.7 or 3.4 <http://www.python.org>`_
  45. * `Numpy 1.8 <http://www.numpy.org>`_
  46. * `Transformations.c 2015.01.29 <http://www.lfd.uci.edu/~gohlke/>`_
  47. (recommended for speedup of some functions)
  48. Notes
  49. -----
  50. The API is not stable yet and is expected to change between revisions.
  51. This Python code is not optimized for speed. Refer to the transformations.c
  52. module for a faster implementation of some functions.
  53. Documentation in HTML format can be generated with epydoc.
  54. Matrices (M) can be inverted using numpy.linalg.inv(M), be concatenated using
  55. numpy.dot(M0, M1), or transform homogeneous coordinate arrays (v) using
  56. numpy.dot(M, v) for shape (4, \*) column vectors, respectively
  57. numpy.dot(v, M.T) for shape (\*, 4) row vectors ("array of points").
  58. This module follows the "column vectors on the right" and "row major storage"
  59. (C contiguous) conventions. The translation components are in the right column
  60. of the transformation matrix, i.e. M[:3, 3].
  61. The transpose of the transformation matrices may have to be used to interface
  62. with other graphics systems, e.g. with OpenGL's glMultMatrixd(). See also [16].
  63. Calculations are carried out with numpy.float64 precision.
  64. Vector, point, quaternion, and matrix function arguments are expected to be
  65. "array like", i.e. tuple, list, or numpy arrays.
  66. Return types are numpy arrays unless specified otherwise.
  67. Angles are in radians unless specified otherwise.
  68. Quaternions w+ix+jy+kz are represented as [w, x, y, z].
  69. A triple of Euler angles can be applied/interpreted in 24 ways, which can
  70. be specified using a 4 character string or encoded 4-tuple:
  71. *Axes 4-string*: e.g. 'sxyz' or 'ryxy'
  72. - first character : rotations are applied to 's'tatic or 'r'otating frame
  73. - remaining characters : successive rotation axis 'x', 'y', or 'z'
  74. *Axes 4-tuple*: e.g. (0, 0, 0, 0) or (1, 1, 1, 1)
  75. - inner axis: code of axis ('x':0, 'y':1, 'z':2) of rightmost matrix.
  76. - parity : even (0) if inner axis 'x' is followed by 'y', 'y' is followed
  77. by 'z', or 'z' is followed by 'x'. Otherwise odd (1).
  78. - repetition : first and last axis are same (1) or different (0).
  79. - frame : rotations are applied to static (0) or rotating (1) frame.
  80. Other Python packages and modules for 3D transformations and quaternions:
  81. * `Transforms3d <https://pypi.python.org/pypi/transforms3d>`_
  82. includes most code of this module.
  83. * `Blender.mathutils <http://www.blender.org/api/blender_python_api>`_
  84. * `numpy-dtypes <https://github.com/numpy/numpy-dtypes>`_
  85. References
  86. ----------
  87. (1) Matrices and transformations. Ronald Goldman.
  88. In "Graphics Gems I", pp 472-475. Morgan Kaufmann, 1990.
  89. (2) More matrices and transformations: shear and pseudo-perspective.
  90. Ronald Goldman. In "Graphics Gems II", pp 320-323. Morgan Kaufmann, 1991.
  91. (3) Decomposing a matrix into simple transformations. Spencer Thomas.
  92. In "Graphics Gems II", pp 320-323. Morgan Kaufmann, 1991.
  93. (4) Recovering the data from the transformation matrix. Ronald Goldman.
  94. In "Graphics Gems II", pp 324-331. Morgan Kaufmann, 1991.
  95. (5) Euler angle conversion. Ken Shoemake.
  96. In "Graphics Gems IV", pp 222-229. Morgan Kaufmann, 1994.
  97. (6) Arcball rotation control. Ken Shoemake.
  98. In "Graphics Gems IV", pp 175-192. Morgan Kaufmann, 1994.
  99. (7) Representing attitude: Euler angles, unit quaternions, and rotation
  100. vectors. James Diebel. 2006.
  101. (8) A discussion of the solution for the best rotation to relate two sets
  102. of vectors. W Kabsch. Acta Cryst. 1978. A34, 827-828.
  103. (9) Closed-form solution of absolute orientation using unit quaternions.
  104. BKP Horn. J Opt Soc Am A. 1987. 4(4):629-642.
  105. (10) Quaternions. Ken Shoemake.
  106. http://www.sfu.ca/~jwa3/cmpt461/files/quatut.pdf
  107. (11) From quaternion to matrix and back. JMP van Waveren. 2005.
  108. http://www.intel.com/cd/ids/developer/asmo-na/eng/293748.htm
  109. (12) Uniform random rotations. Ken Shoemake.
  110. In "Graphics Gems III", pp 124-132. Morgan Kaufmann, 1992.
  111. (13) Quaternion in molecular modeling. CFF Karney.
  112. J Mol Graph Mod, 25(5):595-604
  113. (14) New method for extracting the quaternion from a rotation matrix.
  114. Itzhack Y Bar-Itzhack, J Guid Contr Dynam. 2000. 23(6): 1085-1087.
  115. (15) Multiple View Geometry in Computer Vision. Hartley and Zissermann.
  116. Cambridge University Press; 2nd Ed. 2004. Chapter 4, Algorithm 4.7, p 130.
  117. (16) Column Vectors vs. Row Vectors.
  118. http://steve.hollasch.net/cgindex/math/matrix/column-vec.html
  119. Examples
  120. --------
  121. >>> alpha, beta, gamma = 0.123, -1.234, 2.345
  122. >>> origin, xaxis, yaxis, zaxis = [0, 0, 0], [1, 0, 0], [0, 1, 0], [0, 0, 1]
  123. >>> I = identity_matrix()
  124. >>> Rx = rotation_matrix(alpha, xaxis)
  125. >>> Ry = rotation_matrix(beta, yaxis)
  126. >>> Rz = rotation_matrix(gamma, zaxis)
  127. >>> R = concatenate_matrices(Rx, Ry, Rz)
  128. >>> euler = euler_from_matrix(R, 'rxyz')
  129. >>> numpy.allclose([alpha, beta, gamma], euler)
  130. True
  131. >>> Re = euler_matrix(alpha, beta, gamma, 'rxyz')
  132. >>> is_same_transform(R, Re)
  133. True
  134. >>> al, be, ga = euler_from_matrix(Re, 'rxyz')
  135. >>> is_same_transform(Re, euler_matrix(al, be, ga, 'rxyz'))
  136. True
  137. >>> qx = quaternion_about_axis(alpha, xaxis)
  138. >>> qy = quaternion_about_axis(beta, yaxis)
  139. >>> qz = quaternion_about_axis(gamma, zaxis)
  140. >>> q = quaternion_multiply(qx, qy)
  141. >>> q = quaternion_multiply(q, qz)
  142. >>> Rq = quaternion_matrix(q)
  143. >>> is_same_transform(R, Rq)
  144. True
  145. >>> S = scale_matrix(1.23, origin)
  146. >>> T = translation_matrix([1, 2, 3])
  147. >>> Z = shear_matrix(beta, xaxis, origin, zaxis)
  148. >>> R = random_rotation_matrix(numpy.random.rand(3))
  149. >>> M = concatenate_matrices(T, R, Z, S)
  150. >>> scale, shear, angles, trans, persp = decompose_matrix(M)
  151. >>> numpy.allclose(scale, 1.23)
  152. True
  153. >>> numpy.allclose(trans, [1, 2, 3])
  154. True
  155. >>> numpy.allclose(shear, [0, math.tan(beta), 0])
  156. True
  157. >>> is_same_transform(R, euler_matrix(axes='sxyz', *angles))
  158. True
  159. >>> M1 = compose_matrix(scale, shear, angles, trans, persp)
  160. >>> is_same_transform(M, M1)
  161. True
  162. >>> v0, v1 = random_vector(3), random_vector(3)
  163. >>> M = rotation_matrix(angle_between_vectors(v0, v1), vector_product(v0, v1))
  164. >>> v2 = numpy.dot(v0, M[:3,:3].T)
  165. >>> numpy.allclose(unit_vector(v1), unit_vector(v2))
  166. True
  167. """
  168. import math
  169. import numpy
  170. __version__ = '2015.01.29'
  171. __docformat__ = 'restructuredtext en'
  172. __all__ = ()
  173. def identity_matrix():
  174. """Return 4x4 identity/unit matrix.
  175. >>> I = identity_matrix()
  176. >>> numpy.allclose(I, numpy.dot(I, I))
  177. True
  178. >>> numpy.sum(I), numpy.trace(I)
  179. (4.0, 4.0)
  180. >>> numpy.allclose(I, numpy.identity(4))
  181. True
  182. """
  183. return numpy.identity(4)
  184. def translation_matrix(direction):
  185. """Return matrix to translate by direction vector.
  186. >>> v = numpy.random.random(3) - 0.5
  187. >>> numpy.allclose(v, translation_matrix(v)[:3, 3])
  188. True
  189. """
  190. M = numpy.identity(4)
  191. M[:3, 3] = direction[:3]
  192. return M
  193. def translation_from_matrix(matrix):
  194. """Return translation vector from translation matrix.
  195. >>> v0 = numpy.random.random(3) - 0.5
  196. >>> v1 = translation_from_matrix(translation_matrix(v0))
  197. >>> numpy.allclose(v0, v1)
  198. True
  199. """
  200. return numpy.array(matrix, copy=False)[:3, 3].copy()
  201. def reflection_matrix(point, normal):
  202. """Return matrix to mirror at plane defined by point and normal vector.
  203. >>> v0 = numpy.random.random(4) - 0.5
  204. >>> v0[3] = 1.
  205. >>> v1 = numpy.random.random(3) - 0.5
  206. >>> R = reflection_matrix(v0, v1)
  207. >>> numpy.allclose(2, numpy.trace(R))
  208. True
  209. >>> numpy.allclose(v0, numpy.dot(R, v0))
  210. True
  211. >>> v2 = v0.copy()
  212. >>> v2[:3] += v1
  213. >>> v3 = v0.copy()
  214. >>> v2[:3] -= v1
  215. >>> numpy.allclose(v2, numpy.dot(R, v3))
  216. True
  217. """
  218. normal = unit_vector(normal[:3])
  219. M = numpy.identity(4)
  220. M[:3, :3] -= 2.0 * numpy.outer(normal, normal)
  221. M[:3, 3] = (2.0 * numpy.dot(point[:3], normal)) * normal
  222. return M
  223. def reflection_from_matrix(matrix):
  224. """Return mirror plane point and normal vector from reflection matrix.
  225. >>> v0 = numpy.random.random(3) - 0.5
  226. >>> v1 = numpy.random.random(3) - 0.5
  227. >>> M0 = reflection_matrix(v0, v1)
  228. >>> point, normal = reflection_from_matrix(M0)
  229. >>> M1 = reflection_matrix(point, normal)
  230. >>> is_same_transform(M0, M1)
  231. True
  232. """
  233. M = numpy.array(matrix, dtype=numpy.float64, copy=False)
  234. # normal: unit eigenvector corresponding to eigenvalue -1
  235. w, V = numpy.linalg.eig(M[:3, :3])
  236. i = numpy.where(abs(numpy.real(w) + 1.0) < 1e-8)[0]
  237. if not len(i):
  238. raise ValueError("no unit eigenvector corresponding to eigenvalue -1")
  239. normal = numpy.real(V[:, i[0]]).squeeze()
  240. # point: any unit eigenvector corresponding to eigenvalue 1
  241. w, V = numpy.linalg.eig(M)
  242. i = numpy.where(abs(numpy.real(w) - 1.0) < 1e-8)[0]
  243. if not len(i):
  244. raise ValueError("no unit eigenvector corresponding to eigenvalue 1")
  245. point = numpy.real(V[:, i[-1]]).squeeze()
  246. point /= point[3]
  247. return point, normal
  248. def rotation_matrix(angle, direction, point=None):
  249. """Return matrix to rotate about axis defined by point and direction.
  250. >>> R = rotation_matrix(math.pi/2, [0, 0, 1], [1, 0, 0])
  251. >>> numpy.allclose(numpy.dot(R, [0, 0, 0, 1]), [1, -1, 0, 1])
  252. True
  253. >>> angle = (random.random() - 0.5) * (2*math.pi)
  254. >>> direc = numpy.random.random(3) - 0.5
  255. >>> point = numpy.random.random(3) - 0.5
  256. >>> R0 = rotation_matrix(angle, direc, point)
  257. >>> R1 = rotation_matrix(angle-2*math.pi, direc, point)
  258. >>> is_same_transform(R0, R1)
  259. True
  260. >>> R0 = rotation_matrix(angle, direc, point)
  261. >>> R1 = rotation_matrix(-angle, -direc, point)
  262. >>> is_same_transform(R0, R1)
  263. True
  264. >>> I = numpy.identity(4, numpy.float64)
  265. >>> numpy.allclose(I, rotation_matrix(math.pi*2, direc))
  266. True
  267. >>> numpy.allclose(2, numpy.trace(rotation_matrix(math.pi/2,
  268. ... direc, point)))
  269. True
  270. """
  271. sina = math.sin(angle)
  272. cosa = math.cos(angle)
  273. direction = unit_vector(direction[:3])
  274. # rotation matrix around unit vector
  275. R = numpy.diag([cosa, cosa, cosa])
  276. R += numpy.outer(direction, direction) * (1.0 - cosa)
  277. direction *= sina
  278. R += numpy.array([[ 0.0, -direction[2], direction[1]],
  279. [ direction[2], 0.0, -direction[0]],
  280. [-direction[1], direction[0], 0.0]])
  281. M = numpy.identity(4)
  282. M[:3, :3] = R
  283. if point is not None:
  284. # rotation not around origin
  285. point = numpy.array(point[:3], dtype=numpy.float64, copy=False)
  286. M[:3, 3] = point - numpy.dot(R, point)
  287. return M
  288. def rotation_from_matrix(matrix):
  289. """Return rotation angle and axis from rotation matrix.
  290. >>> angle = (random.random() - 0.5) * (2*math.pi)
  291. >>> direc = numpy.random.random(3) - 0.5
  292. >>> point = numpy.random.random(3) - 0.5
  293. >>> R0 = rotation_matrix(angle, direc, point)
  294. >>> angle, direc, point = rotation_from_matrix(R0)
  295. >>> R1 = rotation_matrix(angle, direc, point)
  296. >>> is_same_transform(R0, R1)
  297. True
  298. """
  299. R = numpy.array(matrix, dtype=numpy.float64, copy=False)
  300. R33 = R[:3, :3]
  301. # direction: unit eigenvector of R33 corresponding to eigenvalue of 1
  302. w, W = numpy.linalg.eig(R33.T)
  303. i = numpy.where(abs(numpy.real(w) - 1.0) < 1e-8)[0]
  304. if not len(i):
  305. raise ValueError("no unit eigenvector corresponding to eigenvalue 1")
  306. direction = numpy.real(W[:, i[-1]]).squeeze()
  307. # point: unit eigenvector of R33 corresponding to eigenvalue of 1
  308. w, Q = numpy.linalg.eig(R)
  309. i = numpy.where(abs(numpy.real(w) - 1.0) < 1e-8)[0]
  310. if not len(i):
  311. raise ValueError("no unit eigenvector corresponding to eigenvalue 1")
  312. point = numpy.real(Q[:, i[-1]]).squeeze()
  313. point /= point[3]
  314. # rotation angle depending on direction
  315. cosa = (numpy.trace(R33) - 1.0) / 2.0
  316. if abs(direction[2]) > 1e-8:
  317. sina = (R[1, 0] + (cosa-1.0)*direction[0]*direction[1]) / direction[2]
  318. elif abs(direction[1]) > 1e-8:
  319. sina = (R[0, 2] + (cosa-1.0)*direction[0]*direction[2]) / direction[1]
  320. else:
  321. sina = (R[2, 1] + (cosa-1.0)*direction[1]*direction[2]) / direction[0]
  322. angle = math.atan2(sina, cosa)
  323. return angle, direction, point
  324. def scale_matrix(factor, origin=None, direction=None):
  325. """Return matrix to scale by factor around origin in direction.
  326. Use factor -1 for point symmetry.
  327. >>> v = (numpy.random.rand(4, 5) - 0.5) * 20
  328. >>> v[3] = 1
  329. >>> S = scale_matrix(-1.234)
  330. >>> numpy.allclose(numpy.dot(S, v)[:3], -1.234*v[:3])
  331. True
  332. >>> factor = random.random() * 10 - 5
  333. >>> origin = numpy.random.random(3) - 0.5
  334. >>> direct = numpy.random.random(3) - 0.5
  335. >>> S = scale_matrix(factor, origin)
  336. >>> S = scale_matrix(factor, origin, direct)
  337. """
  338. if direction is None:
  339. # uniform scaling
  340. M = numpy.diag([factor, factor, factor, 1.0])
  341. if origin is not None:
  342. M[:3, 3] = origin[:3]
  343. M[:3, 3] *= 1.0 - factor
  344. else:
  345. # nonuniform scaling
  346. direction = unit_vector(direction[:3])
  347. factor = 1.0 - factor
  348. M = numpy.identity(4)
  349. M[:3, :3] -= factor * numpy.outer(direction, direction)
  350. if origin is not None:
  351. M[:3, 3] = (factor * numpy.dot(origin[:3], direction)) * direction
  352. return M
  353. def scale_from_matrix(matrix):
  354. """Return scaling factor, origin and direction from scaling matrix.
  355. >>> factor = random.random() * 10 - 5
  356. >>> origin = numpy.random.random(3) - 0.5
  357. >>> direct = numpy.random.random(3) - 0.5
  358. >>> S0 = scale_matrix(factor, origin)
  359. >>> factor, origin, direction = scale_from_matrix(S0)
  360. >>> S1 = scale_matrix(factor, origin, direction)
  361. >>> is_same_transform(S0, S1)
  362. True
  363. >>> S0 = scale_matrix(factor, origin, direct)
  364. >>> factor, origin, direction = scale_from_matrix(S0)
  365. >>> S1 = scale_matrix(factor, origin, direction)
  366. >>> is_same_transform(S0, S1)
  367. True
  368. """
  369. M = numpy.array(matrix, dtype=numpy.float64, copy=False)
  370. M33 = M[:3, :3]
  371. factor = numpy.trace(M33) - 2.0
  372. try:
  373. # direction: unit eigenvector corresponding to eigenvalue factor
  374. w, V = numpy.linalg.eig(M33)
  375. i = numpy.where(abs(numpy.real(w) - factor) < 1e-8)[0][0]
  376. direction = numpy.real(V[:, i]).squeeze()
  377. direction /= vector_norm(direction)
  378. except IndexError:
  379. # uniform scaling
  380. factor = (factor + 2.0) / 3.0
  381. direction = None
  382. # origin: any eigenvector corresponding to eigenvalue 1
  383. w, V = numpy.linalg.eig(M)
  384. i = numpy.where(abs(numpy.real(w) - 1.0) < 1e-8)[0]
  385. if not len(i):
  386. raise ValueError("no eigenvector corresponding to eigenvalue 1")
  387. origin = numpy.real(V[:, i[-1]]).squeeze()
  388. origin /= origin[3]
  389. return factor, origin, direction
  390. def projection_matrix(point, normal, direction=None,
  391. perspective=None, pseudo=False):
  392. """Return matrix to project onto plane defined by point and normal.
  393. Using either perspective point, projection direction, or none of both.
  394. If pseudo is True, perspective projections will preserve relative depth
  395. such that Perspective = dot(Orthogonal, PseudoPerspective).
  396. >>> P = projection_matrix([0, 0, 0], [1, 0, 0])
  397. >>> numpy.allclose(P[1:, 1:], numpy.identity(4)[1:, 1:])
  398. True
  399. >>> point = numpy.random.random(3) - 0.5
  400. >>> normal = numpy.random.random(3) - 0.5
  401. >>> direct = numpy.random.random(3) - 0.5
  402. >>> persp = numpy.random.random(3) - 0.5
  403. >>> P0 = projection_matrix(point, normal)
  404. >>> P1 = projection_matrix(point, normal, direction=direct)
  405. >>> P2 = projection_matrix(point, normal, perspective=persp)
  406. >>> P3 = projection_matrix(point, normal, perspective=persp, pseudo=True)
  407. >>> is_same_transform(P2, numpy.dot(P0, P3))
  408. True
  409. >>> P = projection_matrix([3, 0, 0], [1, 1, 0], [1, 0, 0])
  410. >>> v0 = (numpy.random.rand(4, 5) - 0.5) * 20
  411. >>> v0[3] = 1
  412. >>> v1 = numpy.dot(P, v0)
  413. >>> numpy.allclose(v1[1], v0[1])
  414. True
  415. >>> numpy.allclose(v1[0], 3-v1[1])
  416. True
  417. """
  418. M = numpy.identity(4)
  419. point = numpy.array(point[:3], dtype=numpy.float64, copy=False)
  420. normal = unit_vector(normal[:3])
  421. if perspective is not None:
  422. # perspective projection
  423. perspective = numpy.array(perspective[:3], dtype=numpy.float64,
  424. copy=False)
  425. M[0, 0] = M[1, 1] = M[2, 2] = numpy.dot(perspective-point, normal)
  426. M[:3, :3] -= numpy.outer(perspective, normal)
  427. if pseudo:
  428. # preserve relative depth
  429. M[:3, :3] -= numpy.outer(normal, normal)
  430. M[:3, 3] = numpy.dot(point, normal) * (perspective+normal)
  431. else:
  432. M[:3, 3] = numpy.dot(point, normal) * perspective
  433. M[3, :3] = -normal
  434. M[3, 3] = numpy.dot(perspective, normal)
  435. elif direction is not None:
  436. # parallel projection
  437. direction = numpy.array(direction[:3], dtype=numpy.float64, copy=False)
  438. scale = numpy.dot(direction, normal)
  439. M[:3, :3] -= numpy.outer(direction, normal) / scale
  440. M[:3, 3] = direction * (numpy.dot(point, normal) / scale)
  441. else:
  442. # orthogonal projection
  443. M[:3, :3] -= numpy.outer(normal, normal)
  444. M[:3, 3] = numpy.dot(point, normal) * normal
  445. return M
  446. def projection_from_matrix(matrix, pseudo=False):
  447. """Return projection plane and perspective point from projection matrix.
  448. Return values are same as arguments for projection_matrix function:
  449. point, normal, direction, perspective, and pseudo.
  450. >>> point = numpy.random.random(3) - 0.5
  451. >>> normal = numpy.random.random(3) - 0.5
  452. >>> direct = numpy.random.random(3) - 0.5
  453. >>> persp = numpy.random.random(3) - 0.5
  454. >>> P0 = projection_matrix(point, normal)
  455. >>> result = projection_from_matrix(P0)
  456. >>> P1 = projection_matrix(*result)
  457. >>> is_same_transform(P0, P1)
  458. True
  459. >>> P0 = projection_matrix(point, normal, direct)
  460. >>> result = projection_from_matrix(P0)
  461. >>> P1 = projection_matrix(*result)
  462. >>> is_same_transform(P0, P1)
  463. True
  464. >>> P0 = projection_matrix(point, normal, perspective=persp, pseudo=False)
  465. >>> result = projection_from_matrix(P0, pseudo=False)
  466. >>> P1 = projection_matrix(*result)
  467. >>> is_same_transform(P0, P1)
  468. True
  469. >>> P0 = projection_matrix(point, normal, perspective=persp, pseudo=True)
  470. >>> result = projection_from_matrix(P0, pseudo=True)
  471. >>> P1 = projection_matrix(*result)
  472. >>> is_same_transform(P0, P1)
  473. True
  474. """
  475. M = numpy.array(matrix, dtype=numpy.float64, copy=False)
  476. M33 = M[:3, :3]
  477. w, V = numpy.linalg.eig(M)
  478. i = numpy.where(abs(numpy.real(w) - 1.0) < 1e-8)[0]
  479. if not pseudo and len(i):
  480. # point: any eigenvector corresponding to eigenvalue 1
  481. point = numpy.real(V[:, i[-1]]).squeeze()
  482. point /= point[3]
  483. # direction: unit eigenvector corresponding to eigenvalue 0
  484. w, V = numpy.linalg.eig(M33)
  485. i = numpy.where(abs(numpy.real(w)) < 1e-8)[0]
  486. if not len(i):
  487. raise ValueError("no eigenvector corresponding to eigenvalue 0")
  488. direction = numpy.real(V[:, i[0]]).squeeze()
  489. direction /= vector_norm(direction)
  490. # normal: unit eigenvector of M33.T corresponding to eigenvalue 0
  491. w, V = numpy.linalg.eig(M33.T)
  492. i = numpy.where(abs(numpy.real(w)) < 1e-8)[0]
  493. if len(i):
  494. # parallel projection
  495. normal = numpy.real(V[:, i[0]]).squeeze()
  496. normal /= vector_norm(normal)
  497. return point, normal, direction, None, False
  498. else:
  499. # orthogonal projection, where normal equals direction vector
  500. return point, direction, None, None, False
  501. else:
  502. # perspective projection
  503. i = numpy.where(abs(numpy.real(w)) > 1e-8)[0]
  504. if not len(i):
  505. raise ValueError(
  506. "no eigenvector not corresponding to eigenvalue 0")
  507. point = numpy.real(V[:, i[-1]]).squeeze()
  508. point /= point[3]
  509. normal = - M[3, :3]
  510. perspective = M[:3, 3] / numpy.dot(point[:3], normal)
  511. if pseudo:
  512. perspective -= normal
  513. return point, normal, None, perspective, pseudo
  514. def clip_matrix(left, right, bottom, top, near, far, perspective=False):
  515. """Return matrix to obtain normalized device coordinates from frustum.
  516. The frustum bounds are axis-aligned along x (left, right),
  517. y (bottom, top) and z (near, far).
  518. Normalized device coordinates are in range [-1, 1] if coordinates are
  519. inside the frustum.
  520. If perspective is True the frustum is a truncated pyramid with the
  521. perspective point at origin and direction along z axis, otherwise an
  522. orthographic canonical view volume (a box).
  523. Homogeneous coordinates transformed by the perspective clip matrix
  524. need to be dehomogenized (divided by w coordinate).
  525. >>> frustum = numpy.random.rand(6)
  526. >>> frustum[1] += frustum[0]
  527. >>> frustum[3] += frustum[2]
  528. >>> frustum[5] += frustum[4]
  529. >>> M = clip_matrix(perspective=False, *frustum)
  530. >>> numpy.dot(M, [frustum[0], frustum[2], frustum[4], 1])
  531. array([-1., -1., -1., 1.])
  532. >>> numpy.dot(M, [frustum[1], frustum[3], frustum[5], 1])
  533. array([ 1., 1., 1., 1.])
  534. >>> M = clip_matrix(perspective=True, *frustum)
  535. >>> v = numpy.dot(M, [frustum[0], frustum[2], frustum[4], 1])
  536. >>> v / v[3]
  537. array([-1., -1., -1., 1.])
  538. >>> v = numpy.dot(M, [frustum[1], frustum[3], frustum[4], 1])
  539. >>> v / v[3]
  540. array([ 1., 1., -1., 1.])
  541. """
  542. if left >= right or bottom >= top or near >= far:
  543. raise ValueError("invalid frustum")
  544. if perspective:
  545. if near <= _EPS:
  546. raise ValueError("invalid frustum: near <= 0")
  547. t = 2.0 * near
  548. M = [[t/(left-right), 0.0, (right+left)/(right-left), 0.0],
  549. [0.0, t/(bottom-top), (top+bottom)/(top-bottom), 0.0],
  550. [0.0, 0.0, (far+near)/(near-far), t*far/(far-near)],
  551. [0.0, 0.0, -1.0, 0.0]]
  552. else:
  553. M = [[2.0/(right-left), 0.0, 0.0, (right+left)/(left-right)],
  554. [0.0, 2.0/(top-bottom), 0.0, (top+bottom)/(bottom-top)],
  555. [0.0, 0.0, 2.0/(far-near), (far+near)/(near-far)],
  556. [0.0, 0.0, 0.0, 1.0]]
  557. return numpy.array(M)
  558. def shear_matrix(angle, direction, point, normal):
  559. """Return matrix to shear by angle along direction vector on shear plane.
  560. The shear plane is defined by a point and normal vector. The direction
  561. vector must be orthogonal to the plane's normal vector.
  562. A point P is transformed by the shear matrix into P" such that
  563. the vector P-P" is parallel to the direction vector and its extent is
  564. given by the angle of P-P'-P", where P' is the orthogonal projection
  565. of P onto the shear plane.
  566. >>> angle = (random.random() - 0.5) * 4*math.pi
  567. >>> direct = numpy.random.random(3) - 0.5
  568. >>> point = numpy.random.random(3) - 0.5
  569. >>> normal = numpy.cross(direct, numpy.random.random(3))
  570. >>> S = shear_matrix(angle, direct, point, normal)
  571. >>> numpy.allclose(1, numpy.linalg.det(S))
  572. True
  573. """
  574. normal = unit_vector(normal[:3])
  575. direction = unit_vector(direction[:3])
  576. if abs(numpy.dot(normal, direction)) > 1e-6:
  577. raise ValueError("direction and normal vectors are not orthogonal")
  578. angle = math.tan(angle)
  579. M = numpy.identity(4)
  580. M[:3, :3] += angle * numpy.outer(direction, normal)
  581. M[:3, 3] = -angle * numpy.dot(point[:3], normal) * direction
  582. return M
  583. def shear_from_matrix(matrix):
  584. """Return shear angle, direction and plane from shear matrix.
  585. >>> angle = (random.random() - 0.5) * 4*math.pi
  586. >>> direct = numpy.random.random(3) - 0.5
  587. >>> point = numpy.random.random(3) - 0.5
  588. >>> normal = numpy.cross(direct, numpy.random.random(3))
  589. >>> S0 = shear_matrix(angle, direct, point, normal)
  590. >>> angle, direct, point, normal = shear_from_matrix(S0)
  591. >>> S1 = shear_matrix(angle, direct, point, normal)
  592. >>> is_same_transform(S0, S1)
  593. True
  594. """
  595. M = numpy.array(matrix, dtype=numpy.float64, copy=False)
  596. M33 = M[:3, :3]
  597. # normal: cross independent eigenvectors corresponding to the eigenvalue 1
  598. w, V = numpy.linalg.eig(M33)
  599. i = numpy.where(abs(numpy.real(w) - 1.0) < 1e-4)[0]
  600. if len(i) < 2:
  601. raise ValueError("no two linear independent eigenvectors found %s" % w)
  602. V = numpy.real(V[:, i]).squeeze().T
  603. lenorm = -1.0
  604. for i0, i1 in ((0, 1), (0, 2), (1, 2)):
  605. n = numpy.cross(V[i0], V[i1])
  606. w = vector_norm(n)
  607. if w > lenorm:
  608. lenorm = w
  609. normal = n
  610. normal /= lenorm
  611. # direction and angle
  612. direction = numpy.dot(M33 - numpy.identity(3), normal)
  613. angle = vector_norm(direction)
  614. direction /= angle
  615. angle = math.atan(angle)
  616. # point: eigenvector corresponding to eigenvalue 1
  617. w, V = numpy.linalg.eig(M)
  618. i = numpy.where(abs(numpy.real(w) - 1.0) < 1e-8)[0]
  619. if not len(i):
  620. raise ValueError("no eigenvector corresponding to eigenvalue 1")
  621. point = numpy.real(V[:, i[-1]]).squeeze()
  622. point /= point[3]
  623. return angle, direction, point, normal
  624. def decompose_matrix(matrix):
  625. """Return sequence of transformations from transformation matrix.
  626. matrix : array_like
  627. Non-degenerative homogeneous transformation matrix
  628. Return tuple of:
  629. scale : vector of 3 scaling factors
  630. shear : list of shear factors for x-y, x-z, y-z axes
  631. angles : list of Euler angles about static x, y, z axes
  632. translate : translation vector along x, y, z axes
  633. perspective : perspective partition of matrix
  634. Raise ValueError if matrix is of wrong type or degenerative.
  635. >>> T0 = translation_matrix([1, 2, 3])
  636. >>> scale, shear, angles, trans, persp = decompose_matrix(T0)
  637. >>> T1 = translation_matrix(trans)
  638. >>> numpy.allclose(T0, T1)
  639. True
  640. >>> S = scale_matrix(0.123)
  641. >>> scale, shear, angles, trans, persp = decompose_matrix(S)
  642. >>> scale[0]
  643. 0.123
  644. >>> R0 = euler_matrix(1, 2, 3)
  645. >>> scale, shear, angles, trans, persp = decompose_matrix(R0)
  646. >>> R1 = euler_matrix(*angles)
  647. >>> numpy.allclose(R0, R1)
  648. True
  649. """
  650. M = numpy.array(matrix, dtype=numpy.float64, copy=True).T
  651. if abs(M[3, 3]) < _EPS:
  652. raise ValueError("M[3, 3] is zero")
  653. M /= M[3, 3]
  654. P = M.copy()
  655. P[:, 3] = 0.0, 0.0, 0.0, 1.0
  656. if not numpy.linalg.det(P):
  657. raise ValueError("matrix is singular")
  658. scale = numpy.zeros((3, ))
  659. shear = [0.0, 0.0, 0.0]
  660. angles = [0.0, 0.0, 0.0]
  661. if any(abs(M[:3, 3]) > _EPS):
  662. perspective = numpy.dot(M[:, 3], numpy.linalg.inv(P.T))
  663. M[:, 3] = 0.0, 0.0, 0.0, 1.0
  664. else:
  665. perspective = numpy.array([0.0, 0.0, 0.0, 1.0])
  666. translate = M[3, :3].copy()
  667. M[3, :3] = 0.0
  668. row = M[:3, :3].copy()
  669. scale[0] = vector_norm(row[0])
  670. row[0] /= scale[0]
  671. shear[0] = numpy.dot(row[0], row[1])
  672. row[1] -= row[0] * shear[0]
  673. scale[1] = vector_norm(row[1])
  674. row[1] /= scale[1]
  675. shear[0] /= scale[1]
  676. shear[1] = numpy.dot(row[0], row[2])
  677. row[2] -= row[0] * shear[1]
  678. shear[2] = numpy.dot(row[1], row[2])
  679. row[2] -= row[1] * shear[2]
  680. scale[2] = vector_norm(row[2])
  681. row[2] /= scale[2]
  682. shear[1:] /= scale[2]
  683. if numpy.dot(row[0], numpy.cross(row[1], row[2])) < 0:
  684. numpy.negative(scale, scale)
  685. numpy.negative(row, row)
  686. angles[1] = math.asin(-row[0, 2])
  687. if math.cos(angles[1]):
  688. angles[0] = math.atan2(row[1, 2], row[2, 2])
  689. angles[2] = math.atan2(row[0, 1], row[0, 0])
  690. else:
  691. #angles[0] = math.atan2(row[1, 0], row[1, 1])
  692. angles[0] = math.atan2(-row[2, 1], row[1, 1])
  693. angles[2] = 0.0
  694. return scale, shear, angles, translate, perspective
  695. def compose_matrix(scale=None, shear=None, angles=None, translate=None,
  696. perspective=None):
  697. """Return transformation matrix from sequence of transformations.
  698. This is the inverse of the decompose_matrix function.
  699. Sequence of transformations:
  700. scale : vector of 3 scaling factors
  701. shear : list of shear factors for x-y, x-z, y-z axes
  702. angles : list of Euler angles about static x, y, z axes
  703. translate : translation vector along x, y, z axes
  704. perspective : perspective partition of matrix
  705. >>> scale = numpy.random.random(3) - 0.5
  706. >>> shear = numpy.random.random(3) - 0.5
  707. >>> angles = (numpy.random.random(3) - 0.5) * (2*math.pi)
  708. >>> trans = numpy.random.random(3) - 0.5
  709. >>> persp = numpy.random.random(4) - 0.5
  710. >>> M0 = compose_matrix(scale, shear, angles, trans, persp)
  711. >>> result = decompose_matrix(M0)
  712. >>> M1 = compose_matrix(*result)
  713. >>> is_same_transform(M0, M1)
  714. True
  715. """
  716. M = numpy.identity(4)
  717. if perspective is not None:
  718. P = numpy.identity(4)
  719. P[3, :] = perspective[:4]
  720. M = numpy.dot(M, P)
  721. if translate is not None:
  722. T = numpy.identity(4)
  723. T[:3, 3] = translate[:3]
  724. M = numpy.dot(M, T)
  725. if angles is not None:
  726. R = euler_matrix(angles[0], angles[1], angles[2], 'sxyz')
  727. M = numpy.dot(M, R)
  728. if shear is not None:
  729. Z = numpy.identity(4)
  730. Z[1, 2] = shear[2]
  731. Z[0, 2] = shear[1]
  732. Z[0, 1] = shear[0]
  733. M = numpy.dot(M, Z)
  734. if scale is not None:
  735. S = numpy.identity(4)
  736. S[0, 0] = scale[0]
  737. S[1, 1] = scale[1]
  738. S[2, 2] = scale[2]
  739. M = numpy.dot(M, S)
  740. M /= M[3, 3]
  741. return M
  742. def orthogonalization_matrix(lengths, angles):
  743. """Return orthogonalization matrix for crystallographic cell coordinates.
  744. Angles are expected in degrees.
  745. The de-orthogonalization matrix is the inverse.
  746. >>> O = orthogonalization_matrix([10, 10, 10], [90, 90, 90])
  747. >>> numpy.allclose(O[:3, :3], numpy.identity(3, float) * 10)
  748. True
  749. >>> O = orthogonalization_matrix([9.8, 12.0, 15.5], [87.2, 80.7, 69.7])
  750. >>> numpy.allclose(numpy.sum(O), 43.063229)
  751. True
  752. """
  753. a, b, c = lengths
  754. angles = numpy.radians(angles)
  755. sina, sinb, _ = numpy.sin(angles)
  756. cosa, cosb, cosg = numpy.cos(angles)
  757. co = (cosa * cosb - cosg) / (sina * sinb)
  758. return numpy.array([
  759. [ a*sinb*math.sqrt(1.0-co*co), 0.0, 0.0, 0.0],
  760. [-a*sinb*co, b*sina, 0.0, 0.0],
  761. [ a*cosb, b*cosa, c, 0.0],
  762. [ 0.0, 0.0, 0.0, 1.0]])
  763. def affine_matrix_from_points(v0, v1, shear=True, scale=True, usesvd=True):
  764. """Return affine transform matrix to register two point sets.
  765. v0 and v1 are shape (ndims, \*) arrays of at least ndims non-homogeneous
  766. coordinates, where ndims is the dimensionality of the coordinate space.
  767. If shear is False, a similarity transformation matrix is returned.
  768. If also scale is False, a rigid/Euclidean transformation matrix
  769. is returned.
  770. By default the algorithm by Hartley and Zissermann [15] is used.
  771. If usesvd is True, similarity and Euclidean transformation matrices
  772. are calculated by minimizing the weighted sum of squared deviations
  773. (RMSD) according to the algorithm by Kabsch [8].
  774. Otherwise, and if ndims is 3, the quaternion based algorithm by Horn [9]
  775. is used, which is slower when using this Python implementation.
  776. The returned matrix performs rotation, translation and uniform scaling
  777. (if specified).
  778. >>> v0 = [[0, 1031, 1031, 0], [0, 0, 1600, 1600]]
  779. >>> v1 = [[675, 826, 826, 677], [55, 52, 281, 277]]
  780. >>> affine_matrix_from_points(v0, v1)
  781. array([[ 0.14549, 0.00062, 675.50008],
  782. [ 0.00048, 0.14094, 53.24971],
  783. [ 0. , 0. , 1. ]])
  784. >>> T = translation_matrix(numpy.random.random(3)-0.5)
  785. >>> R = random_rotation_matrix(numpy.random.random(3))
  786. >>> S = scale_matrix(random.random())
  787. >>> M = concatenate_matrices(T, R, S)
  788. >>> v0 = (numpy.random.rand(4, 100) - 0.5) * 20
  789. >>> v0[3] = 1
  790. >>> v1 = numpy.dot(M, v0)
  791. >>> v0[:3] += numpy.random.normal(0, 1e-8, 300).reshape(3, -1)
  792. >>> M = affine_matrix_from_points(v0[:3], v1[:3])
  793. >>> numpy.allclose(v1, numpy.dot(M, v0))
  794. True
  795. More examples in superimposition_matrix()
  796. """
  797. v0 = numpy.array(v0, dtype=numpy.float64, copy=True)
  798. v1 = numpy.array(v1, dtype=numpy.float64, copy=True)
  799. ndims = v0.shape[0]
  800. if ndims < 2 or v0.shape[1] < ndims or v0.shape != v1.shape:
  801. raise ValueError("input arrays are of wrong shape or type")
  802. # move centroids to origin
  803. t0 = -numpy.mean(v0, axis=1)
  804. M0 = numpy.identity(ndims+1)
  805. M0[:ndims, ndims] = t0
  806. v0 += t0.reshape(ndims, 1)
  807. t1 = -numpy.mean(v1, axis=1)
  808. M1 = numpy.identity(ndims+1)
  809. M1[:ndims, ndims] = t1
  810. v1 += t1.reshape(ndims, 1)
  811. if shear:
  812. # Affine transformation
  813. A = numpy.concatenate((v0, v1), axis=0)
  814. u, s, vh = numpy.linalg.svd(A.T)
  815. vh = vh[:ndims].T
  816. B = vh[:ndims]
  817. C = vh[ndims:2*ndims]
  818. t = numpy.dot(C, numpy.linalg.pinv(B))
  819. t = numpy.concatenate((t, numpy.zeros((ndims, 1))), axis=1)
  820. M = numpy.vstack((t, ((0.0,)*ndims) + (1.0,)))
  821. elif usesvd or ndims != 3:
  822. # Rigid transformation via SVD of covariance matrix
  823. u, s, vh = numpy.linalg.svd(numpy.dot(v1, v0.T))
  824. # rotation matrix from SVD orthonormal bases
  825. R = numpy.dot(u, vh)
  826. if numpy.linalg.det(R) < 0.0:
  827. # R does not constitute right handed system
  828. R -= numpy.outer(u[:, ndims-1], vh[ndims-1, :]*2.0)
  829. s[-1] *= -1.0
  830. # homogeneous transformation matrix
  831. M = numpy.identity(ndims+1)
  832. M[:ndims, :ndims] = R
  833. else:
  834. # Rigid transformation matrix via quaternion
  835. # compute symmetric matrix N
  836. xx, yy, zz = numpy.sum(v0 * v1, axis=1)
  837. xy, yz, zx = numpy.sum(v0 * numpy.roll(v1, -1, axis=0), axis=1)
  838. xz, yx, zy = numpy.sum(v0 * numpy.roll(v1, -2, axis=0), axis=1)
  839. N = [[xx+yy+zz, 0.0, 0.0, 0.0],
  840. [yz-zy, xx-yy-zz, 0.0, 0.0],
  841. [zx-xz, xy+yx, yy-xx-zz, 0.0],
  842. [xy-yx, zx+xz, yz+zy, zz-xx-yy]]
  843. # quaternion: eigenvector corresponding to most positive eigenvalue
  844. w, V = numpy.linalg.eigh(N)
  845. q = V[:, numpy.argmax(w)]
  846. q /= vector_norm(q) # unit quaternion
  847. # homogeneous transformation matrix
  848. M = quaternion_matrix(q)
  849. if scale and not shear:
  850. # Affine transformation; scale is ratio of RMS deviations from centroid
  851. v0 *= v0
  852. v1 *= v1
  853. M[:ndims, :ndims] *= math.sqrt(numpy.sum(v1) / numpy.sum(v0))
  854. # move centroids back
  855. M = numpy.dot(numpy.linalg.inv(M1), numpy.dot(M, M0))
  856. M /= M[ndims, ndims]
  857. return M
  858. def superimposition_matrix(v0, v1, scale=False, usesvd=True):
  859. """Return matrix to transform given 3D point set into second point set.
  860. v0 and v1 are shape (3, \*) or (4, \*) arrays of at least 3 points.
  861. The parameters scale and usesvd are explained in the more general
  862. affine_matrix_from_points function.
  863. The returned matrix is a similarity or Euclidean transformation matrix.
  864. This function has a fast C implementation in transformations.c.
  865. >>> v0 = numpy.random.rand(3, 10)
  866. >>> M = superimposition_matrix(v0, v0)
  867. >>> numpy.allclose(M, numpy.identity(4))
  868. True
  869. >>> R = random_rotation_matrix(numpy.random.random(3))
  870. >>> v0 = [[1,0,0], [0,1,0], [0,0,1], [1,1,1]]
  871. >>> v1 = numpy.dot(R, v0)
  872. >>> M = superimposition_matrix(v0, v1)
  873. >>> numpy.allclose(v1, numpy.dot(M, v0))
  874. True
  875. >>> v0 = (numpy.random.rand(4, 100) - 0.5) * 20
  876. >>> v0[3] = 1
  877. >>> v1 = numpy.dot(R, v0)
  878. >>> M = superimposition_matrix(v0, v1)
  879. >>> numpy.allclose(v1, numpy.dot(M, v0))
  880. True
  881. >>> S = scale_matrix(random.random())
  882. >>> T = translation_matrix(numpy.random.random(3)-0.5)
  883. >>> M = concatenate_matrices(T, R, S)
  884. >>> v1 = numpy.dot(M, v0)
  885. >>> v0[:3] += numpy.random.normal(0, 1e-9, 300).reshape(3, -1)
  886. >>> M = superimposition_matrix(v0, v1, scale=True)
  887. >>> numpy.allclose(v1, numpy.dot(M, v0))
  888. True
  889. >>> M = superimposition_matrix(v0, v1, scale=True, usesvd=False)
  890. >>> numpy.allclose(v1, numpy.dot(M, v0))
  891. True
  892. >>> v = numpy.empty((4, 100, 3))
  893. >>> v[:, :, 0] = v0
  894. >>> M = superimposition_matrix(v0, v1, scale=True, usesvd=False)
  895. >>> numpy.allclose(v1, numpy.dot(M, v[:, :, 0]))
  896. True
  897. """
  898. v0 = numpy.array(v0, dtype=numpy.float64, copy=False)[:3]
  899. v1 = numpy.array(v1, dtype=numpy.float64, copy=False)[:3]
  900. return affine_matrix_from_points(v0, v1, shear=False,
  901. scale=scale, usesvd=usesvd)
  902. def euler_matrix(ai, aj, ak, axes='sxyz'):
  903. """Return homogeneous rotation matrix from Euler angles and axis sequence.
  904. ai, aj, ak : Euler's roll, pitch and yaw angles
  905. axes : One of 24 axis sequences as string or encoded tuple
  906. >>> R = euler_matrix(1, 2, 3, 'syxz')
  907. >>> numpy.allclose(numpy.sum(R[0]), -1.34786452)
  908. True
  909. >>> R = euler_matrix(1, 2, 3, (0, 1, 0, 1))
  910. >>> numpy.allclose(numpy.sum(R[0]), -0.383436184)
  911. True
  912. >>> ai, aj, ak = (4*math.pi) * (numpy.random.random(3) - 0.5)
  913. >>> for axes in _AXES2TUPLE.keys():
  914. ... R = euler_matrix(ai, aj, ak, axes)
  915. >>> for axes in _TUPLE2AXES.keys():
  916. ... R = euler_matrix(ai, aj, ak, axes)
  917. """
  918. try:
  919. firstaxis, parity, repetition, frame = _AXES2TUPLE[axes]
  920. except (AttributeError, KeyError):
  921. _TUPLE2AXES[axes] # validation
  922. firstaxis, parity, repetition, frame = axes
  923. i = firstaxis
  924. j = _NEXT_AXIS[i+parity]
  925. k = _NEXT_AXIS[i-parity+1]
  926. if frame:
  927. ai, ak = ak, ai
  928. if parity:
  929. ai, aj, ak = -ai, -aj, -ak
  930. si, sj, sk = math.sin(ai), math.sin(aj), math.sin(ak)
  931. ci, cj, ck = math.cos(ai), math.cos(aj), math.cos(ak)
  932. cc, cs = ci*ck, ci*sk
  933. sc, ss = si*ck, si*sk
  934. M = numpy.identity(4)
  935. if repetition:
  936. M[i, i] = cj
  937. M[i, j] = sj*si
  938. M[i, k] = sj*ci
  939. M[j, i] = sj*sk
  940. M[j, j] = -cj*ss+cc
  941. M[j, k] = -cj*cs-sc
  942. M[k, i] = -sj*ck
  943. M[k, j] = cj*sc+cs
  944. M[k, k] = cj*cc-ss
  945. else:
  946. M[i, i] = cj*ck
  947. M[i, j] = sj*sc-cs
  948. M[i, k] = sj*cc+ss
  949. M[j, i] = cj*sk
  950. M[j, j] = sj*ss+cc
  951. M[j, k] = sj*cs-sc
  952. M[k, i] = -sj
  953. M[k, j] = cj*si
  954. M[k, k] = cj*ci
  955. return M
  956. def euler_from_matrix(matrix, axes='sxyz'):
  957. """Return Euler angles from rotation matrix for specified axis sequence.
  958. axes : One of 24 axis sequences as string or encoded tuple
  959. Note that many Euler angle triplets can describe one matrix.
  960. >>> R0 = euler_matrix(1, 2, 3, 'syxz')
  961. >>> al, be, ga = euler_from_matrix(R0, 'syxz')
  962. >>> R1 = euler_matrix(al, be, ga, 'syxz')
  963. >>> numpy.allclose(R0, R1)
  964. True
  965. >>> angles = (4*math.pi) * (numpy.random.random(3) - 0.5)
  966. >>> for axes in _AXES2TUPLE.keys():
  967. ... R0 = euler_matrix(axes=axes, *angles)
  968. ... R1 = euler_matrix(axes=axes, *euler_from_matrix(R0, axes))
  969. ... if not numpy.allclose(R0, R1): print(axes, "failed")
  970. """
  971. try:
  972. firstaxis, parity, repetition, frame = _AXES2TUPLE[axes.lower()]
  973. except (AttributeError, KeyError):
  974. _TUPLE2AXES[axes] # validation
  975. firstaxis, parity, repetition, frame = axes
  976. i = firstaxis
  977. j = _NEXT_AXIS[i+parity]
  978. k = _NEXT_AXIS[i-parity+1]
  979. M = numpy.array(matrix, dtype=numpy.float64, copy=False)[:3, :3]
  980. if repetition:
  981. sy = math.sqrt(M[i, j]*M[i, j] + M[i, k]*M[i, k])
  982. if sy > _EPS:
  983. ax = math.atan2( M[i, j], M[i, k])
  984. ay = math.atan2( sy, M[i, i])
  985. az = math.atan2( M[j, i], -M[k, i])
  986. else:
  987. ax = math.atan2(-M[j, k], M[j, j])
  988. ay = math.atan2( sy, M[i, i])
  989. az = 0.0
  990. else:
  991. cy = math.sqrt(M[i, i]*M[i, i] + M[j, i]*M[j, i])
  992. if cy > _EPS:
  993. ax = math.atan2( M[k, j], M[k, k])
  994. ay = math.atan2(-M[k, i], cy)
  995. az = math.atan2( M[j, i], M[i, i])
  996. else:
  997. ax = math.atan2(-M[j, k], M[j, j])
  998. ay = math.atan2(-M[k, i], cy)
  999. az = 0.0
  1000. if parity:
  1001. ax, ay, az = -ax, -ay, -az
  1002. if frame:
  1003. ax, az = az, ax
  1004. return ax, ay, az
  1005. def euler_from_quaternion(quaternion, axes='sxyz'):
  1006. """Return Euler angles from quaternion for specified axis sequence.
  1007. >>> angles = euler_from_quaternion([0.99810947, 0.06146124, 0, 0])
  1008. >>> numpy.allclose(angles, [0.123, 0, 0])
  1009. True
  1010. """
  1011. return euler_from_matrix(quaternion_matrix(quaternion), axes)
  1012. def quaternion_from_euler(ai, aj, ak, axes='sxyz'):
  1013. """Return quaternion from Euler angles and axis sequence.
  1014. ai, aj, ak : Euler's roll, pitch and yaw angles
  1015. axes : One of 24 axis sequences as string or encoded tuple
  1016. >>> q = quaternion_from_euler(1, 2, 3, 'ryxz')
  1017. >>> numpy.allclose(q, [0.435953, 0.310622, -0.718287, 0.444435])
  1018. True
  1019. """
  1020. try:
  1021. firstaxis, parity, repetition, frame = _AXES2TUPLE[axes.lower()]
  1022. except (AttributeError, KeyError):
  1023. _TUPLE2AXES[axes] # validation
  1024. firstaxis, parity, repetition, frame = axes
  1025. i = firstaxis + 1
  1026. j = _NEXT_AXIS[i+parity-1] + 1
  1027. k = _NEXT_AXIS[i-parity] + 1
  1028. if frame:
  1029. ai, ak = ak, ai
  1030. if parity:
  1031. aj = -aj
  1032. ai /= 2.0
  1033. aj /= 2.0
  1034. ak /= 2.0
  1035. ci = math.cos(ai)
  1036. si = math.sin(ai)
  1037. cj = math.cos(aj)
  1038. sj = math.sin(aj)
  1039. ck = math.cos(ak)
  1040. sk = math.sin(ak)
  1041. cc = ci*ck
  1042. cs = ci*sk
  1043. sc = si*ck
  1044. ss = si*sk
  1045. q = numpy.empty((4, ))
  1046. if repetition:
  1047. q[0] = cj*(cc - ss)
  1048. q[i] = cj*(cs + sc)
  1049. q[j] = sj*(cc + ss)
  1050. q[k] = sj*(cs - sc)
  1051. else:
  1052. q[0] = cj*cc + sj*ss
  1053. q[i] = cj*sc - sj*cs
  1054. q[j] = cj*ss + sj*cc
  1055. q[k] = cj*cs - sj*sc
  1056. if parity:
  1057. q[j] *= -1.0
  1058. return q
  1059. def quaternion_about_axis(angle, axis):
  1060. """Return quaternion for rotation about axis.
  1061. >>> q = quaternion_about_axis(0.123, [1, 0, 0])
  1062. >>> numpy.allclose(q, [0.99810947, 0.06146124, 0, 0])
  1063. True
  1064. """
  1065. q = numpy.array([0.0, axis[0], axis[1], axis[2]])
  1066. qlen = vector_norm(q)
  1067. if qlen > _EPS:
  1068. q *= math.sin(angle/2.0) / qlen
  1069. q[0] = math.cos(angle/2.0)
  1070. return q
  1071. def quaternion_matrix(quaternion):
  1072. """Return homogeneous rotation matrix from quaternion.
  1073. >>> M = quaternion_matrix([0.99810947, 0.06146124, 0, 0])
  1074. >>> numpy.allclose(M, rotation_matrix(0.123, [1, 0, 0]))
  1075. True
  1076. >>> M = quaternion_matrix([1, 0, 0, 0])
  1077. >>> numpy.allclose(M, numpy.identity(4))
  1078. True
  1079. >>> M = quaternion_matrix([0, 1, 0, 0])
  1080. >>> numpy.allclose(M, numpy.diag([1, -1, -1, 1]))
  1081. True
  1082. """
  1083. q = numpy.array(quaternion, dtype=numpy.float64, copy=True)
  1084. n = numpy.dot(q, q)
  1085. if n < _EPS:
  1086. return numpy.identity(4)
  1087. q *= math.sqrt(2.0 / n)
  1088. q = numpy.outer(q, q)
  1089. return numpy.array([
  1090. [1.0-q[2, 2]-q[3, 3], q[1, 2]-q[3, 0], q[1, 3]+q[2, 0], 0.0],
  1091. [ q[1, 2]+q[3, 0], 1.0-q[1, 1]-q[3, 3], q[2, 3]-q[1, 0], 0.0],
  1092. [ q[1, 3]-q[2, 0], q[2, 3]+q[1, 0], 1.0-q[1, 1]-q[2, 2], 0.0],
  1093. [ 0.0, 0.0, 0.0, 1.0]])
  1094. def quaternion_from_matrix(matrix, isprecise=False):
  1095. """Return quaternion from rotation matrix.
  1096. If isprecise is True, the input matrix is assumed to be a precise rotation
  1097. matrix and a faster algorithm is used.
  1098. >>> q = quaternion_from_matrix(numpy.identity(4), True)
  1099. >>> numpy.allclose(q, [1, 0, 0, 0])
  1100. True
  1101. >>> q = quaternion_from_matrix(numpy.diag([1, -1, -1, 1]))
  1102. >>> numpy.allclose(q, [0, 1, 0, 0]) or numpy.allclose(q, [0, -1, 0, 0])
  1103. True
  1104. >>> R = rotation_matrix(0.123, (1, 2, 3))
  1105. >>> q = quaternion_from_matrix(R, True)
  1106. >>> numpy.allclose(q, [0.9981095, 0.0164262, 0.0328524, 0.0492786])
  1107. True
  1108. >>> R = [[-0.545, 0.797, 0.260, 0], [0.733, 0.603, -0.313, 0],
  1109. ... [-0.407, 0.021, -0.913, 0], [0, 0, 0, 1]]
  1110. >>> q = quaternion_from_matrix(R)
  1111. >>> numpy.allclose(q, [0.19069, 0.43736, 0.87485, -0.083611])
  1112. True
  1113. >>> R = [[0.395, 0.362, 0.843, 0], [-0.626, 0.796, -0.056, 0],
  1114. ... [-0.677, -0.498, 0.529, 0], [0, 0, 0, 1]]
  1115. >>> q = quaternion_from_matrix(R)
  1116. >>> numpy.allclose(q, [0.82336615, -0.13610694, 0.46344705, -0.29792603])
  1117. True
  1118. >>> R = random_rotation_matrix()
  1119. >>> q = quaternion_from_matrix(R)
  1120. >>> is_same_transform(R, quaternion_matrix(q))
  1121. True
  1122. """
  1123. M = numpy.array(matrix, dtype=numpy.float64, copy=False)[:4, :4]
  1124. if isprecise:
  1125. q = numpy.empty((4, ))
  1126. t = numpy.trace(M)
  1127. if t > M[3, 3]:
  1128. q[0] = t
  1129. q[3] = M[1, 0] - M[0, 1]
  1130. q[2] = M[0, 2] - M[2, 0]
  1131. q[1] = M[2, 1] - M[1, 2]
  1132. else:
  1133. i, j, k = 1, 2, 3
  1134. if M[1, 1] > M[0, 0]:
  1135. i, j, k = 2, 3, 1
  1136. if M[2, 2] > M[i, i]:
  1137. i, j, k = 3, 1, 2
  1138. t = M[i, i] - (M[j, j] + M[k, k]) + M[3, 3]
  1139. q[i] = t
  1140. q[j] = M[i, j] + M[j, i]
  1141. q[k] = M[k, i] + M[i, k]
  1142. q[3] = M[k, j] - M[j, k]
  1143. q *= 0.5 / math.sqrt(t * M[3, 3])
  1144. else:
  1145. m00 = M[0, 0]
  1146. m01 = M[0, 1]
  1147. m02 = M[0, 2]
  1148. m10 = M[1, 0]
  1149. m11 = M[1, 1]
  1150. m12 = M[1, 2]
  1151. m20 = M[2, 0]
  1152. m21 = M[2, 1]
  1153. m22 = M[2, 2]
  1154. # symmetric matrix K
  1155. K = numpy.array([[m00-m11-m22, 0.0, 0.0, 0.0],
  1156. [m01+m10, m11-m00-m22, 0.0, 0.0],
  1157. [m02+m20, m12+m21, m22-m00-m11, 0.0],
  1158. [m21-m12, m02-m20, m10-m01, m00+m11+m22]])
  1159. K /= 3.0
  1160. # quaternion is eigenvector of K that corresponds to largest eigenvalue
  1161. w, V = numpy.linalg.eigh(K)
  1162. q = V[[3, 0, 1, 2], numpy.argmax(w)]
  1163. if q[0] < 0.0:
  1164. numpy.negative(q, q)
  1165. return q
  1166. def quaternion_multiply(quaternion1, quaternion0):
  1167. """Return multiplication of two quaternions.
  1168. >>> q = quaternion_multiply([4, 1, -2, 3], [8, -5, 6, 7])
  1169. >>> numpy.allclose(q, [28, -44, -14, 48])
  1170. True
  1171. """
  1172. w0, x0, y0, z0 = quaternion0
  1173. w1, x1, y1, z1 = quaternion1
  1174. return numpy.array([-x1*x0 - y1*y0 - z1*z0 + w1*w0,
  1175. x1*w0 + y1*z0 - z1*y0 + w1*x0,
  1176. -x1*z0 + y1*w0 + z1*x0 + w1*y0,
  1177. x1*y0 - y1*x0 + z1*w0 + w1*z0], dtype=numpy.float64)
  1178. def quaternion_conjugate(quaternion):
  1179. """Return conjugate of quaternion.
  1180. >>> q0 = random_quaternion()
  1181. >>> q1 = quaternion_conjugate(q0)
  1182. >>> q1[0] == q0[0] and all(q1[1:] == -q0[1:])
  1183. True
  1184. """
  1185. q = numpy.array(quaternion, dtype=numpy.float64, copy=True)
  1186. numpy.negative(q[1:], q[1:])
  1187. return q
  1188. def quaternion_inverse(quaternion):
  1189. """Return inverse of quaternion.
  1190. >>> q0 = random_quaternion()
  1191. >>> q1 = quaternion_inverse(q0)
  1192. >>> numpy.allclose(quaternion_multiply(q0, q1), [1, 0, 0, 0])
  1193. True
  1194. """
  1195. q = numpy.array(quaternion, dtype=numpy.float64, copy=True)
  1196. numpy.negative(q[1:], q[1:])
  1197. return q / numpy.dot(q, q)
  1198. def quaternion_real(quaternion):
  1199. """Return real part of quaternion.
  1200. >>> quaternion_real([3, 0, 1, 2])
  1201. 3.0
  1202. """
  1203. return float(quaternion[0])
  1204. def quaternion_imag(quaternion):
  1205. """Return imaginary part of quaternion.
  1206. >>> quaternion_imag([3, 0, 1, 2])
  1207. array([ 0., 1., 2.])
  1208. """
  1209. return numpy.array(quaternion[1:4], dtype=numpy.float64, copy=True)
  1210. def quaternion_slerp(quat0, quat1, fraction, spin=0, shortestpath=True):
  1211. """Return spherical linear interpolation between two quaternions.
  1212. >>> q0 = random_quaternion()
  1213. >>> q1 = random_quaternion()
  1214. >>> q = quaternion_slerp(q0, q1, 0)
  1215. >>> numpy.allclose(q, q0)
  1216. True
  1217. >>> q = quaternion_slerp(q0, q1, 1, 1)
  1218. >>> numpy.allclose(q, q1)
  1219. True
  1220. >>> q = quaternion_slerp(q0, q1, 0.5)
  1221. >>> angle = math.acos(numpy.dot(q0, q))
  1222. >>> numpy.allclose(2, math.acos(numpy.dot(q0, q1)) / angle) or \
  1223. numpy.allclose(2, math.acos(-numpy.dot(q0, q1)) / angle)
  1224. True
  1225. """
  1226. q0 = unit_vector(quat0[:4])
  1227. q1 = unit_vector(quat1[:4])
  1228. if fraction == 0.0:
  1229. return q0
  1230. elif fraction == 1.0:
  1231. return q1
  1232. d = numpy.dot(q0, q1)
  1233. if abs(abs(d) - 1.0) < _EPS:
  1234. return q0
  1235. if shortestpath and d < 0.0:
  1236. # invert rotation
  1237. d = -d
  1238. numpy.negative(q1, q1)
  1239. angle = math.acos(d) + spin * math.pi
  1240. if abs(angle) < _EPS:
  1241. return q0
  1242. isin = 1.0 / math.sin(angle)
  1243. q0 *= math.sin((1.0 - fraction) * angle) * isin
  1244. q1 *= math.sin(fraction * angle) * isin
  1245. q0 += q1
  1246. return q0
  1247. def random_quaternion(rand=None):
  1248. """Return uniform random unit quaternion.
  1249. rand: array like or None
  1250. Three independent random variables that are uniformly distributed
  1251. between 0 and 1.
  1252. >>> q = random_quaternion()
  1253. >>> numpy.allclose(1, vector_norm(q))
  1254. True
  1255. >>> q = random_quaternion(numpy.random.random(3))
  1256. >>> len(q.shape), q.shape[0]==4
  1257. (1, True)
  1258. """
  1259. if rand is None:
  1260. rand = numpy.random.rand(3)
  1261. else:
  1262. assert len(rand) == 3
  1263. r1 = numpy.sqrt(1.0 - rand[0])
  1264. r2 = numpy.sqrt(rand[0])
  1265. pi2 = math.pi * 2.0
  1266. t1 = pi2 * rand[1]
  1267. t2 = pi2 * rand[2]
  1268. return numpy.array([numpy.cos(t2)*r2, numpy.sin(t1)*r1,
  1269. numpy.cos(t1)*r1, numpy.sin(t2)*r2])
  1270. def random_rotation_matrix(rand=None):
  1271. """Return uniform random rotation matrix.
  1272. rand: array like
  1273. Three independent random variables that are uniformly distributed
  1274. between 0 and 1 for each returned quaternion.
  1275. >>> R = random_rotation_matrix()
  1276. >>> numpy.allclose(numpy.dot(R.T, R), numpy.identity(4))
  1277. True
  1278. """
  1279. return quaternion_matrix(random_quaternion(rand))
  1280. class Arcball(object):
  1281. """Virtual Trackball Control.
  1282. >>> ball = Arcball()
  1283. >>> ball = Arcball(initial=numpy.identity(4))
  1284. >>> ball.place([320, 320], 320)
  1285. >>> ball.down([500, 250])
  1286. >>> ball.drag([475, 275])
  1287. >>> R = ball.matrix()
  1288. >>> numpy.allclose(numpy.sum(R), 3.90583455)
  1289. True
  1290. >>> ball = Arcball(initial=[1, 0, 0, 0])
  1291. >>> ball.place([320, 320], 320)
  1292. >>> ball.setaxes([1, 1, 0], [-1, 1, 0])
  1293. >>> ball.constrain = True
  1294. >>> ball.down([400, 200])
  1295. >>> ball.drag([200, 400])
  1296. >>> R = ball.matrix()
  1297. >>> numpy.allclose(numpy.sum(R), 0.2055924)
  1298. True
  1299. >>> ball.next()
  1300. """
  1301. def __init__(self, initial=None):
  1302. """Initialize virtual trackball control.
  1303. initial : quaternion or rotation matrix
  1304. """
  1305. self._axis = None
  1306. self._axes = None
  1307. self._radius = 1.0
  1308. self._center = [0.0, 0.0]
  1309. self._vdown = numpy.array([0.0, 0.0, 1.0])
  1310. self._constrain = False
  1311. if initial is None:
  1312. self._qdown = numpy.array([1.0, 0.0, 0.0, 0.0])
  1313. else:
  1314. initial = numpy.array(initial, dtype=numpy.float64)
  1315. if initial.shape == (4, 4):
  1316. self._qdown = quaternion_from_matrix(initial)
  1317. elif initial.shape == (4, ):
  1318. initial /= vector_norm(initial)
  1319. self._qdown = initial
  1320. else:
  1321. raise ValueError("initial not a quaternion or matrix")
  1322. self._qnow = self._qpre = self._qdown
  1323. def place(self, center, radius):
  1324. """Place Arcball, e.g. when window size changes.
  1325. center : sequence[2]
  1326. Window coordinates of trackball center.
  1327. radius : float
  1328. Radius of trackball in window coordinates.
  1329. """
  1330. self._radius = float(radius)
  1331. self._center[0] = center[0]
  1332. self._center[1] = center[1]
  1333. def setaxes(self, *axes):
  1334. """Set axes to constrain rotations."""
  1335. if axes is None:
  1336. self._axes = None
  1337. else:
  1338. self._axes = [unit_vector(axis) for axis in axes]
  1339. @property
  1340. def constrain(self):
  1341. """Return state of constrain to axis mode."""
  1342. return self._constrain
  1343. @constrain.setter
  1344. def constrain(self, value):
  1345. """Set state of constrain to axis mode."""
  1346. self._constrain = bool(value)
  1347. def down(self, point):
  1348. """Set initial cursor window coordinates and pick constrain-axis."""
  1349. self._vdown = arcball_map_to_sphere(point, self._center, self._radius)
  1350. self._qdown = self._qpre = self._qnow
  1351. if self._constrain and self._axes is not None:
  1352. self._axis = arcball_nearest_axis(self._vdown, self._axes)
  1353. self._vdown = arcball_constrain_to_axis(self._vdown, self._axis)
  1354. else:
  1355. self._axis = None
  1356. def drag(self, point):
  1357. """Update current cursor window coordinates."""
  1358. vnow = arcball_map_to_sphere(point, self._center, self._radius)
  1359. if self._axis is not None:
  1360. vnow = arcball_constrain_to_axis(vnow, self._axis)
  1361. self._qpre = self._qnow
  1362. t = numpy.cross(self._vdown, vnow)
  1363. if numpy.dot(t, t) < _EPS:
  1364. self._qnow = self._qdown
  1365. else:
  1366. q = [numpy.dot(self._vdown, vnow), t[0], t[1], t[2]]
  1367. self._qnow = quaternion_multiply(q, self._qdown)
  1368. def next(self, acceleration=0.0):
  1369. """Continue rotation in direction of last drag."""
  1370. q = quaternion_slerp(self._qpre, self._qnow, 2.0+acceleration, False)
  1371. self._qpre, self._qnow = self._qnow, q
  1372. def matrix(self):
  1373. """Return homogeneous rotation matrix."""
  1374. return quaternion_matrix(self._qnow)
  1375. def arcball_map_to_sphere(point, center, radius):
  1376. """Return unit sphere coordinates from window coordinates."""
  1377. v0 = (point[0] - center[0]) / radius
  1378. v1 = (center[1] - point[1]) / radius
  1379. n = v0*v0 + v1*v1
  1380. if n > 1.0:
  1381. # position outside of sphere
  1382. n = math.sqrt(n)
  1383. return numpy.array([v0/n, v1/n, 0.0])
  1384. else:
  1385. return numpy.array([v0, v1, math.sqrt(1.0 - n)])
  1386. def arcball_constrain_to_axis(point, axis):
  1387. """Return sphere point perpendicular to axis."""
  1388. v = numpy.array(point, dtype=numpy.float64, copy=True)
  1389. a = numpy.array(axis, dtype=numpy.float64, copy=True)
  1390. v -= a * numpy.dot(a, v) # on plane
  1391. n = vector_norm(v)
  1392. if n > _EPS:
  1393. if v[2] < 0.0:
  1394. numpy.negative(v, v)
  1395. v /= n
  1396. return v
  1397. if a[2] == 1.0:
  1398. return numpy.array([1.0, 0.0, 0.0])
  1399. return unit_vector([-a[1], a[0], 0.0])
  1400. def arcball_nearest_axis(point, axes):
  1401. """Return axis, which arc is nearest to point."""
  1402. point = numpy.array(point, dtype=numpy.float64, copy=False)
  1403. nearest = None
  1404. mx = -1.0
  1405. for axis in axes:
  1406. t = numpy.dot(arcball_constrain_to_axis(point, axis), point)
  1407. if t > mx:
  1408. nearest = axis
  1409. mx = t
  1410. return nearest
  1411. # epsilon for testing whether a number is close to zero
  1412. _EPS = numpy.finfo(float).eps * 4.0
  1413. # axis sequences for Euler angles
  1414. _NEXT_AXIS = [1, 2, 0, 1]
  1415. # map axes strings to/from tuples of inner axis, parity, repetition, frame
  1416. _AXES2TUPLE = {
  1417. 'sxyz': (0, 0, 0, 0), 'sxyx': (0, 0, 1, 0), 'sxzy': (0, 1, 0, 0),
  1418. 'sxzx': (0, 1, 1, 0), 'syzx': (1, 0, 0, 0), 'syzy': (1, 0, 1, 0),
  1419. 'syxz': (1, 1, 0, 0), 'syxy': (1, 1, 1, 0), 'szxy': (2, 0, 0, 0),
  1420. 'szxz': (2, 0, 1, 0), 'szyx': (2, 1, 0, 0), 'szyz': (2, 1, 1, 0),
  1421. 'rzyx': (0, 0, 0, 1), 'rxyx': (0, 0, 1, 1), 'ryzx': (0, 1, 0, 1),
  1422. 'rxzx': (0, 1, 1, 1), 'rxzy': (1, 0, 0, 1), 'ryzy': (1, 0, 1, 1),
  1423. 'rzxy': (1, 1, 0, 1), 'ryxy': (1, 1, 1, 1), 'ryxz': (2, 0, 0, 1),
  1424. 'rzxz': (2, 0, 1, 1), 'rxyz': (2, 1, 0, 1), 'rzyz': (2, 1, 1, 1)}
  1425. _TUPLE2AXES = dict((v, k) for k, v in list(_AXES2TUPLE.items()))
  1426. def vector_norm(data, axis=None, out=None):
  1427. """Return length, i.e. Euclidean norm, of ndarray along axis.
  1428. >>> v = numpy.random.random(3)
  1429. >>> n = vector_norm(v)
  1430. >>> numpy.allclose(n, numpy.linalg.norm(v))
  1431. True
  1432. >>> v = numpy.random.rand(6, 5, 3)
  1433. >>> n = vector_norm(v, axis=-1)
  1434. >>> numpy.allclose(n, numpy.sqrt(numpy.sum(v*v, axis=2)))
  1435. True
  1436. >>> n = vector_norm(v, axis=1)
  1437. >>> numpy.allclose(n, numpy.sqrt(numpy.sum(v*v, axis=1)))
  1438. True
  1439. >>> v = numpy.random.rand(5, 4, 3)
  1440. >>> n = numpy.empty((5, 3))
  1441. >>> vector_norm(v, axis=1, out=n)
  1442. >>> numpy.allclose(n, numpy.sqrt(numpy.sum(v*v, axis=1)))
  1443. True
  1444. >>> vector_norm([])
  1445. 0.0
  1446. >>> vector_norm([1])
  1447. 1.0
  1448. """
  1449. data = numpy.array(data, dtype=numpy.float64, copy=True)
  1450. if out is None:
  1451. if data.ndim == 1:
  1452. return math.sqrt(numpy.dot(data, data))
  1453. data *= data
  1454. out = numpy.atleast_1d(numpy.sum(data, axis=axis))
  1455. numpy.sqrt(out, out)
  1456. return out
  1457. else:
  1458. data *= data
  1459. numpy.sum(data, axis=axis, out=out)
  1460. numpy.sqrt(out, out)
  1461. def unit_vector(data, axis=None, out=None):
  1462. """Return ndarray normalized by length, i.e. Euclidean norm, along axis.
  1463. >>> v0 = numpy.random.random(3)
  1464. >>> v1 = unit_vector(v0)
  1465. >>> numpy.allclose(v1, v0 / numpy.linalg.norm(v0))
  1466. True
  1467. >>> v0 = numpy.random.rand(5, 4, 3)
  1468. >>> v1 = unit_vector(v0, axis=-1)
  1469. >>> v2 = v0 / numpy.expand_dims(numpy.sqrt(numpy.sum(v0*v0, axis=2)), 2)
  1470. >>> numpy.allclose(v1, v2)
  1471. True
  1472. >>> v1 = unit_vector(v0, axis=1)
  1473. >>> v2 = v0 / numpy.expand_dims(numpy.sqrt(numpy.sum(v0*v0, axis=1)), 1)
  1474. >>> numpy.allclose(v1, v2)
  1475. True
  1476. >>> v1 = numpy.empty((5, 4, 3))
  1477. >>> unit_vector(v0, axis=1, out=v1)
  1478. >>> numpy.allclose(v1, v2)
  1479. True
  1480. >>> list(unit_vector([]))
  1481. []
  1482. >>> list(unit_vector([1]))
  1483. [1.0]
  1484. """
  1485. if out is None:
  1486. data = numpy.array(data, dtype=numpy.float64, copy=True)
  1487. if data.ndim == 1:
  1488. data /= math.sqrt(numpy.dot(data, data))
  1489. return data
  1490. else:
  1491. if out is not data:
  1492. out[:] = numpy.array(data, copy=False)
  1493. data = out
  1494. length = numpy.atleast_1d(numpy.sum(data*data, axis))
  1495. numpy.sqrt(length, length)
  1496. if axis is not None:
  1497. length = numpy.expand_dims(length, axis)
  1498. data /= length
  1499. if out is None:
  1500. return data
  1501. def random_vector(size):
  1502. """Return array of random doubles in the half-open interval [0.0, 1.0).
  1503. >>> v = random_vector(10000)
  1504. >>> numpy.all(v >= 0) and numpy.all(v < 1)
  1505. True
  1506. >>> v0 = random_vector(10)
  1507. >>> v1 = random_vector(10)
  1508. >>> numpy.any(v0 == v1)
  1509. False
  1510. """
  1511. return numpy.random.random(size)
  1512. def vector_product(v0, v1, axis=0):
  1513. """Return vector perpendicular to vectors.
  1514. >>> v = vector_product([2, 0, 0], [0, 3, 0])
  1515. >>> numpy.allclose(v, [0, 0, 6])
  1516. True
  1517. >>> v0 = [[2, 0, 0, 2], [0, 2, 0, 2], [0, 0, 2, 2]]
  1518. >>> v1 = [[3], [0], [0]]
  1519. >>> v = vector_product(v0, v1)
  1520. >>> numpy.allclose(v, [[0, 0, 0, 0], [0, 0, 6, 6], [0, -6, 0, -6]])
  1521. True
  1522. >>> v0 = [[2, 0, 0], [2, 0, 0], [0, 2, 0], [2, 0, 0]]
  1523. >>> v1 = [[0, 3, 0], [0, 0, 3], [0, 0, 3], [3, 3, 3]]
  1524. >>> v = vector_product(v0, v1, axis=1)
  1525. >>> numpy.allclose(v, [[0, 0, 6], [0, -6, 0], [6, 0, 0], [0, -6, 6]])
  1526. True
  1527. """
  1528. return numpy.cross(v0, v1, axis=axis)
  1529. def angle_between_vectors(v0, v1, directed=True, axis=0):
  1530. """Return angle between vectors.
  1531. If directed is False, the input vectors are interpreted as undirected axes,
  1532. i.e. the maximum angle is pi/2.
  1533. >>> a = angle_between_vectors([1, -2, 3], [-1, 2, -3])
  1534. >>> numpy.allclose(a, math.pi)
  1535. True
  1536. >>> a = angle_between_vectors([1, -2, 3], [-1, 2, -3], directed=False)
  1537. >>> numpy.allclose(a, 0)
  1538. True
  1539. >>> v0 = [[2, 0, 0, 2], [0, 2, 0, 2], [0, 0, 2, 2]]
  1540. >>> v1 = [[3], [0], [0]]
  1541. >>> a = angle_between_vectors(v0, v1)
  1542. >>> numpy.allclose(a, [0, 1.5708, 1.5708, 0.95532])
  1543. True
  1544. >>> v0 = [[2, 0, 0], [2, 0, 0], [0, 2, 0], [2, 0, 0]]
  1545. >>> v1 = [[0, 3, 0], [0, 0, 3], [0, 0, 3], [3, 3, 3]]
  1546. >>> a = angle_between_vectors(v0, v1, axis=1)
  1547. >>> numpy.allclose(a, [1.5708, 1.5708, 1.5708, 0.95532])
  1548. True
  1549. """
  1550. v0 = numpy.array(v0, dtype=numpy.float64, copy=False)
  1551. v1 = numpy.array(v1, dtype=numpy.float64, copy=False)
  1552. dot = numpy.sum(v0 * v1, axis=axis)
  1553. dot /= vector_norm(v0, axis=axis) * vector_norm(v1, axis=axis)
  1554. return numpy.arccos(dot if directed else numpy.fabs(dot))
  1555. def inverse_matrix(matrix):
  1556. """Return inverse of square transformation matrix.
  1557. >>> M0 = random_rotation_matrix()
  1558. >>> M1 = inverse_matrix(M0.T)
  1559. >>> numpy.allclose(M1, numpy.linalg.inv(M0.T))
  1560. True
  1561. >>> for size in range(1, 7):
  1562. ... M0 = numpy.random.rand(size, size)
  1563. ... M1 = inverse_matrix(M0)
  1564. ... if not numpy.allclose(M1, numpy.linalg.inv(M0)): print(size)
  1565. """
  1566. return numpy.linalg.inv(matrix)
  1567. def concatenate_matrices(*matrices):
  1568. """Return concatenation of series of transformation matrices.
  1569. >>> M = numpy.random.rand(16).reshape((4, 4)) - 0.5
  1570. >>> numpy.allclose(M, concatenate_matrices(M))
  1571. True
  1572. >>> numpy.allclose(numpy.dot(M, M.T), concatenate_matrices(M, M.T))
  1573. True
  1574. """
  1575. M = numpy.identity(4)
  1576. for i in matrices:
  1577. M = numpy.dot(M, i)
  1578. return M
  1579. def is_same_transform(matrix0, matrix1):
  1580. """Return True if two matrices perform same transformation.
  1581. >>> is_same_transform(numpy.identity(4), numpy.identity(4))
  1582. True
  1583. >>> is_same_transform(numpy.identity(4), random_rotation_matrix())
  1584. False
  1585. """
  1586. matrix0 = numpy.array(matrix0, dtype=numpy.float64, copy=True)
  1587. matrix0 /= matrix0[3, 3]
  1588. matrix1 = numpy.array(matrix1, dtype=numpy.float64, copy=True)
  1589. matrix1 /= matrix1[3, 3]
  1590. return numpy.allclose(matrix0, matrix1)
  1591. # def _import_module(name, package=None, warn=True, prefix='_py_', ignore='_'):
  1592. # """Try import all public attributes from module into global namespace.
  1593. #
  1594. # Existing attributes with name clashes are renamed with prefix.
  1595. # Attributes starting with underscore are ignored by default.
  1596. #
  1597. # Return True on successful import.
  1598. #
  1599. # """
  1600. # import warnings
  1601. # from importlib import import_module
  1602. # try:
  1603. # if not package:
  1604. # module = import_module(name)
  1605. # else:
  1606. # module = import_module('.' + name, package=package)
  1607. # except ImportError:
  1608. # if warn:
  1609. # warnings.warn("failed to import module %s" % name)
  1610. # else:
  1611. # for attr in dir(module):
  1612. # if ignore and attr.startswith(ignore):
  1613. # continue
  1614. # if prefix:
  1615. # if attr in globals():
  1616. # globals()[prefix + attr] = globals()[attr]
  1617. # elif warn:
  1618. # warnings.warn("no Python implementation of " + attr)
  1619. # globals()[attr] = getattr(module, attr)
  1620. # return True
  1621. #
  1622. #
  1623. # _import_module('_transformations')
  1624. #
  1625. # if __name__ == "__main__":
  1626. # import doctest
  1627. # import random # used in doctests
  1628. # numpy.set_printoptions(suppress=True, precision=5)
  1629. # doctest.testmod()