File size: 6,375 Bytes
dc2106c
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
"""

Discrete Fourier Transforms - helper.py



"""
from numpy.core import integer, empty, arange, asarray, roll
from numpy.core.overrides import array_function_dispatch, set_module

# Created by Pearu Peterson, September 2002

__all__ = ['fftshift', 'ifftshift', 'fftfreq', 'rfftfreq']

integer_types = (int, integer)


def _fftshift_dispatcher(x, axes=None):
    return (x,)


@array_function_dispatch(_fftshift_dispatcher, module='numpy.fft')
def fftshift(x, axes=None):
    """

    Shift the zero-frequency component to the center of the spectrum.



    This function swaps half-spaces for all axes listed (defaults to all).

    Note that ``y[0]`` is the Nyquist component only if ``len(x)`` is even.



    Parameters

    ----------

    x : array_like

        Input array.

    axes : int or shape tuple, optional

        Axes over which to shift.  Default is None, which shifts all axes.



    Returns

    -------

    y : ndarray

        The shifted array.



    See Also

    --------

    ifftshift : The inverse of `fftshift`.



    Examples

    --------

    >>> freqs = np.fft.fftfreq(10, 0.1)

    >>> freqs

    array([ 0.,  1.,  2., ..., -3., -2., -1.])

    >>> np.fft.fftshift(freqs)

    array([-5., -4., -3., -2., -1.,  0.,  1.,  2.,  3.,  4.])



    Shift the zero-frequency component only along the second axis:



    >>> freqs = np.fft.fftfreq(9, d=1./9).reshape(3, 3)

    >>> freqs

    array([[ 0.,  1.,  2.],

           [ 3.,  4., -4.],

           [-3., -2., -1.]])

    >>> np.fft.fftshift(freqs, axes=(1,))

    array([[ 2.,  0.,  1.],

           [-4.,  3.,  4.],

           [-1., -3., -2.]])



    """
    x = asarray(x)
    if axes is None:
        axes = tuple(range(x.ndim))
        shift = [dim // 2 for dim in x.shape]
    elif isinstance(axes, integer_types):
        shift = x.shape[axes] // 2
    else:
        shift = [x.shape[ax] // 2 for ax in axes]

    return roll(x, shift, axes)


@array_function_dispatch(_fftshift_dispatcher, module='numpy.fft')
def ifftshift(x, axes=None):
    """

    The inverse of `fftshift`. Although identical for even-length `x`, the

    functions differ by one sample for odd-length `x`.



    Parameters

    ----------

    x : array_like

        Input array.

    axes : int or shape tuple, optional

        Axes over which to calculate.  Defaults to None, which shifts all axes.



    Returns

    -------

    y : ndarray

        The shifted array.



    See Also

    --------

    fftshift : Shift zero-frequency component to the center of the spectrum.



    Examples

    --------

    >>> freqs = np.fft.fftfreq(9, d=1./9).reshape(3, 3)

    >>> freqs

    array([[ 0.,  1.,  2.],

           [ 3.,  4., -4.],

           [-3., -2., -1.]])

    >>> np.fft.ifftshift(np.fft.fftshift(freqs))

    array([[ 0.,  1.,  2.],

           [ 3.,  4., -4.],

           [-3., -2., -1.]])



    """
    x = asarray(x)
    if axes is None:
        axes = tuple(range(x.ndim))
        shift = [-(dim // 2) for dim in x.shape]
    elif isinstance(axes, integer_types):
        shift = -(x.shape[axes] // 2)
    else:
        shift = [-(x.shape[ax] // 2) for ax in axes]

    return roll(x, shift, axes)


@set_module('numpy.fft')
def fftfreq(n, d=1.0):
    """

    Return the Discrete Fourier Transform sample frequencies.



    The returned float array `f` contains the frequency bin centers in cycles

    per unit of the sample spacing (with zero at the start).  For instance, if

    the sample spacing is in seconds, then the frequency unit is cycles/second.



    Given a window length `n` and a sample spacing `d`::



      f = [0, 1, ...,   n/2-1,     -n/2, ..., -1] / (d*n)   if n is even

      f = [0, 1, ..., (n-1)/2, -(n-1)/2, ..., -1] / (d*n)   if n is odd



    Parameters

    ----------

    n : int

        Window length.

    d : scalar, optional

        Sample spacing (inverse of the sampling rate). Defaults to 1.



    Returns

    -------

    f : ndarray

        Array of length `n` containing the sample frequencies.



    Examples

    --------

    >>> signal = np.array([-2, 8, 6, 4, 1, 0, 3, 5], dtype=float)

    >>> fourier = np.fft.fft(signal)

    >>> n = signal.size

    >>> timestep = 0.1

    >>> freq = np.fft.fftfreq(n, d=timestep)

    >>> freq

    array([ 0.  ,  1.25,  2.5 , ..., -3.75, -2.5 , -1.25])



    """
    if not isinstance(n, integer_types):
        raise ValueError("n should be an integer")
    val = 1.0 / (n * d)
    results = empty(n, int)
    N = (n-1)//2 + 1
    p1 = arange(0, N, dtype=int)
    results[:N] = p1
    p2 = arange(-(n//2), 0, dtype=int)
    results[N:] = p2
    return results * val


@set_module('numpy.fft')
def rfftfreq(n, d=1.0):
    """

    Return the Discrete Fourier Transform sample frequencies

    (for usage with rfft, irfft).



    The returned float array `f` contains the frequency bin centers in cycles

    per unit of the sample spacing (with zero at the start).  For instance, if

    the sample spacing is in seconds, then the frequency unit is cycles/second.



    Given a window length `n` and a sample spacing `d`::



      f = [0, 1, ...,     n/2-1,     n/2] / (d*n)   if n is even

      f = [0, 1, ..., (n-1)/2-1, (n-1)/2] / (d*n)   if n is odd



    Unlike `fftfreq` (but like `scipy.fftpack.rfftfreq`)

    the Nyquist frequency component is considered to be positive.



    Parameters

    ----------

    n : int

        Window length.

    d : scalar, optional

        Sample spacing (inverse of the sampling rate). Defaults to 1.



    Returns

    -------

    f : ndarray

        Array of length ``n//2 + 1`` containing the sample frequencies.



    Examples

    --------

    >>> signal = np.array([-2, 8, 6, 4, 1, 0, 3, 5, -3, 4], dtype=float)

    >>> fourier = np.fft.rfft(signal)

    >>> n = signal.size

    >>> sample_rate = 100

    >>> freq = np.fft.fftfreq(n, d=1./sample_rate)

    >>> freq

    array([  0.,  10.,  20., ..., -30., -20., -10.])

    >>> freq = np.fft.rfftfreq(n, d=1./sample_rate)

    >>> freq

    array([  0.,  10.,  20.,  30.,  40.,  50.])



    """
    if not isinstance(n, integer_types):
        raise ValueError("n should be an integer")
    val = 1.0/(n*d)
    N = n//2 + 1
    results = arange(0, N, dtype=int)
    return results * val