helper.py 6.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224
  1. """
  2. Discrete Fourier Transforms - helper.py
  3. """
  4. from __future__ import division, absolute_import, print_function
  5. from numpy.compat import integer_types
  6. from numpy.core import (
  7. asarray, concatenate, arange, take, integer, empty
  8. )
  9. # Created by Pearu Peterson, September 2002
  10. __all__ = ['fftshift', 'ifftshift', 'fftfreq', 'rfftfreq']
  11. integer_types = integer_types + (integer,)
  12. def fftshift(x, axes=None):
  13. """
  14. Shift the zero-frequency component to the center of the spectrum.
  15. This function swaps half-spaces for all axes listed (defaults to all).
  16. Note that ``y[0]`` is the Nyquist component only if ``len(x)`` is even.
  17. Parameters
  18. ----------
  19. x : array_like
  20. Input array.
  21. axes : int or shape tuple, optional
  22. Axes over which to shift. Default is None, which shifts all axes.
  23. Returns
  24. -------
  25. y : ndarray
  26. The shifted array.
  27. See Also
  28. --------
  29. ifftshift : The inverse of `fftshift`.
  30. Examples
  31. --------
  32. >>> freqs = np.fft.fftfreq(10, 0.1)
  33. >>> freqs
  34. array([ 0., 1., 2., 3., 4., -5., -4., -3., -2., -1.])
  35. >>> np.fft.fftshift(freqs)
  36. array([-5., -4., -3., -2., -1., 0., 1., 2., 3., 4.])
  37. Shift the zero-frequency component only along the second axis:
  38. >>> freqs = np.fft.fftfreq(9, d=1./9).reshape(3, 3)
  39. >>> freqs
  40. array([[ 0., 1., 2.],
  41. [ 3., 4., -4.],
  42. [-3., -2., -1.]])
  43. >>> np.fft.fftshift(freqs, axes=(1,))
  44. array([[ 2., 0., 1.],
  45. [-4., 3., 4.],
  46. [-1., -3., -2.]])
  47. """
  48. tmp = asarray(x)
  49. ndim = len(tmp.shape)
  50. if axes is None:
  51. axes = list(range(ndim))
  52. elif isinstance(axes, integer_types):
  53. axes = (axes,)
  54. y = tmp
  55. for k in axes:
  56. n = tmp.shape[k]
  57. p2 = (n+1)//2
  58. mylist = concatenate((arange(p2, n), arange(p2)))
  59. y = take(y, mylist, k)
  60. return y
  61. def ifftshift(x, axes=None):
  62. """
  63. The inverse of `fftshift`. Although identical for even-length `x`, the
  64. functions differ by one sample for odd-length `x`.
  65. Parameters
  66. ----------
  67. x : array_like
  68. Input array.
  69. axes : int or shape tuple, optional
  70. Axes over which to calculate. Defaults to None, which shifts all axes.
  71. Returns
  72. -------
  73. y : ndarray
  74. The shifted array.
  75. See Also
  76. --------
  77. fftshift : Shift zero-frequency component to the center of the spectrum.
  78. Examples
  79. --------
  80. >>> freqs = np.fft.fftfreq(9, d=1./9).reshape(3, 3)
  81. >>> freqs
  82. array([[ 0., 1., 2.],
  83. [ 3., 4., -4.],
  84. [-3., -2., -1.]])
  85. >>> np.fft.ifftshift(np.fft.fftshift(freqs))
  86. array([[ 0., 1., 2.],
  87. [ 3., 4., -4.],
  88. [-3., -2., -1.]])
  89. """
  90. tmp = asarray(x)
  91. ndim = len(tmp.shape)
  92. if axes is None:
  93. axes = list(range(ndim))
  94. elif isinstance(axes, integer_types):
  95. axes = (axes,)
  96. y = tmp
  97. for k in axes:
  98. n = tmp.shape[k]
  99. p2 = n-(n+1)//2
  100. mylist = concatenate((arange(p2, n), arange(p2)))
  101. y = take(y, mylist, k)
  102. return y
  103. def fftfreq(n, d=1.0):
  104. """
  105. Return the Discrete Fourier Transform sample frequencies.
  106. The returned float array `f` contains the frequency bin centers in cycles
  107. per unit of the sample spacing (with zero at the start). For instance, if
  108. the sample spacing is in seconds, then the frequency unit is cycles/second.
  109. Given a window length `n` and a sample spacing `d`::
  110. f = [0, 1, ..., n/2-1, -n/2, ..., -1] / (d*n) if n is even
  111. f = [0, 1, ..., (n-1)/2, -(n-1)/2, ..., -1] / (d*n) if n is odd
  112. Parameters
  113. ----------
  114. n : int
  115. Window length.
  116. d : scalar, optional
  117. Sample spacing (inverse of the sampling rate). Defaults to 1.
  118. Returns
  119. -------
  120. f : ndarray
  121. Array of length `n` containing the sample frequencies.
  122. Examples
  123. --------
  124. >>> signal = np.array([-2, 8, 6, 4, 1, 0, 3, 5], dtype=float)
  125. >>> fourier = np.fft.fft(signal)
  126. >>> n = signal.size
  127. >>> timestep = 0.1
  128. >>> freq = np.fft.fftfreq(n, d=timestep)
  129. >>> freq
  130. array([ 0. , 1.25, 2.5 , 3.75, -5. , -3.75, -2.5 , -1.25])
  131. """
  132. if not isinstance(n, integer_types):
  133. raise ValueError("n should be an integer")
  134. val = 1.0 / (n * d)
  135. results = empty(n, int)
  136. N = (n-1)//2 + 1
  137. p1 = arange(0, N, dtype=int)
  138. results[:N] = p1
  139. p2 = arange(-(n//2), 0, dtype=int)
  140. results[N:] = p2
  141. return results * val
  142. #return hstack((arange(0,(n-1)/2 + 1), arange(-(n/2),0))) / (n*d)
  143. def rfftfreq(n, d=1.0):
  144. """
  145. Return the Discrete Fourier Transform sample frequencies
  146. (for usage with rfft, irfft).
  147. The returned float array `f` contains the frequency bin centers in cycles
  148. per unit of the sample spacing (with zero at the start). For instance, if
  149. the sample spacing is in seconds, then the frequency unit is cycles/second.
  150. Given a window length `n` and a sample spacing `d`::
  151. f = [0, 1, ..., n/2-1, n/2] / (d*n) if n is even
  152. f = [0, 1, ..., (n-1)/2-1, (n-1)/2] / (d*n) if n is odd
  153. Unlike `fftfreq` (but like `scipy.fftpack.rfftfreq`)
  154. the Nyquist frequency component is considered to be positive.
  155. Parameters
  156. ----------
  157. n : int
  158. Window length.
  159. d : scalar, optional
  160. Sample spacing (inverse of the sampling rate). Defaults to 1.
  161. Returns
  162. -------
  163. f : ndarray
  164. Array of length ``n//2 + 1`` containing the sample frequencies.
  165. Examples
  166. --------
  167. >>> signal = np.array([-2, 8, 6, 4, 1, 0, 3, 5, -3, 4], dtype=float)
  168. >>> fourier = np.fft.rfft(signal)
  169. >>> n = signal.size
  170. >>> sample_rate = 100
  171. >>> freq = np.fft.fftfreq(n, d=1./sample_rate)
  172. >>> freq
  173. array([ 0., 10., 20., 30., 40., -50., -40., -30., -20., -10.])
  174. >>> freq = np.fft.rfftfreq(n, d=1./sample_rate)
  175. >>> freq
  176. array([ 0., 10., 20., 30., 40., 50.])
  177. """
  178. if not isinstance(n, integer_types):
  179. raise ValueError("n should be an integer")
  180. val = 1.0/(n*d)
  181. N = n//2 + 1
  182. results = arange(0, N, dtype=int)
  183. return results * val