_pocketfft.py 47 KB

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