_special_matrices.py 40 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311
  1. import math
  2. import warnings
  3. import numpy as np
  4. from numpy.lib.stride_tricks import as_strided
  5. from scipy._lib._util import _apply_over_batch
  6. from scipy._lib._array_api import array_namespace, xp_capabilities, xp_size
  7. import scipy._lib.array_api_extra as xpx
  8. __all__ = ['toeplitz', 'circulant', 'hankel',
  9. 'hadamard', 'leslie', 'block_diag', 'companion',
  10. 'helmert', 'hilbert', 'invhilbert', 'pascal', 'invpascal', 'dft',
  11. 'fiedler', 'fiedler_companion', 'convolution_matrix']
  12. # -----------------------------------------------------------------------------
  13. # matrix construction functions
  14. # -----------------------------------------------------------------------------
  15. def toeplitz(c, r=None):
  16. r"""
  17. Construct a Toeplitz matrix.
  18. The Toeplitz matrix has constant diagonals, with c as its first column
  19. and r as its first row. If r is not given, ``r == conjugate(c)`` is
  20. assumed.
  21. The documentation is written assuming array arguments are of specified
  22. "core" shapes. However, array argument(s) of this function may have additional
  23. "batch" dimensions prepended to the core shape. In this case, the array is treated
  24. as a batch of lower-dimensional slices; see :ref:`linalg_batch` for details.
  25. Parameters
  26. ----------
  27. c : array_like
  28. First column of the matrix.
  29. r : array_like, optional
  30. First row of the matrix. If None, ``r = conjugate(c)`` is assumed;
  31. in this case, if c[0] is real, the result is a Hermitian matrix.
  32. r[0] is ignored; the first row of the returned matrix is
  33. ``[c[0], r[1:]]``.
  34. Returns
  35. -------
  36. A : (len(c), len(r)) ndarray
  37. The Toeplitz matrix. Dtype is the same as ``(c[0] + r[0]).dtype``.
  38. See Also
  39. --------
  40. circulant : circulant matrix
  41. hankel : Hankel matrix
  42. solve_toeplitz : Solve a Toeplitz system.
  43. Notes
  44. -----
  45. The behavior when `c` or `r` is a scalar, or when `c` is complex and
  46. `r` is None, was changed in version 0.8.0. The behavior in previous
  47. versions was undocumented and is no longer supported.
  48. Examples
  49. --------
  50. >>> from scipy.linalg import toeplitz
  51. >>> toeplitz([1,2,3], [1,4,5,6])
  52. array([[1, 4, 5, 6],
  53. [2, 1, 4, 5],
  54. [3, 2, 1, 4]])
  55. >>> toeplitz([1.0, 2+3j, 4-1j])
  56. array([[ 1.+0.j, 2.-3.j, 4.+1.j],
  57. [ 2.+3.j, 1.+0.j, 2.-3.j],
  58. [ 4.-1.j, 2.+3.j, 1.+0.j]])
  59. """
  60. c = np.atleast_1d(c)
  61. if r is None:
  62. r = c.conjugate()
  63. else:
  64. r = np.atleast_1d(r)
  65. return _toeplitz(c, r)
  66. @_apply_over_batch(("c", 1), ("r", 1))
  67. def _toeplitz(c, r):
  68. # Form a 1-D array containing a reversed c followed by r[1:] that could be
  69. # strided to give us toeplitz matrix.
  70. vals = np.concatenate((c[::-1], r[1:]))
  71. out_shp = len(c), len(r)
  72. n = vals.strides[0]
  73. return as_strided(vals[len(c)-1:], shape=out_shp, strides=(-n, n)).copy()
  74. def circulant(c):
  75. """
  76. Construct a circulant matrix.
  77. Array argument(s) of this function may have additional
  78. "batch" dimensions prepended to the core shape. In this case, the array is treated
  79. as a batch of lower-dimensional slices; see :ref:`linalg_batch` for details.
  80. Parameters
  81. ----------
  82. c : (..., N,) array_like
  83. The first column(s) of the matrix. Multidimensional arrays are treated as a
  84. batch: each slice along the last axis is the first column of an output matrix.
  85. Returns
  86. -------
  87. A : (..., N, N) ndarray
  88. A circulant matrix whose first column is given by `c`. For batch input, each
  89. slice of shape ``(N, N)`` along the last two dimensions of the output
  90. corresponds with a slice of shape ``(N,)`` along the last dimension of the
  91. input.
  92. See Also
  93. --------
  94. toeplitz : Toeplitz matrix
  95. hankel : Hankel matrix
  96. solve_circulant : Solve a circulant system.
  97. Notes
  98. -----
  99. .. versionadded:: 0.8.0
  100. Examples
  101. --------
  102. >>> from scipy.linalg import circulant
  103. >>> circulant([1, 2, 3])
  104. array([[1, 3, 2],
  105. [2, 1, 3],
  106. [3, 2, 1]])
  107. >>> circulant([[1, 2, 3], [4, 5, 6]])
  108. array([[[1, 3, 2],
  109. [2, 1, 3],
  110. [3, 2, 1]],
  111. [[4, 6, 5],
  112. [5, 4, 6],
  113. [6, 5, 4]]])
  114. """
  115. c = np.atleast_1d(c)
  116. batch_shape, N = c.shape[:-1], c.shape[-1]
  117. # Need to use `prod(batch_shape)` instead of `-1` in case array has zero size
  118. c = c.reshape(math.prod(batch_shape), N) if batch_shape else c
  119. # Form an extended array that could be strided to give circulant version
  120. c_ext = np.concatenate((c[..., ::-1], c[..., :0:-1]), axis=-1).ravel()
  121. L = c.shape[-1]
  122. n = c_ext.strides[-1]
  123. if c.ndim == 1:
  124. A = as_strided(c_ext[L-1:], shape=(L, L), strides=(-n, n))
  125. else:
  126. m = c.shape[0]
  127. A = as_strided(c_ext[L-1:], shape=(m, L, L), strides=(n*(2*L-1), -n, n))
  128. return A.reshape(batch_shape + (N, N)).copy()
  129. def hankel(c, r=None):
  130. r"""
  131. Construct a Hankel matrix.
  132. The Hankel matrix has constant anti-diagonals, with `c` as its
  133. first column and `r` as its last row. If the first element of `r`
  134. differs from the last element of `c`, the first element of `r` is
  135. replaced by the last element of `c` to ensure that anti-diagonals
  136. remain constant. If `r` is not given, then `r = zeros_like(c)` is
  137. assumed.
  138. Parameters
  139. ----------
  140. c : array_like
  141. First column of the matrix.
  142. r : array_like, optional
  143. Last row of the matrix. If None, ``r = zeros_like(c)`` is assumed.
  144. r[0] is ignored; the last row of the returned matrix is
  145. ``[c[-1], r[1:]]``.
  146. .. warning::
  147. Beginning in SciPy 1.19, multidimensional input will be treated as a batch,
  148. not ``ravel``\ ed. To preserve the existing behavior, ``ravel`` arguments
  149. before passing them to `toeplitz`.
  150. Returns
  151. -------
  152. A : (len(c), len(r)) ndarray
  153. The Hankel matrix. Dtype is the same as ``(c[0] + r[0]).dtype``.
  154. See Also
  155. --------
  156. toeplitz : Toeplitz matrix
  157. circulant : circulant matrix
  158. Examples
  159. --------
  160. >>> from scipy.linalg import hankel
  161. >>> hankel([1, 17, 99])
  162. array([[ 1, 17, 99],
  163. [17, 99, 0],
  164. [99, 0, 0]])
  165. >>> hankel([1,2,3,4], [4,7,7,8,9])
  166. array([[1, 2, 3, 4, 7],
  167. [2, 3, 4, 7, 7],
  168. [3, 4, 7, 7, 8],
  169. [4, 7, 7, 8, 9]])
  170. """
  171. c = np.asarray(c)
  172. if r is None:
  173. r = np.zeros_like(c)
  174. else:
  175. r = np.asarray(r)
  176. if c.ndim > 1 or r.ndim > 1:
  177. msg = ("Beginning in SciPy 1.19, multidimensional input will be treated as a "
  178. "batch, not `ravel`ed. To preserve the existing behavior and silence "
  179. "this warning, `ravel` arguments before passing them to `hankel`.")
  180. warnings.warn(msg, FutureWarning, stacklevel=2)
  181. c, r = c.ravel(), r.ravel()
  182. # Form a 1-D array of values to be used in the matrix, containing `c`
  183. # followed by r[1:].
  184. vals = np.concatenate((c, r[1:]))
  185. # Stride on concatenated array to get hankel matrix
  186. out_shp = len(c), len(r)
  187. n = vals.strides[0]
  188. return as_strided(vals, shape=out_shp, strides=(n, n)).copy()
  189. def hadamard(n, dtype=int):
  190. """
  191. Construct an Hadamard matrix.
  192. Constructs an n-by-n Hadamard matrix, using Sylvester's
  193. construction. `n` must be a power of 2.
  194. Parameters
  195. ----------
  196. n : int
  197. The order of the matrix. `n` must be a power of 2.
  198. dtype : dtype, optional
  199. The data type of the array to be constructed.
  200. Returns
  201. -------
  202. H : (n, n) ndarray
  203. The Hadamard matrix.
  204. Notes
  205. -----
  206. .. versionadded:: 0.8.0
  207. Examples
  208. --------
  209. >>> from scipy.linalg import hadamard
  210. >>> hadamard(2, dtype=complex)
  211. array([[ 1.+0.j, 1.+0.j],
  212. [ 1.+0.j, -1.-0.j]])
  213. >>> hadamard(4)
  214. array([[ 1, 1, 1, 1],
  215. [ 1, -1, 1, -1],
  216. [ 1, 1, -1, -1],
  217. [ 1, -1, -1, 1]])
  218. """
  219. # This function is a slightly modified version of the
  220. # function contributed by Ivo in ticket #675.
  221. if n < 1:
  222. lg2 = 0
  223. else:
  224. lg2 = int(math.log(n, 2))
  225. if 2 ** lg2 != n:
  226. raise ValueError("n must be an positive integer, and n must be "
  227. "a power of 2")
  228. H = np.array([[1]], dtype=dtype)
  229. # Sylvester's construction
  230. for i in range(0, lg2):
  231. H = np.vstack((np.hstack((H, H)), np.hstack((H, -H))))
  232. return H
  233. @_apply_over_batch(("f", 1), ("s", 1))
  234. def leslie(f, s):
  235. """
  236. Create a Leslie matrix.
  237. Given the length n array of fecundity coefficients `f` and the length
  238. n-1 array of survival coefficients `s`, return the associated Leslie
  239. matrix.
  240. The documentation is written assuming array arguments are of specified
  241. "core" shapes. However, array argument(s) of this function may have additional
  242. "batch" dimensions prepended to the core shape. In this case, the array is treated
  243. as a batch of lower-dimensional slices; see :ref:`linalg_batch` for details.
  244. Parameters
  245. ----------
  246. f : (N,) array_like
  247. The "fecundity" coefficients.
  248. s : (N-1,) array_like
  249. The "survival" coefficients. The length of `s` must be one less
  250. than the length of `f`, and it must be at least 1.
  251. Returns
  252. -------
  253. L : (N, N) ndarray
  254. The array is zero except for the first row,
  255. which is `f`, and the first sub-diagonal, which is `s`.
  256. The data-type of the array will be the data-type of
  257. ``f[0]+s[0]``.
  258. Notes
  259. -----
  260. The Leslie matrix is used to model discrete-time, age-structured
  261. population growth [1]_ [2]_. In a population with `n` age classes, two sets
  262. of parameters define a Leslie matrix: the `n` "fecundity coefficients",
  263. which give the number of offspring per-capita produced by each age
  264. class, and the `n` - 1 "survival coefficients", which give the
  265. per-capita survival rate of each age class.
  266. References
  267. ----------
  268. .. [1] P. H. Leslie, On the use of matrices in certain population
  269. mathematics, Biometrika, Vol. 33, No. 3, 183--212 (Nov. 1945)
  270. .. [2] P. H. Leslie, Some further notes on the use of matrices in
  271. population mathematics, Biometrika, Vol. 35, No. 3/4, 213--245
  272. (Dec. 1948)
  273. Examples
  274. --------
  275. >>> from scipy.linalg import leslie
  276. >>> leslie([0.1, 2.0, 1.0, 0.1], [0.2, 0.8, 0.7])
  277. array([[ 0.1, 2. , 1. , 0.1],
  278. [ 0.2, 0. , 0. , 0. ],
  279. [ 0. , 0.8, 0. , 0. ],
  280. [ 0. , 0. , 0.7, 0. ]])
  281. """
  282. f = np.atleast_1d(f)
  283. s = np.atleast_1d(s)
  284. if f.shape[-1] != s.shape[-1] + 1:
  285. raise ValueError("Incorrect lengths for f and s. The length of s along "
  286. "the last axis must be one less than the length of f.")
  287. if s.shape[-1] == 0:
  288. raise ValueError("The length of s must be at least 1.")
  289. n = f.shape[-1]
  290. tmp = f[0] + s[0]
  291. a = np.zeros((n, n), dtype=tmp.dtype)
  292. a[0] = f
  293. a[list(range(1, n)), list(range(0, n - 1))] = s
  294. return a
  295. @xp_capabilities(jax_jit=False, allow_dask_compute=2)
  296. def block_diag(*arrs):
  297. """
  298. Create a block diagonal array from provided arrays.
  299. For example, given 2-D inputs `A`, `B` and `C`, the output will have these
  300. arrays arranged on the diagonal::
  301. [[A, 0, 0],
  302. [0, B, 0],
  303. [0, 0, C]]
  304. The documentation is written assuming array arguments are of specified
  305. "core" shapes. However, array argument(s) of this function may have additional
  306. "batch" dimensions prepended to the core shape. In this case, the array is treated
  307. as a batch of lower-dimensional slices; see :ref:`linalg_batch` for details.
  308. Parameters
  309. ----------
  310. A, B, C, ... : array_like
  311. Input arrays. A 1-D array or array_like sequence of length ``n`` is
  312. treated as a 2-D array with shape ``(1, n)``.
  313. Returns
  314. -------
  315. D : ndarray
  316. Array with `A`, `B`, `C`, ... on the diagonal of the last two
  317. dimensions. `D` has the same dtype as the result type of the
  318. inputs.
  319. Notes
  320. -----
  321. If all the input arrays are square, the output is known as a
  322. block diagonal matrix.
  323. Empty sequences (i.e., array-likes of zero size) will not be ignored.
  324. Noteworthy, both ``[]`` and ``[[]]`` are treated as matrices with shape
  325. ``(1,0)``.
  326. Examples
  327. --------
  328. >>> import numpy as np
  329. >>> from scipy.linalg import block_diag
  330. >>> A = [[1, 0],
  331. ... [0, 1]]
  332. >>> B = [[3, 4, 5],
  333. ... [6, 7, 8]]
  334. >>> C = [[7]]
  335. >>> P = np.zeros((2, 0), dtype='int32')
  336. >>> block_diag(A, B, C)
  337. array([[1, 0, 0, 0, 0, 0],
  338. [0, 1, 0, 0, 0, 0],
  339. [0, 0, 3, 4, 5, 0],
  340. [0, 0, 6, 7, 8, 0],
  341. [0, 0, 0, 0, 0, 7]])
  342. >>> block_diag(A, P, B, C)
  343. array([[1, 0, 0, 0, 0, 0],
  344. [0, 1, 0, 0, 0, 0],
  345. [0, 0, 0, 0, 0, 0],
  346. [0, 0, 0, 0, 0, 0],
  347. [0, 0, 3, 4, 5, 0],
  348. [0, 0, 6, 7, 8, 0],
  349. [0, 0, 0, 0, 0, 7]])
  350. >>> block_diag(1.0, [2, 3], [[4, 5], [6, 7]])
  351. array([[ 1., 0., 0., 0., 0.],
  352. [ 0., 2., 3., 0., 0.],
  353. [ 0., 0., 0., 4., 5.],
  354. [ 0., 0., 0., 6., 7.]])
  355. """
  356. xp = array_namespace(*arrs)
  357. if arrs == ():
  358. arrs = ([],)
  359. arrs = [xpx.atleast_nd(xp.asarray(a), ndim=2) for a in arrs]
  360. batch_shapes = [a.shape[:-2] for a in arrs]
  361. batch_shape = np.broadcast_shapes(*batch_shapes)
  362. arrs = [xp.broadcast_to(a, batch_shape + a.shape[-2:]) for a in arrs]
  363. out_dtype = xp.result_type(*arrs)
  364. block_shapes = [a.shape[-2:] for a in arrs]
  365. out = xp.zeros(batch_shape +
  366. tuple(map(int, xp.sum(xp.asarray(block_shapes), axis=0))),
  367. dtype=out_dtype)
  368. r, c = 0, 0
  369. for i, (rr, cc) in enumerate(block_shapes):
  370. out = xpx.at(out)[..., r:r+rr, c:c+cc].set(arrs[i])
  371. r += rr
  372. c += cc
  373. return out
  374. def companion(a):
  375. """
  376. Create a companion matrix.
  377. Create the companion matrix [1]_ associated with the polynomial whose
  378. coefficients are given in `a`.
  379. Array argument(s) of this function may have additional
  380. "batch" dimensions prepended to the core shape. In this case, the array is treated
  381. as a batch of lower-dimensional slices; see :ref:`linalg_batch` for details.
  382. Parameters
  383. ----------
  384. a : (..., N) array_like
  385. 1-D array of polynomial coefficients. The length of `a` must be
  386. at least two, and ``a[0]`` must not be zero.
  387. M-dimensional arrays are treated as a batch: each slice along the last
  388. axis is a 1-D array of polynomial coefficients.
  389. Returns
  390. -------
  391. c : (..., N-1, N-1) ndarray
  392. For 1-D input, the first row of `c` is ``-a[1:]/a[0]``, and the first
  393. sub-diagonal is all ones. The data-type of the array is the same
  394. as the data-type of ``1.0*a[0]``.
  395. For batch input, each slice of shape ``(N-1, N-1)`` along the last two
  396. dimensions of the output corresponds with a slice of shape ``(N,)``
  397. along the last dimension of the input.
  398. Raises
  399. ------
  400. ValueError
  401. If any of the following are true: a) ``a.shape[-1] < 2``; b) ``a[..., 0] == 0``.
  402. Notes
  403. -----
  404. .. versionadded:: 0.8.0
  405. References
  406. ----------
  407. .. [1] R. A. Horn & C. R. Johnson, *Matrix Analysis*. Cambridge, UK:
  408. Cambridge University Press, 1999, pp. 146-7.
  409. Examples
  410. --------
  411. >>> from scipy.linalg import companion
  412. >>> companion([1, -10, 31, -30])
  413. array([[ 10., -31., 30.],
  414. [ 1., 0., 0.],
  415. [ 0., 1., 0.]])
  416. """
  417. a = np.atleast_1d(a)
  418. n = a.shape[-1]
  419. if n < 2:
  420. raise ValueError("The length of `a` along the last axis must be at least 2.")
  421. if np.any(a[..., 0] == 0):
  422. raise ValueError("The first coefficient(s) of `a` (i.e. elements "
  423. "of `a[..., 0]`) must not be zero.")
  424. first_row = -a[..., 1:] / (1.0 * a[..., 0:1])
  425. c = np.zeros(a.shape[:-1] + (n - 1, n - 1), dtype=first_row.dtype)
  426. c[..., 0, :] = first_row
  427. c[..., np.arange(1, n - 1), np.arange(0, n - 2)] = 1
  428. return c
  429. def helmert(n, full=False):
  430. """
  431. Create an Helmert matrix of order `n`.
  432. This has applications in statistics, compositional or simplicial analysis,
  433. and in Aitchison geometry.
  434. Parameters
  435. ----------
  436. n : int
  437. The size of the array to create.
  438. full : bool, optional
  439. If True the (n, n) ndarray will be returned.
  440. Otherwise the submatrix that does not include the first
  441. row will be returned.
  442. Default: False.
  443. Returns
  444. -------
  445. M : ndarray
  446. The Helmert matrix.
  447. The shape is (n, n) or (n-1, n) depending on the `full` argument.
  448. Examples
  449. --------
  450. >>> from scipy.linalg import helmert
  451. >>> helmert(5, full=True)
  452. array([[ 0.4472136 , 0.4472136 , 0.4472136 , 0.4472136 , 0.4472136 ],
  453. [ 0.70710678, -0.70710678, 0. , 0. , 0. ],
  454. [ 0.40824829, 0.40824829, -0.81649658, 0. , 0. ],
  455. [ 0.28867513, 0.28867513, 0.28867513, -0.8660254 , 0. ],
  456. [ 0.2236068 , 0.2236068 , 0.2236068 , 0.2236068 , -0.89442719]])
  457. """
  458. H = np.tril(np.ones((n, n)), -1) - np.diag(np.arange(n))
  459. d = np.arange(n) * np.arange(1, n+1)
  460. H[0] = 1
  461. d[0] = n
  462. H_full = H / np.sqrt(d)[:, np.newaxis]
  463. if full:
  464. return H_full
  465. else:
  466. return H_full[1:]
  467. def hilbert(n):
  468. """
  469. Create a Hilbert matrix of order `n`.
  470. Returns the `n` by `n` array with entries `h[i,j] = 1 / (i + j + 1)`.
  471. Parameters
  472. ----------
  473. n : int
  474. The size of the array to create.
  475. Returns
  476. -------
  477. h : (n, n) ndarray
  478. The Hilbert matrix.
  479. See Also
  480. --------
  481. invhilbert : Compute the inverse of a Hilbert matrix.
  482. Notes
  483. -----
  484. .. versionadded:: 0.10.0
  485. Examples
  486. --------
  487. >>> from scipy.linalg import hilbert
  488. >>> hilbert(3)
  489. array([[ 1. , 0.5 , 0.33333333],
  490. [ 0.5 , 0.33333333, 0.25 ],
  491. [ 0.33333333, 0.25 , 0.2 ]])
  492. """
  493. values = 1.0 / (1.0 + np.arange(2 * n - 1))
  494. h = hankel(values[:n], r=values[n - 1:])
  495. return h
  496. def invhilbert(n, exact=False):
  497. """
  498. Compute the inverse of the Hilbert matrix of order `n`.
  499. The entries in the inverse of a Hilbert matrix are integers. When `n`
  500. is greater than 14, some entries in the inverse exceed the upper limit
  501. of 64 bit integers. The `exact` argument provides two options for
  502. dealing with these large integers.
  503. Parameters
  504. ----------
  505. n : int
  506. The order of the Hilbert matrix.
  507. exact : bool, optional
  508. If False, the data type of the array that is returned is np.float64,
  509. and the array is an approximation of the inverse.
  510. If True, the array is the exact integer inverse array. To represent
  511. the exact inverse when n > 14, the returned array is an object array
  512. of long integers. For n <= 14, the exact inverse is returned as an
  513. array with data type np.int64.
  514. Returns
  515. -------
  516. invh : (n, n) ndarray
  517. The data type of the array is np.float64 if `exact` is False.
  518. If `exact` is True, the data type is either np.int64 (for n <= 14)
  519. or object (for n > 14). In the latter case, the objects in the
  520. array will be long integers.
  521. See Also
  522. --------
  523. hilbert : Create a Hilbert matrix.
  524. Notes
  525. -----
  526. .. versionadded:: 0.10.0
  527. Examples
  528. --------
  529. >>> from scipy.linalg import invhilbert
  530. >>> invhilbert(4)
  531. array([[ 16., -120., 240., -140.],
  532. [ -120., 1200., -2700., 1680.],
  533. [ 240., -2700., 6480., -4200.],
  534. [ -140., 1680., -4200., 2800.]])
  535. >>> invhilbert(4, exact=True)
  536. array([[ 16, -120, 240, -140],
  537. [ -120, 1200, -2700, 1680],
  538. [ 240, -2700, 6480, -4200],
  539. [ -140, 1680, -4200, 2800]], dtype=int64)
  540. >>> invhilbert(16)[7,7]
  541. 4.2475099528537506e+19
  542. >>> invhilbert(16, exact=True)[7,7]
  543. 42475099528537378560
  544. """
  545. from scipy.special import comb
  546. if exact:
  547. if n > 14:
  548. dtype = object
  549. else:
  550. dtype = np.int64
  551. else:
  552. dtype = np.float64
  553. invh = np.empty((n, n), dtype=dtype)
  554. for i in range(n):
  555. for j in range(0, i + 1):
  556. s = i + j
  557. invh[i, j] = ((-1) ** s * (s + 1) *
  558. comb(n + i, n - j - 1, exact=exact) *
  559. comb(n + j, n - i - 1, exact=exact) *
  560. comb(s, i, exact=exact) ** 2)
  561. if i != j:
  562. invh[j, i] = invh[i, j]
  563. return invh
  564. def pascal(n, kind='symmetric', exact=True):
  565. """
  566. Returns the n x n Pascal matrix.
  567. The Pascal matrix is a matrix containing the binomial coefficients as
  568. its elements.
  569. Parameters
  570. ----------
  571. n : int
  572. The size of the matrix to create; that is, the result is an n x n
  573. matrix.
  574. kind : str, optional
  575. Must be one of 'symmetric', 'lower', or 'upper'.
  576. Default is 'symmetric'.
  577. exact : bool, optional
  578. If `exact` is True, the result is either an array of type
  579. numpy.uint64 (if n < 35) or an object array of Python long integers.
  580. If `exact` is False, the coefficients in the matrix are computed using
  581. `scipy.special.comb` with ``exact=False``. The result will be a floating
  582. point array, and the values in the array will not be the exact
  583. coefficients, but this version is much faster than ``exact=True``.
  584. Returns
  585. -------
  586. p : (n, n) ndarray
  587. The Pascal matrix.
  588. See Also
  589. --------
  590. invpascal
  591. Notes
  592. -----
  593. See https://en.wikipedia.org/wiki/Pascal_matrix for more information
  594. about Pascal matrices.
  595. .. versionadded:: 0.11.0
  596. Examples
  597. --------
  598. >>> from scipy.linalg import pascal
  599. >>> pascal(4)
  600. array([[ 1, 1, 1, 1],
  601. [ 1, 2, 3, 4],
  602. [ 1, 3, 6, 10],
  603. [ 1, 4, 10, 20]], dtype=uint64)
  604. >>> pascal(4, kind='lower')
  605. array([[1, 0, 0, 0],
  606. [1, 1, 0, 0],
  607. [1, 2, 1, 0],
  608. [1, 3, 3, 1]], dtype=uint64)
  609. >>> pascal(50)[-1, -1]
  610. 25477612258980856902730428600
  611. >>> from scipy.special import comb
  612. >>> comb(98, 49, exact=True)
  613. 25477612258980856902730428600
  614. """
  615. from scipy.special import comb
  616. if kind not in ['symmetric', 'lower', 'upper']:
  617. raise ValueError("kind must be 'symmetric', 'lower', or 'upper'")
  618. if exact:
  619. if n >= 35:
  620. L_n = np.empty((n, n), dtype=object)
  621. L_n.fill(0)
  622. else:
  623. L_n = np.zeros((n, n), dtype=np.uint64)
  624. for i in range(n):
  625. for j in range(i + 1):
  626. L_n[i, j] = comb(i, j, exact=True)
  627. else:
  628. L_n = comb(*np.ogrid[:n, :n])
  629. if kind == 'lower':
  630. p = L_n
  631. elif kind == 'upper':
  632. p = L_n.T
  633. else:
  634. p = np.dot(L_n, L_n.T)
  635. return p
  636. def invpascal(n, kind='symmetric', exact=True):
  637. """
  638. Returns the inverse of the n x n Pascal matrix.
  639. The Pascal matrix is a matrix containing the binomial coefficients as
  640. its elements.
  641. Parameters
  642. ----------
  643. n : int
  644. The size of the matrix to create; that is, the result is an n x n
  645. matrix.
  646. kind : str, optional
  647. Must be one of 'symmetric', 'lower', or 'upper'.
  648. Default is 'symmetric'.
  649. exact : bool, optional
  650. If `exact` is True, the result is either an array of type
  651. ``numpy.int64`` (if `n` <= 35) or an object array of Python integers.
  652. If `exact` is False, the coefficients in the matrix are computed using
  653. `scipy.special.comb` with `exact=False`. The result will be a floating
  654. point array, and for large `n`, the values in the array will not be the
  655. exact coefficients.
  656. Returns
  657. -------
  658. invp : (n, n) ndarray
  659. The inverse of the Pascal matrix.
  660. See Also
  661. --------
  662. pascal
  663. Notes
  664. -----
  665. .. versionadded:: 0.16.0
  666. References
  667. ----------
  668. .. [1] "Pascal matrix", https://en.wikipedia.org/wiki/Pascal_matrix
  669. .. [2] Cohen, A. M., "The inverse of a Pascal matrix", Mathematical
  670. Gazette, 59(408), pp. 111-112, 1975.
  671. Examples
  672. --------
  673. >>> from scipy.linalg import invpascal, pascal
  674. >>> invp = invpascal(5)
  675. >>> invp
  676. array([[ 5, -10, 10, -5, 1],
  677. [-10, 30, -35, 19, -4],
  678. [ 10, -35, 46, -27, 6],
  679. [ -5, 19, -27, 17, -4],
  680. [ 1, -4, 6, -4, 1]])
  681. >>> p = pascal(5)
  682. >>> p.dot(invp)
  683. array([[ 1., 0., 0., 0., 0.],
  684. [ 0., 1., 0., 0., 0.],
  685. [ 0., 0., 1., 0., 0.],
  686. [ 0., 0., 0., 1., 0.],
  687. [ 0., 0., 0., 0., 1.]])
  688. An example of the use of `kind` and `exact`:
  689. >>> invpascal(5, kind='lower', exact=False)
  690. array([[ 1., -0., 0., -0., 0.],
  691. [-1., 1., -0., 0., -0.],
  692. [ 1., -2., 1., -0., 0.],
  693. [-1., 3., -3., 1., -0.],
  694. [ 1., -4., 6., -4., 1.]])
  695. """
  696. from scipy.special import comb
  697. if kind not in ['symmetric', 'lower', 'upper']:
  698. raise ValueError("'kind' must be 'symmetric', 'lower' or 'upper'.")
  699. if kind == 'symmetric':
  700. if exact:
  701. if n > 34:
  702. dt = object
  703. else:
  704. dt = np.int64
  705. else:
  706. dt = np.float64
  707. invp = np.empty((n, n), dtype=dt)
  708. for i in range(n):
  709. for j in range(0, i + 1):
  710. v = 0
  711. for k in range(n - i):
  712. v += comb(i + k, k, exact=exact) * comb(i + k, i + k - j,
  713. exact=exact)
  714. invp[i, j] = (-1)**(i - j) * v
  715. if i != j:
  716. invp[j, i] = invp[i, j]
  717. else:
  718. # For the 'lower' and 'upper' cases, we computer the inverse by
  719. # changing the sign of every other diagonal of the pascal matrix.
  720. invp = pascal(n, kind=kind, exact=exact)
  721. if invp.dtype == np.uint64:
  722. # This cast from np.uint64 to int64 OK, because if `kind` is not
  723. # "symmetric", the values in invp are all much less than 2**63.
  724. invp = invp.view(np.int64)
  725. # The toeplitz matrix has alternating bands of 1 and -1.
  726. invp *= toeplitz((-1)**np.arange(n)).astype(invp.dtype)
  727. return invp
  728. def dft(n, scale=None):
  729. """
  730. Discrete Fourier transform matrix.
  731. Create the matrix that computes the discrete Fourier transform of a
  732. sequence [1]_. The nth primitive root of unity used to generate the
  733. matrix is exp(-2*pi*i/n), where i = sqrt(-1).
  734. Parameters
  735. ----------
  736. n : int
  737. Size the matrix to create.
  738. scale : str, optional
  739. Must be None, 'sqrtn', or 'n'.
  740. If `scale` is 'sqrtn', the matrix is divided by `sqrt(n)`.
  741. If `scale` is 'n', the matrix is divided by `n`.
  742. If `scale` is None (the default), the matrix is not normalized, and the
  743. return value is simply the Vandermonde matrix of the roots of unity.
  744. Returns
  745. -------
  746. m : (n, n) ndarray
  747. The DFT matrix.
  748. Notes
  749. -----
  750. When `scale` is None, multiplying a vector by the matrix returned by
  751. `dft` is mathematically equivalent to (but much less efficient than)
  752. the calculation performed by `scipy.fft.fft`.
  753. .. versionadded:: 0.14.0
  754. References
  755. ----------
  756. .. [1] "DFT matrix", https://en.wikipedia.org/wiki/DFT_matrix
  757. Examples
  758. --------
  759. >>> import numpy as np
  760. >>> from scipy.linalg import dft
  761. >>> np.set_printoptions(precision=2, suppress=True) # for compact output
  762. >>> m = dft(5)
  763. >>> m
  764. array([[ 1. +0.j , 1. +0.j , 1. +0.j , 1. +0.j , 1. +0.j ],
  765. [ 1. +0.j , 0.31-0.95j, -0.81-0.59j, -0.81+0.59j, 0.31+0.95j],
  766. [ 1. +0.j , -0.81-0.59j, 0.31+0.95j, 0.31-0.95j, -0.81+0.59j],
  767. [ 1. +0.j , -0.81+0.59j, 0.31-0.95j, 0.31+0.95j, -0.81-0.59j],
  768. [ 1. +0.j , 0.31+0.95j, -0.81+0.59j, -0.81-0.59j, 0.31-0.95j]])
  769. >>> x = np.array([1, 2, 3, 0, 3])
  770. >>> m @ x # Compute the DFT of x
  771. array([ 9. +0.j , 0.12-0.81j, -2.12+3.44j, -2.12-3.44j, 0.12+0.81j])
  772. Verify that ``m @ x`` is the same as ``fft(x)``.
  773. >>> from scipy.fft import fft
  774. >>> fft(x) # Same result as m @ x
  775. array([ 9. +0.j , 0.12-0.81j, -2.12+3.44j, -2.12-3.44j, 0.12+0.81j])
  776. """
  777. if scale not in [None, 'sqrtn', 'n']:
  778. raise ValueError("scale must be None, 'sqrtn', or 'n'; "
  779. f"{scale!r} is not valid.")
  780. omegas = np.exp(-2j * np.pi * np.arange(n) / n).reshape(-1, 1)
  781. m = omegas ** np.arange(n)
  782. if scale == 'sqrtn':
  783. m /= math.sqrt(n)
  784. elif scale == 'n':
  785. m /= n
  786. return m
  787. @xp_capabilities()
  788. def fiedler(a):
  789. """Returns a symmetric Fiedler matrix
  790. Given an sequence of numbers `a`, Fiedler matrices have the structure
  791. ``F[i, j] = np.abs(a[i] - a[j])``, and hence zero diagonals and nonnegative
  792. entries. A Fiedler matrix has a dominant positive eigenvalue and other
  793. eigenvalues are negative. Although not valid generally, for certain inputs,
  794. the inverse and the determinant can be derived explicitly as given in [1]_.
  795. Array argument(s) of this function may have additional
  796. "batch" dimensions prepended to the core shape. In this case, the array is treated
  797. as a batch of lower-dimensional slices; see :ref:`linalg_batch` for details.
  798. Parameters
  799. ----------
  800. a : (..., n,) array_like
  801. Coefficient array. N-dimensional arrays are treated as a batch:
  802. each slice along the last axis is a 1-D coefficient array.
  803. Returns
  804. -------
  805. F : (..., n, n) ndarray
  806. Fiedler matrix. For batch input, each slice of shape ``(n, n)``
  807. along the last two dimensions of the output corresponds with a
  808. slice of shape ``(n,)`` along the last dimension of the input.
  809. See Also
  810. --------
  811. circulant, toeplitz
  812. Notes
  813. -----
  814. .. versionadded:: 1.3.0
  815. References
  816. ----------
  817. .. [1] J. Todd, "Basic Numerical Mathematics: Vol.2 : Numerical Algebra",
  818. 1977, Birkhauser, :doi:`10.1007/978-3-0348-7286-7`
  819. Examples
  820. --------
  821. >>> import numpy as np
  822. >>> from scipy.linalg import det, inv, fiedler
  823. >>> a = [1, 4, 12, 45, 77]
  824. >>> n = len(a)
  825. >>> A = fiedler(a)
  826. >>> A
  827. array([[ 0, 3, 11, 44, 76],
  828. [ 3, 0, 8, 41, 73],
  829. [11, 8, 0, 33, 65],
  830. [44, 41, 33, 0, 32],
  831. [76, 73, 65, 32, 0]])
  832. The explicit formulas for determinant and inverse seem to hold only for
  833. monotonically increasing/decreasing arrays. Note the tridiagonal structure
  834. and the corners.
  835. >>> Ai = inv(A)
  836. >>> Ai[np.abs(Ai) < 1e-12] = 0. # cleanup the numerical noise for display
  837. >>> Ai
  838. array([[-0.16008772, 0.16666667, 0. , 0. , 0.00657895],
  839. [ 0.16666667, -0.22916667, 0.0625 , 0. , 0. ],
  840. [ 0. , 0.0625 , -0.07765152, 0.01515152, 0. ],
  841. [ 0. , 0. , 0.01515152, -0.03077652, 0.015625 ],
  842. [ 0.00657895, 0. , 0. , 0.015625 , -0.00904605]])
  843. >>> det(A)
  844. 15409151.999999998
  845. >>> (-1)**(n-1) * 2**(n-2) * np.diff(a).prod() * (a[-1] - a[0])
  846. 15409152
  847. """
  848. xp = array_namespace(a)
  849. a = xpx.atleast_nd(xp.asarray(a), ndim=1)
  850. if xp_size(a) == 0:
  851. return xp.asarray([], dtype=xp.float64)
  852. elif xp_size(a) == 1:
  853. return xp.asarray([[0.]])
  854. else:
  855. return xp.abs(a[..., :, xp.newaxis] - a[..., xp.newaxis, :])
  856. def fiedler_companion(a):
  857. """ Returns a Fiedler companion matrix
  858. Given a polynomial coefficient array ``a``, this function forms a
  859. pentadiagonal matrix with a special structure whose eigenvalues coincides
  860. with the roots of ``a``.
  861. Array argument(s) of this function may have additional
  862. "batch" dimensions prepended to the core shape. In this case, the array is treated
  863. as a batch of lower-dimensional slices; see :ref:`linalg_batch` for details.
  864. Parameters
  865. ----------
  866. a : (..., N) array_like
  867. 1-D array of polynomial coefficients in descending order with a nonzero
  868. leading coefficient. For ``N < 2``, an empty array is returned.
  869. N-dimensional arrays are treated as a batch: each slice along the last
  870. axis is a 1-D array of polynomial coefficients.
  871. Returns
  872. -------
  873. c : (..., N-1, N-1) ndarray
  874. Resulting companion matrix. For batch input, each slice of shape
  875. ``(N-1, N-1)`` along the last two dimensions of the output corresponds
  876. with a slice of shape ``(N,)`` along the last dimension of the input.
  877. See Also
  878. --------
  879. companion
  880. Notes
  881. -----
  882. Similar to `companion`, each leading coefficient along the last axis of the
  883. input should be nonzero.
  884. If the leading coefficient is not 1, other coefficients are rescaled before
  885. the array generation. To avoid numerical issues, it is best to provide a
  886. monic polynomial.
  887. .. versionadded:: 1.3.0
  888. References
  889. ----------
  890. .. [1] M. Fiedler, " A note on companion matrices", Linear Algebra and its
  891. Applications, 2003, :doi:`10.1016/S0024-3795(03)00548-2`
  892. Examples
  893. --------
  894. >>> import numpy as np
  895. >>> from scipy.linalg import fiedler_companion, eigvals
  896. >>> p = np.poly(np.arange(1, 9, 2)) # [1., -16., 86., -176., 105.]
  897. >>> fc = fiedler_companion(p)
  898. >>> fc
  899. array([[ 16., -86., 1., 0.],
  900. [ 1., 0., 0., 0.],
  901. [ 0., 176., 0., -105.],
  902. [ 0., 1., 0., 0.]])
  903. >>> eigvals(fc)
  904. array([7.+0.j, 5.+0.j, 3.+0.j, 1.+0.j])
  905. """
  906. a = np.atleast_1d(a)
  907. if a.ndim > 1:
  908. return np.apply_along_axis(fiedler_companion, -1, a)
  909. if a.size <= 2:
  910. if a.size == 2:
  911. return np.array([[-(a/a[0])[-1]]])
  912. return np.array([], dtype=a.dtype)
  913. if a[0] == 0.:
  914. raise ValueError('Leading coefficient is zero.')
  915. a = a/a[0]
  916. n = a.size - 1
  917. c = np.zeros((n, n), dtype=a.dtype)
  918. # subdiagonals
  919. c[range(3, n, 2), range(1, n-2, 2)] = 1.
  920. c[range(2, n, 2), range(1, n-1, 2)] = -a[3::2]
  921. # superdiagonals
  922. c[range(0, n-2, 2), range(2, n, 2)] = 1.
  923. c[range(0, n-1, 2), range(1, n, 2)] = -a[2::2]
  924. c[[0, 1], 0] = [-a[1], 1]
  925. return c
  926. def convolution_matrix(a, n, mode='full'):
  927. """
  928. Construct a convolution matrix.
  929. Constructs the Toeplitz matrix representing one-dimensional
  930. convolution [1]_. See the notes below for details.
  931. Array argument(s) of this function may have additional
  932. "batch" dimensions prepended to the core shape. In this case, the array is treated
  933. as a batch of lower-dimensional slices; see :ref:`linalg_batch` for details.
  934. Parameters
  935. ----------
  936. a : (..., m) array_like
  937. The 1-D array to convolve. N-dimensional arrays are treated as a
  938. batch: each slice along the last axis is a 1-D array to convolve.
  939. n : int
  940. The number of columns in the resulting matrix. It gives the length
  941. of the input to be convolved with `a`. This is analogous to the
  942. length of `v` in ``numpy.convolve(a, v)``.
  943. mode : str
  944. This is analogous to `mode` in ``numpy.convolve(v, a, mode)``.
  945. It must be one of ('full', 'valid', 'same').
  946. See below for how `mode` determines the shape of the result.
  947. Returns
  948. -------
  949. A : (..., k, n) ndarray
  950. The convolution matrix whose row count `k` depends on `mode`::
  951. ======= =========================
  952. mode k
  953. ======= =========================
  954. 'full' m + n -1
  955. 'same' max(m, n)
  956. 'valid' max(m, n) - min(m, n) + 1
  957. ======= =========================
  958. For batch input, each slice of shape ``(k, n)`` along the last two
  959. dimensions of the output corresponds with a slice of shape ``(m,)``
  960. along the last dimension of the input.
  961. See Also
  962. --------
  963. toeplitz : Toeplitz matrix
  964. Notes
  965. -----
  966. The code::
  967. A = convolution_matrix(a, n, mode)
  968. creates a Toeplitz matrix `A` such that ``A @ v`` is equivalent to
  969. using ``convolve(a, v, mode)``. The returned array always has `n`
  970. columns. The number of rows depends on the specified `mode`, as
  971. explained above.
  972. In the default 'full' mode, the entries of `A` are given by::
  973. A[i, j] == (a[i-j] if (0 <= (i-j) < m) else 0)
  974. where ``m = len(a)``. Suppose, for example, the input array is
  975. ``[x, y, z]``. The convolution matrix has the form::
  976. [x, 0, 0, ..., 0, 0]
  977. [y, x, 0, ..., 0, 0]
  978. [z, y, x, ..., 0, 0]
  979. ...
  980. [0, 0, 0, ..., x, 0]
  981. [0, 0, 0, ..., y, x]
  982. [0, 0, 0, ..., z, y]
  983. [0, 0, 0, ..., 0, z]
  984. In 'valid' mode, the entries of `A` are given by::
  985. A[i, j] == (a[i-j+m-1] if (0 <= (i-j+m-1) < m) else 0)
  986. This corresponds to a matrix whose rows are the subset of those from
  987. the 'full' case where all the coefficients in `a` are contained in the
  988. row. For input ``[x, y, z]``, this array looks like::
  989. [z, y, x, 0, 0, ..., 0, 0, 0]
  990. [0, z, y, x, 0, ..., 0, 0, 0]
  991. [0, 0, z, y, x, ..., 0, 0, 0]
  992. ...
  993. [0, 0, 0, 0, 0, ..., x, 0, 0]
  994. [0, 0, 0, 0, 0, ..., y, x, 0]
  995. [0, 0, 0, 0, 0, ..., z, y, x]
  996. In the 'same' mode, the entries of `A` are given by::
  997. d = (m - 1) // 2
  998. A[i, j] == (a[i-j+d] if (0 <= (i-j+d) < m) else 0)
  999. The typical application of the 'same' mode is when one has a signal of
  1000. length `n` (with `n` greater than ``len(a)``), and the desired output
  1001. is a filtered signal that is still of length `n`.
  1002. For input ``[x, y, z]``, this array looks like::
  1003. [y, x, 0, 0, ..., 0, 0, 0]
  1004. [z, y, x, 0, ..., 0, 0, 0]
  1005. [0, z, y, x, ..., 0, 0, 0]
  1006. [0, 0, z, y, ..., 0, 0, 0]
  1007. ...
  1008. [0, 0, 0, 0, ..., y, x, 0]
  1009. [0, 0, 0, 0, ..., z, y, x]
  1010. [0, 0, 0, 0, ..., 0, z, y]
  1011. .. versionadded:: 1.5.0
  1012. References
  1013. ----------
  1014. .. [1] "Convolution", https://en.wikipedia.org/wiki/Convolution
  1015. Examples
  1016. --------
  1017. >>> import numpy as np
  1018. >>> from scipy.linalg import convolution_matrix
  1019. >>> A = convolution_matrix([-1, 4, -2], 5, mode='same')
  1020. >>> A
  1021. array([[ 4, -1, 0, 0, 0],
  1022. [-2, 4, -1, 0, 0],
  1023. [ 0, -2, 4, -1, 0],
  1024. [ 0, 0, -2, 4, -1],
  1025. [ 0, 0, 0, -2, 4]])
  1026. Compare multiplication by `A` with the use of `numpy.convolve`.
  1027. >>> x = np.array([1, 2, 0, -3, 0.5])
  1028. >>> A @ x
  1029. array([ 2. , 6. , -1. , -12.5, 8. ])
  1030. Verify that ``A @ x`` produced the same result as applying the
  1031. convolution function.
  1032. >>> np.convolve([-1, 4, -2], x, mode='same')
  1033. array([ 2. , 6. , -1. , -12.5, 8. ])
  1034. For comparison to the case ``mode='same'`` shown above, here are the
  1035. matrices produced by ``mode='full'`` and ``mode='valid'`` for the
  1036. same coefficients and size.
  1037. >>> convolution_matrix([-1, 4, -2], 5, mode='full')
  1038. array([[-1, 0, 0, 0, 0],
  1039. [ 4, -1, 0, 0, 0],
  1040. [-2, 4, -1, 0, 0],
  1041. [ 0, -2, 4, -1, 0],
  1042. [ 0, 0, -2, 4, -1],
  1043. [ 0, 0, 0, -2, 4],
  1044. [ 0, 0, 0, 0, -2]])
  1045. >>> convolution_matrix([-1, 4, -2], 5, mode='valid')
  1046. array([[-2, 4, -1, 0, 0],
  1047. [ 0, -2, 4, -1, 0],
  1048. [ 0, 0, -2, 4, -1]])
  1049. """
  1050. if n <= 0:
  1051. raise ValueError('n must be a positive integer.')
  1052. a = np.asarray(a)
  1053. if a.size == 0:
  1054. raise ValueError('len(a) must be at least 1.')
  1055. if mode not in ('full', 'valid', 'same'):
  1056. raise ValueError(
  1057. "'mode' argument must be one of ('full', 'valid', 'same')")
  1058. if a.ndim > 1:
  1059. return np.apply_along_axis(lambda a: convolution_matrix(a, n, mode), -1, a)
  1060. # create zero padded versions of the array
  1061. az = np.pad(a, (0, n-1), 'constant')
  1062. raz = np.pad(a[::-1], (0, n-1), 'constant')
  1063. if mode == 'same':
  1064. trim = min(n, len(a)) - 1
  1065. tb = trim//2
  1066. te = trim - tb
  1067. col0 = az[tb:len(az)-te]
  1068. row0 = raz[-n-tb:len(raz)-tb]
  1069. elif mode == 'valid':
  1070. tb = min(n, len(a)) - 1
  1071. te = tb
  1072. col0 = az[tb:len(az)-te]
  1073. row0 = raz[-n-tb:len(raz)-tb]
  1074. else: # 'full'
  1075. col0 = az
  1076. row0 = raz[-n:]
  1077. return toeplitz(col0, row0)