cupyx.scipy.fftpack.get_fft_plan#
- cupyx.scipy.fftpack.get_fft_plan(a, shape=None, axes=None, value_type='C2C')[source]#
Generate a CUDA FFT plan for transforming up to three axes.
- Parameters:
a (cupy.ndarray) – Array to be transform, assumed to be either C- or F- contiguous.
shape (None or tuple of ints) – Shape of the transformed axes of the output. If
shape
is not given, the lengths of the input along the axes specified byaxes
are used.axes (None or int or tuple of int) –
The axes of the array to transform. If None, it is assumed that all axes are transformed.
Currently, for performing N-D transform these must be a set of up to three adjacent axes, and must include either the first or the last axis of the array.
value_type (str) –
The FFT type to perform. Acceptable values are:
’C2C’: complex-to-complex transform (default)
’R2C’: real-to-complex transform
’C2R’: complex-to-real transform
- Returns:
a cuFFT plan for either 1D transform (
cupy.cuda.cufft.Plan1d
) or N-D transform (cupy.cuda.cufft.PlanNd
).
Note
The returned plan can not only be passed as one of the arguments of the functions in
cupyx.scipy.fftpack
, but also be used as a context manager for bothcupy.fft
andcupyx.scipy.fftpack
functions:x = cupy.random.random(16).reshape(4, 4).astype(complex) plan = cupyx.scipy.fftpack.get_fft_plan(x) with plan: y = cupy.fft.fftn(x) # alternatively: y = cupyx.scipy.fftpack.fftn(x) # no explicit plan is given! # alternatively: y = cupyx.scipy.fftpack.fftn(x, plan=plan) # pass plan explicitly
In the first case, no cuFFT plan will be generated automatically, even if
cupy.fft.config.enable_nd_planning = True
is set.Note
If this function is called under the context of
set_cufft_callbacks()
, the generated plan will have callbacks enabled.Warning
This API is a deviation from SciPy’s, is currently experimental, and may be changed in the future version.