special_matrices.py 29 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038
  1. from __future__ import division, print_function, absolute_import
  2. import math
  3. import numpy as np
  4. from scipy._lib.six import xrange
  5. from scipy._lib.six import string_types
  6. from numpy.lib.stride_tricks import as_strided
  7. __all__ = ['tri', 'tril', 'triu', 'toeplitz', 'circulant', 'hankel',
  8. 'hadamard', 'leslie', 'kron', 'block_diag', 'companion',
  9. 'helmert', 'hilbert', 'invhilbert', 'pascal', 'invpascal', 'dft']
  10. #-----------------------------------------------------------------------------
  11. # matrix construction functions
  12. #-----------------------------------------------------------------------------
  13. #
  14. # *Note*: tri{,u,l} is implemented in numpy, but an important bug was fixed in
  15. # 2.0.0.dev-1af2f3, the following tri{,u,l} definitions are here for backwards
  16. # compatibility.
  17. def tri(N, M=None, k=0, dtype=None):
  18. """
  19. Construct (N, M) matrix filled with ones at and below the k-th diagonal.
  20. The matrix has A[i,j] == 1 for i <= j + k
  21. Parameters
  22. ----------
  23. N : int
  24. The size of the first dimension of the matrix.
  25. M : int or None, optional
  26. The size of the second dimension of the matrix. If `M` is None,
  27. `M = N` is assumed.
  28. k : int, optional
  29. Number of subdiagonal below which matrix is filled with ones.
  30. `k` = 0 is the main diagonal, `k` < 0 subdiagonal and `k` > 0
  31. superdiagonal.
  32. dtype : dtype, optional
  33. Data type of the matrix.
  34. Returns
  35. -------
  36. tri : (N, M) ndarray
  37. Tri matrix.
  38. Examples
  39. --------
  40. >>> from scipy.linalg import tri
  41. >>> tri(3, 5, 2, dtype=int)
  42. array([[1, 1, 1, 0, 0],
  43. [1, 1, 1, 1, 0],
  44. [1, 1, 1, 1, 1]])
  45. >>> tri(3, 5, -1, dtype=int)
  46. array([[0, 0, 0, 0, 0],
  47. [1, 0, 0, 0, 0],
  48. [1, 1, 0, 0, 0]])
  49. """
  50. if M is None:
  51. M = N
  52. if isinstance(M, string_types):
  53. #pearu: any objections to remove this feature?
  54. # As tri(N,'d') is equivalent to tri(N,dtype='d')
  55. dtype = M
  56. M = N
  57. m = np.greater_equal.outer(np.arange(k, N+k), np.arange(M))
  58. if dtype is None:
  59. return m
  60. else:
  61. return m.astype(dtype)
  62. def tril(m, k=0):
  63. """
  64. Make a copy of a matrix with elements above the k-th diagonal zeroed.
  65. Parameters
  66. ----------
  67. m : array_like
  68. Matrix whose elements to return
  69. k : int, optional
  70. Diagonal above which to zero elements.
  71. `k` == 0 is the main diagonal, `k` < 0 subdiagonal and
  72. `k` > 0 superdiagonal.
  73. Returns
  74. -------
  75. tril : ndarray
  76. Return is the same shape and type as `m`.
  77. Examples
  78. --------
  79. >>> from scipy.linalg import tril
  80. >>> tril([[1,2,3],[4,5,6],[7,8,9],[10,11,12]], -1)
  81. array([[ 0, 0, 0],
  82. [ 4, 0, 0],
  83. [ 7, 8, 0],
  84. [10, 11, 12]])
  85. """
  86. m = np.asarray(m)
  87. out = tri(m.shape[0], m.shape[1], k=k, dtype=m.dtype.char) * m
  88. return out
  89. def triu(m, k=0):
  90. """
  91. Make a copy of a matrix with elements below the k-th diagonal zeroed.
  92. Parameters
  93. ----------
  94. m : array_like
  95. Matrix whose elements to return
  96. k : int, optional
  97. Diagonal below which to zero elements.
  98. `k` == 0 is the main diagonal, `k` < 0 subdiagonal and
  99. `k` > 0 superdiagonal.
  100. Returns
  101. -------
  102. triu : ndarray
  103. Return matrix with zeroed elements below the k-th diagonal and has
  104. same shape and type as `m`.
  105. Examples
  106. --------
  107. >>> from scipy.linalg import triu
  108. >>> triu([[1,2,3],[4,5,6],[7,8,9],[10,11,12]], -1)
  109. array([[ 1, 2, 3],
  110. [ 4, 5, 6],
  111. [ 0, 8, 9],
  112. [ 0, 0, 12]])
  113. """
  114. m = np.asarray(m)
  115. out = (1 - tri(m.shape[0], m.shape[1], k - 1, m.dtype.char)) * m
  116. return out
  117. def toeplitz(c, r=None):
  118. """
  119. Construct a Toeplitz matrix.
  120. The Toeplitz matrix has constant diagonals, with c as its first column
  121. and r as its first row. If r is not given, ``r == conjugate(c)`` is
  122. assumed.
  123. Parameters
  124. ----------
  125. c : array_like
  126. First column of the matrix. Whatever the actual shape of `c`, it
  127. will be converted to a 1-D array.
  128. r : array_like, optional
  129. First row of the matrix. If None, ``r = conjugate(c)`` is assumed;
  130. in this case, if c[0] is real, the result is a Hermitian matrix.
  131. r[0] is ignored; the first row of the returned matrix is
  132. ``[c[0], r[1:]]``. Whatever the actual shape of `r`, it will be
  133. converted to a 1-D array.
  134. Returns
  135. -------
  136. A : (len(c), len(r)) ndarray
  137. The Toeplitz matrix. Dtype is the same as ``(c[0] + r[0]).dtype``.
  138. See Also
  139. --------
  140. circulant : circulant matrix
  141. hankel : Hankel matrix
  142. solve_toeplitz : Solve a Toeplitz system.
  143. Notes
  144. -----
  145. The behavior when `c` or `r` is a scalar, or when `c` is complex and
  146. `r` is None, was changed in version 0.8.0. The behavior in previous
  147. versions was undocumented and is no longer supported.
  148. Examples
  149. --------
  150. >>> from scipy.linalg import toeplitz
  151. >>> toeplitz([1,2,3], [1,4,5,6])
  152. array([[1, 4, 5, 6],
  153. [2, 1, 4, 5],
  154. [3, 2, 1, 4]])
  155. >>> toeplitz([1.0, 2+3j, 4-1j])
  156. array([[ 1.+0.j, 2.-3.j, 4.+1.j],
  157. [ 2.+3.j, 1.+0.j, 2.-3.j],
  158. [ 4.-1.j, 2.+3.j, 1.+0.j]])
  159. """
  160. c = np.asarray(c).ravel()
  161. if r is None:
  162. r = c.conjugate()
  163. else:
  164. r = np.asarray(r).ravel()
  165. # Form a 1D array containing a reversed c followed by r[1:] that could be
  166. # strided to give us toeplitz matrix.
  167. vals = np.concatenate((c[::-1], r[1:]))
  168. out_shp = len(c), len(r)
  169. n = vals.strides[0]
  170. return as_strided(vals[len(c)-1:], shape=out_shp, strides=(-n, n)).copy()
  171. def circulant(c):
  172. """
  173. Construct a circulant matrix.
  174. Parameters
  175. ----------
  176. c : (N,) array_like
  177. 1-D array, the first column of the matrix.
  178. Returns
  179. -------
  180. A : (N, N) ndarray
  181. A circulant matrix whose first column is `c`.
  182. See Also
  183. --------
  184. toeplitz : Toeplitz matrix
  185. hankel : Hankel matrix
  186. solve_circulant : Solve a circulant system.
  187. Notes
  188. -----
  189. .. versionadded:: 0.8.0
  190. Examples
  191. --------
  192. >>> from scipy.linalg import circulant
  193. >>> circulant([1, 2, 3])
  194. array([[1, 3, 2],
  195. [2, 1, 3],
  196. [3, 2, 1]])
  197. """
  198. c = np.asarray(c).ravel()
  199. # Form an extended array that could be strided to give circulant version
  200. c_ext = np.concatenate((c[::-1], c[:0:-1]))
  201. L = len(c)
  202. n = c_ext.strides[0]
  203. return as_strided(c_ext[L-1:], shape=(L, L), strides=(-n, n)).copy()
  204. def hankel(c, r=None):
  205. """
  206. Construct a Hankel matrix.
  207. The Hankel matrix has constant anti-diagonals, with `c` as its
  208. first column and `r` as its last row. If `r` is not given, then
  209. `r = zeros_like(c)` is assumed.
  210. Parameters
  211. ----------
  212. c : array_like
  213. First column of the matrix. Whatever the actual shape of `c`, it
  214. will be converted to a 1-D array.
  215. r : array_like, optional
  216. Last row of the matrix. If None, ``r = zeros_like(c)`` is assumed.
  217. r[0] is ignored; the last row of the returned matrix is
  218. ``[c[-1], r[1:]]``. Whatever the actual shape of `r`, it will be
  219. converted to a 1-D array.
  220. Returns
  221. -------
  222. A : (len(c), len(r)) ndarray
  223. The Hankel matrix. Dtype is the same as ``(c[0] + r[0]).dtype``.
  224. See Also
  225. --------
  226. toeplitz : Toeplitz matrix
  227. circulant : circulant matrix
  228. Examples
  229. --------
  230. >>> from scipy.linalg import hankel
  231. >>> hankel([1, 17, 99])
  232. array([[ 1, 17, 99],
  233. [17, 99, 0],
  234. [99, 0, 0]])
  235. >>> hankel([1,2,3,4], [4,7,7,8,9])
  236. array([[1, 2, 3, 4, 7],
  237. [2, 3, 4, 7, 7],
  238. [3, 4, 7, 7, 8],
  239. [4, 7, 7, 8, 9]])
  240. """
  241. c = np.asarray(c).ravel()
  242. if r is None:
  243. r = np.zeros_like(c)
  244. else:
  245. r = np.asarray(r).ravel()
  246. # Form a 1D array of values to be used in the matrix, containing `c`
  247. # followed by r[1:].
  248. vals = np.concatenate((c, r[1:]))
  249. # Stride on concatenated array to get hankel matrix
  250. out_shp = len(c), len(r)
  251. n = vals.strides[0]
  252. return as_strided(vals, shape=out_shp, strides=(n, n)).copy()
  253. def hadamard(n, dtype=int):
  254. """
  255. Construct a Hadamard matrix.
  256. Constructs an n-by-n Hadamard matrix, using Sylvester's
  257. construction. `n` must be a power of 2.
  258. Parameters
  259. ----------
  260. n : int
  261. The order of the matrix. `n` must be a power of 2.
  262. dtype : dtype, optional
  263. The data type of the array to be constructed.
  264. Returns
  265. -------
  266. H : (n, n) ndarray
  267. The Hadamard matrix.
  268. Notes
  269. -----
  270. .. versionadded:: 0.8.0
  271. Examples
  272. --------
  273. >>> from scipy.linalg import hadamard
  274. >>> hadamard(2, dtype=complex)
  275. array([[ 1.+0.j, 1.+0.j],
  276. [ 1.+0.j, -1.-0.j]])
  277. >>> hadamard(4)
  278. array([[ 1, 1, 1, 1],
  279. [ 1, -1, 1, -1],
  280. [ 1, 1, -1, -1],
  281. [ 1, -1, -1, 1]])
  282. """
  283. # This function is a slightly modified version of the
  284. # function contributed by Ivo in ticket #675.
  285. if n < 1:
  286. lg2 = 0
  287. else:
  288. lg2 = int(math.log(n, 2))
  289. if 2 ** lg2 != n:
  290. raise ValueError("n must be an positive integer, and n must be "
  291. "a power of 2")
  292. H = np.array([[1]], dtype=dtype)
  293. # Sylvester's construction
  294. for i in range(0, lg2):
  295. H = np.vstack((np.hstack((H, H)), np.hstack((H, -H))))
  296. return H
  297. def leslie(f, s):
  298. """
  299. Create a Leslie matrix.
  300. Given the length n array of fecundity coefficients `f` and the length
  301. n-1 array of survival coefficients `s`, return the associated Leslie matrix.
  302. Parameters
  303. ----------
  304. f : (N,) array_like
  305. The "fecundity" coefficients.
  306. s : (N-1,) array_like
  307. The "survival" coefficients, has to be 1-D. The length of `s`
  308. must be one less than the length of `f`, and it must be at least 1.
  309. Returns
  310. -------
  311. L : (N, N) ndarray
  312. The array is zero except for the first row,
  313. which is `f`, and the first sub-diagonal, which is `s`.
  314. The data-type of the array will be the data-type of ``f[0]+s[0]``.
  315. Notes
  316. -----
  317. .. versionadded:: 0.8.0
  318. The Leslie matrix is used to model discrete-time, age-structured
  319. population growth [1]_ [2]_. In a population with `n` age classes, two sets
  320. of parameters define a Leslie matrix: the `n` "fecundity coefficients",
  321. which give the number of offspring per-capita produced by each age
  322. class, and the `n` - 1 "survival coefficients", which give the
  323. per-capita survival rate of each age class.
  324. References
  325. ----------
  326. .. [1] P. H. Leslie, On the use of matrices in certain population
  327. mathematics, Biometrika, Vol. 33, No. 3, 183--212 (Nov. 1945)
  328. .. [2] P. H. Leslie, Some further notes on the use of matrices in
  329. population mathematics, Biometrika, Vol. 35, No. 3/4, 213--245
  330. (Dec. 1948)
  331. Examples
  332. --------
  333. >>> from scipy.linalg import leslie
  334. >>> leslie([0.1, 2.0, 1.0, 0.1], [0.2, 0.8, 0.7])
  335. array([[ 0.1, 2. , 1. , 0.1],
  336. [ 0.2, 0. , 0. , 0. ],
  337. [ 0. , 0.8, 0. , 0. ],
  338. [ 0. , 0. , 0.7, 0. ]])
  339. """
  340. f = np.atleast_1d(f)
  341. s = np.atleast_1d(s)
  342. if f.ndim != 1:
  343. raise ValueError("Incorrect shape for f. f must be one-dimensional")
  344. if s.ndim != 1:
  345. raise ValueError("Incorrect shape for s. s must be one-dimensional")
  346. if f.size != s.size + 1:
  347. raise ValueError("Incorrect lengths for f and s. The length"
  348. " of s must be one less than the length of f.")
  349. if s.size == 0:
  350. raise ValueError("The length of s must be at least 1.")
  351. tmp = f[0] + s[0]
  352. n = f.size
  353. a = np.zeros((n, n), dtype=tmp.dtype)
  354. a[0] = f
  355. a[list(range(1, n)), list(range(0, n - 1))] = s
  356. return a
  357. def kron(a, b):
  358. """
  359. Kronecker product.
  360. The result is the block matrix::
  361. a[0,0]*b a[0,1]*b ... a[0,-1]*b
  362. a[1,0]*b a[1,1]*b ... a[1,-1]*b
  363. ...
  364. a[-1,0]*b a[-1,1]*b ... a[-1,-1]*b
  365. Parameters
  366. ----------
  367. a : (M, N) ndarray
  368. Input array
  369. b : (P, Q) ndarray
  370. Input array
  371. Returns
  372. -------
  373. A : (M*P, N*Q) ndarray
  374. Kronecker product of `a` and `b`.
  375. Examples
  376. --------
  377. >>> from numpy import array
  378. >>> from scipy.linalg import kron
  379. >>> kron(array([[1,2],[3,4]]), array([[1,1,1]]))
  380. array([[1, 1, 1, 2, 2, 2],
  381. [3, 3, 3, 4, 4, 4]])
  382. """
  383. if not a.flags['CONTIGUOUS']:
  384. a = np.reshape(a, a.shape)
  385. if not b.flags['CONTIGUOUS']:
  386. b = np.reshape(b, b.shape)
  387. o = np.outer(a, b)
  388. o = o.reshape(a.shape + b.shape)
  389. return np.concatenate(np.concatenate(o, axis=1), axis=1)
  390. def block_diag(*arrs):
  391. """
  392. Create a block diagonal matrix from provided arrays.
  393. Given the inputs `A`, `B` and `C`, the output will have these
  394. arrays arranged on the diagonal::
  395. [[A, 0, 0],
  396. [0, B, 0],
  397. [0, 0, C]]
  398. Parameters
  399. ----------
  400. A, B, C, ... : array_like, up to 2-D
  401. Input arrays. A 1-D array or array_like sequence of length `n` is
  402. treated as a 2-D array with shape ``(1,n)``.
  403. Returns
  404. -------
  405. D : ndarray
  406. Array with `A`, `B`, `C`, ... on the diagonal. `D` has the
  407. same dtype as `A`.
  408. Notes
  409. -----
  410. If all the input arrays are square, the output is known as a
  411. block diagonal matrix.
  412. Empty sequences (i.e., array-likes of zero size) will not be ignored.
  413. Noteworthy, both [] and [[]] are treated as matrices with shape ``(1,0)``.
  414. Examples
  415. --------
  416. >>> from scipy.linalg import block_diag
  417. >>> A = [[1, 0],
  418. ... [0, 1]]
  419. >>> B = [[3, 4, 5],
  420. ... [6, 7, 8]]
  421. >>> C = [[7]]
  422. >>> P = np.zeros((2, 0), dtype='int32')
  423. >>> block_diag(A, B, C)
  424. array([[1, 0, 0, 0, 0, 0],
  425. [0, 1, 0, 0, 0, 0],
  426. [0, 0, 3, 4, 5, 0],
  427. [0, 0, 6, 7, 8, 0],
  428. [0, 0, 0, 0, 0, 7]])
  429. >>> block_diag(A, P, B, C)
  430. array([[1, 0, 0, 0, 0, 0],
  431. [0, 1, 0, 0, 0, 0],
  432. [0, 0, 0, 0, 0, 0],
  433. [0, 0, 0, 0, 0, 0],
  434. [0, 0, 3, 4, 5, 0],
  435. [0, 0, 6, 7, 8, 0],
  436. [0, 0, 0, 0, 0, 7]])
  437. >>> block_diag(1.0, [2, 3], [[4, 5], [6, 7]])
  438. array([[ 1., 0., 0., 0., 0.],
  439. [ 0., 2., 3., 0., 0.],
  440. [ 0., 0., 0., 4., 5.],
  441. [ 0., 0., 0., 6., 7.]])
  442. """
  443. if arrs == ():
  444. arrs = ([],)
  445. arrs = [np.atleast_2d(a) for a in arrs]
  446. bad_args = [k for k in range(len(arrs)) if arrs[k].ndim > 2]
  447. if bad_args:
  448. raise ValueError("arguments in the following positions have dimension "
  449. "greater than 2: %s" % bad_args)
  450. shapes = np.array([a.shape for a in arrs])
  451. out_dtype = np.find_common_type([arr.dtype for arr in arrs], [])
  452. out = np.zeros(np.sum(shapes, axis=0), dtype=out_dtype)
  453. r, c = 0, 0
  454. for i, (rr, cc) in enumerate(shapes):
  455. out[r:r + rr, c:c + cc] = arrs[i]
  456. r += rr
  457. c += cc
  458. return out
  459. def companion(a):
  460. """
  461. Create a companion matrix.
  462. Create the companion matrix [1]_ associated with the polynomial whose
  463. coefficients are given in `a`.
  464. Parameters
  465. ----------
  466. a : (N,) array_like
  467. 1-D array of polynomial coefficients. The length of `a` must be
  468. at least two, and ``a[0]`` must not be zero.
  469. Returns
  470. -------
  471. c : (N-1, N-1) ndarray
  472. The first row of `c` is ``-a[1:]/a[0]``, and the first
  473. sub-diagonal is all ones. The data-type of the array is the same
  474. as the data-type of ``1.0*a[0]``.
  475. Raises
  476. ------
  477. ValueError
  478. If any of the following are true: a) ``a.ndim != 1``;
  479. b) ``a.size < 2``; c) ``a[0] == 0``.
  480. Notes
  481. -----
  482. .. versionadded:: 0.8.0
  483. References
  484. ----------
  485. .. [1] R. A. Horn & C. R. Johnson, *Matrix Analysis*. Cambridge, UK:
  486. Cambridge University Press, 1999, pp. 146-7.
  487. Examples
  488. --------
  489. >>> from scipy.linalg import companion
  490. >>> companion([1, -10, 31, -30])
  491. array([[ 10., -31., 30.],
  492. [ 1., 0., 0.],
  493. [ 0., 1., 0.]])
  494. """
  495. a = np.atleast_1d(a)
  496. if a.ndim != 1:
  497. raise ValueError("Incorrect shape for `a`. `a` must be "
  498. "one-dimensional.")
  499. if a.size < 2:
  500. raise ValueError("The length of `a` must be at least 2.")
  501. if a[0] == 0:
  502. raise ValueError("The first coefficient in `a` must not be zero.")
  503. first_row = -a[1:] / (1.0 * a[0])
  504. n = a.size
  505. c = np.zeros((n - 1, n - 1), dtype=first_row.dtype)
  506. c[0] = first_row
  507. c[list(range(1, n - 1)), list(range(0, n - 2))] = 1
  508. return c
  509. def helmert(n, full=False):
  510. """
  511. Create a Helmert matrix of order `n`.
  512. This has applications in statistics, compositional or simplicial analysis,
  513. and in Aitchison geometry.
  514. Parameters
  515. ----------
  516. n : int
  517. The size of the array to create.
  518. full : bool, optional
  519. If True the (n, n) ndarray will be returned.
  520. Otherwise the submatrix that does not include the first
  521. row will be returned.
  522. Default: False.
  523. Returns
  524. -------
  525. M : ndarray
  526. The Helmert matrix.
  527. The shape is (n, n) or (n-1, n) depending on the `full` argument.
  528. Examples
  529. --------
  530. >>> from scipy.linalg import helmert
  531. >>> helmert(5, full=True)
  532. array([[ 0.4472136 , 0.4472136 , 0.4472136 , 0.4472136 , 0.4472136 ],
  533. [ 0.70710678, -0.70710678, 0. , 0. , 0. ],
  534. [ 0.40824829, 0.40824829, -0.81649658, 0. , 0. ],
  535. [ 0.28867513, 0.28867513, 0.28867513, -0.8660254 , 0. ],
  536. [ 0.2236068 , 0.2236068 , 0.2236068 , 0.2236068 , -0.89442719]])
  537. """
  538. H = np.tril(np.ones((n, n)), -1) - np.diag(np.arange(n))
  539. d = np.arange(n) * np.arange(1, n+1)
  540. H[0] = 1
  541. d[0] = n
  542. H_full = H / np.sqrt(d)[:, np.newaxis]
  543. if full:
  544. return H_full
  545. else:
  546. return H_full[1:]
  547. def hilbert(n):
  548. """
  549. Create a Hilbert matrix of order `n`.
  550. Returns the `n` by `n` array with entries `h[i,j] = 1 / (i + j + 1)`.
  551. Parameters
  552. ----------
  553. n : int
  554. The size of the array to create.
  555. Returns
  556. -------
  557. h : (n, n) ndarray
  558. The Hilbert matrix.
  559. See Also
  560. --------
  561. invhilbert : Compute the inverse of a Hilbert matrix.
  562. Notes
  563. -----
  564. .. versionadded:: 0.10.0
  565. Examples
  566. --------
  567. >>> from scipy.linalg import hilbert
  568. >>> hilbert(3)
  569. array([[ 1. , 0.5 , 0.33333333],
  570. [ 0.5 , 0.33333333, 0.25 ],
  571. [ 0.33333333, 0.25 , 0.2 ]])
  572. """
  573. values = 1.0 / (1.0 + np.arange(2 * n - 1))
  574. h = hankel(values[:n], r=values[n - 1:])
  575. return h
  576. def invhilbert(n, exact=False):
  577. """
  578. Compute the inverse of the Hilbert matrix of order `n`.
  579. The entries in the inverse of a Hilbert matrix are integers. When `n`
  580. is greater than 14, some entries in the inverse exceed the upper limit
  581. of 64 bit integers. The `exact` argument provides two options for
  582. dealing with these large integers.
  583. Parameters
  584. ----------
  585. n : int
  586. The order of the Hilbert matrix.
  587. exact : bool, optional
  588. If False, the data type of the array that is returned is np.float64,
  589. and the array is an approximation of the inverse.
  590. If True, the array is the exact integer inverse array. To represent
  591. the exact inverse when n > 14, the returned array is an object array
  592. of long integers. For n <= 14, the exact inverse is returned as an
  593. array with data type np.int64.
  594. Returns
  595. -------
  596. invh : (n, n) ndarray
  597. The data type of the array is np.float64 if `exact` is False.
  598. If `exact` is True, the data type is either np.int64 (for n <= 14)
  599. or object (for n > 14). In the latter case, the objects in the
  600. array will be long integers.
  601. See Also
  602. --------
  603. hilbert : Create a Hilbert matrix.
  604. Notes
  605. -----
  606. .. versionadded:: 0.10.0
  607. Examples
  608. --------
  609. >>> from scipy.linalg import invhilbert
  610. >>> invhilbert(4)
  611. array([[ 16., -120., 240., -140.],
  612. [ -120., 1200., -2700., 1680.],
  613. [ 240., -2700., 6480., -4200.],
  614. [ -140., 1680., -4200., 2800.]])
  615. >>> invhilbert(4, exact=True)
  616. array([[ 16, -120, 240, -140],
  617. [ -120, 1200, -2700, 1680],
  618. [ 240, -2700, 6480, -4200],
  619. [ -140, 1680, -4200, 2800]], dtype=int64)
  620. >>> invhilbert(16)[7,7]
  621. 4.2475099528537506e+19
  622. >>> invhilbert(16, exact=True)[7,7]
  623. 42475099528537378560L
  624. """
  625. from scipy.special import comb
  626. if exact:
  627. if n > 14:
  628. dtype = object
  629. else:
  630. dtype = np.int64
  631. else:
  632. dtype = np.float64
  633. invh = np.empty((n, n), dtype=dtype)
  634. for i in xrange(n):
  635. for j in xrange(0, i + 1):
  636. s = i + j
  637. invh[i, j] = ((-1) ** s * (s + 1) *
  638. comb(n + i, n - j - 1, exact) *
  639. comb(n + j, n - i - 1, exact) *
  640. comb(s, i, exact) ** 2)
  641. if i != j:
  642. invh[j, i] = invh[i, j]
  643. return invh
  644. def pascal(n, kind='symmetric', exact=True):
  645. """
  646. Returns the n x n Pascal matrix.
  647. The Pascal matrix is a matrix containing the binomial coefficients as
  648. its elements.
  649. Parameters
  650. ----------
  651. n : int
  652. The size of the matrix to create; that is, the result is an n x n
  653. matrix.
  654. kind : str, optional
  655. Must be one of 'symmetric', 'lower', or 'upper'.
  656. Default is 'symmetric'.
  657. exact : bool, optional
  658. If `exact` is True, the result is either an array of type
  659. numpy.uint64 (if n < 35) or an object array of Python long integers.
  660. If `exact` is False, the coefficients in the matrix are computed using
  661. `scipy.special.comb` with `exact=False`. The result will be a floating
  662. point array, and the values in the array will not be the exact
  663. coefficients, but this version is much faster than `exact=True`.
  664. Returns
  665. -------
  666. p : (n, n) ndarray
  667. The Pascal matrix.
  668. See Also
  669. --------
  670. invpascal
  671. Notes
  672. -----
  673. See https://en.wikipedia.org/wiki/Pascal_matrix for more information
  674. about Pascal matrices.
  675. .. versionadded:: 0.11.0
  676. Examples
  677. --------
  678. >>> from scipy.linalg import pascal
  679. >>> pascal(4)
  680. array([[ 1, 1, 1, 1],
  681. [ 1, 2, 3, 4],
  682. [ 1, 3, 6, 10],
  683. [ 1, 4, 10, 20]], dtype=uint64)
  684. >>> pascal(4, kind='lower')
  685. array([[1, 0, 0, 0],
  686. [1, 1, 0, 0],
  687. [1, 2, 1, 0],
  688. [1, 3, 3, 1]], dtype=uint64)
  689. >>> pascal(50)[-1, -1]
  690. 25477612258980856902730428600L
  691. >>> from scipy.special import comb
  692. >>> comb(98, 49, exact=True)
  693. 25477612258980856902730428600L
  694. """
  695. from scipy.special import comb
  696. if kind not in ['symmetric', 'lower', 'upper']:
  697. raise ValueError("kind must be 'symmetric', 'lower', or 'upper'")
  698. if exact:
  699. if n >= 35:
  700. L_n = np.empty((n, n), dtype=object)
  701. L_n.fill(0)
  702. else:
  703. L_n = np.zeros((n, n), dtype=np.uint64)
  704. for i in range(n):
  705. for j in range(i + 1):
  706. L_n[i, j] = comb(i, j, exact=True)
  707. else:
  708. L_n = comb(*np.ogrid[:n, :n])
  709. if kind == 'lower':
  710. p = L_n
  711. elif kind == 'upper':
  712. p = L_n.T
  713. else:
  714. p = np.dot(L_n, L_n.T)
  715. return p
  716. def invpascal(n, kind='symmetric', exact=True):
  717. """
  718. Returns the inverse of the n x n Pascal matrix.
  719. The Pascal matrix is a matrix containing the binomial coefficients as
  720. its elements.
  721. Parameters
  722. ----------
  723. n : int
  724. The size of the matrix to create; that is, the result is an n x n
  725. matrix.
  726. kind : str, optional
  727. Must be one of 'symmetric', 'lower', or 'upper'.
  728. Default is 'symmetric'.
  729. exact : bool, optional
  730. If `exact` is True, the result is either an array of type
  731. `numpy.int64` (if `n` <= 35) or an object array of Python integers.
  732. If `exact` is False, the coefficients in the matrix are computed using
  733. `scipy.special.comb` with `exact=False`. The result will be a floating
  734. point array, and for large `n`, the values in the array will not be the
  735. exact coefficients.
  736. Returns
  737. -------
  738. invp : (n, n) ndarray
  739. The inverse of the Pascal matrix.
  740. See Also
  741. --------
  742. pascal
  743. Notes
  744. -----
  745. .. versionadded:: 0.16.0
  746. References
  747. ----------
  748. .. [1] "Pascal matrix", https://en.wikipedia.org/wiki/Pascal_matrix
  749. .. [2] Cohen, A. M., "The inverse of a Pascal matrix", Mathematical
  750. Gazette, 59(408), pp. 111-112, 1975.
  751. Examples
  752. --------
  753. >>> from scipy.linalg import invpascal, pascal
  754. >>> invp = invpascal(5)
  755. >>> invp
  756. array([[ 5, -10, 10, -5, 1],
  757. [-10, 30, -35, 19, -4],
  758. [ 10, -35, 46, -27, 6],
  759. [ -5, 19, -27, 17, -4],
  760. [ 1, -4, 6, -4, 1]])
  761. >>> p = pascal(5)
  762. >>> p.dot(invp)
  763. array([[ 1., 0., 0., 0., 0.],
  764. [ 0., 1., 0., 0., 0.],
  765. [ 0., 0., 1., 0., 0.],
  766. [ 0., 0., 0., 1., 0.],
  767. [ 0., 0., 0., 0., 1.]])
  768. An example of the use of `kind` and `exact`:
  769. >>> invpascal(5, kind='lower', exact=False)
  770. array([[ 1., -0., 0., -0., 0.],
  771. [-1., 1., -0., 0., -0.],
  772. [ 1., -2., 1., -0., 0.],
  773. [-1., 3., -3., 1., -0.],
  774. [ 1., -4., 6., -4., 1.]])
  775. """
  776. from scipy.special import comb
  777. if kind not in ['symmetric', 'lower', 'upper']:
  778. raise ValueError("'kind' must be 'symmetric', 'lower' or 'upper'.")
  779. if kind == 'symmetric':
  780. if exact:
  781. if n > 34:
  782. dt = object
  783. else:
  784. dt = np.int64
  785. else:
  786. dt = np.float64
  787. invp = np.empty((n, n), dtype=dt)
  788. for i in range(n):
  789. for j in range(0, i + 1):
  790. v = 0
  791. for k in range(n - i):
  792. v += comb(i + k, k, exact=exact) * comb(i + k, i + k - j,
  793. exact=exact)
  794. invp[i, j] = (-1)**(i - j) * v
  795. if i != j:
  796. invp[j, i] = invp[i, j]
  797. else:
  798. # For the 'lower' and 'upper' cases, we computer the inverse by
  799. # changing the sign of every other diagonal of the pascal matrix.
  800. invp = pascal(n, kind=kind, exact=exact)
  801. if invp.dtype == np.uint64:
  802. # This cast from np.uint64 to int64 OK, because if `kind` is not
  803. # "symmetric", the values in invp are all much less than 2**63.
  804. invp = invp.view(np.int64)
  805. # The toeplitz matrix has alternating bands of 1 and -1.
  806. invp *= toeplitz((-1)**np.arange(n)).astype(invp.dtype)
  807. return invp
  808. def dft(n, scale=None):
  809. """
  810. Discrete Fourier transform matrix.
  811. Create the matrix that computes the discrete Fourier transform of a
  812. sequence [1]_. The n-th primitive root of unity used to generate the
  813. matrix is exp(-2*pi*i/n), where i = sqrt(-1).
  814. Parameters
  815. ----------
  816. n : int
  817. Size the matrix to create.
  818. scale : str, optional
  819. Must be None, 'sqrtn', or 'n'.
  820. If `scale` is 'sqrtn', the matrix is divided by `sqrt(n)`.
  821. If `scale` is 'n', the matrix is divided by `n`.
  822. If `scale` is None (the default), the matrix is not normalized, and the
  823. return value is simply the Vandermonde matrix of the roots of unity.
  824. Returns
  825. -------
  826. m : (n, n) ndarray
  827. The DFT matrix.
  828. Notes
  829. -----
  830. When `scale` is None, multiplying a vector by the matrix returned by
  831. `dft` is mathematically equivalent to (but much less efficient than)
  832. the calculation performed by `scipy.fftpack.fft`.
  833. .. versionadded:: 0.14.0
  834. References
  835. ----------
  836. .. [1] "DFT matrix", https://en.wikipedia.org/wiki/DFT_matrix
  837. Examples
  838. --------
  839. >>> from scipy.linalg import dft
  840. >>> np.set_printoptions(precision=5, suppress=True)
  841. >>> x = np.array([1, 2, 3, 0, 3, 2, 1, 0])
  842. >>> m = dft(8)
  843. >>> m.dot(x) # Compute the DFT of x
  844. array([ 12.+0.j, -2.-2.j, 0.-4.j, -2.+2.j, 4.+0.j, -2.-2.j,
  845. -0.+4.j, -2.+2.j])
  846. Verify that ``m.dot(x)`` is the same as ``fft(x)``.
  847. >>> from scipy.fftpack import fft
  848. >>> fft(x) # Same result as m.dot(x)
  849. array([ 12.+0.j, -2.-2.j, 0.-4.j, -2.+2.j, 4.+0.j, -2.-2.j,
  850. 0.+4.j, -2.+2.j])
  851. """
  852. if scale not in [None, 'sqrtn', 'n']:
  853. raise ValueError("scale must be None, 'sqrtn', or 'n'; "
  854. "%r is not valid." % (scale,))
  855. omegas = np.exp(-2j * np.pi * np.arange(n) / n).reshape(-1, 1)
  856. m = omegas ** np.arange(n)
  857. if scale == 'sqrtn':
  858. m /= math.sqrt(n)
  859. elif scale == 'n':
  860. m /= n
  861. return m