summaryrefslogtreecommitdiff
path: root/numpy/fft/helper.py
diff options
context:
space:
mode:
authorPauli Virtanen <pav@iki.fi>2009-10-02 19:30:47 +0000
committerPauli Virtanen <pav@iki.fi>2009-10-02 19:30:47 +0000
commitbede419d707fef62166352a46fa7b6b76e1a13e9 (patch)
treec7ec89cbae61cac99c5b1eec7870dee2b0f01eb7 /numpy/fft/helper.py
parent30d4be0af19a433f9cc1f26142d509be3a8a8de5 (diff)
downloadnumpy-bede419d707fef62166352a46fa7b6b76e1a13e9.tar.gz
Docstring update: fft
Diffstat (limited to 'numpy/fft/helper.py')
-rw-r--r--numpy/fft/helper.py73
1 files changed, 57 insertions, 16 deletions
diff --git a/numpy/fft/helper.py b/numpy/fft/helper.py
index e2e36c323..e08f4b905 100644
--- a/numpy/fft/helper.py
+++ b/numpy/fft/helper.py
@@ -11,21 +11,46 @@ import types
def fftshift(x,axes=None):
"""
- Shift zero-frequency component to center of spectrum.
+ Shift the zero-frequency component to the center of the spectrum.
This function swaps half-spaces for all axes listed (defaults to all).
- If len(x) is even then the Nyquist component is y[0].
+ 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.
+ Axes over which to shift. Default is None, which shifts all axes.
+
+ Returns
+ -------
+ y : ndarray
+ The shifted array.
See Also
--------
- ifftshift
+ ifftshift : The inverse of `fftshift`.
+
+ Examples
+ --------
+ >>> freqs = np.fft.fftfreq(10, 0.1)
+ >>> freqs
+ array([ 0., 1., 2., 3., 4., -5., -4., -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.]])
"""
tmp = asarray(x)
@@ -43,18 +68,35 @@ def fftshift(x,axes=None):
def ifftshift(x,axes=None):
"""
- Inverse of fftshift.
+ The inverse of fftshift.
Parameters
----------
x : array_like
Input array.
axes : int or shape tuple, optional
- Axes over which to calculate. Defaults to None which is over all axes.
+ Axes over which to calculate. Defaults to None, which shifts all axes.
+
+ Returns
+ -------
+ y : ndarray
+ The shifted array.
See Also
--------
- fftshift
+ 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.]])
"""
tmp = asarray(x)
@@ -71,15 +113,14 @@ def ifftshift(x,axes=None):
def fftfreq(n,d=1.0):
"""
- Discrete Fourier Transform sample frequencies.
+ Return the Discrete Fourier Transform sample frequencies.
The returned float array contains the frequency bins in
cycles/unit (with zero at the start) given a window length `n` and a
- sample spacing `d`.
- ::
+ 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
+ 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
----------
@@ -90,14 +131,14 @@ def fftfreq(n,d=1.0):
Returns
-------
- out : ndarray, shape(`n`,)
- Sample frequencies.
+ out : ndarray
+ The array of length `n`, containing the sample frequencies.
Examples
--------
- >>> signal = np.array([-2., 8., -6., 4., 1., 0., 3., 5.])
+ >>> signal = np.array([-2, 8, 6, 4, 1, 0, 3, 5], dtype=float)
>>> fourier = np.fft.fft(signal)
- >>> n = len(signal)
+ >>> n = signal.size
>>> timestep = 0.1
>>> freq = np.fft.fftfreq(n, d=timestep)
>>> freq