DiscreteFourierTransformInverse¶
- class odl.trafos.fourier.DiscreteFourierTransformInverse(*args, **kwargs)[source]¶
Bases:
DiscreteFourierTransformBasePlain backward DFT, only evaluating the trigonometric sum.
This operator calculates the inverse DFT:
f[k] = 1/prod(N) * sum_j( f_hat[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
DiscreteFourierTransformFourierTransformInversenumpy.fft.ifftnn-dimensional inverse FFT routine
numpy.fft.irfftnn-dimensional half-complex inverse 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__(range, domain=None, axes=None, sign='+', halfcomplex=False, impl=None)[source]¶
Initialize a new instance.
- Parameters:
- range
DiscretizedSpace Range of the inverse Fourier transform. If its
DiscretizedSpace.exponentis equal to 2.0, this operator has an adjoint which is equal to the inverse.- domain
DiscretizedSpace, optional Domain of the inverse Fourier transform. If not given, the domain is determined from
rangeand the other parameters as auniform_discrwith unit cell size and exponentp / (p - 1)(read as 'inf' for p=1 and 1 for p='inf').- axessequence 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, interpret the last axis inaxesas the negative frequency part of the transform of a real signal and calculate a "half-complex-to-real" inverse FFT. In this case, the domain has by default the shapefloor(N[i]/2) + 1in this axisi. Otherwise, domain and range have the same shape. Ifrangeis a complex space, this option has no effect.- impl{'numpy', 'pyfftw'}, optional
Backend for the FFT implementation. The 'pyfftw' backend is faster but requires the
pyfftwpackage.Noneselects the fastest available backend.
- range
Examples
Complex-to-complex (default) transforms have the same grids in domain and range:
>>> range = odl.uniform_discr([0, 0], [2, 4], (2, 4), ... nodes_on_bdry=True) >>> ifft = DiscreteFourierTransformInverse(range) >>> ifft.domain.shape (2, 4) >>> ifft.range.shape (2, 4)
Complex-to-real transforms have a domain grid with shape
n // 2 + 1in the last tranform axis:>>> range = odl.uniform_discr([0, 0, 0], [2, 3, 4], (2, 3, 4), ... nodes_on_bdry=True) >>> axes = (0, 1) >>> ifft = DiscreteFourierTransformInverse( ... range, halfcomplex=True, axes=axes) >>> ifft.domain.shape # shortened in the second axis (2, 2, 4) >>> ifft.range.shape (2, 3, 4)