DiscreteFourierTransform¶
- class odl.trafos.fourier.DiscreteFourierTransform(*args, **kwargs)[source]¶
Bases:
DiscreteFourierTransformBasePlain forward DFT, only evaluating the trigonometric sum.
This operator calculates the forward DFT:
f_hat[k] = sum_j( f[j] * exp(-+ 1j*2*pi * j*k/N) )
without any further shifting or scaling compensation. See the Numpy FFT documentation, the pyfftw API documentation or What FFTW really computes for further information.
See also
numpy.fft.fftnn-dimensional FFT routine
numpy.fft.rfftnn-dimensional half-complex FFT
odl.trafos.backends.pyfftw_bindings.pyfftw_callapply an FFTW transform
References
- Attributes:
adjointAdjoint transform, equal to the inverse.
axesAxes along the FT is calculated by this operator.
domainSet of objects on which this operator can be evaluated.
halfcomplexReturn
Trueif the last transform axis is halved.implBackend for the FFT implementation.
inverseInverse Fourier transform.
is_functionalTrueif this operator's range is aField.is_linearTrueif this operator is linear.rangeSet in which the result of an evaluation of this operator lies.
signSign of the complex exponent in the transform.
Methods
__call__(x[, out])Return
self(x[, out, **kwargs]).Delete the FFTW plan of this transform.
derivative(point)Return the operator derivative at
point.init_fftw_plan([planning_effort])Initialize the FFTW plan for this transform for later use.
norm([estimate])Return the operator norm of this operator.
- __init__(domain, range=None, axes=None, sign='-', halfcomplex=False, impl=None)[source]¶
Initialize a new instance.
- Parameters:
- domain
DiscretizedSpace Domain of the Fourier transform. If its
DiscretizedSpace.exponentis equal to 2.0, this operator has an adjoint which is equal to the inverse.- range
DiscretizedSpace, optional Range of the Fourier transform. If not given, the range is determined from
domainand the other parameters as auniform_discrwith unit cell size and exponentp / (p - 1)(read as 'inf' for p=1 and 1 for p='inf').- axesint or sequence of ints, optional
Dimensions in which a transform is to be calculated.
Nonemeans all axes.- sign{'-', '+'}, optional
Sign of the complex exponent.
- halfcomplexbool, optional
If
True, calculate only the negative frequency part along the last axis inaxesfor real input. This reduces the size of the range tofloor(N[i]/2) + 1in this axisi, whereNis the shape of the input arrays. Otherwise, calculate the full complex FFT. Ifdom_dtypeis a complex type, this option has no effect.- impl{'numpy', 'pyfftw'}, optional
Backend for the FFT implementation. The
'pyfftw'backend is faster but requires thepyfftwpackage.Noneselects the fastest available backend.
- domain
Examples
Complex-to-complex (default) transforms have the same grids in domain and range:
>>> domain = odl.uniform_discr([0, 0], [2, 4], (2, 4), ... nodes_on_bdry=True) >>> fft = DiscreteFourierTransform(domain) >>> fft.domain.shape (2, 4) >>> fft.range.shape (2, 4)
Real-to-complex transforms have a range grid with shape
n // 2 + 1in the last tranform axis:>>> domain = odl.uniform_discr([0, 0, 0], [2, 3, 4], (2, 3, 4), ... nodes_on_bdry=True) >>> axes = (0, 1) >>> fft = DiscreteFourierTransform( ... domain, halfcomplex=True, axes=axes) >>> fft.range.shape # shortened in the second axis (2, 2, 4) >>> fft.domain.shape (2, 3, 4)