DiscreteFourierTransform¶
-
class
odl.trafos.fourier.
DiscreteFourierTransform
(*args, **kwargs)[source]¶ Bases:
odl.trafos.fourier.DiscreteFourierTransformBase
Plain 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.fftn
n-dimensional FFT routine
numpy.fft.rfftn
n-dimensional half-complex FFT
odl.trafos.backends.pyfftw_bindings.pyfftw_call
apply an FFTW transform
References
- Attributes
adjoint
Adjoint transform, equal to the inverse.
axes
Axes along the FT is calculated by this operator.
domain
Set of objects on which this operator can be evaluated.
halfcomplex
Return
True
if the last transform axis is halved.impl
Backend for the FFT implementation.
inverse
Inverse Fourier transform.
is_functional
True
if this operator’s range is aField
.is_linear
True
if this operator is linear.range
Set in which the result of an evaluation of this operator lies.
sign
Sign of the complex exponent in the transform.
Methods
_call
(self, x, out, \*\*kwargs)Implement
self(x, out[, **kwargs])
.clear_fftw_plan
(self)Delete the FFTW plan of this transform.
derivative
(self, point)Return the operator derivative at
point
.init_fftw_plan
(self[, planning_effort])Initialize the FFTW plan for this transform for later use.
norm
(self[, estimate])Return the operator norm of this operator.
-
__init__
(self, 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.exponent
is 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
domain
and the other parameters as auniform_discr
with 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.
None
means all axes.- sign{‘-‘, ‘+’}, optional
Sign of the complex exponent.
- halfcomplexbool, optional
If
True
, calculate only the negative frequency part along the last axis inaxes
for real input. This reduces the size of the range tofloor(N[i]/2) + 1
in this axisi
, whereN
is the shape of the input arrays. Otherwise, calculate the full complex FFT. Ifdom_dtype
is a complex type, this option has no effect.- impl{‘numpy’, ‘pyfftw’}, optional
Backend for the FFT implementation. The
'pyfftw'
backend is faster but requires thepyfftw
package.None
selects 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 + 1
in 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)