fftpack.py 46 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303
  1. """
  2. Discrete Fourier Transforms
  3. Routines in this module:
  4. fft(a, n=None, axis=-1)
  5. ifft(a, n=None, axis=-1)
  6. rfft(a, n=None, axis=-1)
  7. irfft(a, n=None, axis=-1)
  8. hfft(a, n=None, axis=-1)
  9. ihfft(a, n=None, axis=-1)
  10. fftn(a, s=None, axes=None)
  11. ifftn(a, s=None, axes=None)
  12. rfftn(a, s=None, axes=None)
  13. irfftn(a, s=None, axes=None)
  14. fft2(a, s=None, axes=(-2,-1))
  15. ifft2(a, s=None, axes=(-2, -1))
  16. rfft2(a, s=None, axes=(-2,-1))
  17. irfft2(a, s=None, axes=(-2, -1))
  18. i = inverse transform
  19. r = transform of purely real data
  20. h = Hermite transform
  21. n = n-dimensional transform
  22. 2 = 2-dimensional transform
  23. (Note: 2D routines are just nD routines with different default
  24. behavior.)
  25. The underlying code for these functions is an f2c-translated and modified
  26. version of the FFTPACK routines.
  27. """
  28. from __future__ import division, absolute_import, print_function
  29. __all__ = ['fft', 'ifft', 'rfft', 'irfft', 'hfft', 'ihfft', 'rfftn',
  30. 'irfftn', 'rfft2', 'irfft2', 'fft2', 'ifft2', 'fftn', 'ifftn']
  31. import functools
  32. from numpy.core import (array, asarray, zeros, swapaxes, shape, conjugate,
  33. take, sqrt)
  34. from numpy.core.multiarray import normalize_axis_index
  35. from numpy.core import overrides
  36. from . import fftpack_lite as fftpack
  37. from .helper import _FFTCache
  38. _fft_cache = _FFTCache(max_size_in_mb=100, max_item_count=32)
  39. _real_fft_cache = _FFTCache(max_size_in_mb=100, max_item_count=32)
  40. array_function_dispatch = functools.partial(
  41. overrides.array_function_dispatch, module='numpy.fft')
  42. def _raw_fft(a, n=None, axis=-1, init_function=fftpack.cffti,
  43. work_function=fftpack.cfftf, fft_cache=_fft_cache):
  44. a = asarray(a)
  45. axis = normalize_axis_index(axis, a.ndim)
  46. if n is None:
  47. n = a.shape[axis]
  48. if n < 1:
  49. raise ValueError("Invalid number of FFT data points (%d) specified."
  50. % n)
  51. # We have to ensure that only a single thread can access a wsave array
  52. # at any given time. Thus we remove it from the cache and insert it
  53. # again after it has been used. Multiple threads might create multiple
  54. # copies of the wsave array. This is intentional and a limitation of
  55. # the current C code.
  56. wsave = fft_cache.pop_twiddle_factors(n)
  57. if wsave is None:
  58. wsave = init_function(n)
  59. if a.shape[axis] != n:
  60. s = list(a.shape)
  61. if s[axis] > n:
  62. index = [slice(None)]*len(s)
  63. index[axis] = slice(0, n)
  64. a = a[tuple(index)]
  65. else:
  66. index = [slice(None)]*len(s)
  67. index[axis] = slice(0, s[axis])
  68. s[axis] = n
  69. z = zeros(s, a.dtype.char)
  70. z[tuple(index)] = a
  71. a = z
  72. if axis != a.ndim - 1:
  73. a = swapaxes(a, axis, -1)
  74. r = work_function(a, wsave)
  75. if axis != a.ndim - 1:
  76. r = swapaxes(r, axis, -1)
  77. # As soon as we put wsave back into the cache, another thread could pick it
  78. # up and start using it, so we must not do this until after we're
  79. # completely done using it ourselves.
  80. fft_cache.put_twiddle_factors(n, wsave)
  81. return r
  82. def _unitary(norm):
  83. if norm not in (None, "ortho"):
  84. raise ValueError("Invalid norm value %s, should be None or \"ortho\"."
  85. % norm)
  86. return norm is not None
  87. def _fft_dispatcher(a, n=None, axis=None, norm=None):
  88. return (a,)
  89. @array_function_dispatch(_fft_dispatcher)
  90. def fft(a, n=None, axis=-1, norm=None):
  91. """
  92. Compute the one-dimensional discrete Fourier Transform.
  93. This function computes the one-dimensional *n*-point discrete Fourier
  94. Transform (DFT) with the efficient Fast Fourier Transform (FFT)
  95. algorithm [CT].
  96. Parameters
  97. ----------
  98. a : array_like
  99. Input array, can be complex.
  100. n : int, optional
  101. Length of the transformed axis of the output.
  102. If `n` is smaller than the length of the input, the input is cropped.
  103. If it is larger, the input is padded with zeros. If `n` is not given,
  104. the length of the input along the axis specified by `axis` is used.
  105. axis : int, optional
  106. Axis over which to compute the FFT. If not given, the last axis is
  107. used.
  108. norm : {None, "ortho"}, optional
  109. .. versionadded:: 1.10.0
  110. Normalization mode (see `numpy.fft`). Default is None.
  111. Returns
  112. -------
  113. out : complex ndarray
  114. The truncated or zero-padded input, transformed along the axis
  115. indicated by `axis`, or the last one if `axis` is not specified.
  116. Raises
  117. ------
  118. IndexError
  119. if `axes` is larger than the last axis of `a`.
  120. See Also
  121. --------
  122. numpy.fft : for definition of the DFT and conventions used.
  123. ifft : The inverse of `fft`.
  124. fft2 : The two-dimensional FFT.
  125. fftn : The *n*-dimensional FFT.
  126. rfftn : The *n*-dimensional FFT of real input.
  127. fftfreq : Frequency bins for given FFT parameters.
  128. Notes
  129. -----
  130. FFT (Fast Fourier Transform) refers to a way the discrete Fourier
  131. Transform (DFT) can be calculated efficiently, by using symmetries in the
  132. calculated terms. The symmetry is highest when `n` is a power of 2, and
  133. the transform is therefore most efficient for these sizes.
  134. The DFT is defined, with the conventions used in this implementation, in
  135. the documentation for the `numpy.fft` module.
  136. References
  137. ----------
  138. .. [CT] Cooley, James W., and John W. Tukey, 1965, "An algorithm for the
  139. machine calculation of complex Fourier series," *Math. Comput.*
  140. 19: 297-301.
  141. Examples
  142. --------
  143. >>> np.fft.fft(np.exp(2j * np.pi * np.arange(8) / 8))
  144. array([ -3.44505240e-16 +1.14383329e-17j,
  145. 8.00000000e+00 -5.71092652e-15j,
  146. 2.33482938e-16 +1.22460635e-16j,
  147. 1.64863782e-15 +1.77635684e-15j,
  148. 9.95839695e-17 +2.33482938e-16j,
  149. 0.00000000e+00 +1.66837030e-15j,
  150. 1.14383329e-17 +1.22460635e-16j,
  151. -1.64863782e-15 +1.77635684e-15j])
  152. In this example, real input has an FFT which is Hermitian, i.e., symmetric
  153. in the real part and anti-symmetric in the imaginary part, as described in
  154. the `numpy.fft` documentation:
  155. >>> import matplotlib.pyplot as plt
  156. >>> t = np.arange(256)
  157. >>> sp = np.fft.fft(np.sin(t))
  158. >>> freq = np.fft.fftfreq(t.shape[-1])
  159. >>> plt.plot(freq, sp.real, freq, sp.imag)
  160. [<matplotlib.lines.Line2D object at 0x...>, <matplotlib.lines.Line2D object at 0x...>]
  161. >>> plt.show()
  162. """
  163. a = asarray(a).astype(complex, copy=False)
  164. if n is None:
  165. n = a.shape[axis]
  166. output = _raw_fft(a, n, axis, fftpack.cffti, fftpack.cfftf, _fft_cache)
  167. if _unitary(norm):
  168. output *= 1 / sqrt(n)
  169. return output
  170. @array_function_dispatch(_fft_dispatcher)
  171. def ifft(a, n=None, axis=-1, norm=None):
  172. """
  173. Compute the one-dimensional inverse discrete Fourier Transform.
  174. This function computes the inverse of the one-dimensional *n*-point
  175. discrete Fourier transform computed by `fft`. In other words,
  176. ``ifft(fft(a)) == a`` to within numerical accuracy.
  177. For a general description of the algorithm and definitions,
  178. see `numpy.fft`.
  179. The input should be ordered in the same way as is returned by `fft`,
  180. i.e.,
  181. * ``a[0]`` should contain the zero frequency term,
  182. * ``a[1:n//2]`` should contain the positive-frequency terms,
  183. * ``a[n//2 + 1:]`` should contain the negative-frequency terms, in
  184. increasing order starting from the most negative frequency.
  185. For an even number of input points, ``A[n//2]`` represents the sum of
  186. the values at the positive and negative Nyquist frequencies, as the two
  187. are aliased together. See `numpy.fft` for details.
  188. Parameters
  189. ----------
  190. a : array_like
  191. Input array, can be complex.
  192. n : int, optional
  193. Length of the transformed axis of the output.
  194. If `n` is smaller than the length of the input, the input is cropped.
  195. If it is larger, the input is padded with zeros. If `n` is not given,
  196. the length of the input along the axis specified by `axis` is used.
  197. See notes about padding issues.
  198. axis : int, optional
  199. Axis over which to compute the inverse DFT. If not given, the last
  200. axis is used.
  201. norm : {None, "ortho"}, optional
  202. .. versionadded:: 1.10.0
  203. Normalization mode (see `numpy.fft`). Default is None.
  204. Returns
  205. -------
  206. out : complex ndarray
  207. The truncated or zero-padded input, transformed along the axis
  208. indicated by `axis`, or the last one if `axis` is not specified.
  209. Raises
  210. ------
  211. IndexError
  212. If `axes` is larger than the last axis of `a`.
  213. See Also
  214. --------
  215. numpy.fft : An introduction, with definitions and general explanations.
  216. fft : The one-dimensional (forward) FFT, of which `ifft` is the inverse
  217. ifft2 : The two-dimensional inverse FFT.
  218. ifftn : The n-dimensional inverse FFT.
  219. Notes
  220. -----
  221. If the input parameter `n` is larger than the size of the input, the input
  222. is padded by appending zeros at the end. Even though this is the common
  223. approach, it might lead to surprising results. If a different padding is
  224. desired, it must be performed before calling `ifft`.
  225. Examples
  226. --------
  227. >>> np.fft.ifft([0, 4, 0, 0])
  228. array([ 1.+0.j, 0.+1.j, -1.+0.j, 0.-1.j])
  229. Create and plot a band-limited signal with random phases:
  230. >>> import matplotlib.pyplot as plt
  231. >>> t = np.arange(400)
  232. >>> n = np.zeros((400,), dtype=complex)
  233. >>> n[40:60] = np.exp(1j*np.random.uniform(0, 2*np.pi, (20,)))
  234. >>> s = np.fft.ifft(n)
  235. >>> plt.plot(t, s.real, 'b-', t, s.imag, 'r--')
  236. ...
  237. >>> plt.legend(('real', 'imaginary'))
  238. ...
  239. >>> plt.show()
  240. """
  241. # The copy may be required for multithreading.
  242. a = array(a, copy=True, dtype=complex)
  243. if n is None:
  244. n = a.shape[axis]
  245. unitary = _unitary(norm)
  246. output = _raw_fft(a, n, axis, fftpack.cffti, fftpack.cfftb, _fft_cache)
  247. return output * (1 / (sqrt(n) if unitary else n))
  248. @array_function_dispatch(_fft_dispatcher)
  249. def rfft(a, n=None, axis=-1, norm=None):
  250. """
  251. Compute the one-dimensional discrete Fourier Transform for real input.
  252. This function computes the one-dimensional *n*-point discrete Fourier
  253. Transform (DFT) of a real-valued array by means of an efficient algorithm
  254. called the Fast Fourier Transform (FFT).
  255. Parameters
  256. ----------
  257. a : array_like
  258. Input array
  259. n : int, optional
  260. Number of points along transformation axis in the input to use.
  261. If `n` is smaller than the length of the input, the input is cropped.
  262. If it is larger, the input is padded with zeros. If `n` is not given,
  263. the length of the input along the axis specified by `axis` is used.
  264. axis : int, optional
  265. Axis over which to compute the FFT. If not given, the last axis is
  266. used.
  267. norm : {None, "ortho"}, optional
  268. .. versionadded:: 1.10.0
  269. Normalization mode (see `numpy.fft`). Default is None.
  270. Returns
  271. -------
  272. out : complex ndarray
  273. The truncated or zero-padded input, transformed along the axis
  274. indicated by `axis`, or the last one if `axis` is not specified.
  275. If `n` is even, the length of the transformed axis is ``(n/2)+1``.
  276. If `n` is odd, the length is ``(n+1)/2``.
  277. Raises
  278. ------
  279. IndexError
  280. If `axis` is larger than the last axis of `a`.
  281. See Also
  282. --------
  283. numpy.fft : For definition of the DFT and conventions used.
  284. irfft : The inverse of `rfft`.
  285. fft : The one-dimensional FFT of general (complex) input.
  286. fftn : The *n*-dimensional FFT.
  287. rfftn : The *n*-dimensional FFT of real input.
  288. Notes
  289. -----
  290. When the DFT is computed for purely real input, the output is
  291. Hermitian-symmetric, i.e. the negative frequency terms are just the complex
  292. conjugates of the corresponding positive-frequency terms, and the
  293. negative-frequency terms are therefore redundant. This function does not
  294. compute the negative frequency terms, and the length of the transformed
  295. axis of the output is therefore ``n//2 + 1``.
  296. When ``A = rfft(a)`` and fs is the sampling frequency, ``A[0]`` contains
  297. the zero-frequency term 0*fs, which is real due to Hermitian symmetry.
  298. If `n` is even, ``A[-1]`` contains the term representing both positive
  299. and negative Nyquist frequency (+fs/2 and -fs/2), and must also be purely
  300. real. If `n` is odd, there is no term at fs/2; ``A[-1]`` contains
  301. the largest positive frequency (fs/2*(n-1)/n), and is complex in the
  302. general case.
  303. If the input `a` contains an imaginary part, it is silently discarded.
  304. Examples
  305. --------
  306. >>> np.fft.fft([0, 1, 0, 0])
  307. array([ 1.+0.j, 0.-1.j, -1.+0.j, 0.+1.j])
  308. >>> np.fft.rfft([0, 1, 0, 0])
  309. array([ 1.+0.j, 0.-1.j, -1.+0.j])
  310. Notice how the final element of the `fft` output is the complex conjugate
  311. of the second element, for real input. For `rfft`, this symmetry is
  312. exploited to compute only the non-negative frequency terms.
  313. """
  314. # The copy may be required for multithreading.
  315. a = array(a, copy=True, dtype=float)
  316. output = _raw_fft(a, n, axis, fftpack.rffti, fftpack.rfftf,
  317. _real_fft_cache)
  318. if _unitary(norm):
  319. if n is None:
  320. n = a.shape[axis]
  321. output *= 1 / sqrt(n)
  322. return output
  323. @array_function_dispatch(_fft_dispatcher)
  324. def irfft(a, n=None, axis=-1, norm=None):
  325. """
  326. Compute the inverse of the n-point DFT for real input.
  327. This function computes the inverse of the one-dimensional *n*-point
  328. discrete Fourier Transform of real input computed by `rfft`.
  329. In other words, ``irfft(rfft(a), len(a)) == a`` to within numerical
  330. accuracy. (See Notes below for why ``len(a)`` is necessary here.)
  331. The input is expected to be in the form returned by `rfft`, i.e. the
  332. real zero-frequency term followed by the complex positive frequency terms
  333. in order of increasing frequency. Since the discrete Fourier Transform of
  334. real input is Hermitian-symmetric, the negative frequency terms are taken
  335. to be the complex conjugates of the corresponding positive frequency terms.
  336. Parameters
  337. ----------
  338. a : array_like
  339. The input array.
  340. n : int, optional
  341. Length of the transformed axis of the output.
  342. For `n` output points, ``n//2+1`` input points are necessary. If the
  343. input is longer than this, it is cropped. If it is shorter than this,
  344. it is padded with zeros. If `n` is not given, it is determined from
  345. the length of the input along the axis specified by `axis`.
  346. axis : int, optional
  347. Axis over which to compute the inverse FFT. If not given, the last
  348. axis is used.
  349. norm : {None, "ortho"}, optional
  350. .. versionadded:: 1.10.0
  351. Normalization mode (see `numpy.fft`). Default is None.
  352. Returns
  353. -------
  354. out : ndarray
  355. The truncated or zero-padded input, transformed along the axis
  356. indicated by `axis`, or the last one if `axis` is not specified.
  357. The length of the transformed axis is `n`, or, if `n` is not given,
  358. ``2*(m-1)`` where ``m`` is the length of the transformed axis of the
  359. input. To get an odd number of output points, `n` must be specified.
  360. Raises
  361. ------
  362. IndexError
  363. If `axis` is larger than the last axis of `a`.
  364. See Also
  365. --------
  366. numpy.fft : For definition of the DFT and conventions used.
  367. rfft : The one-dimensional FFT of real input, of which `irfft` is inverse.
  368. fft : The one-dimensional FFT.
  369. irfft2 : The inverse of the two-dimensional FFT of real input.
  370. irfftn : The inverse of the *n*-dimensional FFT of real input.
  371. Notes
  372. -----
  373. Returns the real valued `n`-point inverse discrete Fourier transform
  374. of `a`, where `a` contains the non-negative frequency terms of a
  375. Hermitian-symmetric sequence. `n` is the length of the result, not the
  376. input.
  377. If you specify an `n` such that `a` must be zero-padded or truncated, the
  378. extra/removed values will be added/removed at high frequencies. One can
  379. thus resample a series to `m` points via Fourier interpolation by:
  380. ``a_resamp = irfft(rfft(a), m)``.
  381. Examples
  382. --------
  383. >>> np.fft.ifft([1, -1j, -1, 1j])
  384. array([ 0.+0.j, 1.+0.j, 0.+0.j, 0.+0.j])
  385. >>> np.fft.irfft([1, -1j, -1])
  386. array([ 0., 1., 0., 0.])
  387. Notice how the last term in the input to the ordinary `ifft` is the
  388. complex conjugate of the second term, and the output has zero imaginary
  389. part everywhere. When calling `irfft`, the negative frequencies are not
  390. specified, and the output array is purely real.
  391. """
  392. # The copy may be required for multithreading.
  393. a = array(a, copy=True, dtype=complex)
  394. if n is None:
  395. n = (a.shape[axis] - 1) * 2
  396. unitary = _unitary(norm)
  397. output = _raw_fft(a, n, axis, fftpack.rffti, fftpack.rfftb,
  398. _real_fft_cache)
  399. return output * (1 / (sqrt(n) if unitary else n))
  400. @array_function_dispatch(_fft_dispatcher)
  401. def hfft(a, n=None, axis=-1, norm=None):
  402. """
  403. Compute the FFT of a signal that has Hermitian symmetry, i.e., a real
  404. spectrum.
  405. Parameters
  406. ----------
  407. a : array_like
  408. The input array.
  409. n : int, optional
  410. Length of the transformed axis of the output. For `n` output
  411. points, ``n//2 + 1`` input points are necessary. If the input is
  412. longer than this, it is cropped. If it is shorter than this, it is
  413. padded with zeros. If `n` is not given, it is determined from the
  414. length of the input along the axis specified by `axis`.
  415. axis : int, optional
  416. Axis over which to compute the FFT. If not given, the last
  417. axis is used.
  418. norm : {None, "ortho"}, optional
  419. Normalization mode (see `numpy.fft`). Default is None.
  420. .. versionadded:: 1.10.0
  421. Returns
  422. -------
  423. out : ndarray
  424. The truncated or zero-padded input, transformed along the axis
  425. indicated by `axis`, or the last one if `axis` is not specified.
  426. The length of the transformed axis is `n`, or, if `n` is not given,
  427. ``2*m - 2`` where ``m`` is the length of the transformed axis of
  428. the input. To get an odd number of output points, `n` must be
  429. specified, for instance as ``2*m - 1`` in the typical case,
  430. Raises
  431. ------
  432. IndexError
  433. If `axis` is larger than the last axis of `a`.
  434. See also
  435. --------
  436. rfft : Compute the one-dimensional FFT for real input.
  437. ihfft : The inverse of `hfft`.
  438. Notes
  439. -----
  440. `hfft`/`ihfft` are a pair analogous to `rfft`/`irfft`, but for the
  441. opposite case: here the signal has Hermitian symmetry in the time
  442. domain and is real in the frequency domain. So here it's `hfft` for
  443. which you must supply the length of the result if it is to be odd.
  444. * even: ``ihfft(hfft(a, 2*len(a) - 2) == a``, within roundoff error,
  445. * odd: ``ihfft(hfft(a, 2*len(a) - 1) == a``, within roundoff error.
  446. Examples
  447. --------
  448. >>> signal = np.array([1, 2, 3, 4, 3, 2])
  449. >>> np.fft.fft(signal)
  450. array([ 15.+0.j, -4.+0.j, 0.+0.j, -1.-0.j, 0.+0.j, -4.+0.j])
  451. >>> np.fft.hfft(signal[:4]) # Input first half of signal
  452. array([ 15., -4., 0., -1., 0., -4.])
  453. >>> np.fft.hfft(signal, 6) # Input entire signal and truncate
  454. array([ 15., -4., 0., -1., 0., -4.])
  455. >>> signal = np.array([[1, 1.j], [-1.j, 2]])
  456. >>> np.conj(signal.T) - signal # check Hermitian symmetry
  457. array([[ 0.-0.j, 0.+0.j],
  458. [ 0.+0.j, 0.-0.j]])
  459. >>> freq_spectrum = np.fft.hfft(signal)
  460. >>> freq_spectrum
  461. array([[ 1., 1.],
  462. [ 2., -2.]])
  463. """
  464. # The copy may be required for multithreading.
  465. a = array(a, copy=True, dtype=complex)
  466. if n is None:
  467. n = (a.shape[axis] - 1) * 2
  468. unitary = _unitary(norm)
  469. return irfft(conjugate(a), n, axis) * (sqrt(n) if unitary else n)
  470. @array_function_dispatch(_fft_dispatcher)
  471. def ihfft(a, n=None, axis=-1, norm=None):
  472. """
  473. Compute the inverse FFT of a signal that has Hermitian symmetry.
  474. Parameters
  475. ----------
  476. a : array_like
  477. Input array.
  478. n : int, optional
  479. Length of the inverse FFT, the number of points along
  480. transformation axis in the input to use. If `n` is smaller than
  481. the length of the input, the input is cropped. If it is larger,
  482. the input is padded with zeros. If `n` is not given, the length of
  483. the input along the axis specified by `axis` is used.
  484. axis : int, optional
  485. Axis over which to compute the inverse FFT. If not given, the last
  486. axis is used.
  487. norm : {None, "ortho"}, optional
  488. Normalization mode (see `numpy.fft`). Default is None.
  489. .. versionadded:: 1.10.0
  490. Returns
  491. -------
  492. out : complex ndarray
  493. The truncated or zero-padded input, transformed along the axis
  494. indicated by `axis`, or the last one if `axis` is not specified.
  495. The length of the transformed axis is ``n//2 + 1``.
  496. See also
  497. --------
  498. hfft, irfft
  499. Notes
  500. -----
  501. `hfft`/`ihfft` are a pair analogous to `rfft`/`irfft`, but for the
  502. opposite case: here the signal has Hermitian symmetry in the time
  503. domain and is real in the frequency domain. So here it's `hfft` for
  504. which you must supply the length of the result if it is to be odd:
  505. * even: ``ihfft(hfft(a, 2*len(a) - 2) == a``, within roundoff error,
  506. * odd: ``ihfft(hfft(a, 2*len(a) - 1) == a``, within roundoff error.
  507. Examples
  508. --------
  509. >>> spectrum = np.array([ 15, -4, 0, -1, 0, -4])
  510. >>> np.fft.ifft(spectrum)
  511. array([ 1.+0.j, 2.-0.j, 3.+0.j, 4.+0.j, 3.+0.j, 2.-0.j])
  512. >>> np.fft.ihfft(spectrum)
  513. array([ 1.-0.j, 2.-0.j, 3.-0.j, 4.-0.j])
  514. """
  515. # The copy may be required for multithreading.
  516. a = array(a, copy=True, dtype=float)
  517. if n is None:
  518. n = a.shape[axis]
  519. unitary = _unitary(norm)
  520. output = conjugate(rfft(a, n, axis))
  521. return output * (1 / (sqrt(n) if unitary else n))
  522. def _cook_nd_args(a, s=None, axes=None, invreal=0):
  523. if s is None:
  524. shapeless = 1
  525. if axes is None:
  526. s = list(a.shape)
  527. else:
  528. s = take(a.shape, axes)
  529. else:
  530. shapeless = 0
  531. s = list(s)
  532. if axes is None:
  533. axes = list(range(-len(s), 0))
  534. if len(s) != len(axes):
  535. raise ValueError("Shape and axes have different lengths.")
  536. if invreal and shapeless:
  537. s[-1] = (a.shape[axes[-1]] - 1) * 2
  538. return s, axes
  539. def _raw_fftnd(a, s=None, axes=None, function=fft, norm=None):
  540. a = asarray(a)
  541. s, axes = _cook_nd_args(a, s, axes)
  542. itl = list(range(len(axes)))
  543. itl.reverse()
  544. for ii in itl:
  545. a = function(a, n=s[ii], axis=axes[ii], norm=norm)
  546. return a
  547. def _fftn_dispatcher(a, s=None, axes=None, norm=None):
  548. return (a,)
  549. @array_function_dispatch(_fftn_dispatcher)
  550. def fftn(a, s=None, axes=None, norm=None):
  551. """
  552. Compute the N-dimensional discrete Fourier Transform.
  553. This function computes the *N*-dimensional discrete Fourier Transform over
  554. any number of axes in an *M*-dimensional array by means of the Fast Fourier
  555. Transform (FFT).
  556. Parameters
  557. ----------
  558. a : array_like
  559. Input array, can be complex.
  560. s : sequence of ints, optional
  561. Shape (length of each transformed axis) of the output
  562. (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.).
  563. This corresponds to ``n`` for ``fft(x, n)``.
  564. Along any axis, if the given shape is smaller than that of the input,
  565. the input is cropped. If it is larger, the input is padded with zeros.
  566. if `s` is not given, the shape of the input along the axes specified
  567. by `axes` is used.
  568. axes : sequence of ints, optional
  569. Axes over which to compute the FFT. If not given, the last ``len(s)``
  570. axes are used, or all axes if `s` is also not specified.
  571. Repeated indices in `axes` means that the transform over that axis is
  572. performed multiple times.
  573. norm : {None, "ortho"}, optional
  574. .. versionadded:: 1.10.0
  575. Normalization mode (see `numpy.fft`). Default is None.
  576. Returns
  577. -------
  578. out : complex ndarray
  579. The truncated or zero-padded input, transformed along the axes
  580. indicated by `axes`, or by a combination of `s` and `a`,
  581. as explained in the parameters section above.
  582. Raises
  583. ------
  584. ValueError
  585. If `s` and `axes` have different length.
  586. IndexError
  587. If an element of `axes` is larger than than the number of axes of `a`.
  588. See Also
  589. --------
  590. numpy.fft : Overall view of discrete Fourier transforms, with definitions
  591. and conventions used.
  592. ifftn : The inverse of `fftn`, the inverse *n*-dimensional FFT.
  593. fft : The one-dimensional FFT, with definitions and conventions used.
  594. rfftn : The *n*-dimensional FFT of real input.
  595. fft2 : The two-dimensional FFT.
  596. fftshift : Shifts zero-frequency terms to centre of array
  597. Notes
  598. -----
  599. The output, analogously to `fft`, contains the term for zero frequency in
  600. the low-order corner of all axes, the positive frequency terms in the
  601. first half of all axes, the term for the Nyquist frequency in the middle
  602. of all axes and the negative frequency terms in the second half of all
  603. axes, in order of decreasingly negative frequency.
  604. See `numpy.fft` for details, definitions and conventions used.
  605. Examples
  606. --------
  607. >>> a = np.mgrid[:3, :3, :3][0]
  608. >>> np.fft.fftn(a, axes=(1, 2))
  609. array([[[ 0.+0.j, 0.+0.j, 0.+0.j],
  610. [ 0.+0.j, 0.+0.j, 0.+0.j],
  611. [ 0.+0.j, 0.+0.j, 0.+0.j]],
  612. [[ 9.+0.j, 0.+0.j, 0.+0.j],
  613. [ 0.+0.j, 0.+0.j, 0.+0.j],
  614. [ 0.+0.j, 0.+0.j, 0.+0.j]],
  615. [[ 18.+0.j, 0.+0.j, 0.+0.j],
  616. [ 0.+0.j, 0.+0.j, 0.+0.j],
  617. [ 0.+0.j, 0.+0.j, 0.+0.j]]])
  618. >>> np.fft.fftn(a, (2, 2), axes=(0, 1))
  619. array([[[ 2.+0.j, 2.+0.j, 2.+0.j],
  620. [ 0.+0.j, 0.+0.j, 0.+0.j]],
  621. [[-2.+0.j, -2.+0.j, -2.+0.j],
  622. [ 0.+0.j, 0.+0.j, 0.+0.j]]])
  623. >>> import matplotlib.pyplot as plt
  624. >>> [X, Y] = np.meshgrid(2 * np.pi * np.arange(200) / 12,
  625. ... 2 * np.pi * np.arange(200) / 34)
  626. >>> S = np.sin(X) + np.cos(Y) + np.random.uniform(0, 1, X.shape)
  627. >>> FS = np.fft.fftn(S)
  628. >>> plt.imshow(np.log(np.abs(np.fft.fftshift(FS))**2))
  629. <matplotlib.image.AxesImage object at 0x...>
  630. >>> plt.show()
  631. """
  632. return _raw_fftnd(a, s, axes, fft, norm)
  633. @array_function_dispatch(_fftn_dispatcher)
  634. def ifftn(a, s=None, axes=None, norm=None):
  635. """
  636. Compute the N-dimensional inverse discrete Fourier Transform.
  637. This function computes the inverse of the N-dimensional discrete
  638. Fourier Transform over any number of axes in an M-dimensional array by
  639. means of the Fast Fourier Transform (FFT). In other words,
  640. ``ifftn(fftn(a)) == a`` to within numerical accuracy.
  641. For a description of the definitions and conventions used, see `numpy.fft`.
  642. The input, analogously to `ifft`, should be ordered in the same way as is
  643. returned by `fftn`, i.e. it should have the term for zero frequency
  644. in all axes in the low-order corner, the positive frequency terms in the
  645. first half of all axes, the term for the Nyquist frequency in the middle
  646. of all axes and the negative frequency terms in the second half of all
  647. axes, in order of decreasingly negative frequency.
  648. Parameters
  649. ----------
  650. a : array_like
  651. Input array, can be complex.
  652. s : sequence of ints, optional
  653. Shape (length of each transformed axis) of the output
  654. (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.).
  655. This corresponds to ``n`` for ``ifft(x, n)``.
  656. Along any axis, if the given shape is smaller than that of the input,
  657. the input is cropped. If it is larger, the input is padded with zeros.
  658. if `s` is not given, the shape of the input along the axes specified
  659. by `axes` is used. See notes for issue on `ifft` zero padding.
  660. axes : sequence of ints, optional
  661. Axes over which to compute the IFFT. If not given, the last ``len(s)``
  662. axes are used, or all axes if `s` is also not specified.
  663. Repeated indices in `axes` means that the inverse transform over that
  664. axis is performed multiple times.
  665. norm : {None, "ortho"}, optional
  666. .. versionadded:: 1.10.0
  667. Normalization mode (see `numpy.fft`). Default is None.
  668. Returns
  669. -------
  670. out : complex ndarray
  671. The truncated or zero-padded input, transformed along the axes
  672. indicated by `axes`, or by a combination of `s` or `a`,
  673. as explained in the parameters section above.
  674. Raises
  675. ------
  676. ValueError
  677. If `s` and `axes` have different length.
  678. IndexError
  679. If an element of `axes` is larger than than the number of axes of `a`.
  680. See Also
  681. --------
  682. numpy.fft : Overall view of discrete Fourier transforms, with definitions
  683. and conventions used.
  684. fftn : The forward *n*-dimensional FFT, of which `ifftn` is the inverse.
  685. ifft : The one-dimensional inverse FFT.
  686. ifft2 : The two-dimensional inverse FFT.
  687. ifftshift : Undoes `fftshift`, shifts zero-frequency terms to beginning
  688. of array.
  689. Notes
  690. -----
  691. See `numpy.fft` for definitions and conventions used.
  692. Zero-padding, analogously with `ifft`, is performed by appending zeros to
  693. the input along the specified dimension. Although this is the common
  694. approach, it might lead to surprising results. If another form of zero
  695. padding is desired, it must be performed before `ifftn` is called.
  696. Examples
  697. --------
  698. >>> a = np.eye(4)
  699. >>> np.fft.ifftn(np.fft.fftn(a, axes=(0,)), axes=(1,))
  700. array([[ 1.+0.j, 0.+0.j, 0.+0.j, 0.+0.j],
  701. [ 0.+0.j, 1.+0.j, 0.+0.j, 0.+0.j],
  702. [ 0.+0.j, 0.+0.j, 1.+0.j, 0.+0.j],
  703. [ 0.+0.j, 0.+0.j, 0.+0.j, 1.+0.j]])
  704. Create and plot an image with band-limited frequency content:
  705. >>> import matplotlib.pyplot as plt
  706. >>> n = np.zeros((200,200), dtype=complex)
  707. >>> n[60:80, 20:40] = np.exp(1j*np.random.uniform(0, 2*np.pi, (20, 20)))
  708. >>> im = np.fft.ifftn(n).real
  709. >>> plt.imshow(im)
  710. <matplotlib.image.AxesImage object at 0x...>
  711. >>> plt.show()
  712. """
  713. return _raw_fftnd(a, s, axes, ifft, norm)
  714. @array_function_dispatch(_fftn_dispatcher)
  715. def fft2(a, s=None, axes=(-2, -1), norm=None):
  716. """
  717. Compute the 2-dimensional discrete Fourier Transform
  718. This function computes the *n*-dimensional discrete Fourier Transform
  719. over any axes in an *M*-dimensional array by means of the
  720. Fast Fourier Transform (FFT). By default, the transform is computed over
  721. the last two axes of the input array, i.e., a 2-dimensional FFT.
  722. Parameters
  723. ----------
  724. a : array_like
  725. Input array, can be complex
  726. s : sequence of ints, optional
  727. Shape (length of each transformed axis) of the output
  728. (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.).
  729. This corresponds to ``n`` for ``fft(x, n)``.
  730. Along each axis, if the given shape is smaller than that of the input,
  731. the input is cropped. If it is larger, the input is padded with zeros.
  732. if `s` is not given, the shape of the input along the axes specified
  733. by `axes` is used.
  734. axes : sequence of ints, optional
  735. Axes over which to compute the FFT. If not given, the last two
  736. axes are used. A repeated index in `axes` means the transform over
  737. that axis is performed multiple times. A one-element sequence means
  738. that a one-dimensional FFT is performed.
  739. norm : {None, "ortho"}, optional
  740. .. versionadded:: 1.10.0
  741. Normalization mode (see `numpy.fft`). Default is None.
  742. Returns
  743. -------
  744. out : complex ndarray
  745. The truncated or zero-padded input, transformed along the axes
  746. indicated by `axes`, or the last two axes if `axes` is not given.
  747. Raises
  748. ------
  749. ValueError
  750. If `s` and `axes` have different length, or `axes` not given and
  751. ``len(s) != 2``.
  752. IndexError
  753. If an element of `axes` is larger than than the number of axes of `a`.
  754. See Also
  755. --------
  756. numpy.fft : Overall view of discrete Fourier transforms, with definitions
  757. and conventions used.
  758. ifft2 : The inverse two-dimensional FFT.
  759. fft : The one-dimensional FFT.
  760. fftn : The *n*-dimensional FFT.
  761. fftshift : Shifts zero-frequency terms to the center of the array.
  762. For two-dimensional input, swaps first and third quadrants, and second
  763. and fourth quadrants.
  764. Notes
  765. -----
  766. `fft2` is just `fftn` with a different default for `axes`.
  767. The output, analogously to `fft`, contains the term for zero frequency in
  768. the low-order corner of the transformed axes, the positive frequency terms
  769. in the first half of these axes, the term for the Nyquist frequency in the
  770. middle of the axes and the negative frequency terms in the second half of
  771. the axes, in order of decreasingly negative frequency.
  772. See `fftn` for details and a plotting example, and `numpy.fft` for
  773. definitions and conventions used.
  774. Examples
  775. --------
  776. >>> a = np.mgrid[:5, :5][0]
  777. >>> np.fft.fft2(a)
  778. array([[ 50.0 +0.j , 0.0 +0.j , 0.0 +0.j ,
  779. 0.0 +0.j , 0.0 +0.j ],
  780. [-12.5+17.20477401j, 0.0 +0.j , 0.0 +0.j ,
  781. 0.0 +0.j , 0.0 +0.j ],
  782. [-12.5 +4.0614962j , 0.0 +0.j , 0.0 +0.j ,
  783. 0.0 +0.j , 0.0 +0.j ],
  784. [-12.5 -4.0614962j , 0.0 +0.j , 0.0 +0.j ,
  785. 0.0 +0.j , 0.0 +0.j ],
  786. [-12.5-17.20477401j, 0.0 +0.j , 0.0 +0.j ,
  787. 0.0 +0.j , 0.0 +0.j ]])
  788. """
  789. return _raw_fftnd(a, s, axes, fft, norm)
  790. @array_function_dispatch(_fftn_dispatcher)
  791. def ifft2(a, s=None, axes=(-2, -1), norm=None):
  792. """
  793. Compute the 2-dimensional inverse discrete Fourier Transform.
  794. This function computes the inverse of the 2-dimensional discrete Fourier
  795. Transform over any number of axes in an M-dimensional array by means of
  796. the Fast Fourier Transform (FFT). In other words, ``ifft2(fft2(a)) == a``
  797. to within numerical accuracy. By default, the inverse transform is
  798. computed over the last two axes of the input array.
  799. The input, analogously to `ifft`, should be ordered in the same way as is
  800. returned by `fft2`, i.e. it should have the term for zero frequency
  801. in the low-order corner of the two axes, the positive frequency terms in
  802. the first half of these axes, the term for the Nyquist frequency in the
  803. middle of the axes and the negative frequency terms in the second half of
  804. both axes, in order of decreasingly negative frequency.
  805. Parameters
  806. ----------
  807. a : array_like
  808. Input array, can be complex.
  809. s : sequence of ints, optional
  810. Shape (length of each axis) of the output (``s[0]`` refers to axis 0,
  811. ``s[1]`` to axis 1, etc.). This corresponds to `n` for ``ifft(x, n)``.
  812. Along each axis, if the given shape is smaller than that of the input,
  813. the input is cropped. If it is larger, the input is padded with zeros.
  814. if `s` is not given, the shape of the input along the axes specified
  815. by `axes` is used. See notes for issue on `ifft` zero padding.
  816. axes : sequence of ints, optional
  817. Axes over which to compute the FFT. If not given, the last two
  818. axes are used. A repeated index in `axes` means the transform over
  819. that axis is performed multiple times. A one-element sequence means
  820. that a one-dimensional FFT is performed.
  821. norm : {None, "ortho"}, optional
  822. .. versionadded:: 1.10.0
  823. Normalization mode (see `numpy.fft`). Default is None.
  824. Returns
  825. -------
  826. out : complex ndarray
  827. The truncated or zero-padded input, transformed along the axes
  828. indicated by `axes`, or the last two axes if `axes` is not given.
  829. Raises
  830. ------
  831. ValueError
  832. If `s` and `axes` have different length, or `axes` not given and
  833. ``len(s) != 2``.
  834. IndexError
  835. If an element of `axes` is larger than than the number of axes of `a`.
  836. See Also
  837. --------
  838. numpy.fft : Overall view of discrete Fourier transforms, with definitions
  839. and conventions used.
  840. fft2 : The forward 2-dimensional FFT, of which `ifft2` is the inverse.
  841. ifftn : The inverse of the *n*-dimensional FFT.
  842. fft : The one-dimensional FFT.
  843. ifft : The one-dimensional inverse FFT.
  844. Notes
  845. -----
  846. `ifft2` is just `ifftn` with a different default for `axes`.
  847. See `ifftn` for details and a plotting example, and `numpy.fft` for
  848. definition and conventions used.
  849. Zero-padding, analogously with `ifft`, is performed by appending zeros to
  850. the input along the specified dimension. Although this is the common
  851. approach, it might lead to surprising results. If another form of zero
  852. padding is desired, it must be performed before `ifft2` is called.
  853. Examples
  854. --------
  855. >>> a = 4 * np.eye(4)
  856. >>> np.fft.ifft2(a)
  857. array([[ 1.+0.j, 0.+0.j, 0.+0.j, 0.+0.j],
  858. [ 0.+0.j, 0.+0.j, 0.+0.j, 1.+0.j],
  859. [ 0.+0.j, 0.+0.j, 1.+0.j, 0.+0.j],
  860. [ 0.+0.j, 1.+0.j, 0.+0.j, 0.+0.j]])
  861. """
  862. return _raw_fftnd(a, s, axes, ifft, norm)
  863. @array_function_dispatch(_fftn_dispatcher)
  864. def rfftn(a, s=None, axes=None, norm=None):
  865. """
  866. Compute the N-dimensional discrete Fourier Transform for real input.
  867. This function computes the N-dimensional discrete Fourier Transform over
  868. any number of axes in an M-dimensional real array by means of the Fast
  869. Fourier Transform (FFT). By default, all axes are transformed, with the
  870. real transform performed over the last axis, while the remaining
  871. transforms are complex.
  872. Parameters
  873. ----------
  874. a : array_like
  875. Input array, taken to be real.
  876. s : sequence of ints, optional
  877. Shape (length along each transformed axis) to use from the input.
  878. (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.).
  879. The final element of `s` corresponds to `n` for ``rfft(x, n)``, while
  880. for the remaining axes, it corresponds to `n` for ``fft(x, n)``.
  881. Along any axis, if the given shape is smaller than that of the input,
  882. the input is cropped. If it is larger, the input is padded with zeros.
  883. if `s` is not given, the shape of the input along the axes specified
  884. by `axes` is used.
  885. axes : sequence of ints, optional
  886. Axes over which to compute the FFT. If not given, the last ``len(s)``
  887. axes are used, or all axes if `s` is also not specified.
  888. norm : {None, "ortho"}, optional
  889. .. versionadded:: 1.10.0
  890. Normalization mode (see `numpy.fft`). Default is None.
  891. Returns
  892. -------
  893. out : complex ndarray
  894. The truncated or zero-padded input, transformed along the axes
  895. indicated by `axes`, or by a combination of `s` and `a`,
  896. as explained in the parameters section above.
  897. The length of the last axis transformed will be ``s[-1]//2+1``,
  898. while the remaining transformed axes will have lengths according to
  899. `s`, or unchanged from the input.
  900. Raises
  901. ------
  902. ValueError
  903. If `s` and `axes` have different length.
  904. IndexError
  905. If an element of `axes` is larger than than the number of axes of `a`.
  906. See Also
  907. --------
  908. irfftn : The inverse of `rfftn`, i.e. the inverse of the n-dimensional FFT
  909. of real input.
  910. fft : The one-dimensional FFT, with definitions and conventions used.
  911. rfft : The one-dimensional FFT of real input.
  912. fftn : The n-dimensional FFT.
  913. rfft2 : The two-dimensional FFT of real input.
  914. Notes
  915. -----
  916. The transform for real input is performed over the last transformation
  917. axis, as by `rfft`, then the transform over the remaining axes is
  918. performed as by `fftn`. The order of the output is as for `rfft` for the
  919. final transformation axis, and as for `fftn` for the remaining
  920. transformation axes.
  921. See `fft` for details, definitions and conventions used.
  922. Examples
  923. --------
  924. >>> a = np.ones((2, 2, 2))
  925. >>> np.fft.rfftn(a)
  926. array([[[ 8.+0.j, 0.+0.j],
  927. [ 0.+0.j, 0.+0.j]],
  928. [[ 0.+0.j, 0.+0.j],
  929. [ 0.+0.j, 0.+0.j]]])
  930. >>> np.fft.rfftn(a, axes=(2, 0))
  931. array([[[ 4.+0.j, 0.+0.j],
  932. [ 4.+0.j, 0.+0.j]],
  933. [[ 0.+0.j, 0.+0.j],
  934. [ 0.+0.j, 0.+0.j]]])
  935. """
  936. # The copy may be required for multithreading.
  937. a = array(a, copy=True, dtype=float)
  938. s, axes = _cook_nd_args(a, s, axes)
  939. a = rfft(a, s[-1], axes[-1], norm)
  940. for ii in range(len(axes)-1):
  941. a = fft(a, s[ii], axes[ii], norm)
  942. return a
  943. @array_function_dispatch(_fftn_dispatcher)
  944. def rfft2(a, s=None, axes=(-2, -1), norm=None):
  945. """
  946. Compute the 2-dimensional FFT of a real array.
  947. Parameters
  948. ----------
  949. a : array
  950. Input array, taken to be real.
  951. s : sequence of ints, optional
  952. Shape of the FFT.
  953. axes : sequence of ints, optional
  954. Axes over which to compute the FFT.
  955. norm : {None, "ortho"}, optional
  956. .. versionadded:: 1.10.0
  957. Normalization mode (see `numpy.fft`). Default is None.
  958. Returns
  959. -------
  960. out : ndarray
  961. The result of the real 2-D FFT.
  962. See Also
  963. --------
  964. rfftn : Compute the N-dimensional discrete Fourier Transform for real
  965. input.
  966. Notes
  967. -----
  968. This is really just `rfftn` with different default behavior.
  969. For more details see `rfftn`.
  970. """
  971. return rfftn(a, s, axes, norm)
  972. @array_function_dispatch(_fftn_dispatcher)
  973. def irfftn(a, s=None, axes=None, norm=None):
  974. """
  975. Compute the inverse of the N-dimensional FFT of real input.
  976. This function computes the inverse of the N-dimensional discrete
  977. Fourier Transform for real input over any number of axes in an
  978. M-dimensional array by means of the Fast Fourier Transform (FFT). In
  979. other words, ``irfftn(rfftn(a), a.shape) == a`` to within numerical
  980. accuracy. (The ``a.shape`` is necessary like ``len(a)`` is for `irfft`,
  981. and for the same reason.)
  982. The input should be ordered in the same way as is returned by `rfftn`,
  983. i.e. as for `irfft` for the final transformation axis, and as for `ifftn`
  984. along all the other axes.
  985. Parameters
  986. ----------
  987. a : array_like
  988. Input array.
  989. s : sequence of ints, optional
  990. Shape (length of each transformed axis) of the output
  991. (``s[0]`` refers to axis 0, ``s[1]`` to axis 1, etc.). `s` is also the
  992. number of input points used along this axis, except for the last axis,
  993. where ``s[-1]//2+1`` points of the input are used.
  994. Along any axis, if the shape indicated by `s` is smaller than that of
  995. the input, the input is cropped. If it is larger, the input is padded
  996. with zeros. If `s` is not given, the shape of the input along the
  997. axes specified by `axes` is used.
  998. axes : sequence of ints, optional
  999. Axes over which to compute the inverse FFT. If not given, the last
  1000. `len(s)` axes are used, or all axes if `s` is also not specified.
  1001. Repeated indices in `axes` means that the inverse transform over that
  1002. axis is performed multiple times.
  1003. norm : {None, "ortho"}, optional
  1004. .. versionadded:: 1.10.0
  1005. Normalization mode (see `numpy.fft`). Default is None.
  1006. Returns
  1007. -------
  1008. out : ndarray
  1009. The truncated or zero-padded input, transformed along the axes
  1010. indicated by `axes`, or by a combination of `s` or `a`,
  1011. as explained in the parameters section above.
  1012. The length of each transformed axis is as given by the corresponding
  1013. element of `s`, or the length of the input in every axis except for the
  1014. last one if `s` is not given. In the final transformed axis the length
  1015. of the output when `s` is not given is ``2*(m-1)`` where ``m`` is the
  1016. length of the final transformed axis of the input. To get an odd
  1017. number of output points in the final axis, `s` must be specified.
  1018. Raises
  1019. ------
  1020. ValueError
  1021. If `s` and `axes` have different length.
  1022. IndexError
  1023. If an element of `axes` is larger than than the number of axes of `a`.
  1024. See Also
  1025. --------
  1026. rfftn : The forward n-dimensional FFT of real input,
  1027. of which `ifftn` is the inverse.
  1028. fft : The one-dimensional FFT, with definitions and conventions used.
  1029. irfft : The inverse of the one-dimensional FFT of real input.
  1030. irfft2 : The inverse of the two-dimensional FFT of real input.
  1031. Notes
  1032. -----
  1033. See `fft` for definitions and conventions used.
  1034. See `rfft` for definitions and conventions used for real input.
  1035. Examples
  1036. --------
  1037. >>> a = np.zeros((3, 2, 2))
  1038. >>> a[0, 0, 0] = 3 * 2 * 2
  1039. >>> np.fft.irfftn(a)
  1040. array([[[ 1., 1.],
  1041. [ 1., 1.]],
  1042. [[ 1., 1.],
  1043. [ 1., 1.]],
  1044. [[ 1., 1.],
  1045. [ 1., 1.]]])
  1046. """
  1047. # The copy may be required for multithreading.
  1048. a = array(a, copy=True, dtype=complex)
  1049. s, axes = _cook_nd_args(a, s, axes, invreal=1)
  1050. for ii in range(len(axes)-1):
  1051. a = ifft(a, s[ii], axes[ii], norm)
  1052. a = irfft(a, s[-1], axes[-1], norm)
  1053. return a
  1054. @array_function_dispatch(_fftn_dispatcher)
  1055. def irfft2(a, s=None, axes=(-2, -1), norm=None):
  1056. """
  1057. Compute the 2-dimensional inverse FFT of a real array.
  1058. Parameters
  1059. ----------
  1060. a : array_like
  1061. The input array
  1062. s : sequence of ints, optional
  1063. Shape of the inverse FFT.
  1064. axes : sequence of ints, optional
  1065. The axes over which to compute the inverse fft.
  1066. Default is the last two axes.
  1067. norm : {None, "ortho"}, optional
  1068. .. versionadded:: 1.10.0
  1069. Normalization mode (see `numpy.fft`). Default is None.
  1070. Returns
  1071. -------
  1072. out : ndarray
  1073. The result of the inverse real 2-D FFT.
  1074. See Also
  1075. --------
  1076. irfftn : Compute the inverse of the N-dimensional FFT of real input.
  1077. Notes
  1078. -----
  1079. This is really `irfftn` with different defaults.
  1080. For more details see `irfftn`.
  1081. """
  1082. return irfftn(a, s, axes, norm)