dense.py 30 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094
  1. from __future__ import annotations
  2. import random
  3. from sympy.core.basic import Basic
  4. from sympy.core.singleton import S
  5. from sympy.core.symbol import Symbol
  6. from sympy.core.sympify import sympify
  7. from sympy.functions.elementary.trigonometric import cos, sin
  8. from sympy.utilities.decorator import doctest_depends_on
  9. from sympy.utilities.exceptions import sympy_deprecation_warning
  10. from sympy.utilities.iterables import is_sequence
  11. from .exceptions import ShapeError
  12. from .decompositions import _cholesky, _LDLdecomposition
  13. from .matrixbase import MatrixBase
  14. from .repmatrix import MutableRepMatrix, RepMatrix
  15. from .solvers import _lower_triangular_solve, _upper_triangular_solve
  16. __doctest_requires__ = {('symarray',): ['numpy']}
  17. def _iszero(x):
  18. """Returns True if x is zero."""
  19. return x.is_zero
  20. class DenseMatrix(RepMatrix):
  21. """Matrix implementation based on DomainMatrix as the internal representation"""
  22. #
  23. # DenseMatrix is a superclass for both MutableDenseMatrix and
  24. # ImmutableDenseMatrix. Methods shared by both classes but not for the
  25. # Sparse classes should be implemented here.
  26. #
  27. is_MatrixExpr: bool = False
  28. _op_priority = 10.01
  29. _class_priority = 4
  30. @property
  31. def _mat(self):
  32. sympy_deprecation_warning(
  33. """
  34. The private _mat attribute of Matrix is deprecated. Use the
  35. .flat() method instead.
  36. """,
  37. deprecated_since_version="1.9",
  38. active_deprecations_target="deprecated-private-matrix-attributes"
  39. )
  40. return self.flat()
  41. def _eval_inverse(self, **kwargs):
  42. return self.inv(method=kwargs.get('method', 'GE'),
  43. iszerofunc=kwargs.get('iszerofunc', _iszero),
  44. try_block_diag=kwargs.get('try_block_diag', False))
  45. def as_immutable(self):
  46. """Returns an Immutable version of this Matrix
  47. """
  48. from .immutable import ImmutableDenseMatrix as cls
  49. return cls._fromrep(self._rep.copy())
  50. def as_mutable(self):
  51. """Returns a mutable version of this matrix
  52. Examples
  53. ========
  54. >>> from sympy import ImmutableMatrix
  55. >>> X = ImmutableMatrix([[1, 2], [3, 4]])
  56. >>> Y = X.as_mutable()
  57. >>> Y[1, 1] = 5 # Can set values in Y
  58. >>> Y
  59. Matrix([
  60. [1, 2],
  61. [3, 5]])
  62. """
  63. return Matrix(self)
  64. def cholesky(self, hermitian=True):
  65. return _cholesky(self, hermitian=hermitian)
  66. def LDLdecomposition(self, hermitian=True):
  67. return _LDLdecomposition(self, hermitian=hermitian)
  68. def lower_triangular_solve(self, rhs):
  69. return _lower_triangular_solve(self, rhs)
  70. def upper_triangular_solve(self, rhs):
  71. return _upper_triangular_solve(self, rhs)
  72. cholesky.__doc__ = _cholesky.__doc__
  73. LDLdecomposition.__doc__ = _LDLdecomposition.__doc__
  74. lower_triangular_solve.__doc__ = _lower_triangular_solve.__doc__
  75. upper_triangular_solve.__doc__ = _upper_triangular_solve.__doc__
  76. def _force_mutable(x):
  77. """Return a matrix as a Matrix, otherwise return x."""
  78. if getattr(x, 'is_Matrix', False):
  79. return x.as_mutable()
  80. elif isinstance(x, Basic):
  81. return x
  82. elif hasattr(x, '__array__'):
  83. a = x.__array__()
  84. if len(a.shape) == 0:
  85. return sympify(a)
  86. return Matrix(x)
  87. return x
  88. class MutableDenseMatrix(DenseMatrix, MutableRepMatrix):
  89. def simplify(self, **kwargs):
  90. """Applies simplify to the elements of a matrix in place.
  91. This is a shortcut for M.applyfunc(lambda x: simplify(x, ratio, measure))
  92. See Also
  93. ========
  94. sympy.simplify.simplify.simplify
  95. """
  96. from sympy.simplify.simplify import simplify as _simplify
  97. for (i, j), element in self.todok().items():
  98. self[i, j] = _simplify(element, **kwargs)
  99. MutableMatrix = Matrix = MutableDenseMatrix
  100. ###########
  101. # Numpy Utility Functions:
  102. # list2numpy, matrix2numpy, symmarray
  103. ###########
  104. def list2numpy(l, dtype=object): # pragma: no cover
  105. """Converts Python list of SymPy expressions to a NumPy array.
  106. See Also
  107. ========
  108. matrix2numpy
  109. """
  110. from numpy import empty
  111. a = empty(len(l), dtype)
  112. for i, s in enumerate(l):
  113. a[i] = s
  114. return a
  115. def matrix2numpy(m, dtype=object): # pragma: no cover
  116. """Converts SymPy's matrix to a NumPy array.
  117. See Also
  118. ========
  119. list2numpy
  120. """
  121. from numpy import empty
  122. a = empty(m.shape, dtype)
  123. for i in range(m.rows):
  124. for j in range(m.cols):
  125. a[i, j] = m[i, j]
  126. return a
  127. ###########
  128. # Rotation matrices:
  129. # rot_givens, rot_axis[123], rot_ccw_axis[123]
  130. ###########
  131. def rot_givens(i, j, theta, dim=3):
  132. r"""Returns a a Givens rotation matrix, a a rotation in the
  133. plane spanned by two coordinates axes.
  134. Explanation
  135. ===========
  136. The Givens rotation corresponds to a generalization of rotation
  137. matrices to any number of dimensions, given by:
  138. .. math::
  139. G(i, j, \theta) =
  140. \begin{bmatrix}
  141. 1 & \cdots & 0 & \cdots & 0 & \cdots & 0 \\
  142. \vdots & \ddots & \vdots & & \vdots & & \vdots \\
  143. 0 & \cdots & c & \cdots & -s & \cdots & 0 \\
  144. \vdots & & \vdots & \ddots & \vdots & & \vdots \\
  145. 0 & \cdots & s & \cdots & c & \cdots & 0 \\
  146. \vdots & & \vdots & & \vdots & \ddots & \vdots \\
  147. 0 & \cdots & 0 & \cdots & 0 & \cdots & 1
  148. \end{bmatrix}
  149. Where $c = \cos(\theta)$ and $s = \sin(\theta)$ appear at the intersections
  150. ``i``\th and ``j``\th rows and columns.
  151. For fixed ``i > j``\, the non-zero elements of a Givens matrix are
  152. given by:
  153. - $g_{kk} = 1$ for $k \ne i,\,j$
  154. - $g_{kk} = c$ for $k = i,\,j$
  155. - $g_{ji} = -g_{ij} = -s$
  156. Parameters
  157. ==========
  158. i : int between ``0`` and ``dim - 1``
  159. Represents first axis
  160. j : int between ``0`` and ``dim - 1``
  161. Represents second axis
  162. dim : int bigger than 1
  163. Number of dimensions. Defaults to 3.
  164. Examples
  165. ========
  166. >>> from sympy import pi, rot_givens
  167. A counterclockwise rotation of pi/3 (60 degrees) around
  168. the third axis (z-axis):
  169. >>> rot_givens(1, 0, pi/3)
  170. Matrix([
  171. [ 1/2, -sqrt(3)/2, 0],
  172. [sqrt(3)/2, 1/2, 0],
  173. [ 0, 0, 1]])
  174. If we rotate by pi/2 (90 degrees):
  175. >>> rot_givens(1, 0, pi/2)
  176. Matrix([
  177. [0, -1, 0],
  178. [1, 0, 0],
  179. [0, 0, 1]])
  180. This can be generalized to any number
  181. of dimensions:
  182. >>> rot_givens(1, 0, pi/2, dim=4)
  183. Matrix([
  184. [0, -1, 0, 0],
  185. [1, 0, 0, 0],
  186. [0, 0, 1, 0],
  187. [0, 0, 0, 1]])
  188. References
  189. ==========
  190. .. [1] https://en.wikipedia.org/wiki/Givens_rotation
  191. See Also
  192. ========
  193. rot_axis1: Returns a rotation matrix for a rotation of theta (in radians)
  194. about the 1-axis (clockwise around the x axis)
  195. rot_axis2: Returns a rotation matrix for a rotation of theta (in radians)
  196. about the 2-axis (clockwise around the y axis)
  197. rot_axis3: Returns a rotation matrix for a rotation of theta (in radians)
  198. about the 3-axis (clockwise around the z axis)
  199. rot_ccw_axis1: Returns a rotation matrix for a rotation of theta (in radians)
  200. about the 1-axis (counterclockwise around the x axis)
  201. rot_ccw_axis2: Returns a rotation matrix for a rotation of theta (in radians)
  202. about the 2-axis (counterclockwise around the y axis)
  203. rot_ccw_axis3: Returns a rotation matrix for a rotation of theta (in radians)
  204. about the 3-axis (counterclockwise around the z axis)
  205. """
  206. if not isinstance(dim, int) or dim < 2:
  207. raise ValueError('dim must be an integer biggen than one, '
  208. 'got {}.'.format(dim))
  209. if i == j:
  210. raise ValueError('i and j must be different, '
  211. 'got ({}, {})'.format(i, j))
  212. for ij in [i, j]:
  213. if not isinstance(ij, int) or ij < 0 or ij > dim - 1:
  214. raise ValueError('i and j must be integers between 0 and '
  215. '{}, got i={} and j={}.'.format(dim-1, i, j))
  216. theta = sympify(theta)
  217. c = cos(theta)
  218. s = sin(theta)
  219. M = eye(dim)
  220. M[i, i] = c
  221. M[j, j] = c
  222. M[i, j] = s
  223. M[j, i] = -s
  224. return M
  225. def rot_axis3(theta):
  226. r"""Returns a rotation matrix for a rotation of theta (in radians)
  227. about the 3-axis.
  228. Explanation
  229. ===========
  230. For a right-handed coordinate system, this corresponds to a
  231. clockwise rotation around the `z`-axis, given by:
  232. .. math::
  233. R = \begin{bmatrix}
  234. \cos(\theta) & \sin(\theta) & 0 \\
  235. -\sin(\theta) & \cos(\theta) & 0 \\
  236. 0 & 0 & 1
  237. \end{bmatrix}
  238. Examples
  239. ========
  240. >>> from sympy import pi, rot_axis3
  241. A rotation of pi/3 (60 degrees):
  242. >>> theta = pi/3
  243. >>> rot_axis3(theta)
  244. Matrix([
  245. [ 1/2, sqrt(3)/2, 0],
  246. [-sqrt(3)/2, 1/2, 0],
  247. [ 0, 0, 1]])
  248. If we rotate by pi/2 (90 degrees):
  249. >>> rot_axis3(pi/2)
  250. Matrix([
  251. [ 0, 1, 0],
  252. [-1, 0, 0],
  253. [ 0, 0, 1]])
  254. See Also
  255. ========
  256. rot_givens: Returns a Givens rotation matrix (generalized rotation for
  257. any number of dimensions)
  258. rot_ccw_axis3: Returns a rotation matrix for a rotation of theta (in radians)
  259. about the 3-axis (counterclockwise around the z axis)
  260. rot_axis1: Returns a rotation matrix for a rotation of theta (in radians)
  261. about the 1-axis (clockwise around the x axis)
  262. rot_axis2: Returns a rotation matrix for a rotation of theta (in radians)
  263. about the 2-axis (clockwise around the y axis)
  264. """
  265. return rot_givens(0, 1, theta, dim=3)
  266. def rot_axis2(theta):
  267. r"""Returns a rotation matrix for a rotation of theta (in radians)
  268. about the 2-axis.
  269. Explanation
  270. ===========
  271. For a right-handed coordinate system, this corresponds to a
  272. clockwise rotation around the `y`-axis, given by:
  273. .. math::
  274. R = \begin{bmatrix}
  275. \cos(\theta) & 0 & -\sin(\theta) \\
  276. 0 & 1 & 0 \\
  277. \sin(\theta) & 0 & \cos(\theta)
  278. \end{bmatrix}
  279. Examples
  280. ========
  281. >>> from sympy import pi, rot_axis2
  282. A rotation of pi/3 (60 degrees):
  283. >>> theta = pi/3
  284. >>> rot_axis2(theta)
  285. Matrix([
  286. [ 1/2, 0, -sqrt(3)/2],
  287. [ 0, 1, 0],
  288. [sqrt(3)/2, 0, 1/2]])
  289. If we rotate by pi/2 (90 degrees):
  290. >>> rot_axis2(pi/2)
  291. Matrix([
  292. [0, 0, -1],
  293. [0, 1, 0],
  294. [1, 0, 0]])
  295. See Also
  296. ========
  297. rot_givens: Returns a Givens rotation matrix (generalized rotation for
  298. any number of dimensions)
  299. rot_ccw_axis2: Returns a rotation matrix for a rotation of theta (in radians)
  300. about the 2-axis (clockwise around the y axis)
  301. rot_axis1: Returns a rotation matrix for a rotation of theta (in radians)
  302. about the 1-axis (counterclockwise around the x axis)
  303. rot_axis3: Returns a rotation matrix for a rotation of theta (in radians)
  304. about the 3-axis (counterclockwise around the z axis)
  305. """
  306. return rot_givens(2, 0, theta, dim=3)
  307. def rot_axis1(theta):
  308. r"""Returns a rotation matrix for a rotation of theta (in radians)
  309. about the 1-axis.
  310. Explanation
  311. ===========
  312. For a right-handed coordinate system, this corresponds to a
  313. clockwise rotation around the `x`-axis, given by:
  314. .. math::
  315. R = \begin{bmatrix}
  316. 1 & 0 & 0 \\
  317. 0 & \cos(\theta) & \sin(\theta) \\
  318. 0 & -\sin(\theta) & \cos(\theta)
  319. \end{bmatrix}
  320. Examples
  321. ========
  322. >>> from sympy import pi, rot_axis1
  323. A rotation of pi/3 (60 degrees):
  324. >>> theta = pi/3
  325. >>> rot_axis1(theta)
  326. Matrix([
  327. [1, 0, 0],
  328. [0, 1/2, sqrt(3)/2],
  329. [0, -sqrt(3)/2, 1/2]])
  330. If we rotate by pi/2 (90 degrees):
  331. >>> rot_axis1(pi/2)
  332. Matrix([
  333. [1, 0, 0],
  334. [0, 0, 1],
  335. [0, -1, 0]])
  336. See Also
  337. ========
  338. rot_givens: Returns a Givens rotation matrix (generalized rotation for
  339. any number of dimensions)
  340. rot_ccw_axis1: Returns a rotation matrix for a rotation of theta (in radians)
  341. about the 1-axis (counterclockwise around the x axis)
  342. rot_axis2: Returns a rotation matrix for a rotation of theta (in radians)
  343. about the 2-axis (clockwise around the y axis)
  344. rot_axis3: Returns a rotation matrix for a rotation of theta (in radians)
  345. about the 3-axis (clockwise around the z axis)
  346. """
  347. return rot_givens(1, 2, theta, dim=3)
  348. def rot_ccw_axis3(theta):
  349. r"""Returns a rotation matrix for a rotation of theta (in radians)
  350. about the 3-axis.
  351. Explanation
  352. ===========
  353. For a right-handed coordinate system, this corresponds to a
  354. counterclockwise rotation around the `z`-axis, given by:
  355. .. math::
  356. R = \begin{bmatrix}
  357. \cos(\theta) & -\sin(\theta) & 0 \\
  358. \sin(\theta) & \cos(\theta) & 0 \\
  359. 0 & 0 & 1
  360. \end{bmatrix}
  361. Examples
  362. ========
  363. >>> from sympy import pi, rot_ccw_axis3
  364. A rotation of pi/3 (60 degrees):
  365. >>> theta = pi/3
  366. >>> rot_ccw_axis3(theta)
  367. Matrix([
  368. [ 1/2, -sqrt(3)/2, 0],
  369. [sqrt(3)/2, 1/2, 0],
  370. [ 0, 0, 1]])
  371. If we rotate by pi/2 (90 degrees):
  372. >>> rot_ccw_axis3(pi/2)
  373. Matrix([
  374. [0, -1, 0],
  375. [1, 0, 0],
  376. [0, 0, 1]])
  377. See Also
  378. ========
  379. rot_givens: Returns a Givens rotation matrix (generalized rotation for
  380. any number of dimensions)
  381. rot_axis3: Returns a rotation matrix for a rotation of theta (in radians)
  382. about the 3-axis (clockwise around the z axis)
  383. rot_ccw_axis1: Returns a rotation matrix for a rotation of theta (in radians)
  384. about the 1-axis (counterclockwise around the x axis)
  385. rot_ccw_axis2: Returns a rotation matrix for a rotation of theta (in radians)
  386. about the 2-axis (counterclockwise around the y axis)
  387. """
  388. return rot_givens(1, 0, theta, dim=3)
  389. def rot_ccw_axis2(theta):
  390. r"""Returns a rotation matrix for a rotation of theta (in radians)
  391. about the 2-axis.
  392. Explanation
  393. ===========
  394. For a right-handed coordinate system, this corresponds to a
  395. counterclockwise rotation around the `y`-axis, given by:
  396. .. math::
  397. R = \begin{bmatrix}
  398. \cos(\theta) & 0 & \sin(\theta) \\
  399. 0 & 1 & 0 \\
  400. -\sin(\theta) & 0 & \cos(\theta)
  401. \end{bmatrix}
  402. Examples
  403. ========
  404. >>> from sympy import pi, rot_ccw_axis2
  405. A rotation of pi/3 (60 degrees):
  406. >>> theta = pi/3
  407. >>> rot_ccw_axis2(theta)
  408. Matrix([
  409. [ 1/2, 0, sqrt(3)/2],
  410. [ 0, 1, 0],
  411. [-sqrt(3)/2, 0, 1/2]])
  412. If we rotate by pi/2 (90 degrees):
  413. >>> rot_ccw_axis2(pi/2)
  414. Matrix([
  415. [ 0, 0, 1],
  416. [ 0, 1, 0],
  417. [-1, 0, 0]])
  418. See Also
  419. ========
  420. rot_givens: Returns a Givens rotation matrix (generalized rotation for
  421. any number of dimensions)
  422. rot_axis2: Returns a rotation matrix for a rotation of theta (in radians)
  423. about the 2-axis (clockwise around the y axis)
  424. rot_ccw_axis1: Returns a rotation matrix for a rotation of theta (in radians)
  425. about the 1-axis (counterclockwise around the x axis)
  426. rot_ccw_axis3: Returns a rotation matrix for a rotation of theta (in radians)
  427. about the 3-axis (counterclockwise around the z axis)
  428. """
  429. return rot_givens(0, 2, theta, dim=3)
  430. def rot_ccw_axis1(theta):
  431. r"""Returns a rotation matrix for a rotation of theta (in radians)
  432. about the 1-axis.
  433. Explanation
  434. ===========
  435. For a right-handed coordinate system, this corresponds to a
  436. counterclockwise rotation around the `x`-axis, given by:
  437. .. math::
  438. R = \begin{bmatrix}
  439. 1 & 0 & 0 \\
  440. 0 & \cos(\theta) & -\sin(\theta) \\
  441. 0 & \sin(\theta) & \cos(\theta)
  442. \end{bmatrix}
  443. Examples
  444. ========
  445. >>> from sympy import pi, rot_ccw_axis1
  446. A rotation of pi/3 (60 degrees):
  447. >>> theta = pi/3
  448. >>> rot_ccw_axis1(theta)
  449. Matrix([
  450. [1, 0, 0],
  451. [0, 1/2, -sqrt(3)/2],
  452. [0, sqrt(3)/2, 1/2]])
  453. If we rotate by pi/2 (90 degrees):
  454. >>> rot_ccw_axis1(pi/2)
  455. Matrix([
  456. [1, 0, 0],
  457. [0, 0, -1],
  458. [0, 1, 0]])
  459. See Also
  460. ========
  461. rot_givens: Returns a Givens rotation matrix (generalized rotation for
  462. any number of dimensions)
  463. rot_axis1: Returns a rotation matrix for a rotation of theta (in radians)
  464. about the 1-axis (clockwise around the x axis)
  465. rot_ccw_axis2: Returns a rotation matrix for a rotation of theta (in radians)
  466. about the 2-axis (counterclockwise around the y axis)
  467. rot_ccw_axis3: Returns a rotation matrix for a rotation of theta (in radians)
  468. about the 3-axis (counterclockwise around the z axis)
  469. """
  470. return rot_givens(2, 1, theta, dim=3)
  471. @doctest_depends_on(modules=('numpy',))
  472. def symarray(prefix, shape, **kwargs): # pragma: no cover
  473. r"""Create a numpy ndarray of symbols (as an object array).
  474. The created symbols are named ``prefix_i1_i2_``... You should thus provide a
  475. non-empty prefix if you want your symbols to be unique for different output
  476. arrays, as SymPy symbols with identical names are the same object.
  477. Parameters
  478. ----------
  479. prefix : string
  480. A prefix prepended to the name of every symbol.
  481. shape : int or tuple
  482. Shape of the created array. If an int, the array is one-dimensional; for
  483. more than one dimension the shape must be a tuple.
  484. \*\*kwargs : dict
  485. keyword arguments passed on to Symbol
  486. Examples
  487. ========
  488. These doctests require numpy.
  489. >>> from sympy import symarray
  490. >>> symarray('', 3)
  491. [_0 _1 _2]
  492. If you want multiple symarrays to contain distinct symbols, you *must*
  493. provide unique prefixes:
  494. >>> a = symarray('', 3)
  495. >>> b = symarray('', 3)
  496. >>> a[0] == b[0]
  497. True
  498. >>> a = symarray('a', 3)
  499. >>> b = symarray('b', 3)
  500. >>> a[0] == b[0]
  501. False
  502. Creating symarrays with a prefix:
  503. >>> symarray('a', 3)
  504. [a_0 a_1 a_2]
  505. For more than one dimension, the shape must be given as a tuple:
  506. >>> symarray('a', (2, 3))
  507. [[a_0_0 a_0_1 a_0_2]
  508. [a_1_0 a_1_1 a_1_2]]
  509. >>> symarray('a', (2, 3, 2))
  510. [[[a_0_0_0 a_0_0_1]
  511. [a_0_1_0 a_0_1_1]
  512. [a_0_2_0 a_0_2_1]]
  513. <BLANKLINE>
  514. [[a_1_0_0 a_1_0_1]
  515. [a_1_1_0 a_1_1_1]
  516. [a_1_2_0 a_1_2_1]]]
  517. For setting assumptions of the underlying Symbols:
  518. >>> [s.is_real for s in symarray('a', 2, real=True)]
  519. [True, True]
  520. """
  521. from numpy import empty, ndindex
  522. arr = empty(shape, dtype=object)
  523. for index in ndindex(shape):
  524. arr[index] = Symbol('%s_%s' % (prefix, '_'.join(map(str, index))),
  525. **kwargs)
  526. return arr
  527. ###############
  528. # Functions
  529. ###############
  530. def casoratian(seqs, n, zero=True):
  531. """Given linear difference operator L of order 'k' and homogeneous
  532. equation Ly = 0 we want to compute kernel of L, which is a set
  533. of 'k' sequences: a(n), b(n), ... z(n).
  534. Solutions of L are linearly independent iff their Casoratian,
  535. denoted as C(a, b, ..., z), do not vanish for n = 0.
  536. Casoratian is defined by k x k determinant::
  537. + a(n) b(n) . . . z(n) +
  538. | a(n+1) b(n+1) . . . z(n+1) |
  539. | . . . . |
  540. | . . . . |
  541. | . . . . |
  542. + a(n+k-1) b(n+k-1) . . . z(n+k-1) +
  543. It proves very useful in rsolve_hyper() where it is applied
  544. to a generating set of a recurrence to factor out linearly
  545. dependent solutions and return a basis:
  546. >>> from sympy import Symbol, casoratian, factorial
  547. >>> n = Symbol('n', integer=True)
  548. Exponential and factorial are linearly independent:
  549. >>> casoratian([2**n, factorial(n)], n) != 0
  550. True
  551. """
  552. seqs = list(map(sympify, seqs))
  553. if not zero:
  554. f = lambda i, j: seqs[j].subs(n, n + i)
  555. else:
  556. f = lambda i, j: seqs[j].subs(n, i)
  557. k = len(seqs)
  558. return Matrix(k, k, f).det()
  559. def eye(*args, **kwargs):
  560. """Create square identity matrix n x n
  561. See Also
  562. ========
  563. diag
  564. zeros
  565. ones
  566. """
  567. return Matrix.eye(*args, **kwargs)
  568. def diag(*values, strict=True, unpack=False, **kwargs):
  569. """Returns a matrix with the provided values placed on the
  570. diagonal. If non-square matrices are included, they will
  571. produce a block-diagonal matrix.
  572. Examples
  573. ========
  574. This version of diag is a thin wrapper to Matrix.diag that differs
  575. in that it treats all lists like matrices -- even when a single list
  576. is given. If this is not desired, either put a `*` before the list or
  577. set `unpack=True`.
  578. >>> from sympy import diag
  579. >>> diag([1, 2, 3], unpack=True) # = diag(1,2,3) or diag(*[1,2,3])
  580. Matrix([
  581. [1, 0, 0],
  582. [0, 2, 0],
  583. [0, 0, 3]])
  584. >>> diag([1, 2, 3]) # a column vector
  585. Matrix([
  586. [1],
  587. [2],
  588. [3]])
  589. See Also
  590. ========
  591. .matrixbase.MatrixBase.eye
  592. .matrixbase.MatrixBase.diagonal
  593. .matrixbase.MatrixBase.diag
  594. .expressions.blockmatrix.BlockMatrix
  595. """
  596. return Matrix.diag(*values, strict=strict, unpack=unpack, **kwargs)
  597. def GramSchmidt(vlist, orthonormal=False):
  598. """Apply the Gram-Schmidt process to a set of vectors.
  599. Parameters
  600. ==========
  601. vlist : List of Matrix
  602. Vectors to be orthogonalized for.
  603. orthonormal : Bool, optional
  604. If true, return an orthonormal basis.
  605. Returns
  606. =======
  607. vlist : List of Matrix
  608. Orthogonalized vectors
  609. Notes
  610. =====
  611. This routine is mostly duplicate from ``Matrix.orthogonalize``,
  612. except for some difference that this always raises error when
  613. linearly dependent vectors are found, and the keyword ``normalize``
  614. has been named as ``orthonormal`` in this function.
  615. See Also
  616. ========
  617. .matrixbase.MatrixBase.orthogonalize
  618. References
  619. ==========
  620. .. [1] https://en.wikipedia.org/wiki/Gram%E2%80%93Schmidt_process
  621. """
  622. return MutableDenseMatrix.orthogonalize(
  623. *vlist, normalize=orthonormal, rankcheck=True
  624. )
  625. def hessian(f, varlist, constraints=()):
  626. """Compute Hessian matrix for a function f wrt parameters in varlist
  627. which may be given as a sequence or a row/column vector. A list of
  628. constraints may optionally be given.
  629. Examples
  630. ========
  631. >>> from sympy import Function, hessian, pprint
  632. >>> from sympy.abc import x, y
  633. >>> f = Function('f')(x, y)
  634. >>> g1 = Function('g')(x, y)
  635. >>> g2 = x**2 + 3*y
  636. >>> pprint(hessian(f, (x, y), [g1, g2]))
  637. [ d d ]
  638. [ 0 0 --(g(x, y)) --(g(x, y)) ]
  639. [ dx dy ]
  640. [ ]
  641. [ 0 0 2*x 3 ]
  642. [ ]
  643. [ 2 2 ]
  644. [d d d ]
  645. [--(g(x, y)) 2*x ---(f(x, y)) -----(f(x, y))]
  646. [dx 2 dy dx ]
  647. [ dx ]
  648. [ ]
  649. [ 2 2 ]
  650. [d d d ]
  651. [--(g(x, y)) 3 -----(f(x, y)) ---(f(x, y)) ]
  652. [dy dy dx 2 ]
  653. [ dy ]
  654. References
  655. ==========
  656. .. [1] https://en.wikipedia.org/wiki/Hessian_matrix
  657. See Also
  658. ========
  659. sympy.matrices.matrixbase.MatrixBase.jacobian
  660. wronskian
  661. """
  662. # f is the expression representing a function f, return regular matrix
  663. if isinstance(varlist, MatrixBase):
  664. if 1 not in varlist.shape:
  665. raise ShapeError("`varlist` must be a column or row vector.")
  666. if varlist.cols == 1:
  667. varlist = varlist.T
  668. varlist = varlist.tolist()[0]
  669. if is_sequence(varlist):
  670. n = len(varlist)
  671. if not n:
  672. raise ShapeError("`len(varlist)` must not be zero.")
  673. else:
  674. raise ValueError("Improper variable list in hessian function")
  675. if not getattr(f, 'diff'):
  676. # check differentiability
  677. raise ValueError("Function `f` (%s) is not differentiable" % f)
  678. m = len(constraints)
  679. N = m + n
  680. out = zeros(N)
  681. for k, g in enumerate(constraints):
  682. if not getattr(g, 'diff'):
  683. # check differentiability
  684. raise ValueError("Function `f` (%s) is not differentiable" % f)
  685. for i in range(n):
  686. out[k, i + m] = g.diff(varlist[i])
  687. for i in range(n):
  688. for j in range(i, n):
  689. out[i + m, j + m] = f.diff(varlist[i]).diff(varlist[j])
  690. for i in range(N):
  691. for j in range(i + 1, N):
  692. out[j, i] = out[i, j]
  693. return out
  694. def jordan_cell(eigenval, n):
  695. """
  696. Create a Jordan block:
  697. Examples
  698. ========
  699. >>> from sympy import jordan_cell
  700. >>> from sympy.abc import x
  701. >>> jordan_cell(x, 4)
  702. Matrix([
  703. [x, 1, 0, 0],
  704. [0, x, 1, 0],
  705. [0, 0, x, 1],
  706. [0, 0, 0, x]])
  707. """
  708. return Matrix.jordan_block(size=n, eigenvalue=eigenval)
  709. def matrix_multiply_elementwise(A, B):
  710. """Return the Hadamard product (elementwise product) of A and B
  711. >>> from sympy import Matrix, matrix_multiply_elementwise
  712. >>> A = Matrix([[0, 1, 2], [3, 4, 5]])
  713. >>> B = Matrix([[1, 10, 100], [100, 10, 1]])
  714. >>> matrix_multiply_elementwise(A, B)
  715. Matrix([
  716. [ 0, 10, 200],
  717. [300, 40, 5]])
  718. See Also
  719. ========
  720. sympy.matrices.matrixbase.MatrixBase.__mul__
  721. """
  722. return A.multiply_elementwise(B)
  723. def ones(*args, **kwargs):
  724. """Returns a matrix of ones with ``rows`` rows and ``cols`` columns;
  725. if ``cols`` is omitted a square matrix will be returned.
  726. See Also
  727. ========
  728. zeros
  729. eye
  730. diag
  731. """
  732. if 'c' in kwargs:
  733. kwargs['cols'] = kwargs.pop('c')
  734. return Matrix.ones(*args, **kwargs)
  735. def randMatrix(r, c=None, min=0, max=99, seed=None, symmetric=False,
  736. percent=100, prng=None):
  737. """Create random matrix with dimensions ``r`` x ``c``. If ``c`` is omitted
  738. the matrix will be square. If ``symmetric`` is True the matrix must be
  739. square. If ``percent`` is less than 100 then only approximately the given
  740. percentage of elements will be non-zero.
  741. The pseudo-random number generator used to generate matrix is chosen in the
  742. following way.
  743. * If ``prng`` is supplied, it will be used as random number generator.
  744. It should be an instance of ``random.Random``, or at least have
  745. ``randint`` and ``shuffle`` methods with same signatures.
  746. * if ``prng`` is not supplied but ``seed`` is supplied, then new
  747. ``random.Random`` with given ``seed`` will be created;
  748. * otherwise, a new ``random.Random`` with default seed will be used.
  749. Examples
  750. ========
  751. >>> from sympy import randMatrix
  752. >>> randMatrix(3) # doctest:+SKIP
  753. [25, 45, 27]
  754. [44, 54, 9]
  755. [23, 96, 46]
  756. >>> randMatrix(3, 2) # doctest:+SKIP
  757. [87, 29]
  758. [23, 37]
  759. [90, 26]
  760. >>> randMatrix(3, 3, 0, 2) # doctest:+SKIP
  761. [0, 2, 0]
  762. [2, 0, 1]
  763. [0, 0, 1]
  764. >>> randMatrix(3, symmetric=True) # doctest:+SKIP
  765. [85, 26, 29]
  766. [26, 71, 43]
  767. [29, 43, 57]
  768. >>> A = randMatrix(3, seed=1)
  769. >>> B = randMatrix(3, seed=2)
  770. >>> A == B
  771. False
  772. >>> A == randMatrix(3, seed=1)
  773. True
  774. >>> randMatrix(3, symmetric=True, percent=50) # doctest:+SKIP
  775. [77, 70, 0],
  776. [70, 0, 0],
  777. [ 0, 0, 88]
  778. """
  779. # Note that ``Random()`` is equivalent to ``Random(None)``
  780. prng = prng or random.Random(seed)
  781. if c is None:
  782. c = r
  783. if symmetric and r != c:
  784. raise ValueError('For symmetric matrices, r must equal c, but %i != %i' % (r, c))
  785. ij = range(r * c)
  786. if percent != 100:
  787. ij = prng.sample(ij, int(len(ij)*percent // 100))
  788. m = zeros(r, c)
  789. if not symmetric:
  790. for ijk in ij:
  791. i, j = divmod(ijk, c)
  792. m[i, j] = prng.randint(min, max)
  793. else:
  794. for ijk in ij:
  795. i, j = divmod(ijk, c)
  796. if i <= j:
  797. m[i, j] = m[j, i] = prng.randint(min, max)
  798. return m
  799. def wronskian(functions, var, method='bareiss'):
  800. """
  801. Compute Wronskian for [] of functions
  802. ::
  803. | f1 f2 ... fn |
  804. | f1' f2' ... fn' |
  805. | . . . . |
  806. W(f1, ..., fn) = | . . . . |
  807. | . . . . |
  808. | (n) (n) (n) |
  809. | D (f1) D (f2) ... D (fn) |
  810. see: https://en.wikipedia.org/wiki/Wronskian
  811. See Also
  812. ========
  813. sympy.matrices.matrixbase.MatrixBase.jacobian
  814. hessian
  815. """
  816. functions = [sympify(f) for f in functions]
  817. n = len(functions)
  818. if n == 0:
  819. return S.One
  820. W = Matrix(n, n, lambda i, j: functions[i].diff(var, j))
  821. return W.det(method)
  822. def zeros(*args, **kwargs):
  823. """Returns a matrix of zeros with ``rows`` rows and ``cols`` columns;
  824. if ``cols`` is omitted a square matrix will be returned.
  825. See Also
  826. ========
  827. ones
  828. eye
  829. diag
  830. """
  831. if 'c' in kwargs:
  832. kwargs['cols'] = kwargs.pop('c')
  833. return Matrix.zeros(*args, **kwargs)